fvDOM.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 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  In 1-D the ray directions are bound to one of the X, Y or Z directions. The
53  total number of solid angles is 2. nPhi and nTheta are ignored.
54 
55  In 2-D the ray directions are within one of the X-Y, X-Z or Y-Z planes. The
56  total number of solid angles is 4*nPhi. nTheta is ignored.
57 
58  In 3D the rays span all directions. The total number of solid angles is
59  4*nPhi*nTheta.
60 
61 SourceFiles
62  fvDOM.C
63 
64 \*---------------------------------------------------------------------------*/
65 
66 #ifndef radiation_fvDOM_H
67 #define radiation_fvDOM_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  //- Emitted 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_;
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 tolerance
128  scalar tolerance_;
129 
130  //- Maximum number of iterations
131  label maxIter_;
132 
133  //- Maximum omega weight
134  scalar omegaMax_;
135 
136 
137  // Private Member Functions
138 
139  //- Initialise
140  void initialise();
141 
142  //- Disallow default bitwise copy construct
143  fvDOM(const fvDOM&);
144 
145  //- Disallow default bitwise assignment
146  void operator=(const fvDOM&);
147 
148  //- Update nlack body emission
149  void updateBlackBodyEmission();
150 
151 
152 public:
153 
154  //- Runtime type information
155  TypeName("fvDOM");
156 
157 
158  // Constructors
159 
160  //- Construct from components
161  fvDOM(const volScalarField& T);
162 
163  //- Construct from components
164  fvDOM(const dictionary& dict, const volScalarField& T);
165 
166 
167  //- Destructor
168  virtual ~fvDOM();
169 
170 
171  // Member functions
172 
173  // Edit
174 
175  //- Solve radiation equation(s)
176  void calculate();
177 
178  //- Read radiation properties dictionary
179  bool read();
180 
181  //- Update G and calculate total heat flux on boundary
182  void updateG();
183 
184  //- Set the rayId and lambdaId from by decomposing an intensity
185  // field name
186  void setRayIdLambdaId
187  (
188  const word& name,
189  label& rayId,
190  label& lambdaId
191  ) const;
192 
193  //- Source term component (for power of T^4)
194  virtual tmp<volScalarField> Rp() const;
195 
196  //- Source term component (constant)
197  virtual tmp<volScalarField::Internal> Ru() const;
198 
199 
200  // Access
201 
202  //- Ray intensity for rayI
203  inline const radiativeIntensityRay& IRay(const label rayI) const;
204 
205  //- Ray intensity for rayI and lambda bandwidth
206  inline const volScalarField& IRayLambda
207  (
208  const label rayI,
209  const label lambdaI
210  ) const;
211 
212  //- Number of angles in theta
213  inline label nTheta() const;
214 
215  //- Number of angles in phi
216  inline label nPhi() const;
217 
218  //- Number of rays
219  inline label nRay() const;
220 
221  //- Number of wavelengths
222  inline label nLambda() const;
223 
224  //- Const access to total absorption coefficient
225  inline const volScalarField& a() const;
226 
227  //- Const access to wavelength total absorption coefficient
228  inline const volScalarField& aLambda(const label lambdaI) const;
229 
230  //- Const access to incident radiation field
231  inline const volScalarField& G() const;
232 
233  //- Const access to total radiative heat flux field
234  inline const volScalarField& qr() const;
235 
236  //- Const access to incident radiative heat flux field
237  inline const volScalarField& qin() const;
238 
239  //- Const access to emitted radiative heat flux field
240  inline const volScalarField& qem() const;
241 
242  //- Const access to black body
243  inline const blackBodyEmission& blackBody() const;
244 
245  //- Return omegaMax
246  inline scalar omegaMax() const;
247 };
248 
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
252 #include "fvDOMI.H"
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 } // End namespace radiation
257 } // End namespace Foam
258 
259 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260 
261 #endif
262 
263 // ************************************************************************* //
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:562
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:543
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:389
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:490
label nTheta() const
Number of angles in theta.
Definition: fvDOM.H:45
virtual ~fvDOM()
Destructor.
Definition: fvDOM.C:383
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:409
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:446
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.