Go to the documentation of this file.
35 template<
class CloudType>
40 scalar& newVolumeFraction
45 newVolumeFraction = 0.0;
46 bool validInjection =
false;
52 return validInjection;
56 scalar t0 = timeStep0_ - SOI_;
57 scalar t1 = time - SOI_;
60 newParcels = this->parcelsToInject(t0, t1);
64 this->volumeToInject(t0, t1)
65 /(volumeTotal_ + ROOTVSMALL);
67 if (newVolumeFraction > 0)
72 validInjection =
true;
78 validInjection =
false;
84 validInjection =
false;
87 return validInjection;
91 template<
class CloudType>
103 const vector p0 = position;
139 position += SMALL*(cellCentres[cellI] - position);
162 <<
"Cannot find parcel injection cell. "
163 <<
"Parcel position = " << p0 <<
nl
176 template<
class CloudType>
180 const scalar volumeFraction,
181 const scalar diameter,
186 switch (parcelBasis_)
190 scalar volumep =
pi/6.0*
pow3(diameter);
191 scalar volumeTot = massTotal_/
rho;
193 nP = (volumeFraction*volumeTot + delayedVolume_)/(parcels*volumep);
198 nP = massTotal_/(
rho*volumeTotal_);
203 nP = nParticleFixed_;
210 <<
"Unknown parcelBasis type" <<
nl
219 template<
class CloudType>
222 const label parcelsAdded,
223 const scalar massAdded
228 if (allParcelsAdded > 0)
231 <<
"Cloud: " << this->owner().name()
232 <<
" injector: " << this->modelName() <<
nl
233 <<
" Added " << allParcelsAdded <<
" new parcels" <<
nl <<
endl;
237 parcelsAddedTotal_ += allParcelsAdded;
243 time0_ = this->owner().db().time().value();
252 template<
class CloudType>
259 massFlowRate_(owner.db().time(),
"massFlowRate"),
260 massInjected_(this->template getModelProperty<scalar>(
"massInjected")),
261 nInjections_(this->template getModelProperty<
label>(
"nInjections")),
264 this->template getModelProperty<scalar>(
"parcelsAddedTotal")
266 parcelBasis_(pbNumber),
267 nParticleFixed_(0.0),
269 timeStep0_(this->template getModelProperty<scalar>(
"timeStep0")),
275 template<
class CloudType>
280 const word& modelName,
281 const word& modelType
288 massFlowRate_(owner.db().time(),
"massFlowRate"),
289 massInjected_(this->
template getModelProperty<scalar>(
"massInjected")),
290 nInjections_(this->
template getModelProperty<scalar>(
"nInjections")),
293 this->
template getModelProperty<scalar>(
"parcelsAddedTotal")
295 parcelBasis_(pbNumber),
296 nParticleFixed_(0.0),
297 time0_(owner.db().time().value()),
298 timeStep0_(this->
template getModelProperty<scalar>(
"timeStep0")),
300 injectorID_(this->coeffDict().lookupOrDefault(
"injectorID", -1))
308 if (injectorID_ != -1)
310 Info<<
" injector ID: " << injectorID_ <<
endl;
313 if (owner.solution().transient())
315 this->coeffDict().lookup(
"massTotal") >> massTotal_;
316 this->coeffDict().lookup(
"SOI") >> SOI_;
320 massFlowRate_.reset(this->coeffDict());
321 massTotal_ = massFlowRate_.value(owner.db().time().value());
322 this->coeffDict().readIfPresent(
"SOI", SOI_);
325 SOI_ = owner.db().time().userTimeToTime(SOI_);
327 const word parcelBasisType = this->coeffDict().lookup(
"parcelBasisType");
329 if (parcelBasisType ==
"mass")
331 parcelBasis_ = pbMass;
333 else if (parcelBasisType ==
"number")
335 parcelBasis_ = pbNumber;
337 else if (parcelBasisType ==
"fixed")
339 parcelBasis_ = pbFixed;
341 Info<<
" Choosing nParticle to be a fixed value, massTotal "
342 <<
"variable now does not determine anything."
350 <<
"parcelBasisType must be either 'number', 'mass' or 'fixed'" <<
nl
356 template<
class CloudType>
381 template<
class CloudType>
388 template<
class CloudType>
395 template<
class CloudType>
399 if (this->owner().
solution().
transient())
401 nTotal = parcelsToInject(0.0, timeEnd() - timeStart());
405 nTotal = parcelsToInject(0.0, 1.0);
408 return massTotal_/nTotal;
412 template<
class CloudType>
413 template<
class TrackData>
421 const scalar time = this->owner().db().time().value();
424 label parcelsAdded = 0;
425 scalar massAdded = 0.0;
426 label newParcels = 0;
427 scalar newVolumeFraction = 0.0;
428 scalar delayedVolume = 0;
430 if (prepareForNextTimeStep(time, newParcels, newVolumeFraction))
433 const scalar trackTime = this->owner().solution().trackTime();
435 typename TrackData::cloudType&
cloud = td.
cloud();
438 const scalar deltaT =
439 max(0.0,
min(trackTime,
min(time - SOI_, timeEnd() - time0_)));
442 const scalar padTime =
max(0.0, SOI_ - time0_);
445 for (
label parcelI = 0; parcelI < newParcels; parcelI++)
447 if (validInjection(parcelI))
450 scalar timeInj = time0_ + padTime + deltaT*parcelI/newParcels;
474 const scalar dt = time - timeInj;
484 cloud.setParcelThermoProperties(*pPtr, dt);
487 setProperties(parcelI, newParcels, timeInj, *pPtr);
490 cloud.checkParcelProperties(*pPtr, dt, fullyDescribed());
510 if (pPtr->nParticle() >= 1.0)
513 massAdded += pPtr->nParticle()*pPtr->mass();
515 if (pPtr->move(td, dt))
517 pPtr->typeId() = injectorID_;
518 cloud.addParticle(pPtr);
527 delayedVolume += pPtr->nParticle()*pPtr->volume();
537 postInjectCheck(parcelsAdded, massAdded);
541 template<
class CloudType>
542 template<
class TrackData>
546 const scalar trackTime
554 const scalar time = this->owner().db().time().value();
562 typename TrackData::cloudType&
cloud = td.
cloud();
568 label parcelsAdded = 0;
569 scalar massAdded = 0.0;
572 label newParcels = parcelsToInject(0.0, 1.0);
575 for (
label parcelI = 0; parcelI < newParcels; parcelI++)
578 scalar newVolumeFraction = 1.0/scalar(newParcels);
609 cloud.setParcelThermoProperties(*pPtr, 0.0);
612 setProperties(parcelI, newParcels, 0.0, *pPtr);
615 cloud.checkParcelProperties(*pPtr, 0.0, fullyDescribed());
630 pPtr->typeId() = injectorID_;
633 cloud.addParticle(pPtr);
635 massAdded += pPtr->nParticle()*pPtr->mass();
640 postInjectCheck(parcelsAdded, massAdded);
644 template<
class CloudType>
647 os <<
" Injector " << this->modelName() <<
":" <<
nl
648 <<
" - parcels added = " << parcelsAddedTotal_ <<
nl
649 <<
" - mass introduced = " << massInjected_ <<
nl;
651 if (this->outputTime())
653 this->setModelProperty(
"massInjected", massInjected_);
654 this->setModelProperty(
"nInjections", nInjections_);
655 this->setModelProperty(
"parcelsAddedTotal", parcelsAddedTotal_);
656 this->setModelProperty(
"timeStep0", timeStep0_);
virtual void info(Ostream &os)
Write injection info to stream.
Selector class for relaxation factors, solver type and solution.
A class for handling words, derived from string.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
InjectionModel(CloudType &owner)
Construct null from owner.
virtual void updateMesh()
Update mesh.
Templated injection model class.
Base class for cloud sub-models.
virtual bool prepareForNextTimeStep(const scalar time, label &newParcels, scalar &newVolumeFraction)
Determine properties for next time step/injection interval.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type & value() const
Return const reference to value.
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
scalar massInjected_
Total mass injected to date [kg].
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 delayedVolume_
Volume that should have been injected, but would lead to.
scalar volumeTotal_
Total volume of particles introduced by this injector [m^3].
bool pointInCell(const point &p, label celli, const cellDecomposition=CELL_TETS) const
Test if point p is in the celli.
scalar timeStep0_
Time at start of injection time step [s].
label injectorID_
Optional injector ID.
Mesh consisting of general polyhedral cells.
TimeDataEntry< scalar > massFlowRate_
Mass flow rate profile for steady calculations.
label findNearestCell(const point &location) const
Find the cell with the nearest cell centre to location.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
dimensionedScalar pow3(const dimensionedScalar &ds)
const fvMesh & mesh() const
Return refernce to the mesh.
scalar time0_
Continuous phase time at start of injection time step [s].
Templated base class for dsmc cloud.
virtual void postInjectCheck(const label parcelsAdded, const scalar massAdded)
Post injection checks.
virtual scalar averageParcelMass()
Return the average parcel mass over the injection period.
A list of keyword definitions, which are a keyword followed by any number of values (e....
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
scalar massTotal_
Total mass to inject [kg].
errorManip< error > abort(error &err)
void inject(TrackData &td)
Main injection loop.
label parcelsAddedTotal_
Running counter of total number of parcels added.
errorManipArg< error, int > exit(error &err, const int errNo=1)
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
A cloud is a collection of lagrangian particles.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const volVectorField & C() const
Return cell centres as volVectorField.
virtual bool findCellAtPosition(label &cellI, label &tetFaceI, label &tetPtI, vector &position, bool errorOnNotFound=true)
Find the cell that contains the supplied position.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
CloudType::parcelType parcelType
Convenience typedef for parcelType.
virtual ~InjectionModel()
Destructor.
scalar SOI_
Start of injection [s].
const Time & time() const
Return the top-level database.
label nInjections_
Number of injections counter.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
parcelBasis parcelBasis_
Parcel basis enumeration.
Generic GeometricField class.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
void injectSteadyState(TrackData &td, const scalar trackTime)
Main injection loop - steady-state.
cloud(const cloud &)
Disallow default bitwise copy construct.
scalar nParticleFixed_
nParticle to assign to parcels when the 'fixed' basis
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar pos(const dimensionedScalar &ds)