Go to the documentation of this file.
41 fixedValueFvPatchScalarField(
p, iF),
55 fixedValueFvPatchScalarField(
p, iF),
56 rhoName_(
dict.lookupOrDefault<
word>(
"rhoName",
"rho")),
57 p_(
"p",
dict,
p.size())
59 if (
dict.found(
"value"))
61 fvPatchScalarField::operator=
82 fixedValueFvPatchScalarField(ptf,
p, iF, mapper),
94 fixedValueFvPatchScalarField(ptf),
107 fixedValueFvPatchScalarField(ptf, iF),
120 fixedValueFvPatchScalarField::autoMap(m);
131 fixedValueFvPatchScalarField::rmap(ptf, addr);
134 refCast<const prghPressureFvPatchScalarField>(ptf);
136 p_.rmap(tiptf.
p_, addr);
160 mag(
g.value()) > SMALL
165 operator==(
p_ - rhop*((
g.value() & patch().Cf()) - ghRef.value()));
167 fixedValueFvPatchScalarField::updateCoeffs();
174 if (rhoName_ !=
"rho")
179 p_.writeEntry(
"p", os);
180 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 rhoName_
Name of phase-fraction field.
A class for handling words, derived from string.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const dimensionedVector & g
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
dimensioned< scalar > mag(const dimensioned< Type > &)
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Pre-declare SubField and related Field type.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
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.
Generic dimensioned Type class.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
virtual void write(Ostream &) const
Write.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
volScalarField scalarField(fieldObject, mesh)
scalarField p_
Static pressure.
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
prghPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Foam::fvPatchFieldMapper.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
This boundary condition provides static pressure condition for p_rgh, calculated as:
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...