viewFactorsGen.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-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
27Application
28 viewFactorsGen
29
30Group
31 grpPreProcessingUtilities
32
33Description
34 View factors are calculated based on a face agglomeration array
35 (finalAgglom generated by faceAgglomerate utility).
36
37 Each view factor between the agglomerated faces i and j (Fij) is calculated
38 using a double integral of the sub-areas composing the agglomeration.
39
40 The patches involved in the view factor calculation are taken from the
41 boundary file and should be part on the group viewFactorWall. ie.:
42
43 floor
44 {
45 type wall;
46 inGroups 2(wall viewFactorWall);
47 nFaces 100;
48 startFace 3100;
49 }
50
51\*---------------------------------------------------------------------------*/
52
53#include "argList.H"
54#include "fvMesh.H"
55#include "volFields.H"
56#include "surfaceFields.H"
58#include "meshTools.H"
59
61#include "DynamicField.H"
62#include "unitConversion.H"
63
64#include "scalarMatrices.H"
65#include "labelListIOList.H"
66#include "scalarListIOList.H"
67
68#include "singleCellFvMesh.H"
69
70#include "IOmapDistribute.H"
71
72using namespace Foam;
73
74
75triSurface triangulate
76(
78 const labelHashSet& includePatches,
79 const labelListIOList& finalAgglom,
81 const globalIndex& globalNumbering,
82 const polyBoundaryMesh& coarsePatches
83)
84{
85 const polyMesh& mesh = bMesh.mesh();
86
87 // Storage for surfaceMesh. Size estimate.
88 DynamicList<labelledTri> triangles(mesh.nBoundaryFaces());
89
90 label newPatchI = 0;
91 label localTriFaceI = 0;
92
93 for (const label patchI : includePatches)
94 {
95 const polyPatch& patch = bMesh[patchI];
96 const pointField& points = patch.points();
97
98 label nTriTotal = 0;
99
100 forAll(patch, patchFaceI)
101 {
102 const face& f = patch[patchFaceI];
103
104 faceList triFaces(f.nTriangles(points));
105
106 label nTri = 0;
107
108 f.triangles(points, nTri, triFaces);
109
110 forAll(triFaces, triFaceI)
111 {
112 const face& f = triFaces[triFaceI];
113
114 triangles.append(labelledTri(f[0], f[1], f[2], newPatchI));
115
116 nTriTotal++;
117
118 triSurfaceToAgglom[localTriFaceI++] = globalNumbering.toGlobal
119 (
120 Pstream::myProcNo(),
121 finalAgglom[patchI][patchFaceI]
122 + coarsePatches[patchI].start()
123 );
124 }
125 }
126
127 newPatchI++;
128 }
129
130 //striSurfaceToAgglom.resize(localTriFaceI-1);
131
132 triangles.shrink();
133
134 // Create globally numbered tri surface
135 triSurface rawSurface(triangles, mesh.points());
136
137 // Create locally numbered tri surface
139 (
140 rawSurface.localFaces(),
141 rawSurface.localPoints()
142 );
143
144 // Add patch names to surface
145 surface.patches().setSize(newPatchI);
146
147 newPatchI = 0;
148
149 for (const label patchI : includePatches)
150 {
151 const polyPatch& patch = bMesh[patchI];
152
153 surface.patches()[newPatchI].index() = patchI;
154 surface.patches()[newPatchI].name() = patch.name();
155 surface.patches()[newPatchI].geometricType() = patch.type();
156
157 newPatchI++;
158 }
159
160 return surface;
161}
162
163
164void writeRays
165(
166 const fileName& fName,
167 const pointField& compactCf,
168 const pointField& myFc,
169 const labelListList& visibleFaceFaces
170)
171{
172 OFstream str(fName);
173 label vertI = 0;
174
175 Pout<< "Dumping rays to " << str.name() << endl;
176
177 forAll(myFc, faceI)
178 {
179 const labelList visFaces = visibleFaceFaces[faceI];
180 forAll(visFaces, faceRemote)
181 {
182 label compactI = visFaces[faceRemote];
183 const point& remoteFc = compactCf[compactI];
184
185 meshTools::writeOBJ(str, myFc[faceI]);
186 vertI++;
187 meshTools::writeOBJ(str, remoteFc);
188 vertI++;
189 str << "l " << vertI-1 << ' ' << vertI << nl;
190 }
191 }
192 str.flush();
193}
194
195
196scalar calculateViewFactorFij
197(
198 const vector& i,
199 const vector& j,
200 const vector& dAi,
201 const vector& dAj
202)
203{
204 vector r = i - j;
205 scalar rMag = mag(r);
206
207 if (rMag > SMALL)
208 {
209 scalar dAiMag = mag(dAi);
210 scalar dAjMag = mag(dAj);
211
212 vector ni = dAi/dAiMag;
213 vector nj = dAj/dAjMag;
214 scalar cosThetaJ = mag(nj & r)/rMag;
215 scalar cosThetaI = mag(ni & r)/rMag;
216
217 return
218 (
219 (cosThetaI*cosThetaJ*dAjMag*dAiMag)
220 /(sqr(rMag)*constant::mathematical::pi)
221 );
222 }
223 else
224 {
225 return 0;
226 }
227}
228
229
230void insertMatrixElements
231(
232 const globalIndex& globalNumbering,
233 const label fromProcI,
234 const labelListList& globalFaceFaces,
235 const scalarListList& viewFactors,
236 scalarSquareMatrix& matrix
237)
238{
239 forAll(viewFactors, faceI)
240 {
241 const scalarList& vf = viewFactors[faceI];
242 const labelList& globalFaces = globalFaceFaces[faceI];
243
244 label globalI = globalNumbering.toGlobal(fromProcI, faceI);
245 forAll(globalFaces, i)
246 {
247 matrix[globalI][globalFaces[i]] = vf[i];
248 }
249 }
250}
251
252
253// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254
255int main(int argc, char *argv[])
256{
257 argList::addNote
258 (
259 "Calculate view factors from face agglomeration array."
260 " The finalAgglom generated by faceAgglomerate utility."
261 );
262
263 #include "addRegionOption.H"
264 #include "setRootCase.H"
265 #include "createTime.H"
266 #include "createNamedMesh.H"
267
268 // Read view factor dictionary
269 IOdictionary viewFactorDict
270 (
272 (
273 "viewFactorsDict",
274 runTime.constant(),
275 mesh,
276 IOobject::MUST_READ_IF_MODIFIED,
277 IOobject::NO_WRITE
278 )
279 );
280
281 const word viewFactorWall("viewFactorWall");
282
283 const bool writeViewFactors =
284 viewFactorDict.getOrDefault("writeViewFactorMatrix", false);
285
286 const bool dumpRays =
287 viewFactorDict.getOrDefault("dumpRays", false);
288
289 const label debug = viewFactorDict.getOrDefault<label>("debug", 0);
290
291 // Read agglomeration map
292 labelListIOList finalAgglom
293 (
295 (
296 "finalAgglom",
297 mesh.facesInstance(),
298 mesh,
299 IOobject::MUST_READ,
300 IOobject::NO_WRITE,
301 false
302 )
303 );
304
305 // Create the coarse mesh using agglomeration
306 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
307
308 if (debug)
309 {
310 Pout << "\nCreating single cell mesh..." << endl;
311 }
312
313 singleCellFvMesh coarseMesh
314 (
316 (
317 "coarse:" + mesh.name(),
318 runTime.timeName(),
319 runTime,
320 IOobject::NO_READ,
321 IOobject::NO_WRITE
322 ),
323 mesh,
324 finalAgglom
325 );
326
327 if (debug)
328 {
329 Pout << "\nCreated single cell mesh..." << endl;
330 }
331
332
333 // Calculate total number of fine and coarse faces
334 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
335
336 label nCoarseFaces = 0; //total number of coarse faces
337 label nFineFaces = 0; //total number of fine faces
338
339 const polyBoundaryMesh& patches = mesh.boundaryMesh();
340 const polyBoundaryMesh& coarsePatches = coarseMesh.boundaryMesh();
341
342 labelList viewFactorsPatches(patches.indices(viewFactorWall));
343 for (const label patchi : viewFactorsPatches)
344 {
345 nCoarseFaces += coarsePatches[patchi].size();
346 nFineFaces += patches[patchi].size();
347 }
348
349 // total number of coarse faces
350 label totalNCoarseFaces = nCoarseFaces;
351
352 reduce(totalNCoarseFaces, sumOp<label>());
353
354 if (Pstream::master())
355 {
356 Info << "\nTotal number of coarse faces: "<< totalNCoarseFaces << endl;
357 }
358
359 if (Pstream::master() && debug)
360 {
361 Pout << "\nView factor patches included in the calculation : "
362 << viewFactorsPatches << endl;
363 }
364
365 // Collect local Cf and Sf on coarse mesh
366 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
367
368 DynamicList<point> localCoarseCf(nCoarseFaces);
369 DynamicList<point> localCoarseSf(nCoarseFaces);
370 DynamicList<label> localAgg(nCoarseFaces);
371 labelHashSet includePatches;
372
373 for (const label patchID : viewFactorsPatches)
374 {
375 const polyPatch& pp = patches[patchID];
376 const labelList& agglom = finalAgglom[patchID];
377
378 includePatches.insert(patchID);
379
380 if (agglom.size() > 0)
381 {
382 label nAgglom = max(agglom)+1;
383 labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
384 const labelList& coarsePatchFace =
385 coarseMesh.patchFaceMap()[patchID];
386
387 const pointField& coarseCf =
388 coarseMesh.Cf().boundaryField()[patchID];
389 const pointField& coarseSf =
390 coarseMesh.Sf().boundaryField()[patchID];
391
392 forAll(coarseCf, faceI)
393 {
394 point cf = coarseCf[faceI];
395
396 const label coarseFaceI = coarsePatchFace[faceI];
397 const labelList& fineFaces = coarseToFine[coarseFaceI];
398 const label agglomI =
399 agglom[fineFaces[0]] + coarsePatches[patchID].start();
400
401 // Construct single face
403 (
404 UIndirectList<face>(pp, fineFaces),
405 pp.points()
406 );
407
408 List<point> availablePoints
409 (
410 upp.faceCentres().size()
411 + upp.localPoints().size()
412 );
413
415 (
416 availablePoints,
417 upp.faceCentres().size()
418 ) = upp.faceCentres();
419
421 (
422 availablePoints,
423 upp.localPoints().size(),
424 upp.faceCentres().size()
425 ) = upp.localPoints();
426
427 point cfo = cf;
428 scalar dist = GREAT;
429 forAll(availablePoints, iPoint)
430 {
431 point cfFine = availablePoints[iPoint];
432 if (mag(cfFine-cfo) < dist)
433 {
434 dist = mag(cfFine-cfo);
435 cf = cfFine;
436 }
437 }
438
439 point sf = coarseSf[faceI];
440 localCoarseCf.append(cf);
441 localCoarseSf.append(sf);
442 localAgg.append(agglomI);
443
444 }
445 }
446 }
447
448 // Distribute local coarse Cf and Sf for shooting rays
449 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
450
451 List<pointField> remoteCoarseCf(Pstream::nProcs());
452 List<pointField> remoteCoarseSf(Pstream::nProcs());
453 List<labelField> remoteCoarseAgg(Pstream::nProcs());
454
455 remoteCoarseCf[Pstream::myProcNo()] = localCoarseCf;
456 remoteCoarseSf[Pstream::myProcNo()] = localCoarseSf;
457 remoteCoarseAgg[Pstream::myProcNo()] = localAgg;
458
459 Pstream::allGatherList(remoteCoarseCf);
460 Pstream::allGatherList(remoteCoarseSf);
461 Pstream::allGatherList(remoteCoarseAgg);
462
463
464 globalIndex globalNumbering(nCoarseFaces);
465
466 // Set up searching engine for obstacles
467 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
468 #include "searchingEngine.H"
469
470 // Determine rays between coarse face centres
471 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
472 DynamicList<label> rayStartFace(nCoarseFaces + 0.01*nCoarseFaces);
473
474 DynamicList<label> rayEndFace(rayStartFace.size());
475
476
477 // Return rayStartFace in local index and rayEndFace in global index
478 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
479
480 #include "shootRays.H"
481
482 // Calculate number of visible faces from local index
483 labelList nVisibleFaceFaces(nCoarseFaces, Zero);
484
485 forAll(rayStartFace, i)
486 {
487 nVisibleFaceFaces[rayStartFace[i]]++;
488 }
489
490 labelListList visibleFaceFaces(nCoarseFaces);
491
492 label nViewFactors = 0;
493 forAll(nVisibleFaceFaces, faceI)
494 {
495 visibleFaceFaces[faceI].setSize(nVisibleFaceFaces[faceI]);
496 nViewFactors += nVisibleFaceFaces[faceI];
497 }
498
499 // - Construct compact numbering
500 // - return map from remote to compact indices
501 // (per processor (!= myProcNo) a map from remote index to compact index)
502 // - construct distribute map
503 // - renumber rayEndFace into compact addressing
504
505 List<Map<label>> compactMap(Pstream::nProcs());
506
507 mapDistribute map(globalNumbering, rayEndFace, compactMap);
508
509 // visibleFaceFaces has:
510 // (local face, local viewed face) = compact viewed face
511 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
512
513 nVisibleFaceFaces = 0;
514 forAll(rayStartFace, i)
515 {
516 label faceI = rayStartFace[i];
517 label compactI = rayEndFace[i];
518 visibleFaceFaces[faceI][nVisibleFaceFaces[faceI]++] = compactI;
519 }
520
521
522 // Construct data in compact addressing
523 // I need coarse Sf (Ai), fine Sf (dAi) and fine Cf(r) to calculate Fij
524 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
525
526 pointField compactCoarseCf(map.constructSize(), Zero);
527 pointField compactCoarseSf(map.constructSize(), Zero);
528 List<List<point>> compactFineSf(map.constructSize());
529 List<List<point>> compactFineCf(map.constructSize());
530
531 DynamicList<label> compactPatchId(map.constructSize());
532
533 // Insert my coarse local values
534 SubList<point>(compactCoarseSf, nCoarseFaces) = localCoarseSf;
535 SubList<point>(compactCoarseCf, nCoarseFaces) = localCoarseCf;
536
537 // Insert my fine local values
538 label compactI = 0;
539 forAll(viewFactorsPatches, i)
540 {
541 label patchID = viewFactorsPatches[i];
542
543 const labelList& agglom = finalAgglom[patchID];
544 if (agglom.size() > 0)
545 {
546 label nAgglom = max(agglom)+1;
547 labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
548 const labelList& coarsePatchFace =
549 coarseMesh.patchFaceMap()[patchID];
550
551 forAll(coarseToFine, coarseI)
552 {
553 compactPatchId.append(patchID);
554 List<point>& fineCf = compactFineCf[compactI];
555 List<point>& fineSf = compactFineSf[compactI++];
556
557 const label coarseFaceI = coarsePatchFace[coarseI];
558 const labelList& fineFaces = coarseToFine[coarseFaceI];
559
560 fineCf.setSize(fineFaces.size());
561 fineSf.setSize(fineFaces.size());
562
563 fineCf = UIndirectList<point>
564 (
565 mesh.Cf().boundaryField()[patchID],
566 coarseToFine[coarseFaceI]
567 );
568 fineSf = UIndirectList<point>
569 (
570 mesh.Sf().boundaryField()[patchID],
571 coarseToFine[coarseFaceI]
572 );
573 }
574 }
575 }
576
577 // Do all swapping
578 map.distribute(compactCoarseSf);
579 map.distribute(compactCoarseCf);
580 map.distribute(compactFineCf);
581 map.distribute(compactFineSf);
582
583 map.distribute(compactPatchId);
584
585
586 // Plot all rays between visible faces.
587 if (dumpRays)
588 {
589 writeRays
590 (
591 runTime.path()/"allVisibleFaces.obj",
592 compactCoarseCf,
593 remoteCoarseCf[Pstream::myProcNo()],
594 visibleFaceFaces
595 );
596 }
597
598
599 // Fill local view factor matrix
600 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
601
603 (
605 (
606 "F",
607 mesh.facesInstance(),
608 mesh,
609 IOobject::NO_READ,
610 IOobject::NO_WRITE,
611 false
612 ),
613 nCoarseFaces
614 );
615
616 label totalPatches = coarsePatches.size();
617 reduce(totalPatches, maxOp<label>());
618
619 // Matrix sum in j(Fij) for each i (if enclosure sum = 1)
620 scalarSquareMatrix sumViewFactorPatch
621 (
622 totalPatches,
623 0.0
624 );
625
626 scalarList patchArea(totalPatches, Zero);
627
628 if (Pstream::master())
629 {
630 Info<< "\nCalculating view factors..." << endl;
631 }
632
633 if (mesh.nSolutionD() == 3)
634 {
635 forAll(localCoarseSf, coarseFaceI)
636 {
637 const List<point>& localFineSf = compactFineSf[coarseFaceI];
638 const vector Ai = sum(localFineSf);
639 const List<point>& localFineCf = compactFineCf[coarseFaceI];
640 const label fromPatchId = compactPatchId[coarseFaceI];
641 patchArea[fromPatchId] += mag(Ai);
642
643 const labelList& visCoarseFaces = visibleFaceFaces[coarseFaceI];
644
645 forAll(visCoarseFaces, visCoarseFaceI)
646 {
647 F[coarseFaceI].setSize(visCoarseFaces.size());
648 label compactJ = visCoarseFaces[visCoarseFaceI];
649 const List<point>& remoteFineSj = compactFineSf[compactJ];
650 const List<point>& remoteFineCj = compactFineCf[compactJ];
651
652 const label toPatchId = compactPatchId[compactJ];
653
654 scalar Fij = 0;
655 forAll(localFineSf, i)
656 {
657 const vector& dAi = localFineSf[i];
658 const vector& dCi = localFineCf[i];
659
660 forAll(remoteFineSj, j)
661 {
662 const vector& dAj = remoteFineSj[j];
663 const vector& dCj = remoteFineCj[j];
664
665 scalar dIntFij = calculateViewFactorFij
666 (
667 dCi,
668 dCj,
669 dAi,
670 dAj
671 );
672
673 Fij += dIntFij;
674 }
675 }
676 F[coarseFaceI][visCoarseFaceI] = Fij/mag(Ai);
677 sumViewFactorPatch[fromPatchId][toPatchId] += Fij;
678 }
679 }
680 }
681 else if (mesh.nSolutionD() == 2)
682 {
683 const boundBox& box = mesh.bounds();
684 const Vector<label>& dirs = mesh.geometricD();
685 vector emptyDir = Zero;
686 forAll(dirs, i)
687 {
688 if (dirs[i] == -1)
689 {
690 emptyDir[i] = 1.0;
691 }
692 }
693
694 scalar wideBy2 = (box.span() & emptyDir)*2.0;
695
696 forAll(localCoarseSf, coarseFaceI)
697 {
698 const vector& Ai = localCoarseSf[coarseFaceI];
699 const vector& Ci = localCoarseCf[coarseFaceI];
700 vector Ain = Ai/mag(Ai);
701 vector R1i = Ci + (mag(Ai)/wideBy2)*(Ain ^ emptyDir);
702 vector R2i = Ci - (mag(Ai)/wideBy2)*(Ain ^ emptyDir) ;
703
704 const label fromPatchId = compactPatchId[coarseFaceI];
705 patchArea[fromPatchId] += mag(Ai);
706
707 const labelList& visCoarseFaces = visibleFaceFaces[coarseFaceI];
708 forAll(visCoarseFaces, visCoarseFaceI)
709 {
710 F[coarseFaceI].setSize(visCoarseFaces.size());
711 label compactJ = visCoarseFaces[visCoarseFaceI];
712 const vector& Aj = compactCoarseSf[compactJ];
713 const vector& Cj = compactCoarseCf[compactJ];
714
715 const label toPatchId = compactPatchId[compactJ];
716
717 vector Ajn = Aj/mag(Aj);
718 vector R1j = Cj + (mag(Aj)/wideBy2)*(Ajn ^ emptyDir);
719 vector R2j = Cj - (mag(Aj)/wideBy2)*(Ajn ^ emptyDir);
720
721 scalar d1 = mag(R1i - R2j);
722 scalar d2 = mag(R2i - R1j);
723 scalar s1 = mag(R1i - R1j);
724 scalar s2 = mag(R2i - R2j);
725
726 scalar Fij = mag((d1 + d2) - (s1 + s2))/(4.0*mag(Ai)/wideBy2);
727
728 F[coarseFaceI][visCoarseFaceI] = Fij;
729 sumViewFactorPatch[fromPatchId][toPatchId] += Fij*mag(Ai);
730 }
731 }
732 }
733
734 if (Pstream::master())
735 {
736 Info << "Writing view factor matrix..." << endl;
737 }
738
739 // Write view factors matrix in listlist form
740 F.write();
741
742 reduce(sumViewFactorPatch, sumOp<scalarSquareMatrix>());
743 reduce(patchArea, sumOp<scalarList>());
744
745
746 if (Pstream::master() && debug)
747 {
748 forAll(viewFactorsPatches, i)
749 {
750 label patchI = viewFactorsPatches[i];
751 forAll(viewFactorsPatches, i)
752 {
753 label patchJ = viewFactorsPatches[i];
754 Info << "F" << patchI << patchJ << ": "
755 << sumViewFactorPatch[patchI][patchJ]/patchArea[patchI]
756 << endl;
757 }
758 }
759 }
760
761
762 if (writeViewFactors)
763 {
764 volScalarField viewFactorField
765 (
767 (
768 "viewFactorField",
769 mesh.time().timeName(),
770 mesh,
771 IOobject::NO_READ,
772 IOobject::NO_WRITE
773 ),
774 mesh,
775 dimensionedScalar(dimless, Zero)
776 );
777
778 label compactI = 0;
779
780 volScalarField::Boundary& vfbf = viewFactorField.boundaryFieldRef();
781 forAll(viewFactorsPatches, i)
782 {
783 label patchID = viewFactorsPatches[i];
784 const labelList& agglom = finalAgglom[patchID];
785 if (agglom.size() > 0)
786 {
787 label nAgglom = max(agglom)+1;
788 labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
789 const labelList& coarsePatchFace =
790 coarseMesh.patchFaceMap()[patchID];
791
792 forAll(coarseToFine, coarseI)
793 {
794 const scalar Fij = sum(F[compactI]);
795 const label coarseFaceID = coarsePatchFace[coarseI];
796 const labelList& fineFaces = coarseToFine[coarseFaceID];
797 forAll(fineFaces, fineId)
798 {
799 const label faceID = fineFaces[fineId];
800 vfbf[patchID][faceID] = Fij;
801 }
802 compactI++;
803 }
804 }
805 }
806 viewFactorField.write();
807 }
808
809
810 // Invert compactMap (from processor+localface to compact) to go
811 // from compact to processor+localface (expressed as a globalIndex)
812 // globalIndex globalCoarFaceNum(coarseMesh.nFaces());
813 labelList compactToGlobal(map.constructSize());
814
815 // Local indices first (note: are not in compactMap)
816 for (label i = 0; i < globalNumbering.localSize(); i++)
817 {
818 compactToGlobal[i] = globalNumbering.toGlobal(i);
819 }
820
821
822 forAll(compactMap, procI)
823 {
824 const Map<label>& localToCompactMap = compactMap[procI];
825
826 forAllConstIters(localToCompactMap, iter)
827 {
828 compactToGlobal[*iter] = globalNumbering.toGlobal
829 (
830 procI,
831 iter.key()
832 );
833 }
834 }
835
836
837 labelListList globalFaceFaces(visibleFaceFaces.size());
838
839 // Create globalFaceFaces needed to insert view factors
840 // in F to the global matrix Fmatrix
841 forAll(globalFaceFaces, faceI)
842 {
843 globalFaceFaces[faceI] = renumber
844 (
845 compactToGlobal,
846 visibleFaceFaces[faceI]
847 );
848 }
849
850 labelListIOList IOglobalFaceFaces
851 (
853 (
854 "globalFaceFaces",
855 mesh.facesInstance(),
856 mesh,
857 IOobject::NO_READ,
858 IOobject::NO_WRITE,
859 false
860 ),
861 std::move(globalFaceFaces)
862 );
863 IOglobalFaceFaces.write();
864
865
866 IOmapDistribute IOmapDist
867 (
869 (
870 "mapDist",
871 mesh.facesInstance(),
872 mesh,
873 IOobject::NO_READ,
874 IOobject::AUTO_WRITE
875 ),
876 std::move(map)
877 );
878
879 IOmapDist.write();
880
881 Info<< "End\n" << endl;
882 return 0;
883}
884
885
886// ************************************************************************* //
Y[inertIndex] max(0.0)
reduce(hasMovingMesh, orOp< bool >())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOmapDistribute is derived from mapDistribute and IOobject to give the mapDistribute automatic IO fun...
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
void setSize(const label n)
Alias for resize()
Definition: List.H:218
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
Output to file stream, using an OSstream.
Definition: OFstream.H:57
virtual const fileName & name() const
Get the name of the stream.
Definition: OSstream.H:107
A list of faces which address into the list of points.
const Field< point_type > & points() const noexcept
Return reference to global points.
A List obtained as a section of another List.
Definition: SubList.H:70
A List with indirect addressing. Like IndirectList but does not store addressing.
Definition: IndirectList.H:79
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:127
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
A class for handling file names.
Definition: fileName.H:76
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
label localSize() const
My local size.
Definition: globalIndexI.H:207
label toGlobal(const label i) const
From local to global index.
Definition: globalIndexI.H:260
A triFace with additional (region) index.
Definition: labelledTri.H:60
Class containing processor-to-processor mapping information.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
label start() const noexcept
The start label of boundary faces in the polyMesh face list.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
fvMesh as subset of other mesh. Consists of one cell and all original bounday faces....
Triangulated surface description with patch information.
Definition: triSurface.H:79
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition: word.H:68
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
engineTime & runTime
Required Variables.
const pointField & points
volVectorField F(fluid.F())
int debug
Static debugging option.
const wordList surface
Standard surface field types (scalar, vector, tensor, etc)
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
Information stream (stdout output on master, null elsewhere)
IntListType renumber(const labelUList &oldToNew, const IntListType &input)
Renumber the values (not the indices) of a list.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:114
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
labelList triSurfaceToAgglom(5 *nFineFaces)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
Foam::surfaceFields.
Unit conversion functions.