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 
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 
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 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::Tw_
volScalarField Tw_
Temperature - wall [K].
Definition: thermoSingleLayer.H:108
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hs
virtual const volScalarField & hs() const
Return the film sensible enthalpy [J/kg].
Definition: thermoSingleLayer.C:680
Foam::regionModels::surfaceFilmModels::heatTransferModel
Base class for film heat transfer models.
Definition: heatTransferModel.H:56
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::YPrimary
const PtrList< volScalarField > & YPrimary() const
Specie mass fractions [0-1].
Definition: thermoSingleLayerI.H:130
Foam::SLGThermo
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:64
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::htcs_
autoPtr< heatTransferModel > htcs_
Heat transfer coefficient between film surface and primary.
Definition: thermoSingleLayer.H:172
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::htcs
const heatTransferModel & htcs() const
Return const access to the (surface) heat transfer model.
Definition: thermoSingleLayerI.H:136
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::thermo
const SLGThermo & thermo() const
Return const reference to the SLGThermo object.
Definition: thermoSingleLayerI.H:44
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::phaseChange_
autoPtr< phaseChangeModel > phaseChange_
Phase change.
Definition: thermoSingleLayer.H:178
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::T_
volScalarField T_
Temperature - mean [K].
Definition: thermoSingleLayer.H:102
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::qconvp
tmp< scalarField > qconvp(const label patchi) const
Return the convective heat energy from primary region to film.
Definition: thermoSingleLayerI.H:170
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::radiation
const filmRadiationModel & radiation() const
Return const access to the radiation model.
Definition: thermoSingleLayerI.H:154
thermoSingleLayer.H
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hsSp_
volScalarField hsSp_
Energy [J/m2/s].
Definition: thermoSingleLayer.H:145
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::Field< scalar >
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::thermo_
const SLGThermo & thermo_
Reference to the SLGThermo.
Definition: thermoSingleLayer.H:90
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::T
virtual const volScalarField & T() const
Return the film mean temperature [K].
Definition: thermoSingleLayer.C:662
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::YPrimary_
PtrList< volScalarField > YPrimary_
List of specie mass fractions [0-1].
Definition: thermoSingleLayer.H:162
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::TPrimary_
volScalarField TPrimary_
Temperature [K].
Definition: thermoSingleLayer.H:159
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hsSp
const volScalarField & hsSp() const
Energy [J/m2/s].
Definition: thermoSingleLayerI.H:112
filmRadiationModel.H
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hsSpPrimary_
volScalarField hsSpPrimary_
Energy [J/m2/s].
Definition: thermoSingleLayer.H:152
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hsSpPrimary
const volScalarField & hsSpPrimary() const
Energy [J/m2/s].
Definition: thermoSingleLayerI.H:118
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::phaseChange
const phaseChangeModel & phaseChange() const
Return const access to the phase change model.
Definition: thermoSingleLayerI.H:148
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::htcw
const heatTransferModel & htcw() const
Return const access to the (wall) heat transfer model.
Definition: thermoSingleLayerI.H:142
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::qconvw
tmp< scalarField > qconvw(const label patchi) const
Return the convective heat energy from film to wall.
Definition: thermoSingleLayerI.H:160
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::radiation_
autoPtr< filmRadiationModel > radiation_
Radiation.
Definition: thermoSingleLayer.H:181
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::TPrimary
const volScalarField & TPrimary() const
Temperature [K].
Definition: thermoSingleLayerI.H:124
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::htcw_
autoPtr< heatTransferModel > htcw_
Heat transfer coefficient between wall and film [W/m2/K].
Definition: thermoSingleLayer.H:175
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::regionModels::surfaceFilmModels::filmRadiationModel
Base class for film radiation models.
Definition: filmRadiationModel.H:56
phaseChangeModel.H
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::regionModels::surfaceFilmModels::phaseChangeModel
Base class for surface film phase change models.
Definition: phaseChangeModel.H:57