InjectionModel.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) 2016-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 "InjectionModel.H"
30 #include "mathematicalConstants.H"
31 #include "meshTools.H"
32 #include "volFields.H"
33 
34 using namespace Foam::constant::mathematical;
35 
36 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
37 
38 template<class CloudType>
40 (
41  const scalar time,
42  label& newParcels,
43  scalar& newVolumeFraction
44 )
45 {
46  // Initialise values
47  newParcels = 0;
48  newVolumeFraction = 0.0;
49  bool validInjection = false;
50 
51  // Return if not started injection event
52  if (time < SOI_)
53  {
54  timeStep0_ = time;
55  return validInjection;
56  }
57 
58  // Make times relative to SOI
59  scalar t0 = timeStep0_ - SOI_;
60  scalar t1 = time - SOI_;
61 
62  // Number of parcels to inject
63  newParcels = this->parcelsToInject(t0, t1);
64 
65  // Volume of parcels to inject
66  newVolumeFraction =
67  this->volumeToInject(t0, t1)
68  /(volumeTotal_ + ROOTVSMALL);
69 
70  if (newVolumeFraction > 0)
71  {
72  if (newParcels > 0)
73  {
74  timeStep0_ = time;
75  validInjection = true;
76  }
77  else
78  {
79  // Injection should have started, but not sufficient volume to
80  // produce (at least) 1 parcel - hold value of timeStep0_
81  validInjection = false;
82  }
83  }
84  else
85  {
86  timeStep0_ = time;
87  validInjection = false;
88  }
89 
90  return validInjection;
91 }
92 
93 
94 template<class CloudType>
96 (
97  label& celli,
98  label& tetFacei,
99  label& tetPti,
100  vector& position,
101  bool errorOnNotFound
102 )
103 {
104  const volVectorField& cellCentres = this->owner().mesh().C();
105 
106  const vector p0 = position;
107 
108  this->owner().mesh().findCellFacePt
109  (
110  position,
111  celli,
112  tetFacei,
113  tetPti
114  );
115 
116  label proci = -1;
117 
118  if (celli >= 0)
119  {
120  proci = Pstream::myProcNo();
121  }
122 
123  reduce(proci, maxOp<label>());
124 
125  // Ensure that only one processor attempts to insert this Parcel
126 
127  if (proci != Pstream::myProcNo())
128  {
129  celli = -1;
130  tetFacei = -1;
131  tetPti = -1;
132  }
133 
134  // Last chance - find nearest cell and try that one - the point is
135  // probably on an edge
136  if (proci == -1)
137  {
138  celli = this->owner().mesh().findNearestCell(position);
139 
140  if (celli >= 0)
141  {
142  position += SMALL*(cellCentres[celli] - position);
143 
144  this->owner().mesh().findCellFacePt
145  (
146  position,
147  celli,
148  tetFacei,
149  tetPti
150  );
151 
152  if (celli > 0)
153  {
154  proci = Pstream::myProcNo();
155  }
156  }
157 
158  reduce(proci, maxOp<label>());
159 
160  if (proci != Pstream::myProcNo())
161  {
162  celli = -1;
163  tetFacei = -1;
164  tetPti = -1;
165  }
166  }
167 
168  if (proci == -1)
169  {
170  if (errorOnNotFound)
171  {
173  << "Cannot find parcel injection cell. "
174  << "Parcel position = " << p0 << nl
175  << exit(FatalError);
176  }
177  else
178  {
179  return false;
180  }
181  }
182 
183  return true;
184 }
185 
186 
187 template<class CloudType>
189 (
190  const label parcels,
191  const scalar volumeFraction,
192  const scalar diameter,
193  const scalar rho
194 )
195 {
196  scalar nP = 0.0;
197  switch (parcelBasis_)
198  {
199  case pbMass:
200  {
201  scalar volumep = pi/6.0*pow3(diameter);
202  scalar volumeTot = massTotal_/rho;
203 
204  nP = (volumeFraction*volumeTot + delayedVolume_)/(parcels*volumep);
205  break;
206  }
207  case pbNumber:
208  {
209  nP = massTotal_/(rho*volumeTotal_);
210  break;
211  }
212  case pbFixed:
213  {
214  nP = nParticleFixed_;
215  break;
216  }
217  default:
218  {
219  nP = 0.0;
221  << "Unknown parcelBasis type" << nl
222  << exit(FatalError);
223  }
224  }
225 
226  return nP;
227 }
228 
229 
230 template<class CloudType>
232 (
233  const label parcelsAdded,
234  const scalar massAdded
235 )
236 {
237  const label allParcelsAdded = returnReduce(parcelsAdded, sumOp<label>());
238 
239  if (allParcelsAdded > 0)
240  {
241  Info<< nl
242  << "Cloud: " << this->owner().name()
243  << " injector: " << this->modelName() << nl
244  << " Added " << allParcelsAdded << " new parcels" << nl << endl;
245  }
246 
247  // Increment total number of parcels added
248  parcelsAddedTotal_ += allParcelsAdded;
249 
250  // Increment total mass injected
251  massInjected_ += returnReduce(massAdded, sumOp<scalar>());
252 
253  // Update time for start of next injection
254  time0_ = this->owner().db().time().value();
255 
256  // Increment number of injections
257  ++nInjections_;
258 }
259 
260 
261 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
262 
263 template<class CloudType>
265 :
267  SOI_(0.0),
268  volumeTotal_(this->template getModelProperty<scalar>("volumeTotal")),
269  massTotal_(0),
270  massFlowRate_(owner.db().time(), "massFlowRate"),
271  massInjected_(this->template getModelProperty<scalar>("massInjected")),
272  nInjections_(this->template getModelProperty<label>("nInjections")),
273  parcelsAddedTotal_
274  (
275  this->template getModelProperty<scalar>("parcelsAddedTotal")
276  ),
277  parcelBasis_(pbNumber),
278  nParticleFixed_(0.0),
279  time0_(0.0),
280  timeStep0_(this->template getModelProperty<scalar>("timeStep0")),
281  minParticlesPerParcel_(1),
282  delayedVolume_(0.0),
283  injectorID_(-1),
284  ignoreOutOfBounds_(false)
285 {}
286 
287 
288 template<class CloudType>
290 (
291  const dictionary& dict,
292  CloudType& owner,
293  const word& modelName,
294  const word& modelType
295 )
296 :
297  CloudSubModelBase<CloudType>(modelName, owner, dict, typeName, modelType),
298  SOI_(0.0),
299  volumeTotal_(this->template getModelProperty<scalar>("volumeTotal")),
300  massTotal_(0),
301  massFlowRate_(owner.db().time(), "massFlowRate"),
302  massInjected_(this->template getModelProperty<scalar>("massInjected")),
303  nInjections_(this->template getModelProperty<scalar>("nInjections")),
304  parcelsAddedTotal_
305  (
306  this->template getModelProperty<scalar>("parcelsAddedTotal")
307  ),
308  parcelBasis_(pbNumber),
309  nParticleFixed_(0.0),
310  time0_(owner.db().time().value()),
311  timeStep0_(this->template getModelProperty<scalar>("timeStep0")),
312  minParticlesPerParcel_
313  (
314  this->coeffDict().getOrDefault("minParticlesPerParcel", scalar(1))
315  ),
316  delayedVolume_(0.0),
317  injectorID_(this->coeffDict().getOrDefault("injectorID", -1)),
318  ignoreOutOfBounds_
319  (
320  this->coeffDict().getOrDefault("ignoreOutOfBounds", false)
321  )
322 {
323  // Provide some info
324  // - also serves to initialise mesh dimensions - needed for parallel runs
325  // due to lazy evaluation of valid mesh dimensions
326  Info<< " Constructing " << owner.mesh().nGeometricD() << "-D injection"
327  << endl;
328 
329  if (injectorID_ != -1)
330  {
331  Info<< " injector ID: " << injectorID_ << endl;
332  }
333 
334  if (owner.solution().active())
335  {
336  if (owner.solution().transient())
337  {
338  this->coeffDict().readEntry("massTotal", massTotal_);
339  this->coeffDict().readEntry("SOI", SOI_);
340  }
341  else
342  {
343  massFlowRate_.reset(this->coeffDict());
344  massTotal_ = massFlowRate_.value(owner.db().time().value());
345  this->coeffDict().readIfPresent("SOI", SOI_);
346  }
347  }
348 
349  SOI_ = owner.db().time().userTimeToTime(SOI_);
350 
351  const word parcelBasisType(this->coeffDict().getWord("parcelBasisType"));
352 
353  if (parcelBasisType == "mass")
354  {
355  parcelBasis_ = pbMass;
356  }
357  else if (parcelBasisType == "number")
358  {
359  parcelBasis_ = pbNumber;
360  }
361  else if (parcelBasisType == "fixed")
362  {
363  parcelBasis_ = pbFixed;
364  this->coeffDict().readEntry("nParticle", nParticleFixed_);
365 
366  Info<< " Choosing nParticle to be a fixed value, massTotal "
367  << "variable now does not determine anything."
368  << endl;
369  }
370  else
371  {
373  << "parcelBasisType must be either 'number', 'mass' or 'fixed'"
374  << nl << exit(FatalError);
375  }
376 }
377 
378 
379 template<class CloudType>
381 (
382  const InjectionModel<CloudType>& im
383 )
384 :
386  SOI_(im.SOI_),
387  volumeTotal_(im.volumeTotal_),
388  massTotal_(im.massTotal_),
389  massFlowRate_(im.massFlowRate_),
390  massInjected_(im.massInjected_),
391  nInjections_(im.nInjections_),
392  parcelsAddedTotal_(im.parcelsAddedTotal_),
393  parcelBasis_(im.parcelBasis_),
394  nParticleFixed_(im.nParticleFixed_),
395  time0_(im.time0_),
396  timeStep0_(im.timeStep0_),
397  minParticlesPerParcel_(im.minParticlesPerParcel_),
398  delayedVolume_(im.delayedVolume_),
399  injectorID_(im.injectorID_),
400  ignoreOutOfBounds_(im.ignoreOutOfBounds_)
401 {}
402 
403 
404 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
405 
406 template<class CloudType>
408 {}
409 
410 
411 template<class CloudType>
413 {
414  label nTotal = 0.0;
415  if (this->owner().solution().transient())
416  {
417  nTotal = parcelsToInject(0.0, timeEnd() - timeStart());
418  }
419  else
420  {
421  nTotal = parcelsToInject(0.0, 1.0);
422  }
423 
424  return massTotal_/nTotal;
425 }
426 
427 
428 template<class CloudType>
429 template<class TrackCloudType>
431 (
432  TrackCloudType& cloud,
433  typename CloudType::parcelType::trackingData& td
434 )
435 {
436  if (!this->active())
437  {
438  return;
439  }
440 
441  const scalar time = this->owner().db().time().value();
442 
443  // Prepare for next time step
444  label parcelsAdded = 0;
445  scalar massAdded = 0.0;
446  label newParcels = 0;
447  scalar newVolumeFraction = 0.0;
448  scalar delayedVolume = 0;
449 
450  if (prepareForNextTimeStep(time, newParcels, newVolumeFraction))
451  {
452  const scalar trackTime = this->owner().solution().trackTime();
453  const polyMesh& mesh = this->owner().mesh();
454 
455  // Duration of injection period during this timestep
456  const scalar deltaT =
457  max(0.0, min(trackTime, min(time - SOI_, timeEnd() - time0_)));
458 
459  // Pad injection time if injection starts during this timestep
460  const scalar padTime = max(0.0, SOI_ - time0_);
461 
462  // Introduce new parcels linearly across carrier phase timestep
463  for (label parcelI = 0; parcelI < newParcels; parcelI++)
464  {
465  if (validInjection(parcelI))
466  {
467  // Calculate the pseudo time of injection for parcel 'parcelI'
468  scalar timeInj = time0_ + padTime + deltaT*parcelI/newParcels;
469 
470  // Determine the injection position and owner cell,
471  // tetFace and tetPt
472  label celli = -1;
473  label tetFacei = -1;
474  label tetPti = -1;
475 
476  vector pos = Zero;
477 
478  setPositionAndCell
479  (
480  parcelI,
481  newParcels,
482  timeInj,
483  pos,
484  celli,
485  tetFacei,
486  tetPti
487  );
488 
489  if (celli > -1)
490  {
491  // Lagrangian timestep
492  const scalar dt = time - timeInj;
493 
494  // Apply corrections to position for 2-D cases
496 
497  // Create a new parcel
498  parcelType* pPtr = new parcelType(mesh, pos, celli);
499 
500  // Check/set new parcel thermo properties
501  cloud.setParcelThermoProperties(*pPtr, dt);
502 
503  // Assign new parcel properties in injection model
504  setProperties(parcelI, newParcels, timeInj, *pPtr);
505 
506  // Check/set new parcel injection properties
507  cloud.checkParcelProperties(*pPtr, dt, fullyDescribed());
508 
509  // Apply correction to velocity for 2-D cases
511  (
512  mesh,
513  mesh.solutionD(),
514  pPtr->U()
515  );
516 
517  // Number of particles per parcel
518  pPtr->nParticle() =
519  setNumberOfParticles
520  (
521  newParcels,
522  newVolumeFraction,
523  pPtr->d(),
524  pPtr->rho()
525  );
526 
527  if (pPtr->nParticle() >= minParticlesPerParcel_)
528  {
529  parcelsAdded++;
530  massAdded += pPtr->nParticle()*pPtr->mass();
531 
532  if (pPtr->move(cloud, td, dt))
533  {
534  pPtr->typeId() = injectorID_;
535  cloud.addParticle(pPtr);
536  }
537  else
538  {
539  delete pPtr;
540  }
541  }
542  else
543  {
544  delayedVolume += pPtr->nParticle()*pPtr->volume();
545  delete pPtr;
546  }
547  }
548  }
549  }
550  }
551 
552  delayedVolume_ = returnReduce(delayedVolume, sumOp<scalar>());
553 
554  postInjectCheck(parcelsAdded, massAdded);
555 }
556 
557 
558 template<class CloudType>
559 template<class TrackCloudType>
561 (
562  TrackCloudType& cloud,
563  typename CloudType::parcelType::trackingData& td,
564  const scalar trackTime
565 )
566 {
567  if (!this->active())
568  {
569  return;
570  }
571 
572  const scalar time = this->owner().db().time().value();
573 
574  if (time < SOI_)
575  {
576  return;
577  }
578 
579  const polyMesh& mesh = this->owner().mesh();
580 
581  massTotal_ = massFlowRate_.value(mesh.time().value());
582 
583  // Reset counters
584  time0_ = 0.0;
585  label parcelsAdded = 0;
586  scalar massAdded = 0.0;
587 
588  // Set number of new parcels to inject based on first second of injection
589  label newParcels = parcelsToInject(0.0, 1.0);
590 
591  // Inject new parcels
592  for (label parcelI = 0; parcelI < newParcels; parcelI++)
593  {
594  // Volume to inject is split equally amongst all parcel streams
595  scalar newVolumeFraction = 1.0/scalar(newParcels);
596 
597  // Determine the injection position and owner cell,
598  // tetFace and tetPt
599  label celli = -1;
600  label tetFacei = -1;
601  label tetPti = -1;
602 
603  vector pos = Zero;
604 
605  setPositionAndCell
606  (
607  parcelI,
608  newParcels,
609  0.0,
610  pos,
611  celli,
612  tetFacei,
613  tetPti
614  );
615 
616  if (celli > -1)
617  {
618  // Apply corrections to position for 2-D cases
620 
621  // Create a new parcel
622  parcelType* pPtr = new parcelType(mesh, pos, celli);
623 
624  // Check/set new parcel thermo properties
625  cloud.setParcelThermoProperties(*pPtr, 0.0);
626 
627  // Assign new parcel properties in injection model
628  setProperties(parcelI, newParcels, 0.0, *pPtr);
629 
630  // Check/set new parcel injection properties
631  cloud.checkParcelProperties(*pPtr, 0.0, fullyDescribed());
632 
633  // Apply correction to velocity for 2-D cases
635 
636  // Number of particles per parcel
637  pPtr->nParticle() =
638  setNumberOfParticles
639  (
640  1,
641  newVolumeFraction,
642  pPtr->d(),
643  pPtr->rho()
644  );
645 
646  pPtr->typeId() = injectorID_;
647 
648  // Add the new parcel
649  cloud.addParticle(pPtr);
650 
651  massAdded += pPtr->nParticle()*pPtr->mass();
652  parcelsAdded++;
653  }
654  }
655 
656  postInjectCheck(parcelsAdded, massAdded);
657 }
658 
659 
660 template<class CloudType>
662 {
663  os << " Injector " << this->modelName() << ":" << nl
664  << " - parcels added = " << parcelsAddedTotal_ << nl
665  << " - mass introduced = " << massInjected_ << nl;
666 
667  if (this->writeTime())
668  {
669  this->setModelProperty("volumeTotal", volumeTotal_);
670  this->setModelProperty("massInjected", massInjected_);
671  this->setModelProperty("nInjections", nInjections_);
672  this->setModelProperty("parcelsAddedTotal", parcelsAddedTotal_);
673  this->setModelProperty("timeStep0", timeStep0_);
674  }
675 }
676 
677 
678 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
679 
680 #include "InjectionModelNew.C"
681 
682 // ************************************************************************* //
Foam::meshTools::constrainDirection
void constrainDirection(const polyMesh &mesh, const Vector< label > &dirs, vector &d)
Set the constrained components of directions/velocity to zero.
Definition: meshTools.C:687
volFields.H
Foam::maxOp
Definition: ops.H:223
meshTools.H
mathematicalConstants.H
Foam::InjectionModel::info
virtual void info(Ostream &os)
Write injection info to stream.
Definition: InjectionModel.C:661
Foam::InjectionModel::findCellAtPosition
virtual bool findCellAtPosition(label &celli, label &tetFacei, label &tetPti, vector &position, bool errorOnNotFound=true)
Find the cell that contains the supplied position.
Definition: InjectionModel.C:96
Foam::solution
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:51
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::polyMesh::nGeometricD
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:869
Foam::InjectionModel::InjectionModel
InjectionModel(CloudType &owner)
Construct null from owner.
Definition: InjectionModel.C:264
InjectionModelNew.C
Foam::InjectionModel::updateMesh
virtual void updateMesh()
Update mesh.
Definition: InjectionModel.C:407
Foam::InjectionModel::massFlowRate_
TimeFunction1< scalar > massFlowRate_
Mass flow rate profile for steady calculations.
Definition: InjectionModel.H:112
Foam::InjectionModel
Templated injection model class.
Definition: InjectionModel.H:73
Foam::CloudSubModelBase
Base class for cloud sub-models.
Definition: CloudSubModelBase.H:51
Foam::InjectionModel::prepareForNextTimeStep
virtual bool prepareForNextTimeStep(const scalar time, label &newParcels, scalar &newVolumeFraction)
Determine properties for next time step/injection interval.
Definition: InjectionModel.C:40
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::polyMesh::solutionD
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:875
Foam::InjectionModel::massInjected_
scalar massInjected_
Total mass injected to date [kg].
Definition: InjectionModel.H:115
Foam::InjectionModel::setNumberOfParticles
virtual scalar setNumberOfParticles(const label parcels, const scalar volumeFraction, const scalar diameter, const scalar rho)
Set number of particles to inject given parcel properties.
Definition: InjectionModel.C:189
Foam::InjectionModel::delayedVolume_
scalar delayedVolume_
Definition: InjectionModel.H:148
InjectionModel.H
Foam::InjectionModel::volumeTotal_
scalar volumeTotal_
Definition: InjectionModel.H:106
rho
rho
Definition: readInitialConditions.H:88
Foam::InjectionModel::timeStep0_
scalar timeStep0_
Time at start of injection time step [s].
Definition: InjectionModel.H:140
Foam::InjectionModel::injectorID_
label injectorID_
Optional injector ID.
Definition: InjectionModel.H:151
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::primitiveMesh::findNearestCell
label findNearestCell(const point &location) const
Find the cell with the nearest cell centre to location.
Definition: primitiveMeshFindCell.C:88
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:44
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::InjectionModel::injectSteadyState
void injectSteadyState(TrackCloudType &cloud, typename CloudType::parcelType::trackingData &td, const scalar trackTime)
Main injection loop - steady-state.
Definition: InjectionModel.C:561
Foam::InjectionModel::time0_
scalar time0_
Continuous phase time at start of injection time step [s].
Definition: InjectionModel.H:137
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::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
Foam::InjectionModel::minParticlesPerParcel_
scalar minParticlesPerParcel_
Definition: InjectionModel.H:144
Foam::InjectionModel::postInjectCheck
virtual void postInjectCheck(const label parcelsAdded, const scalar massAdded)
Post injection checks.
Definition: InjectionModel.C:232
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:341
Foam::InjectionModel::averageParcelMass
virtual scalar averageParcelMass()
Return the average parcel mass over the injection period.
Definition: InjectionModel.C:412
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::polyMesh::findCellFacePt
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition: polyMesh.C:1356
Foam::InjectionModel::massTotal_
scalar massTotal_
Total mass to inject [kg].
Definition: InjectionModel.H:109
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::InjectionModel::parcelsAddedTotal_
label parcelsAddedTotal_
Running counter of total number of parcels added.
Definition: InjectionModel.H:124
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:57
Foam::InjectionModel::ignoreOutOfBounds_
bool ignoreOutOfBounds_
Definition: InjectionModel.H:155
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
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::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::InjectionModel::inject
void inject(TrackCloudType &cloud, typename CloudType::parcelType::trackingData &td)
Main injection loop.
Definition: InjectionModel.C:431
Foam::InjectionModel::parcelType
CloudType::parcelType parcelType
Convenience typedef for parcelType.
Definition: InjectionModel.H:80
Foam::Vector< scalar >
Foam::InjectionModel::SOI_
scalar SOI_
Start of injection [s].
Definition: InjectionModel.H:102
Foam::meshTools::constrainToMeshCentre
void constrainToMeshCentre(const polyMesh &mesh, point &pt)
Set the constrained components of position to mesh centre.
Definition: meshTools.C:629
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:275
Foam::InjectionModel::nInjections_
label nInjections_
Number of injections counter.
Definition: InjectionModel.H:121
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
p0
const volScalarField & p0
Definition: EEqn.H:36
Foam::InjectionModel::parcelBasis_
parcelBasis parcelBasis_
Parcel basis enumeration.
Definition: InjectionModel.H:130
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::InjectionModel::nParticleFixed_
scalar nParticleFixed_
Definition: InjectionModel.H:134
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177