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