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-2017 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  }
48 
49  solverFreq 1; // Number of flow iterations per radiation iteration
50  \endverbatim
51 
52  The total number of solid angles is 4*nPhi*nTheta.
53 
54  In 1D the direction of the rays is X (nPhi and nTheta are ignored)
55  In 2D the direction of the rays is on X-Y plane (only nPhi is considered)
56  In 3D (nPhi and nTheta are considered)
57 
58 SourceFiles
59  fvDOM.C
60 
61 \*---------------------------------------------------------------------------*/
62 
63 #ifndef radiationModelfvDOM_H
64 #define radiationModelfvDOM_H
65 
66 #include "radiativeIntensityRay.H"
67 #include "radiationModel.H"
68 #include "fvMatrices.H"
69 
70 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72 namespace Foam
73 {
74 namespace radiation
75 {
76 
77 /*---------------------------------------------------------------------------*\
78  Class fvDOM Declaration
79 \*---------------------------------------------------------------------------*/
80 
81 class fvDOM
82 :
83  public radiationModel
84 {
85  // Private data
86 
87 
88  //- Incident radiation [W/m2]
89  volScalarField G_;
90 
91  //- Total radiative heat flux [W/m2]
92  volScalarField qr_;
93 
94  //- Emmited radiative heat flux [W/m2]
96 
97  //- Incidet radiative heat flux [W/m2]
98  volScalarField qin_;
99 
100  //- Total absorption coefficient [1/m]
102 
103  //- Number of solid angles in theta
104  label nTheta_;
105 
106  //- Number of solid angles in phi
107  label nPhi_ ;
109  //- Total number of rays (1 per direction)
110  label nRay_;
111 
112  //- Number of wavelength bands
113  label nLambda_;
115  //- Wavelength total absorption coefficient [1/m]
116  PtrList<volScalarField> aLambda_;
117 
118  //- Black body
119  blackBodyEmission blackBody_;
120 
121  //- List of pointers to radiative intensity rays
123 
124  //- Convergence tolerance
125  scalar tolerance_;
126 
127  //- Maximum number of iterations
128  label maxIter_;
129 
130  //- Maximum omega weight
131  scalar omegaMax_;
132 
133 
134  // Private Member Functions
135 
136  //- Initialise
137  void initialise();
138 
139  //- Disallow default bitwise copy construct
140  fvDOM(const fvDOM&);
141 
142  //- Disallow default bitwise assignment
143  void operator=(const fvDOM&);
144 
145  //- Update nlack body emission
146  void updateBlackBodyEmission();
147 
148 
149 public:
150 
151  //- Runtime type information
152  TypeName("fvDOM");
153 
154 
155  // Constructors
156 
157  //- Construct from components
158  fvDOM(const volScalarField& T);
159 
160  //- Construct from components
161  fvDOM(const dictionary& dict, const volScalarField& T);
162 
163 
164  //- Destructor
165  virtual ~fvDOM();
166 
167 
168  // Member functions
169 
170  // Edit
171 
172  //- Solve radiation equation(s)
173  void calculate();
174 
175  //- Read radiation properties dictionary
176  bool read();
177 
178  //- Update G and calculate total heat flux on boundary
179  void updateG();
180 
181  //- Set the rayId and lambdaId from by decomposing an intensity
182  // field name
183  void setRayIdLambdaId
184  (
185  const word& name,
186  label& rayId,
187  label& lambdaId
188  ) const;
189 
190  //- Source term component (for power of T^4)
191  virtual tmp<volScalarField> Rp() const;
192 
193  //- Source term component (constant)
194  virtual tmp<volScalarField::Internal> Ru() const;
195 
196 
197  // Access
198 
199  //- Ray intensity for rayI
200  inline const radiativeIntensityRay& IRay(const label rayI) const;
201 
202  //- Ray intensity for rayI and lambda bandwidth
203  inline const volScalarField& IRayLambda
204  (
205  const label rayI,
206  const label lambdaI
207  ) const;
208 
209  //- Number of angles in theta
210  inline label nTheta() const;
211 
212  //- Number of angles in phi
213  inline label nPhi() const;
214 
215  //- Number of rays
216  inline label nRay() const;
217 
218  //- Number of wavelengths
219  inline label nLambda() const;
220 
221  //- Const access to total absorption coefficient
222  inline const volScalarField& a() const;
223 
224  //- Const access to wavelength total absorption coefficient
225  inline const volScalarField& aLambda(const label lambdaI) const;
226 
227  //- Const access to incident radiation field
228  inline const volScalarField& G() const;
229 
230  //- Const access to total radiative heat flux field
231  inline const volScalarField& qr() const;
232 
233  //- Const access to incident radiative heat flux field
234  inline const volScalarField& qin() const;
235 
236  //- Const access to emitted radiative heat flux field
237  inline const volScalarField& qem() const;
238 
239  //- Const access to black body
240  inline const blackBodyEmission& blackBody() const;
241 
242  //- Return omegaMax
243  inline scalar omegaMax() const;
244 };
245 
246 
247 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
248 
249 #include "fvDOMI.H"
250 
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 
253 } // End namespace radiation
254 } // End namespace Foam
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
258 #endif
259 
260 // ************************************************************************* //
const volScalarField & a() const
Const access to total absorption coefficient.
Definition: fvDOM.H:69
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 volScalarField & IRayLambda(const label rayI, const label lambdaI) const
Ray intensity for rayI and lambda bandwidth.
Definition: fvDOM.H:36
const volScalarField & qr() const
Const access to total radiative heat flux field.
Definition: fvDOM.H:90
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const blackBodyEmission & blackBody() const
Const access to black body.
Definition: fvDOM.H:108
const volScalarField & qin() const
Const access to incident radiative heat flux field.
Definition: fvDOM.H:95
void setRayIdLambdaId(const word &name, label &rayId, label &lambdaId) const
Set the rayId and lambdaId from by decomposing an intensity.
Definition: fvDOM.C:524
const volScalarField & qem() const
Const access to emitted radiative heat flux field.
Definition: fvDOM.H:101
label nPhi() const
Number of angles in phi.
Definition: fvDOM.H:51
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:505
label nRay() const
Number of rays.
Definition: fvDOM.H:57
Top level model for radiation modelling.
const word & name() const
Name function is needed to disambiguate those inherited.
const volScalarField & G() const
Const access to incident radiation field.
Definition: fvDOM.H:84
bool read()
Read radiation properties dictionary.
Definition: fvDOM.C:397
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const radiativeIntensityRay & IRay(const label rayI) const
Ray intensity for rayI.
Definition: fvDOM.H:28
Radiation intensity for a ray in a given direction.
virtual tmp< volScalarField::Internal > Ru() const
Source term component (constant)
Definition: fvDOM.C:479
label nTheta() const
Number of angles in theta.
Definition: fvDOM.H:45
virtual ~fvDOM()
Destructor.
Definition: fvDOM.C:391
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
scalar omegaMax() const
Return omegaMax.
Definition: fvDOM.H:114
TypeName("fvDOM")
Runtime type information.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
void calculate()
Solve radiation equation(s)
Definition: fvDOM.C:419
label nLambda() const
Number of wavelengths.
Definition: fvDOM.H:63
const volScalarField & aLambda(const label lambdaI) const
Const access to wavelength total absorption coefficient.
Definition: fvDOM.H:76
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: fvDOM.C:456
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:80
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.