Go to the documentation of this file.
43 fixedValueFvPatchScalarField(
p, iF),
61 fixedValueFvPatchScalarField(
p, iF,
dict,
false),
62 UName_(
dict.lookupOrDefault<
word>(
"U",
"U")),
63 phiName_(
dict.lookupOrDefault<
word>(
"phi",
"phi")),
64 rhoName_(
dict.lookupOrDefault<
word>(
"rho",
"rho")),
65 psiName_(
dict.lookupOrDefault<
word>(
"psi",
"none")),
66 gamma_(psiName_ !=
"none" ?
dict.get<scalar>(
"gamma") : 1),
69 if (
dict.found(
"value"))
78 const scalar t = this->db().time().timeOutputValue();
93 fixedValueFvPatchScalarField(
p, iF),
95 phiName_(ptf.phiName_),
96 rhoName_(ptf.rhoName_),
97 psiName_(ptf.psiName_),
101 patchType() = ptf.patchType();
105 const scalar t = this->db().time().timeOutputValue();
116 fixedValueFvPatchScalarField(ptf),
118 phiName_(ptf.phiName_),
119 rhoName_(ptf.rhoName_),
120 psiName_(ptf.psiName_),
133 fixedValueFvPatchScalarField(ptf, iF),
135 phiName_(ptf.phiName_),
136 rhoName_(ptf.rhoName_),
137 psiName_(ptf.psiName_),
155 scalar
p0 = p0_->value(this->db().time().timeOutputValue());
162 if (psiName_ ==
"none")
180 scalar gM1ByG = (gamma_ - 1)/gamma_;
187 (1.0 + 0.5*psip*gM1ByG*(1.0 -
pos0(phip))*
magSqr(Up)),
207 <<
" Incorrect pressure dimensions " << internalField().dimensions()
210 <<
" for compressible/variable density flow" <<
nl
212 <<
" for incompressible flow," <<
nl
213 <<
" on patch " << this->
patch().name()
214 <<
" of field " << this->internalField().name()
215 <<
" in file " << this->internalField().objectPath()
219 fixedValueFvPatchScalarField::updateCoeffs();
238 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.
const dimensionSet dimPressure
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
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)
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,...
Macros for easy insertion into run-time selection tables.
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const std::string patch
OpenFOAM patch number as a std::string.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
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...