42 jumpCyclicFvPatchScalarField(
p, iF,
dict),
46 ?
dict.lookupOrDefault<
word>(
"rho",
"rho")
57 if (
dict.found(
"jump"))
62 evaluateNoUpdateCoeffs();
74 jumpCyclicFvPatchScalarField(ptf,
p, iF, mapper),
75 rhoName_(ptf.rhoName_),
77 jump_(mapper(ptf.jump_))
87 jumpCyclicFvPatchScalarField(ptf, iF),
88 rhoName_(ptf.rhoName_),
109 jumpCyclicFvPatchScalarField::map(ptf, mapper);
112 refCast<const prghCyclicPressureFvPatchScalarField>(ptf);
114 mapper(jump_, tiptf.jump_);
123 jumpCyclicFvPatchScalarField::reset(ptf);
126 refCast<const prghCyclicPressureFvPatchScalarField>(ptf);
128 jump_.reset(tiptf.jump_);
143 refCast<const prghCyclicPressureFvPatchScalarField>(nbrPatchField());
145 const label patchi = patch().index(), nbrPatchi = nbrPf.patch().index();
150 (cyclicPatch().owner() ? rhoName_ : nbrPf.rhoName_);
166 vf.
mesh().schemes().snGrad(internalField().
name())
169 deltaCoeffsSf->boundaryField()[patch().index()];
173 (rhoBf[
patchi] - (cyclicPatch().owner() ? rhoInf_ : nbrPf.rhoInf_))
174 *(ghfBf[nbrPatchi] - ghfBf[
patchi])
177 + phiHbyABf[nbrPatchi]/rAUfBf[nbrPatchi]/nbrPf.patch().magSf()
178 )*patch().weights()/deltaCoeffsPf;
180 jumpCyclicFvPatchScalarField::updateCoeffs();
191 if (cyclicPatch().owner())
193 writeEntryIfDifferent<word>(os,
"rho",
"rho", rhoName_);
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...
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
Generic GeometricBoundaryField class.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
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 any number of values (e....
Abstract base class for field mapping.
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 finiteVolume patch using a polyPatch and a fvBoundaryMesh.
static tmp< snGradScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
This boundary condition provides a cyclic condition for p_rgh. It applies corrections to the value an...
virtual tmp< scalarField > jump() const
Return the "jump".
virtual void write(Ostream &) const
Write.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
prghCyclicPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the patch pressure gradient field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
A class for managing temporary objects.
A class for handling words, derived from string.
Calculate the snGrad of the given volField.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
word name(const bool)
Return a word representation of a bool.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
const dimensionSet dimDensity
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)