Go to the documentation of this file.
44 inletOutletFvPatchScalarField(
p, iF),
48 this->refValue() = 0.0;
49 this->refGrad() = 0.0;
50 this->valueFraction() = 0.0;
62 inletOutletFvPatchScalarField(ptf,
p, iF, mapper),
63 intensity_(ptf.intensity_),
75 inletOutletFvPatchScalarField(
p, iF),
76 intensity_(
dict.
get<scalar>(
"intensity")),
77 UName_(
dict.getOrDefault<
word>(
"U",
"U"))
80 this->phiName_ =
dict.getOrDefault<
word>(
"phi",
"phi");
82 if (intensity_ < 0 || intensity_ > 1)
85 <<
"Turbulence intensity should be specified as a fraction 0-1 "
86 "of the mean velocity\n"
87 " value given is " << intensity_ <<
nl
88 <<
" on patch " << this->
patch().name()
89 <<
" of field " << this->internalField().name()
90 <<
" in file " << this->internalField().objectPath()
96 this->refValue() = 0.0;
97 this->refGrad() = 0.0;
98 this->valueFraction() = 0.0;
107 inletOutletFvPatchScalarField(ptf),
108 intensity_(ptf.intensity_),
120 inletOutletFvPatchScalarField(ptf, iF),
121 intensity_(ptf.intensity_),
142 this->refValue() = 1.5*
sqr(intensity_)*
magSqr(Up);
143 this->valueFraction() = 1.0 -
pos0(phip);
145 inletOutletFvPatchScalarField::updateCoeffs();
158 writeEntry(
"value",
os);
169 turbulentIntensityKineticEnergyInletFvPatchScalarField
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 > &)
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void updateCoeffs()
A class for handling words, derived from Foam::string.
fvsPatchField< scalar > fvsPatchScalarField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
turbulentIntensityKineticEnergyInletFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
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.
Vector< scalar > vector
A scalar version of the templated Vector.
errorManipArg< error, int > exit(error &err, const int errNo=1)
GeometricField< vector, fvPatchField, volMesh > volVectorField
fvPatchField< vector > fvPatchVectorField
#define FatalErrorInFunction
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & writeEntry(const keyType &key, const T &value)
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual void write(Ostream &) const
This boundary condition provides a turbulent kinetic energy condition, based on user-supplied turbule...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...