36 template<
class CloudType>
39 const word injectionMethod =
40 this->coeffDict().template lookupOrDefault<word>
46 if (injectionMethod ==
"point" || injectionMethod ==
word::null)
48 injectionMethod_ = imPoint;
52 else if (injectionMethod ==
"disc")
54 injectionMethod_ = imDisc;
56 this->coeffDict().lookup(
"dInner") >> dInner_;
57 this->coeffDict().lookup(
"dOuter") >> dOuter_;
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;
82 Umag_.reset(this->coeffDict());
84 else if (flowType ==
"pressureDrivenVelocity")
86 flowType_ = ftPressureDrivenVelocity;
88 Pinj_.reset(this->coeffDict());
90 else if (flowType ==
"flowRateAndDischarge")
92 flowType_ = ftFlowRateAndDischarge;
94 this->coeffDict().lookup(
"dInner") >> dInner_;
95 this->coeffDict().lookup(
"dOuter") >> dOuter_;
97 Cd_.reset(this->coeffDict());
102 <<
"flowType must be either 'constantVelocity', " 103 <<
"'pressureDrivenVelocity' or 'flowRateAndDischarge'" 111 template<
class CloudType>
116 const word& modelName
120 injectionMethod_(imPoint),
121 flowType_(ftConstantVelocity),
142 injectorTetFace_(-1),
144 duration_(this->coeffDict().
template lookup<scalar>(
"duration")),
147 this->coeffDict().
template lookup<scalar>(
"parcelsPerSecond")
180 this->coeffDict().subDict(
"sizeDistribution"), owner.rndGen()
185 Umag_(owner.db().time(),
"Umag"),
186 Cd_(owner.db().time(),
"Cd"),
187 Pinj_(owner.db().time(),
"Pinj")
189 duration_ = owner.db().time().userTimeToTime(duration_);
191 setInjectionMethod();
196 this->volumeTotal_ = flowRateProfile_.integral(0, duration_);
202 template<
class CloudType>
209 injectionMethod_(im.injectionMethod_),
210 flowType_(im.flowType_),
211 position_(im.position_),
212 positionIsConstant_(im.positionIsConstant_),
213 direction_(im.direction_),
214 injectorCell_(im.injectorCell_),
215 injectorTetFace_(im.injectorTetFace_),
216 injectorTetPt_(im.injectorTetPt_),
217 duration_(im.duration_),
218 parcelsPerSecond_(im.parcelsPerSecond_),
219 flowRateProfile_(im.flowRateProfile_),
220 thetaInner_(im.thetaInner_),
221 thetaOuter_(im.thetaOuter_),
222 sizeDistribution_(im.sizeDistribution_().clone().ptr()),
233 template<
class CloudType>
240 template<
class CloudType>
243 if (injectionMethod_ == imPoint && positionIsConstant_)
245 vector position = position_.value(0);
246 this->findCellAtPosition
257 template<
class CloudType>
260 return this->SOI_ + duration_;
264 template<
class CloudType>
271 if (time0 >= 0 && time0 < duration_)
277 return floor(parcelsPerSecond_*time1 - this->parcelsAddedTotal());
286 template<
class CloudType>
293 if (time0 >= 0 && time0 < duration_)
295 return flowRateProfile_.integral(time0, time1);
304 template<
class CloudType>
318 const scalar t = time - this->SOI_;
320 switch (injectionMethod_)
324 position = position_.value(t);
325 if (positionIsConstant_)
327 cellOwner = injectorCell_;
328 tetFacei = injectorTetFace_;
329 tetPti = injectorTetPt_;
333 this->findCellAtPosition
352 const scalar d =
sqrt((1 - frac)*
sqr(dInner_) + frac*
sqr(dOuter_));
353 position = position_.value(t) + d/2*tanVec;
354 this->findCellAtPosition
372 template<
class CloudType>
383 const scalar t = time - this->SOI_;
390 scalar theta = vGreat;
392 switch (injectionMethod_)
397 const scalar frac = rndGen.
scalar01();
401 tanVec = t1*
cos(beta) + t2*
sin(beta);
407 (1 - frac)*
sqr(thetaInner_.value(t))
408 + frac*
sqr(thetaOuter_.value(t))
415 const scalar r =
mag(parcel.position() - position_.value(t));
416 const scalar frac = (2*r - dInner_)/(dOuter_ - dInner_);
417 tanVec =
normalised(parcel.position() - position_.value(t));
421 (1 - frac)*thetaInner_.value(t)
422 + frac*thetaOuter_.value(t)
443 case ftConstantVelocity:
445 parcel.U() = Umag_.value(t)*dirVec;
448 case ftPressureDrivenVelocity:
450 const scalar pAmbient = this->owner().pAmbient();
451 const scalar
rho = parcel.rho();
452 const scalar Umag =
::sqrt(2*(Pinj_.value(t) - pAmbient)/rho);
453 parcel.U() = Umag*dirVec;
456 case ftFlowRateAndDischarge:
458 const scalar A = 0.25*
pi*(
sqr(dOuter_) -
sqr(dInner_));
459 const scalar massFlowRate =
460 this->massTotal()*flowRateProfile_.value(t)/this->volumeTotal();
462 massFlowRate/(parcel.rho()*Cd_.value(t)*A);
463 parcel.U() = Umag*dirVec;
473 parcel.d() = sizeDistribution_->sample();
477 template<
class CloudType>
484 template<
class CloudType>
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
ConeInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
static const Vector< scalar > max
errorManipArg< error, int > exit(error &err, const int errNo=1)
A list of keyword definitions, which are a keyword followed by any number of values (e...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Templated injection model class.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
dimensionedScalar sqrt(const dimensionedScalar &ds)
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.
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
Form normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
dimensionedScalar cos(const dimensionedScalar &ds)
A class for handling words, derived from string.
Light wrapper around Function1 to provide a mechanism to update time-based entries.
scalar timeEnd() const
Return the end-of-injection time.
static const word null
An empty word.
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
scalar globalScalar01()
Advance the state and return a scalar sample from a uniform.
const scalar twoPi(2 *pi)
ParcelType parcelType
Type of parcel the cloud was instantiated for.
dimensionedScalar sin(const dimensionedScalar &ds)
virtual void topoChange()
Set injector locations when mesh is updated.
virtual ~ConeInjection()
Destructor.
Vector< Cmpt > perpendicular(const Vector< Cmpt > &v)
static autoPtr< distributionModel > New(const dictionary &dict, Random &rndGen)
Selector.
Templated function that returns a constant value.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
dimensioned< scalar > mag(const dimensioned< Type > &)
Templated base class for dsmc cloud.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
scalar scalar01()
Advance the state and return a scalar sample from a uniform.
This injector injects particles in a number of cones. The user specifies a position and a direction t...