radiationModel.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 | Copyright (C) 2011-2016 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 Namespace
25  Foam::radiation
26 
27 Description
28  Namespace for radiation modelling
29 
30 Class
31  Foam::radiation::radiationModel
32 
33 Description
34  Top level model for radiation modelling
35 
36 SourceFiles
37  radiationModel.C
38  radiationModelNew.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef radiationModel_H
43 #define radiationModel_H
44 
45 #include "IOdictionary.H"
46 #include "autoPtr.H"
47 #include "runTimeSelectionTables.H"
48 #include "volFieldsFwd.H"
49 #include "DimensionedField.H"
50 #include "fvMatricesFwd.H"
51 #include "Switch.H"
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 namespace Foam
56 {
57 
58 class fluidThermo;
59 class fvMesh;
60 
61 namespace radiation
62 {
63 
64 // Forward declaration of classes
65 class absorptionEmissionModel;
66 class scatterModel;
67 class sootModel;
68 
69 /*---------------------------------------------------------------------------*\
70  Class radiationModel Declaration
71 \*---------------------------------------------------------------------------*/
72 
73 class radiationModel
74 :
75  public IOdictionary
76 {
77 protected:
78 
79  // Protected data
80 
81  //- Reference to the mesh database
82  const fvMesh& mesh_;
83 
84  //- Reference to the time database
85  const Time& time_;
86 
87  //- Reference to the temperature field
88  const volScalarField& T_;
89 
90  //- Radiation model on/off flag
92 
93  //- Radiation model dictionary
95 
96  //- Radiation solver frequency - number flow solver iterations per
97  // radiation solver iteration
99 
100  //- Flag to enable radiation model to be evaluated on first iteration
101  bool firstIter_;
102 
103 
104  // References to the radiation sub-models
105 
106  //- Absorption/emission model
108 
109  //- Scatter model
111 
112  //- Soot model
114 
115 
116 private:
117 
118  // Private Member Functions
119 
120  //- Create IO object if dictionary is present
121  IOobject createIOobject(const fvMesh& mesh) const;
122 
123  //- Initialise
124  void initialise();
125 
126  //- Disallow default bitwise copy construct
128 
129  //- Disallow default bitwise assignment
130  void operator=(const radiationModel&);
131 
132 
133 public:
134 
135  //- Runtime type information
136  TypeName("radiationModel");
137 
138 
139  // Declare runtime constructor selection table
140 
142  (
143  autoPtr,
145  T,
146  (
147  const volScalarField& T
148  ),
149  (T)
150  );
151 
153  (
154  autoPtr,
156  dictionary,
157  (
158  const dictionary& dict,
159  const volScalarField& T
160  ),
161  (dict, T)
162  );
163 
164 
165  // Constructors
166 
167  //- Null constructor
168  radiationModel(const volScalarField& T);
169 
170  //- Construct from components
171  radiationModel(const word& type, const volScalarField& T);
172 
173  //- Construct from components
175  (
176  const word& type,
177  const dictionary& dict,
178  const volScalarField& T
179  );
180 
181 
182  // Selectors
183 
184  //- Return a reference to the selected radiation model
185  static autoPtr<radiationModel> New(const volScalarField& T);
186 
187  //- Return a reference to the selected radiation model
189  (
190  const dictionary& dict,
191  const volScalarField& T
192  );
193 
194 
195  //- Destructor
196  virtual ~radiationModel();
197 
198 
199  // Member Functions
200 
201  // Edit
202 
203  //- Main update/correction routine
204  virtual void correct();
205 
206  //- Solve radiation equation(s)
207  virtual void calculate() = 0;
208 
209  //- Read radiationProperties dictionary
210  virtual bool read() = 0;
211 
212 
213  // Access
214 
215  //- Radiation model on/off flag
216  const Switch radiation() const
217  {
218  return radiation_;
219  }
220 
221  //- Source term component (for power of T^4)
222  virtual tmp<volScalarField> Rp() const = 0;
223 
224  //- Source term component (constant)
225  virtual tmp<DimensionedField<scalar, volMesh>> Ru() const = 0;
226 
227  //- Energy source term
228  virtual tmp<fvScalarMatrix> Sh(fluidThermo& thermo) const;
229 
230  //- Temperature source term
231  virtual tmp<fvScalarMatrix> ST
232  (
233  const dimensionedScalar& rhoCp,
234  volScalarField& T
235  ) const;
236 
237  //- Access to absorptionEmission model
239 
240  //- Access to soot Model
241  const sootModel& soot() const;
242 };
243 
244 
245 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 #define addToRadiationRunTimeSelectionTables(model) \
248  \
249  addToRunTimeSelectionTable \
250  ( \
251  radiationModel, \
252  model, \
253  dictionary \
254  ); \
255  \
256  addToRunTimeSelectionTable \
257  ( \
258  radiationModel, \
259  model, \
260  T \
261  );
262 
263 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 
265 } // End namespace radiation
266 } // End namespace Foam
267 
268 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
269 
270 #endif
271 
272 // ************************************************************************* //
dictionary dict
dictionary coeffs_
Radiation model dictionary.
const sootModel & soot() const
Access to soot Model.
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
autoPtr< scatterModel > scatter_
Scatter model.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const volScalarField & T_
Reference to the temperature field.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
Definition: Switch.H:60
TypeName("radiationModel")
Runtime type information.
virtual void calculate()=0
Solve radiation equation(s)
const absorptionEmissionModel & absorptionEmission() const
Access to absorptionEmission model.
Model to supply absorption and emission coefficients for radiation modelling.
declareRunTimeSelectionTable(autoPtr, radiationModel, T,(const volScalarField &T),(T))
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
virtual tmp< fvScalarMatrix > Sh(fluidThermo &thermo) const
Energy source term.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
static autoPtr< radiationModel > New(const volScalarField &T)
Return a reference to the selected radiation model.
virtual ~radiationModel()
Destructor.
dynamicFvMesh & mesh
Base class for soor models.
Definition: sootModel.H:53
virtual bool read()=0
Read radiationProperties dictionary.
autoPtr< sootModel > soot_
Soot model.
A class for handling words, derived from string.
Definition: word.H:59
bool firstIter_
Flag to enable radiation model to be evaluated on first iteration.
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:49
const fvMesh & mesh_
Reference to the mesh database.
const Switch radiation() const
Radiation model on/off flag.
Top level model for radiation modelling.
virtual void correct()
Main update/correction routine.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual tmp< volScalarField > Rp() const =0
Source term component (for power of T^4)
virtual tmp< fvScalarMatrix > ST(const dimensionedScalar &rhoCp, volScalarField &T) const
Temperature source term.
Forward declarations of fvMatrix specializations.
const Time & time_
Reference to the time database.
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const =0
Source term component (constant)
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
Switch radiation_
Radiation model on/off flag.
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: PtrList.H:54
label solverFreq_
Radiation solver frequency - number flow solver iterations per.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
autoPtr< absorptionEmissionModel > absorptionEmission_
Absorption/emission model.
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.