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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::radiation::fvDOM
29 
30 Group
31  grpRadiationModels
32 
33 Description
34 
35  Finite Volume Discrete Ordinates Method. Solves the RTE equation for n
36  directions in a participating media, not including scatter and reflective
37  walls.
38 
39  Available absorption models:
40 
41  - constantAbsorptionEmission
42  - greyMeanAbsoprtionEmission
43  - wideBandAbsorptionEmission
44  - multiBandAbsorptionEmission
45 
46  This model can handle non-grey participating media using
47  multiBandAbsorptionEmission model. Accordingly the BC for rays should
48  be wideBandDiffussive type
49 
50 Usage
51  \verbatim
52  fvDOMCoeffs
53  {
54  nPhi 4; // azimuthal angles in PI/2 on X-Y.
55  //(from Y to X)
56  nTheta 0; // polar angles in PI (from Z to X-Y plane)
57  tolerance 1e-3; // convergence tolerance for radiation
58  // iteration
59  maxIter 4; // maximum number of iterations
60  meshOrientation (1 1 1); //Mesh orientation used for 2D and 1D
61 
62  useSolarLoad false;
63  useExternalBeam true;
64  spectralDistribution (2 1);
65  }
66 
67  solverFreq 1; // Number of flow iterations per radiation iteration
68  \endverbatim
69 
70  The total number of solid angles is 4*nPhi*nTheta in 3-D.
71 
72  Operating modes:
73  - 1-D:
74  - ray directions are on X, Y or Z
75  - \c nPhi and \c nTheta entries are ignored
76  - \c meshOrientation vector can be used for any other 1-D direction.
77  - 2-D:
78  - ray directions are on X-Y, X-Z or Y-Z planes
79  - only the \c nPhi entry is considered
80  - \c meshOrientation vector can be used for non-aligned planes
81  specifying the plane normal vector.
82  - 3-D:
83  - rays geberated in 3-D using the \c nPhi and \c nTheta entries
84  - \c meshOrientation vector is not applicable.
85 
86  useSolarLoad calculates the primary and diffusive Sun fluxes on walls in
87  addition to the RTE equations
88 
89  useExternalBeam add an external collimated beam to the domain. This option
90  is not available if useSolarLoad is true.
91 
92  spectralDistribution is the energy spectral distribution of the collimated
93  external beam.
94 
95 SourceFiles
96  fvDOM.C
97 
98 \*---------------------------------------------------------------------------*/
99 
100 #ifndef radiation_fvDOM_H
101 #define radiation_fvDOM_H
102 
103 #include "radiativeIntensityRay.H"
104 #include "radiationModel.H"
105 #include "fvMatrices.H"
106 #include "solarLoad.H"
107 
108 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
109 
110 namespace Foam
111 {
112 namespace radiation
113 {
114 
115 /*---------------------------------------------------------------------------*\
116  Class fvDOM Declaration
117 \*---------------------------------------------------------------------------*/
118 
119 class fvDOM
120 :
121  public radiationModel
122 {
123  // Private Data
124 
125  //- Incident radiation [W/m2]
126  volScalarField G_;
127 
128  //- Total radiative heat flux [W/m2]
129  volScalarField qr_;
130 
131  //- Emitted radiative heat flux [W/m2]
132  volScalarField qem_;
133 
134  //- Incident radiative heat flux [W/m2]
135  volScalarField qin_;
136 
137  //- Total absorption coefficient [1/m]
138  volScalarField a_;
139 
140  //- Number of solid angles in theta
141  label nTheta_;
142 
143  //- Number of solid angles in phi
144  label nPhi_ ;
145 
146  //- Total number of rays (1 per direction)
147  label nRay_;
148 
149  //- Number of wavelength bands
150  label nLambda_;
151 
152  //- Wavelength total absorption coefficient [1/m]
153  PtrList<volScalarField> aLambda_;
154 
155  //- Black body
156  blackBodyEmission blackBody_;
157 
158  //- List of pointers to radiative intensity rays
160 
161  //- Convergence tolerance
162  scalar tolerance_;
163 
164  //- Maximum number of iterations
165  label maxIter_;
166 
167  //- Maximum omega weight
168  scalar omegaMax_;
169 
170  //- Use Solar Load model
171  bool useSolarLoad_;
172 
173  //- Solar load radiation model
174  autoPtr<solarLoad> solarLoad_;
175 
176  //- Mesh orientation vector
177  vector meshOrientation_;
178 
179  //- Use external parallel irradiation beam
180  bool useExternalBeam_;
181 
182  //- Spectral energy distribution for the external beam
183  scalarList spectralDistribution_;
184 
185  //- Time-dependent spectral distributions
186  autoPtr<Function1<scalarField>> spectralDistributions_;
187 
188  //- Solar calculator
189  autoPtr<solarCalculator> solarCalculator_;
190 
191  //- Update Sun position index
192  label updateTimeIndex_;
193 
194 
195  // Private Member Functions
196 
197  //- Initialise
198  void initialise();
199 
200  //- No copy construct
201  fvDOM(const fvDOM&) = delete;
202 
203  //- No copy assignment
204  void operator=(const fvDOM&) = delete;
205 
206  //- Update black body emission
207  void updateBlackBodyEmission();
208 
209 
210 public:
211 
212  //- Runtime type information
213  TypeName("fvDOM");
214 
215 
216  // Constructors
217 
218  //- Construct from components
219  fvDOM(const volScalarField& T);
220 
221  //- Construct from components
222  fvDOM(const dictionary& dict, const volScalarField& T);
223 
224 
225  //- Destructor
226  virtual ~fvDOM() = default;
227 
228 
229  // Member Functions
230 
231  // Edit
232 
233  //- Solve radiation equation(s)
234  void calculate();
235 
236  //- Read radiation properties dictionary
237  bool read();
238 
239  //- Update G and calculate total heat flux on boundary
240  void updateG();
241 
242  //- Set the rayId and lambdaId from by decomposing an intensity
243  // field name
244  void setRayIdLambdaId
245  (
246  const word& name,
247  label& rayId,
248  label& lambdaId
249  ) const;
250 
251 
252  //- Rotate rays according to Sun direction
253  void updateRaysDir();
254 
255  //- Rotate rays spheric equator to sunDir
256  void rotateInitialRays(const vector& sunDir);
257 
258  //- Align closest ray to sunDir
259  void alignClosestRayToSun(const vector& sunDir);
260 
261  //- Source term component (for power of T^4)
262  virtual tmp<volScalarField> Rp() const;
263 
264  //- Source term component (constant)
265  virtual tmp<volScalarField::Internal> Ru() const;
266 
267 
268  // Access
269 
270  //- Solar calculator
271  const solarCalculator& solarCalc() const;
272 
273  //- Ray intensity for rayI
274  inline const radiativeIntensityRay& IRay(const label rayI) const;
275 
276  //- Ray intensity for rayI and lambda bandwidth
277  inline const volScalarField& IRayLambda
278  (
279  const label rayI,
280  const label lambdaI
281  ) const;
282 
283  //- Number of angles in theta
284  inline label nTheta() const;
285 
286  //- Number of angles in phi
287  inline label nPhi() const;
288 
289  //- Number of rays
290  inline label nRay() const;
291 
292  //- Number of wavelengths
293  inline label nLambda() const;
294 
295  //- Number of bands
296  inline label nBands() const;
297 
298  //- Const access to total absorption coefficient
299  inline const volScalarField& a() const;
300 
301  //- Const access to wavelength total absorption coefficient
302  inline const volScalarField& aLambda(const label lambdaI) const;
303 
304  //- Const access to incident radiation field
305  inline const volScalarField& G() const;
306 
307  //- Const access to total radiative heat flux field
308  inline const volScalarField& qr() const;
309 
310  //- Const access to incident radiative heat flux field
311  inline const volScalarField& qin() const;
312 
313  //- Const access to emitted radiative heat flux field
314  inline const volScalarField& qem() const;
315 
316  //- Const access to black body
317  inline const blackBodyEmission& blackBody() const;
318 
319  //- Return omegaMax
320  inline scalar omegaMax() const;
321 
322  //- Return meshOrientation
323  inline vector meshOrientation() const;
324 
325  //- Use solar load
326  inline bool useSolarLoad() const;
327 
328  //- Use external beam
329  inline bool useExternalBeam() const;
330 
331  //- Energy spectral distribution for external beam
332  inline const scalarList& spectralDistribution() const;
333 };
334 
335 
336 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
337 
338 #include "fvDOMI.H"
339 
340 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
341 
342 } // End namespace radiation
343 } // End namespace Foam
344 
345 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
346 
347 #endif
348 
349 // ************************************************************************* //
Foam::radiation::fvDOM::rotateInitialRays
void rotateInitialRays(const vector &sunDir)
Rotate rays spheric equator to sunDir.
Definition: fvDOM.C:54
Foam::radiation::fvDOM::nRay
label nRay() const
Number of rays.
Definition: fvDOMI.H:59
solarLoad.H
Foam::radiation::fvDOM::nPhi
label nPhi() const
Number of angles in phi.
Definition: fvDOMI.H:53
fvDOMI.H
Foam::radiation::fvDOM::a
const volScalarField & a() const
Const access to total absorption coefficient.
Definition: fvDOMI.H:76
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::radiation::fvDOM::omegaMax
scalar omegaMax() const
Return omegaMax.
Definition: fvDOMI.H:121
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
radiativeIntensityRay.H
Foam::radiation::fvDOM::updateG
void updateG()
Update G and calculate total heat flux on boundary.
Definition: fvDOM.C:742
Foam::radiation::fvDOM::alignClosestRayToSun
void alignClosestRayToSun(const vector &sunDir)
Align closest ray to sunDir.
Definition: fvDOM.C:68
Foam::radiation::radiativeIntensityRay
Radiation intensity for a ray in a given direction.
Definition: radiativeIntensityRay.H:57
Foam::radiation::fvDOM::IRay
const radiativeIntensityRay & IRay(const label rayI) const
Ray intensity for rayI.
Definition: fvDOMI.H:30
Foam::radiation::fvDOM::qem
const volScalarField & qem() const
Const access to emitted radiative heat flux field.
Definition: fvDOMI.H:108
Foam::radiation::fvDOM::qin
const volScalarField & qin() const
Const access to incident radiative heat flux field.
Definition: fvDOMI.H:102
Foam::baseIOdictionary::name
const word & name() const
Definition: baseIOdictionary.C:85
Foam::radiation::fvDOM::nLambda
label nLambda() const
Number of wavelengths.
Definition: fvDOMI.H:65
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::radiation::fvDOM::IRayLambda
const volScalarField & IRayLambda(const label rayI, const label lambdaI) const
Ray intensity for rayI and lambda bandwidth.
Definition: fvDOMI.H:38
Foam::radiation::fvDOM::useExternalBeam
bool useExternalBeam() const
Use external beam.
Definition: fvDOMI.H:139
Foam::solarCalculator
A solar calculator model providing models for the solar direction and solar loads.
Definition: solarCalculator.H:444
Foam::radiation::fvDOM::G
const volScalarField & G() const
Const access to incident radiation field.
Definition: fvDOMI.H:91
Foam::radiation::fvDOM::aLambda
const volScalarField & aLambda(const label lambdaI) const
Const access to wavelength total absorption coefficient.
Definition: fvDOMI.H:83
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:58
Foam::radiation::fvDOM::Rp
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: fvDOM.C:645
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::radiation::fvDOM::calculate
void calculate()
Solve radiation equation(s)
Definition: fvDOM.C:576
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:47
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::radiation::fvDOM::updateRaysDir
void updateRaysDir()
Rotate rays according to Sun direction.
Definition: fvDOM.C:99
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::radiation::fvDOM::solarCalc
const solarCalculator & solarCalc() const
Solar calculator.
Definition: fvDOM.C:776
Foam::radiation::fvDOM::~fvDOM
virtual ~fvDOM()=default
Destructor.
Foam::radiation::fvDOM::meshOrientation
vector meshOrientation() const
Return meshOrientation.
Definition: fvDOMI.H:127
Foam::radiation::fvDOM::qr
const volScalarField & qr() const
Const access to total radiative heat flux field.
Definition: fvDOMI.H:97
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::radiation::fvDOM::read
bool read()
Read radiation properties dictionary.
Definition: fvDOM.C:558
Foam::radiation::fvDOM::useSolarLoad
bool useSolarLoad() const
Use solar load.
Definition: fvDOMI.H:133
Foam::Vector< scalar >
Foam::List< scalar >
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:75
Foam::radiation::fvDOM::spectralDistribution
const scalarList & spectralDistribution() const
Energy spectral distribution for external beam.
Definition: fvDOMI.H:146
Foam::radiation::blackBodyEmission
Class black body emission.
Definition: blackBodyEmission.H:58
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::radiation::fvDOM::blackBody
const blackBodyEmission & blackBody() const
Const access to black body.
Definition: fvDOMI.H:115
radiationModel.H
Foam::radiation::fvDOM::Ru
virtual tmp< volScalarField::Internal > Ru() const
Source term component (constant)
Definition: fvDOM.C:689
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:761
Foam::radiation::fvDOM::nBands
label nBands() const
Number of bands.
Definition: fvDOMI.H:71
Foam::radiation::fvDOM
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:118