FreeFOAM The Cross-Platform CFD Toolkit
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) 1991-2010 OpenCFD Ltd.
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 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef radiationModel_H
42 #define radiationModel_H
43 
44 #include <OpenFOAM/IOdictionary.H>
45 #include <OpenFOAM/autoPtr.H>
47 #include <finiteVolume/volFields.H>
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 namespace radiation
57 {
58 
59 // Forward declaration of classes
60 class absorptionEmissionModel;
61 class scatterModel;
62 
63 /*---------------------------------------------------------------------------*\
64  Class radiationModel Declaration
65 \*---------------------------------------------------------------------------*/
66 
67 class radiationModel
68 :
69  public IOdictionary
70 {
71 protected:
72 
73  // Protected data
74 
75  //- Reference to the mesh database
76  const fvMesh& mesh_;
77 
78  //- Reference to the time database
79  const Time& time_;
80 
81  //- Reference to the temperature field
82  const volScalarField& T_;
83 
84  //- Model specific dictionary input parameters
86 
87  //- Radiation model dictionary
89 
90  //- Radiation solver frequency - number flow solver iterations per
91  // radiation solver iteration
92  label solverFreq_;
93 
94 
95  // References to the radiation sub-models
96 
97  //- Absorption/emission model
99 
100  //- Scatter model
102 
103 
104 private:
105 
106  // Private Member Functions
107 
108  //- Disallow default bitwise copy construct
110 
111  //- Disallow default bitwise assignment
112  void operator=(const radiationModel&);
113 
114 
115 public:
116 
117  //- Runtime type information
118  TypeName("radiationModel");
119 
120 
121  // Declare runtime constructor selection table
122 
124  (
125  autoPtr,
127  dictionary,
128  (
129  const volScalarField& T
130  ),
131  (T)
132  );
133 
134 
135  // Constructors
136 
137  //- Null constructor
138  radiationModel(const volScalarField& T);
139 
140  //- Construct from components
141  radiationModel(const word& type, const volScalarField& T);
142 
143 
144  // Selectors
145 
146  //- Return a reference to the selected radiation model
147  static autoPtr<radiationModel> New(const volScalarField& T);
148 
149 
150  //- Destructor
151  virtual ~radiationModel();
152 
153 
154  // Member Functions
155 
156  // Edit
157 
158  //- Main update/correction routine
159  virtual void correct();
160 
161  //- Solve radiation equation(s)
162  virtual void calculate() = 0;
163 
164  //- Read radiationProperties dictionary
165  virtual bool read() = 0;
166 
167 
168  // Access
169 
170  //- Source term component (for power of T^4)
171  virtual tmp<volScalarField> Rp() const = 0;
172 
173  //- Source term component (constant)
174  virtual tmp<DimensionedField<scalar, volMesh> > Ru() const = 0;
175 
176  //- Enthalpy source term
177  virtual tmp<fvScalarMatrix> Sh(basicThermo& thermo) const;
178 
179  //- Sensible enthalpy source term
180  virtual tmp<fvScalarMatrix> Shs(basicThermo& thermo) const;
181 };
182 
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 } // End namespace radiation
187 } // End namespace Foam
188 
189 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190 
191 #endif
192 
193 // ************************ vim: set sw=4 sts=4 et: ************************ //