liquidFilmBase.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) 2020-2021 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::regionModels::liquidFilmBase
28
29Description
30 Base class for thermal 2D shells
31
32SourceFiles
33 liquidFilmBase.C
34
35\*---------------------------------------------------------------------------*/
36
37#ifndef liquidFilmBase_H
38#define liquidFilmBase_H
39
41#include "autoPtr.H"
42#include "faCFD.H"
43#include "volFieldsFwd.H"
45#include "regionFaModel.H"
46#include "faOptions.H"
47
48// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49
50namespace Foam
51{
52namespace regionModels
53{
54namespace areaSurfaceFilmModels
55{
56
57/*---------------------------------------------------------------------------*\
58 Class liquidFilmBase Declaration
59\*---------------------------------------------------------------------------*/
62:
63 public regionFaModel
64{
65protected:
66
67 // Protected Data
68
69 // Solution parameters
70
71 //- Momentum predictor
73
74 //- Number of outer correctors
75 label nOuterCorr_;
76
77 //- Number of PISO-like correctors
78 label nCorr_;
79
80 //- Number of film thickness correctors
81 label nFilmCorr_;
82
83 //- Cumulative continuity error
84 scalar cumulativeContErr_;
85
86 //- Smallest numerical thickness
88
89 //- Delta wet for sub-models
91
92
93 //- Name of the velocity field
95
96 //- Name of the pressure field
98
99 //- Reference absolute pressure
100 scalar pRef_;
101
102
103 // Fields
104
105 //- Film hight
107
108 //- Film velocity
110
111 //- Film pressure
113
114 //- Primary region pressure
116
117 //- Film momentum flux
119
120 //- Film height flux
122
123 //- Normal gravity field
125
126 //- Gravity
128
129
130 // Mass exchanage fields from the primary region (lagragian)
131
132 //- Mass
134
135 //- Momentum
137
138 //- Normal pressure by particles
140
141 //- Energy
143
144 //- Total mass added
145 scalar addedMassTotal_;
146
147
148 //- faOptions
150
151
152public:
153
154 //- Runtime type information
155 TypeName("liquidFilmBase");
156
157
158 // Declare runtime constructor selection tables
161 (
162 autoPtr,
165 (
166 const word& modelType,
167 const fvPatch& patch,
168 const dictionary& dict
169 ),
170 (modelType, patch, dict)
171 );
172
173
174 // Constructors
175
176
177 //- Construct from type name and mesh and dict
179 (
180 const word& modelType,
181 const fvPatch& patch,
182 const dictionary& dict
183 );
184
185 //- No copy construct
186 liquidFilmBase(const liquidFilmBase&) = delete;
187
188 //- No copy assignment
189 void operator=(const liquidFilmBase&) = delete;
190
191
192 // Selectors
193
194 //- Return a reference to the selected model using dictionary
196 (
197 const fvPatch& patch,
198 const dictionary& dict
199 );
200
201
202 //- Destructor
203 virtual ~liquidFilmBase();
204
205
206 // Member Functions
207
208 //- Courant number evaluation
209 virtual scalar CourantNumber() const;
210
211
212 // Helper functions
213
214 //- Wall velocity
215 tmp<areaVectorField> Uw() const;
216
217 //- Film surface film velocity
218 tmp<areaVectorField> Us() const;
219
220 //- Primary region velocity at film hight. Assume the film to be
221 // in the viscous sub-layer
222 tmp<areaVectorField> Up() const;
223
224 //- Map primary static pressure
225 tmp<areaScalarField> pg() const;
226
227 //- Wet indicator using h0
229
230
231 // Access functions
232
233 //- Return faOptions
235
236 //- Access const reference Uf
237 const areaVectorField& Uf() const;
238
239 //- Access const reference gn
240 const areaScalarField& gn() const;
241
242 //- Gravity
243 const uniformDimensionedVectorField& g() const;
244
245 //- Access const reference h
246 const areaScalarField& h() const;
247
248 //- Access to momentum flux
249 const edgeScalarField& phif() const;
250
251 //- Access continuity flux
252 const edgeScalarField& phi2s() const;
253
254
255 //- Return h0
256 const dimensionedScalar& h0() const;
257
258 //- Access to this region
259 const regionFaModel& region() const;
260
261 //- Access to pRef
262 scalar pRef();
263
264 //- Name of the U field
265 word UName() const;
266
267
268 // Transfer fields - to the primary region (lagragian injection)
269
270 //- Return mass transfer source - Eulerian phase only
271 //virtual tmp<volScalarField> primaryMassTrans() const = 0;
272
273 //- Return the film mass available for transfer to cloud
274 virtual const volScalarField& cloudMassTrans() const = 0;
275
276 //- Return the parcel diameters originating from film to cloud
277 virtual const volScalarField& cloudDiameterTrans() const = 0;
278
279
280
281 // Evolution
282
283 //- Pre-evolve film
284 virtual void preEvolveRegion();
285
286 //- Post-evolve film
287 virtual void postEvolveRegion();
288
289
290 // Thermo variables
291
292 //- Access const reference mu
293 virtual const areaScalarField& mu() const = 0;
294
295 //- Access const reference rho
296 virtual const areaScalarField& rho() const = 0;
297
298 //- Access const reference sigma
299 virtual const areaScalarField& sigma() const = 0;
300
301 //- Access const reference Tf
302 virtual const areaScalarField& Tf() const = 0;
303
304 //- Access const reference Cp
305 virtual const areaScalarField& Cp() const = 0;
306
307
308 // External hook to add sources and mass exchange variables
309
310
311 //- Add sources
312 virtual void addSources
313 (
314 const label patchi, // patchi on primary region
315 const label facei, // facei of patchi
316 const scalar massSource, // [kg]
317 const vector& momentumSource, // [kg.m/s] (tang'l momentum)
318 const scalar pressureSource, // [kg.m/s] (normal momentum)
319 const scalar energySource = 0 // [J]
320 );
321};
322
323
324// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
325
326} // End namespace areaSurfaceFilmModels
327} // End namespace regionModels
328} // End namespace Foam
329
330// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
331
332#endif
333
334// ************************************************************************* //
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:78
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Finite-area options.
Definition: faOptions.H:58
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
const edgeScalarField & phif() const
Access to momentum flux.
scalar cumulativeContErr_
Cumulative continuity error.
areaScalarField ppf_
Primary region pressure.
virtual const areaScalarField & Cp() const =0
Access const reference Cp.
virtual void addSources(const label patchi, const label facei, const scalar massSource, const vector &momentumSource, const scalar pressureSource, const scalar energySource=0)
Add sources.
virtual const areaScalarField & sigma() const =0
Access const reference sigma.
tmp< areaScalarField > alpha() const
Wet indicator using h0.
const regionFaModel & region() const
Access to this region.
label nFilmCorr_
Number of film thickness correctors.
Foam::fa::options & faOptions()
Return faOptions.
tmp< areaScalarField > pg() const
Map primary static pressure.
volScalarField pnSource_
Normal pressure by particles.
dimensionedScalar h0_
Smallest numerical thickness.
void operator=(const liquidFilmBase &)=delete
No copy assignment.
tmp< areaVectorField > Uw() const
Wall velocity.
const areaVectorField & Uf() const
Access const reference Uf.
virtual const areaScalarField & Tf() const =0
Access const reference Tf.
const areaScalarField & h() const
Access const reference h.
tmp< areaVectorField > Up() const
Primary region velocity at film hight. Assume the film to be.
virtual const areaScalarField & mu() const =0
Access const reference mu.
virtual const areaScalarField & rho() const =0
Access const reference rho.
const edgeScalarField & phi2s() const
Access continuity flux.
TypeName("liquidFilmBase")
Runtime type information.
declareRunTimeSelectionTable(autoPtr, liquidFilmBase, dictionary,(const word &modelType, const fvPatch &patch, const dictionary &dict),(modelType, patch, dict))
const dimensionedScalar & h0() const
Return h0.
static autoPtr< liquidFilmBase > New(const fvPatch &patch, const dictionary &dict)
Return a reference to the selected model using dictionary.
dimensionedScalar deltaWet_
Delta wet for sub-models.
virtual scalar CourantNumber() const
Courant number evaluation.
tmp< areaVectorField > Us() const
Film surface film velocity.
virtual const volScalarField & cloudDiameterTrans() const =0
Return the parcel diameters originating from film to cloud.
const uniformDimensionedVectorField & g() const
Gravity.
virtual const volScalarField & cloudMassTrans() const =0
Return mass transfer source - Eulerian phase only.
const areaScalarField & gn() const
Access const reference gn.
liquidFilmBase(const liquidFilmBase &)=delete
No copy construct.
Base class for area region models.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
Namespace for OpenFOAM.
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes)
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73