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-2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 
29 #include "mathematicalConstants.H"
30 #include "bitSet.H"
32 
33 using namespace Foam::constant;
34 
35 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36 
37 template<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  {
64  procTag[Pstream::myProcNo()].transfer(tag);
65  Pstream::allGatherList(procTag);
66  tag =
67  ListListOps::combine<List<label>>
68  (
69  procTag, accessOp<List<label>>()
70  );
71 
72  List<List<point>> procPosition(Pstream::nProcs());
73  procPosition[Pstream::myProcNo()].transfer(position);
74  Pstream::allGatherList(procPosition);
75  position =
76  ListListOps::combine<List<point>>
77  (
78  procPosition, accessOp<List<point>>()
79  );
80 
82  procU[Pstream::myProcNo()].transfer(U);
84  U =
85  ListListOps::combine<List<vector>>
86  (
87  procU, accessOp<List<vector>>()
88  );
89 
91  procSOI[Pstream::myProcNo()].transfer(soi);
92  Pstream::allGatherList(procSOI);
93  soi =
94  ListListOps::combine<List<scalar>>
95  (
96  procSOI, accessOp<List<scalar>>()
97  );
98 
100  procD[Pstream::myProcNo()].transfer(d);
101  Pstream::allGatherList(procD);
102  d =
103  ListListOps::combine<List<scalar>>
104  (
105  procD, accessOp<List<scalar>>()
106  );
107  }
108 
109  label maxTag = -1;
110  forAll(tag, particlei)
111  {
112  maxTag = max(maxTag, tag[particlei]);
113  }
114 
115  label nInjectors = maxTag + 1;
116  List<scalar> injStartTime(nInjectors, GREAT);
117  List<scalar> injEndTime(nInjectors, -GREAT);
118  List<DynamicList<point>> injPosition(nInjectors);
119  List<DynamicList<vector>> injU(nInjectors);
120  List<DynamicList<scalar>> injDiameter(nInjectors);
121 
122  // Cache the particle information per tag
123  forAll(tag, i)
124  {
125  const label tagi = tag[i];
126  const scalar t = soi[i];
127  injStartTime[tagi] = min(t, injStartTime[tagi]);
128  injEndTime[tagi] = max(t, injEndTime[tagi]);
129  injPosition[tagi].append(position[i]);
130  injU[tagi].append(U[i]);
131  injDiameter[tagi].append(d[i]);
132  }
133 
134  // Remove single particles and injectors where injection interval is 0
135  // - cannot generate a volume flow rate
136  scalar sumVolume = 0;
137  startTime_.setSize(nInjectors, 0);
138  endTime_.setSize(nInjectors, 0);
139  sizeDistribution_.setSize(nInjectors);
140  position_.setSize(nInjectors);
141  U_.setSize(nInjectors);
142  volumeFlowRate_.setSize(nInjectors, 0);
143 
144  scalar minTime = GREAT;
145 
146  // Populate injector properties, filtering out invalid entries
147  Random& rnd = this->owner().rndGen();
148  label injectori = 0;
149  forAll(injDiameter, i)
150  {
151  const DynamicList<scalar>& diameters = injDiameter[i];
152  const label nParticle = diameters.size();
153  const scalar dTime = injEndTime[i] - injStartTime[i];
154 
155  if ((nParticle > 1) && (dTime > ROOTVSMALL))
156  {
157  minTime = min(minTime, injStartTime[i]);
158 
159  startTime_[injectori] = injStartTime[i];
160  endTime_[injectori] = injEndTime[i];
161 
162  // Re-sample the cloud data
163  position_[injectori].setSize(resampleSize_);
164  U_[injectori].setSize(resampleSize_);
165  List<point>& positioni = position_[injectori];
166  List<vector>& Ui = U_[injectori];
167 
168  for (label samplei = 0; samplei < resampleSize_; ++samplei)
169  {
170  label posi = rnd.globalPosition<label>(0, nParticle - 1);
171  positioni[samplei] = injPosition[i][posi] + positionOffset_;
172  Ui[samplei] = injU[i][posi];
173  }
174 
175  // Calculate the volume flow rate
176  scalar sumPow3 = 0;
177  forAll(diameters, particlei)
178  {
179  sumPow3 += pow3(diameters[particlei]);
180  }
181 
182  const scalar volume = sumPow3*mathematical::pi/6.0;
183  sumVolume += volume;
184  volumeFlowRate_[injectori] = volume/dTime;
185 
186  // Create the size distribution using the user-specified bin width
187  sizeDistribution_.set
188  (
189  injectori,
191  (
192  diameters,
193  binWidth_,
194  this->owner().rndGen()
195  )
196  );
197 
198  injectori++;
199  }
200  }
201 
202  // Resize
203  startTime_.setSize(injectori);
204  endTime_.setSize(injectori);
205  position_.setSize(injectori);
206  U_.setSize(injectori);
207  volumeFlowRate_.setSize(injectori);
208  sizeDistribution_.setSize(injectori);
209 
210  // Reset start time to zero
211  forAll(startTime_, injectori)
212  {
213  startTime_[injectori] -= minTime;
214  endTime_[injectori] -= minTime;
215  }
216 
217 
218  // Set the volume of parcels to inject
219  this->volumeTotal_ = sumVolume;
220 
221  // Provide some feedback
222  Info<< " Read " << position_.size() << " injectors with "
223  << tag.size() << " total particles" << endl;
224 }
226 
227 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
228 
229 template<class CloudType>
232 (
233  const dictionary& dict,
234  CloudType& owner,
235  const word& modelName
236 )
237 :
238  InjectionModel<CloudType>(dict, owner, modelName, typeName),
239  cloudName_(this->coeffDict().lookup("cloud")),
240  startTime_(this->template getModelProperty<scalarList>("startTime")),
241  endTime_(this->template getModelProperty<scalarList>("endTime")),
242  position_(this->template getModelProperty<List<vectorList>>("position")),
243  positionOffset_(this->coeffDict().lookup("positionOffset")),
244  volumeFlowRate_
245  (
246  this->template getModelProperty<scalarList>("volumeFlowRate")
247  ),
248  U_(this->template getModelProperty<List<vectorList>>("U")),
249  binWidth_(this->coeffDict().getScalar("binWidth")),
250  sizeDistribution_(),
251  parcelsPerInjector_
252  (
253  ceil(this->coeffDict().getScalar("parcelsPerInjector"))
254  ),
255  resampleSize_
256  (
257  this->coeffDict().getOrDefault("resampleSize", label(100))
258  ),
259  applyDistributionMassTotal_
260  (
261  this->coeffDict().getBool("applyDistributionMassTotal")
262  ),
263  ignoreOutOfBounds_
264  (
265  this->coeffDict().getOrDefault("ignoreOutOfBounds", false)
266  ),
267  nParcelsInjected_(this->parcelsAddedTotal()),
268  nParcelsInjected0_(0),
269  currentInjectori_(0),
270  currentSamplei_(0)
271 {
272  if (startTime_.size())
273  {
274  // Restart
275  sizeDistribution_.setSize(startTime_.size());
277  {
278  const word dictName("distribution" + Foam::name(i));
280  this->getModelDict(dictName, dict);
281 
283  (
284  i,
286  );
287  }
288  }
289  else
290  {
291  // Clean start
292  initialise();
293  }
294 
295  // Set the mass of parcels to inject from distribution if requested
297  {
298  this->massTotal_ = this->volumeTotal_*this->owner().constProps().rho0();
299  Info<< " Set mass to inject from distribution: "
300  << this->massTotal_ << endl;
301  }
302 }
303 
304 
305 template<class CloudType>
308 (
309  const InjectedParticleDistributionInjection<CloudType>& im
310 )
311 :
312  InjectionModel<CloudType>(im),
313  cloudName_(im.cloudName_),
314  startTime_(im.startTime_),
315  endTime_(im.endTime_),
316  position_(im.position_),
317  positionOffset_(im.positionOffset_),
318  volumeFlowRate_(im.volumeFlowRate_),
319  U_(im.U_),
320  binWidth_(im.binWidth_),
321  sizeDistribution_(im.sizeDistribution_.size()),
322  parcelsPerInjector_(im.parcelsPerInjector_),
323  resampleSize_(im.resampleSize_),
324  applyDistributionMassTotal_(im.applyDistributionMassTotal_),
325  ignoreOutOfBounds_(im.ignoreOutOfBounds_),
326  nParcelsInjected_(im.nParcelsInjected_),
327  nParcelsInjected0_(im.nParcelsInjected0_),
328  currentInjectori_(0),
329  currentSamplei_(0)
330 {
331  forAll(sizeDistribution_, injectori)
332  {
333  if (sizeDistribution_.set(injectori))
334  {
336  (
337  injectori,
338  new distributionModels::general
339  (
340  im.sizeDistribution_[injectori]
341  )
342  );
343  }
344  }
345 }
346 
347 
348 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
349 
350 template<class CloudType>
353 {}
354 
355 
356 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
357 
358 template<class CloudType>
360 {}
361 
362 
363 template<class CloudType>
364 Foam::scalar
366 {
367  return max(endTime_);
368 }
369 
370 
371 template<class CloudType>
372 Foam::label
374 (
375  const scalar time0,
376  const scalar time1
377 )
378 {
379  // Ensure all procs know the latest parcel count
380  nParcelsInjected_ += returnReduce(nParcelsInjected0_, sumOp<label>());
381  nParcelsInjected0_ = 0;
382 
383  if (startTime_.empty() || this->volumeTotal_ < ROOTVSMALL)
384  {
385  return 0;
386  }
387 
388  scalar targetVolume = 0;
389  forAll(startTime_, injectori)
390  {
391  if (time1 > startTime_[injectori])
392  {
393  scalar totalDuration =
394  min(time1, endTime_[injectori]) - startTime_[injectori];
395 
396  targetVolume += volumeFlowRate_[injectori]*totalDuration;
397  }
398  }
399 
400  const label targetParcels =
401  round
402  (
403  scalar(startTime_.size()*parcelsPerInjector_)
404  *targetVolume/this->volumeTotal_
405  );
406 
407  const label nParcels = targetParcels - nParcelsInjected_;
408 
409  return nParcels;
410 }
411 
412 
413 template<class CloudType>
414 Foam::scalar
416 (
417  const scalar time0,
418  const scalar time1
419 )
420 {
421  scalar volume = 0;
422  forAll(startTime_, injectori)
423  {
424  if ((time1 > startTime_[injectori]) && (time1 <= endTime_[injectori]))
425  {
426  scalar duration = min(time1, endTime_[injectori]) - time0;
427  volume += volumeFlowRate_[injectori]*duration;
428  }
429  }
431  return volume;
432 }
433 
434 
435 template<class CloudType>
437 (
438  const label parcelI,
439  const label nParcels,
440  const scalar time,
441  vector& position,
442  label& cellOwner,
443  label& tetFaceI,
444  label& tetPtI
445 )
446 {
447  Random& rnd = this->owner().rndGen();
448  currentInjectori_ = rnd.globalPosition<label>(0, position_.size() - 1);
449  currentSamplei_ = rnd.globalPosition<label>(0, resampleSize_ - 1);
450 
451  position = position_[currentInjectori_][currentSamplei_];
452 
453  // Cache all mesh props for each position?
454  this->findCellAtPosition
455  (
456  cellOwner,
457  tetFaceI,
458  tetPtI,
459  position
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  // Set particle velocity
474  parcel.U() = U_[currentInjectori_][currentSamplei_];
475 
476  // Set particle diameter
477  parcel.d() = sizeDistribution_[currentInjectori_].sample();
478 
479  // Increment number of particles injected
480  // Note: local processor only!
481  nParcelsInjected0_++;
482 }
483 
484 
485 template<class CloudType>
486 bool
488 {
489  return false;
490 }
491 
492 
493 template<class CloudType>
495 (
496  const label
497 )
498 {
499  return true;
500 }
501 
502 
503 template<class CloudType>
505 {
507 
508  if (this->writeTime())
509  {
510  this->setModelProperty("startTime", startTime_);
511  this->setModelProperty("endTime", endTime_);
512  this->setModelProperty("position", position_);
513  this->setModelProperty("volumeFlowRate", volumeFlowRate_);
514  this->setModelProperty("U", U_);
515  forAll(sizeDistribution_, i)
516  {
517  const distributionModels::general& dist = sizeDistribution_[i];
518  const word dictName("distribution" + Foam::name(i));
520  this->setModelProperty(dictName, dict);
521  }
522  }
523 }
524 
525 
526 // ************************************************************************* //
Different types of constants.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
scalar massTotal_
Total mass to inject [kg].
DSMCCloud< dsmcParcel > CloudType
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:517
Templated injection model class.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
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:489
Random rndGen
Definition: createFields.H:23
const word dictName("faMeshDefinition")
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
Object access operator or list access operator (default is pass-through)
Definition: UList.H:1004
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
Primarily stores particle properties so that it can be injected at a later time. Note that this store...
scalar timeEnd() const
Return the end-of-injection time.
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
virtual void updateMesh()
Set injector locations when mesh is updated.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
Lookup type of boundary radiation properties.
Definition: lookup.H:57
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
scalarList startTime_
List of start time per injector.
const dictionary & dict() const
Return const access to the cloud dictionary.
Definition: subModelBase.C:104
const CloudType & owner() const
Return const access to the owner cloud.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static void allGatherList(List< T > &values, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Gather data, but keep individual values separate. Uses linear/tree communication. ...
label size() const noexcept
The number of elements in table.
Definition: HashTable.H:342
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1065
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
bool applyDistributionMassTotal_
Flag to apply mass calculated from distribution instead of.
A class for handling words, derived from Foam::string.
Definition: word.H:63
bool getModelDict(const word &entryName, dictionary &dict) const
Retrieve dictionary, return true if set.
Definition: subModelBase.C:177
constexpr scalar pi(M_PI)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
Random number generator.
Definition: Random.H:55
virtual void info()
Write injection info.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:290
const wordList volume
Standard volume field types (scalar, vector, tensor, etc)
Type globalPosition(const Type &start, const Type &end)
Return a sample on the interval [start,end].
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
PtrList< distributionModels::general > sizeDistribution_
List of size distribution model per injector.
dimensionedScalar pow3(const dimensionedScalar &ds)
U
Definition: pEqn.H:72
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Particle-size distribution model wherein random samples are drawn from a given arbitrary probability ...
Definition: general.H:165
scalar volumeTotal_
Total volume of particles introduced by this injector [m^3] Note: scaled to ensure massTotal is achie...
virtual dictionary writeDict(const word &dictName) const
Write data in dictionary format.
Definition: general.C:296
volScalarField & p
InjectedParticleDistributionInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
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.
const List< typename ParcelType::constantProperties > & constProps() const
Return all of the constant properties.
Definition: DSMCCloudI.H:91