solarLoad.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) 2015 OpenFOAM Foundation
9  Copyright (C) 2016-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 "solarLoad.H"
30 #include "surfaceFields.H"
31 #include "vectorList.H"
34 #include "gravityMeshObject.H"
35 #include "cyclicAMIPolyPatch.H"
36 #include "mappedPatchBase.H"
37 #include "wallPolyPatch.H"
38 #include "constants.H"
39 
40 using namespace Foam::constant;
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46  namespace radiation
47  {
48  defineTypeNameAndDebug(solarLoad, 0);
50  }
51 }
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
55 void Foam::radiation::solarLoad::updateReflectedRays
56 (
57  const labelHashSet& includePatches
58 )
59 {
60  if (!reflectedFaces_ && hitFaces_)
61  {
62  reflectedFaces_.reset
63  (
64  new faceReflecting
65  (
66  mesh_,
67  hitFaces_(),
68  solarCalc_,
69  spectralDistribution_,
70  dict_
71  )
72  );
73  }
74 
75  reflectedFaces_->correct();
76 
77  volScalarField::Boundary& qrBf = qr_.boundaryFieldRef();
78  const scalarField& V = mesh_.V();
79  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
80 
81  forAll(qrBf, patchID)
82  {
83  if (includePatches[patchID])
84  {
85  for (label bandI = 0; bandI < nBands_; ++bandI)
86  {
87  qrBf[patchID] +=
88  reflectedFaces_->qreflective(bandI).boundaryField()[patchID];
89  }
90  }
91  else
92  {
93  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
94  const labelUList& cellIs = patches[patchID].faceCells();
95 
96  for (label bandI = 0; bandI < nBands_; ++bandI)
97  {
98  forAll(cellIs, i)
99  {
100  const label cellI = cellIs[i];
101 
102  Ru_[cellI] +=
103  (
104  reflectedFaces_->qreflective(bandI).
105  boundaryField()[patchID][i] * sf[i]
106  )/V[cellI];
107  }
108  }
109  }
110  }
111 }
112 
113 
114 bool Foam::radiation::solarLoad::updateHitFaces()
115 {
116  if (!hitFaces_)
117  {
118  hitFaces_.reset(new faceShading(mesh_, solarCalc_.direction()));
119  return true;
120  }
121  else
122  {
123  switch (solarCalc_.sunDirectionModel())
124  {
126  {
127  return false;
128  break;
129  }
131  {
132  const label updateIndex = label
133  (
134  mesh_.time().value()/solarCalc_.sunTrackingUpdateInterval()
135  );
136 
137  if (updateIndex > updateTimeIndex_)
138  {
139  Info << "Updating Sun position..." << endl;
140  updateTimeIndex_ = updateIndex;
141  solarCalc_.correctSunDirection();
142  hitFaces_->direction() = solarCalc_.direction();
143  hitFaces_->correct();
144  return true;
145  }
146  break;
147  }
148  }
149  }
150 
151  return false;
152 }
153 
154 
155 void Foam::radiation::solarLoad::updateAbsorptivity
156 (
157  const labelHashSet& includePatches
158 )
159 {
160  const boundaryRadiationProperties& boundaryRadiation =
162 
163  for (const label patchID : includePatches)
164  {
165  absorptivity_[patchID].setSize(nBands_);
166  for (label bandI = 0; bandI < nBands_; ++bandI)
167  {
168  absorptivity_[patchID][bandI] =
169  boundaryRadiation.absorptivity(patchID, bandI);
170  }
171  }
172 }
173 
174 
175 void Foam::radiation::solarLoad::updateDirectHitRadiation
176 (
177  const labelList& hitFacesId,
178  const labelHashSet& includeMappedPatchBasePatches
179 )
180 {
181  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
182  const scalarField& V = mesh_.V();
183  volScalarField::Boundary& qrBf = qr_.boundaryFieldRef();
184 
185  // Reset qr and qrPrimary
186  qrBf = 0.0;
187 
188  for (label bandI = 0; bandI < nBands_; ++bandI)
189  {
190  volScalarField::Boundary& qprimaryBf =
191  qprimaryRad_[bandI].boundaryFieldRef();
192 
193  qprimaryBf = 0.0;
194 
195  forAll(hitFacesId, i)
196  {
197  const label faceI = hitFacesId[i];
198  const label patchID = patches.whichPatch(faceI);
199  const polyPatch& pp = patches[patchID];
200  const label localFaceI = faceI - pp.start();
201  const vector qPrim =
202  solarCalc_.directSolarRad()*solarCalc_.direction();
203 
204  const vectorField& n = pp.faceNormals();
205 
206  {
207  qprimaryBf[patchID][localFaceI] +=
208  (qPrim & n[localFaceI])
209  * spectralDistribution_[bandI]
210  * absorptivity_[patchID][bandI]()[localFaceI];
211  }
212 
213  if (includeMappedPatchBasePatches[patchID])
214  {
215  qrBf[patchID][localFaceI] += qprimaryBf[patchID][localFaceI];
216  }
217  else
218  {
219  const vectorField& sf = mesh_.Sf().boundaryField()[patchID];
220  const label cellI = pp.faceCells()[localFaceI];
221 
222  Ru_[cellI] +=
223  (qPrim & sf[localFaceI])
224  * spectralDistribution_[bandI]
225  * absorptivity_[patchID][bandI]()[localFaceI]
226  / V[cellI];
227  }
228  }
229  }
230 }
231 
232 
233 void Foam::radiation::solarLoad::updateSkyDiffusiveRadiation
234 (
235  const labelHashSet& includePatches,
236  const labelHashSet& includeMappedPatchBasePatches
237 )
238 {
239  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
240  const scalarField& V = mesh_.V();
241  volScalarField::Boundary& qrBf = qr_.boundaryFieldRef();
242 
243  switch(solarCalc_.sunLoadModel())
244  {
247  {
248  for (const label patchID : includePatches)
249  {
250  const polyPatch& pp = patches[patchID];
251  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
252 
253  const vectorField n = pp.faceNormals();
254  const labelList& cellIds = pp.faceCells();
255 
256  forAll(n, faceI)
257  {
258  const scalar cosEpsilon(verticalDir_ & -n[faceI]);
259 
260  scalar Ed(0.0);
261  scalar Er(0.0);
262  const scalar cosTheta(solarCalc_.direction() & -n[faceI]);
263 
264  {
265  // Above the horizon
266  if (cosEpsilon == 0.0)
267  {
268  // Vertical walls
269  scalar Y(0);
270 
271  if (cosTheta > -0.2)
272  {
273  Y = 0.55+0.437*cosTheta + 0.313*sqr(cosTheta);
274  }
275  else
276  {
277  Y = 0.45;
278  }
279 
280  Ed = solarCalc_.C()*Y*solarCalc_.directSolarRad();
281  }
282  else
283  {
284  //Other than vertical walls
285  Ed =
286  solarCalc_.C()
287  * solarCalc_.directSolarRad()
288  * (1.0 + cosEpsilon)/2.0;
289  }
290 
291  // Ground reflected
292  Er =
293  solarCalc_.directSolarRad()
294  * (solarCalc_.C() + Foam::sin(solarCalc_.beta()))
295  * solarCalc_.groundReflectivity()
296  * (1.0 - cosEpsilon)/2.0;
297  }
298 
299  const label cellI = cellIds[faceI];
300  if (includeMappedPatchBasePatches[patchID])
301  {
302  for (label bandI = 0; bandI < nBands_; ++bandI)
303  {
304  qrBf[patchID][faceI] +=
305  (Ed + Er)
306  * spectralDistribution_[bandI]
307  * absorptivity_[patchID][bandI]()[faceI];
308  }
309  }
310  else
311  {
312  for (label bandI = 0; bandI < nBands_; ++bandI)
313  {
314  Ru_[cellI] +=
315  (Ed + Er)
316  * spectralDistribution_[bandI]
317  * absorptivity_[patchID][bandI]()[faceI]
318  * sf[faceI]/V[cellI];
319  }
320  }
321  }
322  }
323  }
324  break;
325 
328  {
329  for (const label patchID : includePatches)
330  {
331  const polyPatch& pp = patches[patchID];
332  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
333 
334  const labelList& cellIds = pp.faceCells();
335  forAll(pp, faceI)
336  {
337  const label cellI = cellIds[faceI];
338  if (includeMappedPatchBasePatches[patchID])
339  {
340  for (label bandI = 0; bandI < nBands_; ++bandI)
341  {
342  qrBf[patchID][faceI] +=
343  solarCalc_.diffuseSolarRad()
344  * spectralDistribution_[bandI]
345  * absorptivity_[patchID][bandI]()[faceI];
346  }
347  }
348  else
349  {
350  for (label bandI = 0; bandI < nBands_; ++bandI)
351  {
352  Ru_[cellI] +=
353  (
354  spectralDistribution_[bandI]
355  * absorptivity_[patchID][bandI]()[faceI]
356  * solarCalc_.diffuseSolarRad()
357  )*sf[faceI]/V[cellI];
358  }
359  }
360  }
361  }
362  break;
363  }
364  }
365 }
366 
367 
368 void Foam::radiation::solarLoad::initialise(const dictionary& coeffs)
369 {
370  spectralDistributions_.reset
371  (
372  Function1<scalarField>::New
373  (
374  "spectralDistribution",
375  coeffs,
376  &mesh_
377  )
378  );
379 
380  spectralDistribution_ =
381  spectralDistributions_->value(mesh_.time().timeOutputValue());
382 
383  nBands_ = spectralDistribution_.size();
384 
385  spectralDistribution_ =
386  spectralDistribution_/sum(spectralDistribution_);
387 
388  qprimaryRad_.setSize(nBands_);
389 
390  forAll(qprimaryRad_, bandI)
391  {
392  qprimaryRad_.set
393  (
394  bandI,
395  new volScalarField
396  (
397  IOobject
398  (
399  "qprimaryRad_" + Foam::name(bandI) ,
400  mesh_.time().timeName(),
401  mesh_,
404  ),
405  mesh_,
407  )
408  );
409  }
410 
411  if (coeffs.readIfPresent("gridUp", verticalDir_))
412  {
413  verticalDir_.normalise();
414  }
415  else
416  {
418  meshObjects::gravity::New(mesh_.time());
419  verticalDir_ = (-g/mag(g)).value();
420  }
421 
422  coeffs.readIfPresent("solidCoupled", solidCoupled_);
423  coeffs.readIfPresent("wallCoupled", wallCoupled_);
424  coeffs.readIfPresent("updateAbsorptivity", updateAbsorptivity_);
425  coeffs.readEntry("useReflectedRays", useReflectedRays_);
426 }
427 
428 /*
429 void Foam::radiation::solarLoad::calculateQdiff
430 (
431  const labelHashSet& includePatches,
432  const labelHashSet& includeMappedPatchBasePatches
433 )
434 {
435  scalarListIOList FmyProc
436  (
437  IOobject
438  (
439  "F",
440  mesh_.facesInstance(),
441  mesh_,
442  IOobject::MUST_READ,
443  IOobject::NO_WRITE,
444  false
445  )
446  );
447 
448 
449  if (!coarseMesh_ && !finalAgglom_.empty())
450  {
451  coarseMesh_.reset
452  (
453  new singleCellFvMesh
454  (
455  IOobject
456  (
457  "coarse:" + mesh_.name(),
458  mesh_.polyMesh::instance(),
459  mesh_.time(),
460  IOobject::NO_READ,
461  IOobject::NO_WRITE
462  ),
463  mesh_,
464  finalAgglom_
465  )
466  );
467  }
468 
469  label nLocalVFCoarseFaces = 0;
470  forAll(includePatches_, i)
471  {
472  const label patchI = includePatches_[i];
473  nLocalVFCoarseFaces += coarseMesh_->boundaryMesh()[patchI].size();
474  }
475 
476  label totalFVNCoarseFaces = nLocalVFCoarseFaces;
477  reduce(totalFVNCoarseFaces, sumOp<label>());
478 
479  // Calculate weighted absorptivity on coarse patches
480  List<scalar> localCoarseEave(nLocalVFCoarseFaces);
481  List<scalar> localTave(nLocalVFCoarseFaces);
482  List<scalar> localCoarsePartialArea(nLocalVFCoarseFaces);
483  List<vector> localCoarseNorm(nLocalVFCoarseFaces);
484 
485  scalarField compactCoarseEave(map_->constructSize(), 0.0);
486  scalarField compactCoarseTave(map_->constructSize(), 0.0);
487 
488  scalarField compactCoarsePartialArea(map_->constructSize(), 0.0);
489 
490  vectorList compactCoarseNorm(map_->constructSize(), Zero);
491 
492  const boundaryRadiationProperties& boundaryRadiation =
493  boundaryRadiationProperties::New(mesh_);
494 
495  coarseToFine_.setSize(includePatches_.size());
496 
497  const labelList& hitFacesId = hitFaces_->rayStartFaces();
498 
499  label startI = 0;
500  label compactI = 0;
501  forAll(includePatches_, i)
502  {
503  const label patchID = includePatches_[i];
504  const polyPatch& pp = mesh_.boundaryMesh()[patchID];
505 
506  const polyPatch& cpp = coarseMesh_->boundaryMesh()[patchID];
507 
508  const labelList& agglom = finalAgglom_[patchID];
509 
510  if (agglom.size() > 0)
511  {
512  label nAgglom = max(agglom) + 1;
513  coarseToFine_[i] = invertOneToMany(nAgglom, agglom);
514  }
515 
516  // Weight emissivity by spectral distribution
517  scalarField e(pp.size(), 0.0);
518 
519  for (label bandI = 0; bandI < nBands_; bandI++)
520  {
521  const tmp<scalarField> te =
522  spectralDistribution_[bandI]
523  *boundaryRadiation.diffReflectivity(patchID, bandI);
524  e += te();
525  }
526 
527  scalarList Eave(cpp.size(), 0.0);
528  scalarList Tave(cpp.size(), 0.0);
529 
530  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
531  const scalarField& Tf = T_.boundaryField()[patchID];
532 
533  const labelList& coarsePatchFace=coarseMesh_->patchFaceMap()[patchID];
534 
535  forAll(cpp, coarseI)
536  {
537  const label coarseFaceID = coarsePatchFace[coarseI];
538  const labelList& fineFaces = coarseToFine_[i][coarseFaceID];
539 
540  UIndirectList<scalar> fineSf(sf, fineFaces);
541  scalar fineArea = sum(fineSf());
542 
543  scalar fullArea = 0.0;
544  forAll(fineFaces, j)
545  {
546  label faceI = fineFaces[j];
547  label globalFaceI = faceI + pp.start();
548 
549  if (hitFacesId.found(globalFaceI))
550  {
551  fullArea += sf[faceI];
552  }
553  Eave[coarseI] += (e[faceI]*sf[faceI])/fineArea;
554  Tave[coarseI] += (pow4(Tf[faceI])*sf[faceI])/fineArea;
555  }
556  localCoarsePartialArea[compactI++] = fullArea/fineArea;
557  }
558 
559  SubList<scalar>
560  (
561  localCoarseEave,
562  Eave.size(),
563  startI
564  ) = Eave;
565 
566  SubList<scalar>
567  (
568  localTave,
569  Tave.size(),
570  startI
571  ) = Tave;
572 
573  const vectorList coarseNSf = cpp.faceNormals();
574  SubList<vector>
575  (
576  localCoarseNorm,
577  cpp.size(),
578  startI
579  ) = coarseNSf;
580 
581  startI += cpp.size();
582  }
583 
584 
585  SubList<scalar>(compactCoarsePartialArea, nLocalVFCoarseFaces) =
586  localCoarsePartialArea;
587 
588  SubList<scalar>(compactCoarseEave, nLocalVFCoarseFaces) =
589  localCoarseEave;
590 
591  SubList<scalar>(compactCoarseTave, nLocalVFCoarseFaces) =
592  localTave;
593 
594  SubList<vector>(compactCoarseNorm, nLocalVFCoarseFaces) =
595  localCoarseNorm;
596 
597  map_->distribute(compactCoarsePartialArea);
598  map_->distribute(compactCoarseEave);
599  map_->distribute(compactCoarseTave);
600  map_->distribute(compactCoarseNorm);
601 
602 
603  // Calculate coarse hitFaces and Sun direct hit heat fluxes
604  scalarList localqDiffusive(nLocalVFCoarseFaces, Zero);
605 
606  label locaFaceI = 0;
607  forAll(includePatches_, i)
608  {
609  const label patchID = includePatches_[i];
610  const polyPatch& pp = coarseMesh_->boundaryMesh()[patchID];
611  const polyPatch& ppf = mesh_.boundaryMesh()[patchID];
612 
613  const labelList& coarsePatchFace = coarseMesh_->patchFaceMap()[patchID];
614  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
615 
616 
617  scalarField a(ppf.size(), Zero);
618 
619  for (label bandI = 0; bandI < nBands_; bandI++)
620  {
621  const tmp<scalarField> ta =
622  spectralDistribution_[bandI]
623  * absorptivity_[patchID][bandI]();
624  a += ta();
625  }
626 
627  forAll(pp, coarseI)
628  {
629  const label coarseFaceID = coarsePatchFace[coarseI];
630  const labelList& fineFaces = coarseToFine_[i][coarseFaceID];
631  UIndirectList<scalar> fineSf(sf, fineFaces);
632  scalar fineArea = sum(fineSf());
633 
634  // // Weighting absorptivity per area on secondary diffussive flux
635  scalar aAve = 0.0;
636  forAll(fineFaces, j)
637  {
638  label faceI = fineFaces[j];
639  aAve += (a[faceI]*sf[faceI])/fineArea;
640  }
641 
642  const scalarList& vf = FmyProc[locaFaceI];
643  const labelList& compactFaces = visibleFaceFaces_[locaFaceI];
644 
645  forAll(compactFaces, j)
646  {
647  label compactI = compactFaces[j];
648 
649  scalar qin =
650  (
651  solarCalc_.directSolarRad()*solarCalc_.direction()
652  & compactCoarseNorm[compactI]
653  )*compactCoarsePartialArea[compactI];
654 
655  // q emission
656  scalar qem =
657  compactCoarseEave[compactI]
658  *physicoChemical::sigma.value()
659  *compactCoarseTave[compactI];
660 
661  // compactCoarseEave is the diffussive reflected coeff
662  scalar qDiff = (compactCoarseEave[compactI])*qin;
663 
664  localqDiffusive[locaFaceI] += (qDiff)*aAve*vf[j];
665  }
666  locaFaceI++;
667  }
668  }
669 
670  volScalarField::Boundary& qsBf = qsecondRad_.boundaryFieldRef();
671 
672  // Fill qsecondRad_
673  label compactId = 0;
674  forAll(includePatches_, i)
675  {
676  const label patchID = includePatches_[i];
677  const polyPatch& pp = coarseMesh_->boundaryMesh()[patchID];
678 
679  if (pp.size() > 0)
680  {
681  scalarField& qrp = qsBf[patchID];
682 
683  const labelList& coarsePatchFace =
684  coarseMesh_->patchFaceMap()[patchID];
685 
686  forAll(pp, coarseI)
687  {
688  const label coarseFaceID = coarsePatchFace[coarseI];
689  const labelList& fineFaces = coarseToFine_[i][coarseFaceID];
690  forAll(fineFaces, k)
691  {
692  label faceI = fineFaces[k];
693  qrp[faceI] = localqDiffusive[compactId];
694  }
695  compactId ++;
696  }
697  }
698  }
699 
700  const scalarField& V = mesh_.V();
701  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
702  volScalarField::Boundary& qrBf = qr_.boundaryFieldRef();
703 
704  for (const label patchID : includePatches)
705  {
706  const scalarField& qSecond = qsecondRad_.boundaryField()[patchID];
707  if (includeMappedPatchBasePatches[patchID])
708  {
709  qrBf[patchID] += qSecond;
710  }
711  else
712  {
713  const polyPatch& pp = patches[patchID];
714  const labelList& cellIds = pp.faceCells();
715  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
716  forAll(pp, faceI)
717  {
718  const label cellI = cellIds[faceI];
719  Ru_[cellI] += qSecond[faceI]*sf[faceI]/V[cellI];
720  }
721  }
722  }
723 }
724 */
725 
726 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
727 
728 Foam::radiation::solarLoad::solarLoad(const volScalarField& T)
729 :
730  radiationModel(typeName, T),
731  solarCalc_(coeffs_, mesh_),
732  dict_(coeffs_),
733  qr_
734  (
735  IOobject
736  (
737  "qr",
738  mesh_.time().timeName(),
739  mesh_,
740  IOobject::READ_IF_PRESENT,
741  IOobject::AUTO_WRITE
742  ),
743  mesh_,
745  ),
746  hitFaces_(),
747  reflectedFaces_(),
748  Ru_
749  (
750  IOobject
751  (
752  "Ru",
753  mesh_.time().timeName(),
754  mesh_,
755  IOobject::NO_READ,
756  IOobject::NO_WRITE
757  ),
758  mesh_,
760  ),
761  absorptivity_(mesh_.boundaryMesh().size()),
762  spectralDistribution_(),
763  spectralDistributions_(),
764  qprimaryRad_(0),
765  verticalDir_(Zero),
766  nBands_(0),
767  updateTimeIndex_(0),
768  solidCoupled_(true),
769  wallCoupled_(false),
770  updateAbsorptivity_(false),
771  useReflectedRays_(false),
772  firstIter_(true)
773 {
774  initialise(coeffs_);
775 }
776 
777 
778 Foam::radiation::solarLoad::solarLoad
779 (
780  const dictionary& dict,
781  const volScalarField& T
782 )
783 :
784  radiationModel(typeName, dict, T),
785  solarCalc_(dict, mesh_),
786  dict_(dict),
787  qr_
788  (
789  IOobject
790  (
791  "qr",
792  mesh_.time().timeName(),
793  mesh_,
796  ),
797  mesh_,
799  ),
800  hitFaces_(),
801  reflectedFaces_(),
802  Ru_
803  (
804  IOobject
805  (
806  "Ru",
807  mesh_.time().timeName(),
808  mesh_,
811  ),
812  mesh_,
814  ),
815  absorptivity_(mesh_.boundaryMesh().size()),
816  spectralDistribution_(),
817  spectralDistributions_(),
818  qprimaryRad_(0),
819  verticalDir_(Zero),
820  nBands_(0),
821  updateTimeIndex_(0),
822  solidCoupled_(true),
823  wallCoupled_(false),
824  updateAbsorptivity_(false),
825  useReflectedRays_(false),
826  firstIter_(true)
827 {
828  initialise(dict);
829 }
830 
831 
832 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
833 
835 {
836  if (radiationModel::read())
837  {
838  return true;
839  }
840 
841  return false;
842 }
843 
844 
846 {
847  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
848  labelHashSet includePatches;
849  forAll(patches, patchI)
850  {
851  const polyPatch& pp = patches[patchI];
852  if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
853  {
854  includePatches.insert(patchI);
855  }
856  }
857 
858  labelHashSet includeMappedPatchBasePatches;
859  forAll(patches, patchI)
860  {
861  const polyPatch& pp = patches[patchI];
862  if
863  (
864  (isA<mappedPatchBase>(pp) && solidCoupled_)
865  || (isA<wallPolyPatch>(pp) && wallCoupled_)
866  )
867  {
868  includeMappedPatchBasePatches.insert(patchI);
869  }
870  }
871 
872  if (updateAbsorptivity_ || firstIter_)
873  {
874  updateAbsorptivity(includePatches);
875  }
876 
877  const bool facesChanged = updateHitFaces();
878 
879  const bool timeDependentLoad =
880  solarCalc_.sunLoadModel() == solarCalculator::mSunLoadTimeDependent;
881 
882  if (facesChanged || timeDependentLoad)
883  {
884  // Reset Ru
886 
887  solarCalc_.correctDirectSolarRad();
888  solarCalc_.correctDiffuseSolarRad();
889 
890  spectralDistribution_ =
891  spectralDistributions_->value(mesh_.time().value());
892 
893  spectralDistribution_ =
894  spectralDistribution_/sum(spectralDistribution_);
895 
896  // Add direct hit radiation
897  const labelList& hitFacesId = hitFaces_->rayStartFaces();
898  updateDirectHitRadiation(hitFacesId, includeMappedPatchBasePatches);
899 
900  // Add sky diffusive radiation
901  updateSkyDiffusiveRadiation
902  (
903  includePatches,
904  includeMappedPatchBasePatches
905  );
906 
907  // Add specular reflected radiation
908  if (useReflectedRays_)
909  {
910  updateReflectedRays(includeMappedPatchBasePatches);
911  }
912 
913  firstIter_ = false;
914  }
915 
916  if (debug)
917  {
918  if (mesh_.time().writeTime())
919  {
920  Ru_.write();
921  }
922  }
923 }
924 
925 
927 {
929  (
930  IOobject
931  (
932  "Rp",
933  mesh_.time().timeName(),
934  mesh_,
937  false
938  ),
939  mesh_,
941  (
943  Zero
944  )
945  );
946 }
947 
948 
951 {
952  return Ru_;
953 }
954 
955 
956 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
solarLoad.H
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
vectorList.H
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::radiation::radiationModel::coeffs_
dictionary coeffs_
Radiation model dictionary.
Definition: radiationModel.H:112
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::polyPatch::coupled
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:377
Foam::MeshObject< fvMesh, Foam::GeometricMeshObject, boundaryRadiationProperties >::New
static const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
gravityMeshObject.H
wallPolyPatch.H
Foam::solarCalculator::mSunLoadFairWeatherConditions
Definition: solarCalculator.H:465
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
surfaceFields.H
Foam::surfaceFields.
Foam::HashSet< label, Hash< label > >
Foam::radiation::solarLoad::Ru
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant)
Definition: solarLoad.C:950
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::radiation::radiationModel::read
virtual bool read()=0
Read radiationProperties dictionary.
Definition: radiationModel.C:210
Foam::dimensioned::readIfPresent
bool readIfPresent(const dictionary &dict)
Definition: dimensionedType.C:483
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::uniformDimensionedVectorField
UniformDimensionedField< vector > uniformDimensionedVectorField
Definition: uniformDimensionedFields.H:50
patchID
label patchID
Definition: boundaryProcessorFaPatchPoints.H:5
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:187
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::solarCalculator::mSunLoadTimeDependent
Definition: solarCalculator.H:464
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:26
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::solarCalculator::mSunDirConstant
Definition: solarCalculator.H:453
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::radiation::solarLoad::calculate
void calculate()
Solve radiation equations.
Definition: solarLoad.C:845
cyclicAMIPolyPatch.H
constants.H
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::solarCalculator::mSunLoadTheoreticalMaximum
Definition: solarCalculator.H:466
Foam::List< label >
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:75
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:62
addToRadiationRunTimeSelectionTables
#define addToRadiationRunTimeSelectionTables(model)
Definition: radiationModel.H:288
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
Foam::radiation::solarLoad::read
bool read()
Read radiationProperties dictionary.
Definition: solarLoad.C:834
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
boundaryRadiationProperties.H
Foam::solarCalculator::mSunDirTracking
Definition: solarCalculator.H:454
Foam::meshObjects::gravity::New
static const gravity & New(const Time &runTime)
Construct on Time.
Definition: gravityMeshObject.H:93
Foam::radiation::solarLoad::Rp
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: solarLoad.C:926
mappedPatchBase.H
Foam::solarCalculator::mSunLoadConstant
Definition: solarCalculator.H:463