37template<
class CloudType>
42 scalar& newVolumeFraction
47 newVolumeFraction = 0.0;
59 scalar t1 = time -
SOI_;
69 if (newVolumeFraction > 0)
89 return validInjection;
93template<
class CloudType>
107 this->owner().mesh().findCellFacePt
137 celli = this->owner().mesh().findNearestCell(position);
141 position += SMALL*(cellCentres[celli] - position);
143 this->owner().mesh().findCellFacePt
172 <<
"Cannot find parcel injection cell. "
173 <<
"Parcel position = " <<
p0 <<
nl
186template<
class CloudType>
190 const scalar volumeFraction,
191 const scalar diameter,
196 switch (parcelBasis_)
203 nP = (volumeFraction*volumeTot +
delayedVolume_)/(parcels*volumep);
208 nP = massTotal_/(
rho*volumeTotal_);
220 <<
"Unknown parcelBasis type" <<
nl
229template<
class CloudType>
232 const label parcelsAdded,
233 const scalar massAdded
238 if (allParcelsAdded > 0)
241 <<
"Cloud: " << this->
owner().name()
243 <<
" Added " << allParcelsAdded <<
" new parcels" <<
nl <<
endl;
247 parcelsAddedTotal_ += allParcelsAdded;
253 time0_ = this->owner().db().time().value();
262template<
class CloudType>
267 volumeTotal_(this->template getModelProperty<scalar>(
"volumeTotal")),
269 massFlowRate_(nullptr),
270 massInjected_(this->template getModelProperty<scalar>(
"massInjected")),
271 nInjections_(this->template getModelProperty<label>(
"nInjections")),
287template<
class CloudType>
292 const word& modelName,
293 const word& modelType
298 volumeTotal_(this->template getModelProperty<scalar>(
"volumeTotal")),
300 massFlowRate_(nullptr),
301 massInjected_(this->template getModelProperty<scalar>(
"massInjected")),
302 nInjections_(this->template getModelProperty<scalar>(
"nInjections")),
305 this->template getModelProperty<scalar>(
"parcelsAddedTotal")
307 parcelBasis_(pbNumber),
308 nParticleFixed_(0.0),
309 time0_(owner.db().time().value()),
310 timeStep0_(this->template getModelProperty<scalar>(
"timeStep0")),
311 minParticlesPerParcel_
313 this->coeffDict().getOrDefault(
"minParticlesPerParcel", scalar(1))
316 injectorID_(this->coeffDict().getOrDefault(
"injectorID", -1)),
319 this->coeffDict().getOrDefault(
"ignoreOutOfBounds", false)
325 Info<<
" Constructing " <<
owner.
mesh().nGeometricD() <<
"-D injection"
333 if (
owner.solution().active())
335 if (owner.solution().transient())
337 this->coeffDict().readEntry(
"massTotal", massTotal_);
338 this->coeffDict().readEntry(
"SOI", SOI_);
344 Function1<scalar>::New
351 massFlowRate_->userTimeToTime(owner.db().time());
352 massTotal_ = massFlowRate_->value(owner.db().time().value());
353 this->coeffDict().readIfPresent(
"SOI", SOI_);
357 SOI_ = owner.
db().time().userTimeToTime(SOI_);
359 const word parcelBasisType(this->coeffDict().getWord(
"parcelBasisType"));
361 if (parcelBasisType ==
"mass")
363 parcelBasis_ = pbMass;
365 else if (parcelBasisType ==
"number")
367 parcelBasis_ = pbNumber;
369 else if (parcelBasisType ==
"fixed")
371 parcelBasis_ = pbFixed;
372 this->coeffDict().readEntry(
"nParticle", nParticleFixed_);
374 Info<<
" Choosing nParticle to be a fixed value, massTotal "
375 <<
"variable now does not determine anything."
381 <<
"parcelBasisType must be either 'number', 'mass' or 'fixed'"
387template<
class CloudType>
414template<
class CloudType>
419template<
class CloudType>
423 if (this->owner().
solution().transient())
436template<
class CloudType>
437template<
class TrackCloudType>
440 TrackCloudType&
cloud,
441 typename CloudType::parcelType::trackingData&
td
449 const scalar time = this->owner().
db().
time().
value();
452 label parcelsAdded = 0;
453 scalar massAdded = 0.0;
454 label newParcels = 0;
455 scalar newVolumeFraction = 0.0;
456 scalar delayedVolume = 0;
458 if (prepareForNextTimeStep(time, newParcels, newVolumeFraction))
460 const scalar trackTime = this->owner().solution().trackTime();
464 const scalar deltaT =
465 max(0.0,
min(trackTime,
min(time - SOI_, timeEnd() - time0_)));
471 for (label parcelI = 0; parcelI < newParcels; parcelI++)
476 scalar timeInj =
time0_ + padTime + deltaT*parcelI/newParcels;
500 const scalar dt = time - timeInj;
509 cloud.setParcelThermoProperties(*pPtr, dt);
538 massAdded += pPtr->nParticle()*pPtr->mass();
543 cloud.addParticle(pPtr);
552 delayedVolume += pPtr->nParticle()*pPtr->volume();
566template<
class CloudType>
567template<
class TrackCloudType>
570 TrackCloudType&
cloud,
571 typename CloudType::parcelType::trackingData&
td,
572 const scalar trackTime
580 const scalar time = this->owner().
db().
time().
value();
593 label parcelsAdded = 0;
594 scalar massAdded = 0.0;
597 label newParcels = parcelsToInject(0.0, 1.0);
600 for (label parcelI = 0; parcelI < newParcels; parcelI++)
603 scalar newVolumeFraction = 1.0/scalar(newParcels);
630 parcelType* pPtr =
new parcelType(
mesh,
pos, celli);
633 cloud.setParcelThermoProperties(*pPtr, 0.0);
636 setProperties(parcelI, newParcels, 0.0, *pPtr);
639 cloud.checkParcelProperties(*pPtr, 0.0, fullyDescribed());
654 pPtr->typeId() = injectorID_;
657 cloud.addParticle(pPtr);
659 massAdded += pPtr->nParticle()*pPtr->mass();
664 postInjectCheck(parcelsAdded, massAdded);
668template<
class CloudType>
673 Log_<<
" Injector " << this->modelName() <<
":" <<
nl
674 <<
" - parcels added = " << parcelsAddedTotal_ <<
nl
675 <<
" - mass introduced = " << massInjected_ <<
nl;
677 if (this->writeTime())
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_);
Base class for cloud sub-models.
const CloudType & owner() const
Return const access to the owner cloud.
virtual void info()
Write to info.
CloudSubModelBase(CloudType &owner)
Construct null from owner cloud.
virtual bool writeTime() const
Flag to indicate when to write a property.
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
const fvMesh & mesh() const
Return reference to the mesh.
label typeId() const
Return type id.
const vector & U() const
Return const access to velocity.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
virtual bool fullyDescribed() const=0
virtual bool prepareForNextTimeStep(const scalar time, label &newParcels, scalar &newVolumeFraction)
Determine properties for next time step/injection interval.
scalar timeStep0_
Time at start of injection time step [s].
virtual scalar timeEnd() const =0
Return the end-of-injection time.
void inject(TrackCloudType &cloud, typename CloudType::parcelType::trackingData &td)
Main injection loop.
parcelBasis parcelBasis_
Parcel basis enumeration.
virtual autoPtr< InjectionModel< Foam::KinematicCloud< Foam::DSMCCloud< dsmcParcel > > > > clone() const=0
scalar minParticlesPerParcel_
scalar delayedVolume_
Volume that should have been injected, but would lead to less than minParticlesPerParcel_ particle pe...
virtual scalar volumeToInject(const scalar time0, const scalar time1)=0
Volume of parcels to introduce relative to SOI.
virtual scalar setNumberOfParticles(const label parcels, const scalar volumeFraction, const scalar diameter, const scalar rho)
Set number of particles to inject given parcel properties.
virtual void setPositionAndCell(const label parcelI, const label nParcels, const scalar time, vector &position, label &cellOwner, label &tetFacei, label &tetPti)=0
scalar timeStart() const
Return the start-of-injection time.
virtual label parcelsToInject(const scalar time0, const scalar time1)=0
Number of parcels to introduce relative to SOI.
virtual scalar averageParcelMass()
Return the average parcel mass over the injection period.
scalar nParticleFixed_
nParticle to assign to parcels when the 'fixed' basis is selected
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, parcelType &parcel)=0
void injectSteadyState(TrackCloudType &cloud, typename CloudType::parcelType::trackingData &td, const scalar trackTime)
Main injection loop - steady-state.
virtual void info()
Write injection info.
virtual bool validInjection(const label parcelI)=0
Additional flag to identify whether or not injection of parcelI is.
label nInjections_
Number of injections counter.
scalar volumeTotal_
Total volume of particles introduced by this injector [m^3] Note: scaled to ensure massTotal is achie...
InjectionModel(CloudType &owner)
Construct null from owner.
virtual void updateMesh()
Update mesh.
scalar massTotal_
Total mass to inject [kg].
virtual bool findCellAtPosition(label &celli, label &tetFacei, label &tetPti, vector &position, bool errorOnNotFound=true)
Find the cell that contains the supplied position.
label parcelsAddedTotal_
Running counter of total number of parcels added.
Foam::KinematicCloud< Foam::DSMCCloud< dsmcParcel > >::parcelType parcelType
scalar time0_
Continuous phase time at start of injection time step [s].
scalar SOI_
Start of injection [s].
virtual void postInjectCheck(const label parcelsAdded, const scalar massAdded)
Post injection checks.
autoPtr< Function1< scalar > > massFlowRate_
Mass flow rate profile for steady calculations.
scalar massInjected_
Total mass injected to date [kg].
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
A cloud is a registry collection of lagrangian particles.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
const Type & value() const noexcept
Return const reference to value.
const Time & time() const
Return the top-level database.
const Time & time() const noexcept
Return time registry.
Mesh consisting of general polyhedral cells.
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Selector class for relaxation factors, solver type and solution.
bool getModelProperty(const word &entryName, Type &value) const
Retrieve generic property from the sub-model.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
const word & modelName() const
Return const access to the name of the sub-model.
const word & modelType() const
Return const access to the sub-model type.
virtual bool active() const
Return the model 'active' status - default active = true.
void setModelProperty(const word &entryName, const Type &value)
Add generic property to the sub-model.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
const volScalarField & p0
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
return returnReduce(nRefine-oldNRefine, sumOp< label >())
#define Log_
Report write to Foam::Info if the class log switch is true.
constexpr scalar pi(M_PI)
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
DSMCCloud< dsmcParcel > CloudType
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedScalar pow3(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
constexpr char nl
The newline '\n' character (0x0a).