Go to the documentation of this file.
44 inletOutletFvPatchScalarField(
p, iF),
50 this->refValue() =
Zero;
51 this->refGrad() =
Zero;
52 this->valueFraction() = 0.0;
65 inletOutletFvPatchScalarField(ptf,
p, iF, mapper),
67 psiName_(ptf.psiName_),
81 inletOutletFvPatchScalarField(
p, iF),
82 UName_(
dict.getOrDefault<
word>(
"U",
"U")),
83 psiName_(
dict.getOrDefault<
word>(
"psi",
"thermo:psi")),
84 gamma_(
dict.get<scalar>(
"gamma")),
85 T0_(
"T0",
dict,
p.size())
87 this->patchType() =
dict.getOrDefault<
word>(
"patchType", word::null);
89 this->phiName_ =
dict.getOrDefault<
word>(
"phi",
"phi");
91 this->refValue() =
Zero;
92 if (
dict.found(
"value"))
104 this->refGrad() =
Zero;
105 this->valueFraction() = 0.0;
115 inletOutletFvPatchScalarField(tppsf),
116 UName_(tppsf.UName_),
117 psiName_(tppsf.psiName_),
118 gamma_(tppsf.gamma_),
130 inletOutletFvPatchScalarField(tppsf, iF),
131 UName_(tppsf.UName_),
132 psiName_(tppsf.psiName_),
133 gamma_(tppsf.gamma_),
145 inletOutletFvPatchScalarField::autoMap(m);
156 inletOutletFvPatchScalarField::rmap(ptf, addr);
159 refCast<const inletOutletTotalTemperatureFvPatchScalarField>(ptf);
161 T0_.rmap(tiptf.T0_, addr);
181 scalar gM1ByG = (gamma_ - 1.0)/gamma_;
184 T0_/(1.0 + 0.5*psip*gM1ByG*(1.0 -
pos0(phip))*
magSqr(Up));
185 this->valueFraction() = 1.0 -
pos0(phip);
187 inletOutletFvPatchScalarField::updateCoeffs();
200 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.
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.
This boundary condition provides an outflow condition for total temperature for use with supersonic c...
static constexpr const zero Zero
Global zero (0)
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
dimensionedScalar pos0(const dimensionedScalar &ds)
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
inletOutletTotalTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
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.
GeometricField< vector, fvPatchField, volMesh > volVectorField
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,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...