Go to the documentation of this file.
41 fluctuationScale_(
Zero),
42 referenceField_(
p.size()),
58 fluctuationScale_(
dict.
get<Type>(
"fluctuationScale")),
59 referenceField_(
"referenceField",
dict,
p.size()),
60 alpha_(
dict.getOrDefault<scalar>(
"alpha", 0.1)),
63 if (
dict.found(
"value"))
80 const turbulentInletFvPatchField<Type>& ptf,
82 const DimensionedField<Type, volMesh>& iF,
83 const fvPatchFieldMapper& mapper
86 fixedValueFvPatchField<Type>(ptf,
p, iF, mapper),
88 fluctuationScale_(ptf.fluctuationScale_),
89 referenceField_(ptf.referenceField_, mapper),
102 ranGen_(ptf.ranGen_),
103 fluctuationScale_(ptf.fluctuationScale_),
104 referenceField_(ptf.referenceField_),
118 ranGen_(ptf.ranGen_),
119 fluctuationScale_(ptf.fluctuationScale_),
120 referenceField_(ptf.referenceField_),
135 referenceField_.autoMap(m);
149 refCast<const turbulentInletFvPatchField<Type>>(ptf);
151 referenceField_.
rmap(tiptf.referenceField_, addr);
163 if (curTimeIndex_ != this->db().time().
timeIndex())
165 Field<Type>& patchField = *
this;
167 Field<Type> randomField(this->size());
171 ranGen_.randomise01<Type>(randomField[facei]);
176 scalar rmsCorr =
sqrt(12*(2*alpha_ -
sqr(alpha_)))/alpha_;
179 (1 - alpha_)*patchField
185 randomField - 0.5*pTraits<Type>::one,
187 )*
mag(referenceField_)
190 curTimeIndex_ = this->db().time().timeIndex();
204 this->writeEntry(
"value",
os);
virtual void write(Ostream &) const
virtual void operator==(const fvPatchField< Type > &)
virtual void rmap(const fvPatchField< Type > &, const labelList &)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
static constexpr const zero Zero
turbulentInletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
virtual void write(Ostream &) const
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
OBJstream os(runTime.globalPath()/outputName)
virtual void updateCoeffs()
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual void autoMap(const fvPatchFieldMapper &)
This boundary condition produces spatiotemporal-variant field by summing a set of pseudo-random numbe...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
dimensionedScalar sqrt(const dimensionedScalar &ds)
A traits class, which is primarily used for primitives.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Ostream & writeEntry(const keyType &key, const T &value)
void write(vtk::formatter &fmt, const Type &val, const label n=1)
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,...
virtual void updateCoeffs()
virtual void autoMap(const fvPatchFieldMapper &)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...