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
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 (
120 IOobject
121 (
122 fieldName,
123 mesh.time().timeName(),
124 mesh,
127 ),
128 mesh,
130 );
131
132 ptr->store();
133 }
134
135 return *ptr;
136}
137
138
139// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
140
142(
143 const fvPatch& p,
145)
146:
147 zeroGradientFvPatchScalarField(p, iF),
148 equilibriumModel_(equilibriumModelType::LANGMUIR),
149 kinematicModel_(kineticModelType::PseudoFirstOrder),
150 thicknessPtr_(nullptr),
151 kabs_(1),
152 kl_(0),
153 max_(1),
154 rhoS_(0),
155 pName_("p"),
156 dfldp_(p.size(), 0),
157 mass_(p.size(), 0)
158{}
159
160
162(
163 const fvPatch& p,
165 const dictionary& dict
166)
167:
168 zeroGradientFvPatchScalarField(p, iF, dict),
169 equilibriumModel_(equilibriumModelTypeNames.get("equilibriumModel", dict)),
170 kinematicModel_(kinematicModelTypeNames.get("kinematicModel", dict)),
171 thicknessPtr_(PatchFunction1<scalar>::New(p.patch(), "thickness", dict)),
172 kabs_(dict.getCheck<scalar>("kabs", scalarMinMax::ge(0))),
173 kl_(dict.getCheck<scalar>("kl", scalarMinMax::ge(0))),
174 max_(dict.getCheck<scalar>("max", scalarMinMax::ge(0))),
175 rhoS_(dict.get<scalar>("rhoS")),
176 pName_(dict.getOrDefault<word>("p", "p")),
177 dfldp_
178 (
179 dict.found("dfldp")
180 ? scalarField("dfldp", dict, p.size())
181 : scalarField(p.size(), 0)
182 ),
183 mass_
184 (
185 dict.found("mass")
186 ? scalarField("mass", dict, p.size())
187 : scalarField(p.size(), 0)
188 )
189{
190 if (dict.found("value"))
191 {
193 (
194 scalarField("value", dict, p.size())
195 );
196 }
197 else
198 {
200 }
201}
202
203
205(
207 const fvPatch& p,
209 const fvPatchFieldMapper& mapper
210)
211:
212 zeroGradientFvPatchScalarField(ptf, p, iF, mapper),
213 equilibriumModel_(ptf.equilibriumModel_),
214 kinematicModel_(ptf.kinematicModel_),
215 thicknessPtr_(ptf.thicknessPtr_.clone(patch().patch())),
216 kabs_(ptf.kabs_),
217 kl_(ptf.kl_),
218 max_(ptf.max_),
219 rhoS_(ptf.rhoS_),
220 pName_(ptf.pName_),
221 dfldp_(ptf.dfldp_, mapper),
222 mass_(ptf.mass_, mapper)
223{}
224
225
227(
229)
230:
231 zeroGradientFvPatchScalarField(ptf),
232 equilibriumModel_(ptf.equilibriumModel_),
233 kinematicModel_(ptf.kinematicModel_),
234 thicknessPtr_(ptf.thicknessPtr_.clone(patch().patch())),
235 kabs_(ptf.kabs_),
236 kl_(ptf.kl_),
237 max_(ptf.max_),
238 rhoS_(ptf.rhoS_),
239 pName_(ptf.pName_),
240 dfldp_(ptf.dfldp_),
241 mass_(ptf.mass_)
242{}
243
244
246(
249)
250:
251 zeroGradientFvPatchScalarField(ptf, iF),
252 equilibriumModel_(ptf.equilibriumModel_),
253 kinematicModel_(ptf.kinematicModel_),
254 thicknessPtr_(ptf.thicknessPtr_.clone(patch().patch())),
255 kabs_(ptf.kabs_),
256 kl_(ptf.kl_),
257 max_(ptf.max_),
258 rhoS_(ptf.rhoS_),
259 pName_(ptf.pName_),
260 dfldp_(ptf.dfldp_),
261 mass_(ptf.mass_)
262{}
263
264
265// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
266
268(
269 const fvPatchFieldMapper& m
270)
271{
272 zeroGradientFvPatchScalarField::autoMap(m);
273
274 dfldp_.autoMap(m);
275 mass_.autoMap(m);
276
277 if (thicknessPtr_)
278 {
279 thicknessPtr_->autoMap(m);
280 }
281}
282
283
285(
286 const fvPatchScalarField& ptf,
287 const labelList& addr
288)
289{
290 zeroGradientFvPatchScalarField::rmap(ptf, addr);
291
292 const auto& tiptf = refCast<const speciesSorptionFvPatchScalarField>(ptf);
293
294 dfldp_.rmap(tiptf.dfldp_, addr);
295 mass_.rmap(tiptf.mass_, addr);
296
297 if (thicknessPtr_)
298 {
299 thicknessPtr_->rmap(tiptf.thicknessPtr_(), addr);
300 }
301}
302
303
305patchSource() const
306{
307 const auto& thermo = db().lookupObject<rhoReactionThermo>
308 (
310 );
311
312 const label speciesId =
313 thermo.composition().species()[this->internalField().name()];
314
315 const scalar Wi(thermo.composition().W(speciesId));
316
317 const scalar t = db().time().timeOutputValue();
318
319 const scalarField h(thicknessPtr_->value(t));
320
321 const scalarField AbyV(this->patch().magSf());
322
323 // Solid mass [kg]
324 const scalarField mass(h*AbyV*rhoS_);
325
326 scalarField Vol(this->patch().size());
327
328 forAll(AbyV, facei)
329 {
330 const label faceCelli = this->patch().faceCells()[facei];
331 Vol[facei] = this->internalField().mesh().V()[faceCelli];
332 }
333
334 // The moles absorbed by the solid
335 // dfldp[mol/kg/sec]* mass[kg]* Wi[kg/mol] / Vol[m3]= [kg/sec/m3]
336 const scalarField dfldp(-dfldp_*mass*Wi*1e-3/Vol);
337
338 if (debug)
339 {
340 Info<< " Patch mass rate min/max [kg/m3/sec]: "
341 << gMin(dfldp) << " - " << gMax(dfldp) << endl;
342 }
343
344 return tmp<scalarField>::New(dfldp);
345}
346
347
349mass() const
350{
351 return tmp<scalarField>::New(mass_);
352}
353
354
356{
357 if (updated())
358 {
359 return;
360 }
361
362 // equilibrium in mol/kg
363 scalarField cEq(patch().size(), 0);
364
365 switch (equilibriumModel_)
366 {
367 case equilibriumModelType::LANGMUIR:
368 {
369 // mole fraction
370 tmp<scalarField> tco = calcMoleFractions();
371
372 const fvPatchField<scalar>& pp =
373 patch().lookupPatchField<volScalarField, scalar>(pName_);
374
375 cEq = max_*(kl_*tco()*pp/(1 + kl_*tco()*pp));
376 break;
377 }
378 default:
379 break;
380 }
381
382 // source [mol/kg/sec]
383 dfldp_ = Zero;
384
385 switch (kinematicModel_)
386 {
387 case kineticModelType::PseudoFirstOrder:
388 {
389 dfldp_ = kabs_*(cEq - mass_);
390 }
391 default:
392 break;
393 }
394
395 // mass [mol/kg]
396 const scalar dt = db().time().deltaTValue();
397 mass_ += dfldp_*dt;
398 mass_ = max(mass_, scalar(0));
399
400 scalarField& pMass =
401 field
402 (
403 "absorbedMass" + this->internalField().name(),
405 ).boundaryFieldRef()[patch().index()];
406
407 pMass = mass_;
408
409 if (debug)
410 {
411 Info<< " Absorption rate min/max [mol/kg/sec]: "
412 << gMin(dfldp_) << " - " << gMax(dfldp_) << endl;
413 }
414
415 zeroGradientFvPatchScalarField::updateCoeffs();
416}
417
418
420{
422
424 (
425 "equilibriumModel", equilibriumModelTypeNames[equilibriumModel_]
426 );
428 (
429 "kinematicModel", kinematicModelTypeNames[kinematicModel_]
430 );
431 if (thicknessPtr_)
432 {
433 thicknessPtr_->writeData(os);
434 }
435 os.writeEntry("kabs", kabs_);
436 os.writeEntry("kl", kl_);
437 os.writeEntry("max", max_);
438 os.writeEntry("rhoS", rhoS_);
439
440 dfldp_.writeEntry("dfldp", os);
441 mass_.writeEntry("mass", os);
442 os.writeEntryIfDifferent<word>("p", "p", pName_);
443
444 writeEntry("value", os);
445}
446
447
448// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
449
450namespace Foam
451{
453 (
456 );
457}
458
459// ************************************************************************* //
bool found
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Field< Type > & field() const
Return field.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
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:251
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
label size() const noexcept
The number of elements in the UList.
Definition: UListI.H:420
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static const word dictName
Definition: basicThermo.H:256
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
virtual bool write()
Write the output fields.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
A FieldMapper for finite-volume patch fields.
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:210
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:408
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:368
friend Ostream & operator(Ostream &, const fvPatchField< scalar > &)
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:362
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:113
const Type & lookupObject(const word &name, const bool recursive=false) const
Type * getObjectPtr(const word &name, const bool recursive=false) const
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 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:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
PtrList< volScalarField > & Y
rDeltaTY field()
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
label cellId
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:55
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
Type gMax(const FieldField< Field, Type > &f)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
volScalarField & h
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333