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  Copyright (C) 2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "thermoSingleLayer.H"
30 #include "filmRadiationModel.H"
31 #include "heatTransferModel.H"
32 #include "phaseChangeModel.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace regionModels
39 {
40 namespace surfaceFilmModels
41 {
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 inline const SLGThermo& thermoSingleLayer::thermo() const
46 {
47  return thermo_;
48 }
49 
50 
52 (
53  const scalarField& T,
54  const label patchi
55 ) const
56 {
57  const scalarField& Cp = Cp_.boundaryField()[patchi];
58  return Cp*(T - Tref.value());
59 }
60 
61 
63 (
64  const volScalarField& T
65 ) const
66 {
67  return tmp<volScalarField>
68  (
69  new volScalarField
70  (
71  IOobject
72  (
73  "hs(" + T.name() + ")",
74  regionMesh().time().timeName(),
75  regionMesh().thisDb(),
78  ),
79  Cp_*(T - Tref)
80  )
81  );
82 }
83 
84 
86 (
87  const volScalarField& hs
88 ) const
89 {
91  (
92  new volScalarField
93  (
94  IOobject
95  (
96  "T(" + hs.name() + ")",
97  regionMesh().time().timeName(),
98  regionMesh().thisDb(),
101  ),
102  hs/Cp_ + Tref
103  )
104  );
105 
106  if (limitType::CLAMP_MIN == withTbounds_)
107  {
108  tT.ref().clamp_min(Tbounds_.min());
109  }
110  else if (limitType::CLAMP_MAX == withTbounds_)
111  {
112  tT.ref().clamp_max(Tbounds_.max());
113  }
114  else if (limitType::CLAMP_RANGE == withTbounds_)
115  {
116  tT.ref().clamp_range(Tbounds_);
117  }
118 
119  return tT;
120 }
121 
123 inline const volScalarField& thermoSingleLayer::hsSp() const
124 {
125  return hsSp_;
126 }
127 
130 {
131  return hsSpPrimary_;
132 }
133 
135 inline const volScalarField& thermoSingleLayer::TPrimary() const
136 {
137  return TPrimary_;
138 }
139 
142 {
143  return YPrimary_;
144 }
145 
147 inline const heatTransferModel& thermoSingleLayer::htcs() const
148 {
149  return *htcs_;
150 }
151 
153 inline const heatTransferModel& thermoSingleLayer::htcw() const
154 {
155  return *htcw_;
156 }
157 
160 {
161  return *phaseChange_;
162 }
163 
166 {
167  return *radiation_;
168 }
169 
170 
171 inline tmp<scalarField> thermoSingleLayer::qconvw(const label patchi) const
172 {
173  const scalarField htc(htcw_->h()().boundaryField()[patchi]);
174  const scalarField& Tp = T_.boundaryField()[patchi];
175  const scalarField& Twp = Tw_.boundaryField()[patchi];
176 
177  return htc*(Tp - Twp);
178 }
179 
180 
181 inline tmp<scalarField> thermoSingleLayer::qconvp(const label patchi) const
182 {
183  const scalarField htc(htcs_->h()().boundaryField()[patchi]);
184  const scalarField& Tp = T_.boundaryField()[patchi];
185  const scalarField& Tpp = TPrimary_.boundaryField()[patchi];
186 
187  return htc*(Tp - Tpp);
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 } // End namespace surfaceFilmModels
194 } // End namespace regionModels
195 } // End namespace Foam
196 
197 // ************************************************************************* //
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].
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
autoPtr< heatTransferModel > htcs_
Heat transfer coefficient between film surface and primary.
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
tmp< scalarField > qconvp(const label patchi) const
Return the convective heat energy from primary region to film.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
Base class for surface film phase change models.
const volScalarField & TPrimary() const
Temperature [K].
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
const T & min() const noexcept
The min value (first)
Definition: MinMaxI.H:107
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
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 a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
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].
const T & max() const noexcept
The max value (second)
Definition: MinMaxI.H:121
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:172
scalarMinMax Tbounds_
Temperature limits (optional)
const SLGThermo & thermo() const
Return const reference to the SLGThermo object.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Namespace for OpenFOAM.
Base class for film heat transfer models.
virtual const volScalarField & T() const
Return the film mean temperature [K].