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 volScalarField::New
68  (
69  "hs(" + T.name() + ")",
71  Cp_*(T - Tref)
72  );
73 }
74 
75 
77 (
78  const volScalarField& hs
79 ) const
80 {
81  auto tT = volScalarField::New
82  (
83  "T(" + hs.name() + ")",
85  hs/Cp_ + Tref
86  );
87 
88  if (limitType::CLAMP_MIN == withTbounds_)
89  {
90  tT.ref().clamp_min(Tbounds_.min());
91  }
92  else if (limitType::CLAMP_MAX == withTbounds_)
93  {
94  tT.ref().clamp_max(Tbounds_.max());
95  }
96  else if (limitType::CLAMP_RANGE == withTbounds_)
97  {
98  tT.ref().clamp_range(Tbounds_);
99  }
100 
101  return tT;
102 }
103 
105 inline const volScalarField& thermoSingleLayer::hsSp() const
106 {
107  return hsSp_;
108 }
109 
112 {
113  return hsSpPrimary_;
114 }
115 
117 inline const volScalarField& thermoSingleLayer::TPrimary() const
118 {
119  return TPrimary_;
120 }
121 
124 {
125  return YPrimary_;
126 }
127 
129 inline const heatTransferModel& thermoSingleLayer::htcs() const
130 {
131  return *htcs_;
132 }
133 
135 inline const heatTransferModel& thermoSingleLayer::htcw() const
136 {
137  return *htcw_;
138 }
139 
142 {
143  return *phaseChange_;
144 }
145 
148 {
149  return *radiation_;
150 }
151 
152 
153 inline tmp<scalarField> thermoSingleLayer::qconvw(const label patchi) const
154 {
155  const scalarField htc(htcw_->h()().boundaryField()[patchi]);
156  const scalarField& Tp = T_.boundaryField()[patchi];
157  const scalarField& Twp = Tw_.boundaryField()[patchi];
158 
159  return htc*(Tp - Twp);
160 }
161 
162 
163 inline tmp<scalarField> thermoSingleLayer::qconvp(const label patchi) const
164 {
165  const scalarField htc(htcs_->h()().boundaryField()[patchi]);
166  const scalarField& Tp = T_.boundaryField()[patchi];
167  const scalarField& Tpp = TPrimary_.boundaryField()[patchi];
168 
169  return htc*(Tp - Tpp);
170 }
171 
172 
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 
175 } // End namespace surfaceFilmModels
176 } // End namespace regionModels
177 } // End namespace Foam
178 
179 // ************************************************************************* //
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
autoPtr< heatTransferModel > htcs_
Heat transfer coefficient between film surface and primary.
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].
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 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 tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
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
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.
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.
Do not request registration (bool: false)
Namespace for OpenFOAM.
Base class for film heat transfer models.
virtual const volScalarField & T() const
Return the film mean temperature [K].