Go to the documentation of this file.
46 mixedFvPatchScalarField(
p, iF),
52 valueFraction() =
Zero;
65 mixedFvPatchScalarField(
p, iF),
66 c_(
dict.getOrDefault<scalar>(
"c", 0)),
67 phiName_(
dict.getOrDefault<
word>(
"phi",
"phi"))
73 valueFraction() =
Zero;
87 mixedFvPatchScalarField(ptf,
p, iF, mapper),
89 phiName_(ptf.phiName_)
100 mixedFvPatchScalarField(ptf),
102 phiName_(ptf.phiName_)
114 mixedFvPatchScalarField(ptf, iF),
116 phiName_(ptf.phiName_)
130 const word& YName = internalField().name();
132 const label nbrPatchi = samplePolyPatch().index();
133 const fvPatch& nbrPatch =
patch().boundaryMesh()[nbrPatchi];
140 return c_*
patch().magSf()*(patchInternalField() - nbrYc);
162 valueFraction() = phip/(phip -
patch().deltaCoeffs()*AMuEffp);
163 refGrad() = - phiY()/AMuEffp;
165 mixedFvPatchScalarField::updateCoeffs();
178 writeEntry(
"value",
os);
189 semiPermeableBaffleMassFractionFvPatchScalarField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
fvPatchField< scalar > fvPatchScalarField
virtual void write(Ostream &) const
virtual void operator=(const UList< Type > &)
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
This is a mass-fraction boundary condition for a semi-permeable baffle.
virtual void write(Ostream &) const
A class for handling words, derived from Foam::string.
tmp< scalarField > phiY() const
A class for managing temporary objects.
static constexpr const zero Zero
const mapDistribute & map() const
const polyPatch & samplePolyPatch() const
static const word propertiesName
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Generic templated field type.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
virtual void updateCoeffs()
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
virtual tmp< Field< Type > > patchInternalField() const
semiPermeableBaffleMassFractionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
Abstract base class for turbulence models (RAS, LES and laminar).
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.
label index() const noexcept
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual void write(Ostream &os) const
static tmp< T > New(Args &&... args)
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...