38template<
class CloudType>
45 { injectionMethod::imPoint,
"point" },
46 { injectionMethod::imDisc,
"disc" },
47 { injectionMethod::imDiscSegments,
"discSegments" },
51template<
class CloudType>
58 { flowType::ftConstantVelocity,
"constantVelocity" },
59 { flowType::ftPressureDrivenVelocity,
"pressureDrivenVelocity" },
60 { flowType::ftFlowRateAndDischarge,
"flowRateAndDischarge" },
66template<
class CloudType>
67void Foam::ConeNozzleInjection<CloudType>::setInjectionGeometry()
69 const auto&
mesh = this->owner().mesh();
74 Function1<vector>::New(
"position", this->coeffDict(), &
mesh)
77 positionVsTime_->userTimeToTime(this->owner().time());
79 if (positionVsTime_->constant())
81 position_ = positionVsTime_->value(0);
85 directionVsTime_.reset
87 Function1<vector>::New(
"direction", this->coeffDict(), &
mesh)
90 directionVsTime_->userTimeToTime(this->owner().time());
92 if (directionVsTime_->constant())
94 direction_ = directionVsTime_->value(0);
95 direction_.normalise();
101 scalar magTangent = 0.0;
103 while(magTangent < SMALL)
107 tangent = v - (v & direction_)*direction_;
108 magTangent =
mag(tangent);
111 tanVec1_ = tangent/magTangent;
112 tanVec2_ = direction_^tanVec1_;
117template<
class CloudType>
118void Foam::ConeNozzleInjection<CloudType>::setFlowType()
122 case flowType::ftConstantVelocity:
124 this->coeffDict().readCompat(
"Umag", {{
"UMag", -2506}}, UMag_);
127 case flowType::ftPressureDrivenVelocity:
131 Function1<scalar>::New
135 &this->owner().
mesh()
138 Pinj_->userTimeToTime(this->owner().time());
141 case flowType::ftFlowRateAndDischarge:
145 Function1<scalar>::New
149 &this->owner().
mesh()
152 Cd_->userTimeToTime(this->owner().time());
158 <<
"Unhandled flow type "
159 << flowTypeNames[flowType_]
168template<
class CloudType>
173 const word& modelName
179 injectionMethodNames.get(
"injectionMethod", this->coeffDict())
181 flowType_(flowTypeNames.get(
"flowType", this->coeffDict())),
182 outerDiameter_(this->coeffDict().getScalar(
"outerDiameter")),
183 innerDiameter_(this->coeffDict().getScalar(
"innerDiameter")),
184 duration_(this->coeffDict().getScalar(
"duration")),
185 positionVsTime_(nullptr),
190 directionVsTime_(nullptr),
201 parcelsPerSecond_(this->coeffDict().getScalar(
"parcelsPerSecond")),
233 this->coeffDict().subDict(
"sizeDistribution"),
237 t0_(this->template getModelProperty<scalar>(
"t0")),
240 this->coeffDict().template getOrDefault<scalar>(
"nInjectors", 1)
245 this->coeffDict().template getOrDefault<
vector>
247 "initialInjectorDir",
259 if (innerDiameter_ >= outerDiameter_)
262 <<
"Inner diameter must be less than the outer diameter:" <<
nl
263 <<
" innerDiameter: " << innerDiameter_ <<
nl
264 <<
" outerDiameter: " << outerDiameter_
268 if (nInjectors_ < SMALL)
271 <<
"Number of injectors in angular-segmented disc "
272 <<
"must be positive" <<
nl
273 <<
" nInjectors: " << nInjectors_ <<
nl
278 const Time& time =
owner.db().time();
280 flowRateProfile_->userTimeToTime(time);
281 thetaInner_->userTimeToTime(time);
282 thetaOuter_->userTimeToTime(time);
286 omegaPtr_->userTimeToTime(time);
289 setInjectionGeometry();
301template<
class CloudType>
308 injectionMethod_(im.injectionMethod_),
309 flowType_(im.flowType_),
310 outerDiameter_(im.outerDiameter_),
311 innerDiameter_(im.innerDiameter_),
312 duration_(im.duration_),
313 positionVsTime_(im.positionVsTime_.clone()),
314 position_(im.position_),
315 injectorCell_(im.injectorCell_),
316 tetFacei_(im.tetFacei_),
318 directionVsTime_(im.directionVsTime_.clone()),
319 direction_(im.direction_),
320 omegaPtr_(im.omegaPtr_.clone()),
321 parcelsPerSecond_(im.parcelsPerSecond_),
322 flowRateProfile_(im.flowRateProfile_.clone()),
323 thetaInner_(im.thetaInner_.clone()),
324 thetaOuter_(im.thetaOuter_.clone()),
325 sizeDistribution_(im.sizeDistribution_.clone()),
327 nInjectors_(im.nInjectors_),
328 Uinjector_(im.Uinjector_),
329 initialInjectorDir_(im.initialInjectorDir_),
330 tanVec1_(im.tanVec1_),
331 tanVec2_(im.tanVec2_),
341template<
class CloudType>
345 if (positionVsTime_->constant())
347 position_ = positionVsTime_->value(0);
349 this->findCellAtPosition
360template<
class CloudType>
363 return this->
SOI_ + duration_;
367template<
class CloudType>
374 if ((time0 >= 0.0) && (time0 < duration_))
376 return floor((time1 - time0)*parcelsPerSecond_);
383template<
class CloudType>
390 if ((time0 >= 0.0) && (time0 < duration_))
392 return flowRateProfile_->integrate(time0, time1);
399template<
class CloudType>
412 const scalar t = time - this->SOI_;
414 if (!directionVsTime_->constant())
416 direction_ = directionVsTime_->value(t);
417 direction_.normalise();
421 scalar magTangent = 0.0;
423 while(magTangent < SMALL)
427 tangent = v - (v & direction_)*direction_;
428 magTangent =
mag(tangent);
431 tanVec1_ = tangent/magTangent;
432 tanVec2_ = direction_^tanVec1_;
438 switch (injectionMethod_)
440 case injectionMethod::imPoint:
442 if (positionVsTime_->constant())
444 position = position_;
445 cellOwner = injectorCell_;
446 tetFacei = tetFacei_;
451 position = positionVsTime_->value(t);
454 const vector position0(positionVsTime_->value(t0_));
455 const scalar dt = t - t0_;
458 Uinjector_ = (position - position0)/dt;
461 this->findCellAtPosition
471 case injectionMethod::imDisc:
474 scalar dr = outerDiameter_ - innerDiameter_;
475 scalar r = 0.5*(innerDiameter_ + frac*dr);
477 position = positionVsTime_->value(t) + r*normal_;
480 const vector position0(positionVsTime_->value(t0_) + r*normal_);
481 const scalar dt = t - t0_;
484 Uinjector_ = (position - position0)/dt;
487 this->findCellAtPosition
496 case injectionMethod::imDiscSegments:
502 const label injectorIndex =
506 const scalar angle = injectorIndex*angleIncrement;
521 const scalar dr = outerDiameter_ - innerDiameter_;
522 const scalar radius = 0.5*(innerDiameter_ + fraction*dr);
526 normal_ =
R & initialInjectorDir_;
529 const vector radialOffset(radius*normal_);
530 position = positionVsTime_->value(t) + radialOffset;
533 const scalar dt = t - t0_;
538 positionVsTime_->value(t0_) + radialOffset
540 Uinjector_ = (position - position0)/dt;
543 this->findCellAtPosition
555 <<
"Unhandled injection method "
556 << injectionMethodNames[injectionMethod_]
566template<
class CloudType>
578 scalar t = time - this->SOI_;
579 scalar ti = thetaInner_->value(t);
580 scalar to = thetaOuter_->value(t);
581 scalar coneAngle =
degToRad(
rndGen.sample01<scalar>()*(to - ti) + ti);
584 scalar dcorr =
cos(coneAngle);
587 vector dirVec = dcorr*direction_;
593 case flowType::ftConstantVelocity:
595 parcel.
U() = UMag_*dirVec + Uinjector_;
598 case flowType::ftPressureDrivenVelocity:
600 scalar pAmbient = this->owner().pAmbient();
601 scalar
rho = parcel.rho();
602 scalar UMag =
::sqrt(2.0*(Pinj_->value(t) - pAmbient)/
rho);
603 parcel.
U() = UMag*dirVec + Uinjector_;
606 case flowType::ftFlowRateAndDischarge:
610 scalar massFlowRate =
612 *flowRateProfile_->value(t)
613 /this->volumeTotal();
615 scalar Umag = massFlowRate/(parcel.rho()*Cd_->value(t)*(Ao - Ai));
616 parcel.
U() = Umag*dirVec + Uinjector_;
622 <<
"Unhandled injection method "
623 << flowTypeNames[flowType_]
631 const scalar omega = omegaPtr_->value(t);
634 const vector r(
p0 - (
p0 & direction_)*direction_);
635 const scalar rMag =
mag(r);
639 parcel.
U() += omega*rMag*d;
643 parcel.d() = sizeDistribution_->sample();
647template<
class CloudType>
654template<
class CloudType>
661template<
class CloudType>
666 if (this->writeTime())
668 this->setModelProperty(
"t0", t0_);
#define R(A, B, C, D, E, F, K, M)
const CloudType & owner() const
Return const access to the owner cloud.
virtual bool writeTime() const
Flag to indicate when to write a property.
ConeNozzleInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
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.
flowType
Flow type enumeration.
@ ftPressureDrivenVelocity
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.
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
static const Enum< injectionMethod > injectionMethodNames
injectionMethod
Injection method enumeration.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
static const Enum< flowType > flowTypeNames
virtual void info()
Write injection info.
virtual void updateMesh()
Set injector locations when mesh is updated.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
scalar timeEnd() const
Return the end-of-injection time.
Random & rndGen()
Return reference to the random object.
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.
scalar massTotal() const
Return mass of particles to introduce.
static autoPtr< InjectionModel< CloudType > > New(const dictionary &dict, CloudType &owner)
Selector with lookup from dictionary.
virtual void info()
Write injection info.
scalar volumeTotal() const
Return the total volume to be injected across the event.
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.
virtual bool findCellAtPosition(label &celli, label &tetFacei, label &tetPti, vector &position, bool errorOnNotFound=true)
Find the cell that contains the supplied position.
scalar SOI_
Start of injection [s].
Type globalSample01()
Return a sample whose components lie in the range [0,1].
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.
Vector< Cmpt > & normalise(const scalar tol=ROOTVSMALL)
Inplace normalise the vector by its magnitude.
static tensor rotation(const vector &axis, const scalar angle, bool degrees=false)
The rotation tensor for given axis/angle.
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...
vector position() const
Return current particle position.
bool getModelProperty(const word &entryName, Type &value) const
Retrieve generic property from the sub-model.
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.
void setModelProperty(const word &entryName, const Type &value)
Add generic property to 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.
const volScalarField & p0
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
constexpr scalar pi(M_PI)
constexpr scalar twoPi(2 *M_PI)
Different types of constants.
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
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
dimensionedScalar sin(const dimensionedScalar &ds)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
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...
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Unit conversion functions.