thermoSingleLayer.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) 2011-2017 OpenFOAM Foundation
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::surfaceFilmModels::thermoSingleLayer
28
29Description
30 Thermodynamic form of single-cell layer surface film model
31
32 Note: defining enthalpy as Cp(T - Tstd) - when using liquids from the
33 thermophysical library, their enthalpies are calculated similarly, where
34 Tstd = 298.15K. This is clearly non-conservative unless the heat-capacity
35 is constant and should be rewritten to use the standard thermodynamics
36 packages.
37
38SourceFiles
39 thermoSingleLayer.C
40
41\*---------------------------------------------------------------------------*/
42
43#ifndef thermoSingleLayer_H
44#define thermoSingleLayer_H
45
47#include "SLGThermo.H"
48
49// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50
51namespace Foam
52{
53namespace regionModels
54{
55namespace surfaceFilmModels
56{
57
58// Forward declaration of classes
59class filmViscosityModel;
60class heatTransferModel;
61class phaseChangeModel;
62class filmRadiationModel;
63
64/*---------------------------------------------------------------------------*\
65 Class thermoSingleLayer Declaration
66\*---------------------------------------------------------------------------*/
69:
71{
72 // Private member functions
73
74 //- No copy construct
75 thermoSingleLayer(const thermoSingleLayer&) = delete;
76
77 //- No copy assignment
78 void operator=(const thermoSingleLayer&) = delete;
79
80 //- Return boundary types for sensible enthalpy field
81 wordList hsBoundaryTypes();
82
83
84protected:
85
86 // Protected data
87
88 // Thermo properties
89
90 //- Reference to the SLGThermo
91 const SLGThermo& thermo_;
92
93
94 // Fields
95
96 //- Specific heat capacity [J/kg/K]
98
99 //- Thermal conductivity [W/m/K]
101
102 //- Temperature - mean [K]
104
105 //- Temperature - surface [K]
107
108 //- Temperature - wall [K]
110
111 //- Sensible enthalpy [J/kg]
113
114
115 // Transfer fields - to the primary region
116
117 //- Film energy transfer
119
120
121 //- Threshold film thickness beyond which the film is considered 'wet'
122 scalar deltaWet_;
123
124
125 // Hydrophilic/phobic properties
126
127 //- Activation flag
128 bool hydrophilic_;
129
130 //- Length scale applied to deltaWet_ to determine when a wet
131 // surface becomes dry, typically 0.5
133
134 //- Length scale applied to deltaWet_ to determine when a dry
135 // surface becomes wet, typically 0.001
137
138
139 // Source term fields
140
141 // Film region - registered to the film region mesh
142 // Note: need boundary value mapped from primary region, and then
143 // pushed into the patch internal field
144
145 //- Energy [J/m2/s]
147
148
149 // Primary region - registered to the primary region mesh
150 // Internal use only - not read-in
151
152 //- Energy [J/m2/s]
154
155
156 // Fields mapped from primary region - registered to the film region
157 // Note: need both boundary AND patch internal fields to be mapped
158
159 //- Temperature [K]
161
162 //- List of specie mass fractions [0-1]
164
165
166 // Sub-models
167
168 //- Viscosity model
170
171 //- Heat transfer coefficient between film surface and primary
172 // region [W/m2/K]
174
175 //- Heat transfer coefficient between wall and film [W/m2/K]
177
178 //- Phase change
180
181 //- Radiation
183
184
185 // Limits
186
187 //- Minimum temperature limit (optional)
188 scalar Tmin_;
189
190 //- Maximum temperature limit (optional)
191 scalar Tmax_;
192
193
194 // Protected member functions
195
196 //- Read control parameters from dictionary
197 virtual bool read();
198
199 //- Correct the thermo fields
200 virtual void correctThermoFields();
201
202 //- Correct sensible enthalpy for mapped temperature fields
203 virtual void correctHsForMappedT();
204
205 //- Correct the film surface and wall temperatures
206 virtual void updateSurfaceTemperatures();
207
208 //- Reset source term fields
209 virtual void resetPrimaryRegionSourceTerms();
210
211 //- Transfer thermo fields from the primary region to the film region
213
214 //- Transfer source fields from the primary region to the film region
216
217 //- Correct film coverage field
218 virtual void correctAlpha();
219
220 //- Update the film sub-models
221 virtual void updateSubmodels();
222
223 //- Return the wall/surface heat transfer term for the enthalpy equation
224 virtual tmp<fvScalarMatrix> q(volScalarField& hs) const;
225
226
227 // Equations
228
229 //- Solve energy equation
230 virtual void solveEnergy();
231
232
233public:
234
235 //- Runtime type information
236 TypeName("thermoSingleLayer");
237
238
239 // Constructors
240
241 //- Construct from components
243 (
244 const word& modelType,
245 const fvMesh& mesh,
246 const dimensionedVector& g,
247 const word& regionType,
248 const bool readFields = true
249 );
250
251
252 //- Destructor
253 virtual ~thermoSingleLayer();
254
255
256 // Member Functions
257
258 // Thermo properties
259
260 //- Return const reference to the SLGThermo object
261 inline const SLGThermo& thermo() const;
262
263
264 // Fields
265
266 //- Return the film specific heat capacity [J/kg/K]
267 virtual const volScalarField& Cp() const;
268
269 //- Return the film thermal conductivity [W/m/K]
270 virtual const volScalarField& kappa() const;
271
272 //- Return the film mean temperature [K]
273 virtual const volScalarField& T() const;
274
275 //- Return the film surface temperature [K]
276 virtual const volScalarField& Ts() const;
277
278 //- Return the film wall temperature [K]
279 virtual const volScalarField& Tw() const;
280
281 //- Return the film sensible enthalpy [J/kg]
282 virtual const volScalarField& hs() const;
283
284
285 // Helper functions
286
287 //- Return sensible enthalpy as a function of temperature
288 // for a patch
289 inline tmp<scalarField> hs
290 (
291 const scalarField& T,
292 const label patchi
293 ) const;
294
295 //- Return sensible enthalpy as a function of temperature
297 (
298 const volScalarField& T
299 ) const;
300
301 //- Return temperature as a function of sensible enthalpy
303 (
304 const volScalarField& hs
305 ) const;
306
307
308 // Source fields (read/write access)
309
310 //- External hook to add sources to the film
311 virtual void addSources
312 (
313 const label patchi, // patchi on primary region
314 const label facei, // facei of patchi
315 const scalar massSource, // [kg]
316 const vector& momentumSource, // [kg.m/s] (tangential momentum)
317 const scalar pressureSource, // [kg.m/s] (normal momentum)
318 const scalar energySource // [J]
319 );
320
321
322 // Source term fields
323
324 // Film region
325
326 //- Energy [J/m2/s]
327 inline const volScalarField& hsSp() const;
328
329
330 // Primary region
331
332 //- Energy [J/m2/s]
333 inline const volScalarField& hsSpPrimary() const;
334
335
336 // Fields mapped from the primary region
337
338 //- Temperature [K]
339 inline const volScalarField& TPrimary() const;
340
341 //- Specie mass fractions [0-1]
342 inline const PtrList<volScalarField>& YPrimary() const;
343
344
345
346 // Sub-models
347
348 //- Return const access to the (surface) heat transfer model
349 inline const heatTransferModel& htcs() const;
350
351 //- Return const access to the (wall) heat transfer model
352 inline const heatTransferModel& htcw() const;
353
354 //- Return const access to the phase change model
355 inline const phaseChangeModel& phaseChange() const;
356
357 //- Return const access to the radiation model
358 inline const filmRadiationModel& radiation() const;
359
360
361 // Derived fields (calculated on-the-fly)
362
363 //- Return the convective heat energy from film to wall
364 inline tmp<scalarField> qconvw(const label patchi) const;
365
366 //- Return the convective heat energy from primary region to film
367 inline tmp<scalarField> qconvp(const label patchi) const;
368
369
370 // Evolution
371
372 //- Pre-evolve film hook
373 virtual void preEvolveRegion();
374
375 //- Evolve the film equations
376 virtual void evolveRegion();
377
378
379 // Source fields
380
381 // Mapped into primary region
382
383 //- Return total mass source - Eulerian phase only
384 virtual tmp<volScalarField::Internal> Srho() const;
385
386 //- Return mass source for specie i - Eulerian phase only
388 (
389 const label i
390 ) const;
391
392 //- Return enthalpy source - Eulerian phase only
393 virtual tmp<volScalarField::Internal> Sh() const;
394
395
396 // I-O
397
398 //- Provide some feedback
399 virtual void info();
400};
401
402
403// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
404
405} // End namespace surfaceFilmModels
406} // End namespace regionModels
407} // End namespace Foam
408
409// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
410
411#include "thermoSingleLayerI.H"
412
413// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
414
415#endif
416
417// ************************************************************************* //
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:67
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Base class for film heat transfer models.
Kinematic form of single-cell layer surface film model.
Base class for surface film phase change models.
const dimensionedVector & g() const
Return the acceleration due to gravity.
Thermodynamic form of single-cell layer surface film model.
virtual void transferPrimaryRegionThermoFields()
Transfer thermo fields from the primary region to the film region.
volScalarField kappa_
Thermal conductivity [W/m/K].
autoPtr< heatTransferModel > htcw_
Heat transfer coefficient between wall and film [W/m2/K].
virtual const volScalarField & hs() const
Return the film sensible enthalpy [J/kg].
virtual tmp< fvScalarMatrix > q(volScalarField &hs) const
Return the wall/surface heat transfer term for the enthalpy equation.
virtual void updateSurfaceTemperatures()
Correct the film surface and wall temperatures.
autoPtr< filmRadiationModel > radiation_
Radiation.
virtual void addSources(const label patchi, const label facei, const scalar massSource, const vector &momentumSource, const scalar pressureSource, const scalar energySource)
External hook to add sources to the film.
const SLGThermo & thermo() const
Return const reference to the SLGThermo object.
tmp< scalarField > qconvp(const label patchi) const
Return the convective heat energy from primary region to film.
virtual const volScalarField & kappa() const
Return the film thermal conductivity [W/m/K].
virtual const volScalarField & T() const
Return the film mean temperature [K].
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
volScalarField Cp_
Specific heat capacity [J/kg/K].
const PtrList< volScalarField > & YPrimary() const
Specie mass fractions [0-1].
scalar hydrophilicWetScale_
Length scale applied to deltaWet_ to determine when a dry.
virtual const volScalarField & Tw() const
Return the film wall temperature [K].
const heatTransferModel & htcw() const
Return const access to the (wall) heat transfer model.
const volScalarField & hsSp() const
Energy [J/m2/s].
const phaseChangeModel & phaseChange() const
Return const access to the phase change model.
virtual const volScalarField & Cp() const
Return the film specific heat capacity [J/kg/K].
virtual void resetPrimaryRegionSourceTerms()
Reset source term fields.
scalar deltaWet_
Threshold film thickness beyond which the film is considered 'wet'.
virtual void correctAlpha()
Correct film coverage field.
PtrList< volScalarField > YPrimary_
List of specie mass fractions [0-1].
const heatTransferModel & htcs() const
Return const access to the (surface) heat transfer model.
volScalarField primaryEnergyTrans_
Film energy transfer.
const filmRadiationModel & radiation() const
Return const access to the radiation model.
const volScalarField & TPrimary() const
Temperature [K].
scalar Tmax_
Maximum temperature limit (optional)
scalar Tmin_
Minimum temperature limit (optional)
scalar hydrophilicDryScale_
Length scale applied to deltaWet_ to determine when a wet.
virtual void correctThermoFields()
Correct the thermo fields.
virtual void updateSubmodels()
Update the film sub-models.
virtual void transferPrimaryRegionSourceFields()
Transfer source fields from the primary region to the film region.
TypeName("thermoSingleLayer")
Runtime type information.
const SLGThermo & thermo_
Reference to the SLGThermo.
const volScalarField & hsSpPrimary() const
Energy [J/m2/s].
autoPtr< phaseChangeModel > phaseChange_
Phase change.
tmp< scalarField > qconvw(const label patchi) const
Return the convective heat energy from film to wall.
virtual tmp< volScalarField::Internal > Srho() const
Return total mass source - Eulerian phase only.
autoPtr< filmViscosityModel > viscosity_
Viscosity model.
autoPtr< heatTransferModel > htcs_
Heat transfer coefficient between film surface and primary.
virtual bool read()
Read control parameters from dictionary.
virtual void correctHsForMappedT()
Correct sensible enthalpy for mapped temperature fields.
virtual tmp< volScalarField::Internal > Sh() const
Return enthalpy source - Eulerian phase only.
virtual void evolveRegion()
Evolve the film equations.
volScalarField Ts_
Temperature - surface [K].
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
Namespace for OpenFOAM.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73