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