Go to the documentation of this file.
44 fixedValueFvPatchScalarField(
p, iF),
62 fixedValueFvPatchScalarField(
p, iF,
dict, false),
63 UName_(
dict.getOrDefault<
word>(
"U",
"U")),
64 phiName_(
dict.getOrDefault<
word>(
"phi",
"phi")),
65 rhoName_(
dict.getOrDefault<
word>(
"rho",
"rho")),
66 psiName_(
dict.getOrDefault<
word>(
"psi",
"none")),
67 gamma_(psiName_ !=
"none" ?
dict.
get<scalar>(
"gamma") : 1),
70 if (
dict.found(
"value"))
79 const scalar t = this->db().time().timeOutputValue();
88 const uniformTotalPressureFvPatchScalarField& ptf,
90 const DimensionedField<scalar, volMesh>& iF,
91 const fvPatchFieldMapper& mapper
94 fixedValueFvPatchScalarField(
p, iF),
96 phiName_(ptf.phiName_),
97 rhoName_(ptf.rhoName_),
98 psiName_(ptf.psiName_),
102 patchType() = ptf.patchType();
106 const scalar t = this->db().time().timeOutputValue();
117 fixedValueFvPatchScalarField(ptf),
119 phiName_(ptf.phiName_),
120 rhoName_(ptf.rhoName_),
121 psiName_(ptf.psiName_),
134 fixedValueFvPatchScalarField(ptf, iF),
136 phiName_(ptf.phiName_),
137 rhoName_(ptf.rhoName_),
138 psiName_(ptf.psiName_),
156 scalar
p0 = p0_->value(this->db().time().timeOutputValue());
158 const fvsPatchField<scalar>& phip =
163 if (psiName_ ==
"none")
167 const fvPatchField<scalar>&
rho =
176 const fvPatchField<scalar>& psip =
181 scalar gM1ByG = (gamma_ - 1)/gamma_;
188 (1.0 + 0.5*psip*gM1ByG*(1.0 -
pos0(phip))*
magSqr(Up)),
208 <<
" Incorrect pressure dimensions " << internalField().dimensions()
211 <<
" for compressible/variable density flow" <<
nl
213 <<
" for incompressible flow," <<
nl
214 <<
" on patch " << this->
patch().name()
215 <<
" of field " << this->internalField().name()
216 <<
" in file " << this->internalField().objectPath()
220 fixedValueFvPatchScalarField::updateCoeffs();
226 updateCoeffs(
patch().lookupPatchField<volVectorField, vector>(UName_));
239 writeEntry(
"value",
os);
250 uniformTotalPressureFvPatchScalarField
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
const dimensionSet dimPressure
virtual void operator==(const fvPatchField< Type > &)
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A class for handling words, derived from Foam::string.
const dimensionSet dimDensity
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
dimensionedScalar pos0(const dimensionedScalar &ds)
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Generic templated field type.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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.
errorManipArg< error, int > exit(error &err, const int errNo=1)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
#define FatalErrorInFunction
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,...
const volScalarField & p0
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...