initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / thermophysicalModels / radiation / radiationModel / fvDOM / fvDOM / fvDOM.H
blobffbe5400b6979c5c6c51021d36fe29240b47a378
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::radiation::fvDOM
28 Description
30     Finite Volume Discrete Ordinates Method. Solves the RTE equation for n
31     directions in a participating media, not including scatter.
33     Available absorption models:
34         greyMeanAbsoprtionEmission
35         wideBandAbsorptionEmission
37     i.e. dictionary
38     fvDOMCoeffs
39     {
40         nPhi    1;          // azimuthal angles in PI/2 on X-Y.(from Y to X)
41         nTheta  2;          // polar angles in PI (from Z to X-Y plane)
42         convergence 1e-4;   // convergence criteria for radiation iteration
43     }
45     solverFreq   1; // Number of flow iterations per radiation iteration
47     The total number of solid angles is  4*nPhi*nTheta.
49     In 1D the direction of the rays is X (nPhi and nTheta are ignored)
50     In 2D the direction of the rays is on X-Y plane (only nPhi is considered)
51     In 3D (nPhi and nTheta are considered)
53 SourceFiles
54     fvDOM.C
56 \*---------------------------------------------------------------------------*/
58 #ifndef radiationModelfvDOM_H
59 #define radiationModelfvDOM_H
61 #include "radiativeIntensityRay.H"
62 #include "radiationModel.H"
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 namespace Foam
68 namespace radiation
71 /*---------------------------------------------------------------------------*\
72                            Class fvDOM Declaration
73 \*---------------------------------------------------------------------------*/
75 class fvDOM
77     public radiationModel
79     // Private data
81         //- Incident radiation  [W/m2]
82         volScalarField G_;
84         //- Total radiative heat flux [W/m2]
85         volScalarField Qr_;
87         //- Total absorption coefficient [1/m]
88         volScalarField a_;
90         //- Total emission coefficient [1/m]
91         volScalarField e_;
93         //- Emission contribution [Kg/m/s^3]
94         volScalarField E_;
96         //- Number of solid angles in theta
97         label nTheta_;
99         //- Number of solid angles in phi
100         label nPhi_ ;
102         //- Total number of rays (1 per direction)
103         label nRay_;
105         //- Number of wavelength bands
106         label nLambda_;
108         //- Wavelength total absorption coefficient [1/m]
109         PtrList<volScalarField> aLambda_;
111         //- Black body
112         blackBodyEmission blackBody_;
114         //- List of pointers to radiative intensity rays
115         PtrList<radiativeIntensityRay> IRay_;
117         //- Convergence criterion
118         scalar convergence_;
120         //- Maximum number of iterations
121         scalar maxIter_;
124     // Private member functions
126         //- Disallow default bitwise copy construct
127         fvDOM(const fvDOM&);
129         //- Disallow default bitwise assignment
130         void operator=(const fvDOM&);
132         //- Update nlack body emission
133         void updateBlackBodyEmission();
136 public:
138     //- Runtime type information
139     TypeName("fvDOM");
142     // Constructors
144         //- Construct from components
145         fvDOM(const volScalarField& T);
148     //- Destructor
149     virtual ~fvDOM();
152     // Member functions
154         // Edit
156             //- Solve radiation equation(s)
157             void calculate();
159             //- Read radiation properties dictionary
160             bool read();
162             //- Update G and calculate total heat flux on boundary
163             void updateG();
165             //- Set the rayId and lambdaId from by decomposing an intensity
166             //  field name
167             void setRayIdLambdaId
168             (
169                 const word& name,
170                 label& rayId,
171                 label& lambdaId
172             ) const;
174             //- Source term component (for power of T^4)
175             virtual tmp<volScalarField> Rp() const;
177             //- Source term component (constant)
178             virtual tmp<DimensionedField<scalar, volMesh> > Ru() const;
181         // Access
183             //- Ray intensity for rayI
184             inline const radiativeIntensityRay& IRay(const label rayI) const;
186             //- Ray intensity for rayI and lambda bandwidth
187             inline const volScalarField& IRayLambda
188             (
189                 const label rayI,
190                 const label lambdaI
191             ) const;
193             //- Number of angles in theta
194             inline label nTheta() const;
196             //- Number of angles in phi
197             inline label nPhi() const;
199             //- Number of rays
200             inline label nRay() const;
202             //- Number of wavelengths
203             inline label nLambda() const;
205             //- Const access to total absorption coefficient
206             inline const volScalarField& a() const;
208             //- Const access to wavelength total absorption coefficient
209             inline const volScalarField& aLambda(const label lambdaI) const;
211             //- Const access to incident radiation field
212             inline const volScalarField& G() const;
214             //- Const access to total radiative heat flux field
215             inline const volScalarField& Qr() const;
217             //- Const access to black body
218             inline const blackBodyEmission& blackBody() const;
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224 #include "fvDOMI.H"
226 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 } // End namespace radiation
229 } // End namespace Foam
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 #endif
235 // ************************************************************************* //