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-2016 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 {}
272 
273 
274 template<class CloudType>
276 (
277  const dictionary& dict,
278  CloudType& owner,
279  const word& modelName,
280  const word& modelType
281 )
282 :
284  SOI_(0.0),
285  volumeTotal_(0.0),
286  massTotal_(0.0),
287  massFlowRate_(owner.db().time(), "massFlowRate"),
288  massInjected_(this->template getModelProperty<scalar>("massInjected")),
289  nInjections_(this->template getModelProperty<scalar>("nInjections")),
291  (
292  this->template getModelProperty<scalar>("parcelsAddedTotal")
293  ),
295  nParticleFixed_(0.0),
296  time0_(owner.db().time().value()),
297  timeStep0_(this->template getModelProperty<scalar>("timeStep0")),
298  delayedVolume_(0.0)
299 {
300  // Provide some info
301  // - also serves to initialise mesh dimensions - needed for parallel runs
302  // due to lazy evaluation of valid mesh dimensions
303  Info<< " Constructing " << owner.mesh().nGeometricD() << "-D injection"
304  << endl;
305 
306  if (owner.solution().transient())
307  {
308  this->coeffDict().lookup("massTotal") >> massTotal_;
309  this->coeffDict().lookup("SOI") >> SOI_;
310  SOI_ = owner.db().time().userTimeToTime(SOI_);
311  }
312  else
313  {
314  massFlowRate_.reset(this->coeffDict());
315  massTotal_ = massFlowRate_.value(owner.db().time().value());
316  }
317 
318  const word parcelBasisType = this->coeffDict().lookup("parcelBasisType");
319 
320  if (parcelBasisType == "mass")
321  {
323  }
324  else if (parcelBasisType == "number")
325  {
327  }
328  else if (parcelBasisType == "fixed")
329  {
331 
332  Info<< " Choosing nParticle to be a fixed value, massTotal "
333  << "variable now does not determine anything."
334  << endl;
335 
336  nParticleFixed_ = readScalar(this->coeffDict().lookup("nParticle"));
337  }
338  else
339  {
341  << "parcelBasisType must be either 'number', 'mass' or 'fixed'" << nl
342  << exit(FatalError);
343  }
344 }
345 
346 
347 template<class CloudType>
349 (
350  const InjectionModel<CloudType>& im
351 )
352 :
354  SOI_(im.SOI_),
363  time0_(im.time0_),
366 {}
367 
368 
369 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
370 
371 template<class CloudType>
373 {}
374 
375 
376 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
377 
378 template<class CloudType>
380 {
381  // do nothing
382 }
383 
384 
385 template<class CloudType>
387 {
388  label nTotal = 0.0;
389  if (this->owner().solution().transient())
390  {
391  nTotal = parcelsToInject(0.0, timeEnd() - timeStart());
392  }
393  else
394  {
395  nTotal = parcelsToInject(0.0, 1.0);
396  }
397 
398  return massTotal_/nTotal;
399 }
400 
401 
402 template<class CloudType>
403 template<class TrackData>
405 {
406  if (!this->active())
407  {
408  return;
409  }
410 
411  const scalar time = this->owner().db().time().value();
412 
413  // Prepare for next time step
414  label parcelsAdded = 0;
415  scalar massAdded = 0.0;
416  label newParcels = 0;
417  scalar newVolumeFraction = 0.0;
418 
419  if (prepareForNextTimeStep(time, newParcels, newVolumeFraction))
420  {
421  scalar delayedVolume = 0;
422 
423  const scalar trackTime = this->owner().solution().trackTime();
424  const polyMesh& mesh = this->owner().mesh();
425  typename TrackData::cloudType& cloud = td.cloud();
426 
427  // Duration of injection period during this timestep
428  const scalar deltaT =
429  max(0.0, min(trackTime, min(time - SOI_, timeEnd() - time0_)));
430 
431  // Pad injection time if injection starts during this timestep
432  const scalar padTime = max(0.0, SOI_ - time0_);
433 
434  // Introduce new parcels linearly across carrier phase timestep
435  for (label parcelI = 0; parcelI < newParcels; parcelI++)
436  {
437  if (validInjection(parcelI))
438  {
439  // Calculate the pseudo time of injection for parcel 'parcelI'
440  scalar timeInj = time0_ + padTime + deltaT*parcelI/newParcels;
441 
442  // Determine the injection position and owner cell,
443  // tetFace and tetPt
444  label celli = -1;
445  label tetFacei = -1;
446  label tetPtI = -1;
447 
448  vector pos = Zero;
449 
451  (
452  parcelI,
453  newParcels,
454  timeInj,
455  pos,
456  celli,
457  tetFacei,
458  tetPtI
459  );
460 
461  if (celli > -1)
462  {
463  // Lagrangian timestep
464  const scalar dt = time - timeInj;
465 
466  // Apply corrections to position for 2-D cases
468 
469  // Create a new parcel
470  parcelType* pPtr =
471  new parcelType(mesh, pos, celli, tetFacei, tetPtI);
472 
473  // Check/set new parcel thermo properties
474  cloud.setParcelThermoProperties(*pPtr, dt);
475 
476  // Assign new parcel properties in injection model
477  setProperties(parcelI, newParcels, timeInj, *pPtr);
478 
479  // Check/set new parcel injection properties
480  cloud.checkParcelProperties(*pPtr, dt, fullyDescribed());
481 
482  // Apply correction to velocity for 2-D cases
484  (
485  mesh,
486  mesh.solutionD(),
487  pPtr->U()
488  );
489 
490  // Number of particles per parcel
491  pPtr->nParticle() =
493  (
494  newParcels,
495  newVolumeFraction,
496  pPtr->d(),
497  pPtr->rho()
498  );
499 
500  if (pPtr->nParticle() >= 1.0)
501  {
502  parcelsAdded++;
503  massAdded += pPtr->nParticle()*pPtr->mass();
504 
505  if (pPtr->move(td, dt))
506  {
507  td.cloud().addParticle(pPtr);
508  }
509  else
510  {
511  delete pPtr;
512  }
513  }
514  else
515  {
516  delayedVolume += pPtr->nParticle()*pPtr->volume();
517  delete pPtr;
518  }
519  }
520  }
521  }
522 
523  delayedVolume_ = delayedVolume;
524  }
525 
526  postInjectCheck(parcelsAdded, massAdded);
527 }
528 
529 
530 template<class CloudType>
531 template<class TrackData>
533 (
534  TrackData& td,
535  const scalar trackTime
536 )
537 {
538  if (!this->active())
539  {
540  return;
541  }
542 
543  const polyMesh& mesh = this->owner().mesh();
544  typename TrackData::cloudType& cloud = td.cloud();
545 
547 
548  // Reset counters
549  time0_ = 0.0;
550  label parcelsAdded = 0;
551  scalar massAdded = 0.0;
552 
553  // Set number of new parcels to inject based on first second of injection
554  label newParcels = parcelsToInject(0.0, 1.0);
555 
556  // Inject new parcels
557  for (label parcelI = 0; parcelI < newParcels; parcelI++)
558  {
559  // Volume to inject is split equally amongst all parcel streams
560  scalar newVolumeFraction = 1.0/scalar(newParcels);
561 
562  // Determine the injection position and owner cell,
563  // tetFace and tetPt
564  label celli = -1;
565  label tetFacei = -1;
566  label tetPtI = -1;
567 
568  vector pos = Zero;
569 
571  (
572  parcelI,
573  newParcels,
574  0.0,
575  pos,
576  celli,
577  tetFacei,
578  tetPtI
579  );
580 
581  if (celli > -1)
582  {
583  // Apply corrections to position for 2-D cases
585 
586  // Create a new parcel
587  parcelType* pPtr =
588  new parcelType(mesh, pos, celli, tetFacei, tetPtI);
589 
590  // Check/set new parcel thermo properties
591  cloud.setParcelThermoProperties(*pPtr, 0.0);
592 
593  // Assign new parcel properties in injection model
594  setProperties(parcelI, newParcels, 0.0, *pPtr);
595 
596  // Check/set new parcel injection properties
597  cloud.checkParcelProperties(*pPtr, 0.0, fullyDescribed());
598 
599  // Apply correction to velocity for 2-D cases
600  meshTools::constrainDirection(mesh, mesh.solutionD(), pPtr->U());
601 
602  // Number of particles per parcel
603  pPtr->nParticle() =
605  (
606  1,
607  newVolumeFraction,
608  pPtr->d(),
609  pPtr->rho()
610  );
611 
612  // Add the new parcel
613  td.cloud().addParticle(pPtr);
614 
615  massAdded += pPtr->nParticle()*pPtr->mass();
616  parcelsAdded++;
617  }
618  }
619 
620  postInjectCheck(parcelsAdded, massAdded);
621 }
622 
623 
624 template<class CloudType>
626 {
627  os << " " << this->modelName() << ":" << nl
628  << " number of parcels added = " << parcelsAddedTotal_ << nl
629  << " mass introduced = " << massInjected_ << nl;
630 
631  if (this->writeTime())
632  {
633  this->setModelProperty("massInjected", massInjected_);
634  this->setModelProperty("nInjections", nInjections_);
635  this->setModelProperty("parcelsAddedTotal", parcelsAddedTotal_);
636  this->setModelProperty("timeStep0", timeStep0_);
637  }
638 }
639 
640 
641 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
642 
643 #include "InjectionModelNew.C"
644 
645 // ************************************************************************* //
const Time & time() const
Return time.
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
const dictionary & dict() const
Return const access to the cloud dictionary.
Definition: subModelBase.C:110
void inject(TrackData &td)
Main injection loop.
scalar time0_
Continuous phase time at start of injection time step [s].
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:221
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
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:801
Templated injection model class.
const word & modelName() const
Return const access to the name of the sub-model.
Definition: subModelBase.C:104
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:417
label parcelsAddedTotal_
Running counter of total number of parcels added.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
virtual bool findCellAtPosition(label &celli, label &tetFacei, label &tetPtI, vector &position, bool errorOnNotFound=true)
Find the cell that contains the supplied position.
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.
scalar SOI_
Start of injection [s].
const Type & value() const
Return const reference to value.
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.
const word & modelType() const
Return const access to the sub-model type.
Definition: subModelBase.C:122
virtual ~InjectionModel()
Destructor.
virtual bool writeTime() const
Flag to indicate when to write a property.
parcelBasis parcelBasis_
Parcel basis enumeration.
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.
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
virtual Type value(const scalar x) const
Return value as a function of (scalar) independent variable.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
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.
scalar timeStart() const
Return the start-of-injection 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 bool prepareForNextTimeStep(const scalar time, label &newParcels, scalar &newVolumeFraction)
Determine properties for next time step/injection interval.
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].
const fvMesh & mesh() const
Return refernce to the mesh.
Definition: DSMCCloudI.H:41
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.
const CloudType & owner() const
Return const access to the owner cloud.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
Definition: subModelBase.C:128
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
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:68
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:451