Go to the documentation of this file.
34 template<
class CloudType>
44 phiName_(this->coeffDict().
template lookupOrDefault<word>(
"phi",
"phi")),
45 rhoName_(this->coeffDict().
template lookupOrDefault<word>(
"rho",
"rho")),
64 this->coeffDict().subDict(
"sizeDistribution"),
69 duration_ = owner.db().time().userTimeToTime(duration_);
71 patchInjectionBase::updateMesh(owner.
mesh());
75 this->volumeTotal_ = 0.0;
76 this->massTotal_ = 0.0;
80 template<
class CloudType>
99 template<
class CloudType>
106 template<
class CloudType>
109 patchInjectionBase::updateMesh(this->owner().
mesh());
113 template<
class CloudType>
116 return this->SOI_ + duration_;
120 template<
class CloudType>
130 scalar flowRateIn = 0.0;
133 flowRateIn =
max(0.0, -
sum(phip));
141 flowRateIn =
max(0.0, -
sum(phip/rhop));
150 template<
class CloudType>
157 if ((time0 >= 0.0) && (time0 < duration_))
159 scalar dt = time1 - time0;
161 scalar
c = concentration_.
value(0.5*(time0 + time1));
163 scalar nParcels = parcelConcentration_*
c*flowRate()*dt;
164 label nParcelsToInject = floor(nParcels);
172 nParcels - scalar(nParcelsToInject)
173 > this->owner().
rndGen().position(scalar(0), scalar(1))
180 return nParcelsToInject;
189 template<
class CloudType>
198 if ((time0 >= 0.0) && (time0 < duration_))
200 scalar
c = concentration_.
value(0.5*(time0 + time1));
202 volume =
c*(time1 - time0)*flowRate();
205 this->volumeTotal_ = volume;
206 this->massTotal_ = volume*this->owner().
constProps().rho0();
212 template<
class CloudType>
224 patchInjectionBase::setPositionAndCell
226 this->owner().
mesh(),
236 template<
class CloudType>
246 parcel.U() = this->owner().U()[parcel.cell()];
249 parcel.d() = sizeDistribution_->sample();
253 template<
class CloudType>
260 template<
class CloudType>
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 scalar flowRate() const
Return the total volumetric flow rate across the patch [m3/s].
A class for handling words, derived from string.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
const List< typename ParcelType::constantProperties > & constProps() const
Return all of the constant properties.
const dimensionSet dimVelocity
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
const autoPtr< distributionModels::distributionModel > sizeDistribution_
Parcel size distribution model.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
Templated injection model class.
virtual ~PatchFlowRateInjection()
Destructor.
const Type & value() const
Return const reference to value.
virtual void updateMesh()
Set injector locations when mesh is updated.
Mesh consisting of general polyhedral cells.
Random & rndGen()
Return refernce to the random object.
scalar timeEnd() const
Return the end-of-injection time.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
const dimensionSet dimArea(sqr(dimLength))
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
const word rhoName_
Name of carrier density field.
const fvMesh & mesh() const
Return refernce to the mesh.
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
const word phiName_
Name of carrier (mass or volume) flux field.
Templated base class for dsmc cloud.
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
A list of keyword definitions, which are a keyword followed by any number of values (e....
const scalar parcelConcentration_
Parcels to introduce per unit volume flow rate m3 [n/m3].
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
PatchFlowRateInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
scalar duration_
Injection duration [s].
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const dimensionedScalar c
Speed of light in a vacuum.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
const TimeDataEntry< scalar > concentration_
Concentration profile of particle volume to carrier volume [-].
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
Generic GeometricField class.
cachedRandom rndGen(label(0), -1)
stressControl lookup("compactNormalStress") >> compactNormalStress