thermalShell.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) 2019-2020 OpenCFD Ltd.
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 Class
27  Foam::regionFaModels::thermalShell
28 
29 Description
30  Thermal-shell finite-area model. It solves the energy
31  equation in 2D. The coupling with the 3D region is done through
32  the \c temperatureCoupledBase, plus \c faOption is available to
33  add extra sources on the shell such as \c externalHeatFluxSource etc.
34 
35 Usage
36  Example of the boundary condition specification:
37  \verbatim
38  <patchName>
39  {
40  // Mandatory entries
41  thermalShellModel thermalShell;
42  thermo
43  {
44  // subdictionary entries
45  }
46 
47  // Optional entries
48  qr <word>;
49  thickness <scalar>;
50 
51  // Inherited entries
52  ...
53  nNonOrthCorr <int>; // read from another dict
54  }
55  \endverbatim
56 
57  where the entries mean:
58  \table
59  Property | Description | Type | Reqd | Deflt
60  thermalShellModel | Type name: thermalShell | word | yes | -
61  thermo | Solid thermal properties | dictionary | yes | -
62  qr | Name of radiative heat flux field | word | no | none
63  thickness | Uniform film thickness [m] | scalar | choice | -
64  \endtable
65 
66  The inherited entries are elaborated in:
67  - \link thermalShellModel.H \endlink
68 
69 SourceFiles
70  thermalShell.C
71  thermalShellI.H
72 
73 \*---------------------------------------------------------------------------*/
74 
75 #ifndef Foam_regionModels_thermalShell_H
76 #define Foam_regionModels_thermalShell_H
77 
78 #include "volFieldsFwd.H"
79 #include "thermalShellModel.H"
80 #include "solidProperties.H"
81 #include "faCFD.H"
82 
83 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 
85 namespace Foam
86 {
87 namespace regionModels
88 {
89 
90 /*---------------------------------------------------------------------------*\
91  Class thermalShell Declaration
92 \*---------------------------------------------------------------------------*/
93 
94 class thermalShell
95 :
96  public thermalShellModel
97 {
98  // Private Member Functions
99 
100  //- Initialize thermalShell
101  bool init(const dictionary& dict);
102 
103  //- Return radiative heat flux
105 
106 
107 protected:
108 
109  // Protected Data
110 
111  // Solution parameters
112 
113  //- Number of non orthogonal correctors
114  label nNonOrthCorr_;
115 
116 
117  // Thermo properties
118 
119  //- Solid properties
121 
122 
123  // Source term fields
124 
125  //- External surface energy source [J/m2/s]
127 
128  //- Film thickness [m]
130 
131  //- Name of the primary region radiative flux
132  const word qrName_;
133 
134  //- Uniform film thickness [m]
135  scalar thickness_;
136 
137 
138  // Protected Member Functions
139 
140  // Equations
141 
142  //- Solve energy equation
143  void solveEnergy();
144 
146 public:
147 
148  //- Runtime type information
149  TypeName("thermalShell");
150 
151 
152  // Constructors
154  //- Construct from components and dict
156  (
157  const word& modelType,
158  const fvMesh& mesh,
159  const dictionary& dict
160  );
162  //- No copy construct
163  thermalShell(const thermalShell&) = delete;
164 
165  //- No copy assignment
166  void operator=(const thermalShell&) = delete;
167 
168 
169  //- Destructor
170  virtual ~thermalShell() = default;
172 
173  // Member Functions
174 
175  // Fields
177  //- Return the film specific heat capacity [J/kg/K]
178  const tmp<areaScalarField> Cp() const;
179 
180  //- Return density [kg/m3]
181  const tmp<areaScalarField> rho() const;
182 
183  //- Return thermal conductivity [W/m/K]
184  const tmp<areaScalarField> kappa() const;
185 
186 
187  // Evolution
188 
189  //- Pre-evolve thermal baffle
190  virtual void preEvolveRegion();
191 
192  //- Evolve the thermal baffle
193  virtual void evolveRegion();
194 
195 
196  // IO
197 
198  //- Provide some feedback
199  virtual void info();
200 };
201 
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 } // End namespace regionModels
206 } // End namespace Foam
207 
208 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 
210 
211 #endif
212 
213 // ************************************************************************* //
label nNonOrthCorr_
Number of non orthogonal correctors.
Definition: thermalShell.H:145
dictionary dict
Forwards and collection of common volume field types.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void info()
Provide some feedback.
Definition: thermalShell.C:215
TypeName("thermalShell")
Runtime type information.
const word qrName_
Name of the primary region radiative flux.
Definition: thermalShell.H:171
const tmp< areaScalarField > rho() const
Return density [kg/m3].
Definition: thermalShell.C:185
void operator=(const thermalShell &)=delete
No copy assignment.
areaScalarField qs_
External surface energy source [J/m2/s].
Definition: thermalShell.H:161
solidProperties thermo_
Solid properties.
Definition: thermalShell.H:153
const tmp< areaScalarField > Cp() const
Return the film specific heat capacity [J/kg/K].
Definition: thermalShell.C:172
virtual void preEvolveRegion()
Pre-evolve thermal baffle.
Definition: thermalShell.C:155
dynamicFvMesh & mesh
scalar thickness_
Uniform film thickness [m].
Definition: thermalShell.H:176
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual void evolveRegion()
Evolve the thermal baffle.
Definition: thermalShell.C:159
The thermophysical properties of a solid.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void solveEnergy()
Solve energy equation.
Definition: thermalShell.C:83
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual ~thermalShell()=default
Destructor.
thermalShell(const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components and dict.
Definition: thermalShell.C:112
const tmp< areaScalarField > kappa() const
Return thermal conductivity [W/m/K].
Definition: thermalShell.C:198
areaScalarField h_
Film thickness [m].
Definition: thermalShell.H:166
Namespace for OpenFOAM.