42 { equilibriumModelType::LANGMUIR,
"Langmuir" }
52 { kineticModelType::PseudoFirstOrder,
"PseudoFirstOrder" }
59Foam::speciesSorptionFvPatchScalarField::calcMoleFractions()
const
62 auto& Mole = tMole.ref();
71 const PtrList<volScalarField>&
Y =
thermo.composition().Y();
77 const label speciesId =
83 thermo.composition().W(speciesId)
90 const label
cellId = faceCells[i];
97 <<
"Thermo type is not 'rhoReactionThermo'. " <<
nl
98 <<
"This BC is designed to operate with a rho based thermo."
107Foam::speciesSorptionFvPatchScalarField::field
109 const word& fieldName,
113 const fvMesh&
mesh = this->internalField().mesh();
148 zeroGradientFvPatchScalarField(
p, iF),
149 equilibriumModel_(equilibriumModelType::LANGMUIR),
150 kinematicModel_(kineticModelType::PseudoFirstOrder),
151 thicknessPtr_(nullptr),
169 zeroGradientFvPatchScalarField(
p, iF,
dict),
170 equilibriumModel_(equilibriumModelTypeNames.get(
"equilibriumModel",
dict)),
171 kinematicModel_(kinematicModelTypeNames.get(
"kinematicModel",
dict)),
176 rhoS_(
dict.get<scalar>(
"rhoS")),
177 pName_(
dict.getOrDefault<
word>(
"p",
"p")),
190 const speciesSorptionFvPatchScalarField& ptf,
192 const DimensionedField<scalar, volMesh>& iF,
193 const fvPatchFieldMapper& mapper
196 zeroGradientFvPatchScalarField(ptf,
p, iF, mapper),
197 equilibriumModel_(ptf.equilibriumModel_),
198 kinematicModel_(ptf.kinematicModel_),
199 thicknessPtr_(ptf.thicknessPtr_.clone(patch().patch())),
205 dfldp_(ptf.dfldp_, mapper),
206 mass_(ptf.mass_, mapper)
215 zeroGradientFvPatchScalarField(ptf),
216 equilibriumModel_(ptf.equilibriumModel_),
217 kinematicModel_(ptf.kinematicModel_),
218 thicknessPtr_(ptf.thicknessPtr_.clone(patch().patch())),
235 zeroGradientFvPatchScalarField(ptf, iF),
236 equilibriumModel_(ptf.equilibriumModel_),
237 kinematicModel_(ptf.kinematicModel_),
238 thicknessPtr_(ptf.thicknessPtr_.clone(patch().patch())),
256 zeroGradientFvPatchScalarField::autoMap(m);
263 thicknessPtr_->autoMap(m);
274 zeroGradientFvPatchScalarField::rmap(ptf, addr);
278 dfldp_.rmap(tiptf.dfldp_, addr);
279 mass_.rmap(tiptf.mass_, addr);
283 thicknessPtr_->rmap(tiptf.thicknessPtr_(), addr);
291 const auto&
thermo = db().lookupObject<rhoReactionThermo>
296 const label speciesId =
297 thermo.composition().species()[this->internalField().name()];
299 const scalar Wi(
thermo.composition().W(speciesId));
301 const scalar t = db().time().timeOutputValue();
314 const label faceCelli = this->patch().faceCells()[facei];
315 Vol[facei] = this->internalField().mesh().V()[faceCelli];
326 Info<<
" Patch mass rate min/max [kg/m3/sec]: "
351 switch (equilibriumModel_)
353 case equilibriumModelType::LANGMUIR:
356 tmp<scalarField> tco = calcMoleFractions();
360 cEq = max_*(kl_*tco()*
pp/(1 + kl_*tco()*
pp));
370 switch (kinematicModel_)
372 case kineticModelType::PseudoFirstOrder:
374 dfldp_ = kabs_*(cEq - mass_);
381 const scalar dt = db().time().deltaTValue();
383 mass_ =
max(mass_, scalar(0));
388 "absorbedMass" + this->internalField().
name(),
390 ).boundaryFieldRef()[
patch().index()];
398 Info<<
" Absorption rate min/max [mol/kg/sec]: "
402 zeroGradientFvPatchScalarField::updateCoeffs();
412 "equilibriumModel", equilibriumModelTypeNames[equilibriumModel_]
416 "kinematicModel", kinematicModelTypeNames[kinematicModel_]
420 thicknessPtr_->writeData(
os);
442 speciesSorptionFvPatchScalarField
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...
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
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.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
bool get(const label i) const
label size() const noexcept
static const word dictName
The dictionary name ("thermophysicalProperties").
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Mesh data needed to do the Finite Volume discretisation.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
const Time & time() const
Return the top-level database.
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.
static tmp< fvPatchField< scalar > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< scalar, volMesh > &)
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
const DimensionedField< scalar, volMesh > & internalField() const noexcept
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.
virtual const labelUList & faceCells() const
Return faceCells.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
This boundary condition provides a first-order zero-gradient condition for a given scalar field to mo...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
static const Enum< kineticModelType > kinematicModelTypeNames
Names for kineticModelType.
virtual void write(Ostream &) const
Write.
virtual tmp< fvPatchField< scalar > > clone() const
Return a clone.
speciesSorptionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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.
tmp< scalarField > mass() const
Access to mass.
equilibriumModelType
Options for the equilibrum model.
static const Enum< equilibriumModelType > equilibriumModelTypeNames
Names for equilibriumModelType.
kineticModelType
Options for the kinematic model.
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.
A class for handling words, derived from Foam::string.
PtrList< volScalarField > & Y
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
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.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
List< label > labelList
A List of labels.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
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).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
UList< label > labelUList
A UList of labels.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
psiReactionThermo & thermo
#define forAll(list, i)
Loop across all elements in list.