fvDOM.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-2018 OpenFOAM Foundation
9  Copyright (C) 2019 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 "fvDOM.H"
31 #include "scatterModel.H"
32 #include "constants.H"
33 #include "unitConversion.H"
34 #include "fvm.H"
36 
37 using namespace Foam::constant;
38 using namespace Foam::constant::mathematical;
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44  namespace radiation
45  {
46  defineTypeNameAndDebug(fvDOM, 0);
48  }
49 }
50 
51 
52 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
55 {
56  // Rotate Y spherical cordinates to Sun direction.
57  // Solid angles on the equator are better fit for planar radiation
58  const tensor coordRot = rotationTensor(vector(0, 1, 0), sunDir);
59 
60  forAll(IRay_, rayId)
61  {
62  IRay_[rayId].dAve() = coordRot & IRay_[rayId].dAve();
63  IRay_[rayId].d() = coordRot & IRay_[rayId].d();
64  }
65 }
66 
67 
69 {
70  label SunRayId(-1);
71  scalar maxSunRay = -GREAT;
72 
73  // Looking for the ray closest to the Sun direction
74  forAll(IRay_, rayId)
75  {
76  const vector& iD = IRay_[rayId].d();
77  scalar dir = sunDir & iD;
78  if (dir > maxSunRay)
79  {
80  maxSunRay = dir;
81  SunRayId = rayId;
82  }
83  }
84 
85  // Second rotation to align colimated radiation with the closest ray
86  const tensor coordRot = rotationTensor(IRay_[SunRayId].d(), sunDir);
87 
88  forAll(IRay_, rayId)
89  {
90  IRay_[rayId].dAve() = coordRot & IRay_[rayId].dAve();
91  IRay_[rayId].d() = coordRot & IRay_[rayId].d();
92  }
93 
94  Info << "Sun direction : " << sunDir << nl << endl;
95  Info << "Sun ray ID : " << SunRayId << nl << endl;
96 }
97 
98 
100 {
101  solarCalculator_->correctSunDirection();
102  const vector sunDir = solarCalculator_->direction();
103 
104  // First iteration
105  if (updateTimeIndex_ == 0)
106  {
107  rotateInitialRays(sunDir);
108  alignClosestRayToSun(sunDir);
109  }
110  else if (updateTimeIndex_ > 0)
111  {
112  alignClosestRayToSun(sunDir);
113  }
114 }
115 
116 
117 void Foam::radiation::fvDOM::initialise()
118 {
119  coeffs_.readIfPresent("useExternalBeam", useExternalBeam_);
120 
121  if (useExternalBeam_)
122  {
123  coeffs_.readEntry("spectralDistribution", spectralDistribution_);
124 
125  spectralDistribution_ =
126  spectralDistribution_/sum(spectralDistribution_);
127 
128  const dictionary& solarDict = this->subDict("solarCalculatorCoeffs");
129  solarCalculator_.reset(new solarCalculator(solarDict, mesh_));
130 
131  if (mesh_.nSolutionD() != 3)
132  {
134  << "External beam model only available in 3D meshes "
135  << abort(FatalError);
136  }
137 
138  if (solarCalculator_->diffuseSolarRad() > 0)
139  {
141  << "External beam model does not support Diffuse "
142  << "Solar Radiation. Set diffuseSolarRad to zero"
143  << abort(FatalError);
144  }
145  if (spectralDistribution_.size() != nLambda_)
146  {
148  << "The epectral energy distribution has different bands "
149  << "than the absoprtivity model "
150  << abort(FatalError);
151  }
152  }
153 
154  // 3D
155  if (mesh_.nSolutionD() == 3)
156  {
157  nRay_ = 4*nPhi_*nTheta_;
158 
159  IRay_.setSize(nRay_);
160 
161  const scalar deltaPhi = pi/(2*nPhi_);
162  const scalar deltaTheta = pi/nTheta_;
163 
164  label i = 0;
165 
166  for (label n = 1; n <= nTheta_; n++)
167  {
168  for (label m = 1; m <= 4*nPhi_; m++)
169  {
170  scalar thetai = (2*n - 1)*deltaTheta/2.0;
171  scalar phii = (2*m - 1)*deltaPhi/2.0;
172 
173  IRay_.set
174  (
175  i,
176  new radiativeIntensityRay
177  (
178  *this,
179  mesh_,
180  phii,
181  thetai,
182  deltaPhi,
183  deltaTheta,
184  nLambda_,
185  *absorptionEmission_,
186  blackBody_,
187  i
188  )
189  );
190  i++;
191  }
192  }
193  }
194  // 2D
195  else if (mesh_.nSolutionD() == 2)
196  {
197  const scalar thetai = piByTwo;
198  const scalar deltaTheta = pi;
199  nRay_ = 4*nPhi_;
200  IRay_.setSize(nRay_);
201  const scalar deltaPhi = pi/(2.0*nPhi_);
202  label i = 0;
203  for (label m = 1; m <= 4*nPhi_; m++)
204  {
205  const scalar phii = (2*m - 1)*deltaPhi/2.0;
206  IRay_.set
207  (
208  i,
209  new radiativeIntensityRay
210  (
211  *this,
212  mesh_,
213  phii,
214  thetai,
215  deltaPhi,
216  deltaTheta,
217  nLambda_,
218  *absorptionEmission_,
219  blackBody_,
220  i
221  )
222  );
223  i++;
224  }
225  }
226  // 1D
227  else
228  {
229  const scalar thetai = piByTwo;
230  const scalar deltaTheta = pi;
231  nRay_ = 2;
232  IRay_.setSize(nRay_);
233  const scalar deltaPhi = pi;
234  label i = 0;
235  for (label m = 1; m <= 2; m++)
236  {
237  const scalar phii = (2*m - 1)*deltaPhi/2.0;
238  IRay_.set
239  (
240  i,
241  new radiativeIntensityRay
242  (
243  *this,
244  mesh_,
245  phii,
246  thetai,
247  deltaPhi,
248  deltaTheta,
249  nLambda_,
250  *absorptionEmission_,
251  blackBody_,
252  i
253  )
254  );
255  i++;
256  }
257  }
258 
259 
260  // Construct absorption field for each wavelength
261  forAll(aLambda_, lambdaI)
262  {
263  aLambda_.set
264  (
265  lambdaI,
266  new volScalarField
267  (
268  IOobject
269  (
270  "aLambda_" + Foam::name(lambdaI) ,
271  mesh_.time().timeName(),
272  mesh_,
275  ),
276  a_
277  )
278  );
279  }
280 
281  Info<< "fvDOM : Allocated " << IRay_.size()
282  << " rays with average orientation:" << nl;
283 
284  if (useExternalBeam_)
285  {
286  // Rotate rays for Sun direction
287  updateRaysDir();
288  }
289 
290  scalar totalOmega = 0;
291  forAll(IRay_, rayId)
292  {
293  if (omegaMax_ < IRay_[rayId].omega())
294  {
295  omegaMax_ = IRay_[rayId].omega();
296  }
297  totalOmega += IRay_[rayId].omega();
298  Info<< '\t' << IRay_[rayId].I().name() << " : " << "dAve : "
299  << '\t' << IRay_[rayId].dAve() << " : " << "omega : "
300  << '\t' << IRay_[rayId].omega() << " : " << "d : "
301  << '\t' << IRay_[rayId].d() << nl;
302  }
303 
304  Info << "Total omega : " << totalOmega << endl;
305 
306  Info<< endl;
307 
308  coeffs_.readIfPresent("useSolarLoad", useSolarLoad_);
309 
310  if (useSolarLoad_)
311  {
312  if (useExternalBeam_)
313  {
315  << "External beam with fvDOM can not be used "
316  << "with the solar load model"
317  << abort(FatalError);
318  }
319  const dictionary& solarDict = this->subDict("solarLoadCoeffs");
320  solarLoad_.reset(new solarLoad(solarDict, T_));
321 
322  if (solarLoad_->nBands() != this->nBands())
323  {
325  << "Requested solar radiation with fvDOM. Using "
326  << "different number of bands for the solar load is not allowed"
327  << abort(FatalError);
328  }
329 
330  Info<< "Creating Solar Load Model " << nl;
331  }
332 }
333 
334 
335 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
336 
337 Foam::radiation::fvDOM::fvDOM(const volScalarField& T)
338 :
339  radiationModel(typeName, T),
340  G_
341  (
342  IOobject
343  (
344  "G",
345  mesh_.time().timeName(),
346  mesh_,
347  IOobject::NO_READ,
348  IOobject::AUTO_WRITE
349  ),
350  mesh_,
352  ),
353  qr_
354  (
355  IOobject
356  (
357  "qr",
358  mesh_.time().timeName(),
359  mesh_,
360  IOobject::READ_IF_PRESENT,
361  IOobject::AUTO_WRITE
362  ),
363  mesh_,
365  ),
366  qem_
367  (
368  IOobject
369  (
370  "qem",
371  mesh_.time().timeName(),
372  mesh_,
373  IOobject::NO_READ,
374  IOobject::NO_WRITE
375  ),
376  mesh_,
378  ),
379  qin_
380  (
381  IOobject
382  (
383  "qin",
384  mesh_.time().timeName(),
385  mesh_,
386  IOobject::READ_IF_PRESENT,
387  IOobject::AUTO_WRITE
388  ),
389  mesh_,
391  ),
392  a_
393  (
394  IOobject
395  (
396  "a",
397  mesh_.time().timeName(),
398  mesh_,
399  IOobject::NO_READ,
400  IOobject::NO_WRITE
401  ),
402  mesh_,
404  ),
405  nTheta_(coeffs_.get<label>("nTheta")),
406  nPhi_(coeffs_.get<label>("nPhi")),
407  nRay_(0),
408  nLambda_(absorptionEmission_->nBands()),
409  aLambda_(nLambda_),
410  blackBody_(nLambda_, T),
411  IRay_(0),
412  tolerance_
413  (
414  coeffs_.lookupOrDefaultCompat("tolerance", {{"convergence", 1712}}, 0.0)
415  ),
416  maxIter_(coeffs_.lookupOrDefault<label>("maxIter", 50)),
417  omegaMax_(0),
418  useSolarLoad_(false),
419  solarLoad_(),
420  meshOrientation_
421  (
422  coeffs_.lookupOrDefault<vector>("meshOrientation", Zero)
423  ),
424  useExternalBeam_(false),
425  spectralDistribution_(),
426  solarCalculator_(),
427  updateTimeIndex_(0)
428 {
429  initialise();
430 }
431 
432 
433 Foam::radiation::fvDOM::fvDOM
434 (
435  const dictionary& dict,
436  const volScalarField& T
437 )
438 :
439  radiationModel(typeName, dict, T),
440  G_
441  (
442  IOobject
443  (
444  "G",
445  mesh_.time().timeName(),
446  mesh_,
449  ),
450  mesh_,
452  ),
453  qr_
454  (
455  IOobject
456  (
457  "qr",
458  mesh_.time().timeName(),
459  mesh_,
462  ),
463  mesh_,
465  ),
466  qem_
467  (
468  IOobject
469  (
470  "qem",
471  mesh_.time().timeName(),
472  mesh_,
475  ),
476  mesh_,
478  ),
479  qin_
480  (
481  IOobject
482  (
483  "qin",
484  mesh_.time().timeName(),
485  mesh_,
488  ),
489  mesh_,
491  ),
492  a_
493  (
494  IOobject
495  (
496  "a",
497  mesh_.time().timeName(),
498  mesh_,
501  ),
502  mesh_,
504  ),
505  nTheta_(coeffs_.get<label>("nTheta")),
506  nPhi_(coeffs_.get<label>("nPhi")),
507  nRay_(0),
508  nLambda_(absorptionEmission_->nBands()),
509  aLambda_(nLambda_),
510  blackBody_(nLambda_, T),
511  IRay_(0),
512  tolerance_
513  (
514  coeffs_.lookupOrDefaultCompat("tolerance", {{"convergence", 1712}}, 0.0)
515  ),
516  maxIter_(coeffs_.lookupOrDefault<label>("maxIter", 50)),
517  omegaMax_(0),
518  useSolarLoad_(false),
519  solarLoad_(),
520  meshOrientation_
521  (
522  coeffs_.lookupOrDefault<vector>("meshOrientation", Zero)
523  ),
524  useExternalBeam_(false),
525  spectralDistribution_(),
526  solarCalculator_(),
527  updateTimeIndex_(0)
528 {
529  initialise();
530 }
531 
532 
533 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
534 
536 {}
537 
538 
539 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
540 
542 {
543  if (radiationModel::read())
544  {
545  // Only reading solution parameters - not changing ray geometry
546  coeffs_.readIfPresentCompat
547  (
548  "tolerance", {{"convergence", 1712}}, tolerance_
549  );
550  coeffs_.readIfPresent("maxIter", maxIter_);
551 
552  return true;
553  }
554 
555  return false;
556 }
557 
558 
560 {
561  absorptionEmission_->correct(a_, aLambda_);
562 
563  updateBlackBodyEmission();
564 
565  if (useSolarLoad_)
566  {
567  solarLoad_->calculate();
568  }
569 
570  if (useExternalBeam_)
571  {
572  switch (solarCalculator_->sunDirectionModel())
573  {
575  {
576  break;
577  }
579  {
580  label updateIndex = label
581  (
582  mesh_.time().value()
583  /solarCalculator_->sunTrackingUpdateInterval()
584  );
585 
586  if (updateIndex > updateTimeIndex_)
587  {
588  Info << "Updating Sun position..." << endl;
589  updateTimeIndex_ = updateIndex;
590  updateRaysDir();
591  }
592  break;
593  }
594  }
595  }
596 
597  // Set rays convergence false
598  List<bool> rayIdConv(nRay_, false);
599 
600  scalar maxResidual = 0;
601  label radIter = 0;
602  do
603  {
604  Info<< "Radiation solver iter: " << radIter << endl;
605 
606  radIter++;
607  maxResidual = 0;
608  forAll(IRay_, rayI)
609  {
610  if (!rayIdConv[rayI])
611  {
612  scalar maxBandResidual = IRay_[rayI].correct();
613  maxResidual = max(maxBandResidual, maxResidual);
614 
615  if (maxBandResidual < tolerance_)
616  {
617  rayIdConv[rayI] = true;
618  }
619  }
620  }
621 
622  } while (maxResidual > tolerance_ && radIter < maxIter_);
623 
624  updateG();
625 }
626 
627 
629 {
630  // Construct using contribution from first frequency band
632  (
633  new volScalarField
634  (
635  IOobject
636  (
637  "Rp",
638  mesh_.time().timeName(),
639  mesh_,
642  false
643  ),
644  (
645  4
647  *(aLambda_[0] - absorptionEmission_->aDisp(0)())
648  *blackBody_.deltaLambdaT(T_, absorptionEmission_->bands(0))
649  )
650  )
651  );
652 
653  volScalarField& Rp=tRp.ref();
654 
655  // Add contributions over remaining frequency bands
656  for (label j=1; j < nLambda_; j++)
657  {
658  Rp +=
659  (
660  4
662  *(aLambda_[j] - absorptionEmission_->aDisp(j)())
663  *blackBody_.deltaLambdaT(T_, absorptionEmission_->bands(j))
664  );
665  }
666 
667  return tRp;
668 }
669 
670 
673 {
675  (
677  (
678  IOobject
679  (
680  "Ru",
681  mesh_.time().timeName(),
682  mesh_,
685  false
686  ),
687  mesh_,
688  dimensionedScalar(dimensionSet(1, -1, -3, 0, 0), Zero)
689  )
690  );
691 
692  DimensionedField<scalar, volMesh>& Ru=tRu.ref();
693 
694  // Sum contributions over all frequency bands
695  for (label j=0; j < nLambda_; j++)
696  {
697  // Compute total incident radiation within frequency band
699  (
700  IRay_[0].ILambda(j)()*IRay_[0].omega()
701  );
702 
703  for (label rayI=1; rayI < nRay_; rayI++)
704  {
705  Gj.ref() += IRay_[rayI].ILambda(j)()*IRay_[rayI].omega();
706  }
707 
708  Ru += (aLambda_[j]() - absorptionEmission_->aDisp(j)()())*Gj
709  - absorptionEmission_->ECont(j)()();
710  }
711 
712  return tRu;
713 }
714 
715 
716 void Foam::radiation::fvDOM::updateBlackBodyEmission()
717 {
718  for (label j=0; j < nLambda_; j++)
719  {
720  blackBody_.correct(j, absorptionEmission_->bands(j));
721  }
722 }
723 
724 
726 {
731 
732  forAll(IRay_, rayI)
733  {
734  IRay_[rayI].addIntensity();
735  G_ += IRay_[rayI].I()*IRay_[rayI].omega();
736  qr_.boundaryFieldRef() += IRay_[rayI].qr().boundaryField();
737  qem_.boundaryFieldRef() += IRay_[rayI].qem().boundaryField();
738  qin_.boundaryFieldRef() += IRay_[rayI].qin().boundaryField();
739  }
740 }
741 
742 
744 (
745  const word& name,
746  label& rayId,
747  label& lambdaId
748 ) const
749 {
750  // Assuming name is in the form: CHARS_rayId_lambdaId
751  const auto i1 = name.find('_');
752  const auto i2 = name.find('_', i1+1);
753 
754  rayId = readLabel(name.substr(i1+1, i2-i1-1));
755  lambdaId = readLabel(name.substr(i2+1));
756 }
757 
758 
760 {
761  return solarCalculator_();
762 }
763 
764 
765 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::radiation::fvDOM::rotateInitialRays
void rotateInitialRays(const vector &sunDir)
Rotate rays spheric equator to sunDir.
Definition: fvDOM.C:54
Foam::Tensor< scalar >
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::radiation::fvDOM::~fvDOM
virtual ~fvDOM()
Destructor.
Definition: fvDOM.C:535
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
scatterModel.H
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::radiation::fvDOM::updateG
void updateG()
Update G and calculate total heat flux on boundary.
Definition: fvDOM.C:725
Foam::radiation::fvDOM::alignClosestRayToSun
void alignClosestRayToSun(const vector &sunDir)
Align closest ray to sunDir.
Definition: fvDOM.C:68
unitConversion.H
Unit conversion functions.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:65
fvDOM.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::solarCalculator
The solar calculator model provides information about the Sun direction and Sun load model....
Definition: solarCalculator.H:105
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::radiation::radiationModel::read
virtual bool read()=0
Read radiationProperties dictionary.
Definition: radiationModel.C:210
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::label
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:62
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
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:628
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::radiation::fvDOM::calculate
void calculate()
Solve radiation equation(s)
Definition: fvDOM.C:559
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
fvm.H
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
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::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::solarCalculator::mSunDirConstant
Definition: solarCalculator.H:114
Foam::radiation::fvDOM::solarCalc
const solarCalculator & solarCalc() const
Solar calculator.
Definition: fvDOM.C:759
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
constants.H
Foam::constant::mathematical::piByTwo
constexpr scalar piByTwo(0.5 *M_PI)
Foam::radiation::fvDOM::read
bool read()
Read radiation properties dictionary.
Definition: fvDOM.C:541
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::constant::mathematical
Mathematical constants.
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:372
absorptionEmissionModel.H
Foam::Vector< scalar >
Foam::readLabel
label readLabel(const char *buf)
Parse entire buffer as a label, skipping leading/trailing whitespace.
Definition: label.H:70
Foam::List< bool >
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:75
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
addToRadiationRunTimeSelectionTables
#define addToRadiationRunTimeSelectionTables(model)
Definition: radiationModel.H:288
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Foam::rotationTensor
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:49
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::radiation::fvDOM::Ru
virtual tmp< volScalarField::Internal > Ru() const
Source term component (constant)
Definition: fvDOM.C:672
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:744
Foam::solarCalculator::mSunDirTracking
Definition: solarCalculator.H:115
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54