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-2015 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  meshOrientation (1 1 1); //Mesh ortientation used for 2D and 1D
51  }
52 
53  solverFreq 1; // Number of flow iterations per radiation iteration
54  \endverbatim
55 
56  The total number of solid angles is 4*nPhi*nTheta in 3-D.
57 
58  In 1-D the ray directions are on X, Y or Z (nPhi and nTheta are ignored).
59  'meshOrientation' vector can be used for any other 1D direction.
60 
61  In 2-D the ray directions are on X-Y, X-Z or Y-Z planes.
62  (only nPhi is considered). 'meshOrientation' vector can be used for
63  not-aligned planes specifying the plane normal vector.
64 
65  In 3D (nPhi and nTheta are considered). 'meshOrientation' vector is not
66  considered.
67 
68 SourceFiles
69  fvDOM.C
70 
71 \*---------------------------------------------------------------------------*/
72 
73 #ifndef radiationModelfvDOM_H
74 #define radiationModelfvDOM_H
75 
76 #include "radiativeIntensityRay.H"
77 #include "radiationModel.H"
78 #include "fvMatrices.H"
79 #include "solarLoad.H"
80 
81 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82 
83 namespace Foam
84 {
85 namespace radiation
86 {
87 
88 /*---------------------------------------------------------------------------*\
89  Class fvDOM Declaration
90 \*---------------------------------------------------------------------------*/
91 
92 class fvDOM
93 :
94  public radiationModel
95 {
96  // Private data
97 
98 
99  //- Incident radiation [W/m2]
101 
102  //- Total radiative heat flux [W/m2]
104 
105  //- Emmited radiative heat flux [W/m2]
107 
108  //- Incidet radiative heat flux [W/m2]
110 
111  //- Total absorption coefficient [1/m]
113 
114  //- Number of solid angles in theta
115  label nTheta_;
116 
117  //- Number of solid angles in phi
118  label nPhi_ ;
119 
120  //- Total number of rays (1 per direction)
121  label nRay_;
122 
123  //- Number of wavelength bands
124  label nLambda_;
125 
126  //- Wavelength total absorption coefficient [1/m]
128 
129  //- Black body
131 
132  //- List of pointers to radiative intensity rays
134 
135  //- Convergence criterion
136  scalar convergence_;
137 
138  //- Maximum number of iterations
139  scalar maxIter_;
140 
141  //- List of cached fvMatrices for rays
143 
144  //- Cache convection div matrix
145  bool cacheDiv_;
146 
147  //- Maximum omega weight
148  scalar omegaMax_;
149 
150  //- Use Solar Load model
151  bool useSolarLoad_;
152 
153  //- Solar load radiation model
155 
156  //- Mesh orientation vector
158 
159 
160  // Private Member Functions
161 
162  //- Initialise
163  void initialise();
164 
165  //- Disallow default bitwise copy construct
166  fvDOM(const fvDOM&);
167 
168  //- Disallow default bitwise assignment
169  void operator=(const fvDOM&);
170 
171  //- Update nlack body emission
173 
174 
175 public:
176 
177  //- Runtime type information
178  TypeName("fvDOM");
179 
180 
181  // Constructors
182 
183  //- Construct from components
184  fvDOM(const volScalarField& T);
185 
186  //- Construct from components
187  fvDOM(const dictionary& dict, const volScalarField& T);
188 
189 
190  //- Destructor
191  virtual ~fvDOM();
192 
193 
194  // Member functions
195 
196  // Edit
197 
198  //- Solve radiation equation(s)
199  void calculate();
200 
201  //- Read radiation properties dictionary
202  bool read();
203 
204  //- Update G and calculate total heat flux on boundary
205  void updateG();
206 
207  //- Set the rayId and lambdaId from by decomposing an intensity
208  // field name
209  void setRayIdLambdaId
210  (
211  const word& name,
212  label& rayId,
213  label& lambdaId
214  ) const;
215 
216  //- Source term component (for power of T^4)
217  virtual tmp<volScalarField> Rp() const;
218 
219  //- Source term component (constant)
220  virtual tmp<DimensionedField<scalar, volMesh> > Ru() const;
221 
222 
223  // Access
224 
225  //- Ray intensity for rayI
226  inline const radiativeIntensityRay& IRay(const label rayI) const;
227 
228  //- Ray intensity for rayI and lambda bandwidth
229  inline const volScalarField& IRayLambda
230  (
231  const label rayI,
232  const label lambdaI
233  ) const;
234 
235  //- Number of angles in theta
236  inline label nTheta() const;
237 
238  //- Number of angles in phi
239  inline label nPhi() const;
240 
241  //- Number of rays
242  inline label nRay() const;
243 
244  //- Number of wavelengths
245  inline label nLambda() const;
246 
247  //- Const access to total absorption coefficient
248  inline const volScalarField& a() const;
249 
250  //- Const access to wavelength total absorption coefficient
251  inline const volScalarField& aLambda(const label lambdaI) const;
252 
253  //- Const access to incident radiation field
254  inline const volScalarField& G() const;
255 
256  //- Const access to total radiative heat flux field
257  inline const volScalarField& Qr() const;
258 
259  //- Const access to incident radiative heat flux field
260  inline const volScalarField& Qin() const;
261 
262  //- Const access to emitted radiative heat flux field
263  inline const volScalarField& Qem() const;
264 
265  //- Const access to black body
266  inline const blackBodyEmission& blackBody() const;
267 
268  //- Const access to cached fvMatrix
269  inline const fvScalarMatrix& fvRayDiv
270  (
271  const label lambdaI,
272  const label rayId
273  ) const;
274 
275  //- Caching div(Ji, Ilamda)
276  inline bool cacheDiv() const;
277 
278  //- Return omegaMax
279  inline scalar omegaMax() const;
280 
281  //- Return meshOrientation
282  inline vector meshOrientation() const;
283 };
284 
285 
286 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287 
288 #include "fvDOMI.H"
289 
290 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
291 
292 } // End namespace radiation
293 } // End namespace Foam
294 
295 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
296 
297 #endif
298 
299 // ************************************************************************* //
Foam::radiation::fvDOM::nRay
label nRay() const
Number of rays.
Definition: fvDOMI.H:56
solarLoad.H
Foam::radiation::fvDOM::nPhi
label nPhi() const
Number of angles in phi.
Definition: fvDOMI.H:50
Foam::radiation::fvDOM::fvRayDiv
const fvScalarMatrix & fvRayDiv(const label lambdaI, const label rayId) const
Const access to cached fvMatrix.
Definition: fvDOMI.H:114
fvDOMI.H
Foam::radiation::fvDOM::cacheDiv_
bool cacheDiv_
Cache convection div matrix.
Definition: fvDOM.H:144
Foam::radiation::fvDOM::a
const volScalarField & a() const
Const access to total absorption coefficient.
Definition: fvDOMI.H:68
Foam::radiation::fvDOM::convergence_
scalar convergence_
Convergence criterion.
Definition: fvDOM.H:135
Foam::radiation::fvDOM::~fvDOM
virtual ~fvDOM()
Destructor.
Definition: fvDOM.C:430
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::radiation::fvDOM::omegaMax
scalar omegaMax() const
Return omegaMax.
Definition: fvDOMI.H:129
Foam::radiation::fvDOM::a_
volScalarField a_
Total absorption coefficient [1/m].
Definition: fvDOM.H:111
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::radiation::fvDOM::Qin
const volScalarField & Qin() const
Const access to incident radiative heat flux field.
Definition: fvDOMI.H:94
Foam::radiation::fvDOM::maxIter_
scalar maxIter_
Maximum number of iterations.
Definition: fvDOM.H:138
Foam::radiation::fvDOM::fvDOM
fvDOM(const fvDOM &)
Disallow default bitwise copy construct.
radiativeIntensityRay.H
Foam::radiation::fvDOM::updateG
void updateG()
Update G and calculate total heat flux on boundary.
Definition: fvDOM.C:540
Foam::radiation::radiativeIntensityRay
Radiation intensity for a ray in a given direction.
Definition: radiativeIntensityRay.H:55
Foam::radiation::fvDOM::IRay
const radiativeIntensityRay & IRay(const label rayI) const
Ray intensity for rayI.
Definition: fvDOMI.H:27
Foam::radiation::fvDOM::IRay_
PtrList< radiativeIntensityRay > IRay_
List of pointers to radiative intensity rays.
Definition: fvDOM.H:132
Foam::radiation::fvDOM::nLambda
label nLambda() const
Number of wavelengths.
Definition: fvDOMI.H:62
Foam::radiation::fvDOM::Qr_
volScalarField Qr_
Total radiative heat flux [W/m2].
Definition: fvDOM.H:102
Foam::radiation::fvDOM::aLambda_
PtrList< volScalarField > aLambda_
Wavelength total absorption coefficient [1/m].
Definition: fvDOM.H:126
Foam::radiation::fvDOM::nLambda_
label nLambda_
Number of wavelength bands.
Definition: fvDOM.H:123
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::radiation::fvDOM::initialise
void initialise()
Initialise.
Definition: fvDOM.C:50
Foam::radiation::fvDOM::IRayLambda
const volScalarField & IRayLambda(const label rayI, const label lambdaI) const
Ray intensity for rayI and lambda bandwidth.
Definition: fvDOMI.H:35
Foam::label
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
Foam::radiation::fvDOM::G
const volScalarField & G() const
Const access to incident radiation field.
Definition: fvDOMI.H:83
Foam::radiation::fvDOM::Qem
const volScalarField & Qem() const
Const access to emitted radiative heat flux field.
Definition: fvDOMI.H:100
Foam::radiation::fvDOM::aLambda
const volScalarField & aLambda(const label lambdaI) const
Const access to wavelength total absorption coefficient.
Definition: fvDOMI.H:75
Foam::radiation::fvDOM::nPhi_
label nPhi_
Number of solid angles in phi.
Definition: fvDOM.H:117
Foam::radiation::fvDOM::updateBlackBodyEmission
void updateBlackBodyEmission()
Update nlack body emission.
Definition: fvDOM.C:531
Foam::radiation::fvDOM::TypeName
TypeName("fvDOM")
Runtime type information.
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::radiation::fvDOM::Rp
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: fvDOM.C:495
Foam::radiation::fvDOM::useSolarLoad_
bool useSolarLoad_
Use Solar Load model.
Definition: fvDOM.H:150
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::radiation::fvDOM::calculate
void calculate()
Solve radiation equation(s)
Definition: fvDOM.C:453
Foam::radiation::fvDOM::Qem_
volScalarField Qem_
Emmited radiative heat flux [W/m2].
Definition: fvDOM.H:105
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::radiation::fvDOM::nTheta
label nTheta() const
Number of angles in theta.
Definition: fvDOMI.H:44
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::radiation::fvDOM::nRay_
label nRay_
Total number of rays (1 per direction)
Definition: fvDOM.H:120
Foam::radiation::fvDOM::Qr
const volScalarField & Qr() const
Const access to total radiative heat flux field.
Definition: fvDOMI.H:89
Foam::radiation::fvDOM::meshOrientation
vector meshOrientation() const
Return meshOrientation.
Definition: fvDOMI.H:135
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::radiation::fvDOM::meshOrientation_
vector meshOrientation_
Mesh orientation vector.
Definition: fvDOM.H:156
Foam::radiation::fvDOM::read
bool read()
Read radiation properties dictionary.
Definition: fvDOM.C:436
Foam::IOdictionary::name
const word & name() const
Name function is needed to disambiguate those inherited.
Definition: IOdictionary.C:181
Foam::radiation::fvDOM::Qin_
volScalarField Qin_
Incidet radiative heat flux [W/m2].
Definition: fvDOM.H:108
Foam::radiation::fvDOM::omegaMax_
scalar omegaMax_
Maximum omega weight.
Definition: fvDOM.H:147
Foam::Vector< scalar >
Foam::List
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
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:73
Foam::radiation::fvDOM::nTheta_
label nTheta_
Number of solid angles in theta.
Definition: fvDOM.H:114
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::radiation::blackBodyEmission
Class black body emission.
Definition: blackBodyEmission.H:56
Foam::radiation::fvDOM::operator=
void operator=(const fvDOM &)
Disallow default bitwise assignment.
Foam::radiation::fvDOM::G_
volScalarField G_
Incident radiation [W/m2].
Definition: fvDOM.H:99
Foam::radiation::fvDOM::blackBody_
blackBodyEmission blackBody_
Black body.
Definition: fvDOM.H:129
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::radiation::fvDOM::blackBody
const blackBodyEmission & blackBody() const
Const access to black body.
Definition: fvDOMI.H:107
Foam::radiation::fvDOM::solarLoad_
autoPtr< solarLoad > solarLoad_
Solar load radiation model.
Definition: fvDOM.H:153
radiationModel.H
Foam::radiation::fvDOM::Ru
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant)
Definition: fvDOM.C:517
Foam::radiation::fvDOM::fvRayDiv_
List< PtrList< fvScalarMatrix > > fvRayDiv_
List of cached fvMatrices for rays.
Definition: fvDOM.H:141
Foam::radiation::fvDOM::setRayIdLambdaId
void setRayIdLambdaId(const word &name, label &rayId, label &lambdaId) const
Set the rayId and lambdaId from by decomposing an intensity.
Definition: fvDOM.C:559
Foam::radiation::fvDOM::cacheDiv
bool cacheDiv() const
Caching div(Ji, Ilamda)
Definition: fvDOMI.H:123
Foam::radiation::fvDOM
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:91