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