radiativeIntensityRay.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-2019 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::radiationModels::radiativeIntensityRay
26 
27 Description
28  Radiation intensity for a ray in a given direction
29 
30 SourceFiles
31  radiativeIntensityRay.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef radiativeIntensityRay_H
36 #define radiativeIntensityRay_H
37 
39 #include "blackBodyEmission.H"
40 
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 namespace radiationModels
47 {
48 
49 // Forward declaration of classes
50 class fvDOM;
51 
52 /*---------------------------------------------------------------------------*\
53  Class radiativeIntensityRay Declaration
54 \*---------------------------------------------------------------------------*/
55 
57 {
58 public:
59 
60  static const word intensityPrefix;
61 
62 
63 private:
64 
65  // Private Data
66 
67  //- Reference to the owner fvDOM object
68  const fvDOM& dom_;
69 
70  //- Reference to the mesh
71  const fvMesh& mesh_;
72 
73  //- Absorption/emission model
74  const absorptionEmissionModel& absorptionEmission_;
75 
76  //- Black body
77  const blackBodyEmission& blackBody_;
78 
79  //- Total radiative intensity / [W/m^2]
80  volScalarField I_;
81 
82  //- Total radiative heat flux on boundary
83  volScalarField qr_;
84 
85  //- Incident radiative heat flux on boundary
86  volScalarField qin_;
87 
88  //- Emitted radiative heat flux on boundary
89  volScalarField qem_;
90 
91  //- Direction
92  vector d_;
93 
94  //- Average direction vector inside the solid angle
95  vector dAve_;
96 
97  //- Theta angle
98  scalar theta_;
99 
100  //- Phi angle
101  scalar phi_;
102 
103  //- Solid angle
104  scalar omega_;
105 
106  //- Number of wavelengths/bands
107  label nLambda_;
108 
109  //- List of pointers to radiative intensity fields for given wavelengths
110  PtrList<volScalarField> ILambda_;
111 
112  //- Global ray id - incremented in constructor
113  static label rayId;
114 
115  //- My ray Id
116  label myRayId_;
117 
118 
119 public:
120 
121  // Constructors
122 
123  //- Construct form components
125  (
126  const fvDOM& dom,
127  const fvMesh& mesh,
128  const scalar phi,
129  const scalar theta,
130  const scalar deltaPhi,
131  const scalar deltaTheta,
132  const label lambda,
133  const absorptionEmissionModel& absEmmModel_,
134  const blackBodyEmission& blackBody,
135  const label rayId
136  );
137 
138  //- Disallow default bitwise copy construction
140 
141 
142  //- Destructor
144 
145 
146  // Member Functions
147 
148  // Edit
149 
150  //- Update radiative intensity on i direction
151  scalar correct();
152 
153  //- Initialise the ray in i direction
154  void init
155  (
156  const scalar phi,
157  const scalar theta,
158  const scalar deltaPhi,
159  const scalar deltaTheta,
160  const scalar lambda
161  );
162 
163  //- Add radiative intensities from all the bands
164  void addIntensity();
165 
166 
167  // Access
168 
169  //- Return intensity
170  inline const volScalarField& I() const;
171 
172  //- Return const access to the boundary heat flux
173  inline const volScalarField& qr() const;
174 
175  //- Return non-const access to the boundary heat flux
176  inline volScalarField& qr();
177 
178  //- Return non-const access to the boundary incident heat flux
179  inline volScalarField& qin();
180 
181  //- Return non-const access to the boundary emitted heat flux
182  inline volScalarField& qem();
183 
184  //- Return const access to the boundary incident heat flux
185  inline const volScalarField& qin() const;
186 
187  //- Return const access to the boundary emitted heat flux
188  inline const volScalarField& qem() const;
189 
190  //- Return direction
191  inline const vector& d() const;
192 
193  //- Return the average vector inside the solid angle
194  inline const vector& dAve() const;
195 
196  //- Return the number of bands
197  inline scalar nLambda() const;
198 
199  //- Return the phi angle
200  inline scalar phi() const;
201 
202  //- Return the theta angle
203  inline scalar theta() const;
204 
205  //- Return the solid angle
206  inline scalar omega() const;
207 
208  //- Return the radiative intensity for a given wavelength
209  inline const volScalarField& ILambda(const label lambdaI) const;
210 
211 
212  // Member Operators
213 
214  //- Disallow default bitwise assignment
215  void operator=(const radiativeIntensityRay&) = delete;
216 };
217 
218 
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 
221 } // End namespace radiationModels
222 } // End namespace Foam
223 
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
225 
226 #include "radiativeIntensityRayI.H"
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 #endif
231 
232 // ************************************************************************* //
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:74
scalar nLambda() const
Return the number of bands.
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
void addIntensity()
Add radiative intensities from all the bands.
scalar theta() const
Return the theta angle.
const vector & d() const
Return direction.
scalar correct()
Update radiative intensity on i direction.
radiativeIntensityRay(const fvDOM &dom, const fvMesh &mesh, const scalar phi, const scalar theta, const scalar deltaPhi, const scalar deltaTheta, const label lambda, const absorptionEmissionModel &absEmmModel_, const blackBodyEmission &blackBody, const label rayId)
Construct form components.
const volScalarField & qr() const
Return const access to the boundary heat flux.
const vector & dAve() const
Return the average vector inside the solid angle.
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
volScalarField & qem()
Return non-const access to the boundary emitted heat flux.
void init(const scalar phi, const scalar theta, const scalar deltaPhi, const scalar deltaTheta, const scalar lambda)
Initialise the ray in i direction.
void operator=(const radiativeIntensityRay &)=delete
Disallow default bitwise assignment.
volScalarField & qin()
Return non-const access to the boundary incident heat flux.
Radiation intensity for a ray in a given direction.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
scalar phi() const
Return the phi angle.
const volScalarField & I() const
Return intensity.
Model to supply absorption and emission coefficients for radiation modelling.
const volScalarField & ILambda(const label lambdaI) const
Return the radiative intensity for a given wavelength.
scalar omega() const
Return the solid angle.
Namespace for OpenFOAM.