Go to the documentation of this file.
41 fixedValueFvPatchScalarField(
p, iF),
59 fixedValueFvPatchScalarField(
p, iF),
60 UName_(
dict.lookupOrDefault<
word>(
"U",
"U")),
61 phiName_(
dict.lookupOrDefault<
word>(
"phi",
"phi")),
62 rhoName_(
dict.lookupOrDefault<
word>(
"rho",
"none")),
63 psiName_(
dict.lookupOrDefault<
word>(
"psi",
"none")),
67 if (
dict.found(
"value"))
76 const scalar t = this->db().time().timeOutputValue();
91 fixedValueFvPatchScalarField(
p, iF),
99 patchType() = ptf.patchType();
103 const scalar t = this->db().time().timeOutputValue();
114 fixedValueFvPatchScalarField(ptf),
131 fixedValueFvPatchScalarField(ptf, iF),
153 scalar p0 = pressure_->value(this->db().time().timeOutputValue());
158 if (psiName_ ==
"none" && rhoName_ ==
"none")
162 else if (rhoName_ ==
"none")
169 scalar gM1ByG = (gamma_ - 1.0)/gamma_;
176 (1.0 + 0.5*psip*gM1ByG*(1.0 -
pos(phip))*
magSqr(Up)),
186 else if (psiName_ ==
"none")
196 <<
" rho or psi set inconsitently, rho = " << rhoName_
197 <<
", psi = " << psiName_ <<
".\n"
198 <<
" Set either rho or psi or neither depending on the "
199 "definition of total pressure.\n"
200 <<
" Set the unused variables to 'none'.\n"
201 <<
" on patch " << this->patch().name()
207 fixedValueFvPatchScalarField::updateCoeffs();
220 writeEntryIfDifferent<word>(os,
"U",
"U", UName_);
221 writeEntryIfDifferent<word>(os,
"phi",
"phi", phiName_);
225 pressure_->writeData(os);
226 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.
A class for handling words, derived from string.
rDeltaT dimensionedInternalField()
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Pre-declare SubField and related Field 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 any number of values (e....
Macros for easy insertion into run-time selection tables.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
volScalarField scalarField(fieldObject, mesh)
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensionedScalar pos(const dimensionedScalar &ds)