33 template<
class CloudType>
43 phiName_(this->coeffDict().template lookupOrDefault<
word>(
"phi",
"phi")),
44 rhoName_(this->coeffDict().template lookupOrDefault<
word>(
"rho",
"rho")),
45 duration_(this->readDuration(
dict, owner)),
51 this->owner().db().time().userUnits(),
58 this->coeffDict().template lookup<scalar>(
"parcelConcentration")
65 this->coeffDict().subDict(
"sizeDistribution"),
73 template<
class CloudType>
81 phiName_(im.phiName_),
82 rhoName_(im.rhoName_),
83 duration_(im.duration_),
84 concentration_(im.concentration_, false),
85 parcelConcentration_(im.parcelConcentration_),
86 sizeDistribution_(im.sizeDistribution_, false)
92 template<
class CloudType>
99 template<
class CloudType>
106 template<
class CloudType>
109 return this->SOI_ + duration_;
113 template<
class CloudType>
116 const polyMesh& mesh = this->owner().mesh();
123 scalar flowRateIn = 0.0;
126 flowRateIn =
max(0.0, -
sum(phip));
134 flowRateIn =
max(0.0, -
sum(phip/rhop));
143 template<
class CloudType>
150 if (time0 >= 0 && time0 < duration_)
152 scalar dt = time1 - time0;
154 scalar
c = concentration_->
value(0.5*(time0 + time1));
156 scalar nParcels = parcelConcentration_*
c*flowRate()*dt;
160 label nParcelsToInject = floor(nParcels);
164 if (nParcels - scalar(nParcelsToInject) > this->globalScalar01(
rndGen))
169 return nParcelsToInject;
178 template<
class CloudType>
187 if (time0 >= 0 && time0 < duration_)
189 scalar
c = concentration_->
value(0.5*(time0 + time1));
191 volume =
c*(time1 - time0)*flowRate();
194 return volume*this->owner().constProps().rho0();
198 template<
class CloudType>
213 this->owner().mesh(),
224 template<
class CloudType>
230 typename CloudType::parcelType::trackingData& td,
235 parcel.U() = this->owner().U()[parcel.cell()];
238 parcel.d() = sizeDistribution_->sample();
242 template<
class CloudType>
Templated base class for dsmc cloud.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
const dimensionSet & dimensions() const
Return dimensions.
Run-time selectable general function of one variable.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Templated injection model class.
Patch injection, by using patch flow rate to determine concentration and velocity.
virtual void topoChange()
Set injector locations when mesh is updated.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType::trackingData &td, typename CloudType::parcelType &parcel)
virtual ~PatchFlowRateInjection()
Destructor.
PatchFlowRateInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
scalar flowRate() const
Return the total volumetric flow rate across the patch [m^3/s].
virtual label nParcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
virtual void setPositionAndCell(const fvMesh &mesh, randomGenerator &rndGen, barycentric &coordinates, label &celli, label &tetFacei, label &tetPti, label &facei)
Inherit setPositionAndCell from patchInjectionBase.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
virtual scalar timeEnd() const
Return the end-of-injection time.
virtual scalar massToInject(const scalar time0, const scalar time1)
Parcel mass to introduce relative to SOI.
A list of keyword definitions, which are a keyword followed by any number of values (e....
const Type & value() const
Return const reference to value.
Base class for statistical distributions.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
Base class for patch-based injection models.
virtual void topoChange(const polyMesh &mesh)
Update patch geometry and derived info for injection locations.
virtual void setPositionAndCell(const fvMesh &mesh, randomGenerator &rndGen, barycentric &coordinates, label &celli, label &tetFacei, label &tetPti, label &facei)
Set the injection position and owner cell, tetFace and tetPt.
Mesh consisting of general polyhedral cells.
A class for handling words, derived from string.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar c
Speed of light in a vacuum.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const dimensionSet dimVolumetricFlux
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const dimensionSet dimless
const dimensionSet dimLength
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
randomGenerator rndGen(653213)