Go to the documentation of this file.
41 Foam::semiPermeableBaffleVelocityFvPatchVectorField::composition()
const
45 if (db().foundObject<psiReactionThermo>(
name))
47 return db().lookupObject<psiReactionThermo>(
name).composition();
49 else if (db().foundObject<rhoReactionThermo>(
name))
51 return db().lookupObject<rhoReactionThermo>(
name).composition();
56 <<
"Could not find a multi-component thermodynamic model."
59 return NullObjectRef<basicSpecieMixture>();
73 fixedValueFvPatchVectorField(
p, iF),
86 fixedValueFvPatchVectorField(
p, iF),
102 fixedValueFvPatchVectorField(ptf,
p, iF, mapper),
103 rhoName_(ptf.rhoName_)
113 fixedValueFvPatchVectorField(ptf),
114 rhoName_(ptf.rhoName_)
125 fixedValueFvPatchVectorField(ptf, iF),
126 rhoName_(ptf.rhoName_)
151 if (!isA<YBCType>(Yp))
154 <<
"The mass-fraction condition on patch " <<
patch().name()
155 <<
" is not of type " << YBCType::typeName <<
"."
159 phip += refCast<const YBCType>(Yp).phiY();
164 fixedValueFvPatchVectorField::updateCoeffs();
175 writeEntry(
"value", os);
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.
This is a mass-fraction boundary condition for a semi-permeable baffle.
A class for handling words, derived from Foam::string.
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
static constexpr const zero Zero
Global zero (0)
virtual void operator==(const fvPatchField< vector > &)
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
basicSpecieMixture & composition
#define forAll(list, i)
Loop across all elements in list.
Field< vector > vectorField
Specialisation of Field<T> for vector.
This is a velocity boundary condition for a semi-permeable baffle.
semiPermeableBaffleVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
word name(const complex &c)
Return string representation of complex.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
GeometricField< scalar, fvPatchField, volMesh > volScalarField
word dictName() const
The local dictionary name (final part of scoped name)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
PtrList< volScalarField > & Y
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
virtual void write(Ostream &) const
Write.
const std::string patch
OpenFOAM patch number as a std::string.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...