Go to the documentation of this file.
44 fixedValueFvPatchVectorField(
p, iF),
56 fixedValueFvPatchVectorField(
p, iF,
dict,
false),
71 fixedValueFvPatchVectorField(ptf,
p, iF, mapper),
81 fixedValueFvPatchVectorField(ptf),
92 fixedValueFvPatchVectorField(ptf, iF),
111 internalField().
group()
119 vector tauHat = tau0_/(
mag(tau0_) + ROOTVSMALL);
123 operator==(tauHat*(tauHat & (tau0_*(1.0/(ry*nuEff)) + Uc)));
125 fixedValueFvPatchVectorField::updateCoeffs();
133 writeEntry(
"value", os);
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
fixedShearStressFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
constexpr const char *const group
Group name for atomic constants.
Set a constant shear stress as tau0 = -nuEff dU/dn.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
static const word propertiesName
Default name of the turbulence properties dictionary.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
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,...
Macros for easy insertion into run-time selection tables.
virtual void write(Ostream &) const
Write.
const std::string patch
OpenFOAM patch number as a std::string.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Foam::fvPatchFieldMapper.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...