thermoSingleLayer.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::regionModels::surfaceFilmModels::thermoSingleLayer
26 
27 Description
28  Thermodynamic form of single-cell layer surface film model
29 
30  Note: defining enthalpy as Cp(T - Tstd) - when using liquids from the
31  thermophysical library, their enthalpies are calculated similarly, where
32  Tstd = 298.15K. This is clearly non-conservative unless the heat-capacity
33  is constant and should be rewritten to use the standard thermodynamics
34  packages.
35 
36 SourceFiles
37  thermoSingleLayer.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef thermoSingleLayer_H
42 #define thermoSingleLayer_H
43 
44 #include "kinematicSingleLayer.H"
45 #include "SLGThermo.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 namespace regionModels
52 {
53 namespace surfaceFilmModels
54 {
55 
56 // Forward declaration of classes
57 class filmViscosityModel;
58 class heatTransferModel;
59 class phaseChangeModel;
60 class filmRadiationModel;
61 
62 /*---------------------------------------------------------------------------*\
63  Class thermoSingleLayer Declaration
64 \*---------------------------------------------------------------------------*/
65 
67 :
69 {
70  // Private Member Functions
71 
72  //- Return boundary types for sensible enthalpy field
73  wordList hsBoundaryTypes();
74 
75 
76 protected:
77 
78  // Protected data
79 
80  // Thermo properties
81 
82  //- Reference to the SLGThermo
83  const SLGThermo& thermo_;
84 
85 
86  // Fields
87 
88  //- Specific heat capacity [J/kg/K]
90 
91  //- Thermal conductivity [W/m/K]
93 
94  //- Temperature - mean [K]
96 
97  //- Temperature - surface [K]
99 
100  //- Temperature - wall [K]
102 
103  //- Sensible enthalpy [J/kg]
105 
106 
107  // Transfer fields - to the primary region
108 
109  //- Film energy transfer
111 
112 
113  //- Threshold film thickness beyond which the film is considered 'wet'
114  scalar deltaWet_;
115 
116 
117  // Hyprophilic/phobic properties
118 
119  //- Activation flag
120  bool hydrophilic_;
121 
122  //- Length scale applied to deltaWet_ to determine when a wet
123  // surface becomes dry, typically 0.5
124  scalar hydrophilicDryScale_;
125 
126  //- Length scale applied to deltaWet_ to determine when a dry
127  // surface becomes wet, typically 0.001
128  scalar hydrophilicWetScale_;
129 
130 
131  // Source term fields
132 
133  // Film region - registered to the film region mesh
134  // Note: need boundary value mapped from primary region, and then
135  // pushed into the patch internal field
136 
137  //- Energy [J/m2/s]
139 
140 
141  // Primary region - registered to the primary region mesh
142  // Internal use only - not read-in
143 
144  //- Energy [J/m2/s]
146 
147 
148  // Fields mapped from primary region - registered to the film region
149  // Note: need both boundary AND patch internal fields to be mapped
150 
151  //- Temperature [K]
153 
154  //- List of specie mass fractions [0-1]
156 
157 
158  // Sub-models
159 
160  //- Viscosity model
162 
163  //- Heat transfer coefficient between film surface and primary
164  // region [W/m^2/K]
166 
167  //- Heat transfer coefficient between wall and film [W/m^2/K]
169 
170  //- Phase change
172 
173  //- Radiation
175 
176 
177  // Limits
178 
179  //- Minimum temperature limit (optional)
180  scalar Tmin_;
181 
182  //- Maximum temperature limit (optional)
183  scalar Tmax_;
184 
185 
186  // Protected member functions
187 
188  //- Read control parameters from dictionary
189  virtual bool read();
190 
191  //- Correct the thermo fields
192  virtual void correctThermoFields();
193 
194  //- Correct sensible enthalpy for mapped temperature fields
195  virtual void correctHsForMappedT();
196 
197  //- Correct the film surface and wall temperatures
198  virtual void updateSurfaceTemperatures();
199 
200  //- Reset source term fields
201  virtual void resetPrimaryRegionSourceTerms();
202 
203  //- Transfer thermo fields from the primary region to the film region
204  virtual void transferPrimaryRegionThermoFields();
205 
206  //- Transfer source fields from the primary region to the film region
207  virtual void transferPrimaryRegionSourceFields();
208 
209  //- Correct film coverage field
210  virtual void correctAlpha();
211 
212  //- Update the film sub-models
213  virtual void updateSubmodels();
214 
215  //- Return the wall/surface heat transfer term for the enthalpy equation
216  virtual tmp<fvScalarMatrix> q(volScalarField& hs) const;
217 
218 
219  // Equations
220 
221  //- Solve energy equation
222  virtual void solveEnergy();
223 
224 
225 public:
226 
227  //- Runtime type information
228  TypeName("thermoSingleLayer");
229 
230 
231  // Constructors
232 
233  //- Construct from components
235  (
236  const word& modelType,
237  const fvMesh& mesh,
238  const dimensionedVector& g,
239  const word& regionType,
240  const bool readFields = true
241  );
242 
243  //- Disallow default bitwise copy construction
244  thermoSingleLayer(const thermoSingleLayer&) = delete;
245 
246 
247  //- Destructor
248  virtual ~thermoSingleLayer();
249 
250 
251  // Member Functions
252 
253  // Thermo properties
254 
255  //- Return const reference to the SLGThermo object
256  inline const SLGThermo& thermo() const;
257 
258 
259  // Fields
260 
261  //- Return the film specific heat capacity [J/kg/K]
262  virtual const volScalarField& Cp() const;
263 
264  //- Return the film thermal conductivity [W/m/K]
265  virtual const volScalarField& kappa() const;
266 
267  //- Return the film mean temperature [K]
268  virtual const volScalarField& T() const;
269 
270  //- Return the film surface temperature [K]
271  virtual const volScalarField& Ts() const;
272 
273  //- Return the film wall temperature [K]
274  virtual const volScalarField& Tw() const;
275 
276  //- Return the film sensible enthalpy [J/kg]
277  virtual const volScalarField& hs() const;
278 
279 
280  // Helper functions
281 
282  //- Return sensible enthalpy as a function of temperature
283  // for a patch
284  inline tmp<scalarField> hs
285  (
286  const scalarField& T,
287  const label patchi
288  ) const;
289 
290  //- Return sensible enthalpy as a function of temperature
291  inline tmp<volScalarField> hs
292  (
293  const volScalarField& T
294  ) const;
295 
296  //- Return temperature as a function of sensible enthalpy
297  inline tmp<volScalarField> T
298  (
299  const volScalarField& hs
300  ) const;
301 
302 
303  // Source fields (read/write access)
304 
305  //- External hook to add sources to the film
306  virtual void addSources
307  (
308  const label patchi, // patchi on primary region
309  const label facei, // facei of patchi
310  const scalar massSource, // [kg]
311  const vector& momentumSource, // [kg m/s] (tangential momentum)
312  const scalar pressureSource, // [kg m/s] (normal momentum)
313  const scalar energySource // [J]
314  );
315 
316 
317  // Source term fields
318 
319  // Film region
320 
321  //- Energy [J/m2/s]
322  inline const volScalarField& hsSp() const;
323 
324 
325  // Primary region
326 
327  //- Energy [J/m2/s]
328  inline const volScalarField& hsSpPrimary() const;
329 
330 
331  // Fields mapped from the primary region
332 
333  //- Temperature [K]
334  inline const volScalarField& TPrimary() const;
335 
336  //- Specie mass fractions [0-1]
337  inline const PtrList<volScalarField>& YPrimary() const;
338 
339 
340 
341  // Sub-models
342 
343  //- Return const access to the (surface) heat transfer model
344  inline const heatTransferModel& htcs() const;
345 
346  //- Return const access to the (wall) heat transfer model
347  inline const heatTransferModel& htcw() const;
348 
349  //- Return const access to the phase change model
350  inline const phaseChangeModel& phaseChange() const;
351 
352  //- Return const access to the radiation model
353  inline const filmRadiationModel& radiation() const;
354 
355 
356  // Derived fields (calculated on-the-fly)
357 
358  //- Return the convective heat energy from film to wall
359  inline tmp<scalarField> Qconvw(const label patchi) const;
360 
361  //- Return the convective heat energy from primary region to film
362  inline tmp<scalarField> Qconvp(const label patchi) const;
363 
364 
365  // Evolution
366 
367  //- Pre-evolve film hook
368  virtual void preEvolveRegion();
369 
370  //- Evolve the film equations
371  virtual void evolveRegion();
372 
373 
374  // Source fields
375 
376  // Mapped into primary region
377 
378  //- Return total mass source - Eulerian phase only
379  virtual tmp<volScalarField::Internal> Srho() const;
380 
381  //- Return mass source for specie i - Eulerian phase only
383  (
384  const label i
385  ) const;
386 
387  //- Return enthalpy source - Eulerian phase only
388  virtual tmp<volScalarField::Internal> Sh() const;
389 
390 
391  // I-O
392 
393  //- Provide some feedback
394  virtual void info();
395 
396 
397  // Member Operators
398 
399  //- Disallow default bitwise assignment
400  void operator=(const thermoSingleLayer&) = delete;
401 };
402 
403 
404 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
405 
406 } // End namespace surfaceFilmModels
407 } // End namespace regionModels
408 } // End namespace Foam
409 
410 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
411 
412 #include "thermoSingleLayerI.H"
413 
414 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
415 
416 #endif
417 
418 // ************************************************************************* //
volScalarField Ts_
Temperature - surface [K].
autoPtr< heatTransferModel > htcw_
Heat transfer coefficient between wall and film [W/m^2/K].
virtual void correctAlpha()
Correct film coverage field.
const heatTransferModel & htcs() const
Return const access to the (surface) heat transfer model.
const volScalarField & hsSpPrimary() const
Energy [J/m2/s].
autoPtr< phaseChangeModel > phaseChange_
Phase change.
const volScalarField & hsSp() const
Energy [J/m2/s].
virtual const volScalarField & Cp() const
Return the film specific heat capacity [J/kg/K].
thermoSingleLayer(const word &modelType, const fvMesh &mesh, const dimensionedVector &g, const word &regionType, const bool readFields=true)
Construct from components.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual void resetPrimaryRegionSourceTerms()
Reset source term fields.
Kinematic form of single-cell layer surface film model.
scalar Tmax_
Maximum temperature limit (optional)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
virtual void transferPrimaryRegionSourceFields()
Transfer source fields from the primary region to the film region.
virtual const volScalarField & kappa() const
Return the film thermal conductivity [W/m/K].
virtual void transferPrimaryRegionThermoFields()
Transfer thermo fields from the primary region to the film region.
tmp< scalarField > Qconvp(const label patchi) const
Return the convective heat energy from primary region to film.
autoPtr< heatTransferModel > htcs_
Heat transfer coefficient between film surface and primary.
virtual void updateSurfaceTemperatures()
Correct the film surface and wall temperatures.
virtual tmp< volScalarField::Internal > Srho() const
Return total mass source - Eulerian phase only.
Base class for surface film phase change models.
const volScalarField & TPrimary() const
Temperature [K].
void operator=(const thermoSingleLayer &)=delete
Disallow default bitwise assignment.
virtual tmp< fvScalarMatrix > q(volScalarField &hs) const
Return the wall/surface heat transfer term for the enthalpy equation.
virtual void updateSubmodels()
Update the film sub-models.
dynamicFvMesh & mesh
TypeName("thermoSingleLayer")
Runtime type information.
scalar hydrophilicDryScale_
Length scale applied to deltaWet_ to determine when a wet.
A class for handling words, derived from string.
Definition: word.H:59
volScalarField kappa_
Thermal conductivity [W/m/K].
const PtrList< volScalarField > & YPrimary() const
Specie mass fractions [0-1].
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
Definition: SLGThermo.H:62
scalar hydrophilicWetScale_
Length scale applied to deltaWet_ to determine when a dry.
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
scalar Tmin_
Minimum temperature limit (optional)
const phaseChangeModel & phaseChange() const
Return const access to the phase change model.
const SLGThermo & thermo_
Reference to the SLGThermo.
const dimensionedVector & g() const
Return the acceleration due to gravity.
virtual void evolveRegion()
Evolve the film equations.
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.
volScalarField primaryEnergyTrans_
Film energy transfer.
virtual void correctThermoFields()
Correct the thermo fields.
autoPtr< filmRadiationModel > radiation_
Radiation.
virtual bool read()
Read control parameters from dictionary.
const heatTransferModel & htcw() const
Return const access to the (wall) heat transfer model.
label patchi
scalar deltaWet_
Threshold film thickness beyond which the film is considered &#39;wet&#39;.
volScalarField Cp_
Specific heat capacity [J/kg/K].
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual void preEvolveRegion()
Pre-evolve film hook.
PtrList< volScalarField > YPrimary_
List of specie mass fractions [0-1].
virtual const volScalarField & hs() const
Return the film sensible enthalpy [J/kg].
autoPtr< filmViscosityModel > viscosity_
Viscosity model.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
const filmRadiationModel & radiation() const
Return const access to the radiation model.
const SLGThermo & thermo() const
Return const reference to the SLGThermo object.
virtual void correctHsForMappedT()
Correct sensible enthalpy for mapped temperature fields.
volScalarField hs_
Sensible enthalpy [J/kg].
virtual tmp< volScalarField::Internal > Sh() const
Return enthalpy source - Eulerian phase only.
tmp< scalarField > Qconvw(const label patchi) const
Return the convective heat energy from film to wall.
Thermodynamic form of single-cell layer surface film model.
Namespace for OpenFOAM.
Base class for film heat transfer models.
virtual const volScalarField & T() const
Return the film mean temperature [K].
virtual const volScalarField & Tw() const
Return the film wall temperature [K].