radiativeIntensityRay.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::radiation::radiativeIntensityRay
28 
29 Description
30  Radiation intensity for a ray in a given direction
31 
32 SourceFiles
33  radiativeIntensityRay.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef radiativeIntensityRay_H
38 #define radiativeIntensityRay_H
39 
41 #include "blackBodyEmission.H"
42 
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 namespace radiation
49 {
50 
51 // Forward declaration of classes
52 class fvDOM;
53 
54 /*---------------------------------------------------------------------------*\
55  Class radiativeIntensityRay Declaration
56 \*---------------------------------------------------------------------------*/
57 
59 {
60 public:
61 
62  static const word intensityPrefix;
63 
64 
65 private:
66 
67  // Private data
68 
69  //- Reference to the owner fvDOM object
70  const fvDOM& dom_;
71 
72  //- Reference to the mesh
73  const fvMesh& mesh_;
74 
75  //- Absorption/emission model
76  const absorptionEmissionModel& absorptionEmission_;
77 
78  //- Black body
79  const blackBodyEmission& blackBody_;
80 
81  //- Total radiative intensity / [W/m2]
82  volScalarField I_;
83 
84  //- Total radiative heat flux on boundary
85  volScalarField qr_;
86 
87  //- Incident radiative heat flux on boundary
88  volScalarField qin_;
89 
90  //- Emitted radiative heat flux on boundary
91  volScalarField qem_;
92 
93  //- Direction
94  vector d_;
95 
96  //- Average direction vector inside the solid angle
97  vector dAve_;
98 
99  //- Theta angle
100  scalar theta_;
101 
102  //- Phi angle
103  scalar phi_;
104 
105  //- Solid angle
106  scalar omega_;
107 
108  //- Number of wavelengths/bands
109  label nLambda_;
110 
111  //- List of pointers to radiative intensity fields for given wavelengths
112  PtrList<volScalarField> ILambda_;
113 
114  //- Global ray id - incremented in constructor
115  static label rayId;
116 
117  //- My ray Id
118  label myRayId_;
119 
120 
121  // Private Member Functions
122 
123  //- No copy construct
125 
126  //- No copy assignment
127  void operator=(const radiativeIntensityRay&) = delete;
128 
129 
130 public:
131 
132  // Constructors
133 
134  //- Construct form components
136  (
137  const fvDOM& dom,
138  const fvMesh& mesh,
139  const scalar phi,
140  const scalar theta,
141  const scalar deltaPhi,
142  const scalar deltaTheta,
143  const label lambda,
144  const absorptionEmissionModel& absEmmModel_,
145  const blackBodyEmission& blackBody,
146  const label rayId
147  );
148 
149 
150  //- Destructor
152 
153 
154  // Member functions
155 
156  // Edit
157 
158  //- Update radiative intensity on i direction
159  scalar correct();
160 
161  //- Initialise the ray in i direction
162  void init
163  (
164  const scalar phi,
165  const scalar theta,
166  const scalar deltaPhi,
167  const scalar deltaTheta,
168  const scalar lambda
169  );
170 
171  //- Add radiative intensities from all the bands
172  void addIntensity();
173 
174 
175  // Access
176 
177  //- Return intensity
178  inline const volScalarField& I() const;
179 
180  //- Return const access to the boundary heat flux
181  inline const volScalarField& qr() const;
182 
183  //- Return non-const access to the boundary heat flux
184  inline volScalarField& qr();
185 
186  //- Return non-const access to the boundary incident heat flux
187  inline volScalarField& qin();
188 
189  //- Return non-const access to the boundary emitted heat flux
190  inline volScalarField& qem();
191 
192  //- Return const access to the boundary incident heat flux
193  inline const volScalarField& qin() const;
194 
195  //- Return const access to the boundary emitted heat flux
196  inline const volScalarField& qem() const;
197 
198  //- Return direction
199  inline const vector& d() const;
200 
201  //- Return the average vector inside the solid angle
202  inline const vector& dAve() const;
203 
204  //- Return direction
205  inline vector& d();
206 
207  //- Return the average vector inside the solid angle
208  inline vector& dAve();
209 
210  //- Return the number of bands
211  inline scalar nLambda() const;
212 
213  //- Return the phi angle
214  inline scalar phi() const;
215 
216  //- Return the theta angle
217  inline scalar theta() const;
218 
219  //- Return the solid angle
220  inline scalar omega() const;
221 
222  //- Return the radiative intensity for a given wavelength
223  inline const volScalarField& ILambda(const label lambdaI) const;
224 
225 };
226 
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 } // End namespace radiation
231 } // End namespace Foam
232 
233 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 
235 #include "radiativeIntensityRayI.H"
236 
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238 
239 #endif
240 
241 // ************************************************************************* //
Foam::radiation::radiativeIntensityRay::d
const vector & d() const
Return direction.
Definition: radiativeIntensityRayI.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::radiation::radiativeIntensityRay::~radiativeIntensityRay
~radiativeIntensityRay()
Destructor.
Definition: radiativeIntensityRay.C:263
Foam::radiation::radiativeIntensityRay::qem
volScalarField & qem()
Return non-const access to the boundary emitted heat flux.
Definition: radiativeIntensityRayI.H:67
Foam::radiation::radiativeIntensityRay::intensityPrefix
static const word intensityPrefix
Definition: radiativeIntensityRay.H:61
Foam::radiation::radiativeIntensityRay::addIntensity
void addIntensity()
Add radiative intensities from all the bands.
Definition: radiativeIntensityRay.C:316
Foam::radiation::radiativeIntensityRay
Radiation intensity for a ray in a given direction.
Definition: radiativeIntensityRay.H:57
Foam::radiation::radiativeIntensityRay::init
void init(const scalar phi, const scalar theta, const scalar deltaPhi, const scalar deltaTheta, const scalar lambda)
Initialise the ray in i direction.
Foam::radiation::radiativeIntensityRay::qr
const volScalarField & qr() const
Return const access to the boundary heat flux.
Definition: radiativeIntensityRayI.H:36
Foam::radiation::radiativeIntensityRay::correct
scalar correct()
Update radiative intensity on i direction.
Definition: radiativeIntensityRay.C:269
blackBodyEmission.H
Foam::radiation::radiativeIntensityRay::qin
volScalarField & qin()
Return non-const access to the boundary incident heat flux.
Definition: radiativeIntensityRayI.H:54
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::radiation::radiativeIntensityRay::phi
scalar phi() const
Return the phi angle.
Definition: radiativeIntensityRayI.H:103
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
lambda
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::radiation::radiativeIntensityRay::nLambda
scalar nLambda() const
Return the number of bands.
Definition: radiativeIntensityRayI.H:97
Foam::radiation::radiativeIntensityRay::dAve
const vector & dAve() const
Return the average vector inside the solid angle.
Definition: radiativeIntensityRayI.H:79
Foam::radiation::radiativeIntensityRay::I
const volScalarField & I() const
Return intensity.
Definition: radiativeIntensityRayI.H:29
absorptionEmissionModel.H
Foam::Vector< scalar >
Foam::radiation::radiativeIntensityRay::theta
scalar theta() const
Return the theta angle.
Definition: radiativeIntensityRayI.H:109
Foam::radiation::radiativeIntensityRay::omega
scalar omega() const
Return the solid angle.
Definition: radiativeIntensityRayI.H:115
radiativeIntensityRayI.H
Foam::radiation::absorptionEmissionModel
Model to supply absorption and emission coefficients for radiation modelling.
Definition: absorptionEmissionModel.H:54
Foam::radiation::blackBodyEmission
Class black body emission.
Definition: blackBodyEmission.H:58
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::radiation::radiativeIntensityRay::ILambda
const volScalarField & ILambda(const label lambdaI) const
Return the radiative intensity for a given wavelength.
Definition: radiativeIntensityRayI.H:123
Foam::radiation::fvDOM
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:118