FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

radiativeIntensityRay.H

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 2008-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 Class
00025     Foam::radiation::radiativeIntensityRay
00026 
00027 Description
00028     Radiation intensity for a ray in a given direction
00029 
00030 SourceFiles
00031     radiativeIntensityRay.C
00032 
00033 \*---------------------------------------------------------------------------*/
00034 
00035 #ifndef radiativeIntensityRay_H
00036 #define radiativeIntensityRay_H
00037 
00038 #include <radiation/absorptionEmissionModel.H>
00039 #include <radiation/blackBodyEmission.H>
00040 
00041 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00042 
00043 namespace Foam
00044 {
00045 namespace radiation
00046 {
00047 
00048 // Forward declaration of classes
00049 class fvDOM;
00050 
00051 /*---------------------------------------------------------------------------*\
00052                     Class radiativeIntensityRay Declaration
00053 \*---------------------------------------------------------------------------*/
00054 
00055 class radiativeIntensityRay
00056 {
00057 public:
00058 
00059     static const word intensityPrefix;
00060 
00061 
00062 private:
00063 
00064     // Private data
00065 
00066         //- Refence to the owner fvDOM object
00067         const fvDOM& dom_;
00068 
00069         //- Reference to the mesh
00070         const fvMesh& mesh_;
00071 
00072         //- Absorption/emission model
00073         const absorptionEmissionModel& absorptionEmission_;
00074 
00075         //- Black body
00076         const blackBodyEmission& blackBody_;
00077 
00078         //- Total radiative intensity / [W/m2]
00079         volScalarField I_;
00080 
00081         //- Total radiative heat flux on boundary
00082         volScalarField Qr_;
00083 
00084         //- Incident radiative heat flux on boundary
00085         volScalarField Qin_;
00086 
00087         //- Emitted radiative heat flux on boundary
00088         volScalarField Qem_;
00089 
00090         //- Direction
00091         vector d_;
00092 
00093         //- Average direction vector inside the solid angle
00094         vector dAve_;
00095 
00096         //- Theta angle
00097         scalar theta_;
00098 
00099         //- Phi angle
00100         scalar phi_;
00101 
00102         //- Solid angle
00103         scalar omega_;
00104 
00105         //- Number of wavelengths/bands
00106         label nLambda_;
00107 
00108         //- List of pointers to radiative intensity fields for given wavelengths
00109         PtrList<volScalarField> ILambda_;
00110 
00111 
00112     // Private member functions
00113 
00114         //- Disallow default bitwise copy construct
00115         radiativeIntensityRay(const radiativeIntensityRay&);
00116 
00117         //- Disallow default bitwise assignment
00118         void operator=(const radiativeIntensityRay&);
00119 
00120 
00121 public:
00122 
00123     // Constructors
00124 
00125         //- Construct form components
00126         radiativeIntensityRay
00127         (
00128             const fvDOM& dom,
00129             const fvMesh& mesh,
00130             const scalar phi,
00131             const scalar theta,
00132             const scalar deltaPhi,
00133             const scalar deltaTheta,
00134             const label lambda,
00135             const absorptionEmissionModel& absEmmModel_,
00136             const blackBodyEmission& blackBody,
00137             const label rayId
00138         );
00139 
00140 
00141     // Destructor
00142     ~radiativeIntensityRay();
00143 
00144 
00145     // Member functions
00146 
00147         // Edit
00148 
00149             //- Update radiative intensity on i direction
00150             scalar correct();
00151 
00152             //- Initialise the ray in i direction
00153             void init
00154             (
00155                 const scalar phi,
00156                 const scalar theta,
00157                 const scalar deltaPhi,
00158                 const scalar deltaTheta,
00159                 const scalar lambda
00160             );
00161 
00162             //- Add radiative intensities from all the bands
00163             void addIntensity();
00164 
00165 
00166         // Access
00167 
00168             //- Return intensity
00169             inline const volScalarField& I() const;
00170 
00171             //- Return const access to the boundary heat flux
00172             inline const volScalarField& Qr() const;
00173 
00174             //- Return non-const access to the boundary heat flux
00175             inline volScalarField& Qr();
00176 
00177             //- Return non-const access to the boundary incident heat flux
00178             inline volScalarField& Qin();
00179 
00180             //- Return non-const access to the boundary emmited heat flux
00181             inline volScalarField& Qem();
00182 
00183             //- Return const access to the boundary incident heat flux
00184             inline const volScalarField& Qin() const;
00185 
00186             //- Return const access to the boundary emmited heat flux
00187             inline const volScalarField& Qem() const;
00188 
00189             //- Return direction
00190             inline const vector& d() const;
00191 
00192             //- Return the average vector inside the solid angle
00193             inline const vector& dAve() const;
00194 
00195             //- Return the number of bands
00196             inline scalar nLambda() const;
00197 
00198             //- Return the phi angle
00199             inline scalar phi() const;
00200 
00201             //- Return the theta angle
00202             inline scalar theta() const;
00203 
00204             //- Return the solid angle
00205             inline scalar omega() const;
00206 
00207             //- Return the radiative intensity for a given wavelength
00208             inline const volScalarField& ILambda(const label lambdaI) const;
00209 };
00210 
00211 
00212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00213 
00214 } // End namespace radiation
00215 } // End namespace Foam
00216 
00217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00218 
00219 #include "radiativeIntensityRayI.H"
00220 
00221 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00222 
00223 #endif
00224 
00225 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines