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 -------------------------------------------------------------------------------
10 License
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 
36 const Foam::Enum
37 <
39 >
41 ({
42  { equilibriumModelType::LANGMUIR, "Langmuir" }
43 });
44 
45 
46 const Foam::Enum
47 <
49 >
51 ({
52  { kineticModelType::PseudoFirstOrder, "PseudoFirstOrder" }
53 });
54 
55 
56 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
57 
59 Foam::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 
107 Foam::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  (
120  IOobject
121  (
122  fieldName,
123  mesh.time().timeName(),
124  mesh.thisDb(),
128  ),
129  mesh,
130  dimensionedScalar(dim, Zero)
131  );
132 
133  ptr->store();
134  }
135 
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 
276  const auto& tiptf = refCast<const speciesSorptionFvPatchScalarField>(ptf);
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 
289 patchSource() 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  Info<< " Patch mass rate min/max [kg/m3/sec]: "
325  << gMin(dfldp) << " - " << gMax(dfldp) << endl;
326  }
327 
328  return tmp<scalarField>::New(dfldp);
329 }
330 
331 
333 mass() const
334 {
335  return tmp<scalarField>::New(mass_);
336 }
337 
338 
340 {
341  if (updated())
342  {
343  return;
344  }
345 
346  // equilibrium in mol/kg
347  scalarField cEq(patch().size(), 0);
348 
349  switch (equilibriumModel_)
350  {
351  case equilibriumModelType::LANGMUIR:
352  {
353  // mole fraction
354  tmp<scalarField> tco = calcMoleFractions();
355 
356  const auto& pp = patch().lookupPatchField<volScalarField>(pName_);
357 
358  cEq = max_*(kl_*tco()*pp/(1 + kl_*tco()*pp));
359  break;
360  }
361  default:
362  break;
363  }
364 
365  // source [mol/kg/sec]
366  dfldp_ = Zero;
367 
368  switch (kinematicModel_)
369  {
370  case kineticModelType::PseudoFirstOrder:
371  {
372  dfldp_ = kabs_*(cEq - mass_);
373  }
374  default:
375  break;
376  }
377 
378  // mass [mol/kg]
379  const scalar dt = db().time().deltaTValue();
380  mass_ += dfldp_*dt;
381  mass_ = max(mass_, scalar(0));
382 
383  scalarField& pMass =
384  field
385  (
386  "absorbedMass" + this->internalField().name(),
387  dimensionSet(dimMoles/dimMass)
388  ).boundaryFieldRef()[patch().index()];
389 
390  pMass = mass_;
391 
392  if (debug)
393  {
394  Info<< " Absorption rate min/max [mol/kg/sec]: "
395  << gMin(dfldp_) << " - " << gMax(dfldp_) << endl;
396  }
397 
398  zeroGradientFvPatchScalarField::updateCoeffs();
399 }
400 
401 
403 {
405 
406  os.writeEntry
407  (
408  "equilibriumModel", equilibriumModelTypeNames[equilibriumModel_]
409  );
410  os.writeEntry
411  (
412  "kinematicModel", kinematicModelTypeNames[kinematicModel_]
413  );
414  if (thicknessPtr_)
415  {
416  thicknessPtr_->writeData(os);
417  }
418  os.writeEntry("kabs", kabs_);
419  os.writeEntry("kl", kl_);
420  os.writeEntry("max", max_);
421  os.writeEntry("rhoS", rhoS_);
422 
423  dfldp_.writeEntry("dfldp", os);
424  mass_.writeEntry("mass", os);
425  os.writeEntryIfDifferent<word>("p", "p", pName_);
426 
428 }
429 
430 
431 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
432 
433 namespace Foam
434 {
436  (
438  speciesSorptionFvPatchScalarField
439  );
440 }
441 
442 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
Definition: fvPatchField.C:32
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
dictionary dict
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...
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
const objectRegistry & db() const
The associated objectRegistry.
rDeltaTY field()
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:269
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Type gMin(const FieldField< Field, Type > &f)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
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.
const speciesTable & species() const
Return the table of species.
bool store()
Register object with its registry and transfer ownership to the registry.
Definition: regIOobjectI.H:36
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
speciesSorptionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
This boundary condition provides a first-order zero-gradient condition for a given scalar field to mo...
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
psiReactionThermo & thermo
Definition: createFields.H:28
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
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...
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:107
A FieldMapper for finite-volume patch fields.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
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:336
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
label size() const noexcept
The number of elements in the container.
Definition: UList.H:671
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:53
static const Enum< equilibriumModelType > equilibriumModelTypeNames
Names for equilibriumModelType.
tmp< scalarField > mass() const
Access to mass.
const dimensionedScalar h
Planck constant.
static const Enum< kineticModelType > kinematicModelTypeNames
Names for kineticModelType.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< volScalarField > & Y
label cellId
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:391
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
Nothing to be read.
Automatically write from objectRegistry::writeObject()
const DimensionedField< scalar, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
Definition: fvPatchField.H:662
const std::string patch
OpenFOAM patch number as a std::string.
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual basicSpecieMixture & composition()=0
Return the composition of the multi-component mixture.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual scalar W(const label speciei) const =0
Molecular weight of the given specie [kg/kmol].
virtual tmp< scalarField > patchSource() const
Source of cells next to the patch.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Request registration (bool: true)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127