34 template<
class CloudType>
37 const word injectionMethod =
38 this->coeffDict().template lookupOrDefault<word>
44 if (injectionMethod ==
"point" || injectionMethod ==
word::null)
46 injectionMethod_ = imPoint;
50 else if (injectionMethod ==
"disc")
52 injectionMethod_ = imDisc;
55 this->coeffDict().template lookup<scalar>(
"dInner",
dimLength);
57 this->coeffDict().template lookup<scalar>(
"dOuter",
dimLength);
62 <<
"injectionMethod must be either 'point' or 'disc'"
68 template<
class CloudType>
72 this->coeffDict().template lookupOrDefault<word>
78 if (flowType ==
"constantVelocity" || flowType ==
word::null)
80 flowType_ = ftConstantVelocity;
87 this->owner().db().time().userUnits(),
93 else if (flowType ==
"pressureDrivenVelocity")
95 flowType_ = ftPressureDrivenVelocity;
102 this->owner().db().time().userUnits(),
108 else if (flowType ==
"flowRateAndDischarge")
110 flowType_ = ftFlowRateAndDischarge;
113 this->coeffDict().template lookup<scalar>(
"dInner",
dimLength);
115 this->coeffDict().template lookup<scalar>(
"dOuter",
dimLength);
122 this->owner().db().time().userUnits(),
131 <<
"flowType must be either 'constantVelocity', "
132 <<
"'pressureDrivenVelocity' or 'flowRateAndDischarge'"
140 template<
class CloudType>
145 const word& modelName
149 injectionMethod_(imPoint),
150 flowType_(ftConstantVelocity),
156 this->owner().db().time().userUnits(),
166 this->owner().db().time().userUnits(),
173 injectorTetFace_(-1),
175 duration_(this->readDuration(
dict, owner)),
176 massFlowRate_(this->readMassFlowRate(
dict, owner, duration_)),
177 parcelsPerSecond_(this->readParcelsPerSecond(
dict, owner)),
183 this->owner().db().time().userUnits(),
193 this->owner().db().time().userUnits(),
203 this->coeffDict().subDict(
"sizeDistribution"),
205 owner.
rndGen().generator()
214 setInjectionMethod();
222 template<
class CloudType>
229 injectionMethod_(im.injectionMethod_),
230 flowType_(im.flowType_),
231 position_(im.position_, false),
232 direction_(im.direction_, false),
233 injectorCoordinates_(im.injectorCoordinates_),
234 injectorCell_(im.injectorCell_),
235 injectorTetFace_(im.injectorTetFace_),
236 injectorTetPt_(im.injectorTetPt_),
237 duration_(im.duration_),
238 massFlowRate_(im.massFlowRate_, false),
239 parcelsPerSecond_(im.parcelsPerSecond_, false),
240 thetaInner_(im.thetaInner_, false),
241 thetaOuter_(im.thetaOuter_, false),
242 sizeDistribution_(im.sizeDistribution_, false),
245 Umag_(im.Umag_, false),
247 Pinj_(im.Pinj_, false)
253 template<
class CloudType>
260 template<
class CloudType>
263 if (injectionMethod_ == imPoint && position_->constant())
265 vector position = position_->value(0);
266 this->findCellAtPosition
269 injectorCoordinates_,
278 template<
class CloudType>
281 return this->SOI_ + duration_;
285 template<
class CloudType>
292 if (time0 >= 0 && time0 < duration_)
301 parcelsPerSecond_->integral(0, time1)
302 - this->parcelsAddedTotal()
312 template<
class CloudType>
319 if (time0 >= 0 && time0 < duration_)
321 return massFlowRate_->integral(time0, time1);
330 template<
class CloudType>
345 const scalar t = time - this->SOI_;
347 switch (injectionMethod_)
352 if (position_->constant())
354 coordinates = injectorCoordinates_;
355 celli = injectorCell_;
356 tetFacei = injectorTetFace_;
357 tetPti = injectorTetPt_;
361 this->findCellAtPosition
375 const scalar beta =
twoPi*this->globalScalar01(
rndGen);
376 const scalar frac = this->globalScalar01(
rndGen);
381 const scalar d =
sqrt((1 - frac)*
sqr(dInner_) + frac*
sqr(dOuter_));
383 this->findCellAtPosition
402 template<
class CloudType>
408 typename CloudType::parcelType::trackingData& td,
412 const polyMesh& mesh = this->owner().mesh();
416 const scalar t = time - this->SOI_;
423 scalar theta = vGreat;
425 switch (injectionMethod_)
434 tanVec = t1*
cos(beta) + t2*
sin(beta);
438 (1 - frac)*
sqr(thetaInner_->value(t))
439 + frac*
sqr(thetaOuter_->value(t))
445 const scalar r =
mag(parcel.position(mesh) - position_->value(t));
446 const scalar frac = (2*r - dInner_)/(dOuter_ - dInner_);
447 tanVec =
normalised(parcel.position(mesh) - position_->value(t));
449 (1 - frac)*thetaInner_->value(t)
450 + frac*thetaOuter_->value(t);
470 case ftConstantVelocity:
472 parcel.U() = Umag_->value(t)*dirVec;
475 case ftPressureDrivenVelocity:
477 const scalar pAmbient = this->owner().pAmbient();
478 const scalar
rho = parcel.rho();
479 const scalar Umag =
::sqrt(2*(Pinj_->value(t) - pAmbient)/
rho);
480 parcel.U() = Umag*dirVec;
483 case ftFlowRateAndDischarge:
485 const scalar
A = 0.25*
pi*(
sqr(dOuter_) -
sqr(dInner_));
487 massFlowRate_->value(t)/(parcel.rho()*Cd_->value(t)*
A);
488 parcel.U() = Umag*dirVec;
498 parcel.d() = sizeDistribution_->sample();
502 template<
class CloudType>
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
This injector injects particles in a number of cones. The user specifies a position and a direction t...
virtual ~ConeInjection()
Destructor.
virtual void topoChange()
Set injector locations when mesh is updated.
ConeInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType::trackingData &td, typename CloudType::parcelType &parcel)
Set the parcel properties.
virtual void setPositionAndCell(const label parcelI, const label nParcels, const scalar time, barycentric &coordinates, label &celli, label &tetFacei, label &tetPti, label &facei)
Set the injection position and owner cell, tetFace and tetPt.
virtual label nParcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
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.
Templated base class for dsmc cloud.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Run-time selectable general function of one variable.
static autoPtr< Function1< scalar > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Templated injection model class.
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.
Mesh consisting of general polyhedral cells.
scalar scalar01()
Return a scalar uniformly distributed between zero and one.
A class for handling words, derived from string.
static const word null
An empty word.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const scalar twoPi(2 *pi)
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensionedScalar pos(const dimensionedScalar &ds)
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 dimPressure
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimless
dimensionedScalar sin(const dimensionedScalar &ds)
const dimensionSet dimLength
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionSet normalised(const dimensionSet &)
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
const dimensionSet dimVelocity
dimensionSet perpendicular(const dimensionSet &)
dimensionedScalar cos(const dimensionedScalar &ds)
const unitConversion unitDegrees
randomGenerator rndGen(653213)