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-2021 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 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 // Forward declaration of classes
52 class fluidThermo;
53 
54 namespace regionModels
55 {
56 namespace surfaceFilmModels
57 {
58 
59 // Forward declaration of classes
60 class viscosityModel;
61 class heatTransferModel;
62 class phaseChangeModel;
63 class radiationModel;
64 
65 /*---------------------------------------------------------------------------*\
66  Class thermoSingleLayer Declaration
67 \*---------------------------------------------------------------------------*/
68 
70 :
72 {
73 
74 protected:
75 
76  // Protected data
77 
78  const word phaseName_;
79 
80  // Thermo properties
81 
82  //- Reference to the primary region thermo
84 
85  // Transfer fields - to the primary region
86 
87  //- Film energy transfer
89 
90 
91  //- Threshold film thickness beyond which the film is considered 'wet'
92  scalar deltaWet_;
93 
94 
95  // Hydrophilic/phobic properties
96 
97  //- Activation flag
98  bool hydrophilic_;
99 
100  //- Length scale applied to deltaWet_ to determine when a wet
101  // surface becomes dry, typically 0.5
102  scalar hydrophilicDryScale_;
103 
104  //- Length scale applied to deltaWet_ to determine when a dry
105  // surface becomes wet, typically 0.001
106  scalar hydrophilicWetScale_;
107 
108 
109  // Source term fields
110 
111  // Film region - registered to the film region mesh
112  // Note: need boundary value mapped from primary region, and then
113  // pushed into the patch internal field
114 
115  //- Energy [J/m2/s]
117 
118 
119  // Primary region - registered to the primary region mesh
120  // Internal use only - not read-in
121 
122  //- Energy [J/m2/s]
124 
125 
126  // Fields mapped from primary region - registered to the film region
127  // Note: need both boundary AND patch internal fields to be mapped
128 
129  //- Temperature [K]
131 
132  //- List of specie mass fractions [0-1]
134 
135 
136  // Sub-models
137 
138  //- Heat transfer coefficient between film surface and primary
139  // region [W/m^2/K]
141 
142  //- Heat transfer coefficient between wall and film [W/m^2/K]
144 
145  //- Phase change
147 
148  //- Radiation
150 
151 
152  // Limits
153 
154  //- Minimum temperature limit (optional)
155  scalar Tmin_;
156 
157  //- Maximum temperature limit (optional)
158  scalar Tmax_;
159 
160 
161  // Protected member functions
162 
163  //- Read control parameters from dictionary
164  virtual bool read();
165 
166  //- Correct sensible enthalpy for mapped temperature fields
167  virtual void correctHforMappedT();
168 
169  //- Reset source term fields
170  virtual void resetPrimaryRegionSourceTerms();
171 
172  //- Transfer thermo fields from the primary region to the film region
173  virtual void transferPrimaryRegionThermoFields();
174 
175  //- Transfer source fields from the primary region to the film region
176  virtual void transferPrimaryRegionSourceFields();
177 
178  //- Correct film coverage field
179  virtual void correctCoverage();
180 
181  //- Update the film sub-models
182  virtual void updateSubmodels();
183 
184  //- Return the wall/surface heat transfer term for the enthalpy equation
185  virtual tmp<fvScalarMatrix> q(volScalarField& h) const;
186 
187 
188  // Equations
189 
190  //- Solve energy equation
191  virtual void solveEnergy();
192 
193 
194 public:
195 
196  //- Runtime type information
197  TypeName("thermoSingleLayer");
198 
199 
200  // Constructors
201 
202  //- Construct from components
204  (
205  const word& modelType,
206  const fvMesh& mesh,
207  const dimensionedVector& g,
208  const word& regionType,
209  const bool readFields = true
210  );
211 
212  //- Disallow default bitwise copy construction
213  thermoSingleLayer(const thermoSingleLayer&) = delete;
214 
215 
216  //- Destructor
217  virtual ~thermoSingleLayer();
218 
219 
220  // Member Functions
221 
222  // Thermo properties
223 
224  //- Return const reference to the primary region thermo object
225  inline const fluidThermo& primaryThermo() const;
226 
227 
228  // Derived Fields
229 
230  //- Return the film surface temperature [K]
231  // Currently this is assumed to be equal to
232  // the film mean temperature
233  virtual tmp<volScalarField::Internal> Ts() const;
234 
235  //- Return the film wall temperature [K]
236  virtual tmp<volScalarField::Internal> Tw() const;
237 
238 
239  // Source fields (read/write access)
240 
241  //- External hook to add sources to the film
242  virtual void addSources
243  (
244  const label patchi, // patchi on primary region
245  const label facei, // facei of patchi
246  const scalar massSource, // [kg]
247  const vector& momentumSource, // [kg m/s] (tangential momentum)
248  const scalar pressureSource, // [kg m/s] (normal momentum)
249  const scalar energySource // [J]
250  );
251 
252 
253  // Source term fields
254 
255  // Film region
256 
257  //- Energy [J/m2/s]
258  inline const volScalarField::Internal& hSp() const;
259 
260 
261  // Primary region
262 
263  //- Energy [J/m2/s]
264  inline const volScalarField& hSpPrimary() const;
265 
266 
267  // Fields mapped from the primary region
268 
269  //- Temperature [K]
270  inline const volScalarField& TPrimary() const;
271 
272  //- Specie mass fractions [0-1]
273  inline const PtrList<volScalarField>& YPrimary() const;
274 
275 
276 
277  // Sub-models
278 
279  //- Return const access to the (surface) heat transfer model
280  inline const heatTransferModel& htcs() const;
281 
282  //- Return const access to the (wall) heat transfer model
283  inline const heatTransferModel& htcw() const;
284 
285  //- Return const access to the phase change model
286  inline const phaseChangeModel& phaseChange() const;
287 
288  //- Return const access to the radiation model
289  inline const radiationModel& radiation() const;
290 
291 
292  // Evolution
293 
294  //- Pre-evolve film hook
295  virtual void preEvolveRegion();
296 
297  //- Evolve the film equations
298  virtual void evolveRegion();
299 
300 
301  // Source fields
302 
303  // Mapped into primary region
304 
305  //- Return mass source for specie i - Eulerian phase only
307  (
308  const label i
309  ) const;
310 
311  //- Return enthalpy source - Eulerian phase only
312  virtual tmp<volScalarField::Internal> Sh() const;
313 
314 
315  // I-O
316 
317  //- Provide some feedback
318  virtual void info();
319 
320 
321  // Member Operators
322 
323  //- Disallow default bitwise assignment
324  void operator=(const thermoSingleLayer&) = delete;
325 };
326 
327 
328 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
329 
330 } // End namespace surfaceFilmModels
331 } // End namespace regionModels
332 } // End namespace Foam
333 
334 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
335 
336 #include "thermoSingleLayerI.H"
337 
338 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
339 
340 #endif
341 
342 // ************************************************************************* //
virtual tmp< volScalarField::Internal > Ts() const
Return the film surface temperature [K].
autoPtr< heatTransferModel > htcw_
Heat transfer coefficient between wall and film [W/m^2/K].
const heatTransferModel & htcs() const
Return const access to the (surface) heat transfer model.
autoPtr< phaseChangeModel > phaseChange_
Phase change.
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.
virtual void correctCoverage()
Correct film coverage field.
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.
const volScalarField & hSpPrimary() const
Energy [J/m2/s].
virtual tmp< fvScalarMatrix > q(volScalarField &h) const
Return the wall/surface heat transfer term for the enthalpy equation.
volScalarField::Internal hSp_
Energy [J/m2/s].
virtual void transferPrimaryRegionThermoFields()
Transfer thermo fields from the primary region to the film region.
autoPtr< heatTransferModel > htcs_
Heat transfer coefficient between film surface and primary.
const dimensionedScalar h
Planck constant.
Base class for surface film phase change models.
const volScalarField & TPrimary() const
Temperature [K].
void operator=(const thermoSingleLayer &)=delete
Disallow default bitwise assignment.
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
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:53
const fluidThermo & primaryThermo_
Reference to the primary region thermo.
const PtrList< volScalarField > & YPrimary() const
Specie mass fractions [0-1].
const radiationModel & radiation() const
Return const access to the radiation model.
virtual tmp< volScalarField::Internal > SYi(const label i) const
Return mass source for specie i - Eulerian phase only.
scalar hydrophilicWetScale_
Length scale applied to deltaWet_ to determine when a dry.
const volScalarField::Internal & hSp() const
Energy [J/m2/s].
scalar Tmin_
Minimum temperature limit (optional)
const phaseChangeModel & phaseChange() const
Return const access to the phase change model.
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 bool read()
Read control parameters from dictionary.
virtual void correctHforMappedT()
Correct sensible enthalpy for mapped temperature fields.
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;.
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.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
PtrList< volScalarField > YPrimary_
List of specie mass fractions [0-1].
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual tmp< volScalarField::Internal > Tw() const
Return the film wall temperature [K].
A class for managing temporary objects.
Definition: PtrList.H:53
const fluidThermo & primaryThermo() const
Return const reference to the primary region thermo object.
virtual tmp< volScalarField::Internal > Sh() const
Return enthalpy source - Eulerian phase only.
Thermodynamic form of single-cell layer surface film model.
Namespace for OpenFOAM.
Base class for film heat transfer models.