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-2018 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 "volFields.H"
49 #include "fvMatricesFwd.H"
50 #include "Switch.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 
57 class basicThermo;
58 class fvMesh;
59 
60 namespace radiation
61 {
62 
63 // Forward declaration of classes
64 class absorptionEmissionModel;
65 class scatterModel;
66 class sootModel;
67 
68 /*---------------------------------------------------------------------------*\
69  Class radiationModel Declaration
70 \*---------------------------------------------------------------------------*/
71 
72 class radiationModel
73 :
74  public IOdictionary
75 {
76 protected:
77 
78  // Protected data
79 
80  //- Reference to the mesh database
81  const fvMesh& mesh_;
82 
83  //- Reference to the time database
84  const Time& time_;
85 
86  //- Reference to the temperature field
87  const volScalarField& T_;
88 
89  //- Radiation model dictionary
91 
92  //- Radiation solver frequency - number flow solver iterations per
93  // radiation solver iteration
95 
96  //- Flag to enable radiation model to be evaluated on first iteration
97  bool firstIter_;
98 
99 
100  // References to the radiation sub-models
101 
102  //- Absorption/emission model
104 
105  //- Scatter model
107 
108  //- Soot model
110 
111 
112 private:
113 
114  // Private Member Functions
115 
116  //- Create IO object if dictionary is present
117  IOobject createIOobject(const fvMesh& mesh) const;
118 
119  //- Initialise
120  void initialise();
121 
122  //- Disallow default bitwise copy construct
124 
125  //- Disallow default bitwise assignment
126  void operator=(const radiationModel&);
127 
128 
129 public:
130 
131  //- Runtime type information
132  TypeName("radiationModel");
133 
134 
135  // Declare runtime constructor selection table
136 
138  (
139  autoPtr,
141  T,
142  (
143  const volScalarField& T
144  ),
145  (T)
146  );
147 
149  (
150  autoPtr,
152  dictionary,
153  (
154  const dictionary& dict,
155  const volScalarField& T
156  ),
157  (dict, T)
158  );
159 
160 
161  // Constructors
162 
163  //- Null constructor
164  radiationModel(const volScalarField& T);
165 
166  //- Construct from components
167  radiationModel(const word& type, const volScalarField& T);
168 
169  //- Construct from components
171  (
172  const word& type,
173  const dictionary& dict,
174  const volScalarField& T
175  );
176 
177 
178  // Selectors
179 
180  //- Return a reference to the selected radiation model
181  static autoPtr<radiationModel> New(const volScalarField& T);
182 
183  //- Return a reference to the selected radiation model
185  (
186  const dictionary& dict,
187  const volScalarField& T
188  );
189 
190 
191  //- Destructor
192  virtual ~radiationModel();
193 
194 
195  // Member Functions
196 
197  // Edit
198 
199  //- Main update/correction routine
200  virtual void correct();
201 
202  //- Solve radiation equation(s)
203  virtual void calculate() = 0;
204 
205  //- Read radiationProperties dictionary
206  virtual bool read() = 0;
207 
208 
209  // Access
210 
211  //- Source term component (for power of T^4)
212  virtual tmp<volScalarField> Rp() const = 0;
213 
214  //- Source term component (constant)
215  virtual tmp<volScalarField::Internal> Ru() const = 0;
216 
217  //- Energy source term
218  virtual tmp<fvScalarMatrix> Sh
219  (
220  const basicThermo& thermo,
221  const volScalarField& he
222  ) const;
223 
224  //- Temperature source term
225  virtual tmp<fvScalarMatrix> ST
226  (
227  const dimensionedScalar& rhoCp,
228  volScalarField& T
229  ) const;
230 
231  //- Access to absorptionEmission model
233 
234  //- Access to soot Model
235  const sootModel& soot() const;
236 };
237 
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 #define addToRadiationRunTimeSelectionTables(model) \
242  \
243  addToRunTimeSelectionTable \
244  ( \
245  radiationModel, \
246  model, \
247  dictionary \
248  ); \
249  \
250  addToRunTimeSelectionTable \
251  ( \
252  radiationModel, \
253  model, \
254  T \
255  );
256 
257 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258 
259 } // End namespace radiation
260 } // End namespace Foam
261 
262 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
263 
264 #endif
265 
266 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:51
dictionary dict
dictionary coeffs_
Radiation model dictionary.
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.
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:52
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.
TypeName("radiationModel")
Runtime type information.
const sootModel & soot() const
Access to soot Model.
virtual void calculate()=0
Solve radiation equation(s)
rhoReactionThermo & thermo
Definition: createFields.H:28
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
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
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.
const fvMesh & mesh_
Reference to the mesh database.
Top level model for radiation modelling.
virtual void correct()
Main update/correction routine.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
virtual tmp< volScalarField::Internal > Ru() const =0
Source term component (constant)
virtual tmp< volScalarField > Rp() const =0
Source term component (for power of T^4)
const absorptionEmissionModel & absorptionEmission() const
Access to absorptionEmission model.
virtual tmp< fvScalarMatrix > Sh(const basicThermo &thermo, const volScalarField &he) const
Energy source term.
Forward declarations of fvMatrix specializations.
const Time & time_
Reference to the time database.
virtual tmp< fvScalarMatrix > ST(const dimensionedScalar &rhoCp, volScalarField &T) const
Temperature source term.
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:52
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: PtrList.H:53
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:92
autoPtr< absorptionEmissionModel > absorptionEmission_
Absorption/emission model.
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.