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-------------------------------------------------------------------------------
10License
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
35namespace Foam
36{
37namespace regionModels
38{
39namespace surfaceFilmModels
40{
41
42// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
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{
67 (
69 (
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 (
92 (
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
111
113{
114 return hsSp_;
115}
116
117
119{
120 return hsSpPrimary_;
121}
122
123
125{
126 return TPrimary_;
127}
128
129
131{
132 return YPrimary_;
133}
134
135
137{
138 return *htcs_;
139}
140
141
143{
144 return *htcw_;
145}
146
147
149{
150 return *phaseChange_;
151}
152
153
155{
156 return *radiation_;
157}
158
159
160inline 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
170inline 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// ************************************************************************* //
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:67
const Type & value() const
Return const reference to value.
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:40
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
static const dimensionedScalar Tref
Reference temperature for enthalpy.
Base class for film heat transfer models.
Base class for surface film phase change models.
autoPtr< heatTransferModel > htcw_
Heat transfer coefficient between wall and film [W/m2/K].
virtual const volScalarField & hs() const
Return the film sensible enthalpy [J/kg].
autoPtr< filmRadiationModel > radiation_
Radiation.
const SLGThermo & thermo() const
Return const reference to the SLGThermo object.
tmp< scalarField > qconvp(const label patchi) const
Return the convective heat energy from primary region to film.
virtual const volScalarField & T() const
Return the film mean temperature [K].
volScalarField Cp_
Specific heat capacity [J/kg/K].
const PtrList< volScalarField > & YPrimary() const
Specie mass fractions [0-1].
const heatTransferModel & htcw() const
Return const access to the (wall) heat transfer model.
const volScalarField & hsSp() const
Energy [J/m2/s].
const phaseChangeModel & phaseChange() const
Return const access to the phase change model.
virtual const volScalarField & Cp() const
Return the film specific heat capacity [J/kg/K].
PtrList< volScalarField > YPrimary_
List of specie mass fractions [0-1].
const heatTransferModel & htcs() const
Return const access to the (surface) heat transfer model.
const filmRadiationModel & radiation() const
Return const access to the radiation model.
const volScalarField & TPrimary() const
Temperature [K].
scalar Tmax_
Maximum temperature limit (optional)
scalar Tmin_
Minimum temperature limit (optional)
const SLGThermo & thermo_
Reference to the SLGThermo.
const volScalarField & hsSpPrimary() const
Energy [J/m2/s].
autoPtr< phaseChangeModel > phaseChange_
Phase change.
tmp< scalarField > qconvw(const label patchi) const
Return the convective heat energy from film to wall.
autoPtr< heatTransferModel > htcs_
Heat transfer coefficient between film surface and primary.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
const volScalarField & T
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.