thermoSingleLayerI.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-2018 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 \*---------------------------------------------------------------------------*/
25 
26 #include "thermoSingleLayer.H"
27 #include "filmRadiationModel.H"
28 #include "heatTransferModel.H"
29 #include "phaseChangeModel.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace regionModels
36 {
37 namespace surfaceFilmModels
38 {
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 inline const SLGThermo& thermoSingleLayer::thermo() const
43 {
44  return thermo_;
45 }
46 
47 
49 (
50  const scalarField& T,
51  const label patchi
52 ) const
53 {
55  return Cp*(T - Tref.value());
56 }
57 
58 
60 (
61  const volScalarField& T
62 ) const
63 {
64  return volScalarField::New
65  (
66  "hs(" + T.name() + ")",
67  Cp_*(T - Tref)
68  );
69 }
70 
71 
73 (
74  const volScalarField& hs
75 ) const
76 {
78  (
79  volScalarField::New("T(" + hs.name() + ")", hs/Cp_ + Tref)
80  );
81 
82  tT.ref().min(Tmax_);
83  tT.ref().max(Tmin_);
84 
85  return tT;
86 }
87 
88 
90 {
91  return hsSp_;
92 }
93 
94 
96 {
97  return hsSpPrimary_;
98 }
99 
100 
102 {
103  return TPrimary_;
104 }
105 
106 
108 {
109  return YPrimary_;
110 }
111 
112 
114 {
115  return htcs_();
116 }
117 
118 
120 {
121  return htcw_();
122 }
123 
124 
126 {
127  return phaseChange_();
128 }
129 
130 
132 {
133  return radiation_();
134 }
135 
136 
138 {
139  const scalarField htc(htcw_->h()().boundaryField()[patchi]);
140  const scalarField& Tp = T_.boundaryField()[patchi];
141  const scalarField& Twp = Tw_.boundaryField()[patchi];
142 
143  return htc*(Tp - Twp);
144 }
145 
146 
148 {
149  const scalarField htc(htcs_->h()().boundaryField()[patchi]);
150  const scalarField& Tp = T_.boundaryField()[patchi];
151  const scalarField& Tpp = TPrimary_.boundaryField()[patchi];
152 
153  return htc*(Tp - Tpp);
154 }
155 
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 } // End namespace surfaceFilmModels
160 } // End namespace regionModels
161 } // End namespace Foam
162 
163 // ************************************************************************* //
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.
const volScalarField & hsSpPrimary() const
Energy [J/m2/s].
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].
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
const word & name() const
Return name.
Definition: IOobject.H:295
scalar Tmax_
Maximum temperature limit (optional)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
tmp< scalarField > Qconvp(const label patchi) const
Return the convective heat energy from primary region to film.
autoPtr< heatTransferModel > htcs_
Heat transfer coefficient between film surface and primary.
Base class for surface film phase change models.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
const volScalarField & TPrimary() const
Temperature [K].
const Type & value() const
Return const reference to value.
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:62
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.
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.
label patchi
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:70
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: PtrList.H:53
const filmRadiationModel & radiation() const
Return const access to the radiation model.
const SLGThermo & thermo() const
Return const reference to the SLGThermo object.
tmp< scalarField > Qconvw(const label patchi) const
Return the convective heat energy from film to wall.
Namespace for OpenFOAM.
Base class for film heat transfer models.
virtual const volScalarField & T() const
Return the film mean temperature [K].