InjectedParticleDistributionInjection.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) 2015-2020 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
30#include "bitSet.H"
32
33using namespace Foam::constant;
34
35// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36
37template<class CloudType>
39{
40 injectedParticleCloud ipCloud(this->owner().mesh(), cloudName_);
41
42 List<label> tag(ipCloud.size());
43 List<point> position(ipCloud.size());
44 List<vector> U(ipCloud.size());
45 List<scalar> soi(ipCloud.size());
46 List<scalar> d(ipCloud.size());
47
48 // Flatten all data
49 label particlei = 0;
50 for (const injectedParticle& p : ipCloud)
51 {
52 tag[particlei] = p.tag();
53 position[particlei] = p.position();
54 U[particlei] = p.U();
55 soi[particlei] = p.soi();
56 d[particlei] = p.d();
57 particlei++;
58 }
59
60 // Combine all proc data
61 if (Pstream::parRun())
62 {
63 List<List<label>> procTag(Pstream::nProcs());
64 procTag[Pstream::myProcNo()].transfer(tag);
65 Pstream::gatherList(procTag);
66 Pstream::scatterList(procTag);
67 tag =
68 ListListOps::combine<List<label>>
69 (
70 procTag, accessOp<List<label>>()
71 );
72
73 List<List<point>> procPosition(Pstream::nProcs());
74 procPosition[Pstream::myProcNo()].transfer(position);
75 Pstream::gatherList(procPosition);
76 Pstream::scatterList(procPosition);
77 position =
78 ListListOps::combine<List<point>>
79 (
80 procPosition, accessOp<List<point>>()
81 );
82
83 List<List<vector>> procU(Pstream::nProcs());
84 procU[Pstream::myProcNo()].transfer(U);
85 Pstream::gatherList(procU);
86 Pstream::scatterList(procU);
87 U =
88 ListListOps::combine<List<vector>>
89 (
90 procU, accessOp<List<vector>>()
91 );
92
93 List<List<scalar>> procSOI(Pstream::nProcs());
94 procSOI[Pstream::myProcNo()].transfer(soi);
95 Pstream::gatherList(procSOI);
96 Pstream::scatterList(procSOI);
97 soi =
98 ListListOps::combine<List<scalar>>
99 (
100 procSOI, accessOp<List<scalar>>()
101 );
102
103 List<List<scalar>> procD(Pstream::nProcs());
104 procD[Pstream::myProcNo()].transfer(d);
105 Pstream::gatherList(procD);
106 Pstream::scatterList(procD);
107 d =
108 ListListOps::combine<List<scalar>>
109 (
110 procD, accessOp<List<scalar>>()
111 );
112 }
113
114 label maxTag = -1;
115 forAll(tag, particlei)
116 {
117 maxTag = max(maxTag, tag[particlei]);
118 }
119
120 label nInjectors = maxTag + 1;
121 List<scalar> injStartTime(nInjectors, GREAT);
122 List<scalar> injEndTime(nInjectors, -GREAT);
123 List<DynamicList<point>> injPosition(nInjectors);
124 List<DynamicList<vector>> injU(nInjectors);
125 List<DynamicList<scalar>> injDiameter(nInjectors);
126
127 // Cache the particle information per tag
128 forAll(tag, i)
129 {
130 const label tagi = tag[i];
131 const scalar t = soi[i];
132 injStartTime[tagi] = min(t, injStartTime[tagi]);
133 injEndTime[tagi] = max(t, injEndTime[tagi]);
134 injPosition[tagi].append(position[i]);
135 injU[tagi].append(U[i]);
136 injDiameter[tagi].append(d[i]);
137 }
138
139 // Remove single particles and injectors where injection interval is 0
140 // - cannot generate a volume flow rate
141 scalar sumVolume = 0;
142 startTime_.setSize(nInjectors, 0);
143 endTime_.setSize(nInjectors, 0);
144 sizeDistribution_.setSize(nInjectors);
145 position_.setSize(nInjectors);
146 U_.setSize(nInjectors);
147 volumeFlowRate_.setSize(nInjectors, 0);
148
149 scalar minTime = GREAT;
150
151 // Populate injector properties, filtering out invalid entries
152 Random& rnd = this->owner().rndGen();
153 label injectori = 0;
154 forAll(injDiameter, i)
155 {
156 const DynamicList<scalar>& diameters = injDiameter[i];
157 const label nParticle = diameters.size();
158 const scalar dTime = injEndTime[i] - injStartTime[i];
159
160 if ((nParticle > 1) && (dTime > ROOTVSMALL))
161 {
162 minTime = min(minTime, injStartTime[i]);
163
164 startTime_[injectori] = injStartTime[i];
165 endTime_[injectori] = injEndTime[i];
166
167 // Re-sample the cloud data
168 position_[injectori].setSize(resampleSize_);
169 U_[injectori].setSize(resampleSize_);
170 List<point>& positioni = position_[injectori];
171 List<vector>& Ui = U_[injectori];
172
173 for (label samplei = 0; samplei < resampleSize_; ++samplei)
174 {
175 label posi = rnd.globalPosition<label>(0, nParticle - 1);
176 positioni[samplei] = injPosition[i][posi] + positionOffset_;
177 Ui[samplei] = injU[i][posi];
178 }
179
180 // Calculate the volume flow rate
181 scalar sumPow3 = 0;
182 forAll(diameters, particlei)
183 {
184 sumPow3 += pow3(diameters[particlei]);
185 }
186
187 const scalar volume = sumPow3*mathematical::pi/16.0;
188 sumVolume += volume;
189 volumeFlowRate_[injectori] = volume/dTime;
190
191 // Create the size distribution using the user-specified bin width
192 sizeDistribution_.set
193 (
194 injectori,
196 (
197 diameters,
198 binWidth_,
199 this->owner().rndGen()
200 )
201 );
202
203 injectori++;
204 }
205 }
206
207 // Resize
208 startTime_.setSize(injectori);
209 endTime_.setSize(injectori);
210 position_.setSize(injectori);
211 U_.setSize(injectori);
212 volumeFlowRate_.setSize(injectori);
213 sizeDistribution_.setSize(injectori);
214
215 // Reset start time to zero
216 forAll(startTime_, injectori)
217 {
218 startTime_[injectori] -= minTime;
219 endTime_[injectori] -= minTime;
220 }
221
222
223 // Set the volume of parcels to inject
224 this->volumeTotal_ = sumVolume;
225
226 // Provide some feedback
227 Info<< " Read " << position_.size() << " injectors with "
228 << tag.size() << " total particles" << endl;
229}
230
231
232// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
233
234template<class CloudType>
237(
238 const dictionary& dict,
239 CloudType& owner,
240 const word& modelName
241)
242:
243 InjectionModel<CloudType>(dict, owner, modelName, typeName),
244 cloudName_(this->coeffDict().lookup("cloud")),
245 startTime_(this->template getModelProperty<scalarList>("startTime")),
246 endTime_(this->template getModelProperty<scalarList>("endTime")),
247 position_(this->template getModelProperty<List<vectorList>>("position")),
248 positionOffset_(this->coeffDict().lookup("positionOffset")),
249 volumeFlowRate_
250 (
251 this->template getModelProperty<scalarList>("volumeFlowRate")
252 ),
253 U_(this->template getModelProperty<List<vectorList>>("U")),
254 binWidth_(this->coeffDict().getScalar("binWidth")),
255 sizeDistribution_(),
256 parcelsPerInjector_
257 (
258 ceil(this->coeffDict().getScalar("parcelsPerInjector"))
259 ),
260 resampleSize_
261 (
262 this->coeffDict().getOrDefault("resampleSize", label(100))
263 ),
264 applyDistributionMassTotal_
265 (
266 this->coeffDict().getBool("applyDistributionMassTotal")
267 ),
268 ignoreOutOfBounds_
269 (
270 this->coeffDict().getOrDefault("ignoreOutOfBounds", false)
271 ),
272 nParcelsInjected_(this->parcelsAddedTotal()),
273 nParcelsInjected0_(0),
274 currentInjectori_(0),
275 currentSamplei_(0)
276{
277 if (startTime_.size())
278 {
279 // Restart
282 {
283 const word dictName("distribution" + Foam::name(i));
285 this->getModelDict(dictName, dict);
286
288 (
289 i,
291 );
292 }
293 }
294 else
295 {
296 // Clean start
297 initialise();
298 }
299
300 // Set the mass of parcels to inject from distribution if requested
302 {
303 this->massTotal_ = this->volumeTotal_*this->owner().constProps().rho0();
304 Info<< " Set mass to inject from distribution: "
305 << this->massTotal_ << endl;
306 }
307}
308
309
310template<class CloudType>
313(
315)
316:
318 cloudName_(im.cloudName_),
319 startTime_(im.startTime_),
320 endTime_(im.endTime_),
321 position_(im.position_),
322 positionOffset_(im.positionOffset_),
323 volumeFlowRate_(im.volumeFlowRate_),
324 U_(im.U_),
325 binWidth_(im.binWidth_),
326 sizeDistribution_(im.sizeDistribution_.size()),
327 parcelsPerInjector_(im.parcelsPerInjector_),
328 resampleSize_(im.resampleSize_),
329 applyDistributionMassTotal_(im.applyDistributionMassTotal_),
330 ignoreOutOfBounds_(im.ignoreOutOfBounds_),
331 nParcelsInjected_(im.nParcelsInjected_),
332 nParcelsInjected0_(im.nParcelsInjected0_),
333 currentInjectori_(0),
334 currentSamplei_(0)
335{
336 forAll(sizeDistribution_, injectori)
337 {
338 if (sizeDistribution_.set(injectori))
339 {
341 (
342 injectori,
344 (
345 im.sizeDistribution_[injectori]
346 )
347 );
348 }
349 }
350}
351
352
353// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
354
355template<class CloudType>
358{}
359
360
361// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
362
363template<class CloudType>
365{}
366
367
368template<class CloudType>
369Foam::scalar
371{
372 return max(endTime_);
373}
374
375
376template<class CloudType>
377Foam::label
379(
380 const scalar time0,
381 const scalar time1
382)
383{
384 // Ensure all procs know the latest parcel count
385 nParcelsInjected_ += returnReduce(nParcelsInjected0_, sumOp<label>());
386 nParcelsInjected0_ = 0;
387
388 if (startTime_.empty() || this->volumeTotal_ < ROOTVSMALL)
389 {
390 return 0;
391 }
392
393 scalar targetVolume = 0;
394 forAll(startTime_, injectori)
395 {
396 if (time1 > startTime_[injectori])
397 {
398 scalar totalDuration =
399 min(time1, endTime_[injectori]) - startTime_[injectori];
400
401 targetVolume += volumeFlowRate_[injectori]*totalDuration;
402 }
403 }
404
405 const label targetParcels =
406 round
407 (
408 scalar(startTime_.size()*parcelsPerInjector_)
409 *targetVolume/this->volumeTotal_
410 );
411
412 const label nParcels = targetParcels - nParcelsInjected_;
413
414 return nParcels;
415}
416
417
418template<class CloudType>
419Foam::scalar
421(
422 const scalar time0,
423 const scalar time1
424)
425{
426 scalar volume = 0;
427 forAll(startTime_, injectori)
428 {
429 if ((time1 > startTime_[injectori]) && (time1 <= endTime_[injectori]))
430 {
431 scalar duration = min(time1, endTime_[injectori]) - time0;
432 volume += volumeFlowRate_[injectori]*duration;
433 }
434 }
435
436 return volume;
437}
438
439
440template<class CloudType>
442(
443 const label parcelI,
444 const label nParcels,
445 const scalar time,
446 vector& position,
447 label& cellOwner,
448 label& tetFaceI,
449 label& tetPtI
450)
451{
452 Random& rnd = this->owner().rndGen();
453 currentInjectori_ = rnd.globalPosition<label>(0, position_.size() - 1);
454 currentSamplei_ = rnd.globalPosition<label>(0, resampleSize_ - 1);
455
456 position = position_[currentInjectori_][currentSamplei_];
457
458 // Cache all mesh props for each position?
459 this->findCellAtPosition
460 (
461 cellOwner,
462 tetFaceI,
463 tetPtI,
464 position
465 );
466}
467
468
469template<class CloudType>
471(
472 const label parcelI,
473 const label,
474 const scalar,
475 typename CloudType::parcelType& parcel
476)
477{
478 // Set particle velocity
479 parcel.U() = U_[currentInjectori_][currentSamplei_];
480
481 // Set particle diameter
482 parcel.d() = sizeDistribution_[currentInjectori_].sample();
483
484 // Increment number of particles injected
485 // Note: local processor only!
486 nParcelsInjected0_++;
487}
488
489
490template<class CloudType>
491bool
493{
494 return false;
495}
496
497
498template<class CloudType>
500(
501 const label
502)
503{
504 return true;
505}
506
507
508template<class CloudType>
510{
512
513 if (this->writeTime())
514 {
515 this->setModelProperty("startTime", startTime_);
516 this->setModelProperty("endTime", endTime_);
517 this->setModelProperty("position", position_);
518 this->setModelProperty("volumeFlowRate", volumeFlowRate_);
519 this->setModelProperty("U", U_);
520 forAll(sizeDistribution_, i)
521 {
522 const distributionModels::general& dist = sizeDistribution_[i];
523 const word dictName("distribution" + Foam::name(i));
525 this->setModelProperty(dictName, dict);
526 }
527 }
528}
529
530
531// ************************************************************************* //
Y[inertIndex] max(0.0)
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
const List< typename ParcelType::constantProperties > & constProps() const
Return all of the constant properties.
Definition: DSMCCloudI.H:98
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
Interrogates an injectedParticleCloud to convert the raw particle data into a set of 'binned' injecto...
bool applyDistributionMassTotal_
Flag to apply mass calculated from distribution instead of.
PtrList< distributionModels::general > sizeDistribution_
List of size distribution model per injector.
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 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 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 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
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:330
void transfer(List< T > &list)
Definition: List.C:447
void setSize(const label n)
Alias for resize()
Definition: List.H:218
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Random number generator.
Definition: Random.H:60
Type globalPosition(const Type &start, const Type &end)
Return a sample on the interval [start,end].
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Particle-size distribution model wherein random samples are drawn from a given arbitrary probability ...
Definition: general.H:176
virtual dictionary writeDict(const word &dictName) const
Write data in dictionary format.
Definition: general.C:303
InfoProxy< ensightCells > info() const
Return info proxy.
Definition: ensightCells.H:254
Primarily stores particle properties so that it can be injected at a later time. Note that this store...
Lookup type of boundary radiation properties.
Definition: lookup.H:66
bool getModelDict(const word &entryName, dictionary &dict) const
Retrieve dictionary, return true if set.
Definition: subModelBase.C:185
const dictionary & dict() const
Return const access to the cloud dictionary.
Definition: subModelBase.C:113
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
const word dictName("faMeshDefinition")
constexpr scalar pi(M_PI)
Different types of constants.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedScalar pow3(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Random rndGen
Definition: createFields.H:23