Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37const Foam::Enum
38<
39 Foam::enthalpySorptionFvPatchScalarField::enthalpyModelType
40>
41Foam::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, &iF.db());
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)
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);
188
189 dhdt_.rmap(tiptf.dhdt_, addr);
190}
191
192
194patchSource() const
195{
196 const auto& Yp =
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 auto limits = gMinMax(dhdt);
238
239 Info<< " Patch enthalpy rate min/max [J/m3/sec]: "
240 << limits.min() << " - " << limits.max() << endl;
241 }
242
243 return tmp<scalarField>::New(dhdt);
244}
245
246
248{
249 if (updated())
250 {
251 return;
252 }
253
254 const auto& Yp =
256 (
257 patch().lookupPatchField<volScalarField>(speciesName_)
258 );
259
260 switch (enthalpyModel_)
261 {
262 case enthalpyModelType::estimated:
263 {
264 dhdt_ = -C_*Hvap_;
265 break;
266 }
267 case enthalpyModelType::calculated:
268 {
269 // mass [mol/kg]
270 tmp<scalarField> tmassb = Yp.mass();
271 const scalarField& massb = tmassb.ref();
272
273 forAll(massb, facei)
274 {
275 const scalar mfacei = massb[facei];
276
277 dhdt_[facei] = enthalpyMassLoadPtr_->value(mfacei);
278 }
279 break;
280 }
281 default:
282 break;
283 }
284
285 if (debug)
286 {
287 auto limits = gMinMax(dhdt_);
288
289 Info<< " Enthalpy change min/max [J/kg]: "
290 << limits.min() << " - " << limits.max() << endl;
291 }
292
293 zeroGradientFvPatchScalarField::updateCoeffs();
294}
295
296
298{
300
301 os.writeEntry("enthalpyModel", enthalpyModelTypeNames[enthalpyModel_]);
302
303 if (enthalpyMassLoadPtr_)
304 {
305 enthalpyMassLoadPtr_->writeData(os);
306 }
307
308 os.writeEntry("species", speciesName_);
309
310 os.writeEntryIfDifferent<bool>("includeHs", true, includeHs_);
311 os.writeEntryIfDifferent<scalar>("C", scalar(0), C_);
312 os.writeEntryIfDifferent<scalar>("Hvap", scalar(0), Hvap_);
313 os.writeEntryIfDifferent<word>("p", "p", pName_);
314 os.writeEntryIfDifferent<word>("T", "T", TName_);
315
316 dhdt_.writeEntry("dhdt", os);
317
320
321
322// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
323
324namespace Foam
325{
327 (
329 enthalpySorptionFvPatchScalarField
330 );
331}
332
333// ************************************************************************* //
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...
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
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
bool get(const label i) const
Definition UList.H:868
void size(const label n)
Definition UList.H:118
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
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
This is a temperature boundary condition which works in conjunction with the speciesSorption conditio...
enthalpySorptionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual tmp< fvPatchField< scalar > > clone() const
Return a clone.
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.
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.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
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
Foam::rhoReactionThermo.
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
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A class for handling words, derived from Foam::string.
Definition word.H:66
basicSpecieMixture & composition
volScalarField & p
auto limits
Definition setRDeltaT.H:186
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
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
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
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
fvPatchField< scalar > fvPatchScalarField
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299