Loading...
Searching...
No Matches
viewFactor.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-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2025 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 "viewFactor.H"
30#include "surfaceFields.H"
31#include "constants.H"
33#include "typeInfo.H"
37
38
39
40using namespace Foam::constant;
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
44namespace Foam
45{
46 namespace radiation
47 {
50 }
52
54 = "viewFactorWall";
55
56// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
57
59{
61
63 {
65 for (label patchi=0; patchi < patches.size(); patchi++)
66 {
67 finalAgglom_[patchi] = identity(patches[patchi].size());
68 }
69 }
70
71 coarseMesh_.reset
72 (
74 (
76 (
77 "coarse:" + mesh_.name(),
78 mesh_.polyMesh::instance(),
79 mesh_.time(),
83 ),
84 mesh_,
86 )
87 );
88
89 const polyBoundaryMesh& coarsePatches = coarseMesh_->boundaryMesh();
90
92
93 for (const label patchi : selectedPatches_)
94 {
95 nLocalCoarseFaces_ += coarsePatches[patchi].size();
96 }
97
98 if (debug)
99 {
100 Pout<< "radiation::viewFactor::initialise() Selected patches:"
102 Pout<< "radiation::viewFactor::initialise() Number of coarse faces:"
104 }
105
107
109 << "Total number of clusters : " << totalNCoarseFaces_ << endl;
110
111 useDirect_ = coeffs_.getOrDefault<bool>("useDirectSolver", true);
112
113
114
115 map_.reset
116 (
117 new IOmapDistribute
118 (
119 IOobject
120 (
121 "mapDist",
122 mesh_.facesInstance(),
123 mesh_,
127 )
128 )
129 );
130
131 FmyProc_.reset
132 (
134 (
135 IOobject
136 (
137 "F",
138 mesh_.facesInstance(),
139 mesh_,
143 )
144 )
145 );
146
147 globalFaceFaces_.reset
148 (
150 (
151 IOobject
152 (
153 "globalFaceFaces",
154 mesh_.facesInstance(),
155 mesh_,
159 )
160 )
161 );
162
163
164 globalIndex globalNumbering(nLocalCoarseFaces_);
165
166 // Size coarse qr
167 qrBandI_.setSize(nBands_);
168 forAll (qrBandI_, bandI)
169 {
170 qrBandI_[bandI].setSize(nLocalCoarseFaces_, 0.0);
171 }
172
173 if (!useDirect_)
174 {
175 DynamicList<label> dfaceJ;
176
177 // Per processor to owner (local)/neighbour (remote)
178 List<DynamicList<label>> procOwner(Pstream::nProcs());
179 List<DynamicList<label>> dynProcNeighbour(Pstream::nProcs());
180
181 forAll (globalFaceFaces_(), iFace)
182 {
183 const labelList& visFaces = globalFaceFaces_()[iFace];
184 forAll (visFaces, j)
185 {
186 label gFacej = visFaces[j];
187 label proci = globalNumbering.whichProcID(gFacej);
188 label faceJ = globalNumbering.toLocal(proci, gFacej);
189
190 if (Pstream::myProcNo() == proci)
191 {
192 edge e(iFace, faceJ);
193 if (rays_.insert(e))
194 {
195 dfaceJ.append(j);
196 }
197 }
198 else
199 {
200 label gFaceI =
201 globalNumbering.toGlobal(Pstream::myProcNo(), iFace);
202
203 label proci = globalNumbering.whichProcID(gFacej);
204
205 label facei =
206 globalNumbering.toLocal(Pstream::myProcNo(), gFaceI);
207
208 label remoteFacei = globalNumbering.toLocal(proci, gFacej);
209
210 procOwner[proci].append(facei);
211 dynProcNeighbour[proci].append(remoteFacei);
212 }
213 }
214 }
215
216 mapRayToFmy_.transfer(dfaceJ);
217
218 labelList upper(rays_.size(), -1);
219 labelList lower(rays_.size(), -1);
220
221 const edgeList raysLst(rays_.sortedToc());
222 label rayI = 0;
223 for (const auto& e : raysLst)
224 {
225 label faceI = e.start();
226 label faceJ = e.end();
227 upper[rayI] = max(faceI, faceJ);
228 lower[rayI] = min(faceI, faceJ);
229 rayI++;
230 }
231
232 labelListList procNeighbour(dynProcNeighbour.size());
233 forAll(procNeighbour, i)
234 {
235 procNeighbour[i] = std::move(dynProcNeighbour[i]);
236 }
237 labelListList mySendCells;
238 Pstream::exchange<labelList, label>(procNeighbour, mySendCells);
239
240 label nbri = 0;
241 forAll(procOwner, proci)
242 {
243 if (procOwner[proci].size())
244 {
245 nbri++;
246 }
247 if (mySendCells[proci].size())
248 {
249 nbri++;
250 }
251 }
252
253 DebugInFunction<< "Number of procBound : " << nbri << endl;
254
255 PtrList<const lduPrimitiveProcessorInterface> primitiveInterfaces(nbri);
256 internalCoeffs_.clear();
257 boundaryCoeffs_.resize_null(nbri);
258
259 nbri = 0;
260
261 procToInterface_.setSize(Pstream::nProcs(), -1);
262
263 forAll(procOwner, proci)
264 {
265 if (proci < Pstream::myProcNo() && procOwner[proci].size())
266 {
267 if (debug)
268 {
269 Pout<< "Adding interface " << nbri
270 << " to receive my " << procOwner[proci].size()
271 << " from " << proci << endl;
272 }
273 procToInterface_[proci] = nbri;
274 primitiveInterfaces.set
275 (
276 nbri++,
277 new lduPrimitiveProcessorInterface
278 (
279 procOwner[proci],
281 proci,
282 tensorField(0),
284 )
285 );
286 }
287 else if (proci > Pstream::myProcNo() && mySendCells[proci].size())
288 {
289 if (debug)
290 {
291 Pout<< "Adding interface " << nbri
292 << " to send my " << mySendCells[proci].size()
293 << " to " << proci << endl;
294 }
295 primitiveInterfaces.set
296 (
297 nbri++,
298 new lduPrimitiveProcessorInterface
299 (
300 mySendCells[proci],
302 proci,
303 tensorField(0),
305 )
306 );
307 }
308 }
309 forAll(procOwner, proci)
310 {
311 if (proci > Pstream::myProcNo() && procOwner[proci].size())
312 {
313 if (debug)
314 {
315 Pout<< "Adding interface " << nbri
316 << " to receive my " << procOwner[proci].size()
317 << " from " << proci << endl;
318 }
319 procToInterface_[proci] = nbri;
320 primitiveInterfaces.set
321 (
322 nbri++,
323 new lduPrimitiveProcessorInterface
324 (
325 procOwner[proci],
327 proci,
328 tensorField(0),
330 )
331 );
332 }
333 else if (proci < Pstream::myProcNo() && mySendCells[proci].size())
334 {
335 if (debug)
336 {
337 Pout<< "Adding interface " << nbri
338 << " to send my " << mySendCells[proci].size()
339 << " to " << proci << endl;
340 }
341 primitiveInterfaces.set
342 (
343 nbri++,
344 new lduPrimitiveProcessorInterface
345 (
346 mySendCells[proci],
348 proci,
349 tensorField(0),
351 )
352 );
353 }
354 }
355
356 forAll (boundaryCoeffs_, proci)
357 {
359 (
360 proci,
361 new Field<scalar>
362 (
363 primitiveInterfaces[proci].faceCells().size(),
364 Zero
365 )
366 );
367 }
368
369 lduInterfacePtrsList allInterfaces(primitiveInterfaces.size());
370
371 forAll(primitiveInterfaces, i)
372 {
373 const lduPrimitiveProcessorInterface& pp = primitiveInterfaces[i];
374
375 allInterfaces.set(i, &pp);
376 }
377
378 const lduSchedule ps
379 (
381 <lduPrimitiveProcessorInterface>(allInterfaces)
382 );
383
384 PtrList<const lduInterface> allInterfacesPtr(allInterfaces.size());
385 forAll (allInterfacesPtr, i)
386 {
387 const lduPrimitiveProcessorInterface& pp = primitiveInterfaces[i];
388
389 allInterfacesPtr.set
390 (
391 i,
392 new lduPrimitiveProcessorInterface(pp)
393 );
394 }
395
396 lduPtr_.reset
397 (
398 new lduPrimitiveMesh
399 (
400 globalFaceFaces_().size(), //rays_.size(),
401 lower,
402 upper,
403 allInterfacesPtr,
404 ps,
406 )
407 );
408
409 // Set size for local lduMatrix
410 matrixPtr_.reset(new lduMatrix(lduPtr_()));
411
412 scalarListList& myF = FmyProc_();
413
414 if (coeffs_.get<bool>("smoothing"))
415 {
416 scalar maxDelta = 0;
417 scalar totalDelta = 0;
418
419 if (myF.size())
420 {
421 forAll (myF, i)
422 {
423 scalar sumF = 0.0;
424 scalarList& myFij = myF[i];
425 forAll (myFij, j)
426 {
427 sumF += myFij[j];
428 }
429 const scalar delta = sumF - 1.0;
430 forAll (myFij, j)
431 {
432 myFij[j] *= (1.0 - delta/(sumF + 0.001));
433 }
434 totalDelta += delta;
435 if (delta > maxDelta)
436 {
437 maxDelta = delta;
438 }
439 }
440 totalDelta /= myF.size();
441 }
442 reduce(totalDelta, sumOp<scalar>());
443 reduce(maxDelta, maxOp<scalar>());
444 Info << "Smoothing average delta : " << totalDelta << endl;
445 Info << "Smoothing maximum delta : " << maxDelta << nl << endl;
446 }
447 }
448
449 if (useDirect_)
450 {
451 List<labelListList> globalFaceFacesProc(Pstream::nProcs());
452 globalFaceFacesProc[Pstream::myProcNo()] = globalFaceFaces_();
453 Pstream::gatherList(globalFaceFacesProc);
454
455 List<scalarListList> F(Pstream::nProcs());
458
459 if (Pstream::master())
460 {
461 Fmatrix_.reset
462 (
464 );
465
467 << "Insert elements in the matrix..." << endl;
468
469 for (const int procI : Pstream::allProcs())
470 {
472 (
473 globalNumbering,
474 procI,
475 globalFaceFacesProc[procI],
476 F[procI],
477 Fmatrix_()
478 );
479 }
480
481 if (coeffs_.get<bool>("smoothing"))
482 {
483 DebugInFunction << "Smoothing the matrix..." << endl;
484
485 for (label i=0; i<totalNCoarseFaces_; i++)
486 {
487 scalar sumF = 0.0;
488 for (label j=0; j<totalNCoarseFaces_; j++)
489 {
490 sumF += Fmatrix_()(i, j);
491 }
492
493 const scalar delta = sumF - 1.0;
494 for (label j=0; j<totalNCoarseFaces_; j++)
495 {
496 Fmatrix_()(i, j) *= (1.0 - delta/(sumF + 0.001));
497 }
498 }
499 }
500
501 coeffs_.readEntry("constantEmissivity", constEmissivity_);
503 {
504 CLU_.reset
505 (
507 );
508
509 pivotIndices_.setSize(CLU_().m());
510 }
511 }
512 }
513
514 coeffs_.readIfPresent("useSolarLoad", useSolarLoad_);
515
516 if (useSolarLoad_)
517 {
518 const dictionary& solarDict = this->subDict("solarLoadCoeffs");
519 solarLoad_.reset(new solarLoad(solarDict, T_));
520
521 if (solarLoad_->nBands() != nBands_)
522 {
524 << "Solar radiation and view factor band numbers "
525 << "are different"
526 << abort(FatalError);
527 }
528
529 Info<< "Creating Solar Load Model " << nl;
530 }
531}
532
533
534// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
535
536Foam::radiation::viewFactor::viewFactor(const volScalarField& T)
537:
539 finalAgglom_
540 (
542 (
543 "finalAgglom",
544 mesh_.facesInstance(),
545 mesh_,
546 IOobject::READ_IF_PRESENT,
547 IOobject::NO_WRITE,
548 IOobject::NO_REGISTER
549 )
550 ),
551 map_(),
552 coarseMesh_(),
553// (
554// IOobject
555// (
556// "coarse:" + mesh_.name(),
557// mesh_.polyMesh::instance(),
558// mesh_.time(),
559// IOobject::NO_READ,
560// IOobject::NO_WRITE,
561// IOobject::NO_REGISTER
562// ),
563// mesh_,
564// finalAgglom_
565// ),
566 qr_
567 (
569 (
570 "qr",
571 mesh_.time().timeName(),
572 mesh_,
573 IOobject::MUST_READ,
574 IOobject::AUTO_WRITE
575 ),
576 mesh_
577 ),
578 Fmatrix_(),
579 CLU_(),
580 selectedPatches_(mesh_.boundary().size(), -1),
581 totalNCoarseFaces_(0),
582 nLocalCoarseFaces_(0),
583 constEmissivity_(false),
584 iterCounter_(0),
585 pivotIndices_(0),
586 useSolarLoad_(false),
587 solarLoad_(),
588 nBands_(coeffs_.getOrDefault<label>("nBands", 1)),
589 FmyProc_()
590{
591 initialise();
592}
593
594
595Foam::radiation::viewFactor::viewFactor
596(
597 const dictionary& dict,
598 const volScalarField& T
599)
600:
602 finalAgglom_
603 (
605 (
606 "finalAgglom",
607 mesh_.facesInstance(),
608 mesh_,
609 IOobject::READ_IF_PRESENT,
610 IOobject::NO_WRITE,
611 IOobject::NO_REGISTER
612 )
613 ),
614 map_(),
615 coarseMesh_(),
616// (
617// IOobject
618// (
619// "coarse:" + mesh_.name(),
620// mesh_.polyMesh::instance(),
621// mesh_.time(),
622// IOobject::NO_READ,
623// IOobject::NO_WRITE,
624// IOobject::NO_REGISTER
625// ),
626// mesh_,
627// finalAgglom_
628// ),
629 qr_
630 (
632 (
633 "qr",
634 mesh_.time().timeName(),
635 mesh_,
636 IOobject::MUST_READ,
637 IOobject::AUTO_WRITE
638 ),
639 mesh_
640 ),
641 Fmatrix_(),
642 CLU_(),
643 selectedPatches_(mesh_.boundary().size(), -1),
644 totalNCoarseFaces_(0),
645 nLocalCoarseFaces_(0),
646 constEmissivity_(false),
647 iterCounter_(0),
648 pivotIndices_(0),
649 useSolarLoad_(false),
650 solarLoad_(),
651 nBands_(coeffs_.getOrDefault<label>("nBands", 1)),
652 FmyProc_()
654 initialise();
655}
656
657
658// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
659
661{
663 {
664 return true;
665 }
666
667 return false;
668}
669
670
672(
673 const globalIndex& globalNumbering,
674 const label procI,
675 const labelListList& globalFaceFaces,
676 const scalarListList& viewFactors,
677 scalarSquareMatrix& Fmatrix
678)
679{
680 forAll(viewFactors, faceI)
681 {
682 const scalarList& vf = viewFactors[faceI];
683 const labelList& globalFaces = globalFaceFaces[faceI];
684
685 label globalI = globalNumbering.toGlobal(procI, faceI);
686 forAll(globalFaces, i)
688 Fmatrix[globalI][globalFaces[i]] = vf[i];
689 }
690 }
691}
692
693
695{
696 // Store previous iteration
697 qr_.storePrevIter();
698
699 if (useSolarLoad_)
700 {
701 solarLoad_->calculate();
702 }
703
704 // Global net radiation
705 scalarField qNet(totalNCoarseFaces_, 0.0);
706
707 // Local net radiation
708 scalarField qTotalCoarse(nLocalCoarseFaces_, 0.0);
709
710 // Referen to fvMesh qr field
711 volScalarField::Boundary& qrBf = qr_.boundaryFieldRef();
712
713 globalIndex globalNumbering(nLocalCoarseFaces_);
714
715 const boundaryRadiationProperties& boundaryRadiation =
717
718 for (label bandI = 0; bandI < nBands_; bandI++)
719 {
720 // Global bandI radiation
721 scalarField qBandI(totalNCoarseFaces_, 0.0);
722
723 scalarField compactCoarseT4(map_->constructSize(), 0.0);
724 scalarField compactCoarseE(map_->constructSize(), 0.0);
725 scalarField compactCoarseHo(map_->constructSize(), 0.0);
726
727 // Fill local averaged(T), emissivity(E) and external heatFlux(Ho)
728 DynamicList<scalar> localCoarseT4ave(nLocalCoarseFaces_);
729 DynamicList<scalar> localCoarseEave(nLocalCoarseFaces_);
730 DynamicList<scalar> localCoarseHoave(nLocalCoarseFaces_);
731
732 forAll(selectedPatches_, i)
733 {
734 label patchID = selectedPatches_[i];
735
736 const scalarField& Tp = T_.boundaryField()[patchID];
737 const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
738
739 fvPatchScalarField& qrPatch = qrBf[patchID];
740
741 auto& qrp =
742 refCast
743 <
745 >(qrPatch);
746
747 const tmp<scalarField> teb =
748 boundaryRadiation.emissivity(patchID, bandI, nullptr, &Tp);
749
750 const scalarField& eb = teb();
751
752 const tmp<scalarField> tHoi = qrp.qro(bandI);
753 const scalarField& Hoi = tHoi();
754
755 const polyPatch& pp = coarseMesh_->boundaryMesh()[patchID];
756
757 const labelList& coarsePatchFace =
758 coarseMesh_->patchFaceMap()[patchID];
759
760 scalarList T4ave(pp.size(), 0.0);
761 scalarList Eave(pp.size(), 0.0);
762 scalarList Hoiave(pp.size(), 0.0);
763
764 if (pp.size() > 0)
765 {
766 const labelList& agglom = finalAgglom_[patchID];
767 label nAgglom = max(agglom) + 1;
768
769 labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
770
771 forAll(coarseToFine, coarseI)
772 {
773 const label coarseFaceID = coarsePatchFace[coarseI];
774 const labelList& fineFaces = coarseToFine[coarseFaceID];
775 UIndirectList<scalar> fineSf
776 (
777 sf,
778 fineFaces
779 );
780
781 const scalar area = sum(fineSf());
782
783 // Temperature, emissivity and external flux area weighting
784 forAll(fineFaces, j)
785 {
786 label facei = fineFaces[j];
787 T4ave[coarseI] += (pow4(Tp[facei])*sf[facei])/area;
788 Eave[coarseI] += (eb[facei]*sf[facei])/area;
789 Hoiave[coarseI] += (Hoi[facei]*sf[facei])/area;
790 }
791 }
792 }
793
794 localCoarseT4ave.append(T4ave);
795 localCoarseEave.append(Eave);
796 localCoarseHoave.append(Hoiave);
797
798 }
799
800 // Fill the local values to distribute
801 SubList<scalar>(compactCoarseT4, nLocalCoarseFaces_) =
802 localCoarseT4ave;
803 SubList<scalar>(compactCoarseE, nLocalCoarseFaces_) = localCoarseEave;
804 SubList<scalar>(compactCoarseHo, nLocalCoarseFaces_) =
805 localCoarseHoave;
806
807
808 // Distribute data
809 map_->distribute(compactCoarseT4);
810 map_->distribute(compactCoarseE);
811 map_->distribute(compactCoarseHo);
812
813 // Distribute local global ID
814 labelList compactGlobalIds(map_->constructSize(), Zero);
815
816 SubList<label>
817 (
818 compactGlobalIds,
819 nLocalCoarseFaces_
820 ) = identity
821 (
822 globalNumbering.localSize(),
823 globalNumbering.localStart()
824 );
825
826 map_->distribute(compactGlobalIds);
827
828 if (!useDirect_)
829 {
830 const labelList globalToCompact
831 (
832 invert(totalNCoarseFaces_, compactGlobalIds)
833 );
834
835 scalarField& diag = matrixPtr_->diag(localCoarseEave.size());
836 scalarField& upper = matrixPtr_->upper(rays_.size());
837 scalarField& lower = matrixPtr_->lower(rays_.size());
838
839 scalarField source(nLocalCoarseFaces_, 0);
840
841
842 // Local diag and source
843 forAll(source, i)
844 {
845 const scalar sigmaT4 =
846 physicoChemical::sigma.value()*localCoarseT4ave[i];
847
848 diag[i] = 1/localCoarseEave[i];
849
850 source[i] += -sigmaT4 + localCoarseHoave[i];
851 }
852
853 // Local matrix coefficients
854 if (!constEmissivity_ || iterCounter_ == 0)
855 {
856 const edgeList raysLst(rays_.sortedToc());
857
858 label rayI = 0;
859 for (const auto& e : raysLst)
860 {
861 label facelJ = e.end();
862 label faceI = e.start();
863
864 label faceJ = mapRayToFmy_[rayI];
865
866 lower[rayI] =
867 (1-1/localCoarseEave[faceI])*FmyProc_()[faceI][faceJ];
868
869 upper[rayI] =
870 (1-1/localCoarseEave[facelJ])*FmyProc_()[faceI][faceJ];
871
872 rayI++;
873 }
874 iterCounter_++;
875 }
876
877 // Extra local contribution to the source and extra matrix coefs
878 // from other procs (boundaryCoeffs)
879 label nInterfaces = lduPtr_().interfaces().size();
880 labelList boundCoeffI(nInterfaces, Zero);
881 forAll (globalFaceFaces_(), iFace)
882 {
883 const labelList& visFaces = globalFaceFaces_()[iFace];
884 forAll (visFaces, jFace)
885 {
886 label gFacej = visFaces[jFace];
887 label proci = globalNumbering.whichProcID(gFacej);
888 if (Pstream::myProcNo() == proci)
889 {
890 label lFacej =
891 globalNumbering.toLocal
892 (
894 gFacej
895 );
896
897 const scalar sigmaT4 =
899 *localCoarseT4ave[lFacej];
900
901 source[iFace] += FmyProc_()[iFace][jFace]*sigmaT4;
902 }
903 else
904 {
905 label compactFaceJ = globalToCompact[gFacej];
906 const scalar sigmaT4 =
908 *compactCoarseT4[compactFaceJ];
909
910 source[iFace] += FmyProc_()[iFace][jFace]*sigmaT4;
911
912 label interfaceI = procToInterface_[proci];
913
914 boundaryCoeffs_
915 [interfaceI][boundCoeffI[interfaceI]++] =
916 -(1-1/compactCoarseE[compactFaceJ])
917 *FmyProc_()[iFace][jFace];
918 }
919 }
920 }
921
922 PtrList<const lduInterfaceField> interfaces(nInterfaces);
923 for(label i = 0; i < interfaces.size(); i++)
924 {
925 interfaces.set
926 (
927 i,
928 new lduCalculatedProcessorField<scalar>
929 (
930 lduPtr_().interfaces()[i]
931 )
932 );
933 }
934
935 const dictionary& solverControls =
936 qr_.mesh().solverDict
937 (
938 qr_.select(qr_.mesh().data().isFinalIteration())
939 );
940
941 // Solver call
943 (
944 "qr",
945 matrixPtr_(),
946 boundaryCoeffs_,
947 internalCoeffs_,
948 interfaces,
949 solverControls
950 )->solve(qrBandI_[bandI], source);
951
952 solverPerf.print(Info.masterStream(qr_.mesh().comm()));
953
954 qTotalCoarse += qrBandI_[bandI];
955 }
956
957 if (useDirect_)
958 {
959 // Create global size vectors
960 scalarField T4(totalNCoarseFaces_, 0.0);
961 scalarField E(totalNCoarseFaces_, 0.0);
962 scalarField qrExt(totalNCoarseFaces_, 0.0);
963
964 // Fill lists from compact to global indexes.
965 forAll(compactCoarseT4, i)
966 {
967 T4[compactGlobalIds[i]] = compactCoarseT4[i];
968 E[compactGlobalIds[i]] = compactCoarseE[i];
969 qrExt[compactGlobalIds[i]] = compactCoarseHo[i];
970 }
971
972 Pstream::listReduce(T4, maxOp<scalar>());
973 Pstream::listReduce(E, maxOp<scalar>());
974 Pstream::listReduce(qrExt, maxOp<scalar>());
975
976 if (Pstream::master())
977 {
978 // Variable emissivity
979 if (!constEmissivity_)
980 {
981 scalarSquareMatrix C(totalNCoarseFaces_, 0.0);
982
983 for (label i=0; i<totalNCoarseFaces_; i++)
984 {
985 for (label j=0; j<totalNCoarseFaces_; j++)
986 {
987 const scalar invEj = 1.0/E[j];
988 const scalar sigmaT4 =
990
991 if (i==j)
992 {
993 C(i, j) = invEj;
994 qBandI[i] += -sigmaT4 + qrExt[j];
995 }
996 else
997 {
998 C(i, j) = (1.0 - invEj)*Fmatrix_()(i, j);
999 qBandI[i] += Fmatrix_()(i, j)*sigmaT4;
1000 }
1001 }
1002 }
1003
1004 Info<< "Solving view factor equations for band :"
1005 << bandI << endl;
1006
1007 // Negative coming into the fluid
1008 LUsolve(C, qBandI);
1009 }
1010 else //Constant emissivity
1011 {
1012 // Initial iter calculates CLU and caches it
1013 if (iterCounter_ == 0)
1014 {
1015 for (label i=0; i<totalNCoarseFaces_; i++)
1016 {
1017 for (label j=0; j<totalNCoarseFaces_; j++)
1018 {
1019 const scalar invEj = 1.0/E[j];
1020 if (i==j)
1021 {
1022 CLU_()(i, j) = invEj;
1023 }
1024 else
1025 {
1026 CLU_()(i, j) =
1027 (1.0-invEj)*Fmatrix_()(i, j);
1028 }
1029 }
1030 }
1031
1032 DebugInFunction << "\nDecomposing C matrix..." << endl;
1033
1034 LUDecompose(CLU_(), pivotIndices_);
1035 }
1036
1037 for (label i=0; i<totalNCoarseFaces_; i++)
1038 {
1039 for (label j=0; j<totalNCoarseFaces_; j++)
1040 {
1041 const scalar sigmaT4 =
1043
1044 if (i==j)
1045 {
1046 qBandI[i] += -sigmaT4 + qrExt[j];
1047 }
1048 else
1049 {
1050 qBandI[i] += Fmatrix_()(i, j)*sigmaT4;
1051 }
1052 }
1053 }
1054
1055 Info<< "Solving view factor equations for band : "
1056 << bandI << endl;
1057
1058 LUBacksubstitute(CLU_(), pivotIndices_, qBandI);
1059 iterCounter_ ++;
1060 }
1061 }
1062 // Broadcast qBandI and fill qr
1063 Pstream::broadcast(qBandI);
1064
1065 qNet += qBandI;
1066 }
1067 }
1068
1069 label globCoarseId = 0;
1070 for (const label patchID : selectedPatches_)
1071 {
1072 const polyPatch& pp = coarseMesh_->boundaryMesh()[patchID];
1073
1074 if (pp.size() > 0)
1075 {
1076 scalarField& qrp = qrBf[patchID];
1077 //const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
1078 const labelList& agglom = finalAgglom_[patchID];
1079 label nAgglom = max(agglom)+1;
1080
1081 labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
1082
1083 const labelList& coarsePatchFace =
1084 coarseMesh_->patchFaceMap()[patchID];
1085
1086 //scalar heatFlux = 0.0;
1087 forAll(coarseToFine, coarseI)
1088 {
1089 label globalCoarse = globalNumbering.toGlobal
1090 (
1092 globCoarseId
1093 );
1094
1095 const label coarseFaceID = coarsePatchFace[coarseI];
1096 const labelList& fineFaces = coarseToFine[coarseFaceID];
1097
1098 for (const label facei : fineFaces)
1099 {
1100 if (useDirect_)
1101 {
1102 qrp[facei] = qNet[globalCoarse];
1103 }
1104 else
1105 {
1106 qrp[facei] = qTotalCoarse[globCoarseId];
1107 }
1108 //heatFlux += qrp[facei]*sf[facei];
1109 }
1110 globCoarseId++;
1111 }
1112 }
1113 }
1114
1115 if (debug)
1116 {
1117 for (const label patchID : selectedPatches_)
1118 {
1119 const scalarField& qrp = qrBf[patchID];
1120 const scalarField& magSf = mesh_.magSf().boundaryField()[patchID];
1121 const scalar heatFlux = gWeightedSum(magSf, qrp);
1122
1123 Info<< "Total heat transfer rate at patch: "
1124 << patchID << " "
1125 << heatFlux << endl;
1126 }
1128
1129 // Relax qr if necessary
1130 qr_.relax();
1131}
1132
1133
1135{
1136 return volScalarField::New
1137 (
1138 "Rp",
1140 mesh_,
1142 (
1153 (
1154 "Ru",
1156 mesh_,
1158 );
1159}
1160
1161
1162// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
scalar delta
Macros for easy insertion into run-time selection tables.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Graphite solid properties.
Definition C.H:49
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const Field< Type > &iField)
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions, copy of internal field....
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
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())
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
IOmapDistribute is derived from mapDistribute and IOobject to give the mapDistribute automatic IO fun...
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ 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
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info. A void type suppresses trait and t...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
void setSize(label n)
Alias for resize().
Definition List.H:536
static FOAM_NO_DANGLING_REFERENCE const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
static void exchange(const UList< Container > &sendBufs, const labelUList &recvSizes, List< Container > &recvBufs, const int tag=UPstream::msgType(), const int comm=UPstream::worldComm, const bool wait=true)
Helper: exchange contiguous data. Sends sendBufs, receives into recvBufs using predetermined receive ...
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
void print(Ostream &os) const
Print summary of solver performance to the given stream.
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition UPstream.H:1069
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition UPstream.H:1857
@ gatherList
gatherList [manual algorithm]
Definition UPstream.H:194
@ broadcast
broadcast [MPI]
Definition UPstream.H:189
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition UPtrList.H:366
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
dictionary()
Default construct, a top-level empty dictionary.
Definition dictionary.C:68
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
const Type & value() const noexcept
Return const reference to value.
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
const word & name() const
Return reference to name.
Definition fvMesh.H:387
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
label toLocal(const label proci, const label i) const
From global to local on proci.
label whichProcID(const label proci, const label i) const
Which processor does global id come from? Checks proci first (assumed to occur reasonably frequently)...
label toGlobal(const label proci, const label i) const
From local to global on proci.
label localStart(const label proci) const
Start of proci data.
label localSize(const label proci) const
Size of proci data.
A lduProcessorField type bypassing coupledFvPatchField.
static autoPtr< solver > New(const word &solverName, const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Return a new solver of given type.
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition lduMatrix.H:81
Simplest concrete lduMesh that stores the addressing needed by lduMatrix.
static lduSchedule nonBlockingSchedule(const lduInterfacePtrsList &)
Get non-scheduled send/receive schedule.
Concrete implementation of processor interface. Used to temporarily store settings.
OSstream & masterStream(int communicator)
Return OSstream for output operations on the master process only, Snull on other processes.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
labelList indices(const wordRe &matcher, const bool useGroups=true) const
The (sorted) patch indices for all matches, optionally matching patch groups.
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
tmp< scalarField > emissivity(const label patchI, const label bandI=0, const vectorField *incomingDirection=nullptr, const scalarField *T=nullptr) const
Access boundary emissivity on patch.
This boundary condition provides a grey-diffuse condition for radiative heat flux,...
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.
const volScalarField & T_
Reference to the temperature field.
The solarLoad radiation model includes Sun primary hits, their reflective fluxes and diffusive sky ra...
Definition solarLoad.H:163
View factor radiation model. The system solved is: C q = b where: Cij = deltaij/Ej - (1/Ej - 1)Fij q ...
Definition viewFactor.H:76
void initialise()
Initialise.
Definition viewFactor.C:51
autoPtr< scalarSquareMatrix > Fmatrix_
View factor matrix.
Definition viewFactor.H:132
autoPtr< singleCellFvMesh > coarseMesh_
Coarse mesh.
Definition viewFactor.H:117
autoPtr< lduMatrix > matrixPtr_
Matrix formed from view factors.
Definition viewFactor.H:197
virtual tmp< volScalarField::Internal > Ru() const
Source term component (constant).
autoPtr< scalarSquareMatrix > CLU_
Inverse of C matrix.
Definition viewFactor.H:137
label nLocalCoarseFaces_
Total local coarse faces.
Definition viewFactor.H:157
static const word viewFactorWalls
Static name for view factor walls.
Definition viewFactor.H:99
label nBands_
Number of bands.
Definition viewFactor.H:187
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4).
autoPtr< solarLoad > solarLoad_
Solar load radiation model.
Definition viewFactor.H:182
autoPtr< scalarListIOList > FmyProc_
Local view factors.
Definition viewFactor.H:224
volScalarField qr_
Net radiative heat flux [W/m2].
Definition viewFactor.H:122
autoPtr< IOmapDistribute > map_
Map distributed.
Definition viewFactor.H:112
FieldField< Field, scalar > boundaryCoeffs_
Boundary scalar field containing pseudo-matrix coeffs for boundary cells.
Definition viewFactor.H:209
labelList procToInterface_
Map from proc to interafce.
Definition viewFactor.H:229
List< label > mapRayToFmy_
Map local-ray to j-column for F.
Definition viewFactor.H:219
FieldField< Field, scalar > internalCoeffs_
Boundary scalar field containing pseudo-matrix coeffs for internal cells.
Definition viewFactor.H:203
labelList pivotIndices_
Pivot Indices for LU decomposition.
Definition viewFactor.H:172
bool useDirect_
Use direct or iterative solver.
Definition viewFactor.H:234
autoPtr< labelListIOList > globalFaceFaces_
Visible global faces.
Definition viewFactor.H:142
labelListIOList finalAgglom_
Agglomeration List.
Definition viewFactor.H:107
List< scalarField > qrBandI_
Coarse radiative heat flux.
Definition viewFactor.H:127
label totalNCoarseFaces_
Total global coarse faces.
Definition viewFactor.H:152
edgeHashSet rays_
Rays on local proc.
Definition viewFactor.H:214
bool constEmissivity_
Constant emissivity.
Definition viewFactor.H:162
autoPtr< lduPrimitiveMesh > lduPtr_
Primitive addressing for lduMatrix.
Definition viewFactor.H:192
labelList selectedPatches_
Selected patches.
Definition viewFactor.H:147
bool useSolarLoad_
Use Solar Load model.
Definition viewFactor.H:177
void insertMatrixElements(const globalIndex &index, const label fromProci, const labelListList &globalFaceFaces, const scalarListList &viewFactors, scalarSquareMatrix &matrix)
Insert view factors into main matrix.
Definition viewFactor.C:665
label iterCounter_
Iterations Counter.
Definition viewFactor.H:167
bool read()
Read radiation properties dictionary.
Definition viewFactor.C:653
void calculate()
Solve system of equation(s).
Definition viewFactor.C:687
fvMesh as subset of other mesh. Consists of one cell and all original boundary faces....
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const polyBoundaryMesh & patches
faceListList boundary
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
word timeName
Definition getTimeIndex.H:3
#define DebugInFunction
Report an information message using Foam::Info.
volVectorField F(fluid.F())
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Different types of constants.
Namespace for handling debugging switches.
Definition debug.C:45
const wordList area
Standard area field types (scalar, vector, tensor, etc).
Namespace for radiation modelling.
string upper(const std::string &s)
Return string copy transformed with std::toupper on each character.
string lower(const std::string &s)
Return string copy transformed with std::tolower on each character.
Namespace for OpenFOAM.
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
List< edge > edgeList
List of edge.
Definition edgeList.H:32
IOList< scalarList > scalarListIOList
IO for a List of scalarList.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
void LUBacksubstitute(const scalarSquareMatrix &luMmatrix, const labelList &pivotIndices, List< Type > &source)
LU back-substitution with given source, returning the solution in the source.
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< lduScheduleEntry > lduSchedule
A List of lduSchedule entries.
Definition lduSchedule.H:46
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Type gWeightedSum(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted sum (integral) of a field, using the mag() of the weights.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
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)
IOList< labelList > labelListIOList
IO for a List of labelList.
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
dimensionedScalar pow4(const dimensionedScalar &ds)
errorManip< error > abort(error &err)
Definition errorManip.H:139
SquareMatrix< scalar > scalarSquareMatrix
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition ListOps.C:28
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
Solve the matrix using LU decomposition with pivoting returning the LU form of the matrix and the sol...
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition ListOps.C:125
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
SolverPerformance< scalar > solverPerformance
SolverPerformance instantiated for a scalar.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define addToRadiationRunTimeSelectionTables(model)
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.
Basic run-time type information using word as the type's name. Used to enhance the standard RTTI to c...