Go to the documentation of this file.
44 #ifndef KinematicParcel_H
45 #define KinematicParcel_H
60 template<
class ParcelType>
65 template<
class ParcelType>
76 template<
class ParcelType>
144 inline scalar
rhoMin()
const;
147 inline scalar
rho0()
const;
154 template<
class CloudType>
220 inline const vector&
g()
const;
282 template<
class TrackData>
317 +
" (UTurbx UTurby UTurbz)"
330 const label tetFaceI,
340 const label tetFaceI,
343 const scalar nParticle0,
345 const scalar dTarget0,
404 inline bool active()
const;
413 inline scalar
d()
const;
419 inline const vector&
U()
const;
422 inline scalar
rho()
const;
425 inline scalar
age()
const;
428 inline scalar
tTurb()
const;
434 inline scalar
rhoc()
const;
440 inline scalar
muc()
const;
464 inline scalar&
rho();
467 inline scalar&
age();
470 inline scalar&
tTurb();
485 inline scalar
mass()
const;
491 inline scalar
volume()
const;
494 inline static scalar
volume(
const scalar
d);
497 inline scalar
areaP()
const;
500 inline static scalar
areaP(
const scalar
d);
503 inline scalar
areaS()
const;
506 inline static scalar
areaS(
const scalar
d);
538 template<
class TrackData>
547 template<
class TrackData>
556 template<
class TrackData>
568 template<
class TrackData>
569 bool move(TrackData& td,
const scalar trackTime);
579 template<
class TrackData>
584 template<
class TrackData>
590 const scalar trackFraction,
596 template<
class TrackData>
604 template<
class TrackData>
613 template<
class TrackData>
636 template<
class CloudType>
640 template<
class CloudType>
646 friend Ostream& operator<< <ParcelType>
static const std::size_t sizeofFields_
Size in bytes of the fields.
Factory class to read-construct particles used for.
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
label parcelTypeId() const
Return const access to the parcel type id.
label typeId() const
Return const access to type id.
scalar tTurb_
Time spent in turbulent eddy [s].
const vector & UTurb() const
Return const access to turbulent velocity fluctuation.
AddToPropertyList(ParcelType, " active"+" typeId"+" nParticle"+" d"+" dTarget "+" (Ux Uy Uz)"+" rho"+" age"+" tTurb"+" (UTurbx UTurby UTurbz)")
String representation of properties.
void cellValueSourceCorrection(TrackData &td, const scalar dt, const label cellI)
Correct cell values using latest transfer information.
const dimensionedScalar mu
Atomic mass unit.
scalar nParticle() const
Return const access to number of particles.
scalar muc_
Viscosity [Pa.s].
scalar dTarget_
Target diameter [m].
demandDrivenEntry< scalar > minParcelMass_
Minimum parcel mass [kg].
const vector & U() const
Return const access to velocity.
scalar tTurb() const
Return const access to time spent in turbulent eddy.
virtual scalar wallImpactDistance(const vector &n) const
The nearest distance to a wall that the particle can be.
scalar dTarget() const
Return const access to target diameter.
TrackingData(CloudType &cloud, trackPart part=tpLinearTrack)
Construct from components.
scalar rho0() const
Return const access to the particle density.
const interpolation< scalar > & muInterp() const
Return conat access to the interpolator for continuous.
demandDrivenEntry< label > parcelTypeId_
Parcel type id - used for post-processing to flag the type.
bool cp(const fileName &src, const fileName &dst)
Copy, recursively if necessary, the source to the destination.
const vector calcVelocity(TrackData &td, const scalar dt, const label cellI, const scalar Re, const scalar mu, const scalar mass, const vector &Su, vector &dUTrans, scalar &Spu) const
Calculate new particle velocity.
This function object reads fields from the time directories and adds them to the mesh database for fu...
virtual autoPtr< particle > clone(const polyMesh &mesh) const
Construct and return a (basic particle) clone.
bool move(TrackData &td, const scalar trackTime)
Move the parcel.
Mesh consisting of general polyhedral cells.
scalar age() const
Return const access to the age.
autoPtr< KinematicParcel< ParcelType > > operator()(Istream &is) const
scalar muc() const
Return const access to carrier viscosity [Pa.s].
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
scalar rhoc() const
Return const access to carrier density [kg/m3].
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
A patch is a list of labels that address the faces in the global face list.
const interpolation< vector > & UInterp() const
Return conat access to the interpolator for continuous.
vector UTurb_
Turbulent velocity fluctuation [m/s].
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
trackPart part() const
Return the part of the tracking operation taking place.
vector Uc_
Velocity [m/s].
scalar momentOfInertia() const
Particle moment of inertia around diameter axis.
void hitProcessorPatch(const processorPolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a.
Neighbour processor patch.
scalar volume() const
Particle volume.
const vector & g_
Local gravitational or other body-force acceleration.
scalar minParcelMass() const
Return const access to the minimum parcel mass.
scalar areaP() const
Particle projected area.
Templated base class for dsmc cloud.
label faceInterpolation() const
Return the index of the face used in the interpolation routine.
vector U_
Velocity of Parcel [m/s].
scalar Re(const vector &U, const scalar d, const scalar rhoc, const scalar muc) const
Reynolds number.
A list of keyword definitions, which are a keyword followed by any number of values (e....
scalar massCell(const label cellI) const
Cell owner mass.
scalar nParticle_
Number of particles in Parcel.
scalar rho_
Density [kg/m3].
bool active_
Active flag - tracking inactive when active = false.
void hitFace(int &td)
Overridable function to handle the particle hitting a face.
bool active() const
Return const access to active flag.
autoPtr< interpolation< scalar > > rhoInterp_
Density interpolator.
label typeId_
Parcel type id.
scalar rhoc_
Density [kg/m3].
scalar areaS() const
Particle surface area.
autoPtr< interpolation< vector > > UInterp_
Velocity interpolator.
bool hitPatch(const polyPatch &p, TrackData &td, const label patchI, const scalar trackFraction, const tetIndices &tetIs)
Overridable function to handle the particle hitting a patch.
scalar rhoMin() const
Return const access to the minimum density.
scalar We(const vector &U, const scalar d, const scalar rhoc, const scalar sigma) const
Weber number.
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
A cloud is a collection of lagrangian particles.
demandDrivenEntry< scalar > rho0_
Particle density [kg/m3] (constant)
scalar Eo(const vector &a, const scalar d, const scalar sigma) const
Eotvos number.
void setCellValues(TrackData &td, const scalar dt, const label cellI)
Set cell values.
static label maxTrackAttempts
Number of particle tracking attempts before we assume that it stalls.
Kinematic parcel class with rotational motion (as spherical particles only) and one/two-way coupling ...
TypeName("KinematicParcel")
Runtime type information.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
const dictionary & dict() const
Return const access to the constant properties dictionary.
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
const dimensionedScalar c
Speed of light in a vacuum.
const interpolation< scalar > & rhoInterp() const
Return conat access to the interpolator for continuous.
const vector & Uc() const
Return const access to carrier velocity [m/s].
scalar d() const
Return const access to diameter.
scalar rho() const
Return const access to density.
KinematicParcel(const polyMesh &mesh, const vector &position, const label cellI, const label tetFaceI, const label tetPtI)
Construct from owner, position, and cloud owner.
void calc(TrackData &td, const scalar dt, const label cellI)
Update parcel properties over the time interval.
static void writeFields(const CloudType &c)
Write.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
autoPtr< interpolation< scalar > > muInterp_
Dynamic viscosity interpolator.
const dictionary dict_
Constant properties dictionary.
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
demandDrivenEntry< scalar > rhoMin_
Minimum density [kg/m3].
constantProperties()
Null constructor.
scalar mass() const
Particle mass.
static void readFields(CloudType &c)
Read.
Class to hold kinematic particle constant properties.
void hitWallPatch(const wallPolyPatch &, TrackData &td, const tetIndices &)
Overridable function to handle the particle hitting a wallPatch.
iNew(const polyMesh &mesh)