Loading...
Searching...
No Matches
enthalpySorptionFvPatchScalarField.H
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
26Class
27 Foam::enthalpySorptionFvPatchScalarField
28
29Group
30 grpGenericBoundaryConditions
31
32Description
33 This is a temperature boundary condition which works
34 in conjunction with the \c speciesSorption condition for species.
35
36 This boundary condition substracts or adds enthalpy associated with the
37 adsorption provided by the \c speciesSorption condition.
38
39 It can handle two enthalpy models:
40
41 1) Estimate
42 2) Calculated
43
44 On top of this, the enthalpy associated with the sensible enthalpy
45 corresponding with the species transport can be added using \c includeHs.
46
47Usage
48 Example of the boundary condition specification:
49 \verbatim
50 <patchName>
51 {
52 // Mandatory entries
53 type enthalpySorption;
54 enthalpyModel <word>;
55 species <word>;
56
57 // Conditional mandatory entries
58
59 // when enthalpyModel == calculated
60 enthalpyTable <Function1<scalar>>
61
62 // enthalpyTable
63 // {
64 // type table;
65 // values ((0 0)(1 50));
66 // }
67
68 // Optional entries
69 includeHs <bool>;
70 C <scalar>;
71 Hvap <scalar>;
72 dhdt <scalarField>;
73 p <word>;
74 T <word>;
75 }
76 \endverbatim
77
78 where the entries mean:
79 \table
80 Property | Description | Type | Reqd | Deflt
81 type | Type name: enthalpySorption | word | yes | -
82 enthalpyModel | Adsorption enthalpy model | word | yes | -
83 species | Name of associated species | word | yes | -
84 enthalpyTable | Calculated enthalpy model table <!--
85 --> | Function1<scalar> | no | -
86 includeHs | Include sensible enthalpy | bool | no | true
87 C | Estimate enthalpy model constant | scalar | no | 0
88 Hvap | Evaporation heat for species | scalar | no | 0
89 p | Name of operand pressure field | word | no | p
90 T | Name of operand temperature field | word | no | T
91 dhdt | Enthalpy change on cells next to patch [J/kg] <!--
92 --> | scalarField | no | Zero
93 \endtable
94
95 Options for the \c enthalpyModel entry:
96 \verbatim
97 estimated | Enthalpy is estimated
98 calculated | Enthalpy is calculated based on enthalpyTable
99 \endverbatim
100
101 The inherited entries are elaborated in:
102 - \link zeroGradientFvPatchFields.H \endlink
103 - \link Function1.H \endlink
104
105SourceFiles
106 enthalpySorptionFvPatchScalarField.C
107
108\*---------------------------------------------------------------------------*/
109
110#ifndef Foam_enthalpySorptionFvPatchScalarField_H
111#define Foam_enthalpySorptionFvPatchScalarField_H
112
113#include "boundarySourcePatch.H"
115#include "Function1.H"
116
117// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
118
119namespace Foam
120{
121
122/*---------------------------------------------------------------------------*\
123 Class enthalpySorptionFvPatchScalarField Declaration
124\*---------------------------------------------------------------------------*/
125
127:
128 public zeroGradientFvPatchField<scalar>,
130{
131 // Private Enumeration
132
133 //- Options for the enthalpy model
134 enum enthalpyModelType : char
135 {
136 estimated,
137 calculated
138 };
139
140 //- Names for enthalpyModelType
141 static const Enum<enthalpyModelType> enthalpyModelTypeNames;
142
143
144 // Private Data
145
146 //- Enthalpy model
147 enum enthalpyModelType enthalpyModel_;
148
149 //- Include sensible enthalpy of the species
150 bool includeHs_;
151
152 //- Load enthalpy table for calculated model
153 autoPtr<Function1<scalar>> enthalpyMassLoadPtr_;
154
155 //- Estimated enthalpy model constant
156 scalar C_;
157
158 //- Heat of evaporation of species
159 scalar Hvap_;
160
161 //- Name of operand species field
162 word speciesName_;
163
164 //- Name of operand pressure field
165 word pName_;
166
167 //- Name of operand temperature field
168 word TName_;
169
170 //- Enthalpy change on cells next to patch [J/kg]
171 scalarField dhdt_;
172
173
174public:
175
176 //- Runtime type information
177 TypeName("enthalpySorption");
178
179
180 // Constructors
181
182 //- Construct from patch and internal field
184 (
185 const fvPatch&,
186 const DimensionedField<scalar, volMesh>&
187 );
188
189 //- Construct from patch, internal field and dictionary
192 const fvPatch&,
194 const dictionary&
195 );
196
197 //- Construct by mapping given
198 //- enthalpySorptionFvPatchScalarField onto a new patch
200 (
202 const fvPatch&,
204 const fvPatchFieldMapper&
205 );
206
207 //- Construct as copy
209 (
211 );
212
213 //- Construct as copy setting internal field reference
215 (
218 );
219
220 //- Return a clone
221 virtual tmp<fvPatchField<scalar>> clone() const
222 {
223 return fvPatchField<scalar>::Clone(*this);
224 }
225
226 //- Clone with an internal field reference
228 (
230 ) const
231 {
232 return fvPatchField<scalar>::Clone(*this, iF);
233 }
234
235
236 // Member Functions
237
238 // Mapping
239
240 //- Map (and resize as needed) from self given a mapping object
241 virtual void autoMap
242 (
243 const fvPatchFieldMapper&
244 );
245
246 //- Reverse map the given fvPatchField onto this fvPatchField
247 virtual void rmap
248 (
249 const fvPatchScalarField&,
250 const labelList&
251 );
252
253
254 // Evaluation
255
256 //- Source of cells next to the patch
257 virtual tmp<scalarField> patchSource() const;
258
259 //- Update the coefficients associated with the patch field
260 virtual void updateCoeffs();
261
262
263 // I-O
264
265 //- Write
266 virtual void write(Ostream&) const;
267};
268
269
270// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
271
272} // End namespace Foam
273
274// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
275
276#endif
277
278// ************************************************************************* //
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
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Pure virtual class for sources on cells next to patches.
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.
TypeName("enthalpySorption")
Runtime type information.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual tmp< fvPatchField< scalar > > clone(const DimensionedField< scalar, volMesh > &iF) const
Clone with an internal field reference.
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.
A FieldMapper for finite-volume patch fields.
static tmp< fvPatchField< Type > > Clone(const DerivedPatchField &pf, Args &&... args)
Clone a patch field, optionally with internal field reference etc.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
zeroGradientFvPatchField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
fvPatchField< scalar > fvPatchScalarField
runTime write()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68