Go to the documentation of this file.
34 template<
class CloudType>
42 turbulenceModel::propertiesName,
54 <<
"Turbulence model not found in mesh database" <<
nl
55 <<
"Database objects include: " << obr.
sortedToc()
62 template<
class CloudType>
70 turbulenceModel::propertiesName,
78 return turb->epsilon();
82 <<
"Turbulence model not found in mesh database" <<
nl
83 <<
"Database objects include: " << obr.
sortedToc()
92 template<
class CloudType>
102 epsilonPtr_(
nullptr),
107 template<
class CloudType>
126 template<
class CloudType>
135 template<
class CloudType>
153 if (tepsilon.
isTmp())
155 epsilonPtr_ = tepsilon.
ptr();
160 epsilonPtr_ = &tepsilon();
171 if (ownEpsilon_ && epsilonPtr_)
180 template<
class CloudType>
185 os.writeEntry(
"ownK", ownK_);
186 os.writeEntry(
"ownEpsilon", ownEpsilon_);
A class for handling words, derived from Foam::string.
constexpr const char *const group
Group name for atomic constants.
A class for managing temporary objects.
Base class for dispersion modelling.
Template functions to aid in the implementation of demand driven data.
const volScalarField * epsilonPtr_
Turbulence epsilon.
virtual void write(Ostream &os) const
Write.
DispersionRASModel(const dictionary &dict, CloudType &owner)
Construct from components.
bool isTmp() const noexcept
Identical to is_pointer()
const volScalarField * kPtr_
Turbulence k.
Registry of regIOobjects.
void deleteDemandDrivenData(DataPtr &dataPtr)
virtual void cacheFields(const bool store)
Cache carrier fields.
const fvMesh & mesh() const
Return reference to the mesh.
Templated base class for dsmc cloud.
Abstract base class for turbulence models (RAS, LES and laminar).
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
OBJstream os(runTime.globalPath()/outputName)
tmp< volScalarField > kModel() const
Return the k field from the turbulence model.
T * ptr() const
Return managed pointer for reuse, or clone() the object reference.
errorManip< error > abort(error &err)
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
bool ownEpsilon_
Take ownership of the epsilon field.
virtual ~DispersionRASModel()
Destructor.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
bool ownK_
Take ownership of the k field.
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
compressible::turbulenceModel & turb
Base class for particle dispersion models based on RAS turbulence.
tmp< volScalarField > epsilonModel() const
Return the epsilon field from the turbulence model.