Go to the documentation of this file.
50 inletOutletFvPatchScalarField(
p, iF),
55 this->refValue() = 0.0;
56 this->refGrad() = 0.0;
57 this->valueFraction() = 0.0;
70 inletOutletFvPatchScalarField(ptf,
p, iF, mapper),
72 mixingLength_(ptf.mixingLength_),
85 inletOutletFvPatchScalarField(
p, iF),
86 kName_(
dict.getOrDefault<
word>(
"k",
"k")),
93 this->phiName_ =
dict.getOrDefault<
word>(
"phi",
"phi");
97 this->refValue() = 0.0;
98 this->refGrad() = 0.0;
99 this->valueFraction() = 0.0;
109 inletOutletFvPatchScalarField(ptf),
111 mixingLength_(ptf.mixingLength_),
123 inletOutletFvPatchScalarField(ptf, iF),
125 mixingLength_(ptf.mixingLength_),
145 internalField().
group()
149 Cmu_ = turbModel.coeffDict().getOrDefault<scalar>(
"Cmu", Cmu_);
151 const scalar Cmu75 =
pow(Cmu_, 0.75);
159 this->refValue() = (Cmu75/mixingLength_)*
pow(kp, 1.5);
160 this->valueFraction() = 1.0 -
pos0(phip);
162 inletOutletFvPatchScalarField::updateCoeffs();
172 os.writeEntry(
"mixingLength", mixingLength_);
173 os.writeEntry(
"phi", this->phiName_);
174 os.writeEntry(
"k", kName_);
175 writeEntry(
"value",
os);
184 turbulentMixingLengthDissipationRateInletFvPatchScalarField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
fvPatchField< scalar > fvPatchScalarField
virtual void write(Ostream &) const
virtual void operator=(const UList< Type > &)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A class for handling words, derived from Foam::string.
constexpr const char *const group
fvsPatchField< scalar > fvsPatchScalarField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
static const word propertiesName
virtual void write(Ostream &) const
dimensionedScalar pos0(const dimensionedScalar &ds)
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
virtual const dictionary & coeffDict() const =0
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
T getCheckOrDefault(const word &keyword, const T &deflt, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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)
Macros for easy insertion into run-time selection tables.
This boundary condition provides an inlet condition for turbulent kinetic energy dissipation rate,...
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
static MinMax< T > ge(const T &minVal)
T getCheck(const word &keyword, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
Ostream & writeEntry(const keyType &key, const T &value)
Foam::fvPatchFieldMapper.
static word groupName(StringType base, const word &group)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
A min/max value pair with additional methods. In addition to conveniently storing values,...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
turbulentMixingLengthDissipationRateInletFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
virtual void updateCoeffs()