43 fixedValueFvPatchScalarField(
p, iF),
46 psiName_(
"thermo:psi"),
60 fixedValueFvPatchScalarField(ptf,
p, iF, mapper),
62 phiName_(ptf.phiName_),
63 psiName_(ptf.psiName_),
77 UName_(
dict.getOrDefault<
word>(
"U",
"U")),
78 phiName_(
dict.getOrDefault<
word>(
"phi",
"phi")),
79 psiName_(
dict.getOrDefault<
word>(
"psi",
"thermo:psi")),
80 gamma_(
dict.get<scalar>(
"gamma")),
81 T0_(
"T0",
dict,
p.size())
83 if (!this->readValueEntry(
dict))
92 const totalTemperatureFvPatchScalarField& tppsf
95 fixedValueFvPatchScalarField(tppsf),
97 phiName_(tppsf.phiName_),
98 psiName_(tppsf.psiName_),
110 fixedValueFvPatchScalarField(tppsf, iF),
111 UName_(tppsf.UName_),
112 phiName_(tppsf.phiName_),
113 psiName_(tppsf.psiName_),
114 gamma_(tppsf.gamma_),
126 fixedValueFvPatchScalarField::autoMap(m);
137 fixedValueFvPatchScalarField::rmap(ptf, addr);
142 T0_.rmap(tiptf.T0_, addr);
162 scalar gM1ByG = (gamma_ - 1.0)/gamma_;
169 fixedValueFvPatchScalarField::updateCoeffs();
176 os.writeEntryIfDifferent<
word>(
"U",
"U", UName_);
177 os.writeEntryIfDifferent<
word>(
"phi",
"phi", phiName_);
178 os.writeEntryIfDifferent<
word>(
"psi",
"thermo:psi", psiName_);
179 os.writeEntry(
"gamma", gamma_);
180 T0_.writeEntry(
"T0",
os);
192 totalTemperatureFvPatchScalarField
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual void operator=(const UList< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
This boundary condition provides a total temperature condition.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
totalTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void write(Ostream &) const
Write.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from Foam::string.
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensionedScalar neg(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
fvPatchField< scalar > fvPatchScalarField