34template<
class CloudType>
38 { vtFixedValue,
"fixedValue" },
40 { vtZeroGradient,
"zeroGradient" },
45template<
class CloudType>
54 patchInjectionBase(owner.
mesh(), this->coeffDict().getWord(
"patch")),
55 duration_(this->coeffDict().getScalar(
"duration")),
58 this->coeffDict().getScalar(
"parcelsPerSecond")
62 velocityTypeNames_.getOrDefault
71 velocityType_ == vtFixedValue
72 ? this->coeffDict().template
get<
vector>(
"U0")
77 Function1<scalar>::
New
86 distributionModel::
New
88 this->coeffDict().subDict(
"sizeDistribution"),
96 const Time& time =
owner.db().time();
97 duration_ = time.userTimeToTime(duration_);
98 flowRateProfile_->userTimeToTime(time);
103 this->
volumeTotal_ = flowRateProfile_->integrate(0.0, duration_);
107template<
class CloudType>
115 duration_(im.duration_),
116 parcelsPerSecond_(im.parcelsPerSecond_),
117 velocityType_(im.velocityType_),
119 flowRateProfile_(im.flowRateProfile_.clone()),
120 sizeDistribution_(im.sizeDistribution_.clone()),
121 currentParceli_(im.currentParceli_),
122 currentFacei_(im.currentFacei_)
128template<
class CloudType>
135template<
class CloudType>
138 return this->
SOI_ + duration_;
142template<
class CloudType>
149 if ((time0 >= 0.0) && (time0 < duration_))
151 scalar nParcels = (time1 - time0)*parcelsPerSecond_;
152 Random& rnd = this->owner().rndGen();
154 label nParcelsToInject = floor(nParcels);
161 && (nParcels - scalar(nParcelsToInject) > rndPos)
167 return nParcelsToInject;
174template<
class CloudType>
181 if ((time0 >= 0.0) && (time0 < duration_))
183 return flowRateProfile_->integrate(time0, time1);
190template<
class CloudType>
194 const label nParcels,
202 currentParceli_ = parcelI;
206 this->owner().
mesh(),
216template<
class CloudType>
220 const label nParcels,
226 switch (velocityType_)
235 if (parcelI != currentParceli_)
238 <<
"Synchronisation problem: "
239 <<
"attempting to set injected parcel " << parcelI
240 <<
" properties using cached parcel " << currentParceli_
241 <<
" properties" <<
endl;
244 const label patchFacei = currentFacei_;
249 <<
"Unable to set parcel velocity using patch value "
250 <<
"due to missing face index: patchFacei=" << patchFacei
255 const label patchi = this->patchId_;
256 parcel.
U() =
U.boundaryField()[patchi][patchFacei];
261 const label celli = parcel.
cell();
266 <<
"Unable to set parcel velocity using zeroGradient "
267 <<
"due to missing cell index"
272 parcel.
U() =
U[celli];
278 <<
"Unhandled velocityType "
279 << velocityTypeNames_[velocityType_]
280 <<
". Available options are:"
281 << velocityTypeNames_.sortedToc()
287 parcel.d() = sizeDistribution_->sample();
291template<
class CloudType>
298template<
class CloudType>
const CloudType & owner() const
Return const access to the owner cloud.
const vector & U() const
Return const access to velocity.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Templated injection model class.
static autoPtr< InjectionModel< CloudType > > New(const dictionary &dict, CloudType &owner)
Selector with lookup from dictionary.
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.
scalar SOI_
Start of injection [s].
@ vtFixedValue
User supplied fixed value.
@ vtPatchValue
Patch face values.
@ vtZeroGradient
Patch internal cell values.
virtual autoPtr< InjectionModel< CloudType > > clone() const
Construct and return a clone.
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
virtual void setPositionAndCell(const label parcelI, const label nParcels, const scalar time, vector &position, label &cellOwner, label &tetFacei, label &tetPti)
Set the injection position and owner cell, tetFace and tetPt.
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
PatchInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
static const Enum< velocityType > velocityTypeNames_
Velocity type names.
virtual void updateMesh()
Set injector locations when mesh is updated.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
virtual scalar timeEnd() const
Return the end-of-injection time.
Type globalPosition(const Type &start, const Type &end)
Return a sample on the interval [start,end].
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (e.g. CA deg) to real-time (s).
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
label cell() const noexcept
Return current cell particle is in.
const label patchId_
Patch ID.
patchInjectionBase(const polyMesh &mesh, const word &patchName)
Construct from mesh and patch name.
virtual void updateMesh(const polyMesh &mesh)
Update patch geometry and derived info for injection locations.
label setPositionAndCell(const fvMesh &mesh, const scalar fraction01, Random &rnd, vector &position, label &cellOwner, label &tetFacei, label &tetPti)
Set the injection position and owner cell, tetFace and tetPt.
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.
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.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define WarningInFunction
Report a warning using Foam::Warning.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
DSMCCloud< dsmcParcel > CloudType
GeometricField< vector, fvPatchField, volMesh > volVectorField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
errorManip< error > abort(error &err)
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...