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_
86  (
87  dict.found("dhdt")
88  ? scalarField("dhdt", dict, p.size())
89  : scalarField(p.size(), 0)
90  )
91 {
92  switch (enthalpyModel_)
93  {
94  case enthalpyModelType::calculated:
95  {
96  enthalpyMassLoadPtr_ =
97  Function1<scalar>::New("enthalpyTable", dict);
98  break;
99  }
100  case enthalpyModelType::estimated:
101  {
102  break;
103  }
104  }
105 
106  if (dict.found("value"))
107  {
109  (
110  scalarField("value", dict, p.size())
111  );
112  }
113  else
114  {
116  }
117 }
118 
119 
121 (
122  const enthalpySorptionFvPatchScalarField& ptf,
123  const fvPatch& p,
124  const DimensionedField<scalar, volMesh>& iF,
125  const fvPatchFieldMapper& mapper
126 )
127 :
128  zeroGradientFvPatchScalarField(ptf, p, iF, mapper),
129  enthalpyModel_(ptf.enthalpyModel_),
130  includeHs_(ptf.includeHs_),
131  enthalpyMassLoadPtr_(ptf.enthalpyMassLoadPtr_.clone()),
132  C_(ptf.C_),
133  Hvap_(ptf.Hvap_),
134  speciesName_(ptf.speciesName_),
135  pName_(ptf.pName_),
136  TName_(ptf.TName_),
137  dhdt_(ptf.dhdt_, mapper)
138 {}
139 
140 
142 (
144 )
145 :
146  zeroGradientFvPatchScalarField(ptf),
147  enthalpyModel_(ptf.enthalpyModel_),
148  includeHs_(ptf.includeHs_),
149  enthalpyMassLoadPtr_(ptf.enthalpyMassLoadPtr_.clone()),
150  C_(ptf.C_),
151  Hvap_(ptf.Hvap_),
152  speciesName_(ptf.speciesName_),
153  pName_(ptf.pName_),
154  TName_(ptf.TName_),
155  dhdt_(ptf.dhdt_)
156 {}
157 
158 
160 (
163 )
164 :
165  zeroGradientFvPatchScalarField(ptf, iF),
166  enthalpyModel_(ptf.enthalpyModel_),
167  includeHs_(ptf.includeHs_),
168  enthalpyMassLoadPtr_(ptf.enthalpyMassLoadPtr_.clone()),
169  C_(ptf.C_),
170  Hvap_(ptf.Hvap_),
171  speciesName_(ptf.speciesName_),
172  pName_(ptf.pName_),
173  TName_(ptf.TName_),
174  dhdt_(ptf.dhdt_)
175 {}
176 
177 
178 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
179 
181 (
182  const fvPatchFieldMapper& m
183 )
184 {
185  zeroGradientFvPatchScalarField::autoMap(m);
186 
187  dhdt_.autoMap(m);
188 }
189 
190 
192 (
193  const fvPatchScalarField& ptf,
194  const labelList& addr
195 )
196 {
197  zeroGradientFvPatchScalarField::rmap(ptf, addr);
198 
199  const auto& tiptf = refCast<const enthalpySorptionFvPatchScalarField>(ptf);
200 
201  dhdt_.rmap(tiptf.dhdt_, addr);
202 }
203 
204 
206 patchSource() const
207 {
208  const auto& Yp =
209  refCast<const speciesSorptionFvPatchScalarField>
210  (
211  patch().lookupPatchField<volScalarField, scalar>
212  (
213  speciesName_
214  )
215  );
216 
217  // mass rate [kg/sec/m3]
218  tmp<scalarField> tmassb = Yp.patchSource();
219  const scalarField& massb = tmassb();
220 
221  // The moles absorbed by the solid
222  // dhdt[J/kg] * kg/sec/m3 = [J/m3/s]
223  scalarField dhdt(dhdt_*massb);
224 
225  if (includeHs_)
226  {
227  const fvPatchField<scalar>& pp =
228  patch().lookupPatchField<volScalarField, scalar>(pName_);
229 
230  const fvPatchField<scalar>& Tp =
231  patch().lookupPatchField<volScalarField, scalar>(TName_);
232 
233  const auto& thermo = db().lookupObject<rhoReactionThermo>
234  (
236  );
237 
238  const basicSpecieMixture& composition = thermo.composition();
239 
240  const label speciesId =
241  thermo.composition().species()[speciesName_];
242 
243  scalarField hsp(this->patch().size(), 0);
244 
245  forAll(pp, facei)
246  {
247  hsp[facei] = composition.Hs(speciesId, pp[facei], Tp[facei]);
248  }
249 
250  dhdt += hsp*massb;
251  }
252 
253  if (debug)
254  {
255  Info<< " Patch enthalpy rate min/max [J/m3/sec]: "
256  << gMin(dhdt) << " - " << gMax(dhdt) << endl;
257  }
258 
259  return tmp<scalarField>::New(dhdt);
260 }
261 
262 
264 {
265  if (updated())
266  {
267  return;
268  }
269 
270  const auto& Yp =
271  refCast<const speciesSorptionFvPatchScalarField>
272  (
273  patch().lookupPatchField<volScalarField, scalar>
274  (
275  speciesName_
276  )
277  );
278 
279  switch (enthalpyModel_)
280  {
281  case enthalpyModelType::estimated:
282  {
283  dhdt_ = -C_*Hvap_;
284  break;
285  }
286  case enthalpyModelType::calculated:
287  {
288  // mass [mol/kg]
289  tmp<scalarField> tmassb = Yp.mass();
290  const scalarField& massb = tmassb.ref();
291 
292  forAll(massb, facei)
293  {
294  const scalar mfacei = massb[facei];
295 
296  dhdt_[facei] = enthalpyMassLoadPtr_->value(mfacei);
297  }
298  break;
299  }
300  default:
301  break;
302  }
303 
304  if (debug)
305  {
306  Info<< " Enthalpy change min/max [J/kg]: "
307  << gMin(dhdt_) << " - " << gMax(dhdt_) << endl;
308  }
309 
310  zeroGradientFvPatchScalarField::updateCoeffs();
311 }
312 
313 
315 {
317 
318  os.writeEntry("enthalpyModel", enthalpyModelTypeNames[enthalpyModel_]);
319 
320  if (enthalpyMassLoadPtr_)
321  {
322  enthalpyMassLoadPtr_->writeData(os);
323  }
324 
325  os.writeEntry("species", speciesName_);
326 
327  os.writeEntryIfDifferent<bool>("includeHs", true, includeHs_);
328  os.writeEntryIfDifferent<scalar>("C", scalar(0), C_);
329  os.writeEntryIfDifferent<scalar>("Hvap", scalar(0), Hvap_);
330  os.writeEntryIfDifferent<word>("p", "p", pName_);
331  os.writeEntryIfDifferent<word>("T", "T", TName_);
332 
333  dhdt_.writeEntry("dhdt", os);
334 
335  writeEntry("value", os);
336 }
337 
338 
339 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
340 
341 namespace Foam
342 {
344  (
346  enthalpySorptionFvPatchScalarField
347  );
348 }
349 
350 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
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:120
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.
friend Ostream & operator(Ostream &, const fvPatchField< scalar > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:68
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:312
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:333
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:413
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
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:203
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:327
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...
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:352
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.
bool found
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157