radiativeIntensityRay.C
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) 2018-2021 OpenCFD Ltd
10-------------------------------------------------------------------------------
11License
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\*---------------------------------------------------------------------------*/
28
30#include "fvm.H"
31#include "fvDOM.H"
32#include "constants.H"
33
34using namespace Foam::constant;
35
36const Foam::word
38
39
40// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41
43(
44 const fvDOM& dom,
45 const fvMesh& mesh,
46 const scalar phi,
47 const scalar theta,
48 const scalar deltaPhi,
49 const scalar deltaTheta,
50 const label nLambda,
51 const absorptionEmissionModel& absorptionEmission,
52 const blackBodyEmission& blackBody,
53 const label rayId
54)
55:
56 dom_(dom),
57 mesh_(mesh),
58 absorptionEmission_(absorptionEmission),
59 blackBody_(blackBody),
60 I_
61 (
63 (
64 "I" + name(rayId),
65 mesh_.time().timeName(),
66 mesh_,
67 IOobject::NO_READ,
68 IOobject::NO_WRITE
69 ),
70 mesh_,
72 ),
73 qr_
74 (
76 (
77 "qr" + name(rayId),
78 mesh_.time().timeName(),
79 mesh_,
80 IOobject::NO_READ,
81 IOobject::NO_WRITE
82 ),
83 mesh_,
85 ),
86 qin_
87 (
89 (
90 "qin" + name(rayId),
91 mesh_.time().timeName(),
92 mesh_,
93 IOobject::NO_READ,
94 IOobject::NO_WRITE
95 ),
96 mesh_,
98 ),
99 qem_
100 (
102 (
103 "qem" + name(rayId),
104 mesh_.time().timeName(),
105 mesh_,
106 IOobject::NO_READ,
107 IOobject::NO_WRITE
108 ),
109 mesh_,
111 ),
112 d_(Zero),
113 dAve_(Zero),
114 theta_(theta),
115 phi_(phi),
116 omega_(0.0),
117 nLambda_(nLambda),
118 ILambda_(nLambda),
119 myRayId_(rayId)
120{
121 scalar sinTheta = Foam::sin(theta);
122 scalar cosTheta = Foam::cos(theta);
123 scalar sinPhi = Foam::sin(phi);
124 scalar cosPhi = Foam::cos(phi);
125
126 omega_ = 2.0*sinTheta*Foam::sin(deltaTheta/2.0)*deltaPhi;
127 d_ = vector(sinTheta*sinPhi, sinTheta*cosPhi, cosTheta);
128 dAve_ = vector
129 (
130 sinPhi
131 *Foam::sin(0.5*deltaPhi)
132 *(deltaTheta - Foam::cos(2.0*theta)
133 *Foam::sin(deltaTheta)),
134 cosPhi
135 *Foam::sin(0.5*deltaPhi)
136 *(deltaTheta - Foam::cos(2.0*theta)
137 *Foam::sin(deltaTheta)),
138 0.5*deltaPhi*Foam::sin(2.0*theta)*Foam::sin(deltaTheta)
139 );
140
141 if (mesh_.nSolutionD() == 2)
142 {
143 // Omega for 2D
144 omega_ = deltaPhi;
145
146 // dAve for 2D
147 dAve_ = vector
148 (
149 2*sinPhi*Foam::sin(0.5*deltaPhi),
150 2*cosPhi*Foam::sin(0.5*deltaPhi),
151 0
152 );
153
154 vector meshDir(Zero);
155 if (dom_.meshOrientation() != vector::zero)
156 {
157 meshDir = dom_.meshOrientation();
158 }
159 else
160 {
161 for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
162 {
163 if (mesh_.geometricD()[cmpt] == -1)
164 {
165 meshDir[cmpt] = 1;
166 }
167 }
168 }
169 const vector normal(vector(0, 0, 1));
170
171 const tensor coordRot = rotationTensor(normal, meshDir);
172
173 dAve_ = coordRot & dAve_;
174 d_ = coordRot & d_;
175
176 }
177 else if (mesh_.nSolutionD() == 1)
178 {
179 vector meshDir(Zero);
180 if (dom_.meshOrientation() != vector::zero)
181 {
182 meshDir = dom_.meshOrientation();
183 }
184 else
185 {
186 for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
187 {
188 if (mesh_.geometricD()[cmpt] == 1)
189 {
190 meshDir[cmpt] = 1;
191 }
192 }
193 }
194 const vector normal(vector(1, 0, 0));
195
196 dAve_ = (dAve_ & normal)*meshDir;
197 d_ = (d_ & normal)*meshDir;
198
199 // Omega normalization for 1D
200 omega_ /= 2;
201 }
202
203 autoPtr<volScalarField> IDefaultPtr;
204
205 forAll(ILambda_, lambdaI)
206 {
207 IOobject IHeader
208 (
209 intensityPrefix + "_" + name(rayId) + "_" + name(lambdaI),
210 mesh_.time().timeName(),
211 mesh_,
214 );
215
216 // Check if field exists and can be read
217 if (IHeader.typeHeaderOk<volScalarField>(true))
218 {
219 ILambda_.set
220 (
221 lambdaI,
222 new volScalarField(IHeader, mesh_)
223 );
224 }
225 else
226 {
227 // Demand driven load the IDefault field
228 if (!IDefaultPtr)
229 {
230 IDefaultPtr.reset
231 (
233 (
235 (
236 "IDefault",
237 mesh_.time().timeName(),
238 mesh_,
241 ),
242 mesh_
243 )
244 );
245 }
246
247 // Reset the MUST_READ flag
248 IOobject noReadHeader(IHeader);
249 noReadHeader.readOpt(IOobject::NO_READ);
250
251 ILambda_.set
252 (
253 lambdaI,
254 new volScalarField(noReadHeader, IDefaultPtr())
255 );
256 }
257 }
258}
259
260
261// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
262
264{}
265
266
267// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
268
270{
271 // Reset boundary heat flux to zero
272 qr_.boundaryFieldRef() = 0.0;
273 qem_.boundaryFieldRef() = 0.0;
274 qin_.boundaryFieldRef() = 0.0;
275
276 scalar maxResidual = -GREAT;
277
278 forAll(ILambda_, lambdaI)
279 {
280 const volScalarField& k = dom_.aLambda(lambdaI);
281
282 const surfaceScalarField Ji(dAve_ & mesh_.Sf());
283
284 fvScalarMatrix IiEq
285 (
286 fvm::div(Ji, ILambda_[lambdaI], "div(Ji,Ii_h)")
287 + fvm::Sp(k*omega_, ILambda_[lambdaI])
288 ==
290 *(
291 (k - absorptionEmission_.aDisp(lambdaI))
292 *blackBody_.bLambda(lambdaI)
293
294 + absorptionEmission_.E(lambdaI)/4
295 )
296 );
297
298 IiEq.relax();
299
300 const solverPerformance ILambdaSol = solve
301 (
302 IiEq,
303 mesh_.solver("Ii")
304 );
305
306 const scalar initialRes =
307 ILambdaSol.initialResidual()*omega_/dom_.omegaMax();
308
309 maxResidual = max(initialRes, maxResidual);
310 }
311
312 return maxResidual;
313}
314
315
317{
319
320 forAll(ILambda_, lambdaI)
321 {
322 I_ += ILambda_[lambdaI];
323 }
324}
325
326
327// ************************************************************************* //
label k
surfaceScalarField & phi
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
readOption readOpt() const noexcept
The read option.
Definition: IOobjectI.H:164
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (uses typeFilePath to find file) and check its info.
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
const Type & initialResidual() const noexcept
Return initial residual.
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1092
autoPtr< fvSolver > solver(const dictionary &)
Construct and return the solver.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
static constexpr direction nComponents
Number of components in bool is 1.
Definition: bool.H:98
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:900
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:872
Model to supply absorption and emission coefficients for radiation modelling.
Class black body emission.
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:121
vector meshOrientation() const
Return meshOrientation.
Definition: fvDOMI.H:127
Radiation intensity for a ray in a given direction.
scalar theta() const
Return the theta angle.
scalar correct()
Update radiative intensity on i direction.
scalar phi() const
Return the phi angle.
void addIntensity()
Add radiative intensities from all the bands.
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
word timeName
Definition: getTimeIndex.H:3
constexpr scalar pi(M_PI)
Different types of constants.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:51
dimensionedScalar sin(const dimensionedScalar &ds)
uint8_t direction
Definition: direction.H:56
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Vector< scalar > vector
Definition: vector.H:61
dimensionedScalar cos(const dimensionedScalar &ds)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
CEqn solve()
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
A non-counting (dummy) refCount.
Definition: refCount.H:59