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