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 | Copyright (C) 2011-2016 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::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.15 K
33 
34 SourceFiles
35  thermoSingleLayer.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef thermoSingleLayer_H
40 #define thermoSingleLayer_H
41 
42 #include "kinematicSingleLayer.H"
43 #include "SLGThermo.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 namespace regionModels
50 {
51 namespace surfaceFilmModels
52 {
53 
54 // Forward declaration of classes
55 class filmViscosityModel;
56 class heatTransferModel;
57 class phaseChangeModel;
58 class filmRadiationModel;
59 
60 /*---------------------------------------------------------------------------*\
61  Class thermoSingleLayer Declaration
62 \*---------------------------------------------------------------------------*/
63 
65 :
67 {
68 private:
69 
70  // Private member functions
71 
72  //- Disallow default bitwise copy construct
74 
75  //- Disallow default bitwise assignment
76  void operator=(const thermoSingleLayer&);
77 
78  //- Return boundary types for sensible enthalpy field
79  wordList hsBoundaryTypes();
80 
81 
82 protected:
83 
84  // Protected data
85 
86  // Thermo properties
87 
88  //- Reference to the SLGThermo
89  const SLGThermo& thermo_;
90 
91 
92  // Fields
93 
94  //- Specific heat capacity / [J/kg/K]
96 
97  //- Thermal conductivity / [W/m/K]
99 
100  //- Temperature - mean / [K]
102 
103  //- Temperature - surface / [K]
105 
106  //- Temperature - wall / [K]
108 
109  //- Sensible enthalpy / [J/kg]
111 
112 
113  // Transfer fields - to the primary region
114 
115  //- Film mass evolved via phase change
117 
118  //- Film energy evolved via phase change
120 
121 
122  //- Threshold film thickness beyond which the film is considered 'wet'
123  scalar deltaWet_;
124 
125 
126  // Hyprophilic/phobic properties
127 
128  //- Activation flag
129  bool hydrophilic_;
130 
131  //- Length scale applied to deltaWet_ to determine when a wet
132  // surface becomes dry, typically 0.5
133  scalar hydrophilicDryScale_;
134 
135  //- Length scale applied to deltaWet_ to determine when a dry
136  // surface becomes wet, typically 0.001
137  scalar hydrophilicWetScale_;
138 
139 
140  // Source term fields
141 
142  // Film region - registered to the film region mesh
143  // Note: need boundary value mapped from primary region, and then
144  // pushed into the patch internal field
145 
146  //- Energy / [J/m2/s]
148 
149 
150  // Primary region - registered to the primary region mesh
151  // Internal use only - not read-in
152 
153  //- Energy / [J/m2/s]
155 
156 
157  // Fields mapped from primary region - registered to the film region
158  // Note: need both boundary AND patch internal fields to be mapped
159 
160  //- Temperature / [K]
162 
163  //- List of specie mass fractions / [0-1]
165 
166 
167  // Sub-models
168 
169  //- Viscosity model
171 
172  //- Heat transfer coefficient bewteen film surface and primary
173  // region [W/m2/K]
175 
176  //- Heat transfer coefficient bewteen wall and film [W/m2/K]
178 
179  //- Phase change
181 
182  //- Radiation
184 
185 
186  // Limits
187 
188  //- Minimum temperature limit (optional)
189  scalar Tmin_;
190 
191  //- Maximum temperature limit (optional)
192  scalar Tmax_;
193 
194 
195  // Protected member functions
196 
197  //- Read control parameters from dictionary
198  virtual bool read();
199 
200  //- Correct the thermo fields
201  virtual void correctThermoFields();
202 
203  //- Correct sensible enthalpy for mapped temperature fields
204  virtual void correctHsForMappedT();
205 
206  //- Correct the film surface and wall temperatures
207  virtual void updateSurfaceTemperatures();
208 
209  //- Reset source term fields
210  virtual void resetPrimaryRegionSourceTerms();
211 
212  //- Transfer thermo fields from the primary region to the film region
213  virtual void transferPrimaryRegionThermoFields();
214 
215  //- Transfer source fields from the primary region to the film region
216  virtual void transferPrimaryRegionSourceFields();
217 
218  //- Correct film coverage field
219  virtual void correctAlpha();
220 
221  //- Update the film sub-models
222  virtual void updateSubmodels();
223 
224  //- Return the wall/surface heat transfer term for the enthalpy equation
225  virtual tmp<fvScalarMatrix> q(volScalarField& hs) const;
226 
227 
228  // Equations
229 
230  //- Solve energy equation
231  virtual void solveEnergy();
232 
233 
234 public:
235 
236  //- Runtime type information
237  TypeName("thermoSingleLayer");
238 
239 
240  // Constructors
241 
242  //- Construct from components
244  (
245  const word& modelType,
246  const fvMesh& mesh,
247  const dimensionedVector& g,
248  const word& regionType,
249  const bool readFields = true
250  );
251 
252 
253  //- Destructor
254  virtual ~thermoSingleLayer();
255 
256 
257  // Member Functions
258 
259  // Thermo properties
260 
261  //- Return const reference to the SLGThermo object
262  inline const SLGThermo& thermo() const;
263 
264 
265  // Fields
266 
267  //- Return the film specific heat capacity [J/kg/K]
268  virtual const volScalarField& Cp() const;
269 
270  //- Return the film thermal conductivity [W/m/K]
271  virtual const volScalarField& kappa() const;
272 
273  //- Return the film mean temperature [K]
274  virtual const volScalarField& T() const;
275 
276  //- Return the film surface temperature [K]
277  virtual const volScalarField& Ts() const;
278 
279  //- Return the film wall temperature [K]
280  virtual const volScalarField& Tw() const;
281 
282  //- Return the film sensible enthalpy [J/kg]
283  virtual const volScalarField& hs() const;
284 
285 
286 
287  // Transfer fields - to the primary region
288 
289  //- Return mass transfer source - Eulerian phase only
290  virtual tmp<volScalarField> primaryMassTrans() const;
291 
292 
293  // Helper functions
294 
295  //- Return sensible enthalpy as a function of temperature
296  // for a patch
297  inline tmp<scalarField> hs
298  (
299  const scalarField& T,
300  const label patchi
301  ) const;
302 
303  //- Return sensible enthalpy as a function of temperature
304  inline tmp<volScalarField> hs
305  (
306  const volScalarField& T
307  ) const;
308 
309  //- Return temperature as a function of sensible enthalpy
310  inline tmp<volScalarField> T
311  (
312  const volScalarField& hs
313  ) const;
314 
315 
316  // Source fields (read/write access)
317 
318  //- External hook to add sources to the film
319  virtual void addSources
320  (
321  const label patchi, // patchi on primary region
322  const label facei, // facei of patchi
323  const scalar massSource, // [kg]
324  const vector& momentumSource, // [kg.m/s] (tangential momentum)
325  const scalar pressureSource, // [kg.m/s] (normal momentum)
326  const scalar energySource // [J]
327  );
328 
329 
330  // Source term fields
331 
332  // Film region
333 
334  //- Energy / [J/m2/s]
335  inline const volScalarField& hsSp() const;
336 
337 
338  // Primary region
339 
340  //- Energy / [J/m2/s]
341  inline const volScalarField& hsSpPrimary() const;
342 
343 
344  // Fields mapped from the primary region
345 
346  //- Temperature / [K]
347  inline const volScalarField& TPrimary() const;
348 
349  //- Specie mass fractions / [0-1]
350  inline const PtrList<volScalarField>& YPrimary() const;
351 
352 
353 
354  // Sub-models
355 
356  //- Return const access to the (surface) heat transfer model
357  inline const heatTransferModel& htcs() const;
358 
359  //- Return const access to the (wall) heat transfer model
360  inline const heatTransferModel& htcw() const;
361 
362  //- Return const access to the phase change model
363  inline const phaseChangeModel& phaseChange() const;
364 
365  //- Return const access to the radiation model
366  inline const filmRadiationModel& radiation() const;
367 
368 
369  // Derived fields (calculated on-the-fly)
370 
371  //- Return the convective heat energy from film to wall
372  inline tmp<scalarField> Qconvw(const label patchi) const;
373 
374  //- Return the convective heat energy from primary region to film
375  inline tmp<scalarField> Qconvp(const label patchi) const;
376 
377 
378  // Evolution
379 
380  //- Pre-evolve film hook
381  virtual void preEvolveRegion();
382 
383  //- Evolve the film equations
384  virtual void evolveRegion();
385 
386 
387  // Source fields
388 
389  // Mapped into primary region
390 
391  //- Return total mass source - Eulerian phase only
393 
394  //- Return mass source for specie i - Eulerian phase only
396  (
397  const label i
398  ) const;
399 
400  //- Return enthalpy source - Eulerian phase only
401  virtual tmp<DimensionedField<scalar, volMesh>> Sh() const;
402 
403 
404  // I-O
405 
406  //- Provide some feedback
407  virtual void info();
408 };
409 
410 
411 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
412 
413 } // End namespace surfaceFilmModels
414 } // End namespace regionModels
415 } // End namespace Foam
416 
417 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
418 
419 #include "thermoSingleLayerI.H"
420 
421 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
422 
423 #endif
424 
425 // ************************************************************************* //
volScalarField Ts_
Temperature - surface / [K].
volScalarField Tw_
Temperature - wall / [K].
autoPtr< heatTransferModel > htcw_
Heat transfer coefficient bewteen wall and film [W/m2/K].
virtual void correctAlpha()
Correct film coverage field.
virtual const volScalarField & Cp() const
Return the film specific heat capacity [J/kg/K].
virtual const volScalarField & kappa() const
Return the film thermal conductivity [W/m/K].
autoPtr< phaseChangeModel > phaseChange_
Phase change.
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.
scalar Tmax_
Maximum temperature limit (optional)
virtual void transferPrimaryRegionSourceFields()
Transfer source fields from the primary region to the film region.
const volScalarField & TPrimary() const
Temperature / [K].
virtual tmp< DimensionedField< scalar, volMesh > > Sh() const
Return enthalpy source - Eulerian phase only.
const dimensionedVector & g() const
Return the accleration due to gravity.
const volScalarField & hsSp() const
Energy / [J/m2/s].
virtual void transferPrimaryRegionThermoFields()
Transfer thermo fields from the primary region to the film region.
const volScalarField & hsSpPrimary() const
Energy / [J/m2/s].
autoPtr< heatTransferModel > htcs_
Heat transfer coefficient bewteen film surface and primary.
virtual tmp< DimensionedField< scalar, volMesh > > Srho() const
Return total mass source - Eulerian phase only.
virtual void updateSurfaceTemperatures()
Correct the film surface and wall temperatures.
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
volScalarField primaryEnergyPCTrans_
Film energy evolved via phase change.
virtual void updateSubmodels()
Update the film sub-models.
tmp< scalarField > Qconvp(const label patchi) const
Return the convective heat energy from primary region to film.
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
virtual tmp< fvScalarMatrix > q(volScalarField &hs) const
Return the wall/surface heat transfer term for the enthalpy equation.
volScalarField kappa_
Thermal conductivity / [W/m/K].
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.
const SLGThermo & thermo() const
Return const reference to the SLGThermo object.
scalar Tmin_
Minimum temperature limit (optional)
virtual tmp< volScalarField > primaryMassTrans() const
Return mass transfer source - Eulerian phase only.
const SLGThermo & thermo_
Reference to the SLGThermo.
virtual const volScalarField & T() const
Return the film mean temperature [K].
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.
const PtrList< volScalarField > & YPrimary() const
Specie mass fractions / [0-1].
const filmRadiationModel & radiation() const
Return const access to the radiation model.
virtual void correctThermoFields()
Correct the thermo fields.
autoPtr< filmRadiationModel > radiation_
Radiation.
virtual bool read()
Read control parameters from dictionary.
virtual const volScalarField & Tw() const
Return the film wall temperature [K].
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:62
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual void preEvolveRegion()
Pre-evolve film hook.
const heatTransferModel & htcs() const
Return const access to the (surface) heat transfer model.
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
PtrList< volScalarField > YPrimary_
List of specie mass fractions / [0-1].
const heatTransferModel & htcw() const
Return const access to the (wall) heat transfer model.
virtual const volScalarField & hs() const
Return the film sensible enthalpy [J/kg].
autoPtr< filmViscosityModel > viscosity_
Viscosity model.
volScalarField primaryMassPCTrans_
Film mass evolved via phase change.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
A class for managing temporary objects.
Definition: PtrList.H:54
virtual void correctHsForMappedT()
Correct sensible enthalpy for mapped temperature fields.
volScalarField hs_
Sensible enthalpy / [J/kg].
tmp< scalarField > Qconvw(const label patchi) const
Return the convective heat energy from film to wall.
const phaseChangeModel & phaseChange() const
Return const access to the phase change model.
Namespace for OpenFOAM.