Loading...
Searching...
No Matches
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 liquid-film models.
31
32Usage
33 Example of the boundary condition specification:
34 \verbatim
35 <patchName>
36 {
37 // Mandatory entries
38 U <word>;
39 pRef <scalar>;
40
41 // Optional entries
42 h0 <scalar>;
43 deltaWet <scalar>;
44 p <word>;
45
46 // Inherited entries
47 ...
48 momentumPredictor <int>;
49 nOuterCorr <label>;
50 nCorr <label>;
51 nFilmCorr <label>;
52 }
53 \endverbatim
54
55 where the entries mean:
56 \table
57 Property | Description | Type | Reqd | Deflt
58 U | Name of velocity field | word | yes | -
59 pRef | Reference absolute pressure | scalar | yes | -
60 h0 | Smallest film thickness | scalar | no | 1e-7
61 deltaWet | Film thickness beyond which face is assumed to be wet <!--
62 --> | scalar | no | 1e-4
63 p | Name of pressure field | word | no | null
64 \endtable
65
66 The inherited entries are elaborated in:
67 - \link regionFaModel.H \endlink
68
69SourceFiles
70 liquidFilmBase.C
71
72\*---------------------------------------------------------------------------*/
73
74#ifndef Foam_liquidFilmBase_H
75#define Foam_liquidFilmBase_H
76
77#include "faCFD.H"
78#include "volFieldsFwd.H"
80#include "regionFaModel.H"
81#include "faOptions.H"
82
83// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84
85namespace Foam
86{
87namespace regionModels
88{
89namespace areaSurfaceFilmModels
90{
91
92/*---------------------------------------------------------------------------*\
93 Class liquidFilmBase Declaration
94\*---------------------------------------------------------------------------*/
95
97:
98 public regionFaModel
99{
100protected:
101
102 // Protected Data
103
104 // Solution parameters
105
106 //- Flag to enable momentum predictor
108
109 //- Number of outer correctors
110 label nOuterCorr_;
111
112 //- Number of PISO-like inner correctors
113 label nCorr_;
114
115 //- Number of film thickness correctors
116 label nFilmCorr_;
117
118 //- Cumulative continuity error
119 scalar cumulativeContErr_;
120
121 //- Smallest film thickness
123
124 //- Film thickness beyond which face is assumed to be wet
126
127
128 //- Name of the velocity field
130
131 //- Name of the pressure field
132 word pName_;
133
134 //- Reference absolute pressure
135 scalar pRef_;
136
137
138 // Fields
139
140 //- Film height
143 //- Film velocity
145
146 //- Film pressure
148
149 //- Primary region pressure
151
152 //- Film momentum flux
154
155 //- Film height flux
158 //- Normal gravity field
160
161 //- Gravity
163
164
165 // Mass exchanage fields from the primary region (lagragian)
166
167 //- Mass
169
170 //- Momentum
173 //- Normal pressure by particles
175
176 //- Total mass added
177 scalar addedMassTotal_;
179
180 //- faOptions
182
184public:
185
186 //- Runtime type information
187 TypeName("liquidFilmBase");
189
190 // Declare runtime constructor selection tables
191
193 (
194 autoPtr,
197 (
198 const word& modelType,
199 const fvMesh& mesh,
200 const dictionary& dict
201 ),
202 (modelType, mesh, dict)
203 );
204
205
206 // Constructors
207
208 //- Construct from type name and mesh and dict
210 (
211 const word& modelType,
212 const fvMesh& mesh,
213 const dictionary& dict
214 );
215
216 //- No copy construct
217 liquidFilmBase(const liquidFilmBase&) = delete;
218
219 //- No copy assignment
220 void operator=(const liquidFilmBase&) = delete;
222
223 // Selectors
224
225 //- Return a reference to the selected model using dictionary
227 (
228 const fvMesh& mesh,
229 const dictionary& dict
230 );
232
233 //- Destructor
234 virtual ~liquidFilmBase();
235
236
237 // Member Functions
238
239 //- Courant number evaluation
240 virtual scalar CourantNumber() const;
241
242
243 // Helper functions
245 //- Wall velocity
246 tmp<areaVectorField> Uw() const;
247
248 //- Film surface film velocity
250
251 //- Primary region velocity at film hight. Assume the film to be
252 // in the viscous sub-layer
253 tmp<areaVectorField> Up() const;
255 //- Map primary static pressure
256 tmp<areaScalarField> pg() const;
257
258 //- Wet indicator using h0
261
262 // Access
263
264 //- Return faOptions
266
267 //- Access const reference Uf
268 const areaVectorField& Uf() const noexcept { return Uf_; }
269
270 //- Access const reference gn
271 const areaScalarField& gn() const noexcept { return gn_; }
272
273 //- Gravity
274 const uniformDimensionedVectorField& g() const noexcept { return g_; }
275
276 //- Access const reference h
277 const areaScalarField& h() const noexcept { return h_; }
278
279 //- Access to momentum flux
280 const edgeScalarField& phif() const noexcept { return phif_; }
281
282 //- Access continuity flux
283 const edgeScalarField& phi2s() const noexcept { return phi2s_; }
284
285 //- Return h0
286 const dimensionedScalar& h0() const noexcept { return h0_; }
287
288 //- Access to this region
289 const regionFaModel& region() const noexcept { return *this; }
290
291 //- Access to pRef
292 scalar pRef() const { return pRef_; }
293
294 //- Name of the U field
295 word UName() const { return UName_; }
296
297
298 // Transfer fields - to the primary region (lagragian injection)
299
300 //- Return mass transfer source - Eulerian phase only
301 //virtual tmp<volScalarField> primaryMassTrans() const = 0;
303 //- Return the film mass available for transfer to cloud
304 virtual const volScalarField& cloudMassTrans() const = 0;
305
306 //- Return the parcel diameters originating from film to cloud
307 virtual const volScalarField& cloudDiameterTrans() const = 0;
308
309
310
311 // Evolution
312
313 //- Pre-evolve film
314 virtual void preEvolveRegion();
315
316 //- Post-evolve film
317 virtual void postEvolveRegion();
318
319
320 // Thermo variables
321
322 //- Access const reference mu
323 virtual const areaScalarField& mu() const = 0;
324
325 //- Access const reference rho
326 virtual const areaScalarField& rho() const = 0;
327
328 //- Access const reference sigma
329 virtual const areaScalarField& sigma() const = 0;
330
331 //- Access const reference Tf
332 virtual const areaScalarField& Tf() const = 0;
333
334 //- Access const reference Cp
335 virtual const areaScalarField& Cp() const = 0;
336
337
338 // External hook to add sources and mass exchange variables
339
340
341 //- Add sources
342 virtual void addSources
343 (
344 const label patchi, // patchi on primary region
345 const label facei, // facei of patchi
346 const scalar massSource, // [kg]
347 const vector& momentumSource, // [kg.m/s] (tang'l momentum)
348 const scalar pressureSource, // [kg.m/s] (normal momentum)
349 const scalar energySource = 0 // [J]
350 );
351};
352
353
354// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
355
356} // End namespace areaSurfaceFilmModels
357} // End namespace regionModels
358} // End namespace Foam
359
360// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
361
362#endif
363
364// ************************************************************************* //
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
dictionary()
Default construct, a top-level empty dictionary.
Definition dictionary.C:68
Finite-area options, which is an IOdictionary of values and a fa::optionList.
Definition faOptions.H:73
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static autoPtr< liquidFilmBase > New(const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected model using dictionary.
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.
label nFilmCorr_
Number of film thickness correctors.
const dimensionedScalar & h0() const noexcept
Return h0.
liquidFilmBase(const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from type name and mesh and dict.
Switch momentumPredictor_
Flag to enable momentum predictor.
tmp< areaScalarField > pg() const
Map primary static pressure.
const areaScalarField & gn() const noexcept
Access const reference gn.
const edgeScalarField & phif() const noexcept
Access to momentum flux.
volScalarField pnSource_
Normal pressure by particles.
const uniformDimensionedVectorField & g() const noexcept
Gravity.
dimensionedScalar h0_
Smallest film thickness.
const regionFaModel & region() const noexcept
Access to this region.
void operator=(const liquidFilmBase &)=delete
No copy assignment.
tmp< areaVectorField > Uw() const
Wall velocity.
virtual const areaScalarField & Tf() const =0
Access const reference Tf.
Foam::fa::options & faOptions() noexcept
Return faOptions.
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.
label nCorr_
Number of PISO-like inner correctors.
virtual const areaScalarField & rho() const =0
Access const reference rho.
TypeName("liquidFilmBase")
Runtime type information.
const areaVectorField & Uf() const noexcept
Access const reference Uf.
declareRunTimeSelectionTable(autoPtr, liquidFilmBase, dictionary,(const word &modelType, const fvMesh &mesh, const dictionary &dict),(modelType, mesh, dict))
dimensionedScalar deltaWet_
Film thickness beyond which face is assumed to be wet.
const areaScalarField & h() const noexcept
Access const reference h.
virtual scalar CourantNumber() const
Courant number evaluation.
const edgeScalarField & phi2s() const noexcept
Access continuity flux.
tmp< areaVectorField > Us() const
Film surface film velocity.
virtual const volScalarField & cloudDiameterTrans() const =0
Return the parcel diameters originating from film to cloud.
virtual const volScalarField & cloudMassTrans() const =0
Return mass transfer source - Eulerian phase only.
liquidFilmBase(const liquidFilmBase &)=delete
No copy construct.
Base class for liquid-film models.
regionFaModel(const fvMesh &mesh, const word &regionType, const word &modelName, const dictionary &dict, bool readFields=true)
Construct from mesh and name and dict.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
dynamicFvMesh & mesh
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< vector, faPatchField, areaMesh > areaVectorField
GeometricField< scalar, faePatchField, edgeMesh > edgeScalarField
UniformDimensionedField< vector > uniformDimensionedVectorField
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
const direction noexcept
Definition scalarImpl.H:265
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
#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:68
Various UniformDimensionedField types.
Forwards and collection of common volume field types.