thermoSingleLayerI.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 -------------------------------------------------------------------------------
10 License
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 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "thermoSingleLayer.H"
29 #include "filmRadiationModel.H"
30 #include "heatTransferModel.H"
31 #include "phaseChangeModel.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace regionModels
38 {
39 namespace surfaceFilmModels
40 {
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 inline const SLGThermo& thermoSingleLayer::thermo() const
45 {
46  return thermo_;
47 }
48 
49 
51 (
52  const scalarField& T,
53  const label patchi
54 ) const
55 {
56  const scalarField& Cp = Cp_.boundaryField()[patchi];
57  return Cp*(T - Tref.value());
58 }
59 
60 
62 (
63  const volScalarField& T
64 ) const
65 {
66  return tmp<volScalarField>
67  (
68  new volScalarField
69  (
70  IOobject
71  (
72  "hs(" + T.name() + ")",
73  time().timeName(),
74  regionMesh(),
77  ),
78  Cp_*(T - Tref)
79  )
80  );
81 }
82 
83 
85 (
86  const volScalarField& hs
87 ) const
88 {
90  (
91  new volScalarField
92  (
93  IOobject
94  (
95  "T(" + hs.name() + ")",
96  time().timeName(),
97  regionMesh(),
100  ),
101  hs/Cp_ + Tref
102  )
103  );
104 
105  tT.ref().min(Tmax_);
106  tT.ref().max(Tmin_);
107 
108  return tT;
109 }
110 
112 inline const volScalarField& thermoSingleLayer::hsSp() const
113 {
114  return hsSp_;
115 }
116 
119 {
120  return hsSpPrimary_;
121 }
122 
124 inline const volScalarField& thermoSingleLayer::TPrimary() const
125 {
126  return TPrimary_;
127 }
128 
131 {
132  return YPrimary_;
133 }
134 
136 inline const heatTransferModel& thermoSingleLayer::htcs() const
137 {
138  return *htcs_;
139 }
140 
142 inline const heatTransferModel& thermoSingleLayer::htcw() const
143 {
144  return *htcw_;
145 }
146 
149 {
150  return *phaseChange_;
151 }
152 
155 {
156  return *radiation_;
157 }
158 
159 
160 inline tmp<scalarField> thermoSingleLayer::qconvw(const label patchi) const
161 {
162  const scalarField htc(htcw_->h()().boundaryField()[patchi]);
163  const scalarField& Tp = T_.boundaryField()[patchi];
164  const scalarField& Twp = Tw_.boundaryField()[patchi];
165 
166  return htc*(Tp - Twp);
167 }
168 
169 
170 inline tmp<scalarField> thermoSingleLayer::qconvp(const label patchi) const
171 {
172  const scalarField htc(htcs_->h()().boundaryField()[patchi]);
173  const scalarField& Tp = T_.boundaryField()[patchi];
174  const scalarField& Tpp = TPrimary_.boundaryField()[patchi];
175 
176  return htc*(Tp - Tpp);
177 }
178 
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 } // End namespace surfaceFilmModels
183 } // End namespace regionModels
184 } // End namespace Foam
185 
186 // ************************************************************************* //
autoPtr< heatTransferModel > htcw_
Heat transfer coefficient between wall and film [W/m2/K].
const heatTransferModel & htcs() const
Return const access to the (surface) heat transfer model.
const volScalarField & hsSpPrimary() const
Energy [J/m2/s].
const Type & value() const noexcept
Return const reference to value.
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].
scalar Tmax_
Maximum temperature limit (optional)
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:150
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:210
autoPtr< heatTransferModel > htcs_
Heat transfer coefficient between film surface and primary.
Ignore writing from objectRegistry::writeObject()
tmp< scalarField > qconvp(const label patchi) const
Return the convective heat energy from primary region to film.
Base class for surface film phase change models.
const volScalarField & TPrimary() const
Temperature [K].
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
tmp< scalarField > qconvw(const label patchi) const
Return the convective heat energy from film to wall.
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:26
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:60
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.
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
autoPtr< filmRadiationModel > radiation_
Radiation.
static const dimensionedScalar Tref
Reference temperature for enthalpy.
const heatTransferModel & htcw() const
Return const access to the (wall) heat transfer model.
volScalarField Cp_
Specific heat capacity [J/kg/K].
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Nothing to be read.
PtrList< volScalarField > YPrimary_
List of specie mass fractions [0-1].
virtual const volScalarField & hs() const
Return the film sensible enthalpy [J/kg].
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const filmRadiationModel & radiation() const
Return const access to the radiation model.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
const SLGThermo & thermo() const
Return const reference to the SLGThermo object.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
const Time & time() const noexcept
Return the reference to the time database.
Definition: regionModel.H:244
Namespace for OpenFOAM.
Base class for film heat transfer models.
virtual const volScalarField & T() const
Return the film mean temperature [K].