Go to the documentation of this file.
43 fixedValueFvPatchScalarField(
p, iF),
60 fixedValueFvPatchScalarField(
p, iF,
dict,
false),
61 UName_(
dict.lookupOrDefault<
word>(
"U",
"U")),
62 phiName_(
dict.lookupOrDefault<
word>(
"phi",
"phi")),
63 rhoName_(
dict.lookupOrDefault<
word>(
"rho",
"rho")),
64 psiName_(
dict.lookupOrDefault<
word>(
"psi",
"none")),
65 gamma_(psiName_ !=
"none" ?
dict.get<scalar>(
"gamma") : 1),
66 p0_(
"p0",
dict,
p.size())
68 if (
dict.found(
"value"))
90 fixedValueFvPatchScalarField(ptf,
p, iF, mapper),
92 phiName_(ptf.phiName_),
93 rhoName_(ptf.rhoName_),
94 psiName_(ptf.psiName_),
105 fixedValueFvPatchScalarField(tppsf),
106 UName_(tppsf.UName_),
107 phiName_(tppsf.phiName_),
108 rhoName_(tppsf.rhoName_),
109 psiName_(tppsf.psiName_),
110 gamma_(tppsf.gamma_),
121 fixedValueFvPatchScalarField(tppsf, iF),
122 UName_(tppsf.UName_),
123 phiName_(tppsf.phiName_),
124 rhoName_(tppsf.rhoName_),
125 psiName_(tppsf.psiName_),
126 gamma_(tppsf.gamma_),
138 fixedValueFvPatchScalarField::autoMap(m);
149 fixedValueFvPatchScalarField::rmap(ptf, addr);
152 refCast<const totalPressureFvPatchScalarField>(ptf);
154 p0_.rmap(tiptf.p0_, addr);
174 if (psiName_ ==
"none")
192 scalar gM1ByG = (gamma_ - 1)/gamma_;
199 (1.0 + 0.5*psip*gM1ByG*(1.0 -
pos0(phip))*
magSqr(Up)),
219 <<
" Incorrect pressure dimensions " << internalField().dimensions()
222 <<
" for compressible/variable density flow" <<
nl
224 <<
" for incompressible flow," <<
nl
225 <<
" on patch " << this->
patch().name()
226 <<
" of field " << this->internalField().name()
227 <<
" in file " << this->internalField().objectPath()
231 fixedValueFvPatchScalarField::updateCoeffs();
240 patch().lookupPatchField<volVectorField, vector>(
UName())
254 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.
static constexpr const zero Zero
Global zero.
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)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
This boundary condition provides a total pressure condition. Four variants are possible:
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
const word & UName() const
Return the name of the velocity field.
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.
const scalarField & p0() const
Return the total pressure.
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.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
virtual void write(Ostream &) const
Write.
Foam::fvPatchFieldMapper.
totalPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.