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-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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"
31#include "bitSet.H"
32#include "cellSet.H"
33#include "ListListOps.H"
34
35using namespace Foam::constant::mathematical;
36
37// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38
39template<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 (
55 Function1<scalar>::New
56 (
57 "flowRateProfile",
58 this->coeffDict(),
59 &owner.mesh()
60 )
61 ),
62 growthRate_
63 (
64 Function1<scalar>::New
65 (
66 "growthRate",
67 this->coeffDict(),
68 &owner.mesh()
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 // Convert from user time to reduce the number of time conversion calls
86 const Time& time = owner.db().time();
87 duration_ = time.userTimeToTime(duration_);
88 flowRateProfile_->userTimeToTime(time);
89 growthRate_->userTimeToTime(time);
90
91 if (selfSeed_)
92 {
93 this->coeffDict().readEntry("dSeed", dSeed_);
94 }
95
96 cellSet generationCells(this->owner().mesh(), generationSetName_);
97
98 generationCells_ = generationCells.toc();
99
100 cellSet inflationCells(this->owner().mesh(), inflationSetName_);
101
102 // Union of cellSets
103 inflationCells |= generationCells;
104
105 inflationCells_ = inflationCells.toc();
106
107 if (Pstream::parRun())
108 {
109 scalar generationVolume = 0.0;
110
111 forAll(generationCells_, gCI)
112 {
113 label cI = generationCells_[gCI];
114
115 generationVolume += this->owner().mesh().cellVolumes()[cI];
116 }
117
118 scalar totalGenerationVolume = generationVolume;
119
120 reduce(totalGenerationVolume, sumOp<scalar>());
121
122 fraction_ = generationVolume/totalGenerationVolume;
123 }
124
125 // Set total volume/mass to inject
126 this->volumeTotal_ = fraction_*flowRateProfile_->integrate(0.0, duration_);
127 this->massTotal_ *= fraction_;
128}
129
130
131template<class CloudType>
133(
135)
136:
138 generationSetName_(im.generationSetName_),
139 inflationSetName_(im.inflationSetName_),
140 generationCells_(im.generationCells_),
141 inflationCells_(im.inflationCells_),
142 duration_(im.duration_),
143 flowRateProfile_(im.flowRateProfile_.clone()),
144 growthRate_(im.growthRate_.clone()),
145 newParticles_(im.newParticles_),
146 volumeAccumulator_(im.volumeAccumulator_),
147 fraction_(im.fraction_),
148 selfSeed_(im.selfSeed_),
149 dSeed_(im.dSeed_),
150 sizeDistribution_(im.sizeDistribution_.clone())
151{}
152
153
154// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
155
156template<class CloudType>
158{}
159
160
161// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162
163template<class CloudType>
165{}
166
167
168template<class CloudType>
170{
171 return this->SOI_ + duration_;
172}
173
174
175template<class CloudType>
177(
178 const scalar time0,
179 const scalar time1
180)
181{
182 const polyMesh& mesh = this->owner().mesh();
183
185 this->owner().cellOccupancy();
186
187 scalar gR = growthRate_->value(time1);
188
189 scalar dT = time1 - time0;
190
191 // Inflate existing particles
192
193 forAll(inflationCells_, iCI)
194 {
195 label cI = inflationCells_[iCI];
196
197 typename CloudType::parcelType* pPtr = nullptr;
198
199 forAll(cellOccupancy[cI], cPI)
200 {
201 pPtr = cellOccupancy[cI][cPI];
202
203 scalar dTarget = pPtr->dTarget();
204
205 pPtr->d() = min(dTarget, pPtr->d() + gR*dT);
206 }
207 }
208
209 // Generate new particles
210
211 newParticles_.clear();
212
213 Random& rnd = this->owner().rndGen();
214
215 // Diameter factor, when splitting particles into 4, this is the
216 // factor that modifies the diameter.
217 scalar dFact = sqrt(2.0)/(sqrt(3.0) + sqrt(2.0));
218
219 if ((time0 >= 0.0) && (time0 < duration_))
220 {
221 volumeAccumulator_ +=
222 fraction_*flowRateProfile_->integrate(time0, time1);
223 }
224
225 labelHashSet cellCentresUsed;
226
227 // Loop escape counter
228 label maxIterations = max
229 (
230 1,
231 (10*volumeAccumulator_)
232 /CloudType::parcelType::volume(sizeDistribution_().minValue())
233 );
234
235 label iterationNo = 0;
236
237 // Info<< "Accumulated volume to inject: "
238 // << returnReduce(volumeAccumulator_, sumOp<scalar>()) << endl;
239
240 while (!generationCells_.empty() && volumeAccumulator_ > 0)
241 {
242 if (iterationNo > maxIterations)
243 {
245 << "Maximum particle split iterations ("
246 << maxIterations << ") exceeded" << endl;
247
248 break;
249 }
250
251 label cI = generationCells_
252 [
253 rnd.position
254 (
255 label(0),
256 generationCells_.size() - 1
257 )
258 ];
259
260 // Pick a particle at random from the cell - if there are
261 // none, insert one at the cell centre. Otherwise, split an
262 // existing particle into four new ones.
263
264 if (cellOccupancy[cI].empty())
265 {
266 if (selfSeed_ && !cellCentresUsed.found(cI))
267 {
268 scalar dNew = sizeDistribution_().sample();
269
270 newParticles_.append
271 (
273 (
275 Pair<scalar>(dSeed_, dNew)
276 )
277 );
278
279 volumeAccumulator_ -= CloudType::parcelType::volume(dNew);
280
281 cellCentresUsed.insert(cI);
282 }
283 }
284 else
285 {
286 label cPI = rnd.position(label(0), cellOccupancy[cI].size() - 1);
287
288 // This has to be a reference to the pointer so that it
289 // can be set to nullptr when the particle is deleted.
290 typename CloudType::parcelType*& pPtr = cellOccupancy[cI][cPI];
291
292 if (pPtr != nullptr)
293 {
294 scalar pD = pPtr->d();
295
296 // Select bigger particles by preference
297 if ((pD/pPtr->dTarget()) < rnd.sample01<scalar>())
298 {
299 continue;
300 }
301
302 const point pP(pPtr->position());
303 const vector& pU = pPtr->U();
304
305 // Generate a tetrahedron of new positions with the
306 // four new spheres fitting inside the old one, where
307 // a is the diameter of the new spheres, and is
308 // related to the diameter of the enclosing sphere, A,
309 // by a = sqrt(2)*A/(sqrt(3) + sqrt(2));
310
311 // Positions around the origin, which is the
312 // tetrahedron centroid (centre of old sphere).
313
314 // x = a/sqrt(3)
315 // r = a/(2*sqrt(6))
316 // R = sqrt(3)*a/(2*sqrt(2))
317 // d = a/(2*sqrt(3))
318
319 // p0(x, 0, -r)
320 // p1(-d, a/2, -r)
321 // p2(-d, -a/2, -r)
322 // p3(0, 0, R)
323
324 scalar a = pD*dFact;
325
326 scalar x = a/sqrt(3.0);
327 scalar r = a/(2.0*sqrt(6.0));
328 scalar R = sqrt(3.0)*a/(2.0*sqrt(2.0));
329 scalar d = a/(2.0*sqrt(3.0));
330
331 scalar dNew = sizeDistribution_().sample();
332 scalar volNew = CloudType::parcelType::volume(dNew);
333
334 newParticles_.append
335 (
337 (
338 Pair<vector>(vector(x, 0, -r) + pP, pU),
339 Pair<scalar>(a, dNew)
340 )
341 );
342 volumeAccumulator_ -= volNew;
343
344 dNew = sizeDistribution_().sample();
345 newParticles_.append
346 (
348 (
349 Pair<vector>(vector(-d, a/2, -r) + pP, pU),
350 Pair<scalar>(a, dNew)
351 )
352 );
353 volumeAccumulator_ -= volNew;
354
355 dNew = sizeDistribution_().sample();
356 newParticles_.append
357 (
359 (
360 Pair<vector>(vector(-d, -a/2, -r) + pP, pU),
361 Pair<scalar>(a, dNew)
362 )
363 );
364 volumeAccumulator_ -= volNew;
365
366 dNew = sizeDistribution_().sample();
367 newParticles_.append
368 (
370 (
371 Pair<vector>(vector(0, 0, R) + pP, pU),
372 Pair<scalar>(a, dNew)
373 )
374 );
375 volumeAccumulator_ -= volNew;
376
377 // Account for the lost volume of the particle which
378 // is to be deleted
379 volumeAccumulator_ += CloudType::parcelType::volume
380 (
381 pPtr->dTarget()
382 );
383
384 this->owner().deleteParticle(*pPtr);
385
386 pPtr = nullptr;
387 }
388 }
389
390 iterationNo++;
391 }
392
393 if (Pstream::parRun())
394 {
395 List<List<vectorPairScalarPair>> gatheredNewParticles
396 (
398 );
399
400 gatheredNewParticles[Pstream::myProcNo()] = newParticles_;
401
402 // Gather data onto master
403 Pstream::gatherList(gatheredNewParticles);
404
405 // Combine
406 List<vectorPairScalarPair> combinedNewParticles
407 (
409 (
410 gatheredNewParticles,
412 )
413 );
414
415 if (Pstream::master())
416 {
417 newParticles_ = combinedNewParticles;
418 }
419
420 Pstream::scatter(newParticles_);
421 }
422
423 return newParticles_.size();
424}
425
426
427template<class CloudType>
429(
430 const scalar time0,
431 const scalar time1
432)
433{
434 if ((time0 >= 0.0) && (time0 < duration_))
435 {
436 return fraction_*flowRateProfile_->integrate(time0, time1);
437 }
438
439 return 0.0;
440}
441
442
443template<class CloudType>
445(
446 const label parcelI,
447 const label,
448 const scalar,
449 vector& position,
450 label& cellOwner,
451 label& tetFacei,
452 label& tetPti
453)
454{
455 position = newParticles_[parcelI].first().first();
456
457 this->findCellAtPosition
458 (
459 cellOwner,
460 tetFacei,
461 tetPti,
462 position,
463 false
464 );
465}
466
467
468template<class CloudType>
470(
471 const label parcelI,
472 const label,
473 const scalar,
474 typename CloudType::parcelType& parcel
475)
476{
477 parcel.U() = newParticles_[parcelI].first().second();
478
479 parcel.d() = newParticles_[parcelI].second().first();
480
481 parcel.dTarget() = newParticles_[parcelI].second().second();
482}
483
484
485template<class CloudType>
487{
488 return false;
489}
490
491
492template<class CloudType>
494{
495 return true;
496}
497
498
499// ************************************************************************* //
scalar minValue
#define R(A, B, C, D, E, F, K, M)
const List< DynamicList< molecule * > > & cellOccupancy
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:44
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:96
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition: HashTable.C:122
bool found(const Key &key) const
Return true if hashed entry is found in table.
Definition: HashTableI.H:100
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:500
Inflation injection - creates new particles by splitting existing particles within in a set of genera...
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
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.
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
virtual ~InflationInjection()
Destructor.
virtual void updateMesh()
Set injector locations when mesh is updated.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
scalar timeEnd() const
Return the end-of-injection time.
Templated injection model class.
scalar massTotal_
Total mass to inject [kg].
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:69
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
static void gatherList(const List< commsStruct > &comms, List< T > &values, const int tag, const label comm)
static void scatter(const List< commsStruct > &comms, T &value, const int tag, const label comm)
Broadcast data: Distribute without modification.
Random number generator.
Definition: Random.H:60
Type position(const Type &start, const Type &end)
Return a sample on the interval [start,end].
Type sample01()
Return a sample whose components lie in the range [0,1].
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (e.g. CA deg) to real-time (s).
Definition: TimeState.C:49
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: Tuple2.H:58
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
A collection of cell labels.
Definition: cellSet.H:54
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
const Time & time() const noexcept
Return time registry.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const scalarField & cellVolumes() const
const vectorField & cellCentres() const
int myProcNo() const noexcept
Return processor number.
Lookup type of boundary radiation properties.
Definition: lookup.H:66
splitCell * master() const
Definition: splitCell.H:113
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
Definition: subModelBase.C:131
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
#define WarningInFunction
Report a warning using Foam::Warning.
AccessType combine(const UList< T > &lists, AccessOp aop=accessOp< T >())
Combines sub-lists into a single list.
Definition: ListListOps.C:69
Mathematical constants.
const wordList volume
Standard volume field types (scalar, vector, tensor, etc)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Vector< scalar > vector
Definition: vector.H:61
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Random rndGen
Definition: createFields.H:23