fvDOM.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 Class
25  Foam::radiation::fvDOM
26 
27 Description
28 
29  Finite Volume Discrete Ordinates Method. Solves the RTE equation for n
30  directions in a participating media, not including scatter.
31 
32  Available absorption models:
33  constantAbsorptionEmission
34  greyMeanAbsoprtionEmission
35  wideBandAbsorptionEmission
36 
37  i.e. dictionary
38  \verbatim
39  fvDOMCoeffs
40  {
41  nPhi 4; // azimuthal angles in PI/2 on X-Y.
42  //(from Y to X)
43  nTheta 0; // polar angles in PI (from Z to X-Y plane)
44  convergence 1e-3; // convergence criteria for radiation
45  //iteration
46  maxIter 4; // maximum number of iterations
47  cacheDiv true; // cache the div of the RTE equation.
48  //NOTE: Caching div is "only" accurate if the upwind scheme is used
49  //in div(Ji,Ii_h)
50  }
51 
52  solverFreq 1; // Number of flow iterations per radiation iteration
53  \endverbatim
54 
55  The total number of solid angles is 4*nPhi*nTheta.
56 
57  In 1D the direction of the rays is X (nPhi and nTheta are ignored)
58  In 2D the direction of the rays is on X-Y plane (only nPhi is considered)
59  In 3D (nPhi and nTheta are considered)
60 
61 SourceFiles
62  fvDOM.C
63 
64 \*---------------------------------------------------------------------------*/
65 
66 #ifndef radiationModelfvDOM_H
67 #define radiationModelfvDOM_H
68 
70 #include "radiationModel.H"
71 #include "fvMatrices.H"
72 
73 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
74 
75 namespace Foam
76 {
77 namespace radiation
78 {
79 
80 /*---------------------------------------------------------------------------*\
81  Class fvDOM Declaration
82 \*---------------------------------------------------------------------------*/
83 
84 class fvDOM
85 :
86  public radiationModel
87 {
88  // Private data
89 
90 
91  //- Incident radiation [W/m2]
92  volScalarField G_;
93 
94  //- Total radiative heat flux [W/m2]
96 
97  //- Emmited radiative heat flux [W/m2]
98  volScalarField Qem_;
99 
100  //- Incidet radiative heat flux [W/m2]
102 
103  //- Total absorption coefficient [1/m]
104  volScalarField a_;
105 
106  //- Number of solid angles in theta
107  label nTheta_;
109  //- Number of solid angles in phi
110  label nPhi_ ;
111 
112  //- Total number of rays (1 per direction)
113  label nRay_;
114 
115  //- Number of wavelength bands
116  label nLambda_;
117 
118  //- Wavelength total absorption coefficient [1/m]
119  PtrList<volScalarField> aLambda_;
120 
121  //- Black body
122  blackBodyEmission blackBody_;
123 
124  //- List of pointers to radiative intensity rays
126 
127  //- Convergence criterion
128  scalar convergence_;
129 
130  //- Maximum number of iterations
131  label maxIter_;
132 
133  //- List of cached fvMatrices for rays
135 
136  //- Cache convection div matrix
137  bool cacheDiv_;
138 
139  //- Maximum omega weight
140  scalar omegaMax_;
141 
142 
143  // Private Member Functions
144 
145  //- Initialise
146  void initialise();
147 
148  //- Disallow default bitwise copy construct
149  fvDOM(const fvDOM&);
150 
151  //- Disallow default bitwise assignment
152  void operator=(const fvDOM&);
153 
154  //- Update nlack body emission
155  void updateBlackBodyEmission();
156 
157 
158 public:
159 
160  //- Runtime type information
161  TypeName("fvDOM");
162 
163 
164  // Constructors
165 
166  //- Construct from components
167  fvDOM(const volScalarField& T);
168 
169  //- Construct from components
170  fvDOM(const dictionary& dict, const volScalarField& T);
171 
172 
173  //- Destructor
174  virtual ~fvDOM();
175 
176 
177  // Member functions
178 
179  // Edit
180 
181  //- Solve radiation equation(s)
182  void calculate();
183 
184  //- Read radiation properties dictionary
185  bool read();
186 
187  //- Update G and calculate total heat flux on boundary
188  void updateG();
189 
190  //- Set the rayId and lambdaId from by decomposing an intensity
191  // field name
192  void setRayIdLambdaId
193  (
194  const word& name,
195  label& rayId,
196  label& lambdaId
197  ) const;
198 
199  //- Source term component (for power of T^4)
200  virtual tmp<volScalarField> Rp() const;
201 
202  //- Source term component (constant)
203  virtual tmp<DimensionedField<scalar, volMesh>> Ru() const;
204 
205 
206  // Access
207 
208  //- Ray intensity for rayI
209  inline const radiativeIntensityRay& IRay(const label rayI) const;
210 
211  //- Ray intensity for rayI and lambda bandwidth
212  inline const volScalarField& IRayLambda
213  (
214  const label rayI,
215  const label lambdaI
216  ) const;
217 
218  //- Number of angles in theta
219  inline label nTheta() const;
220 
221  //- Number of angles in phi
222  inline label nPhi() const;
223 
224  //- Number of rays
225  inline label nRay() const;
226 
227  //- Number of wavelengths
228  inline label nLambda() const;
229 
230  //- Const access to total absorption coefficient
231  inline const volScalarField& a() const;
232 
233  //- Const access to wavelength total absorption coefficient
234  inline const volScalarField& aLambda(const label lambdaI) const;
235 
236  //- Const access to incident radiation field
237  inline const volScalarField& G() const;
238 
239  //- Const access to total radiative heat flux field
240  inline const volScalarField& Qr() const;
241 
242  //- Const access to incident radiative heat flux field
243  inline const volScalarField& Qin() const;
244 
245  //- Const access to emitted radiative heat flux field
246  inline const volScalarField& Qem() const;
247 
248  //- Const access to black body
249  inline const blackBodyEmission& blackBody() const;
250 
251  //- Const access to cached fvMatrix
252  inline const fvScalarMatrix& fvRayDiv
253  (
254  const label lambdaI,
255  const label rayId
256  ) const;
257 
258  //- Caching div(Ji, Ilamda)
259  inline bool cacheDiv() const;
260 
261  //- Return omegaMax
262  inline scalar omegaMax() const;
263 };
264 
265 
266 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
267 
268 #include "fvDOMI.H"
269 
270 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
271 
272 } // End namespace radiation
273 } // End namespace Foam
274 
275 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276 
277 #endif
278 
279 // ************************************************************************* //
Class black body emission.
dictionary dict
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
const radiativeIntensityRay & IRay(const label rayI) const
Ray intensity for rayI.
Definition: fvDOM.H:28
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
label nPhi() const
Number of angles in phi.
Definition: fvDOM.H:51
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
const volScalarField & G() const
Const access to incident radiation field.
Definition: fvDOM.H:84
void setRayIdLambdaId(const word &name, label &rayId, label &lambdaId) const
Set the rayId and lambdaId from by decomposing an intensity.
Definition: fvDOM.C:538
scalar omegaMax() const
Return omegaMax.
Definition: fvDOM.H:130
const volScalarField & aLambda(const label lambdaI) const
Const access to wavelength total absorption coefficient.
Definition: fvDOM.H:76
const volScalarField & IRayLambda(const label rayI, const label lambdaI) const
Ray intensity for rayI and lambda bandwidth.
Definition: fvDOM.H:36
bool cacheDiv() const
Caching div(Ji, Ilamda)
Definition: fvDOM.H:124
label nLambda() const
Number of wavelengths.
Definition: fvDOM.H:63
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant)
Definition: fvDOM.C:493
const volScalarField & Qin() const
Const access to incident radiative heat flux field.
Definition: fvDOM.H:95
const blackBodyEmission & blackBody() const
Const access to black body.
Definition: fvDOM.H:108
A class for handling words, derived from string.
Definition: word.H:59
void updateG()
Update G and calculate total heat flux on boundary.
Definition: fvDOM.C:519
const volScalarField & Qr() const
Const access to total radiative heat flux field.
Definition: fvDOM.H:90
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:71
Top level model for radiation modelling.
const word & name() const
Name function is needed to disambiguate those inherited.
Definition: IOdictionary.C:181
label nTheta() const
Number of angles in theta.
Definition: fvDOM.H:45
const volScalarField & Qem() const
Const access to emitted radiative heat flux field.
Definition: fvDOM.H:101
bool read()
Read radiation properties dictionary.
Definition: fvDOM.C:416
const fvScalarMatrix & fvRayDiv(const label lambdaI, const label rayId) const
Const access to cached fvMatrix.
Definition: fvDOM.H:115
const volScalarField & a() const
Const access to total absorption coefficient.
Definition: fvDOM.H:69
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Radiation intensity for a ray in a given direction.
virtual ~fvDOM()
Destructor.
Definition: fvDOM.C:410
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
TypeName("fvDOM")
Runtime type information.
label nRay() const
Number of rays.
Definition: fvDOM.H:57
A special matrix type and solver, designed for finite volume solutions of scalar equations.
void calculate()
Solve radiation equation(s)
Definition: fvDOM.C:433
A class for managing temporary objects.
Definition: PtrList.H:54
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:83
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: fvDOM.C:470