50 const volScalarField::Internal&
k,
51 const volScalarField::Internal&
epsilon
55 tmp<volScalarField::Internal> L_(
pow(Cmu_, 0.75)*
pow(
k, 1.5)/
epsilon);
59 return (C2_ - C1_)*
pow(L_/Lmax_, n_);
73 tmp<volScalarField::Internal> L_(
sqrt(
k)/(
pow025(Cmu_)*omega));
84 const word& sourceName,
85 const word& modelType,
92 rhoName_(coeffs_.getOrDefault<
word>(
"rho",
"rho")),
98 coeffs_.getCheckOrDefault<scalar>
102 [&](const scalar Lmax){
return Lmax > SMALL; }
111 coeffs_.getCheckOrDefault<scalar>
115 [&](
const scalar
n){
return n > SMALL; }
124 const auto* turbPtr =
133 <<
"Unable to find a turbulence model."
137 fieldNames_.resize(1);
139 tmp<volScalarField> tepsilon = turbPtr->epsilon();
140 tmp<volScalarField> tomega = turbPtr->omega();
142 if (!tepsilon.isTmp())
144 fieldNames_[0] = tepsilon().name();
146 const dictionary& turbDict = turbPtr->coeffDict();
147 Cmu_.read(
"Cmu", turbDict);
148 C1_.read(
"C1", turbDict);
149 C2_.read(
"C2", turbDict);
150 C3_.read(
"C3", turbDict);
152 else if (!tomega.isTmp())
155 fieldNames_[0] = tomega().name();
157 const dictionary& turbDict = turbPtr->coeffDict();
158 Cmu_.read(
"betaStar", turbDict);
163 <<
"Unable to find neither epsilon nor omega field." <<
nl
164 <<
"atmLengthScaleTurbSource needs either epsilon or omega field."
170 Log <<
" Applying atmLengthScaleTurbSource to: " << fieldNames_[0]
185 atmLengthScaleTurbSourceEpsilon
195 atmLengthScaleTurbSourceOmega
234 atmLengthScaleTurbSourceEpsilon(
alpha,
rho, eqn, fieldi);
238 atmLengthScaleTurbSourceOmega(
alpha,
rho, eqn, fieldi);
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
DimensionedField< scalar, volMesh > Internal
The internal field type from which this GeometricField is derived.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Mesh data needed to do the Finite Volume discretisation.
Applies sources on either epsilon or omega to correct mixing-length scale estimations for atmospheric...
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Intermediate abstract class for handling cell-set options for the derived fvOptions.
Base abstract class for handling finite volume options (i.e. fvOption).
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
A class for managing temporary objects.
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
errorManip< error > abort(error &err)
static constexpr const zero Zero
Global zero (0)
dimensionedScalar pow025(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a)
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)