fvDOM.C
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-2024 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 \*---------------------------------------------------------------------------*/
25 
26 #include "fvDOM.H"
28 #include "scatterModel.H"
29 #include "constants.H"
30 #include "fvm.H"
31 #include "wedgePolyPatch.H"
32 #include "cyclicTransform.H"
34 
35 using namespace Foam::constant;
36 using namespace Foam::constant::mathematical;
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace radiationModels
43 {
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 void Foam::radiationModels::fvDOM::initialise()
53 {
54  forAll(mesh_.boundaryMesh(), patchi)
55  {
56  const polyPatch& pp = mesh_.boundaryMesh()[patchi];
57  if
58  (
59  (
60  isA<cyclicTransform>(pp)
61  && refCast<const cyclicTransform>(pp).transform().rotates()
62  )
63  || isA<wedgePolyPatch>(pp)
64  )
65  {
67  << " radiation model does not currently support"
68  " rotationally transforming patches: cyclic and wedge."
69  << exit(FatalError);
70  }
71  }
72 
73  // 3D
74  if (mesh_.nSolutionD() == 3)
75  {
76  nRay_ = 4*nPhi_*nTheta_;
77  IRay_.setSize(nRay_);
78  const scalar deltaPhi = pi/(2.0*nPhi_);
79  const scalar deltaTheta = pi/nTheta_;
80  label i = 0;
81  for (label n = 1; n <= nTheta_; n++)
82  {
83  for (label m = 1; m <= 4*nPhi_; m++)
84  {
85  const scalar thetai = (2*n - 1)*deltaTheta/2.0;
86  const scalar phii = (2*m - 1)*deltaPhi/2.0;
87  IRay_.set
88  (
89  i,
90  new radiativeIntensityRay
91  (
92  *this,
93  mesh_,
94  phii,
95  thetai,
96  deltaPhi,
97  deltaTheta,
98  nLambda_,
99  absorptionEmission_,
100  blackBody_,
101  i
102  )
103  );
104  i++;
105  }
106  }
107  }
108  // 2D
109  else if (mesh_.nSolutionD() == 2)
110  {
111  const scalar thetai = piByTwo;
112  const scalar deltaTheta = pi;
113  nRay_ = 4*nPhi_;
114  IRay_.setSize(nRay_);
115  const scalar deltaPhi = pi/(2.0*nPhi_);
116  label i = 0;
117  for (label m = 1; m <= 4*nPhi_; m++)
118  {
119  const scalar phii = (2*m - 1)*deltaPhi/2.0;
120  IRay_.set
121  (
122  i,
123  new radiativeIntensityRay
124  (
125  *this,
126  mesh_,
127  phii,
128  thetai,
129  deltaPhi,
130  deltaTheta,
131  nLambda_,
132  absorptionEmission_,
133  blackBody_,
134  i
135  )
136  );
137  i++;
138  }
139  }
140  // 1D
141  else
142  {
143  const scalar thetai = piByTwo;
144  const scalar deltaTheta = pi;
145  nRay_ = 2;
146  IRay_.setSize(nRay_);
147  const scalar deltaPhi = pi;
148  label i = 0;
149  for (label m = 1; m <= 2; m++)
150  {
151  const scalar phii = (2*m - 1)*deltaPhi/2.0;
152  IRay_.set
153  (
154  i,
155  new radiativeIntensityRay
156  (
157  *this,
158  mesh_,
159  phii,
160  thetai,
161  deltaPhi,
162  deltaTheta,
163  nLambda_,
164  absorptionEmission_,
165  blackBody_,
166  i
167  )
168  );
169  i++;
170  }
171  }
172 
173 
174  // Construct absorption field for each wavelength
175  forAll(aLambda_, lambdaI)
176  {
177  aLambda_.set
178  (
179  lambdaI,
180  new volScalarField
181  (
182  IOobject
183  (
184  "aLambda_" + Foam::name(lambdaI) ,
185  mesh_.time().name(),
186  mesh_,
189  ),
190  a_
191  )
192  );
193  }
194 
195 
196  // Calculate the maximum solid angle
197  forAll(IRay_, rayId)
198  {
199  if (omegaMax_ < IRay_[rayId].omega())
200  {
201  omegaMax_ = IRay_[rayId].omega();
202  }
203  }
204 
205 
206  Info<< typeName << ": Created " << IRay_.size() << " rays with average "
207  << "directions (dAve) and solid angles (omega)" << endl;
208  Info<< incrIndent;
209  forAll(IRay_, rayId)
210  {
211  Info<< indent
212  << "Ray " << IRay_[rayId].I().name() << ": "
213  << "dAve = " << IRay_[rayId].dAve() << ", "
214  << "omega = " << IRay_[rayId].omega() << endl;
215  }
216  Info<< decrIndent << endl;
217 }
218 
219 
220 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
221 
223 :
224  radiationModel(typeName, T),
225  G_
226  (
227  IOobject
228  (
229  "G",
230  mesh_.time().name(),
231  mesh_,
232  IOobject::NO_READ,
233  IOobject::AUTO_WRITE
234  ),
235  mesh_,
237  ),
238  qr_
239  (
240  IOobject
241  (
242  "qr",
243  mesh_.time().name(),
244  mesh_,
245  IOobject::READ_IF_PRESENT,
246  IOobject::AUTO_WRITE
247  ),
248  mesh_,
250  ),
251  qem_
252  (
253  IOobject
254  (
255  "qem",
256  mesh_.time().name(),
257  mesh_,
258  IOobject::READ_IF_PRESENT,
259  IOobject::AUTO_WRITE
260  ),
261  mesh_,
263  ),
264  qin_
265  (
266  IOobject
267  (
268  "qin",
269  mesh_.time().name(),
270  mesh_,
271  IOobject::READ_IF_PRESENT,
272  IOobject::AUTO_WRITE
273  ),
274  mesh_,
276  ),
277  a_
278  (
279  IOobject
280  (
281  "a",
282  mesh_.time().name(),
283  mesh_,
284  IOobject::NO_READ,
285  IOobject::AUTO_WRITE
286  ),
287  mesh_,
289  ),
290  nTheta_(coeffs_.lookup<label>("nTheta")),
291  nPhi_(coeffs_.lookup<label>("nPhi")),
292  nRay_(0),
293  nLambda_(absorptionEmission_->nBands()),
294  aLambda_(nLambda_),
295  blackBody_(nLambda_, T),
296  IRay_(0),
297  tolerance_
298  (
299  coeffs_.lookupOrDefaultBackwardsCompatible<scalar>
300  (
301  {"tolerance", "convergence"},
302  0
303  )
304  ),
305  maxIter_(coeffs_.lookupOrDefault<label>("maxIter", 50)),
306  omegaMax_(0)
307 {
308  initialise();
309 }
310 
311 
313 (
314  const dictionary& dict,
315  const volScalarField& T
316 )
317 :
318  radiationModel(typeName, dict, T),
319  G_
320  (
321  IOobject
322  (
323  "G",
324  mesh_.time().name(),
325  mesh_,
326  IOobject::READ_IF_PRESENT,
327  IOobject::AUTO_WRITE
328  ),
329  mesh_,
331  ),
332  qr_
333  (
334  IOobject
335  (
336  "qr",
337  mesh_.time().name(),
338  mesh_,
339  IOobject::READ_IF_PRESENT,
340  IOobject::AUTO_WRITE
341  ),
342  mesh_,
344  ),
345  qem_
346  (
347  IOobject
348  (
349  "qem",
350  mesh_.time().name(),
351  mesh_,
352  IOobject::NO_READ,
353  IOobject::NO_WRITE
354  ),
355  mesh_,
357  ),
358  qin_
359  (
360  IOobject
361  (
362  "qin",
363  mesh_.time().name(),
364  mesh_,
365  IOobject::READ_IF_PRESENT,
366  IOobject::AUTO_WRITE
367  ),
368  mesh_,
370  ),
371  a_
372  (
373  IOobject
374  (
375  "a",
376  mesh_.time().name(),
377  mesh_,
378  IOobject::NO_READ,
379  IOobject::NO_WRITE
380  ),
381  mesh_,
383  ),
384  nTheta_(coeffs_.lookup<label>("nTheta")),
385  nPhi_(coeffs_.lookup<label>("nPhi")),
386  nRay_(0),
387  nLambda_(absorptionEmission_->nBands()),
388  aLambda_(nLambda_),
389  blackBody_(nLambda_, T),
390  IRay_(0),
391  tolerance_
392  (
393  coeffs_.lookupOrDefaultBackwardsCompatible<scalar>
394  (
395  {"tolerance", "convergence"},
396  0
397  )
398  ),
399  maxIter_(coeffs_.lookupOrDefault<label>("maxIter", 50)),
400  omegaMax_(0)
401 {
402  initialise();
403 }
404 
405 
406 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
407 
409 {}
410 
411 
412 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
413 
415 {
416  if (radiationModel::read())
417  {
418  // Only reading solution parameters - not changing ray geometry
419 
420  // For backward-compatibility
421  coeffs_.readIfPresent("convergence", tolerance_);
422  coeffs_.readIfPresent("tolerance", tolerance_);
423  coeffs_.readIfPresent("maxIter", maxIter_);
424 
425  return true;
426  }
427  else
428  {
429  return false;
430  }
431 }
432 
433 
435 {
436  absorptionEmission_->correct(a_, aLambda_);
437 
438  updateBlackBodyEmission();
439 
440  // Set rays converged false
441  List<bool> rayIdConv(nRay_, false);
442 
443  scalar maxResidual = 0;
444  label radIter = 0;
445  do
446  {
447  Info<< "Radiation solver iter: " << radIter << endl;
448 
449  radIter++;
450  maxResidual = 0;
451  forAll(IRay_, rayI)
452  {
453  if (!rayIdConv[rayI])
454  {
455  scalar maxBandResidual = IRay_[rayI].correct();
456  maxResidual = max(maxBandResidual, maxResidual);
457 
458  if (maxBandResidual < tolerance_)
459  {
460  rayIdConv[rayI] = true;
461  }
462  }
463  }
464 
465  } while (maxResidual > tolerance_ && radIter < maxIter_);
466 
467  updateG();
468 }
469 
470 
472 {
473  // Construct using contribution from first frequency band
475  (
477  (
478  "Rp",
479  4
481  *(aLambda_[0] - absorptionEmission_->aDisp(0)())
482  *blackBody_.deltaLambdaT(T_, absorptionEmission_->bands(0))
483  )
484  );
485 
486  volScalarField& Rp=tRp.ref();
487 
488  // Add contributions over remaining frequency bands
489  for (label j=1; j < nLambda_; j++)
490  {
491  Rp +=
492  (
493  4
495  *(aLambda_[j] - absorptionEmission_->aDisp(j)())
496  *blackBody_.deltaLambdaT(T_, absorptionEmission_->bands(j))
497  );
498  }
499 
500  return tRp;
501 }
502 
503 
506 {
508  (
510  (
511  IOobject
512  (
513  "Ru",
514  mesh_.time().name(),
515  mesh_,
518  false
519  ),
520  mesh_,
521  dimensionedScalar(dimensionSet(1, -1, -3, 0, 0), 0)
522  )
523  );
524 
526 
527  // Sum contributions over all frequency bands
528  for (label j=0; j < nLambda_; j++)
529  {
530  // Compute total incident radiation within frequency band
532  (
533  IRay_[0].ILambda(j)()*IRay_[0].omega()
534  );
535 
536  for (label rayI=1; rayI < nRay_; rayI++)
537  {
538  Gj.ref() += IRay_[rayI].ILambda(j)()*IRay_[rayI].omega();
539  }
540 
541  Ru += (aLambda_[j]() - absorptionEmission_->aDisp(j)()())*Gj
542  - absorptionEmission_->ECont(j)()();
543  }
544 
545  return tRu;
546 }
547 
548 
549 void Foam::radiationModels::fvDOM::updateBlackBodyEmission()
550 {
551  for (label j=0; j < nLambda_; j++)
552  {
553  blackBody_.correct(j, absorptionEmission_->bands(j));
554  }
555 }
556 
557 
559 {
560  G_ = dimensionedScalar("zero",dimMass/pow3(dimTime), 0);
561  qr_ = dimensionedScalar("zero",dimMass/pow3(dimTime), 0);
564 
565  forAll(IRay_, rayI)
566  {
567  IRay_[rayI].addIntensity();
568  G_ += IRay_[rayI].I()*IRay_[rayI].omega();
569  qr_.boundaryFieldRef() += IRay_[rayI].qr().boundaryField();
570  qem_.boundaryFieldRef() += IRay_[rayI].qem().boundaryField();
571  qin_.boundaryFieldRef() += IRay_[rayI].qin().boundaryField();
572  }
573 }
574 
575 
577 (
578  const word& name,
579  label& rayId,
580  label& lambdaId
581 ) const
582 {
583  // Assume name is in the form: <name>_<rayId>_<lambdaId>
584  const size_type i1 = name.find_first_of("_");
585  const size_type i2 = name.find_last_of("_");
586 
587  rayId = readLabel(IStringStream(name.substr(i1 + 1, i2 - 1))());
588  lambdaId = readLabel(IStringStream(name.substr(i2 + 1, name.size() - 1))());
589 }
590 
591 
592 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:294
Input from memory buffer stream.
Definition: IStringStream.H:52
label size_type
The type that can represent the size of a DLList.
Definition: UILList.H:177
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
Top level model for radiation modelling.
virtual bool read()=0
Read radiationProperties dictionary.
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:77
fvDOM(const volScalarField &T)
Construct from components.
Definition: fvDOM.C:222
virtual tmp< volScalarField::Internal > Ru() const
Source term component (constant)
Definition: fvDOM.C:505
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: fvDOM.C:471
void updateG()
Update G and calculate total heat flux on boundary.
Definition: fvDOM.C:558
void setRayIdLambdaId(const word &name, label &rayId, label &lambdaId) const
Set the rayId and lambdaId from by decomposing an intensity.
Definition: fvDOM.C:577
virtual ~fvDOM()
Destructor.
Definition: fvDOM.C:408
bool read()
Read radiation properties dictionary.
Definition: fvDOM.C:414
void calculate()
Solve radiation equation(s)
Definition: fvDOM.C:434
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const scalar omega
label patchi
mathematical constants.
const scalar piByTwo(0.5 *pi)
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
Collection of constants.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:241
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:234
const dimensionSet dimTime
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:504
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimMass
label readLabel(Istream &is)
Definition: label.H:64
error FatalError
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:227
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
#define addToRadiationRunTimeSelectionTables(model)
dictionary dict