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-2017 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  this->owner().mesh().findCellFacePt
142  (
143  position,
144  celli,
145  tetFacei,
146  tetPti
147  );
148 
149  if (celli > 0)
150  {
151  proci = Pstream::myProcNo();
152  }
153  }
154 
155  reduce(proci, maxOp<label>());
156 
157  if (proci != Pstream::myProcNo())
158  {
159  celli = -1;
160  tetFacei = -1;
161  tetPti = -1;
162  }
163  }
164 
165  if (proci == -1)
166  {
167  if (errorOnNotFound)
168  {
170  << "Cannot find parcel injection cell. "
171  << "Parcel position = " << p0 << nl
172  << exit(FatalError);
173  }
174  else
175  {
176  return false;
177  }
178  }
179 
180  return true;
181 }
182 
183 
184 template<class CloudType>
186 (
187  const label parcels,
188  const scalar volumeFraction,
189  const scalar diameter,
190  const scalar rho
191 )
192 {
193  scalar nP = 0.0;
194  switch (parcelBasis_)
195  {
196  case pbMass:
197  {
198  scalar volumep = pi/6.0*pow3(diameter);
199  scalar volumeTot = massTotal_/rho;
200 
201  nP = (volumeFraction*volumeTot + delayedVolume_)/(parcels*volumep);
202  break;
203  }
204  case pbNumber:
205  {
206  nP = massTotal_/(rho*volumeTotal_);
207  break;
208  }
209  case pbFixed:
210  {
211  nP = nParticleFixed_;
212  break;
213  }
214  default:
215  {
216  nP = 0.0;
218  << "Unknown parcelBasis type" << nl
219  << exit(FatalError);
220  }
221  }
222 
223  return nP;
224 }
225 
226 
227 template<class CloudType>
229 (
230  const label parcelsAdded,
231  const scalar massAdded
232 )
233 {
234  const label allParcelsAdded = returnReduce(parcelsAdded, sumOp<label>());
235 
236  if (allParcelsAdded > 0)
237  {
238  Info<< nl
239  << "Cloud: " << this->owner().name()
240  << " injector: " << this->modelName() << nl
241  << " Added " << allParcelsAdded << " new parcels" << nl << endl;
242  }
243 
244  // Increment total number of parcels added
245  parcelsAddedTotal_ += allParcelsAdded;
246 
247  // Increment total mass injected
248  massInjected_ += returnReduce(massAdded, sumOp<scalar>());
249 
250  // Update time for start of next injection
251  time0_ = this->owner().db().time().value();
252 
253  // Increment number of injections
254  nInjections_++;
255 }
256 
257 
258 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
259 
260 template<class CloudType>
262 :
264  SOI_(0.0),
265  volumeTotal_(0.0),
266  massTotal_(0.0),
267  massFlowRate_(owner.db().time(), "massFlowRate"),
268  massInjected_(this->template getModelProperty<scalar>("massInjected")),
269  nInjections_(this->template getModelProperty<label>("nInjections")),
270  parcelsAddedTotal_
271  (
272  this->template getModelProperty<scalar>("parcelsAddedTotal")
273  ),
274  parcelBasis_(pbNumber),
275  nParticleFixed_(0.0),
276  time0_(0.0),
277  timeStep0_(this->template getModelProperty<scalar>("timeStep0")),
278  delayedVolume_(0.0)
279 {}
280 
281 
282 template<class CloudType>
284 (
285  const dictionary& dict,
286  CloudType& owner,
287  const word& modelName,
288  const word& modelType
289 )
290 :
292  SOI_(0.0),
293  volumeTotal_(0.0),
294  massTotal_(0.0),
295  massFlowRate_(owner.db().time(), "massFlowRate"),
296  massInjected_(this->template getModelProperty<scalar>("massInjected")),
297  nInjections_(this->template getModelProperty<scalar>("nInjections")),
299  (
300  this->template getModelProperty<scalar>("parcelsAddedTotal")
301  ),
303  nParticleFixed_(0.0),
304  time0_(owner.db().time().value()),
305  timeStep0_(this->template getModelProperty<scalar>("timeStep0")),
306  delayedVolume_(0.0)
307 {
308  // Provide some info
309  // - also serves to initialise mesh dimensions - needed for parallel runs
310  // due to lazy evaluation of valid mesh dimensions
311  Info<< " Constructing " << owner.mesh().nGeometricD() << "-D injection"
312  << endl;
313 
314  if (owner.solution().transient())
315  {
316  this->coeffDict().lookup("massTotal") >> massTotal_;
317  this->coeffDict().lookup("SOI") >> SOI_;
318  SOI_ = owner.db().time().userTimeToTime(SOI_);
319  }
320  else
321  {
322  massFlowRate_.reset(this->coeffDict());
323  massTotal_ = massFlowRate_.value(owner.db().time().value());
324  }
325 
326  const word parcelBasisType = this->coeffDict().lookup("parcelBasisType");
327 
328  if (parcelBasisType == "mass")
329  {
331  }
332  else if (parcelBasisType == "number")
333  {
335  }
336  else if (parcelBasisType == "fixed")
337  {
339 
340  Info<< " Choosing nParticle to be a fixed value, massTotal "
341  << "variable now does not determine anything."
342  << endl;
343 
344  nParticleFixed_ = readScalar(this->coeffDict().lookup("nParticle"));
345  }
346  else
347  {
349  << "parcelBasisType must be either 'number', 'mass' or 'fixed'"
350  << nl << exit(FatalError);
351  }
352 }
353 
354 
355 template<class CloudType>
357 (
358  const InjectionModel<CloudType>& im
359 )
360 :
362  SOI_(im.SOI_),
371  time0_(im.time0_),
374 {}
375 
376 
377 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
378 
379 template<class CloudType>
381 {}
382 
383 
384 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
385 
386 template<class CloudType>
388 {}
389 
390 
391 template<class CloudType>
393 {
394  label nTotal = 0.0;
395  if (this->owner().solution().transient())
396  {
397  nTotal = parcelsToInject(0.0, timeEnd() - timeStart());
398  }
399  else
400  {
401  nTotal = parcelsToInject(0.0, 1.0);
402  }
403 
404  return massTotal_/nTotal;
405 }
406 
407 
408 template<class CloudType>
409 template<class TrackData>
411 {
412  if (!this->active())
413  {
414  return;
415  }
416 
417  const scalar time = this->owner().db().time().value();
418 
419  // Prepare for next time step
420  label parcelsAdded = 0;
421  scalar massAdded = 0.0;
422  label newParcels = 0;
423  scalar newVolumeFraction = 0.0;
424 
425  if (prepareForNextTimeStep(time, newParcels, newVolumeFraction))
426  {
427  scalar delayedVolume = 0;
428 
429  const scalar trackTime = this->owner().solution().trackTime();
430  const polyMesh& mesh = this->owner().mesh();
431  typename TrackData::cloudType& cloud = td.cloud();
432 
433  // Duration of injection period during this timestep
434  const scalar deltaT =
435  max(0.0, min(trackTime, min(time - SOI_, timeEnd() - time0_)));
436 
437  // Pad injection time if injection starts during this timestep
438  const scalar padTime = max(0.0, SOI_ - time0_);
439 
440  // Introduce new parcels linearly across carrier phase timestep
441  for (label parcelI = 0; parcelI < newParcels; parcelI++)
442  {
443  if (validInjection(parcelI))
444  {
445  // Calculate the pseudo time of injection for parcel 'parcelI'
446  scalar timeInj = time0_ + padTime + deltaT*parcelI/newParcels;
447 
448  // Determine the injection position and owner cell,
449  // tetFace and tetPt
450  label celli = -1;
451  label tetFacei = -1;
452  label tetPti = -1;
453 
454  vector pos = Zero;
455 
457  (
458  parcelI,
459  newParcels,
460  timeInj,
461  pos,
462  celli,
463  tetFacei,
464  tetPti
465  );
466 
467  if (celli > -1)
468  {
469  // Lagrangian timestep
470  const scalar dt = time - timeInj;
471 
472  // Apply corrections to position for 2-D cases
474 
475  // Create a new parcel
476  parcelType* pPtr = new parcelType(mesh, pos, celli);
477 
478  // Check/set new parcel thermo properties
479  cloud.setParcelThermoProperties(*pPtr, dt);
480 
481  // Assign new parcel properties in injection model
482  setProperties(parcelI, newParcels, timeInj, *pPtr);
483 
484  // Check/set new parcel injection properties
485  cloud.checkParcelProperties(*pPtr, dt, fullyDescribed());
486 
487  // Apply correction to velocity for 2-D cases
489  (
490  mesh,
491  mesh.solutionD(),
492  pPtr->U()
493  );
494 
495  // Number of particles per parcel
496  pPtr->nParticle() =
498  (
499  newParcels,
500  newVolumeFraction,
501  pPtr->d(),
502  pPtr->rho()
503  );
504 
505  if (pPtr->nParticle() >= 1.0)
506  {
507  parcelsAdded++;
508  massAdded += pPtr->nParticle()*pPtr->mass();
509 
510  if (pPtr->move(td, dt))
511  {
512  td.cloud().addParticle(pPtr);
513  }
514  else
515  {
516  delete pPtr;
517  }
518  }
519  else
520  {
521  delayedVolume += pPtr->nParticle()*pPtr->volume();
522  delete pPtr;
523  }
524  }
525  }
526  }
527 
528  delayedVolume_ = delayedVolume;
529  }
530 
531  postInjectCheck(parcelsAdded, massAdded);
532 }
533 
534 
535 template<class CloudType>
536 template<class TrackData>
538 (
539  TrackData& td,
540  const scalar trackTime
541 )
542 {
543  if (!this->active())
544  {
545  return;
546  }
547 
548  const polyMesh& mesh = this->owner().mesh();
549  typename TrackData::cloudType& cloud = td.cloud();
550 
552 
553  // Reset counters
554  time0_ = 0.0;
555  label parcelsAdded = 0;
556  scalar massAdded = 0.0;
557 
558  // Set number of new parcels to inject based on first second of injection
559  label newParcels = parcelsToInject(0.0, 1.0);
560 
561  // Inject new parcels
562  for (label parcelI = 0; parcelI < newParcels; parcelI++)
563  {
564  // Volume to inject is split equally amongst all parcel streams
565  scalar newVolumeFraction = 1.0/scalar(newParcels);
566 
567  // Determine the injection position and owner cell,
568  // tetFace and tetPt
569  label celli = -1;
570  label tetFacei = -1;
571  label tetPti = -1;
572 
573  vector pos = Zero;
574 
576  (
577  parcelI,
578  newParcels,
579  0.0,
580  pos,
581  celli,
582  tetFacei,
583  tetPti
584  );
585 
586  if (celli > -1)
587  {
588  // Apply corrections to position for 2-D cases
590 
591  // Create a new parcel
592  parcelType* pPtr = new parcelType(mesh, pos, celli);
593 
594  // Check/set new parcel thermo properties
595  cloud.setParcelThermoProperties(*pPtr, 0.0);
596 
597  // Assign new parcel properties in injection model
598  setProperties(parcelI, newParcels, 0.0, *pPtr);
599 
600  // Check/set new parcel injection properties
601  cloud.checkParcelProperties(*pPtr, 0.0, fullyDescribed());
602 
603  // Apply correction to velocity for 2-D cases
604  meshTools::constrainDirection(mesh, mesh.solutionD(), pPtr->U());
605 
606  // Number of particles per parcel
607  pPtr->nParticle() =
609  (
610  1,
611  newVolumeFraction,
612  pPtr->d(),
613  pPtr->rho()
614  );
615 
616  // Add the new parcel
617  td.cloud().addParticle(pPtr);
618 
619  massAdded += pPtr->nParticle()*pPtr->mass();
620  parcelsAdded++;
621  }
622  }
623 
624  postInjectCheck(parcelsAdded, massAdded);
625 }
626 
627 
628 template<class CloudType>
630 {
631  os << " " << this->modelName() << ":" << nl
632  << " number of parcels added = " << parcelsAddedTotal_ << nl
633  << " mass introduced = " << massInjected_ << nl;
634 
635  if (this->writeTime())
636  {
637  this->setModelProperty("massInjected", massInjected_);
638  this->setModelProperty("nInjections", nInjections_);
639  this->setModelProperty("parcelsAddedTotal", parcelsAddedTotal_);
640  this->setModelProperty("timeStep0", timeStep0_);
641  }
642 }
643 
644 
645 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
646 
647 #include "InjectionModelNew.C"
648 
649 // ************************************************************************* //
scalar massTotal_
Total mass to inject [kg].
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
void inject(TrackData &td)
Main injection loop.
scalar time0_
Continuous phase time at start of injection time step [s].
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Templated injection model class.
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:824
CloudType::parcelType parcelType
Convenience typedef for parcelType.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:418
label parcelsAddedTotal_
Running counter of total number of parcels added.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
scalar SOI_
Start of injection [s].
Base class for cloud sub-models.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, parcelType &parcel)=0
Set the parcel properties.
virtual scalar timeEnd() const =0
Return the end-of-injection time.
virtual bool writeTime() const
Flag to indicate when to write a property.
const word & modelName() const
Return const access to the name of the sub-model.
Definition: subModelBase.C:104
const dictionary & dict() const
Return const access to the cloud dictionary.
Definition: subModelBase.C:110
virtual ~InjectionModel()
Destructor.
parcelBasis parcelBasis_
Parcel basis enumeration.
const CloudType & owner() const
Return const access to the owner cloud.
void injectSteadyState(TrackData &td, const scalar trackTime)
Main injection loop - steady-state.
dimensionedScalar pos(const dimensionedScalar &ds)
scalar nParticleFixed_
nParticle to assign to parcels when the &#39;fixed&#39; basis
stressControl lookup("compactNormalStress") >> compactNormalStress
label nInjections_
Number of injections counter.
dynamicFvMesh & mesh
void setModelProperty(const word &entryName, const Type &value)
Add generic property to the sub-model.
mathematical constants.
A class for handling words, derived from string.
Definition: word.H:59
virtual void info(Ostream &os)
Write injection info to stream.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
virtual bool validInjection(const label parcelI)=0
Additional flag to identify whether or not injection of parcelI is.
const Type & value() const
Return const reference to value.
static const zero Zero
Definition: zero.H:91
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
void reset(const dictionary &dict)
Reset entry by re-reading from dictionary.
Definition: TimeFunction1.C:82
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
Definition: subModelBase.C:128
virtual bool findCellAtPosition(label &celli, label &tetFacei, label &tetPti, vector &position, bool errorOnNotFound=true)
Find the cell that contains the supplied position.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
InjectionModel(CloudType &owner)
Construct null from owner.
static const char nl
Definition: Ostream.H:262
virtual void updateMesh()
Update mesh.
const word & modelType() const
Return const access to the sub-model type.
Definition: subModelBase.C:122
virtual Type value(const scalar x) const
Return value as a function of (scalar) independent variable.
const Time & time() const
Return time.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
virtual void setPositionAndCell(const label parcelI, const label nParcels, const scalar time, vector &position, label &cellOwner, label &tetFacei, label &tetPti)=0
Set the injection position and owner cell, tetFace and tetPt.
virtual bool prepareForNextTimeStep(const scalar time, label &newParcels, scalar &newVolumeFraction)
Determine properties for next time step/injection interval.
const fvMesh & mesh() const
Return refernce to the mesh.
Definition: DSMCCloudI.H:41
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.
scalar massInjected_
Total mass injected to date [kg].
scalar timeStep0_
Time at start of injection time step [s].
TimeFunction1< scalar > massFlowRate_
Mass flow rate profile for steady calculations.
messageStream Info
scalar delayedVolume_
Volume that should have been injected, but would lead to.
virtual bool fullyDescribed() const =0
Flag to identify whether model fully describes the parcel.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
virtual label parcelsToInject(const scalar time0, const scalar time1)=0
Number of parcels to introduce relative to SOI.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
scalar volumeTotal_
Total volume of particles introduced by this injector [m^3].
void constrainDirection(const polyMesh &mesh, const Vector< label > &dirs, vector &d)
Set the constrained components of directions/velocity to zero.
Definition: meshTools.C:671
scalar timeStart() const
Return the start-of-injection time.
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:331
virtual bool active() const
Return the model &#39;active&#39; status - default active = true.
Definition: subModelBase.C:154
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
void constrainToMeshCentre(const polyMesh &mesh, point &pt)
Set the constrained components of position to mesh centre.
Definition: meshTools.C:613
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576