47Foam::resolutionIndexModels::CelikEtaIndex::eta()
const
60Foam::resolutionIndexModels::CelikEtaIndex::epsilon()
const
62 const auto& kSgs = getOrReadField<volScalarField>(kName_);
63 const auto& Delta = getOrReadField<volScalarField>(deltaName_);
64 tmp<volScalarField> tnuEff = nuEff();
67 return tnuEff*kSgs/(Ck_*
sqr(Delta));
71Foam::tmp<Foam::volScalarField>
72Foam::resolutionIndexModels::CelikEtaIndex::nuEff()
const
74 const auto&
nu = getOrReadField<volScalarField>(nuName_);
75 const auto& nuSgs = getOrReadField<volScalarField>(nutName_);
76 tmp<volScalarField> tnuNum = nuNum();
79 return tnuNum + nuSgs +
nu;
83Foam::tmp<Foam::volScalarField>
84Foam::resolutionIndexModels::CelikEtaIndex::nuNum()
const
86 const auto& Delta = getOrReadField<volScalarField>(deltaName_);
88 tmp<volScalarField> tkNum = kNum();
91 return sign(tkNum.cref())*Cnu_*Delta*
sqrt(
mag(tkNum.cref()));
95Foam::tmp<Foam::volScalarField>
96Foam::resolutionIndexModels::CelikEtaIndex::kNum()
const
98 const auto& kSgs = getOrReadField<volScalarField>(kName_);
99 const auto& Delta = getOrReadField<volScalarField>(deltaName_);
101 tmp<volScalarField> th =
cbrt(
V());
104 return Cn_*
sqr(th/Delta)*kSgs;
142 alphaEta_ =
dict.getOrDefault<scalar>(
"alphaEta", 0.05);
143 m_ =
dict.getOrDefault<scalar>(
"m", 0.5);
144 Cnu_ =
dict.getOrDefault<scalar>(
"Cnu", 0.1);
145 Cn_ =
dict.getOrDefault<scalar>(
"Cn", 1.0);
146 Ck_ =
dict.getOrDefault<scalar>(
"Ck", 0.0376);
147 kName_ =
dict.getOrDefault<
word>(
"k",
"k");
148 deltaName_ =
dict.getOrDefault<
word>(
"delta",
"delta");
150 nutName_ =
dict.getOrDefault<
word>(
"nut",
"nut");
163 auto& index = getOrReadField<volScalarField>(resultName());
166 index = 1.0/(1.0 + alphaEta_*
pow(th/teta, m_));
167 index.correctBoundaryConditions();
175 const auto& index = getOrReadField<volScalarField>(resultName());
177 Info<<
tab <<
"writing field:" << index.name() <<
endl;
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Mesh data needed to do the Finite Volume discretisation.
A base class for resolutionIndex models.
const fvMesh & mesh() const noexcept
Return const reference to the mesh.
tmp< volScalarField > V() const
Return cell volume field.
resolutionIndexModel(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from components.
GeoFieldType & getOrReadField(const word &fieldName) const
Return requested field from the object registry or read+register the field to the object registry.
const word & resultName() const noexcept
Return const reference to the result name.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
Computes a single-mesh resolution index according to Celik et al.'s index using Kolmogorov length sca...
CelikEtaIndex(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from components.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
virtual bool execute()
Execute the function-object operations.
virtual bool write()
Write the function-object results.
A class for managing temporary objects.
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
const expr V(m.psi().mesh().V())
A namespace for various resolutionIndex model implementations.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
bool read(const char *buf, int32_t &val)
Same as readInt32.
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pow025(const dimensionedScalar &ds)
constexpr char tab
The tab '\t' character(0x09).