39 Foam::enthalpySorptionFvPatchScalarField::enthalpyModelType
41Foam::enthalpySorptionFvPatchScalarField::enthalpyModelTypeNames
43 { enthalpyModelType::estimated,
"estimated" },
44 { enthalpyModelType::calculated,
"calculated" }
53 const DimensionedField<scalar, volMesh>& iF
56 zeroGradientFvPatchScalarField(
p, iF),
57 enthalpyModel_(enthalpyModelType::estimated),
59 enthalpyMassLoadPtr_(nullptr),
76 zeroGradientFvPatchScalarField(
p, iF,
dict),
77 enthalpyModel_(enthalpyModelTypeNames.get(
"enthalpyModel",
dict)),
78 includeHs_(
dict.getOrDefault<bool>(
"includeHs", true)),
79 enthalpyMassLoadPtr_(nullptr),
82 speciesName_(
dict.get<
word>(
"species")),
83 pName_(
dict.getOrDefault<
word>(
"p",
"p")),
84 TName_(
dict.getOrDefault<
word>(
"T",
"T")),
87 switch (enthalpyModel_)
89 case enthalpyModelType::calculated:
91 enthalpyMassLoadPtr_ =
92 Function1<scalar>::New(
"enthalpyTable",
dict, &iF.
db());
95 case enthalpyModelType::estimated:
110 const enthalpySorptionFvPatchScalarField& ptf,
112 const DimensionedField<scalar, volMesh>& iF,
113 const fvPatchFieldMapper& mapper
116 zeroGradientFvPatchScalarField(ptf,
p, iF, mapper),
117 enthalpyModel_(ptf.enthalpyModel_),
118 includeHs_(ptf.includeHs_),
119 enthalpyMassLoadPtr_(ptf.enthalpyMassLoadPtr_.clone()),
122 speciesName_(ptf.speciesName_),
125 dhdt_(ptf.dhdt_, mapper)
134 zeroGradientFvPatchScalarField(ptf),
135 enthalpyModel_(ptf.enthalpyModel_),
136 includeHs_(ptf.includeHs_),
137 enthalpyMassLoadPtr_(ptf.enthalpyMassLoadPtr_.clone()),
140 speciesName_(ptf.speciesName_),
153 zeroGradientFvPatchScalarField(ptf, iF),
154 enthalpyModel_(ptf.enthalpyModel_),
155 includeHs_(ptf.includeHs_),
156 enthalpyMassLoadPtr_(ptf.enthalpyMassLoadPtr_.clone()),
159 speciesName_(ptf.speciesName_),
173 zeroGradientFvPatchScalarField::autoMap(m);
185 zeroGradientFvPatchScalarField::rmap(ptf, addr);
189 dhdt_.rmap(tiptf.dhdt_, addr);
199 patch().lookupPatchField<volScalarField>(speciesName_)
213 const auto& Tp = patch().lookupPatchField<
volScalarField>(TName_);
222 const label speciesId =
223 thermo.composition().species()[speciesName_];
229 hsp[facei] =
composition.Hs(speciesId,
pp[facei], Tp[facei]);
239 Info<<
" Patch enthalpy rate min/max [J/m3/sec]: "
257 patch().lookupPatchField<volScalarField>(speciesName_)
260 switch (enthalpyModel_)
262 case enthalpyModelType::estimated:
267 case enthalpyModelType::calculated:
275 const scalar mfacei = massb[facei];
277 dhdt_[facei] = enthalpyMassLoadPtr_->value(mfacei);
289 Info<<
" Enthalpy change min/max [J/kg]: "
293 zeroGradientFvPatchScalarField::updateCoeffs();
301 os.writeEntry(
"enthalpyModel", enthalpyModelTypeNames[enthalpyModel_]);
303 if (enthalpyMassLoadPtr_)
305 enthalpyMassLoadPtr_->writeData(
os);
316 dhdt_.writeEntry(
"dhdt",
os);
329 enthalpySorptionFvPatchScalarField
Macros for easy insertion into run-time selection tables.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
const objectRegistry & db() const noexcept
Return the local objectRegistry.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
bool get(const label i) const
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
static const word dictName
The dictionary name ("thermophysicalProperties").
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
This is a temperature boundary condition which works in conjunction with the speciesSorption conditio...
enthalpySorptionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void write(Ostream &) const
Write.
virtual tmp< fvPatchField< scalar > > clone() const
Return a clone.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual tmp< scalarField > patchSource() const
Source of cells next to the patch.
const objectRegistry & db() const
The associated objectRegistry.
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual void operator=(const UList< Type > &)
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
A class for handling words, derived from Foam::string.
basicSpecieMixture & composition
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Namespace for handling debugging switches.
const std::string patch
OpenFOAM patch number as a std::string.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
List< label > labelList
A List of labels.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
MinMax< scalar > scalarMinMax
A scalar min/max range.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
fvPatchField< scalar > fvPatchScalarField
#define forAll(list, i)
Loop across all elements in list.