Go to the documentation of this file.
40 mixedFvPatchVectorField(
p, iF),
43 psiName_(
"thermo:psi"),
49 refValue() = patchInternalField();
50 refGrad() = vector::zero;
64 mixedFvPatchVectorField(ptf,
p, iF, mapper),
83 mixedFvPatchVectorField(
p, iF),
84 TName_(
dict.lookupOrDefault<
word>(
"T",
"T")),
85 pName_(
dict.lookupOrDefault<
word>(
"p",
"p")),
86 psiName_(
dict.lookupOrDefault<
word>(
"psi",
"thermo:psi")),
87 UInf_(
dict.lookup(
"UInf")),
92 if (
dict.found(
"value"))
105 refGrad() = vector::zero;
113 ) <<
" unphysical pInf specified (pInf <= 0.0)"
114 <<
"\n on patch " << this->patch().name()
128 mixedFvPatchVectorField(sfspvf),
146 mixedFvPatchVectorField(sfspvf, iF),
161 if (!size() || updated())
177 scalar
R = 1.0/(ppsi[0]*pT[0]);
184 <<
"\n on patch " << this->patch().name()
229 if (pp[facei] >=
pInf_)
238 Up[facei] = Ut[facei] + fpp*nHatInf[facei];
247 Up[facei] =
U[facei];
248 valueFraction()[facei] = 0;
274 scalar fpp = (nuMachInf - nuMachf)*
mag(Ut[facei]);
276 Up[facei] = Ut[facei] + fpp*nHatInf[facei];
281 <<
"unphysical subsonic inflow has been generated"
282 <<
"\n on patch " << this->patch().name()
291 mixedFvPatchVectorField::updateCoeffs();
298 writeEntryIfDifferent<word>(os,
"T",
"T", TName_);
299 writeEntryIfDifferent<word>(os,
"p",
"p", pName_);
300 writeEntryIfDifferent<word>(os,
"psi",
"thermo:psi", psiName_);
305 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.
word psiName_
Name of compressibility field field, default = "thermo:psi".
A class for handling words, derived from string.
virtual void write(Ostream &) const
Write.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
This boundary condition provides a supersonic free-stream condition.
rDeltaT dimensionedInternalField()
scalar pInf_
Pressure of the free stream.
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar TInf_
Temperature of the free stream.
#define R(A, B, C, D, E, F, K, M)
Pre-declare SubField and related Field type.
scalar gamma_
Heat capacity ratio.
word pName_
Name of pressure field, default = "p".
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
volVectorField vectorField(fieldObject, mesh)
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....
word TName_
Name of temperature field, default = "T".
dimensionedScalar log(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
errorManipArg< error, int > exit(error &err, const int errNo=1)
vector UInf_
Velocity of the free stream.
supersonicFreestreamFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
dimensionedScalar atan(const dimensionedScalar &ds)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...