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