Go to the documentation of this file.
44 fixedValueFvPatchScalarField(
p, iF),
61 fixedValueFvPatchScalarField(
p, iF,
dict, false),
62 UName_(
dict.getOrDefault<
word>(
"U",
"U")),
63 phiName_(
dict.getOrDefault<
word>(
"phi",
"phi")),
64 rhoName_(
dict.getOrDefault<
word>(
"rho",
"rho")),
65 psiName_(
dict.getOrDefault<
word>(
"psi",
"none")),
66 gamma_(psiName_ !=
"none" ?
dict.
get<scalar>(
"gamma") : 1),
67 p0_(
"p0",
dict,
p.size())
69 if (
dict.found(
"value"))
85 const totalPressureFvPatchScalarField& ptf,
87 const DimensionedField<scalar, volMesh>& iF,
88 const fvPatchFieldMapper& mapper
91 fixedValueFvPatchScalarField(ptf,
p, iF, mapper),
93 phiName_(ptf.phiName_),
94 rhoName_(ptf.rhoName_),
95 psiName_(ptf.psiName_),
106 fixedValueFvPatchScalarField(tppsf),
107 UName_(tppsf.UName_),
108 phiName_(tppsf.phiName_),
109 rhoName_(tppsf.rhoName_),
110 psiName_(tppsf.psiName_),
111 gamma_(tppsf.gamma_),
122 fixedValueFvPatchScalarField(tppsf, iF),
123 UName_(tppsf.UName_),
124 phiName_(tppsf.phiName_),
125 rhoName_(tppsf.rhoName_),
126 psiName_(tppsf.psiName_),
127 gamma_(tppsf.gamma_),
139 fixedValueFvPatchScalarField::autoMap(m);
150 fixedValueFvPatchScalarField::rmap(ptf, addr);
153 refCast<const totalPressureFvPatchScalarField>(ptf);
155 p0_.rmap(tiptf.p0_, addr);
170 const fvsPatchField<scalar>& phip =
175 if (psiName_ ==
"none")
179 const fvPatchField<scalar>&
rho =
188 const fvPatchField<scalar>& psip =
193 scalar gM1ByG = (gamma_ - 1)/gamma_;
200 (1.0 + 0.5*psip*gM1ByG*(1.0 -
pos0(phip))*
magSqr(Up)),
220 <<
" Incorrect pressure dimensions " << internalField().dimensions()
223 <<
" for compressible/variable density flow" <<
nl
225 <<
" for incompressible flow," <<
nl
226 <<
" on patch " << this->
patch().name()
227 <<
" of field " << this->internalField().name()
228 <<
" in file " << this->internalField().objectPath()
232 fixedValueFvPatchScalarField::updateCoeffs();
241 patch().lookupPatchField<volVectorField, vector>(
UName())
255 writeEntry(
"value",
os);
266 totalPressureFvPatchScalarField
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
virtual void operator=(const UList< Type > &)
const dimensionSet dimPressure
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.
static constexpr const zero Zero
const dimensionSet dimDensity
const word UName(propsDict.getOrDefault< word >("U", "U"))
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:
Generic templated field type.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
virtual void autoMap(const fvPatchFieldMapper &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
const word & UName() const
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.
const scalarField & p0() const
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
virtual void updateCoeffs()
Ostream & writeEntry(const keyType &key, const T &value)
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
virtual void write(Ostream &) const
Foam::fvPatchFieldMapper.
totalPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
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...
virtual void rmap(const fvPatchScalarField &, const labelList &)