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 | Copyright (C) 2015 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 "solarLoad.H"
27 #include "surfaceFields.H"
28 #include "vectorList.H"
32 #include "cyclicAMIPolyPatch.H"
33 #include "mappedPatchBase.H"
34 #include "wallPolyPatch.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  namespace radiation
41  {
42  defineTypeNameAndDebug(solarLoad, 0);
44  }
45 }
46 
48  = "viewFactorWall";
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
53 {
54  if (hitFaces_.empty())
55  {
57  return true;
58  }
59  else
60  {
61  switch (solarCalc_.sunDirectionModel())
62  {
64  {
65  return false;
66  break;
67  }
69  {
70  label updateIndex = label
71  (
73  );
74 
75  if (updateIndex > updateTimeIndex_)
76  {
77  Info << "Updating Sun position..." << endl;
78  updateTimeIndex_ = updateIndex;
80  hitFaces_->direction() = solarCalc_.direction();
81  hitFaces_->correct();
82  return true;
83  break;
84  }
85  }
86  }
87  }
88 
89  return false;
90 }
91 
92 
94 (
95  const labelHashSet& includePatches
96 )
97 {
98  const boundaryRadiationProperties& boundaryRadiation =
100 
101  forAllConstIter(labelHashSet, includePatches, iter)
102  {
103  const label patchID = iter.key();
104  absorptivity_[patchID].setSize(nBands_);
105  for (label bandI = 0; bandI < nBands_; bandI++)
106  {
107  absorptivity_[patchID][bandI] =
108  boundaryRadiation.absorptivity(patchID, bandI);
109  }
110  }
111 }
112 
113 
115 (
116  const labelList& hitFacesId,
117  const labelHashSet& includeMappedPatchBasePatches
118 )
119 {
120  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
121  const scalarField& V = mesh_.V();
122 
123  forAll (hitFacesId, i)
124  {
125  const label faceI = hitFacesId[i];
126  label patchID = patches.whichPatch(faceI);
127  const polyPatch& pp = patches[patchID];
128  const label localFaceI = faceI - pp.start();
129  const vector qPrim =
130  solarCalc_.directSolarRad()
131  * solarCalc_.direction();
132 
133  if (includeMappedPatchBasePatches[patchID])
134  {
135  const vectorField n = pp.faceNormals();
136 
137  for (label bandI = 0; bandI < nBands_; bandI++)
138  {
139  Qr_.boundaryField()[patchID][localFaceI] +=
140  (qPrim & n[localFaceI])
141  * spectralDistribution_[bandI]
142  * absorptivity_[patchID][bandI]()[localFaceI];
143  }
144  }
145  else
146  {
147  const vectorField& sf = mesh_.Sf().boundaryField()[patchID];
148  const label cellI = pp.faceCells()[localFaceI];
149 
150  for (label bandI = 0; bandI < nBands_; bandI++)
151  {
152  Ru_[cellI] +=
153  (qPrim & sf[localFaceI])
154  * spectralDistribution_[bandI]
155  * absorptivity_[patchID][bandI]()[localFaceI]
156  / V[cellI];
157  }
158  }
159  }
160 }
161 
163 (
164  const labelHashSet& includePatches,
165  const labelHashSet& includeMappedPatchBasePatches
166 )
167 {
168  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
169  const scalarField& V = mesh_.V();
170 
171  switch(solarCalc_.sunLoadModel())
172  {
175  {
176  forAllConstIter(labelHashSet, includePatches, iter)
177  {
178  const label patchID = iter.key();
179  const polyPatch& pp = patches[patchID];
180  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
181 
182  const vectorField n = pp.faceNormals();
183  const labelList& cellIds = pp.faceCells();
184 
185  forAll (n, faceI)
186  {
187  const scalar cosEpsilon(verticalDir_ & -n[faceI]);
188 
189  scalar Ed(0.0);
190  scalar Er(0.0);
191  const scalar cosTheta(solarCalc_.direction() & -n[faceI]);
192 
193  {
194  // Above the horizon
195  if (cosEpsilon == 0.0)
196  {
197  // Vertical walls
198  scalar Y(0);
199 
200  if (cosTheta > -0.2)
201  {
202  Y = 0.55+0.437*cosTheta + 0.313*sqr(cosTheta);
203  }
204  else
205  {
206  Y = 0.45;
207  }
208 
209  Ed = solarCalc_.C()*Y*solarCalc_.directSolarRad();
210  }
211  else
212  {
213  //Other than vertical walls
214  Ed =
215  solarCalc_.C()
216  * solarCalc_.directSolarRad()
217  * (1.0 + cosEpsilon)/2.0;
218  }
219 
220  // Ground reflected
221  Er =
222  solarCalc_.directSolarRad()
223  * (solarCalc_.C() + Foam::sin(solarCalc_.beta()))
224  * solarCalc_.groundReflectivity()
225  * (1.0 - cosEpsilon)/2.0;
226  }
227 
228  const label cellI = cellIds[faceI];
229  if (includeMappedPatchBasePatches[patchID])
230  {
231  for (label bandI = 0; bandI < nBands_; bandI++)
232  {
233  Qr_.boundaryField()[patchID][faceI] +=
234  (Ed + Er)
235  * spectralDistribution_[bandI]
236  * absorptivity_[patchID][bandI]()[faceI];
237  }
238  }
239  else
240  {
241  for (label bandI = 0; bandI < nBands_; bandI++)
242  {
243  Ru_[cellI] +=
244  (Ed + Er)
245  * spectralDistribution_[bandI]
246  * absorptivity_[patchID][bandI]()[faceI]
247  * sf[faceI]/V[cellI];
248  }
249  }
250  }
251  }
252  }
253  break;
254 
256  {
257  forAllConstIter(labelHashSet, includePatches, iter)
258  {
259  const label patchID = iter.key();
260  const polyPatch& pp = patches[patchID];
261  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
262 
263  const labelList& cellIds = pp.faceCells();
264  forAll (pp, faceI)
265  {
266  const label cellI = cellIds[faceI];
267  if (includeMappedPatchBasePatches[patchID])
268  {
269  for (label bandI = 0; bandI < nBands_; bandI++)
270  {
271  Qr_.boundaryField()[patchID][faceI] +=
272  solarCalc_.diffuseSolarRad()
273  * spectralDistribution_[bandI]
274  * absorptivity_[patchID][bandI]()[faceI];
275  }
276  }
277  else
278  {
279  for (label bandI = 0; bandI < nBands_; bandI++)
280  {
281  Ru_[cellI] +=
282  (
283  spectralDistribution_[bandI]
284  * absorptivity_[patchID][bandI]()[faceI]
285  * solarCalc_.diffuseSolarRad()
286  )*sf[faceI]/V[cellI];
287  }
288  }
289  }
290  }
291  break;
292  }
293  }
294 }
295 
296 
298 {
299 
300  if (coeffs.found("gridUp"))
301  {
302  coeffs.lookup("gridUp") >> verticalDir_;
303  verticalDir_ /= mag(verticalDir_);
304  }
305  else if (mesh_.foundObject<uniformDimensionedVectorField>("g"))
306  {
308  mesh_.lookupObject<uniformDimensionedVectorField>("g");
309  verticalDir_ = (-g/mag(g)).value();
310  }
311 
312  includePatches_ = mesh_.boundaryMesh().findIndices(viewFactorWalls);
313 
314  coeffs.lookup("useVFbeamToDiffuse") >> useVFbeamToDiffuse_;
315 
316  coeffs.lookup("spectralDistribution") >> spectralDistribution_;
317  spectralDistribution_ =
318  spectralDistribution_/sum(spectralDistribution_);
319 
320  nBands_ = spectralDistribution_.size();
321 
322  if (useVFbeamToDiffuse_)
323  {
324  map_.reset
325  (
326  new IOmapDistribute
327  (
328  IOobject
329  (
330  "mapDist",
331  mesh_.facesInstance(),
332  mesh_,
335  )
336  )
337  );
338  }
339 
340  if (coeffs.found("solidCoupled"))
341  {
342  coeffs.lookup("solidCoupled") >> solidCoupled_;
343  }
344 
345  if (coeffs.found("wallCoupled"))
346  {
347  coeffs.lookup("wallCoupled") >> wallCoupled_;
348  }
349 
350  if (coeffs.found("updateAbsorptivity"))
351  {
352  coeffs.lookup("updateAbsorptivity") >> updateAbsorptivity_;
353  }
354 }
355 
356 
358 (
359  const labelHashSet& includePatches,
360  const labelHashSet& includeMappedPatchBasePatches
361 )
362 {
363 
364  scalarListIOList FmyProc
365  (
366  IOobject
367  (
368  "F",
369  mesh_.facesInstance(),
370  mesh_,
373  false
374  )
375  );
376 
377 
378  if (finalAgglom_.size() > 0 && coarseMesh_.empty())
379  {
380  coarseMesh_.reset
381  (
382  new singleCellFvMesh
383  (
384  IOobject
385  (
386  "coarse." + mesh_.name(),
387  mesh_.polyMesh::instance(),
388  mesh_.time(),
391  ),
392  mesh_,
393  finalAgglom_
394  )
395  );
396  }
397 
398  label nLocalVFCoarseFaces = 0;
399  forAll(includePatches_, i)
400  {
401  const label patchI = includePatches_[i];
402  nLocalVFCoarseFaces +=
403  coarseMesh_->boundaryMesh()[patchI].size();
404  }
405 
406  label totalFVNCoarseFaces = nLocalVFCoarseFaces;
407  reduce(totalFVNCoarseFaces, sumOp<label>());
408 
409  // Calculate weighted absorptivity on coarse patches
410  List<scalar> localCoarseRave(nLocalVFCoarseFaces);
411  List<scalar> localCoarsePartialArea(nLocalVFCoarseFaces);
412  List<vector> localCoarseNorm(nLocalVFCoarseFaces);
413 
414  scalarField compactCoarseRave(map_->constructSize(), 0.0);
415  scalarField compactCoarsePartialArea(map_->constructSize(), 0.0);
416  vectorList compactCoarseNorm(map_->constructSize(), vector::zero);
417 
418  const boundaryRadiationProperties& boundaryRadiation =
420 
421  coarseToFine_.setSize(includePatches_.size());
422 
423  const labelList& hitFacesId = hitFaces_->rayStartFaces();
424 
425  label startI = 0;
426  label compactI = 0;
427  forAll (includePatches_, i)
428  {
429  const label patchID = includePatches_[i];
430  const polyPatch& pp = mesh_.boundaryMesh()[patchID];
431 
432  const polyPatch& cpp = coarseMesh_->boundaryMesh()[patchID];
433 
434  const labelList& agglom = finalAgglom_[patchID];
435  //if (pp.size() > 0)
436  if (agglom.size() > 0)
437  {
438  label nAgglom = max(agglom) + 1;
439  coarseToFine_[i] = invertOneToMany(nAgglom, agglom);
440  }
441 
442  scalarField r(pp.size(), 0.0);
443  for (label bandI = 0; bandI < nBands_; bandI++)
444  {
445  const tmp<scalarField> tr =
446  spectralDistribution_[bandI]
447  *boundaryRadiation.reflectivity(patchID, bandI);
448  r += tr();
449  }
450 
451  scalarList Rave(cpp.size(), 0.0);
452  scalarList area(cpp.size(), 0.0);
453 
454  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
455 
456  const labelList& coarsePatchFace =
457  coarseMesh_->patchFaceMap()[patchID];
458 
459  forAll (cpp, coarseI)
460  {
461  const label coarseFaceID = coarsePatchFace[coarseI];
462  const labelList& fineFaces =
463  coarseToFine_[i][coarseFaceID];
464 
465  UIndirectList<scalar> fineSf
466  (
467  sf,
468  fineFaces
469  );
470  scalar fineArea = sum(fineSf());
471 
472  scalar fullArea = 0.0;
473  forAll(fineFaces, j)
474  {
475  label faceI = fineFaces[j];
476  label globalFaceI = faceI + pp.start();
477 
478  if (findIndex(hitFacesId, globalFaceI) != -1)
479  {
480  fullArea += sf[faceI];
481  }
482  Rave[coarseI] += (r[faceI]*sf[faceI])/fineArea;
483  }
484  localCoarsePartialArea[compactI++] = fullArea/fineArea;
485  }
486 
488  (
489  localCoarseRave,
490  Rave.size(),
491  startI
492  ).assign(Rave);
493 
494 
495  const vectorList coarseNSf = cpp.faceNormals();
497  (
498  localCoarseNorm,
499  cpp.size(),
500  startI
501  ).assign(coarseNSf);
502  startI += cpp.size();
503  }
504 
505 
506  SubList<scalar>(compactCoarsePartialArea, nLocalVFCoarseFaces).assign
507  (
508  localCoarsePartialArea
509  );
510 
511  SubList<scalar>(compactCoarseRave, nLocalVFCoarseFaces).assign
512  (
513  localCoarseRave
514  );
515 
516  SubList<vector>(compactCoarseNorm, nLocalVFCoarseFaces).assign
517  (
518  localCoarseNorm
519  );
520 
521  map_->distribute(compactCoarsePartialArea);
522  map_->distribute(compactCoarseRave);
523  map_->distribute(compactCoarseNorm);
524 
525 
526  // Calculate coarse hitFaces and Sun direct hit heat fluxes
527  scalarList localqDiffusive(nLocalVFCoarseFaces, 0.0);
528 
529  label locaFaceI = 0;
530  forAll (includePatches_, i)
531  {
532  const label patchID = includePatches_[i];
533  const polyPatch& pp = coarseMesh_->boundaryMesh()[patchID];
534  const polyPatch& ppf = mesh_.boundaryMesh()[patchID];
535 
536  const labelList& coarsePatchFace =
537  coarseMesh_->patchFaceMap()[patchID];
538  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
539 
540  scalarField a(ppf.size(), 0.0);
541  for (label bandI = 0; bandI < nBands_; bandI++)
542  {
543  const tmp<scalarField> ta =
544  spectralDistribution_[bandI]
545  * absorptivity_[patchID][bandI]();
546  a += ta();
547  }
548 
549  forAll (pp, coarseI)
550  {
551  const label coarseFaceID = coarsePatchFace[coarseI];
552  const labelList& fineFaces = coarseToFine_[i][coarseFaceID];
553  UIndirectList<scalar> fineSf
554  (
555  sf,
556  fineFaces
557  );
558  scalar fineArea = sum(fineSf());
559 
560  scalar aAve = 0.0;
561  forAll(fineFaces, j)
562  {
563  label faceI = fineFaces[j];
564  aAve += (a[faceI]*sf[faceI])/fineArea;
565  }
566 
567  const scalarList& vf = FmyProc[locaFaceI];
568 
569  const labelList& compactFaces = visibleFaceFaces_[locaFaceI];
570 
571 
572  forAll (compactFaces, j)
573  {
574  label compactI = compactFaces[j];
575 
576  localqDiffusive[locaFaceI] +=
577  compactCoarsePartialArea[compactI]
578  * aAve
579  * (solarCalc_.directSolarRad()*solarCalc_.direction())
580  & compactCoarseNorm[compactI]
581  * vf[j]
582  * compactCoarseRave[compactI];
583 
584  }
585  locaFaceI++;
586  }
587  }
588 
589  // Fill QsecondRad_
590  label compactId = 0;
591  forAll (includePatches_, i)
592  {
593  const label patchID = includePatches_[i];
594  const polyPatch& pp = coarseMesh_->boundaryMesh()[patchID];
595 
596  if (pp.size() > 0)
597  {
598  scalarField& Qrp = QsecondRad_.boundaryField()[patchID];
599 
600  const labelList& coarsePatchFace =
601  coarseMesh_->patchFaceMap()[patchID];
602 
603  forAll(pp, coarseI)
604  {
605  const label coarseFaceID = coarsePatchFace[coarseI];
606  const labelList& fineFaces = coarseToFine_[i][coarseFaceID];
607  forAll(fineFaces, k)
608  {
609  label faceI = fineFaces[k];
610  Qrp[faceI] = localqDiffusive[compactId];
611  }
612  compactId ++;
613  }
614  }
615  }
616 
617  const scalarField& V = mesh_.V();
618  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
619 
620  forAllConstIter(labelHashSet, includePatches, iter)
621  {
622  const label patchID = iter.key();
623  const scalarField& qSecond = QsecondRad_.boundaryField()[patchID];
624  if (includeMappedPatchBasePatches[patchID])
625  {
626  Qr_.boundaryField()[patchID] += qSecond;
627  }
628  else
629  {
630  const polyPatch& pp = patches[patchID];
631  const labelList& cellIds = pp.faceCells();
632  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
633  forAll (pp, faceI)
634  {
635  const label cellI = cellIds[faceI];
636  Ru_[cellI] += qSecond[faceI]*sf[faceI]/V[cellI];
637  }
638  }
639  }
640 }
641 
642 
643 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
644 
646 :
647  radiationModel(typeName, T),
648  finalAgglom_
649  (
650  IOobject
651  (
652  "finalAgglom",
653  mesh_.facesInstance(),
654  mesh_,
655  IOobject::READ_IF_PRESENT,
656  IOobject::NO_WRITE,
657  false
658  )
659  ),
660  coarseMesh_(),
661  Qr_
662  (
663  IOobject
664  (
665  "Qr",
666  mesh_.time().timeName(),
667  mesh_,
668  IOobject::READ_IF_PRESENT,
669  IOobject::AUTO_WRITE
670  ),
671  mesh_,
673  ),
674  QsecondRad_
675  (
676  IOobject
677  (
678  "QsecondRad",
679  mesh_.time().timeName(),
680  mesh_,
681  IOobject::NO_READ,
682  IOobject::AUTO_WRITE
683  ),
684  mesh_,
685  dimensionedScalar("QsecondRad", dimMass/pow3(dimTime), 0.0)
686  ),
687  hitFaces_(),
688  Ru_
689  (
690  IOobject
691  (
692  "Ru",
693  mesh_.time().timeName(),
694  mesh_,
695  IOobject::NO_READ,
696  IOobject::NO_WRITE
697  ),
698  mesh_,
700  ),
701  solarCalc_(this->subDict(typeName + "Coeffs"), mesh_),
702  verticalDir_(vector::zero),
703  useVFbeamToDiffuse_(false),
704  includePatches_(mesh_.boundary().size(), -1),
705  coarseToFine_(),
706  nBands_(1),
707  spectralDistribution_(nBands_),
708  map_(),
709  visibleFaceFaces_
710  (
711  IOobject
712  (
713  "visibleFaceFaces",
714  mesh_.facesInstance(),
715  mesh_,
716  IOobject::READ_IF_PRESENT,
717  IOobject::NO_WRITE,
718  false
719  )
720  ),
721  solidCoupled_(true),
722  absorptivity_(mesh_.boundaryMesh().size()),
723  updateAbsorptivity_(false),
724  firstIter_(true),
725  updateTimeIndex_(0)
726 {
728 }
729 
730 
732 (
733  const dictionary& dict,
734  const volScalarField& T
735 )
736 :
737  radiationModel(typeName, dict, T),
738  finalAgglom_
739  (
740  IOobject
741  (
742  "finalAgglom",
743  mesh_.facesInstance(),
744  mesh_,
747  false
748  )
749  ),
750  coarseMesh_(),
751  Qr_
752  (
753  IOobject
754  (
755  "Qr",
756  mesh_.time().timeName(),
757  mesh_,
760  ),
761  mesh_,
763  ),
764  QsecondRad_
765  (
766  IOobject
767  (
768  "QsecondRad",
769  mesh_.time().timeName(),
770  mesh_,
773  ),
774  mesh_,
775  dimensionedScalar("QsecondRad", dimMass/pow3(dimTime), 0.0)
776  ),
777  hitFaces_(),
778  Ru_
779  (
780  IOobject
781  (
782  "Ru",
783  mesh_.time().timeName(),
784  mesh_,
787  ),
788  mesh_,
790  ),
791  solarCalc_(coeffs_, mesh_),
792  verticalDir_(vector::zero),
793  useVFbeamToDiffuse_(false),
794  includePatches_(mesh_.boundary().size(), -1),
795  coarseToFine_(),
796  nBands_(1),
797  spectralDistribution_(nBands_),
798  map_(),
799  visibleFaceFaces_
800  (
801  IOobject
802  (
803  "visibleFaceFaces",
804  mesh_.facesInstance(),
805  mesh_,
808  false
809  )
810  ),
811  solidCoupled_(true),
812  wallCoupled_(false),
813  absorptivity_(mesh_.boundaryMesh().size()),
814  updateAbsorptivity_(false),
815  firstIter_(true),
816  updateTimeIndex_(0)
817 {
818  initialise(coeffs_);
819 }
820 
821 
823 (
824  const dictionary& dict,
825  const volScalarField& T,
826  const word radWallFieldName
827 )
828 :
829  radiationModel("none", T),
830  finalAgglom_
831  (
832  IOobject
833  (
834  "finalAgglom",
835  mesh_.facesInstance(),
836  mesh_,
839  false
840  )
841  ),
842  coarseMesh_(),
843  Qr_
844  (
845  IOobject
846  (
847  radWallFieldName,
848  mesh_.time().timeName(),
849  mesh_,
852  ),
853  mesh_,
855  ),
856  QsecondRad_
857  (
858  IOobject
859  (
860  "QsecondRad",
861  mesh_.time().timeName(),
862  mesh_,
865  ),
866  mesh_,
867  dimensionedScalar("QsecondRad", dimMass/pow3(dimTime), 0.0)
868  ),
869  hitFaces_(),
870  Ru_
871  (
872  IOobject
873  (
874  "Ru",
875  mesh_.time().timeName(),
876  mesh_,
879  ),
880  mesh_,
882  ),
883  solarCalc_(dict, mesh_),
884  verticalDir_(vector::zero),
885  useVFbeamToDiffuse_(false),
886  includePatches_(mesh_.boundary().size(), -1),
887  coarseToFine_(),
888  nBands_(1),
889  spectralDistribution_(nBands_),
890  map_(),
891  visibleFaceFaces_
892  (
893  IOobject
894  (
895  "visibleFaceFaces",
896  mesh_.facesInstance(),
897  mesh_,
900  false
901  )
902  ),
903  solidCoupled_(true),
904  wallCoupled_(false),
905  absorptivity_(mesh_.boundaryMesh().size()),
906  updateAbsorptivity_(false),
907  firstIter_(true)
908 {
909  initialise(dict);
910 }
911 
912 
913 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
914 
916 {}
917 
918 
919 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
920 
922 {
923  return nBands_;
924 }
925 
926 
928 {
929  if (radiationModel::read())
930  {
931  return true;
932  }
933  else
934  {
935  return false;
936  }
937 }
938 
939 
941 {
942  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
943  labelHashSet includePatches;
944  forAll(patches, patchI)
945  {
946  const polyPatch& pp = patches[patchI];
947  if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
948  {
949  includePatches.insert(patchI);
950  }
951  }
952 
953  labelHashSet includeMappedPatchBasePatches;
954  forAll(patches, patchI)
955  {
956  const polyPatch& pp = patches[patchI];
957  if
958  (
959  (isA<mappedPatchBase>(pp) && solidCoupled_)
960  || (isA<wallPolyPatch>(pp) && wallCoupled_)
961  )
962  {
963  includeMappedPatchBasePatches.insert(patchI);
964  }
965  }
966 
967  if (updateAbsorptivity_ || firstIter_)
968  {
969  updateAbsorptivity(includePatches);
970  }
971 
972  bool facesChanged = updateHitFaces();
973 
974  if (facesChanged)
975  {
976  // Reset Ru and Qr
977  Ru_ = dimensionedScalar("Ru", dimMass/dimLength/pow3(dimTime), 0.0);
978  Qr_.boundaryField() = 0.0;
979 
980  // Add direct hit radation
981  const labelList& hitFacesId = hitFaces_->rayStartFaces();
982  updateDirectHitRadiation(hitFacesId, includeMappedPatchBasePatches);
983 
984  // Add sky diffusive radiation
985  updateSkyDiffusiveRadiation
986  (
987  includePatches,
988  includeMappedPatchBasePatches
989  );
990 
991  // Add indirect diffusive radiation
992  if (useVFbeamToDiffuse_)
993  {
994  calculateQdiff(includePatches, includeMappedPatchBasePatches);
995  }
996 
997  firstIter_ = false;
998  }
999 
1000  if (debug)
1001  {
1002  if (mesh_.time().outputTime())
1003  {
1004  Ru_.write();
1005  }
1006  }
1007 }
1008 
1009 
1011 {
1012  return tmp<volScalarField>
1013  (
1014  new volScalarField
1015  (
1016  IOobject
1017  (
1018  "Rp",
1019  mesh_.time().timeName(),
1020  mesh_,
1023  false
1024  ),
1025  mesh_,
1027  (
1028  "zero",
1030  0.0
1031  )
1032  )
1033  );
1034 }
1035 
1036 
1039 {
1040  return Ru_;
1041 }
1042 
1043 // ************************************************************************* //
Foam::radiation::solarLoad::nBands
label nBands() const
Number of bands.
Definition: solarLoad.C:921
solarLoad.H
Foam::radiation::solarLoad::updateAbsorptivity
void updateAbsorptivity(const labelHashSet &includePatches)
Update absorptivity.
Definition: solarLoad.C:94
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::solarCalculator::direction
const vector direction() const
const acess to direction
Definition: solarCalculator.H:231
Foam::radiation::boundaryRadiationProperties::absorptivity
tmp< scalarField > absorptivity(const label patchId, const label bandI=0) const
Access boundary absorptivity on patch.
Definition: boundaryRadiationProperties.C:133
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::solarCalculator::sunTrackingUpdateInterval
scalar sunTrackingUpdateInterval()
Return sunTrackingUpdateInterval.
Definition: solarCalculator.H:297
vectorList.H
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::solarCalculator::mSunDirTraking
@ mSunDirTraking
Definition: solarCalculator.H:111
Foam::radiation::radiationModel::coeffs_
dictionary coeffs_
Radiation model dictionary.
Definition: radiationModel.H:103
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::findIndex
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
Foam::invertOneToMany
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
Foam::IOmapDistribute
IOmapDistribute is derived from mapDistribute and IOobject to give the mapDistribute automatic IO fun...
Definition: IOmapDistribute.H:51
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:255
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::polyPatch::coupled
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:322
boundary
faceListList boundary(nPatches)
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
wallPolyPatch.H
Foam::radiation::solarLoad::updateHitFaces
bool updateHitFaces()
Update hit faces.
Definition: solarLoad.C:52
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::solarCalculator::mSunLoadFairWeatherConditions
@ mSunLoadFairWeatherConditions
Definition: solarCalculator.H:118
g
const dimensionedVector & g
Definition: setRegionFluidFields.H:33
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::radiation::solarLoad::solarCalc_
solarCalculator solarCalc_
Solar calculator.
Definition: solarLoad.H:116
Foam::HashSet< label, Hash< label > >
Foam::radiation::solarLoad::Ru
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant)
Definition: solarLoad.C:1038
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::MeshObject< fvMesh, Foam::GeometricMeshObject, boundaryRadiationProperties >::New
static const boundaryRadiationProperties & New(const fvMesh &mesh)
Definition: MeshObject.C:44
Foam::radiation::solarLoad::calculateQdiff
void calculateQdiff(const labelHashSet &, const labelHashSet &)
Calculate diffusive heat flux.
Definition: solarLoad.C:358
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::radiation::boundaryRadiationProperties::reflectivity
tmp< scalarField > reflectivity(const label patchId, const label bandI=0) const
Access boundary reflectivity on patch.
Definition: boundaryRadiationProperties.C:179
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::singleCellFvMesh
fvMesh as subset of other mesh. Consists of one cell and all original bounday faces....
Definition: singleCellFvMesh.H:53
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:98
Foam::radiation::radiationModel::read
virtual bool read()=0
Read radiationProperties dictionary.
Definition: radiationModel.C:197
Foam::UniformDimensionedField< vector >
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:59
Foam::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::polyPatch::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:302
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::UList::assign
void assign(const UList< T > &)
Assign elements to those from UList.
Definition: UList.C:37
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::radiation::solarLoad::updateSkyDiffusiveRadiation
void updateSkyDiffusiveRadiation(const labelHashSet &, const labelHashSet &)
Update Sky diffusive radiation.
Definition: solarLoad.C:163
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::faceShading
faceShading uses the transmissivity value in the boundaryRadiationProperties in order to evaluate whi...
Definition: faceShading.H:58
Foam::solarCalculator::sunDirectionModel
sunDirModel sunDirectionModel() const
Return Sun direction model.
Definition: solarCalculator.H:273
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Foam::radiation::boundaryRadiationProperties
Definition: boundaryRadiationProperties.H:56
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::solarCalculator::mSunDirConstant
@ mSunDirConstant
Definition: solarCalculator.H:110
Foam::polyPatch::faceCells
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:340
sf
volScalarField sf(fieldObject, mesh)
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::radiation::solarLoad::calculate
void calculate()
Solve.
Definition: solarLoad.C:940
Foam::radiation::radiationModel::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: radiationModel.H:91
cyclicAMIPolyPatch.H
Foam::radiation::defineTypeNameAndDebug
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
Foam::radiation::solarLoad::solarLoad
solarLoad(const solarLoad &)
Disallow default bitwise copy construct.
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
uniformDimensionedFields.H
Foam::radiation::addToRadiationRunTimeSelectionTables
addToRadiationRunTimeSelectionTables(fvDOM)
Foam::solarCalculator::correctSunDirection
void correctSunDirection()
Recalculate.
Definition: solarCalculator.C:298
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::sumOp
Definition: ops.H:162
Foam::PrimitivePatch::faceNormals
const Field< PointType > & faceNormals() const
Return face normals for patch.
Definition: PrimitivePatchTemplate.C:520
Foam::solarCalculator::mSunLoadTheoreticalMaximum
@ mSunLoadTheoreticalMaximum
Definition: solarCalculator.H:119
Foam::Vector< scalar >
Foam::radiation::solarLoad::~solarLoad
virtual ~solarLoad()
Destructor.
Definition: solarLoad.C:915
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::radiation::solarLoad::updateTimeIndex_
label updateTimeIndex_
Update Sun position index.
Definition: solarLoad.H:158
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:73
Foam::radiation::radiationModel::initialise
void initialise()
Initialise.
Definition: radiationModel.C:77
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::IOList
A List of objects of type <T> with automated input and output.
Definition: IOList.H:50
Foam::radiation::solarLoad::hitFaces_
autoPtr< faceShading > hitFaces_
Direct hit faces Ids.
Definition: solarLoad.H:110
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:59
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::radiation::solarLoad::viewFactorWalls
static const word viewFactorWalls
Static name for view factor walls.
Definition: solarLoad.H:89
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:49
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::radiation::solarLoad::read
bool read()
Read radiation properties dictionary.
Definition: solarLoad.C:927
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::IOobject::READ_IF_PRESENT
@ READ_IF_PRESENT
Definition: IOobject.H:110
Foam::radiation::solarLoad::updateDirectHitRadiation
void updateDirectHitRadiation(const labelList &, const labelHashSet &)
Update direct hit faces radiation.
Definition: solarLoad.C:115
boundaryRadiationProperties.H
Y
PtrList< volScalarField > & Y
Definition: createFields.H:36
Foam::radiation::solarLoad::Rp
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: solarLoad.C:1010
Foam::zero
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:47
mappedPatchBase.H
Foam::solarCalculator::mSunLoadConstant
@ mSunLoadConstant
Definition: solarCalculator.H:117