InflationInjection.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "InflationInjection.H"
27 #include "mathematicalConstants.H"
28 #include "PackedBoolList.H"
29 #include "cellSet.H"
30 #include "ListListOps.H"
31 
32 using namespace Foam::constant::mathematical;
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 template<class CloudType>
38 (
39  const dictionary& dict,
40  CloudType& owner,
41  const word& modelName
42 )
43 :
44  InjectionModel<CloudType>(dict, owner, modelName, typeName),
45  generationSetName_(this->coeffDict().lookup("generationCellSet")),
46  inflationSetName_(this->coeffDict().lookup("inflationCellSet")),
47  generationCells_(),
48  inflationCells_(),
49  duration_(readScalar(this->coeffDict().lookup("duration"))),
50  flowRateProfile_
51  (
53  (
54  owner.db().time(),
55  "flowRateProfile",
56  this->coeffDict()
57  )
58  ),
59  growthRate_
60  (
62  (
63  owner.db().time(),
64  "growthRate",
65  this->coeffDict()
66  )
67  ),
68  newParticles_(),
69  volumeAccumulator_(0.0),
70  fraction_(1.0),
71  selfSeed_(this->coeffDict().lookupOrDefault("selfSeed", false)),
72  dSeed_(SMALL),
73  sizeDistribution_
74  (
76  (
77  this->coeffDict().subDict("sizeDistribution"),
78  owner.rndGen()
79  )
80  )
81 {
82  duration_ = owner.db().time().userTimeToTime(duration_);
83 
84  if (selfSeed_)
85  {
86  dSeed_ = readScalar(this->coeffDict().lookup("dSeed"));
87  }
88 
89  cellSet generationCells(this->owner().mesh(), generationSetName_);
90 
91  generationCells_ = generationCells.toc();
92 
93  cellSet inflationCells(this->owner().mesh(), inflationSetName_);
94 
95  // Union of cellSets
96  inflationCells |= generationCells;
97 
98  inflationCells_ = inflationCells.toc();
99 
100  if (Pstream::parRun())
101  {
102  scalar generationVolume = 0.0;
103 
104  forAll(generationCells_, gCI)
105  {
106  label cI = generationCells_[gCI];
107 
108  generationVolume += this->owner().mesh().cellVolumes()[cI];
109  }
110 
111  scalar totalGenerationVolume = generationVolume;
112 
113  reduce(totalGenerationVolume, sumOp<scalar>());
114 
115  fraction_ = generationVolume/totalGenerationVolume;
116  }
117 
118  // Set total volume/mass to inject
119  this->volumeTotal_ = fraction_*flowRateProfile_.integrate(0.0, duration_);
120  this->massTotal_ *= fraction_;
121 }
122 
123 
124 template<class CloudType>
126 (
128 )
129 :
131  generationSetName_(im.generationSetName_),
132  inflationSetName_(im.inflationSetName_),
133  generationCells_(im.generationCells_),
134  inflationCells_(im.inflationCells_),
135  duration_(im.duration_),
136  flowRateProfile_(im.flowRateProfile_),
137  growthRate_(im.growthRate_),
138  newParticles_(im.newParticles_),
139  volumeAccumulator_(im.volumeAccumulator_),
140  fraction_(im.fraction_),
141  selfSeed_(im.selfSeed_),
142  dSeed_(im.dSeed_),
143  sizeDistribution_(im.sizeDistribution_, false)
144 {}
145 
146 
147 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
148 
149 template<class CloudType>
151 {}
152 
153 
154 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
155 
156 template<class CloudType>
158 {
159  // do nothing
160 }
161 
162 
163 template<class CloudType>
165 {
166  return this->SOI_ + duration_;
167 }
168 
169 
170 template<class CloudType>
172 (
173  const scalar time0,
174  const scalar time1
175 )
176 {
177  const polyMesh& mesh = this->owner().mesh();
178 
180  this->owner().cellOccupancy();
181 
182  scalar gR = growthRate_.value(time1);
183 
184  scalar dT = time1 - time0;
185 
186  // Inflate existing particles
187 
188  forAll(inflationCells_, iCI)
189  {
190  label cI = inflationCells_[iCI];
191 
192  typename CloudType::parcelType* pPtr = NULL;
193 
194  forAll(cellOccupancy[cI], cPI)
195  {
196  pPtr = cellOccupancy[cI][cPI];
197 
198  scalar dTarget = pPtr->dTarget();
199 
200  pPtr->d() = min(dTarget, pPtr->d() + gR*dT);
201  }
202  }
203 
204  // Generate new particles
205 
206  newParticles_.clear();
207 
208  cachedRandom& rnd = this->owner().rndGen();
209 
210  // Diameter factor, when splitting particles into 4, this is the
211  // factor that modifies the diameter.
212  scalar dFact = sqrt(2.0)/(sqrt(3.0) + sqrt(2.0));
213 
214  if ((time0 >= 0.0) && (time0 < duration_))
215  {
216  volumeAccumulator_ +=
217  fraction_*flowRateProfile_.integrate(time0, time1);
218  }
219 
220  labelHashSet cellCentresUsed;
221 
222  // Loop escape counter
223  label maxIterations = max
224  (
225  1,
226  (10*volumeAccumulator_)
227  /CloudType::parcelType::volume(sizeDistribution_().minValue())
228  );
229 
230  label iterationNo = 0;
231 
232  // Info<< "Accumulated volume to inject: "
233  // << returnReduce(volumeAccumulator_, sumOp<scalar>()) << endl;
234 
235  while (!generationCells_.empty() && volumeAccumulator_ > 0)
236  {
237  if (iterationNo > maxIterations)
238  {
240  << "Maximum particle split iterations ("
241  << maxIterations << ") exceeded" << endl;
242 
243  break;
244  }
245 
246  label cI = generationCells_
247  [
248  rnd.position
249  (
250  label(0),
251  generationCells_.size() - 1
252  )
253  ];
254 
255  // Pick a particle at random from the cell - if there are
256  // none, insert one at the cell centre. Otherwise, split an
257  // existing particle into four new ones.
258 
259  if (cellOccupancy[cI].empty())
260  {
261  if (selfSeed_ && !cellCentresUsed.found(cI))
262  {
263  scalar dNew = sizeDistribution_().sample();
264 
265  newParticles_.append
266  (
268  (
269  Pair<vector>(mesh.cellCentres()[cI], vector::zero),
270  Pair<scalar>(dSeed_, dNew)
271  )
272  );
273 
274  volumeAccumulator_ -= CloudType::parcelType::volume(dNew);
275 
276  cellCentresUsed.insert(cI);
277  }
278  }
279  else
280  {
281  label cPI = rnd.position(label(0), cellOccupancy[cI].size() - 1);
282 
283  // This has to be a reference to the pointer so that it
284  // can be set to NULL when the particle is deleted.
285  typename CloudType::parcelType*& pPtr = cellOccupancy[cI][cPI];
286 
287  if (pPtr != NULL)
288  {
289  scalar pD = pPtr->d();
290 
291  // Select bigger particles by preference
292  if ((pD/pPtr->dTarget()) < rnd.sample01<scalar>())
293  {
294  continue;
295  }
296 
297  const point& pP = pPtr->position();
298  const vector& pU = pPtr->U();
299 
300  // Generate a tetrahedron of new positions with the
301  // four new spheres fitting inside the old one, where
302  // a is the diameter of the new spheres, and is
303  // related to the diameter of the enclosing sphere, A,
304  // by a = sqrt(2)*A/(sqrt(3) + sqrt(2));
305 
306  // Positions around the origin, which is the
307  // tetrahedron centroid (centre of old sphere).
308 
309  // x = a/sqrt(3)
310  // r = a/(2*sqrt(6))
311  // R = sqrt(3)*a/(2*sqrt(2))
312  // d = a/(2*sqrt(3))
313 
314  // p0(x, 0, -r)
315  // p1(-d, a/2, -r)
316  // p2(-d, -a/2, -r)
317  // p3(0, 0, R)
318 
319  scalar a = pD*dFact;
320 
321  scalar x = a/sqrt(3.0);
322  scalar r = a/(2.0*sqrt(6.0));
323  scalar R = sqrt(3.0)*a/(2.0*sqrt(2.0));
324  scalar d = a/(2.0*sqrt(3.0));
325 
326  scalar dNew = sizeDistribution_().sample();
327  scalar volNew = CloudType::parcelType::volume(dNew);
328 
329  newParticles_.append
330  (
332  (
333  Pair<vector>(vector(x, 0, -r) + pP, pU),
334  Pair<scalar>(a, dNew)
335  )
336  );
337  volumeAccumulator_ -= volNew;
338 
339  dNew = sizeDistribution_().sample();
340  newParticles_.append
341  (
343  (
344  Pair<vector>(vector(-d, a/2, -r) + pP, pU),
345  Pair<scalar>(a, dNew)
346  )
347  );
348  volumeAccumulator_ -= volNew;
349 
350  dNew = sizeDistribution_().sample();
351  newParticles_.append
352  (
354  (
355  Pair<vector>(vector(-d, -a/2, -r) + pP, pU),
356  Pair<scalar>(a, dNew)
357  )
358  );
359  volumeAccumulator_ -= volNew;
360 
361  dNew = sizeDistribution_().sample();
362  newParticles_.append
363  (
365  (
366  Pair<vector>(vector(0, 0, R) + pP, pU),
367  Pair<scalar>(a, dNew)
368  )
369  );
370  volumeAccumulator_ -= volNew;
371 
372  // Account for the lost volume of the particle which
373  // is to be deleted
374  volumeAccumulator_ += CloudType::parcelType::volume
375  (
376  pPtr->dTarget()
377  );
378 
379  this->owner().deleteParticle(*pPtr);
380 
381  pPtr = NULL;
382  }
383  }
384 
385  iterationNo++;
386  }
387 
388  if (Pstream::parRun())
389  {
390  List<List<vectorPairScalarPair> > gatheredNewParticles
391  (
393  );
394 
395  gatheredNewParticles[Pstream::myProcNo()] = newParticles_;
396 
397  // Gather data onto master
398  Pstream::gatherList(gatheredNewParticles);
399 
400  // Combine
401  List<vectorPairScalarPair> combinedNewParticles
402  (
404  (
405  gatheredNewParticles,
407  )
408  );
409 
410  if (Pstream::master())
411  {
412  newParticles_ = combinedNewParticles;
413  }
414 
415  Pstream::scatter(newParticles_);
416  }
417 
418  return newParticles_.size();
419 }
420 
421 
422 template<class CloudType>
424 (
425  const scalar time0,
426  const scalar time1
427 )
428 {
429  if ((time0 >= 0.0) && (time0 < duration_))
430  {
431  return fraction_*flowRateProfile_.integrate(time0, time1);
432  }
433  else
434  {
435  return 0.0;
436  }
437 }
438 
439 
440 template<class CloudType>
442 (
443  const label parcelI,
444  const label,
445  const scalar,
446  vector& position,
447  label& cellOwner,
448  label& tetFaceI,
449  label& tetPtI
450 )
451 {
452  position = newParticles_[parcelI].first().first();
453 
454  this->findCellAtPosition
455  (
456  cellOwner,
457  tetFaceI,
458  tetPtI,
459  position,
460  false
461  );
462 }
463 
464 
465 template<class CloudType>
467 (
468  const label parcelI,
469  const label,
470  const scalar,
471  typename CloudType::parcelType& parcel
472 )
473 {
474  parcel.U() = newParticles_[parcelI].first().second();
475 
476  parcel.d() = newParticles_[parcelI].second().first();
477 
478  parcel.dTarget() = newParticles_[parcelI].second().second();
479 }
480 
481 
482 template<class CloudType>
484 {
485  return false;
486 }
487 
488 
489 template<class CloudType>
491 {
492  return true;
493 }
494 
495 
496 // ************************************************************************* //
Foam::InflationInjection::inflationSetName_
word inflationSetName_
Name of cellSet for inflating new particles.
Definition: InflationInjection.H:73
mathematicalConstants.H
Foam::cachedRandom::sample01
Type sample01()
Return a sample whose components lie in the range 0-1.
Definition: cachedRandomTemplates.C:32
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::InflationInjection::InflationInjection
InflationInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: InflationInjection.C:38
Foam::InflationInjection::parcelsToInject
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
Definition: InflationInjection.C:172
Foam::HashTable::toc
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:201
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
ListListOps.H
Foam::InflationInjection::sizeDistribution_
const autoPtr< distributionModels::distributionModel > sizeDistribution_
Parcel size distribution model.
Definition: InflationInjection.H:110
Foam::ListListOps::combine
AccessType combine(const List< T > &, AccessOp aop=accessOp< T >())
Combines sublists into one big list.
Definition: ListListOps.C:34
Foam::UPstream::nProcs
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:387
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::accessOp
Definition: ListListOps.H:98
Foam::InjectionModel
Templated injection model class.
Definition: InjectionModel.H:67
Foam::Pstream::scatter
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
Definition: gatherScatter.C:147
Foam::distributionModels::distributionModel::New
static autoPtr< distributionModel > New(const dictionary &dict, cachedRandom &rndGen)
Selector.
Definition: distributionModelNew.C:32
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
InflationInjection.H
Foam::InflationInjection::flowRateProfile_
TimeDataEntry< scalar > flowRateProfile_
Flow rate profile relative to SOI [m3/s].
Definition: InflationInjection.H:86
Foam::InflationInjection::fullyDescribed
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
Definition: InflationInjection.C:483
Foam::HashSet< label, Hash< label > >
Foam::InflationInjection
Inflation injection - creates new particles by splitting existing particles within in a set of genera...
Definition: InflationInjection.H:63
Foam::TimeDataEntry< scalar >
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::DSMCCloud::rndGen
Random & rndGen()
Return refernce to the random object.
Definition: DSMCCloudI.H:121
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
R
#define R(A, B, C, D, E, F, K, M)
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
cellOccupancy
const List< DynamicList< molecule * > > & cellOccupancy
Definition: calculateMDFields.H:1
Foam::InflationInjection::generationCells_
labelList generationCells_
Set of cells to generate particles in.
Definition: InflationInjection.H:76
Foam::cachedRandom
Random number generator.
Definition: cachedRandom.H:63
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return refernce to the mesh.
Definition: DSMCCloudI.H:41
Foam::cachedRandom::position
Type position(const Type &start, const Type &end)
Return a sample between start and end.
Definition: cachedRandomTemplates.C:58
Foam::Cloud::deleteParticle
void deleteParticle(ParticleType &)
Remove particle from cloud and delete.
Definition: Cloud.C:169
PackedBoolList.H
Foam::InflationInjection::validInjection
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
Definition: InflationInjection.C:490
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
Foam::DSMCCloud::cellOccupancy
const List< DynamicList< ParcelType * > > & cellOccupancy() const
Return the cell occupancy addressing.
Definition: DSMCCloudI.H:72
Foam::InflationInjection::~InflationInjection
virtual ~InflationInjection()
Destructor.
Definition: InflationInjection.C:150
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::HashTable::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam::primitiveMesh::cellVolumes
const scalarField & cellVolumes() const
Definition: primitiveMeshCellCentresAndVols.C:222
Foam::cellSet
A collection of cell labels.
Definition: cellSet.H:48
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Foam::InflationInjection::duration_
scalar duration_
Injection duration [s].
Definition: InflationInjection.H:83
Foam::InflationInjection::selfSeed_
Switch selfSeed_
Switch to control whether or not the injector is allowed.
Definition: InflationInjection.H:104
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:399
Foam::InflationInjection::volumeToInject
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
Definition: InflationInjection.C:424
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::InflationInjection::setPositionAndCell
virtual void setPositionAndCell(const label parcelI, const label nParcels, const scalar time, vector &position, label &cellOwner, label &tetFaceI, label &tetPtI)
Set the injection position and owner cell, tetFace and tetPt.
Definition: InflationInjection.C:442
Foam::InflationInjection::dSeed_
scalar dSeed_
Diameter with which to create new seed particles.
Definition: InflationInjection.H:107
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:49
Foam::constant::mathematical
mathematical constants.
Foam::sumOp
Definition: ops.H:162
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::InflationInjection::volumeAccumulator_
scalar volumeAccumulator_
Accumulation variable to carry over volume from one injection.
Definition: InflationInjection.H:97
Foam::InflationInjection::setProperties
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
Definition: InflationInjection.C:467
Foam::InflationInjection::growthRate_
TimeDataEntry< scalar > growthRate_
Growth rate of particle diameters towards target [m/s].
Definition: InflationInjection.H:89
Foam::InflationInjection::generationSetName_
word generationSetName_
Name of cellSet for generating new particles.
Definition: InflationInjection.H:70
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::InflationInjection::timeEnd
scalar timeEnd() const
Return the end-of-injection time.
Definition: InflationInjection.C:164
Foam::InflationInjection::fraction_
scalar fraction_
Fraction of injection controlled by this processor.
Definition: InflationInjection.H:100
minValue
scalar minValue
Definition: LISASMDCalcMethod2.H:12
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217
Foam::InflationInjection::newParticles_
DynamicList< vectorPairScalarPair > newParticles_
Positions, velocities, diameters and target diameters of.
Definition: InflationInjection.H:93
cellSet.H
Foam::Tuple2
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:47
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::InflationInjection::inflationCells_
labelList inflationCells_
Set of cells to inflate particles in, includes all.
Definition: InflationInjection.H:80
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress
Foam::InflationInjection::updateMesh
virtual void updateMesh()
Set injector locations when mesh is updated.
Definition: InflationInjection.C:157