InjectionModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
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 "InjectionModel.H"
27 #include "meshTools.H"
28 #include "volFields.H"
29 #include "Scale.H"
30 
31 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
32 
33 template<class CloudType>
35 (
36  const dictionary& dict,
37  CloudType& owner
38 )
39 {
40  if (dict.found("nParticle"))
41  {
42  if (dict.found("massTotal"))
43  {
45  << "If nParticle is specified then the massTotal "
46  << "setting has no effect " << endl;
47  }
48 
49  return NaN;
50  }
51 
52  if (owner.solution().steadyState())
53  {
55  << "The " << type() << " injection model is not compatible with "
56  << "steady state solution"
57  << exit(FatalError);
58 
59  return NaN;
60  }
61 
62  return dict.lookup<scalar>("massTotal", dimMass);
63 }
64 
65 
66 template<class CloudType>
68 (
69  const dictionary& dict,
70  CloudType& owner
71 )
72 {
73  if (owner.solution().steadyState())
74  {
75  return vGreat;
76  }
77 
78  return dict.lookup<scalar>("duration", owner.db().time().userUnits());
79 }
80 
81 
82 template<class CloudType>
85 (
86  const dictionary& dict,
87  CloudType& owner,
88  const scalar duration
89 )
90 {
91  const bool haveMassFlowRate = dict.found("massFlowRate");
92  const bool haveMassTotal = dict.found("massTotal");
93 
94  if (dict.found("nParticle"))
95  {
96  if (haveMassFlowRate || haveMassTotal)
97  {
99  << "If nParticle is specified then massFlowRate and massTotal "
100  << "settings have no effect " << endl;
101  }
102 
103  return
105  (
106  new Function1s::Constant<scalar>("NaN", NaN)
107  );
108  }
109 
110  if (owner.solution().steadyState() && haveMassTotal)
111  {
113  << "Cannot specify the massTotal of a steady injection. Use "
114  << "massFlowRate instead." << exit(FatalIOError);
115  }
116 
117  if (haveMassFlowRate && haveMassTotal)
118  {
120  << "Cannot specify both massFlowRate and massTotal. Use one or "
121  << "the other." << exit(FatalIOError);
122  }
123 
124  if (owner.solution().steadyState() || haveMassFlowRate)
125  {
126  return
128  (
129  "massFlowRate",
130  this->owner().db().time().userUnits(),
132  dict
133  );
134  }
135 
136  const scalar massTotal = dict.lookup<scalar>("massTotal", dimMass);
137 
138  if (!dict.found("flowRateProfile"))
139  {
140  return
142  (
144  (
145  "massFlowRate",
146  massTotal/duration
147  )
148  );
149  }
150 
151  autoPtr<Function1<scalar>> flowRateProfile
152  (
154  (
155  "flowRateProfile",
156  this->owner().db().time().userUnits(),
157  dimless,
158  dict
159  )
160  );
161 
162  const scalar sumFlowRateProfile = flowRateProfile->integral(0, duration);
163 
164  return
166  (
168  (
169  "massFlowRate",
170  Function1s::Constant<scalar>("m", massTotal/sumFlowRateProfile),
171  Function1s::Constant<scalar>("one", scalar(1)),
172  flowRateProfile()
173  )
174  );
175 }
176 
177 
178 template<class CloudType>
181 (
182  const dictionary& dict,
183  CloudType& owner
184 )
185 {
186  return
188  (
189  "parcelsPerSecond",
190  this->owner().db().time().userUnits(),
192  dict
193  );
194 }
195 
196 
197 template<class CloudType>
199 {
200  forAll(this->owner().injectors(), i)
201  {
202  if (this->owner().injectors()(i) == this)
203  {
204  return i;
205  }
206  }
207 
208  return -1;
209 }
210 
211 
212 template<class CloudType>
214 (
215  const point& position,
216  barycentric& coordinates,
217  label& celli,
218  label& tetFacei,
219  label& tetPti,
220  bool errorOnNotFound
221 )
222 {
223  // Subroutine for finding the cell
224  auto findProcAndCell = [this](const point& pos)
225  {
226  // Find the containing cell
227  label celli = this->owner().mesh().findCell(pos);
228 
229  // Synchronise so only a single processor finds this position
230  label proci = celli >= 0 ? Pstream::myProcNo() : -1;
231  reduce(proci, maxOp<label>());
232  if (proci != Pstream::myProcNo())
233  {
234  celli = -1;
235  }
236 
237  return labelPair(proci, celli);
238  };
239 
240  point pos = position;
241 
242  // Try and find the cell at the given position
243  const labelPair procAndCelli = findProcAndCell(pos);
244  label proci = procAndCelli.first();
245  celli = procAndCelli.second();
246 
247  // Didn't find it. The point may be awkwardly on an edge or face. Try
248  // again, but move the point into its nearest cell a little bit.
249  if (proci == -1)
250  {
251  pos += small*(this->owner().mesh().C()[celli] - pos);
252  const labelPair procAndCelli = findProcAndCell(pos);
253  proci = procAndCelli.first();
254  celli = procAndCelli.second();
255  }
256 
257  // Didn't find it. Error or return false.
258  if (proci == -1)
259  {
260  if (errorOnNotFound)
261  {
263  << "Cannot find parcel injection cell. "
264  << "Parcel position = " << position << nl
265  << exit(FatalError);
266  }
267 
268  return false;
269  }
270 
271  // Found it. Construct the barycentric coordinates.
272  if (proci == Pstream::myProcNo())
273  {
274  label nLocateBoundaryHits = 0;
275  particle p(this->owner().mesh(), pos, celli, nLocateBoundaryHits);
276 
277  if (nLocateBoundaryHits != 0)
278  {
280  << "Injection model " << this->modelName()
281  << " for cloud " << this->owner().name()
282  << " did not accurately locate the position "
283  << pos << " within the mesh" << endl;
284  }
285 
286  coordinates = p.coordinates();
287  celli = p.cell();
288  tetFacei = p.tetFace();
289  tetPti = p.tetPt();
290  }
291 
292  return true;
293 }
294 
295 
296 template<class CloudType>
298 (
299  typename CloudType::parcelType::trackingData& td,
300  typename CloudType::parcelType& parcel
301 )
302 {
303  const vector d = parcel.deviationFromMeshCentre(td.mesh);
304 
305  if (d == vector::zero)
306  {
307  return;
308  }
309 
310  const label facei = parcel.face();
311 
312  // If the parcel is not on a face, then just track it to the mesh centre
313  if (facei == -1)
314  {
315  parcel.track(td.mesh, - d, 0);
316  }
317 
318  // If the parcel is on a face, then track in two steps, going slightly into
319  // the current cell. This prevents a boundary hit from ending the track
320  // prematurely.
321  if (facei != -1)
322  {
323  const vector pc =
324  td.mesh.cellCentres()[parcel.cell()] - parcel.position(td.mesh);
325 
326  parcel.track(td.mesh, - d/2 + rootSmall*pc, 0);
327  parcel.track(td.mesh, - d/2 - rootSmall*pc, 0);
328  }
329 
330  // Restore any face-association that got changed during tracking
331  parcel.face() = facei;
332 }
333 
334 
335 template<class CloudType>
337 {
338  switch (uniformParcelSize_)
339  {
340  case uniformParcelSize::nParticle:
341  return 0;
342  case uniformParcelSize::surfaceArea:
343  return 2;
344  case uniformParcelSize::volume:
345  return 3;
346  }
347 
348  return -labelMax;
349 }
350 
351 
352 template<class CloudType>
354 (
355  PtrList<parcelType>& parcelPtrs,
356  const scalar mass
357 ) const
358 {
359  auto size = [&](const parcelType& p)
360  {
361  switch (uniformParcelSize_)
362  {
363  case uniformParcelSize::nParticle:
364  return scalar(1);
365  case uniformParcelSize::surfaceArea:
366  return p.areaS();
367  case uniformParcelSize::volume:
368  return p.volume();
369  }
370  return NaN;
371  };
372 
373  // Determine the total mass and size of all created particles
374  scalar sumMassBySize = 0;
375  forAll(parcelPtrs, parceli)
376  {
377  if (parcelPtrs.set(parceli))
378  {
379  const parcelType& p = parcelPtrs[parceli];
380  sumMassBySize += p.mass()/size(p);
381  }
382  }
383 
384  reduce(sumMassBySize, sumOp<scalar>());
385 
386  // Set the numbers of particles on each parcel
387  forAll(parcelPtrs, parceli)
388  {
389  if (parcelPtrs.set(parceli))
390  {
391  parcelType& p = parcelPtrs[parceli];
392  p.nParticle() = mass/size(p)/sumMassBySize;
393  }
394  }
395 
396  // Check that the constraints are correct
397  if (debug)
398  {
399  scalar massN = 0, minSizeN = vGreat, maxSizeN = -vGreat;
400  forAll(parcelPtrs, parceli)
401  {
402  if (parcelPtrs.set(parceli))
403  {
404  const parcelType& p = parcelPtrs[parceli];
405  massN += p.nParticle()*p.mass();
406  minSizeN = min(minSizeN, p.nParticle()*size(p));
407  maxSizeN = max(minSizeN, p.nParticle()*size(p));
408  }
409  }
410 
411  reduce(massN, sumOp<scalar>());
412 
413  if (mag(massN - mass) > rootSmall*(massN + mass)/2)
414  {
416  << "Parcels do not have the required mass"
417  << exit(FatalError);
418  }
419 
420  reduce(minSizeN, minOp<scalar>());
421  reduce(maxSizeN, maxOp<scalar>());
422 
423  if (maxSizeN - minSizeN > rootSmall*(maxSizeN + minSizeN)/2)
424  {
426  << "Parcel sizes are not uniform"
427  << exit(FatalError);
428  }
429  }
430 }
431 
432 
433 template<class CloudType>
435 (
436  typename parcelType::trackingData& td
437 )
438 {}
439 
440 
441 template<class CloudType>
443 (
444  const label nParcelsAdded,
445  const scalar massAdded,
446  typename parcelType::trackingData& td
447 )
448 {
449  const label allNParcelsAdded = returnReduce(nParcelsAdded, sumOp<label>());
450 
451  if (allNParcelsAdded > 0)
452  {
453  Info<< nl
454  << "Cloud: " << this->owner().name()
455  << " injector: " << this->modelName() << nl
456  << " Added " << allNParcelsAdded << " new parcels" << nl << endl;
457  }
458 
459  // Increment total number of parcels added
460  parcelsAddedTotal_ += allNParcelsAdded;
461 
462  // Increment total mass injected
463  massInjected_ += returnReduce(massAdded, sumOp<scalar>());
464 
465  // Update time for start of next injection
466  time0_ = this->owner().db().time().value();
467 
468  // Increment number of injections
469  nInjections_++;
470 }
471 
472 
473 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
474 
475 template<class CloudType>
477 :
479  SOI_(0),
480  massInjected_(this->template getModelProperty<scalar>("massInjected")),
481  nInjections_(this->template getModelProperty<label>("nInjections")),
482  parcelsAddedTotal_
483  (
484  this->template getModelProperty<scalar>("parcelsAddedTotal")
485  ),
486  nParticleFixed_(-vGreat),
487  uniformParcelSize_(uniformParcelSize::nParticle),
488  time0_(0),
489  timeStep0_(this->template getModelProperty<scalar>("timeStep0"))
490 {}
491 
492 
493 template<class CloudType>
495 (
496  const dictionary& dict,
497  CloudType& owner,
498  const word& modelName,
499  const word& modelType
500 )
501 :
502  CloudSubModelBase<CloudType>(modelName, owner, dict, typeName, modelType),
503  SOI_(0),
504  massInjected_(this->template getModelProperty<scalar>("massInjected")),
505  nInjections_(this->template getModelProperty<scalar>("nInjections")),
506  parcelsAddedTotal_
507  (
508  this->template getModelProperty<scalar>("parcelsAddedTotal")
509  ),
510  nParticleFixed_(dict.lookupOrDefault<scalar>("nParticle", -vGreat)),
511  uniformParcelSize_
512  (
513  uniformParcelSizeNames_
514  [
515  !dict.found("parcelBasisType") && nParticleFixed_ > 0
516  ? dict.lookupOrDefault<word>
517  (
518  "uniformParcelSize",
519  uniformParcelSizeNames_[uniformParcelSize::nParticle]
520  )
521  : dict.lookup<word>("uniformParcelSize")
522  ]
523  ),
524  time0_(owner.db().time().value()),
525  timeStep0_(this->template getModelProperty<scalar>("timeStep0"))
526 {
527  // Provide some info. Also serves to initialise mesh dimensions. This may
528  // be needed for parallel runs due to lazy evaluation of valid mesh
529  // dimensions.
530  Info<< " Constructing " << owner.mesh().nGeometricD() << "-D injection"
531  << endl;
532 
533  if
534  (
535  nParticleFixed_ > 0
537  )
538  {
540  << "If nParticle is specified then the uniformParcelSize must be "
542  << exit(FatalIOError);
543  }
544 
545  if (owner.solution().transient())
546  {
547  SOI_ = dict.lookup<scalar>("SOI", owner.db().time().userUnits());
548  }
549 }
550 
551 
552 template<class CloudType>
554 (
555  const InjectionModel<CloudType>& im
556 )
557 :
559  SOI_(im.SOI_),
560  massInjected_(im.massInjected_),
561  nInjections_(im.nInjections_),
562  parcelsAddedTotal_(im.parcelsAddedTotal_),
563  nParticleFixed_(im.nParticleFixed_),
564  uniformParcelSize_(im.uniformParcelSize_),
565  time0_(im.time0_),
566  timeStep0_(im.timeStep0_)
567 {}
568 
569 
570 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
571 
572 template<class CloudType>
574 {}
575 
576 
577 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
578 
579 template<class CloudType>
581 {}
582 
583 
584 template<class CloudType>
586 {
587  const scalar deltaT =
588  this->owner().solution().transient() ? timeEnd() - timeStart() : 1;
589 
590  return massToInject(0, deltaT)/nParcelsToInject(0, deltaT);
591 }
592 
593 
594 template<class CloudType>
595 template<class TrackCloudType>
597 (
598  TrackCloudType& cloud,
599  typename CloudType::parcelType::trackingData& td
600 )
601 {
602  const polyMesh& mesh = this->owner().mesh();
603 
604  const scalar time = this->owner().db().time().value();
605 
606  preInject(td);
607 
608  // Reset counters
609  label nParcelsAdded = 0;
610  scalar massAdded = 0;
611 
612  // Get amounts to inject
613  label nParcels;
614  scalar mass;
615  bool inject;
616  if (time < SOI_)
617  {
618  // Injection has not started yet
619  nParcels = 0;
620  mass = 0;
621  inject = false;
622  timeStep0_ = time;
623  }
624  else
625  {
626  // Injection has started. Get amounts between times.
627  const scalar t0 = timeStep0_ - SOI_, t1 = time - SOI_;
628  nParcels = nParcelsToInject(t0, t1);
629  mass = nParticleFixed_ < 0 ? massToInject(t0, t1) : NaN;
630 
631  if (nParcels > 0 && (nParticleFixed_ > 0 || mass > 0))
632  {
633  // Injection is valid
634  inject = true;
635  timeStep0_ = time;
636  }
637  else if (nParcels == 0 && (nParticleFixed_ < 0 && mass > 0))
638  {
639  // Injection should have started, but there is not a sufficient
640  // amount to inject a single parcel. Do not inject, but hold
641  // advancement of the old-time so that the mass gets added to a
642  // subsequent injection.
643  inject = false;
644  }
645  else
646  {
647  // Nothing to inject
648  inject = false;
649  timeStep0_ = time;
650  }
651  }
652 
653  // Do injection
654  if (inject)
655  {
656  // Duration of injection period during this timestep
657  const scalar deltaT =
658  max
659  (
660  scalar(0),
661  min
662  (
663  td.trackTime(),
664  min
665  (
666  time - SOI_,
667  timeEnd() - time0_
668  )
669  )
670  );
671 
672  // Pad injection time if injection starts during this timestep
673  const scalar padTime = max(scalar(0), SOI_ - time0_);
674 
675  // Create new parcels linearly across carrier phase timestep
676  PtrList<parcelType> parcelPtrs(nParcels);
677  forAll(parcelPtrs, parceli)
678  {
679  // Calculate the pseudo time of injection for parcel 'parceli'
680  scalar timeInj = time0_ + padTime + deltaT*parceli/nParcels;
681 
682  // Determine the injection coordinates and owner cell,
683  // tetFace and tetPt
684  barycentric coordinates = barycentric::uniform(NaN);
685  label celli = -1, tetFacei = -1, tetPti = -1, facei = -1;
686  setPositionAndCell
687  (
688  parceli,
689  nParcels,
690  timeInj,
691  coordinates,
692  celli,
693  tetFacei,
694  tetPti,
695  facei
696  );
697 
698  if (celli > -1)
699  {
700  // Lagrangian timestep
701  const scalar dt = timeInj - time0_;
702 
703  // Create a new parcel
704  parcelPtrs.set
705  (
706  parceli,
707  new parcelType
708  (
709  mesh,
710  coordinates,
711  celli,
712  tetFacei,
713  tetPti,
714  facei
715  )
716  );
717  parcelType& p = parcelPtrs[parceli];
718 
719  // Correct the position for reduced-dimension cases
720  constrainPosition(td, p);
721 
722  // Check/set new parcel thermo properties
723  cloud.setParcelThermoProperties(p);
724 
725  // Assign new parcel properties in injection model
726  setProperties(parceli, nParcels, timeInj, td, p);
727 
728  // Check/set new parcel injection properties
729  cloud.checkParcelProperties(p, index());
730 
731  // Apply correction to velocity for 2-D cases
733  (
734  mesh,
735  mesh.solutionD(),
736  p.U()
737  );
738 
739  // Modify the step fraction so that the particles are
740  // injected continually through the time-step
741  p.stepFraction() = dt/td.trackTime();
742 
743  // Set the number of particles. If not fixed, this will set
744  // a junk value, which will get corrected below.
745  p.nParticle() = nParticleFixed_;
746  }
747  }
748 
749  // Set the number of particles so that the introduced mass is correct
750  // and the uniform size is as specified
751  if (nParticleFixed_ < 0)
752  {
753  setNumberOfParticles(parcelPtrs, mass);
754  }
755 
756  // Add the new parcels
757  forAll(parcelPtrs, parceli)
758  {
759  if (parcelPtrs.set(parceli))
760  {
761  parcelType& p = parcelPtrs[parceli];
762  nParcelsAdded ++;
763  massAdded += p.nParticle()*p.mass();
764  cloud.addParticle(parcelPtrs.set(parceli, nullptr).ptr());
765  }
766  }
767  }
768 
769  postInject(nParcelsAdded, massAdded, td);
770 }
771 
772 
773 template<class CloudType>
774 template<class TrackCloudType>
776 (
777  TrackCloudType& cloud,
778  typename CloudType::parcelType::trackingData& td
779 )
780 {
781  const polyMesh& mesh = this->owner().mesh();
782 
783  preInject(td);
784 
785  // Reset counters
786  label nParcelsAdded = 0;
787  scalar massAdded = 0;
788 
789  // Get amounts to inject based on first second of injection
790  const label nParcels = nParcelsToInject(0, 1);
791  const scalar mass = nParticleFixed_ < 0 ? massToInject(0, 1) : NaN;
792 
793  // Do injection
794  if (nParcels > 0)
795  {
796  PtrList<parcelType> parcelPtrs(nParcels);
797  forAll(parcelPtrs, parceli)
798  {
799  // Determine the injection coordinates and owner cell,
800  // tetFace and tetPt
801  barycentric coordinates = barycentric::uniform(NaN);
802  label celli = -1, tetFacei = -1, tetPti = -1, facei = -1;
803  setPositionAndCell
804  (
805  parceli,
806  nParcels,
807  0,
808  coordinates,
809  celli,
810  tetFacei,
811  tetPti,
812  facei
813  );
814 
815  if (celli > -1)
816  {
817  // Create a new parcel
818  parcelPtrs.set
819  (
820  parceli,
821  new parcelType
822  (
823  mesh,
824  coordinates,
825  celli,
826  tetFacei,
827  tetPti,
828  facei
829  )
830  );
831  parcelType& p = parcelPtrs[parceli];
832 
833  // Correct the position for reduced-dimension cases
834  constrainPosition(td, p);
835 
836  // Check/set new parcel thermo properties
837  cloud.setParcelThermoProperties(p);
838 
839  // Assign new parcel properties in injection model
840  setProperties(parceli, nParcels, 0, td, p);
841 
842  // Check/set new parcel injection properties
843  cloud.checkParcelProperties(p, index());
844 
845  // Apply correction to velocity for 2-D cases
846  meshTools::constrainDirection(mesh, mesh.solutionD(), p.U());
847 
848  // Initial step fraction
849  p.stepFraction() = 0;
850 
851  // Set the number of particles. If not fixed, this will set
852  // a junk value, which will get corrected below.
853  p.nParticle() = nParticleFixed_;
854  }
855  }
856 
857  // Set the number of particles so that the introduced mass is correct
858  // and the uniform size is as specified
859  if (nParticleFixed_ < 0)
860  {
861  setNumberOfParticles(parcelPtrs, mass);
862  }
863 
864  // Add the new parcels
865  forAll(parcelPtrs, parceli)
866  {
867  if (parcelPtrs.set(parceli))
868  {
869  parcelType& p = parcelPtrs[parceli];
870  nParcelsAdded ++;
871  massAdded += p.nParticle()*p.mass();
872  cloud.addParticle(parcelPtrs.set(parceli, nullptr).ptr());
873  }
874  }
875  }
876 
877  postInject(nParcelsAdded, massAdded, td);
878 }
879 
880 
881 template<class CloudType>
883 {
884  os << " " << this->modelName() << ":" << nl
885  << " number of parcels added = " << parcelsAddedTotal_ << nl
886  << " mass introduced = " << massInjected_ << nl;
887 
888  if (this->writeTime())
889  {
890  this->setModelProperty("massInjected", massInjected_);
891  this->setModelProperty("nInjections", nInjections_);
892  this->setModelProperty("parcelsAddedTotal", parcelsAddedTotal_);
893  this->setModelProperty("timeStep0", timeStep0_);
894  }
895 }
896 
897 
898 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
899 
900 #include "InjectionModelNew.C"
901 
902 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Base class for cloud sub-models.
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:225
const fvMesh & mesh() const
Return references to the mesh.
Definition: DSMCCloudI.H:41
Run-time selectable general function of one variable.
Definition: Function1.H:125
Templated function that returns a constant value.
Definition: Constant.H:60
Function1 which scales a given 'value' function by a 'scale' scalar function and scales the 'x' argum...
Definition: Scale.H:155
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
Templated injection model class.
virtual void topoChange()
Update mesh.
virtual ~InjectionModel()
Destructor.
void inject(TrackCloudType &cloud, typename CloudType::parcelType::trackingData &td)
Main injection loop.
autoPtr< Function1< scalar > > readParcelsPerSecond(const dictionary &dict, CloudType &owner)
Read the number of parcels injected per second for continuous.
void constrainPosition(typename CloudType::parcelType::trackingData &td, typename CloudType::parcelType &parcel)
Constrain a parcel's position appropriately to the geometric.
virtual void postInject(const label parcelsAdded, const scalar massAdded, typename parcelType::trackingData &td)
Post injection hook.
virtual void preInject(typename parcelType::trackingData &td)
Pre injection hook.
bool findCellAtPosition(const point &position, barycentric &coordinates, label &celli, label &tetFacei, label &tetPti, bool errorOnNotFound=true)
Find the cell that contains the supplied position.
scalar readMassTotal(const dictionary &dict, CloudType &owner)
Read the total mass value for instantaneous injections.
virtual void info(Ostream &os)
Write injection info to stream.
void setNumberOfParticles(PtrList< parcelType > &parcelPtrs, const scalar mass) const
Set number of particles to inject given parcel properties.
autoPtr< Function1< scalar > > readMassFlowRate(const dictionary &dict, CloudType &owner, const scalar duration)
Read the mass flow rate function for continuous injections.
uniformParcelSize uniformParcelSize_
Size uniform to all parcels.
label sizeSampleQ() const
Return the sampling moment to be used by the size distribution.
scalar averageParcelMass()
Return the average injected parcel mass.
scalar nParticleFixed_
Fixed nParticle to assign to parcels. Only valid if.
label index() const
Get the index of this injector.
void injectSteadyState(TrackCloudType &cloud, typename CloudType::parcelType::trackingData &td)
Main injection loop - steady-state.
scalar readDuration(const dictionary &dict, CloudType &owner)
Read the duration for continuous injections.
InjectionModel(CloudType &owner)
Construct null from owner.
CloudType::parcelType parcelType
Convenience typedef for parcelType.
scalar SOI_
Start of injection [s].
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
const Type & second() const
Return second.
Definition: Pair.H:110
const Type & first() const
Return first.
Definition: Pair.H:98
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
const unitConversion & userUnits() const
Return the user-time unit conversion.
Definition: Time.C:869
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
const Type & value() const
Return const reference to value.
uniformParcelSize
Enumeration for the parcels' uniform size.
static const NamedEnum< uniformParcelSize, 3 > uniformParcelSizeNames_
Names of the parcels' uniform size.
const Time & time() const
Return time.
Base particle class.
Definition: particle.H:83
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:1000
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:1006
const dictionary & dict() const
Return const access to the cloud dictionary.
Definition: subModelBase.C:110
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
#define WarningInFunction
Report a warning using Foam::Warning.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
void constrainDirection(const polyMesh &mesh, const Vector< label > &dirs, vector &d)
Set the constrained components of directions/velocity to zero.
Definition: meshTools.C:671
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedScalar pos(const dimensionedScalar &ds)
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
const dimensionSet dimless
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
const dimensionSet dimMass
error FatalError
static const label labelMax
Definition: label.H:62
static const char nl
Definition: Ostream.H:266
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict
volScalarField & p