32 template<
class ParcelType>
37 parcelTypeId_(dict_, -1),
40 minParcelMass_(dict_, 0.0)
44 template<
class ParcelType>
51 parcelTypeId_(cp.parcelTypeId_),
54 minParcelMass_(cp.minParcelMass_)
58 template<
class ParcelType>
65 parcelTypeId_(
dict_,
"parcelTypeId", -1),
66 rhoMin_(
dict_,
"rhoMin", 1
e-15),
68 minParcelMass_(
dict_,
"minParcelMass", 1
e-15)
72 template<
class ParcelType>
82 ParcelType(owner, coordinates, celli, tetFacei, tetPti),
99 template<
class ParcelType>
107 ParcelType(owner, position, celli),
124 template<
class ParcelType>
130 const label tetFacei,
133 const scalar nParticle0,
135 const scalar dTarget0,
140 ParcelType(owner, coordinates, celli, tetFacei, tetPti),
159 template<
class ParcelType>
167 template<
class ParcelType>
171 return parcelTypeId_.
value();
175 template<
class ParcelType>
179 return rhoMin_.
value();
183 template<
class ParcelType>
187 return rho0_.
value();
191 template<
class ParcelType>
195 return minParcelMass_.
value();
201 template<
class ParcelType>
208 template<
class ParcelType>
215 template<
class ParcelType>
222 template<
class ParcelType>
229 template<
class ParcelType>
236 template<
class ParcelType>
243 template<
class ParcelType>
250 template<
class ParcelType>
257 template<
class ParcelType>
264 template<
class ParcelType>
271 template<
class ParcelType>
278 template<
class ParcelType>
285 template<
class ParcelType>
292 template<
class ParcelType>
299 template<
class ParcelType>
306 template<
class ParcelType>
313 template<
class ParcelType>
320 template<
class ParcelType>
327 template<
class ParcelType>
334 template<
class ParcelType>
341 template<
class ParcelType>
348 template<
class ParcelType>
355 template<
class ParcelType>
362 template<
class ParcelType>
368 return rhoc_*this->
mesh().cellVolumes()[celli];
372 template<
class ParcelType>
379 template<
class ParcelType>
386 template<
class ParcelType>
393 template<
class ParcelType>
400 template<
class ParcelType>
407 template<
class ParcelType>
410 return 0.25*
areaS(d);
414 template<
class ParcelType>
421 template<
class ParcelType>
428 template<
class ParcelType>
437 return rhoc*
mag(U -
Uc_)*d/(muc + ROOTVSMALL);
441 template<
class ParcelType>
450 return rhoc*
magSqr(U -
Uc_)*d/(sigma + ROOTVSMALL);
454 template<
class ParcelType>
const vector & U() const
Return const access to velocity.
scalar tTurb() const
Return const access to time spent in turbulent eddy.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
scalar rhoMin() const
Return const access to the minimum density.
const double e
Elementary charge.
const vector & UTurb() const
Return const access to turbulent velocity fluctuation.
A list of keyword definitions, which are a keyword followed by any number of values (e...
scalar muc() const
Return const access to carrier viscosity [Pa.s].
scalar nParticle() const
Return const access to number of particles.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Class to hold kinematic particle constant properties.
scalar mass() const
Particle mass.
scalar d() const
Return const access to diameter.
scalar tTurb_
Time spent in turbulent eddy [s].
scalar rho() const
Return const access to density.
scalar rho0() const
Return const access to the particle density.
constantProperties()
Null constructor.
scalar areaP() const
Particle projected area.
KinematicParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
scalar dTarget() const
Return const access to target diameter.
vector Uc_
Velocity [m/s].
const Type & value() const
Return the value.
label parcelTypeId() const
Return const access to the parcel type id.
scalar dTarget_
Target diameter [m].
scalar age() const
Return const access to the age.
scalar Eo(const vector &a, const scalar d, const scalar sigma) const
Eotvos number.
scalar rho_
Density [kg/m3].
scalar momentOfInertia() const
Particle moment of inertia around diameter axis.
scalar muc_
Viscosity [Pa.s].
bool active_
Active flag - tracking inactive when active = false.
const dictionary dict_
Constant properties dictionary.
scalar massCell(const label celli) const
Cell owner mass.
label typeId_
Parcel type id.
vector UTurb_
Turbulent velocity fluctuation [m/s].
bool active() const
Return const access to active flag.
const vector & Uc() const
Return const access to carrier velocity [m/s].
dimensioned< scalar > magSqr(const dimensioned< Type > &)
scalar nParticle_
Number of particles in Parcel.
scalar minParcelMass() const
Return const access to the minimum parcel mass.
scalar rhoc() const
Return const access to carrier density [kg/m3].
dimensionedScalar pow3(const dimensionedScalar &ds)
scalar rhoc_
Density [kg/m3].
vector U_
Velocity of Parcel [m/s].
scalar volume() const
Particle volume.
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar We(const vector &U, const scalar d, const scalar rhoc, const scalar sigma) const
Weber number.
Mesh consisting of general polyhedral cells.
label typeId() const
Return const access to type id.
scalar Re(const vector &U, const scalar d, const scalar rhoc, const scalar muc) const
Reynolds number.
scalar areaS() const
Particle surface area.
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
const dictionary & dict() const
Return const access to the constant properties dictionary.