thermalBaffleModel.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-2013 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 Class
27  Foam::regionModels::thermalBaffleModels::thermalBaffleModel
28 
29 Description
30 
31 SourceFiles
32  thermalBaffleModel.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef thermalBaffleModel_H
37 #define thermalBaffleModel_H
38 
39 #include "runTimeSelectionTables.H"
40 #include "scalarIOField.H"
41 #include "autoPtr.H"
42 #include "volFieldsFwd.H"
43 #include "solidThermo.H"
44 #include "regionModel1D.H"
45 #include "radiationModel.H"
46 
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 namespace regionModels
53 {
54 namespace thermalBaffleModels
55 {
56 
57 /*---------------------------------------------------------------------------*\
58  Class thermalBaffleModel Declaration
59 \*---------------------------------------------------------------------------*/
60 
62 :
63  public regionModel1D
64 {
65 private:
66 
67  // Private Member Functions
68 
69  //- No copy construct
70  thermalBaffleModel(const thermalBaffleModel&) = delete;
71 
72  //- No copy assignment
73  void operator=(const thermalBaffleModel&) = delete;
74 
75  //- Initialize thermal Baffle
76  void init();
77 
78 
79 protected:
80 
81  // Protected Data
82 
83  //- Baffle physical thickness
85 
86  //- Baffle mesh thickness
88 
89  //- Is it one dimension
90  bool oneD_;
91 
92  //- Is thickness constant
93  bool constantThickness_;
94 
95 
96  // Protected Member Functions
97 
98  //- Read control parameters from IO dictionary
99  virtual bool read();
100 
101  //- Read control parameters from dictionary
102  virtual bool read(const dictionary&);
103 
104 
105 public:
106 
107  //- Runtime type information
108  TypeName("thermalBaffleModel");
109 
110 
111  // Declare runtime constructor selection tables
112 
114  (
115  autoPtr,
117  mesh,
118  (
119  const word& modelType,
120  const fvMesh& mesh
121  ),
122  (modelType, mesh)
123  );
124 
126  (
127  autoPtr,
129  dictionary,
130  (
131  const word& modelType,
132  const fvMesh& mesh,
133  const dictionary& dict
134  ),
135  (modelType, mesh, dict)
136  );
137 
138 
139  // Constructors
140 
141  //- Construct null from mesh
143 
144  //- Construct from type name and mesh
145  thermalBaffleModel(const word& modelType, const fvMesh& mesh);
146 
147  //- Construct from type name and mesh and dict
149  (
150  const word& modelType,
151  const fvMesh& mesh,
152  const dictionary& dict
153  );
154 
155 
156  // Selectors
157 
158  //- Return a reference to the selected model
160 
161  //- Return a reference to the selected model using dictionary
163  (
164  const fvMesh& mesh,
165  const dictionary& dict
166  );
167 
168 
169  //- Destructor
170  virtual ~thermalBaffleModel();
171 
172 
173  // Member Functions
174 
175  // Access
176 
177  //- Return solid thermo
178  virtual const solidThermo& thermo() const = 0;
179 
180  //- Return thickness
181  const scalarField& thickness() const
182  {
183  return thickness_;
184  }
185 
186  //- Return geometrical thickness
187  const dimensionedScalar& delta() const
188  {
189  return delta_;
190  }
191 
192  //- Return if region is one dimensional
193  bool oneD() const
194  {
195  return oneD_;
196  }
197 
198  //- Return if region has constant thickness
199  bool constantThickness() const
200  {
201  return constantThickness_;
202  }
203 
204 
205  // Fields
206 
207  //- Return density [kg/m3]
208  virtual const volScalarField& rho() const = 0;
209 
210  //- Return const temperature [K]
211  virtual const volScalarField& T() const = 0;
212 
213  //- Return specific heat capacity [J/kg/K]
214  virtual const tmp<volScalarField> Cp() const = 0;
215 
216  //- Return the region absorptivity [1/m]
217  virtual const volScalarField& kappaRad() const = 0;
218 
219  //- Return the region thermal conductivity [W/m/k]
220  virtual const volScalarField& kappa() const = 0;
221 
222 
223  // Evolution
224 
225  //- Pre-evolve region
226  virtual void preEvolveRegion();
227 };
228 
229 
230 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
231 
232 } // End namespace thermalBaffleModels
233 } // End namespace regionModels
234 } // End namespace Foam
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 #endif
239 
240 // ************************************************************************* //
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::oneD_
bool oneD_
Is it one dimension.
Definition: thermalBaffleModel.H:89
volFieldsFwd.H
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::kappa
virtual const volScalarField & kappa() const =0
Return the region thermal conductivity [W/m/k].
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::Cp
virtual const tmp< volScalarField > Cp() const =0
Return specific heat capacity [J/kg/K].
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::oneD
bool oneD() const
Return if region is one dimensional.
Definition: thermalBaffleModel.H:192
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::delta_
dimensionedScalar delta_
Baffle mesh thickness.
Definition: thermalBaffleModel.H:86
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::constantThickness_
bool constantThickness_
Is thickness constant.
Definition: thermalBaffleModel.H:92
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::thermo
virtual const solidThermo & thermo() const =0
Return solid thermo.
scalarIOField.H
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::delta
const dimensionedScalar & delta() const
Return geometrical thickness.
Definition: thermalBaffleModel.H:186
regionModel1D.H
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::thickness_
scalarField thickness_
Baffle physical thickness.
Definition: thermalBaffleModel.H:83
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::New
static autoPtr< thermalBaffleModel > New(const fvMesh &mesh)
Return a reference to the selected model.
Definition: thermalBaffleModelNew.C:42
solidThermo.H
Foam::solidThermo
Fundamental solid thermodynamic properties.
Definition: solidThermo.H:52
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::thickness
const scalarField & thickness() const
Return thickness.
Definition: thermalBaffleModel.H:180
Foam::regionModels::thermalBaffleModels::thermalBaffleModel
Definition: thermalBaffleModel.H:60
Foam::Field< scalar >
init
mesh init(true)
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::~thermalBaffleModel
virtual ~thermalBaffleModel()
Destructor.
Definition: thermalBaffleModel.C:231
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, thermalBaffleModel, mesh,(const word &modelType, const fvMesh &mesh),(modelType, mesh))
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned< scalar >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::TypeName
TypeName("thermalBaffleModel")
Runtime type information.
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::T
virtual const volScalarField & T() const =0
Return const temperature [K].
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::kappaRad
virtual const volScalarField & kappaRad() const =0
Return the region absorptivity [1/m].
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
runTimeSelectionTables.H
Macros to ease declaration of run-time selection tables.
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::constantThickness
bool constantThickness() const
Return if region has constant thickness.
Definition: thermalBaffleModel.H:198
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve region.
Definition: thermalBaffleModel.C:237
Foam::regionModels::regionModel1D
Base class for 1-D region models.
Definition: regionModel1D.H:54
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::rho
virtual const volScalarField & rho() const =0
Return density [kg/m3].
radiationModel.H
Foam::regionModels::thermalBaffleModels::thermalBaffleModel::read
virtual bool read()
Read control parameters from IO dictionary.
Definition: thermalBaffleModel.C:52
autoPtr.H