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 -------------------------------------------------------------------------------
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 \*---------------------------------------------------------------------------*/
28 
29 #include "radiativeIntensityRay.H"
30 #include "fvm.H"
31 #include "fvDOM.H"
32 #include "constants.H"
33 
34 using namespace Foam::constant;
35 
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 Foam::radiation::radiativeIntensityRay::radiativeIntensityRay
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  (
62  IOobject
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  (
75  IOobject
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  (
88  IOobject
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  (
101  IOobject
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  vector meshDir(Zero);
144  if (dom_.meshOrientation() != vector::zero)
145  {
146  meshDir = dom_.meshOrientation();
147  }
148  else
149  {
150  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
151  {
152  if (mesh_.geometricD()[cmpt] == -1)
153  {
154  meshDir[cmpt] = 1;
155  }
156  }
157  }
158  const vector normal(vector(0, 0, 1));
159 
160  const tensor coordRot = rotationTensor(normal, meshDir);
161 
162  dAve_ = coordRot & dAve_;
163  d_ = coordRot & d_;
164  }
165  else if (mesh_.nSolutionD() == 1)
166  {
167  vector meshDir(Zero);
168  if (dom_.meshOrientation() != vector::zero)
169  {
170  meshDir = dom_.meshOrientation();
171  }
172  else
173  {
174  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
175  {
176  if (mesh_.geometricD()[cmpt] == 1)
177  {
178  meshDir[cmpt] = 1;
179  }
180  }
181  }
182  const vector normal(vector(1, 0, 0));
183 
184  dAve_ = (dAve_ & normal)*meshDir;
185  d_ = (d_ & normal)*meshDir;
186  }
187 
188  autoPtr<volScalarField> IDefaultPtr;
189 
190  forAll(ILambda_, lambdaI)
191  {
192  IOobject IHeader
193  (
194  intensityPrefix + "_" + name(rayId) + "_" + name(lambdaI),
195  mesh_.time().timeName(),
196  mesh_,
199  );
200 
201  // Check if field exists and can be read
202  if (IHeader.typeHeaderOk<volScalarField>(true))
203  {
204  ILambda_.set
205  (
206  lambdaI,
207  new volScalarField(IHeader, mesh_)
208  );
209  }
210  else
211  {
212  // Demand driven load the IDefault field
213  if (!IDefaultPtr)
214  {
215  IDefaultPtr.reset
216  (
217  new volScalarField
218  (
219  IOobject
220  (
221  "IDefault",
222  mesh_.time().timeName(),
223  mesh_,
226  ),
227  mesh_
228  )
229  );
230  }
231 
232  // Reset the MUST_READ flag
233  IOobject noReadHeader(IHeader);
234  noReadHeader.readOpt(IOobject::NO_READ);
235 
236  ILambda_.set
237  (
238  lambdaI,
239  new volScalarField(noReadHeader, IDefaultPtr())
240  );
241  }
242  }
243 }
244 
245 
246 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
249 {}
250 
251 
252 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
253 
255 {
256  // Reset boundary heat flux to zero
257  qr_.boundaryFieldRef() = 0.0;
258  qem_.boundaryFieldRef() = 0.0;
259  qin_.boundaryFieldRef() = 0.0;
260 
261  scalar maxResidual = -GREAT;
262 
263  forAll(ILambda_, lambdaI)
264  {
265  const volScalarField& k = dom_.aLambda(lambdaI);
266 
267  const surfaceScalarField Ji(dAve_ & mesh_.Sf());
268 
269  fvScalarMatrix IiEq
270  (
271  fvm::div(Ji, ILambda_[lambdaI], "div(Ji,Ii_h)")
272  + fvm::Sp(k*omega_, ILambda_[lambdaI])
273  ==
274  1.0/constant::mathematical::pi*omega_
275  *(
276  (k - absorptionEmission_.aDisp(lambdaI))
277  *blackBody_.bLambda(lambdaI)
278 
279  + absorptionEmission_.E(lambdaI)/4
280  )
281  );
282 
283  IiEq.relax();
284 
285  const solverPerformance ILambdaSol = solve
286  (
287  IiEq,
288  mesh_.solver("Ii")
289  );
290 
291  const scalar initialRes =
292  ILambdaSol.initialResidual()*omega_/dom_.omegaMax();
293 
294  maxResidual = max(initialRes, maxResidual);
295  }
296 
297  return maxResidual;
298 }
299 
300 
302 {
304 
305  forAll(ILambda_, lambdaI)
306  {
307  I_ += ILambda_[lambdaI];
308  }
309 }
310 
311 
312 // ************************************************************************* //
Different types of constants.
Class black body emission.
uint8_t direction
Definition: direction.H:46
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Tensor< scalar > tensor
Definition: symmTensor.H:57
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
label k
Boltzmann constant.
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:47
Model to supply absorption and emission coefficients for radiation modelling.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:865
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
T cosPhi(const Vector< T > &a, const Vector< T > &b, const T &tolerance=SMALL)
Calculate angle between a and b in radians.
Definition: vectorTools.H:115
scalar correct()
Update radiative intensity on i direction.
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
scalar phi() const
Return the phi angle.
A class for handling words, derived from Foam::string.
Definition: word.H:63
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
constexpr scalar pi(M_PI)
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
Vector< scalar > vector
Definition: vector.H:57
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
dimensionedScalar sin(const dimensionedScalar &ds)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
void addIntensity()
Add radiative intensities from all the bands.
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:893
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
Automatically write from objectRegistry::writeObject()
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
vector meshOrientation() const
Return meshOrientation.
Definition: fvDOM.H:121
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:114
scalar theta() const
Return the theta angle.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const Type & initialResidual() const noexcept
Return initial residual.