Loading...
Searching...
No Matches
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-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
40using namespace Foam::constant;
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
44namespace Foam
45{
46 namespace radiation
47 {
50 }
51}
52
53// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54
55void Foam::radiation::solarLoad::updateReflectedRays
56(
57 const labelHashSet& includePatches
58)
59{
60 if (!reflectedFaces_ && hitFaces_)
61 {
62 reflectedFaces_.reset
63 (
65 (
66 mesh_,
67 hitFaces_(),
68 solarCalc_,
69 spectralDistribution_,
70 dict_
71 )
72 );
73 }
74
75 reflectedFaces_->correct();
76
78 const scalarField& V = mesh_.V();
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
114bool Foam::radiation::solarLoad::updateHitFaces()
115{
116 if (!hitFaces_)
117 {
118 const boundaryRadiationProperties& boundaryRadiation =
120
121 const labelList activeFaceZoneIds(boundaryRadiation.faceZoneIds());
122
123 if (!activeFaceZoneIds.size())
124 {
125 hitFaces_.reset(new faceShading(mesh_, solarCalc_.direction()));
126 }
127 else
128 {
129 hitFaces_.reset
130 (
131 new faceShading
132 (
133 mesh_,
135 activeFaceZoneIds,
136 solarCalc_.direction()
137 )
138 );
139 }
140
141 return true;
142 }
143 else
144 {
145 switch (solarCalc_.sunDirectionModel())
146 {
148 {
149 return false;
150 break;
151 }
153 {
154 const label updateIndex = label
155 (
156 mesh_.time().value()/solarCalc_.sunTrackingUpdateInterval()
157 );
158
159 if (updateIndex > updateTimeIndex_)
160 {
161 Info << "Updating Sun position..." << endl;
162 updateTimeIndex_ = updateIndex;
163 solarCalc_.correctSunDirection();
164 hitFaces_->direction() = solarCalc_.direction();
165 hitFaces_->correct();
166 return true;
167 }
168 break;
169 }
170 }
171 }
172
173 return false;
174}
175
176
177void Foam::radiation::solarLoad::updateAbsorptivity
178(
179 const labelHashSet& includePatches
180)
181{
182 const boundaryRadiationProperties& boundaryRadiation =
184
185 for (const label patchID : includePatches)
186 {
187 absorptivity_[patchID].setSize(nBands_);
188 for (label bandi = 0; bandi < nBands_; ++bandi)
189 {
190 absorptivity_[patchID][bandi] =
191 boundaryRadiation.absorptivity(patchID, bandi);
192 }
193 }
194}
195
196
197void Foam::radiation::solarLoad::updateDirectHitRadiation
198(
199 const labelList& hitFacesId,
200 const labelHashSet& includeMappedPatchBasePatches
201)
202{
203 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
204 const scalarField& V = mesh_.V();
205 volScalarField::Boundary& qrBf = qr_.boundaryFieldRef();
206
207 // Reset qr and qrPrimary
208 qrBf = 0.0;
209
210 for (label bandi = 0; bandi < nBands_; ++bandi)
211 {
212 volScalarField::Boundary& qprimaryBf =
213 qprimaryRad_[bandi].boundaryFieldRef();
214
215 qprimaryBf = 0.0;
216
217 forAll(hitFacesId, i)
218 {
219 const label facei = hitFacesId[i];
220 const label patchID = patches.whichPatch(facei);
221
222 if (patchID == -1)
223 {
224 continue;
225 }
226
227 const polyPatch& pp = patches[patchID];
228 const label localFaceI = facei - pp.start();
229 const vector qPrim
230 (
231 solarCalc_.directSolarRad()*solarCalc_.direction()
232 );
233
234 const vectorField& n = pp.faceNormals();
235
236 {
237 qprimaryBf[patchID][localFaceI] +=
238 (qPrim & n[localFaceI])
239 * spectralDistribution_[bandi]
240 * absorptivity_[patchID][bandi]()[localFaceI];
241 }
242
243 if (includeMappedPatchBasePatches[patchID])
244 {
245 qrBf[patchID][localFaceI] += qprimaryBf[patchID][localFaceI];
246 }
247 else
248 {
249 const vectorField& sf = mesh_.Sf().boundaryField()[patchID];
250 const label celli = pp.faceCells()[localFaceI];
251
252 Ru_[celli] +=
253 (qPrim & sf[localFaceI])
254 * spectralDistribution_[bandi]
255 * absorptivity_[patchID][bandi]()[localFaceI]
256 / V[celli];
257 }
258 }
259 }
260}
261
262
263void Foam::radiation::solarLoad::updateSkyDiffusiveRadiation
264(
265 const labelHashSet& includePatches,
266 const labelHashSet& includeMappedPatchBasePatches
267)
268{
269 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
270 const scalarField& V = mesh_.V();
271 volScalarField::Boundary& qrBf = qr_.boundaryFieldRef();
272
273 switch(solarCalc_.sunLoadModel())
274 {
277 {
278 for (const label patchID : includePatches)
279 {
280 const polyPatch& pp = patches[patchID];
281 const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
282
283 const vectorField n = pp.faceNormals();
284 const labelUList& cellids = pp.faceCells();
285
286 // Calculate diffusive radiance
287 // contribution from sky and ground
288 tmp<scalarField> tdiffuseSolarRad =
289 solarCalc_.diffuseSolarRad(n);
290 const scalarField& diffuseSolarRad = tdiffuseSolarRad.cref();
291
292 forAll(n, facei)
293 {
294 const label celli = cellids[facei];
295 if (includeMappedPatchBasePatches[patchID])
296 {
297 for (label bandi = 0; bandi < nBands_; ++bandi)
298 {
299 qrBf[patchID][facei] +=
300 diffuseSolarRad[facei]
301 * spectralDistribution_[bandi]
302 * absorptivity_[patchID][bandi]()[facei];
303 }
304 }
305 else
306 {
307 for (label bandi = 0; bandi < nBands_; ++bandi)
308 {
309 Ru_[celli] +=
310 diffuseSolarRad[facei]
311 * spectralDistribution_[bandi]
312 * absorptivity_[patchID][bandi]()[facei]
313 * sf[facei]/V[celli];
314 }
315 }
316 }
317 }
318 }
319 break;
320
323 {
324 for (const label patchID : includePatches)
325 {
326 const polyPatch& pp = patches[patchID];
327 const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
328
329 const labelUList& cellids = pp.faceCells();
330 forAll(pp, facei)
331 {
332 const label celli = cellids[facei];
333 if (includeMappedPatchBasePatches[patchID])
334 {
335 for (label bandi = 0; bandi < nBands_; ++bandi)
336 {
337 qrBf[patchID][facei] +=
338 solarCalc_.diffuseSolarRad()
339 * spectralDistribution_[bandi]
340 * absorptivity_[patchID][bandi]()[facei];
341 }
342 }
343 else
344 {
345 for (label bandi = 0; bandi < nBands_; ++bandi)
346 {
347 Ru_[celli] +=
348 (
349 spectralDistribution_[bandi]
350 * absorptivity_[patchID][bandi]()[facei]
351 * solarCalc_.diffuseSolarRad()
352 )*sf[facei]/V[celli];
353 }
354 }
355 }
356 }
357 break;
358 }
359 }
360}
361
362
363void Foam::radiation::solarLoad::initialise(const dictionary& coeffs)
364{
365 spectralDistributions_.reset
366 (
367 Function1<scalarField>::New
368 (
369 "spectralDistribution",
370 coeffs,
371 &mesh_
372 )
373 );
374
375 spectralDistribution_ =
376 spectralDistributions_->value(mesh_.time().timeOutputValue());
377
378 nBands_ = spectralDistribution_.size();
379
380 spectralDistribution_ =
381 spectralDistribution_/sum(spectralDistribution_);
382
383 qprimaryRad_.setSize(nBands_);
384
385 forAll(qprimaryRad_, bandi)
386 {
387 qprimaryRad_.set
388 (
389 bandi,
391 (
392 IOobject
393 (
394 "qprimaryRad_" + Foam::name(bandi) ,
395 mesh_.time().timeName(),
396 mesh_,
399 ),
400 mesh_,
402 )
403 );
404 }
405
406 coeffs.readIfPresent("solidCoupled", solidCoupled_);
407 coeffs.readIfPresent("wallCoupled", wallCoupled_);
408 coeffs.readIfPresent("updateAbsorptivity", updateAbsorptivity_);
409 coeffs.readEntry("useReflectedRays", useReflectedRays_);
410
411 (void) updateHitFaces();
412}
413
414/*
415void Foam::radiation::solarLoad::calculateQdiff
416(
417 const labelHashSet& includePatches,
418 const labelHashSet& includeMappedPatchBasePatches
419)
420{
421 scalarListIOList FmyProc
422 (
423 IOobject
424 (
425 "F",
426 mesh_.facesInstance(),
427 mesh_,
428 IOobject::MUST_READ,
429 IOobject::NO_WRITE,
430 IOobject::NO_REGISTER
431 )
432 );
433
434
435 if (!coarseMesh_ && !finalAgglom_.empty())
436 {
437 coarseMesh_.reset
438 (
439 new singleCellFvMesh
440 (
441 IOobject
442 (
443 "coarse:" + mesh_.name(),
444 mesh_.polyMesh::instance(),
445 mesh_.time(),
446 IOobject::NO_READ,
447 IOobject::NO_WRITE
448 ),
449 mesh_,
450 finalAgglom_
451 )
452 );
453 }
454
455 label nLocalVFCoarseFaces = 0;
456 forAll(includePatches_, i)
457 {
458 const label patchI = includePatches_[i];
459 nLocalVFCoarseFaces += coarseMesh_->boundaryMesh()[patchI].size();
460 }
461
462 const label totalFVNCoarseFaces =
463 returnReduce(nLocalVFCoarseFaces, sumOp<label>());
464
465 // Calculate weighted absorptivity on coarse patches
466 List<scalar> localCoarseEave(nLocalVFCoarseFaces);
467 List<scalar> localTave(nLocalVFCoarseFaces);
468 List<scalar> localCoarsePartialArea(nLocalVFCoarseFaces);
469 List<vector> localCoarseNorm(nLocalVFCoarseFaces);
470
471 scalarField compactCoarseEave(map_->constructSize(), 0.0);
472 scalarField compactCoarseTave(map_->constructSize(), 0.0);
473
474 scalarField compactCoarsePartialArea(map_->constructSize(), 0.0);
475
476 vectorList compactCoarseNorm(map_->constructSize(), Zero);
477
478 const boundaryRadiationProperties& boundaryRadiation =
479 boundaryRadiationProperties::New(mesh_);
480
481 coarseToFine_.setSize(includePatches_.size());
482
483 const labelList& hitFacesId = hitFaces_->rayStartFaces();
484
485 label startI = 0;
486 label compactI = 0;
487 forAll(includePatches_, i)
488 {
489 const label patchID = includePatches_[i];
490 const polyPatch& pp = mesh_.boundaryMesh()[patchID];
491
492 const polyPatch& cpp = coarseMesh_->boundaryMesh()[patchID];
493
494 const labelList& agglom = finalAgglom_[patchID];
495
496 if (agglom.size() > 0)
497 {
498 label nAgglom = max(agglom) + 1;
499 coarseToFine_[i] = invertOneToMany(nAgglom, agglom);
500 }
501
502 // Weight emissivity by spectral distribution
503 scalarField e(pp.size(), 0.0);
504
505 for (label bandi = 0; bandi < nBands_; bandi++)
506 {
507 const tmp<scalarField> te =
508 spectralDistribution_[bandi]
509 *boundaryRadiation.diffReflectivity(patchID, bandi);
510 e += te();
511 }
512
513 scalarList Eave(cpp.size(), 0.0);
514 scalarList Tave(cpp.size(), 0.0);
515
516 const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
517 const scalarField& Tf = T_.boundaryField()[patchID];
518
519 const labelList& coarsePatchFace=coarseMesh_->patchFaceMap()[patchID];
520
521 forAll(cpp, coarseI)
522 {
523 const label coarseFaceID = coarsePatchFace[coarseI];
524 const labelList& fineFaces = coarseToFine_[i][coarseFaceID];
525
526 UIndirectList<scalar> fineSf(sf, fineFaces);
527 scalar fineArea = sum(fineSf());
528
529 scalar fullArea = 0.0;
530 forAll(fineFaces, j)
531 {
532 label facei = fineFaces[j];
533 label globalFaceI = facei + pp.start();
534
535 if (hitFacesId.found(globalFaceI))
536 {
537 fullArea += sf[facei];
538 }
539 Eave[coarseI] += (e[facei]*sf[facei])/fineArea;
540 Tave[coarseI] += (pow4(Tf[facei])*sf[facei])/fineArea;
541 }
542 localCoarsePartialArea[compactI++] = fullArea/fineArea;
543 }
544
545 SubList<scalar>
546 (
547 localCoarseEave,
548 Eave.size(),
549 startI
550 ) = Eave;
551
552 SubList<scalar>
553 (
554 localTave,
555 Tave.size(),
556 startI
557 ) = Tave;
558
559 const vectorList coarseNSf = cpp.faceNormals();
560 SubList<vector>
561 (
562 localCoarseNorm,
563 cpp.size(),
564 startI
565 ) = coarseNSf;
566
567 startI += cpp.size();
568 }
569
570
571 SubList<scalar>(compactCoarsePartialArea, nLocalVFCoarseFaces) =
572 localCoarsePartialArea;
573
574 SubList<scalar>(compactCoarseEave, nLocalVFCoarseFaces) =
575 localCoarseEave;
576
577 SubList<scalar>(compactCoarseTave, nLocalVFCoarseFaces) =
578 localTave;
579
580 SubList<vector>(compactCoarseNorm, nLocalVFCoarseFaces) =
581 localCoarseNorm;
582
583 map_->distribute(compactCoarsePartialArea);
584 map_->distribute(compactCoarseEave);
585 map_->distribute(compactCoarseTave);
586 map_->distribute(compactCoarseNorm);
587
588
589 // Calculate coarse hitFaces and Sun direct hit heat fluxes
590 scalarList localqDiffusive(nLocalVFCoarseFaces, Zero);
591
592 label locaFaceI = 0;
593 forAll(includePatches_, i)
594 {
595 const label patchID = includePatches_[i];
596 const polyPatch& pp = coarseMesh_->boundaryMesh()[patchID];
597 const polyPatch& ppf = mesh_.boundaryMesh()[patchID];
598
599 const labelList& coarsePatchFace = coarseMesh_->patchFaceMap()[patchID];
600 const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
601
602
603 scalarField a(ppf.size(), Zero);
604
605 for (label bandi = 0; bandi < nBands_; bandi++)
606 {
607 const tmp<scalarField> ta =
608 spectralDistribution_[bandi]
609 * absorptivity_[patchID][bandi]();
610 a += ta();
611 }
612
613 forAll(pp, coarseI)
614 {
615 const label coarseFaceID = coarsePatchFace[coarseI];
616 const labelList& fineFaces = coarseToFine_[i][coarseFaceID];
617 UIndirectList<scalar> fineSf(sf, fineFaces);
618 scalar fineArea = sum(fineSf());
619
620 // // Weighting absorptivity per area on secondary diffussive flux
621 scalar aAve = 0.0;
622 forAll(fineFaces, j)
623 {
624 label facei = fineFaces[j];
625 aAve += (a[facei]*sf[facei])/fineArea;
626 }
627
628 const scalarList& vf = FmyProc[locaFaceI];
629 const labelList& compactFaces = visibleFaceFaces_[locaFaceI];
630
631 forAll(compactFaces, j)
632 {
633 label compactI = compactFaces[j];
634
635 scalar qin =
636 (
637 solarCalc_.directSolarRad()*solarCalc_.direction()
638 & compactCoarseNorm[compactI]
639 )*compactCoarsePartialArea[compactI];
640
641 // q emission
642 scalar qem =
643 compactCoarseEave[compactI]
644 *physicoChemical::sigma.value()
645 *compactCoarseTave[compactI];
646
647 // compactCoarseEave is the diffussive reflected coeff
648 scalar qDiff = (compactCoarseEave[compactI])*qin;
649
650 localqDiffusive[locaFaceI] += (qDiff)*aAve*vf[j];
651 }
652 locaFaceI++;
653 }
654 }
655
656 volScalarField::Boundary& qsBf = qsecondRad_.boundaryFieldRef();
657
658 // Fill qsecondRad_
659 label compactId = 0;
660 forAll(includePatches_, i)
661 {
662 const label patchID = includePatches_[i];
663 const polyPatch& pp = coarseMesh_->boundaryMesh()[patchID];
664
665 if (pp.size() > 0)
666 {
667 scalarField& qrp = qsBf[patchID];
668
669 const labelList& coarsePatchFace =
670 coarseMesh_->patchFaceMap()[patchID];
671
672 forAll(pp, coarseI)
673 {
674 const label coarseFaceID = coarsePatchFace[coarseI];
675 const labelList& fineFaces = coarseToFine_[i][coarseFaceID];
676 forAll(fineFaces, k)
677 {
678 label facei = fineFaces[k];
679 qrp[facei] = localqDiffusive[compactId];
680 }
681 compactId ++;
682 }
683 }
684 }
685
686 const scalarField& V = mesh_.V();
687 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
688 volScalarField::Boundary& qrBf = qr_.boundaryFieldRef();
689
690 for (const label patchID : includePatches)
691 {
692 const scalarField& qSecond = qsecondRad_.boundaryField()[patchID];
693 if (includeMappedPatchBasePatches[patchID])
694 {
695 qrBf[patchID] += qSecond;
696 }
697 else
698 {
699 const polyPatch& pp = patches[patchID];
700 const labelUList& cellids = pp.faceCells();
701 const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
702 forAll(pp, facei)
703 {
704 const label celli = cellids[facei];
705 Ru_[celli] += qSecond[facei]*sf[facei]/V[celli];
706 }
708 }
709}
710*/
711
712// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
713
714Foam::radiation::solarLoad::solarLoad(const volScalarField& T)
715:
716 radiationModel(typeName, T),
718 solarCalc_(coeffs_, mesh_),
719 dict_(coeffs_),
720 qr_
721 (
722 IOobject
723 (
724 "qr",
725 mesh_.time().timeName(),
726 mesh_,
727 IOobject::READ_IF_PRESENT,
728 IOobject::AUTO_WRITE
729 ),
730 mesh_,
732 ),
733 hitFaces_(),
734 reflectedFaces_(),
735 Ru_
736 (
737 IOobject
738 (
739 "Ru",
740 mesh_.time().timeName(),
741 mesh_,
742 IOobject::NO_READ,
743 IOobject::NO_WRITE
744 ),
745 mesh_,
747 ),
748 absorptivity_(mesh_.boundaryMesh().size()),
749 spectralDistribution_(),
750 spectralDistributions_(),
751 qprimaryRad_(0),
752 nBands_(0),
753 updateTimeIndex_(0),
754 solidCoupled_(true),
755 wallCoupled_(false),
756 updateAbsorptivity_(false),
757 useReflectedRays_(false),
758 firstIter_(true)
759{
760 initialise(coeffs_);
761}
762
763
764Foam::radiation::solarLoad::solarLoad
765(
766 const dictionary& dict,
767 const volScalarField& T
768)
769:
771 solarLoadBase(mesh_),
772 solarCalc_(dict, mesh_),
773 dict_(dict),
774 qr_
775 (
777 (
778 "qr",
779 mesh_.time().timeName(),
780 mesh_,
781 IOobject::READ_IF_PRESENT,
782 IOobject::AUTO_WRITE
783 ),
784 mesh_,
786 ),
787 hitFaces_(),
788 reflectedFaces_(),
789 Ru_
790 (
792 (
793 "Ru",
794 mesh_.time().timeName(),
795 mesh_,
796 IOobject::NO_READ,
797 IOobject::NO_WRITE
798 ),
799 mesh_,
801 ),
802 absorptivity_(mesh_.boundaryMesh().size()),
803 spectralDistribution_(),
804 spectralDistributions_(),
805 qprimaryRad_(0),
806 nBands_(0),
807 updateTimeIndex_(0),
808 solidCoupled_(true),
809 wallCoupled_(false),
810 updateAbsorptivity_(false),
811 useReflectedRays_(false),
812 firstIter_(true)
814 initialise(dict);
815}
816
817
818// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
819
821{
823 {
824 return true;
825 }
826
827 return false;
828}
829
830
832{
833 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
834 labelHashSet includePatches;
835 forAll(patches, patchI)
836 {
837 const polyPatch& pp = patches[patchI];
838 if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
839 {
840 includePatches.insert(patchI);
841 }
842 }
843
844 labelHashSet includeMappedPatchBasePatches;
845 forAll(patches, patchI)
846 {
847 const polyPatch& pp = patches[patchI];
848 if
849 (
850 (isA<mappedPatchBase>(pp) && solidCoupled_)
851 || (isA<wallPolyPatch>(pp) && wallCoupled_)
852 )
853 {
854 includeMappedPatchBasePatches.insert(patchI);
855 }
856 }
857
858 if (updateAbsorptivity_ || firstIter_)
859 {
860 updateAbsorptivity(includePatches);
861 }
862
863 const bool facesChanged = updateHitFaces();
864
865 const bool timeDependentLoad =
866 solarCalc_.sunLoadModel() == solarCalculator::mSunLoadTimeDependent;
867
868 if (firstIter_ || facesChanged || timeDependentLoad)
869 {
870 // Reset Ru
871 Ru_ = Zero;
872
873 solarCalc_.correctDirectSolarRad();
874 solarCalc_.correctDiffuseSolarRad();
875
876 spectralDistribution_ =
877 spectralDistributions_->value(mesh_.time().value());
878
879 spectralDistribution_ =
880 spectralDistribution_/sum(spectralDistribution_);
881
882 // Add direct hit radiation
883 const labelList& hitFacesId = hitFaces_->rayStartFaces();
884 updateDirectHitRadiation(hitFacesId, includeMappedPatchBasePatches);
885
886 // Add sky diffusive radiation
887 updateSkyDiffusiveRadiation
888 (
889 includePatches,
890 includeMappedPatchBasePatches
891 );
892
893 // Add specular reflected radiation
894 if (useReflectedRays_)
895 {
896 updateReflectedRays(includeMappedPatchBasePatches);
897 }
898
899 firstIter_ = false;
900 }
901
902 if (debug)
903 {
904 if (mesh_.time().writeTime())
906 Ru_.write();
907 }
908 }
909}
910
911
913{
915 (
916 "Rp",
918 mesh_,
920 (
923 )
924 );
925}
926
927
930{
931 return Ru_;
932}
933
934
937{
938 return solarCalc_;
939}
940
941
944{
945 return hitFaces_();
946}
947
948
949// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
static FOAM_NO_DANGLING_REFERENCE const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
Calculates the reflecting faces from specular surfaces. It only takes into account the first reflecti...
Helper class to calculate visible faces for global, sun-like illumination.
Definition faceShading.H:58
static labelList nonCoupledPatches(const polyMesh &mesh)
Helper: return all uncoupled patches.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
const labelList faceZoneIds() const
Return identifiers of face zones activated for boundary radiation.
tmp< scalarField > absorptivity(const label patchI, const label bandI=0, const vectorField *incomingDirection=nullptr, const scalarField *T=nullptr) const
Access boundary absorptivity on patch.
Top level model for radiation modelling.
const fvMesh & mesh_
Reference to the mesh database.
virtual bool read()=0
Read radiationProperties dictionary.
dictionary coeffs_
Radiation model dictionary.
const volScalarField & T() const noexcept
Return access to the temperature field.
Base class for solarLoad models.
solarLoadBase(const fvMesh &mesh)
Construct.
The solarLoad radiation model includes Sun primary hits, their reflective fluxes and diffusive sky ra...
Definition solarLoad.H:163
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant).
Definition solarLoad.C:922
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4).
Definition solarLoad.C:905
virtual const solarCalculator & solarCalculatorRef() const
Return const reference to the solar calculator.
Definition solarLoad.C:929
virtual const faceShading & faceShadingRef() const
Return const reference to the face shading calculator.
Definition solarLoad.C:936
bool read()
Read radiationProperties dictionary.
Definition solarLoad.C:813
void calculate()
Solve radiation equations.
Definition solarLoad.C:824
A solar calculator model providing models for the solar direction and solar loads.
A class for managing temporary objects.
Definition tmp.H:75
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition tmpI.H:221
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const polyBoundaryMesh & patches
word timeName
Definition getTimeIndex.H:3
const expr V(m.psi().mesh().V())
Different types of constants.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for radiation modelling.
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensionedScalar pow4(const dimensionedScalar &ds)
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
#define addToRadiationRunTimeSelectionTables(model)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.