enthalpySorptionFvPatchScalarField.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"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 const Foam::Enum
38 <
39  Foam::enthalpySorptionFvPatchScalarField::enthalpyModelType
40 >
41 Foam::enthalpySorptionFvPatchScalarField::enthalpyModelTypeNames
42 ({
43  { enthalpyModelType::estimated, "estimated" },
44  { enthalpyModelType::calculated, "calculated" }
45 });
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const fvPatch& p,
53  const DimensionedField<scalar, volMesh>& iF
54 )
55 :
56  zeroGradientFvPatchScalarField(p, iF),
57  enthalpyModel_(enthalpyModelType::estimated),
58  includeHs_(false),
59  enthalpyMassLoadPtr_(nullptr),
60  C_(0),
61  Hvap_(0),
62  speciesName_("none"),
63  pName_("p"),
64  TName_("T"),
65  dhdt_(p.size(), 0)
66 {}
67 
68 
70 (
71  const fvPatch& p,
73  const dictionary& dict
74 )
75 :
76  zeroGradientFvPatchScalarField(p, iF, dict),
77  enthalpyModel_(enthalpyModelTypeNames.get("enthalpyModel", dict)),
78  includeHs_(dict.getOrDefault<bool>("includeHs", true)),
79  enthalpyMassLoadPtr_(nullptr),
80  C_(dict.getCheckOrDefault<scalar>("C", 0, scalarMinMax::ge(0))),
81  Hvap_(dict.getCheckOrDefault<scalar>("Hvap", 0, scalarMinMax::ge(0))),
82  speciesName_(dict.get<word>("species")),
83  pName_(dict.getOrDefault<word>("p", "p")),
84  TName_(dict.getOrDefault<word>("T", "T")),
85  dhdt_("dhdt", dict, p.size(), IOobjectOption::LAZY_READ)
86 {
87  switch (enthalpyModel_)
88  {
89  case enthalpyModelType::calculated:
90  {
91  enthalpyMassLoadPtr_ =
92  Function1<scalar>::New("enthalpyTable", dict);
93  break;
94  }
95  case enthalpyModelType::estimated:
96  {
97  break;
98  }
99  }
100 
101  if (!this->readValueEntry(dict))
102  {
104  }
105 }
106 
107 
109 (
110  const enthalpySorptionFvPatchScalarField& ptf,
111  const fvPatch& p,
112  const DimensionedField<scalar, volMesh>& iF,
113  const fvPatchFieldMapper& mapper
114 )
115 :
116  zeroGradientFvPatchScalarField(ptf, p, iF, mapper),
117  enthalpyModel_(ptf.enthalpyModel_),
118  includeHs_(ptf.includeHs_),
119  enthalpyMassLoadPtr_(ptf.enthalpyMassLoadPtr_.clone()),
120  C_(ptf.C_),
121  Hvap_(ptf.Hvap_),
122  speciesName_(ptf.speciesName_),
123  pName_(ptf.pName_),
124  TName_(ptf.TName_),
125  dhdt_(ptf.dhdt_, mapper)
126 {}
127 
128 
130 (
132 )
133 :
134  zeroGradientFvPatchScalarField(ptf),
135  enthalpyModel_(ptf.enthalpyModel_),
136  includeHs_(ptf.includeHs_),
137  enthalpyMassLoadPtr_(ptf.enthalpyMassLoadPtr_.clone()),
138  C_(ptf.C_),
139  Hvap_(ptf.Hvap_),
140  speciesName_(ptf.speciesName_),
141  pName_(ptf.pName_),
142  TName_(ptf.TName_),
143  dhdt_(ptf.dhdt_)
144 {}
145 
146 
148 (
151 )
152 :
153  zeroGradientFvPatchScalarField(ptf, iF),
154  enthalpyModel_(ptf.enthalpyModel_),
155  includeHs_(ptf.includeHs_),
156  enthalpyMassLoadPtr_(ptf.enthalpyMassLoadPtr_.clone()),
157  C_(ptf.C_),
158  Hvap_(ptf.Hvap_),
159  speciesName_(ptf.speciesName_),
160  pName_(ptf.pName_),
161  TName_(ptf.TName_),
162  dhdt_(ptf.dhdt_)
163 {}
164 
165 
166 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
167 
169 (
170  const fvPatchFieldMapper& m
171 )
172 {
173  zeroGradientFvPatchScalarField::autoMap(m);
174 
175  dhdt_.autoMap(m);
176 }
177 
178 
180 (
181  const fvPatchScalarField& ptf,
182  const labelList& addr
183 )
184 {
185  zeroGradientFvPatchScalarField::rmap(ptf, addr);
186 
187  const auto& tiptf = refCast<const enthalpySorptionFvPatchScalarField>(ptf);
188 
189  dhdt_.rmap(tiptf.dhdt_, addr);
190 }
191 
192 
194 patchSource() const
195 {
196  const auto& Yp =
197  refCast<const speciesSorptionFvPatchScalarField>
198  (
199  patch().lookupPatchField<volScalarField>(speciesName_)
200  );
201 
202  // mass rate [kg/sec/m3]
203  tmp<scalarField> tmassb = Yp.patchSource();
204  const scalarField& massb = tmassb();
205 
206  // The moles absorbed by the solid
207  // dhdt[J/kg] * kg/sec/m3 = [J/m3/s]
208  scalarField dhdt(dhdt_*massb);
209 
210  if (includeHs_)
211  {
212  const auto& pp = patch().lookupPatchField<volScalarField>(pName_);
213  const auto& Tp = patch().lookupPatchField<volScalarField>(TName_);
214 
215  const auto& thermo = db().lookupObject<rhoReactionThermo>
216  (
218  );
219 
220  const basicSpecieMixture& composition = thermo.composition();
221 
222  const label speciesId =
223  thermo.composition().species()[speciesName_];
224 
225  scalarField hsp(this->patch().size(), 0);
226 
227  forAll(pp, facei)
228  {
229  hsp[facei] = composition.Hs(speciesId, pp[facei], Tp[facei]);
230  }
231 
232  dhdt += hsp*massb;
233  }
234 
235  if (debug)
236  {
237  Info<< " Patch enthalpy rate min/max [J/m3/sec]: "
238  << gMin(dhdt) << " - " << gMax(dhdt) << endl;
239  }
240 
241  return tmp<scalarField>::New(dhdt);
242 }
243 
244 
246 {
247  if (updated())
248  {
249  return;
250  }
251 
252  const auto& Yp =
253  refCast<const speciesSorptionFvPatchScalarField>
254  (
255  patch().lookupPatchField<volScalarField>(speciesName_)
256  );
257 
258  switch (enthalpyModel_)
259  {
260  case enthalpyModelType::estimated:
261  {
262  dhdt_ = -C_*Hvap_;
263  break;
264  }
265  case enthalpyModelType::calculated:
266  {
267  // mass [mol/kg]
268  tmp<scalarField> tmassb = Yp.mass();
269  const scalarField& massb = tmassb.ref();
270 
271  forAll(massb, facei)
272  {
273  const scalar mfacei = massb[facei];
274 
275  dhdt_[facei] = enthalpyMassLoadPtr_->value(mfacei);
276  }
277  break;
278  }
279  default:
280  break;
281  }
282 
283  if (debug)
284  {
285  Info<< " Enthalpy change min/max [J/kg]: "
286  << gMin(dhdt_) << " - " << gMax(dhdt_) << endl;
287  }
288 
289  zeroGradientFvPatchScalarField::updateCoeffs();
290 }
291 
292 
294 {
296 
297  os.writeEntry("enthalpyModel", enthalpyModelTypeNames[enthalpyModel_]);
298 
299  if (enthalpyMassLoadPtr_)
300  {
301  enthalpyMassLoadPtr_->writeData(os);
302  }
303 
304  os.writeEntry("species", speciesName_);
305 
306  os.writeEntryIfDifferent<bool>("includeHs", true, includeHs_);
307  os.writeEntryIfDifferent<scalar>("C", scalar(0), C_);
308  os.writeEntryIfDifferent<scalar>("Hvap", scalar(0), Hvap_);
309  os.writeEntryIfDifferent<word>("p", "p", pName_);
310  os.writeEntryIfDifferent<word>("T", "T", TName_);
311 
312  dhdt_.writeEntry("dhdt", os);
313 
315 }
316 
317 
318 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319 
320 namespace Foam
321 {
323  (
325  enthalpySorptionFvPatchScalarField
326  );
327 }
328 
329 // ************************************************************************* //
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
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
dictionary dict
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
This is a temperature boundary condition which works in conjunction with the speciesSorption conditio...
Type gMin(const FieldField< Field, Type > &f)
basicSpecieMixture & composition
enthalpySorptionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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
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
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Foam::rhoReactionThermo.
A FieldMapper for finite-volume patch fields.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
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
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
virtual scalar Hs(const label speciei, const scalar p, const scalar T) const =0
Sensible enthalpy [J/kg].
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
const std::string patch
OpenFOAM patch number as a std::string.
messageStream Info
Information stream (stdout output on master, null elsewhere)
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual tmp< scalarField > patchSource() const
Source of cells next to the patch.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127