Loading...
Searching...
No Matches
speciesSorptionFvPatchScalarField.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
30#include "fvPatchFieldMapper.H"
31#include "volFields.H"
32#include "rhoReactionThermo.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36const Foam::Enum
37<
39>
41({
42 { equilibriumModelType::LANGMUIR, "Langmuir" }
43});
44
45
46const Foam::Enum
47<
49>
51({
52 { kineticModelType::PseudoFirstOrder, "PseudoFirstOrder" }
53});
54
55
56// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
57
59Foam::speciesSorptionFvPatchScalarField::calcMoleFractions() const
60{
61 auto tMole = tmp<scalarField>::New(patch().size(), Zero);
62 auto& Mole = tMole.ref();
63
64 if (db().foundObject<rhoReactionThermo>(basicThermo::dictName))
65 {
66 const auto& thermo = db().lookupObject<rhoReactionThermo>
67 (
69 );
70
71 const PtrList<volScalarField>& Y = thermo.composition().Y();
72
73 const volScalarField W(thermo.W());
74
75 const labelUList& faceCells = patch().faceCells();
76
77 const label speciesId =
78 thermo.composition().species()[this->internalField().name()];
79
80 const dimensionedScalar Wi
81 (
83 thermo.composition().W(speciesId)
84 );
85
86 const volScalarField X(W*Y[speciesId]/Wi);
87
88 forAll(faceCells, i)
89 {
90 const label cellId = faceCells[i];
91 Mole[i] = X[cellId];
92 }
93 }
94 else
95 {
97 << "Thermo type is not 'rhoReactionThermo'. " << nl
98 << "This BC is designed to operate with a rho based thermo."
99 << exit(FatalError);
100 }
101
102 return tMole;
103}
104
105
107Foam::speciesSorptionFvPatchScalarField::field
108(
109 const word& fieldName,
110 const dimensionSet& dim
111) const
112{
113 const fvMesh& mesh = this->internalField().mesh();
114 auto* ptr = mesh.getObjectPtr<volScalarField>(fieldName);
115
116 if (!ptr)
117 {
118 ptr = new volScalarField
119 (
121 (
122 fieldName,
123 mesh.time().timeName(),
124 mesh.thisDb(),
128 ),
129 mesh,
131 );
132
133 ptr->store();
134 }
136 return *ptr;
137}
138
139
140// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
141
143(
144 const fvPatch& p,
146)
147:
148 zeroGradientFvPatchScalarField(p, iF),
149 equilibriumModel_(equilibriumModelType::LANGMUIR),
150 kinematicModel_(kineticModelType::PseudoFirstOrder),
151 thicknessPtr_(nullptr),
152 kabs_(1),
153 kl_(0),
154 max_(1),
155 rhoS_(0),
156 pName_("p"),
157 dfldp_(p.size(), 0),
158 mass_(p.size(), 0)
159{}
160
161
163(
164 const fvPatch& p,
166 const dictionary& dict
167)
168:
169 zeroGradientFvPatchScalarField(p, iF, dict),
170 equilibriumModel_(equilibriumModelTypeNames.get("equilibriumModel", dict)),
171 kinematicModel_(kinematicModelTypeNames.get("kinematicModel", dict)),
172 thicknessPtr_(PatchFunction1<scalar>::New(p.patch(), "thickness", dict)),
173 kabs_(dict.getCheck<scalar>("kabs", scalarMinMax::ge(0))),
174 kl_(dict.getCheck<scalar>("kl", scalarMinMax::ge(0))),
175 max_(dict.getCheck<scalar>("max", scalarMinMax::ge(0))),
176 rhoS_(dict.get<scalar>("rhoS")),
177 pName_(dict.getOrDefault<word>("p", "p")),
178 dfldp_("dfldp", dict, p.size(), IOobjectOption::LAZY_READ),
179 mass_("mass", dict, p.size(), IOobjectOption::LAZY_READ)
180{
181 if (!this->readValueEntry(dict))
182 {
184 }
185}
186
187
189(
190 const speciesSorptionFvPatchScalarField& ptf,
191 const fvPatch& p,
192 const DimensionedField<scalar, volMesh>& iF,
193 const fvPatchFieldMapper& mapper
194)
195:
196 zeroGradientFvPatchScalarField(ptf, p, iF, mapper),
197 equilibriumModel_(ptf.equilibriumModel_),
198 kinematicModel_(ptf.kinematicModel_),
199 thicknessPtr_(ptf.thicknessPtr_.clone(patch().patch())),
200 kabs_(ptf.kabs_),
201 kl_(ptf.kl_),
202 max_(ptf.max_),
203 rhoS_(ptf.rhoS_),
204 pName_(ptf.pName_),
205 dfldp_(ptf.dfldp_, mapper),
206 mass_(ptf.mass_, mapper)
207{}
208
209
211(
213)
214:
215 zeroGradientFvPatchScalarField(ptf),
216 equilibriumModel_(ptf.equilibriumModel_),
217 kinematicModel_(ptf.kinematicModel_),
218 thicknessPtr_(ptf.thicknessPtr_.clone(patch().patch())),
219 kabs_(ptf.kabs_),
220 kl_(ptf.kl_),
221 max_(ptf.max_),
222 rhoS_(ptf.rhoS_),
223 pName_(ptf.pName_),
224 dfldp_(ptf.dfldp_),
225 mass_(ptf.mass_)
226{}
227
228
230(
233)
234:
235 zeroGradientFvPatchScalarField(ptf, iF),
236 equilibriumModel_(ptf.equilibriumModel_),
237 kinematicModel_(ptf.kinematicModel_),
238 thicknessPtr_(ptf.thicknessPtr_.clone(patch().patch())),
239 kabs_(ptf.kabs_),
240 kl_(ptf.kl_),
241 max_(ptf.max_),
242 rhoS_(ptf.rhoS_),
243 pName_(ptf.pName_),
244 dfldp_(ptf.dfldp_),
245 mass_(ptf.mass_)
246{}
247
248
249// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
250
252(
253 const fvPatchFieldMapper& m
254)
255{
256 zeroGradientFvPatchScalarField::autoMap(m);
257
258 dfldp_.autoMap(m);
259 mass_.autoMap(m);
260
261 if (thicknessPtr_)
262 {
263 thicknessPtr_->autoMap(m);
264 }
265}
266
267
269(
270 const fvPatchScalarField& ptf,
271 const labelList& addr
272)
273{
274 zeroGradientFvPatchScalarField::rmap(ptf, addr);
275
277
278 dfldp_.rmap(tiptf.dfldp_, addr);
279 mass_.rmap(tiptf.mass_, addr);
280
281 if (thicknessPtr_)
282 {
283 thicknessPtr_->rmap(tiptf.thicknessPtr_(), addr);
284 }
285}
286
287
289patchSource() const
290{
291 const auto& thermo = db().lookupObject<rhoReactionThermo>
292 (
294 );
295
296 const label speciesId =
297 thermo.composition().species()[this->internalField().name()];
298
299 const scalar Wi(thermo.composition().W(speciesId));
300
301 const scalar t = db().time().timeOutputValue();
302
303 const scalarField h(thicknessPtr_->value(t));
304
305 const scalarField AbyV(this->patch().magSf());
306
307 // Solid mass [kg]
308 const scalarField mass(h*AbyV*rhoS_);
309
310 scalarField Vol(this->patch().size());
311
312 forAll(AbyV, facei)
313 {
314 const label faceCelli = this->patch().faceCells()[facei];
315 Vol[facei] = this->internalField().mesh().V()[faceCelli];
316 }
317
318 // The moles absorbed by the solid
319 // dfldp[mol/kg/sec]* mass[kg]* Wi[kg/mol] / Vol[m3]= [kg/sec/m3]
320 const scalarField dfldp(-dfldp_*mass*Wi*1e-3/Vol);
321
322 if (debug)
323 {
324 auto limits = gMinMax(dfldp);
325
326 Info<< " Patch mass rate min/max [kg/m3/sec]: "
327 << limits.min() << " - " << limits.max() << endl;
328 }
329
330 return tmp<scalarField>::New(dfldp);
331}
332
333
335mass() const
336{
337 return tmp<scalarField>::New(mass_);
338}
339
340
342{
343 if (updated())
344 {
345 return;
346 }
347
348 // equilibrium in mol/kg
349 scalarField cEq(patch().size(), 0);
350
351 switch (equilibriumModel_)
352 {
353 case equilibriumModelType::LANGMUIR:
354 {
355 // mole fraction
356 tmp<scalarField> tco = calcMoleFractions();
357
358 const auto& pp = patch().lookupPatchField<volScalarField>(pName_);
359
360 cEq = max_*(kl_*tco()*pp/(1 + kl_*tco()*pp));
361 break;
362 }
363 default:
364 break;
365 }
366
367 // source [mol/kg/sec]
368 dfldp_ = Zero;
369
370 switch (kinematicModel_)
371 {
372 case kineticModelType::PseudoFirstOrder:
373 {
374 dfldp_ = kabs_*(cEq - mass_);
375 }
376 default:
377 break;
378 }
379
380 // mass [mol/kg]
381 const scalar dt = db().time().deltaTValue();
382 mass_ += dfldp_*dt;
383 mass_ = max(mass_, scalar(0));
384
385 scalarField& pMass =
386 field
387 (
388 "absorbedMass" + this->internalField().name(),
390 ).boundaryFieldRef()[patch().index()];
391
392 pMass = mass_;
393
394 if (debug)
395 {
396 auto limits = gMinMax(dfldp_);
397
398 Info<< " Absorption rate min/max [mol/kg/sec]: "
399 << limits.min() << " - " << limits.max() << endl;
400 }
401
402 zeroGradientFvPatchScalarField::updateCoeffs();
403}
404
405
407{
409
410 os.writeEntry
411 (
412 "equilibriumModel", equilibriumModelTypeNames[equilibriumModel_]
413 );
414 os.writeEntry
415 (
416 "kinematicModel", kinematicModelTypeNames[kinematicModel_]
417 );
418 if (thicknessPtr_)
419 {
420 thicknessPtr_->writeData(os);
421 }
422 os.writeEntry("kabs", kabs_);
423 os.writeEntry("kl", kl_);
424 os.writeEntry("max", max_);
425 os.writeEntry("rhoS", rhoS_);
426
427 dfldp_.writeEntry("dfldp", os);
428 mass_.writeEntry("mass", os);
429 os.writeEntryIfDifferent<word>("p", "p", pName_);
430
433
434
435// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
436
437namespace Foam
438{
440 (
442 speciesSorptionFvPatchScalarField
443 );
444}
445
446// ************************************************************************* //
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...
Definition Enum.H:57
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,...
Definition IOobject.H:191
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition Ostream.H:331
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition Ostream.H:346
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.
Definition Time.C:714
bool get(const label i) const
Definition UList.H:868
label size() const noexcept
Definition UList.H:706
void size(const label n)
Definition UList.H:118
static const word dictName
The dictionary name ("thermophysicalProperties").
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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.
Definition fvMesh.H:85
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition fvMesh.H:376
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
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.
Definition fvPatch.H:71
virtual const labelUList & faceCells() const
Return faceCells.
Definition fvPatch.C:107
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...
Foam::rhoReactionThermo.
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 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.
static const Enum< equilibriumModelType > equilibriumModelTypeNames
Names for equilibriumModelType.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
PtrList< volScalarField > & Y
rDeltaTY field()
auto limits
Definition setRDeltaT.H:186
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
auto & name
label cellId
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
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.
Definition List.H:62
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.
Definition MinMax.H:97
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
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)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
volScalarField & h
psiReactionThermo & thermo
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299