-
Notifications
You must be signed in to change notification settings - Fork 1
/
surfaceFilmModel.H
255 lines (174 loc) · 7.22 KB
/
surfaceFilmModel.H
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::regionModels::surfaceFilmModels::surfaceFilmModel
Description
Base class for surface film models
SourceFiles
surfaceFilmModelI.H
surfaceFilmModel.C
surfaceFilmModelNew.C
\*---------------------------------------------------------------------------*/
#ifndef surfaceFilmModel_H
#define surfaceFilmModel_H
#include "singleLayerRegion.H"
#include "dimensionedVector.H"
#include "runTimeSelectionTables.H"
#include "volFieldsFwd.H"
#include "DimensionedField.H"
#include "labelList.H"
#include "NamedEnum.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace regionModels
{
namespace surfaceFilmModels
{
/*---------------------------------------------------------------------------*\
Class surfaceFilmModel Declaration
\*---------------------------------------------------------------------------*/
class surfaceFilmModel
:
public singleLayerRegion
{
private:
// Private Member Functions
//- Disallow default bitwise copy construct
surfaceFilmModel(const surfaceFilmModel&);
//- Disallow default bitwise assignment
void operator=(const surfaceFilmModel&);
protected:
// Protected data
//- Acceleration due to gravity [m/s2]
const dimensionedVector& g_;
// Protected member functions
//- Read control parameters from dictionary
virtual bool read();
public:
//- Runtime type information
TypeName("surfaceFilmModel");
//- Reference temperature for enthalpy
static const dimensionedScalar Tref;
// Declare runtime constructor selection table
declareRunTimeSelectionTable
(
autoPtr,
surfaceFilmModel,
mesh,
(
const word& modelType,
const fvMesh& mesh,
const dimensionedVector& g,
const word& regionType
),
(modelType, mesh, g, regionType)
);
// Constructors
//- Construct from type name, mesh and gravity vector
surfaceFilmModel
(
const word& modelType,
const fvMesh& mesh,
const dimensionedVector& g,
const word& regionType
);
// Selectors
//- Return a reference to the selected surface film model
static autoPtr<surfaceFilmModel> New
(
const fvMesh& mesh,
const dimensionedVector& g,
const word& regionType = "surfaceFilm"
);
//- Destructor
virtual ~surfaceFilmModel();
// Member Functions
// Access
//- Return the accleration due to gravity
inline const dimensionedVector& g() const;
//- External hook to add sources to the film
virtual void addSources
(
const label patchi,
const label facei,
const scalar massSource,
const vector& momentumSource,
const scalar pressureSource,
const scalar energySource
) = 0;
// Solution parameters
//- Courant number evaluation
virtual scalar CourantNumber() const;
// Fields
//- Return the film thickness [m]
virtual const volScalarField& delta() const = 0;
//- Return the film coverage, 1 = covered, 0 = uncovered / []
virtual const volScalarField& alpha() const = 0;
//- Return the film velocity [m/s]
virtual const volVectorField& U() const = 0;
//- Return the film surface velocity [m/s]
virtual const volVectorField& Us() const = 0;
//- Return the film wall velocity [m/s]
virtual const volVectorField& Uw() const = 0;
//- Return the film density [kg/m3]
virtual const volScalarField& rho() const = 0;
//- Return the film mean temperature [K]
virtual const volScalarField& T() const = 0;
//- Return the film surface temperature [K]
virtual const volScalarField& Ts() const = 0;
//- Return the film wall temperature [K]
virtual const volScalarField& Tw() const = 0;
//- Return the film surface temperature [J/kg]
virtual const volScalarField& hs() const = 0;
//- Return the film specific heat capacity [J/kg/K]
virtual const volScalarField& Cp() const = 0;
//- Return the film thermal conductivity [W/m/K]
virtual const volScalarField& kappa() const = 0;
//- Return the film surface tension [N/m]
virtual const volScalarField& sigma() const = 0;
// Transfer fields - to the primary region
//- Return mass transfer source - Eulerian phase only
virtual tmp<volScalarField> primaryMassTrans() const = 0;
//- Return the film mass available for transfer
virtual const volScalarField& cloudMassTrans() const = 0;
//- Return the parcel diameters originating from film
virtual const volScalarField& cloudDiameterTrans() const = 0;
// Source fields
// Mapped into primary region
//- Return total mass source - Eulerian phase only
virtual tmp<volScalarField::Internal> Srho() const;
//- Return mass source for specie i - Eulerian phase only
virtual tmp<volScalarField::Internal> Srho
(
const label i
) const;
//- Return enthalpy source - Eulerian phase only
virtual tmp<volScalarField::Internal> Sh() const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace surfaceFilmModels
} // End namespace regionModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#include "surfaceFilmModelI.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //