Loading...
Searching...
No Matches
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-2024 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 This view factors generation application uses a combined approach of
35 double area integral (2AI) and double linear integral (2LI). 2AI is used
36 when the two surfaces are 'far' apart and 2LI when they are 'close'.
37 2LI is integrated along edges using Gaussian quadrature.
38 The distance between faces is calculating a ratio between averaged areas
39 and the distance between face centres.
40
41 The input from viewFactorsDict are:
42 \verbatim
43 GaussQuadTol 0.1; // GaussQuad error
44 distTol 8; // R/Average(rm)
45 alpha 0.22; // Use for common edges for 2LI
46 \endverbatim
47
48
49 For debugging purposes, the following entries can be set in viewFactorsDict:
50 \verbatim
51 writeViewFactorMatrix true;
52 writeFacesAgglomeration false;
53 dumpRays false;
54
55 writeViewFactorMatrix writes the sum of the VF on each face.
56 writeFacesAgglomeration writes the agglomeration
57 dumpRays dumps rays
58 \endverbatim
59
60 The participating patches in the VF calculation have to be in the
61 'viewFactorWall' patch group (in the polyMesh/boundary file), e.g.
62
63 \verbatim
64 floor
65 {
66 type wall;
67 inGroups 2(wall viewFactorWall);
68 nFaces 100;
69 startFace 3100;
70 }
71 \endverbatim
72
73 Compile with -DNO_CGAL only if no CGAL present - CGAL AABB tree performs
74 better than the built-in octree.
75
76\*---------------------------------------------------------------------------*/
77
78#include "argList.H"
79#include "fvMesh.H"
80#include "volFields.H"
81#include "surfaceFields.H"
83#include "meshTools.H"
84#include "constants.H"
85
87#include "DynamicField.H"
88
89#include "scalarMatrices.H"
90#include "labelListIOList.H"
91#include "scalarListIOList.H"
92
93#include "singleCellFvMesh.H"
94#include "IOmapDistribute.H"
95
96#ifndef NO_CGAL
97
98// Silence boost bind deprecation warnings (before CGAL-5.2.1)
99#include "CGAL/version.h"
100#if defined(CGAL_VERSION_NR) && (CGAL_VERSION_NR < 1050211000)
101#define BOOST_BIND_GLOBAL_PLACEHOLDERS
102#endif
103#pragma clang diagnostic ignored "-Wbitwise-instead-of-logical"
104#pragma clang diagnostic ignored "-Wdeprecated-builtins"
105#pragma clang diagnostic ignored "-Wdeprecated-declarations"
106
107#include <vector>
108#include <CGAL/Simple_cartesian.h>
109#include <CGAL/AABB_tree.h>
110#if defined(CGAL_VERSION_NR) && (CGAL_VERSION_NR < 1060011000)
111#include <CGAL/AABB_traits.h>
112#include <CGAL/AABB_triangle_primitive.h>
113#else
114#include <CGAL/AABB_traits_3.h>
115#include <CGAL/AABB_triangle_primitive_3.h>
116#endif
117#include <CGAL/Surface_mesh.h>
118
119typedef CGAL::Simple_cartesian<double> K;
120typedef K::Point_3 Point;
121typedef K::Direction_3 Vector3C;
122typedef K::Triangle_3 Triangle;
123typedef K::Segment_3 Segment;
124
125typedef std::vector<Triangle>::iterator Iterator;
126#if defined(CGAL_VERSION_NR) && (CGAL_VERSION_NR < 1060011000)
127typedef CGAL::AABB_triangle_primitive<K, Iterator> Primitive;
128typedef CGAL::AABB_traits<K, Primitive> AABB_triangle_traits;
129#else
130typedef CGAL::AABB_triangle_primitive_3<K, Iterator> Primitive;
131typedef CGAL::AABB_traits_3<K, Primitive> AABB_triangle_traits;
132#endif
133typedef CGAL::AABB_tree<AABB_triangle_traits> Tree;
134
135// Used boost::optional prior to CGAL-6.0
136#if (CGAL_VERSION_NR < 1060000910)
137typedef boost::optional
138#else
139typedef std::optional
140#endif
141<
142 Tree::Intersection_and_primitive_id<Segment>::Type
143> Segment_intersection;
144
145#endif // NO_CGAL
146
147using namespace Foam;
148using namespace Foam::constant;
149using namespace Foam::constant::mathematical;
150
151triSurface triangulate
152(
153 const polyBoundaryMesh& bMesh,
154 const labelHashSet& includePatches,
155 const labelListIOList& finalAgglom,
157 const globalIndex& globalNumbering,
158 const polyBoundaryMesh& coarsePatches
159)
160{
161 const polyMesh& mesh = bMesh.mesh();
162
163 // Storage for surfaceMesh. Size estimate.
164 DynamicList<labelledTri> triangles(mesh.nBoundaryFaces());
165
166 label newPatchI = 0;
167 label localTriFaceI = 0;
168
169 for (const label patchI : includePatches)
170 {
171 const polyPatch& patch = bMesh[patchI];
172 const pointField& points = patch.points();
173
174 //label nTriTotal = 0;
175
176 forAll(patch, patchFaceI)
177 {
178 const face& f = patch[patchFaceI];
179
180 faceList triFaces(f.nTriangles(points));
181
182 label nTri = 0;
183
184 f.triangles(points, nTri, triFaces);
185
186 forAll(triFaces, triFaceI)
187 {
188 const face& f = triFaces[triFaceI];
189
190 triangles.append(labelledTri(f[0], f[1], f[2], newPatchI));
191
192 //nTriTotal++;
193
194 triSurfaceToAgglom[localTriFaceI++] = globalNumbering.toGlobal
195 (
197 finalAgglom[patchI][patchFaceI]
198 + coarsePatches[patchI].start()
199 );
200 }
201 }
202
203 newPatchI++;
204 }
205
206 //triSurfaceToAgglom.resize(localTriFaceI-1);
207
208 triangles.shrink();
209 triSurface surface(triangles, mesh.points());
210 surface.compactPoints();
211
212
213#ifndef NO_CGAL
214
215 // CGAL : every processor has whole surface
216
217 const globalIndex globalFaceIdx
218 (
220 surface.size()
221 );
222 const globalIndex globalPointIdx
223 (
225 surface.points().size()
226 );
227
228 List<labelledTri> globalSurfaceTris(globalFaceIdx.gather(surface));
229 pointField globalSurfacePoints(globalPointIdx.gather(surface.points()));
230
231 //label offset = 0;
232 for (const label proci : globalPointIdx.allProcs())
233 {
234 const label offset = globalPointIdx.localStart(proci);
235
236 if (offset)
237 {
238 for
239 (
240 labelledTri& tri
241 : globalSurfaceTris.slice(globalFaceIdx.range(proci))
242 )
243 {
244 tri[0] += offset;
245 tri[1] += offset;
246 tri[2] += offset;
247 }
248 }
249 }
250
251 surface =
253 (
254 std::move(globalSurfaceTris),
255 std::move(globalSurfacePoints)
256 );
257
258 Pstream::broadcast(surface);
259#endif
260
261 // Add patch names to surface
262 surface.patches().setSize(newPatchI);
263
264 newPatchI = 0;
265
266 for (const label patchI : includePatches)
267 {
268 const polyPatch& patch = bMesh[patchI];
269
270 surface.patches()[newPatchI].index() = patchI;
271 surface.patches()[newPatchI].name() = patch.name();
272 surface.patches()[newPatchI].geometricType() = patch.type();
273
274 newPatchI++;
275 }
276
277 return surface;
278}
279
280
281void writeRays
282(
283 const fileName& fName,
284 const pointField& compactCf,
285 const pointField& myFc,
286 const labelListList& visibleFaceFaces
287)
288{
289 OFstream str(fName);
290 label vertI = 0;
291
292 Pout<< "Dumping rays to " << str.name() << endl;
293
294 forAll(myFc, faceI)
295 {
296 const labelList visFaces = visibleFaceFaces[faceI];
297 forAll(visFaces, faceRemote)
298 {
299 label compactI = visFaces[faceRemote];
300 const point& remoteFc = compactCf[compactI];
301
302 meshTools::writeOBJ(str, myFc[faceI]);
303 vertI++;
304 meshTools::writeOBJ(str, remoteFc);
305 vertI++;
306 str << "l " << vertI-1 << ' ' << vertI << nl;
307 }
308 }
309 str.flush();
310}
311
312
313scalar calculateViewFactorFij2AI
314(
315 const vector& i,
316 const vector& j,
317 const vector& dAi,
318 const vector& dAj
319)
320{
321 vector r = i - j;
322 scalar rMag = mag(r);
323
324 if (rMag > SMALL)
325 {
326 scalar dAiMag = mag(dAi);
327 scalar dAjMag = mag(dAj);
328
329 vector ni = dAi/dAiMag;
330 vector nj = dAj/dAjMag;
331 scalar cosThetaJ = mag(nj & r)/rMag;
332 scalar cosThetaI = mag(ni & r)/rMag;
333
334 return
335 (
336 (cosThetaI*cosThetaJ*dAjMag*dAiMag)
338 );
339 }
340 else
341 {
342 return 0;
343 }
344}
345
346
347void insertMatrixElements
348(
349 const globalIndex& globalNumbering,
350 const label fromProcI,
351 const labelListList& globalFaceFaces,
352 const scalarListList& viewFactors,
353 scalarSquareMatrix& matrix
354)
355{
356 forAll(viewFactors, faceI)
357 {
358 const scalarList& vf = viewFactors[faceI];
359 const labelList& globalFaces = globalFaceFaces[faceI];
360
361 label globalI = globalNumbering.toGlobal(fromProcI, faceI);
362 forAll(globalFaces, i)
363 {
364 matrix[globalI][globalFaces[i]] = vf[i];
365 }
366 }
367}
368
369
370scalar GaussQuad
371(
372
373 const scalarList& w,
374 const scalarList& p,
375 const scalar& magSi,
376 const scalar& magSj,
377 const vector& di,
378 const vector& dj,
379 const vector& ci,
380 const vector& cj,
381 const scalar cosij,
382 const scalar alpha,
383 label gi
384)
385{
386 scalar dIntFij = 0;
387 if (gi == 0)
388 {
389 vector r(ci - cj);
390 if (mag(r) < SMALL)
391 {
392 r = (alpha*magSi)*di;
393 }
394 dIntFij = max(cosij*Foam::log(r&r)*magSi*magSj, 0);
395 }
396 else
397 {
398 List<vector> pi(w.size());
399 forAll (pi, i)
400 {
401 pi[i] = ci + p[i]*(magSi/2)*di;
402 }
403
404 List<vector> pj(w.size());
405 forAll (pj, i)
406 {
407 pj[i] = cj + p[i]*(magSj/2)*dj;
408 }
409
410 forAll (w, i)
411 {
412 forAll (w, j)
413 {
414 vector r(pi[i] - pj[j]);
415
416 if (mag(r) < SMALL)
417 {
418 r = (alpha*magSi)*di;
419 dIntFij +=
420 cosij*w[i]*w[j]*Foam::log(r&r);
421 }
422 else
423 {
424 dIntFij +=
425 cosij*w[i]*w[j]*Foam::log(r&r);
426 }
427
428 }
429 }
430
431 dIntFij *= (magSi/2) * (magSj/2);
432
433 }
434 return dIntFij;
435}
436
437
438// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
439
440int main(int argc, char *argv[])
441{
443 (
444 "Calculate view factors from face agglomeration array."
445 " The finalAgglom generated by faceAgglomerate utility."
446 );
447
448 #include "addRegionOption.H"
449 #include "setRootCase.H"
450 #include "createTime.H"
451 #include "createNamedMesh.H"
452
453 // Read view factor dictionary
454 IOdictionary viewFactorDict
455 (
457 (
458 "viewFactorsDict",
459 runTime.constant(),
460 mesh,
463 )
464 );
465
466 const word viewFactorWall("viewFactorWall");
467
468 const bool writeViewFactors =
469 viewFactorDict.getOrDefault("writeViewFactorMatrix", false);
470
471 const bool dumpRays =
472 viewFactorDict.getOrDefault("dumpRays", false);
473
474 const label debug = viewFactorDict.getOrDefault<label>("debug", 0);
475
476 const scalar GaussQuadTol =
477 viewFactorDict.getOrDefault<scalar>("GaussQuadTol", 0.01);
478
479 const scalar distTol =
480 viewFactorDict.getOrDefault<scalar>("distTol", 8);
481
482 const scalar alpha =
483 viewFactorDict.getOrDefault<scalar>("alpha", 0.21);
484
485 const scalar intTol =
486 viewFactorDict.getOrDefault<scalar>("intTol", 1e-2);
487
488 bool useAgglomeration(true);
489
490 const polyBoundaryMesh& patches = mesh.boundaryMesh();
491 const labelList viewFactorsPatches(patches.indices(viewFactorWall));
492
493 // Read agglomeration map
494 labelListIOList finalAgglom
495 (
497 (
498 "finalAgglom",
499 mesh.facesInstance(),
500 mesh,
504 )
505 );
506
507 if (!finalAgglom.typeHeaderOk<labelListIOList>())
508 {
509 finalAgglom.setSize(patches.size());
510 for (label patchi=0; patchi < patches.size(); patchi++)
511 {
512 finalAgglom[patchi] = identity(patches[patchi].size());
513 }
514 useAgglomeration = false;
515 }
516
517 // Create the coarse mesh using agglomeration
518 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
519
520 if (debug)
521 {
522 Pout << "\nCreating single cell mesh..." << endl;
523 }
524
525 singleCellFvMesh coarseMesh
526 (
528 (
529 "coarse:" + mesh.name(),
530 runTime.timeName(),
531 runTime,
534 ),
535 mesh,
536 finalAgglom
537 );
538
539 if (debug)
540 {
541 Pout << "\nCreated single cell mesh..." << endl;
542 }
543
544
545 // Calculate total number of fine and coarse faces
546 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
547
548 label nCoarseFaces = 0; //total number of coarse faces
549 label nFineFaces = 0; //total number of fine faces
550
551 const polyBoundaryMesh& coarsePatches = coarseMesh.boundaryMesh();
552
553 for (const label patchi : viewFactorsPatches)
554 {
555 nCoarseFaces += coarsePatches[patchi].size();
556 nFineFaces += patches[patchi].size();
557 }
558
559
560 Info<< "\nTotal number of coarse faces: "
561 << returnReduce(nCoarseFaces, sumOp<label>())
562 << endl;
563
564 if (Pstream::master() && debug)
565 {
566 Pout << "\nView factor patches included in the calculation : "
567 << viewFactorsPatches << endl;
568 }
569
570 // Collect local Cf and Sf on coarse mesh
571 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
572
573 DynamicList<point> localCoarseCf(nCoarseFaces);
574 DynamicList<point> localCoarseSf(nCoarseFaces);
575 DynamicList<label> localAgg(nCoarseFaces);
576 labelHashSet includePatches;
577
578 for (const label patchID : viewFactorsPatches)
579 {
580 const polyPatch& pp = patches[patchID];
581 const labelList& agglom = finalAgglom[patchID];
582
583 includePatches.insert(patchID);
584
585 if (agglom.size() > 0)
586 {
587 label nAgglom = max(agglom)+1;
588 labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
589 const labelList& coarsePatchFace =
590 coarseMesh.patchFaceMap()[patchID];
591
592 const pointField& coarseCf =
593 coarseMesh.Cf().boundaryField()[patchID];
594 const pointField& coarseSf =
595 coarseMesh.Sf().boundaryField()[patchID];
596
597 forAll(coarseCf, faceI)
598 {
599 point cf = coarseCf[faceI];
600
601 const label coarseFaceI = coarsePatchFace[faceI];
602 const labelList& fineFaces = coarseToFine[coarseFaceI];
603 const label agglomI =
604 agglom[fineFaces[0]] + coarsePatches[patchID].start();
605
606 // Construct single face
608 (
609 UIndirectList<face>(pp, fineFaces),
610 pp.points()
611 );
612
613 List<point> availablePoints
614 (
615 upp.faceCentres().size()
616 + upp.localPoints().size()
617 );
618
620 (
621 availablePoints,
622 upp.faceCentres().size()
623 ) = upp.faceCentres();
624
626 (
627 availablePoints,
628 upp.localPoints().size(),
629 upp.faceCentres().size()
630 ) = upp.localPoints();
631
632 point cfo = cf;
633 scalar dist = GREAT;
634 forAll(availablePoints, iPoint)
635 {
636 point cfFine = availablePoints[iPoint];
637 if (mag(cfFine-cfo) < dist)
638 {
639 dist = mag(cfFine-cfo);
640 cf = cfFine;
641 }
642 }
643
644 point sf = coarseSf[faceI];
645 localCoarseCf.append(cf);
646 localCoarseSf.append(sf);
647 localAgg.append(agglomI);
648
649 }
650 }
651 }
652
653 // Distribute local coarse Cf and Sf for shooting rays
654 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
655
656 List<pointField> remoteCoarseCf(Pstream::nProcs());
657 List<pointField> remoteCoarseSf(Pstream::nProcs());
658 List<labelField> remoteCoarseAgg(Pstream::nProcs());
659
660 remoteCoarseCf[Pstream::myProcNo()] = localCoarseCf;
661 remoteCoarseSf[Pstream::myProcNo()] = localCoarseSf;
662 remoteCoarseAgg[Pstream::myProcNo()] = localAgg;
663
664 Pstream::allGatherList(remoteCoarseCf);
665 Pstream::allGatherList(remoteCoarseSf);
666 Pstream::allGatherList(remoteCoarseAgg);
667
668
669 globalIndex globalNumbering(nCoarseFaces);
670
671 // Set up searching engine for obstacles
672 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
673#ifdef NO_CGAL
674 // Using octree
675 #include "searchingEngine.H"
676#else
677 // Using CGAL aabbtree (faster, more robust)
678 #include "searchingEngine_CGAL.H"
679#endif
680
681 // Determine rays between coarse face centres
682 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
683 DynamicList<label> rayStartFace(nCoarseFaces + 0.01*nCoarseFaces);
684
685 DynamicList<label> rayEndFace(rayStartFace.size());
686
687 // Return rayStartFace in local index and rayEndFace in global index
688 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
689#ifdef NO_CGAL
690 // Using octree, distributedTriSurfaceMesh
691 #include "shootRays.H"
692#else
693 // Using CGAL aabbtree (faster, more robust)
694 #include "shootRays_CGAL.H"
695#endif
696
697 // Calculate number of visible faces from local index
698 labelList nVisibleFaceFaces(nCoarseFaces, Zero);
699
700 forAll(rayStartFace, i)
701 {
702 nVisibleFaceFaces[rayStartFace[i]]++;
703 }
704
705 labelListList visibleFaceFaces(nCoarseFaces);
706
707 //label nViewFactors = 0;
708 forAll(nVisibleFaceFaces, faceI)
709 {
710 visibleFaceFaces[faceI].setSize(nVisibleFaceFaces[faceI]);
711 //nViewFactors += nVisibleFaceFaces[faceI];
712 }
713
714 // - Construct compact numbering
715 // - return map from remote to compact indices
716 // (per processor (!= myProcNo) a map from remote index to compact index)
717 // - construct distribute map
718 // - renumber rayEndFace into compact addressing
719
720 List<Map<label>> compactMap(Pstream::nProcs());
721
722 mapDistribute map(globalNumbering, rayEndFace, compactMap);
723
724 // visibleFaceFaces has:
725 // (local face, local viewed face) = compact viewed face
726 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
727
728 nVisibleFaceFaces = 0;
729 forAll(rayStartFace, i)
730 {
731 label faceI = rayStartFace[i];
732 label compactI = rayEndFace[i];
733 visibleFaceFaces[faceI][nVisibleFaceFaces[faceI]++] = compactI;
734 }
735
736 // Construct data in compact addressing
737 // (2AA) need coarse (Ai), fine Sf (dAi) and fine Cf(r) to calculate Fij
738 // (2LI) need edges (li)
739 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
740
741 pointField compactCoarseCf(map.constructSize(), Zero);
742 pointField compactCoarseSf(map.constructSize(), Zero);
743 List<List<point>> compactFineSf(map.constructSize());
744 List<List<point>> compactFineCf(map.constructSize());
745
746 DynamicList<List<point>> compactPoints(map.constructSize());
747
748 DynamicList<label> compactPatchId(map.constructSize());
749
750 // Insert my coarse local values
751 SubList<point>(compactCoarseSf, nCoarseFaces) = localCoarseSf;
752 SubList<point>(compactCoarseCf, nCoarseFaces) = localCoarseCf;
753
754 const faceList& faces = mesh.faces();
755
756 // Insert my fine local values
757 label compactI = 0;
758 forAll(viewFactorsPatches, i)
759 {
760 label patchID = viewFactorsPatches[i];
761
762 const labelList& agglom = finalAgglom[patchID];
763 if (agglom.size() > 0)
764 {
765 label nAgglom = max(agglom)+1;
766 labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
767 const labelList& coarsePatchFace =
768 coarseMesh.patchFaceMap()[patchID];
769
770 const polyPatch& pp = patches[patchID];
771
772 forAll(coarseToFine, coarseI)
773 {
774 compactPatchId.append(patchID);
775 List<point>& fineCf = compactFineCf[compactI];
776 List<point>& fineSf = compactFineSf[compactI];
777
778 label startFace = pp.start();
779
780 const vectorField locPoints
781 (
782 mesh.points(),
783 faces[coarseI + startFace]
784 );
785
786 const label coarseFaceI = coarsePatchFace[coarseI];
787 const labelList& fineFaces = coarseToFine[coarseFaceI];
788
789 fineCf.setSize(fineFaces.size());
790 fineSf.setSize(fineFaces.size());
791
792 compactPoints.append(locPoints);
793
794 fineCf = UIndirectList<point>
795 (
796 mesh.Cf().boundaryField()[patchID],
797 coarseToFine[coarseFaceI]
798 );
799 fineSf = UIndirectList<point>
800 (
801 mesh.Sf().boundaryField()[patchID],
802 coarseToFine[coarseFaceI]
803 );
804
805 compactI++;
806 }
807 }
808 }
809
810 if (Pstream::master() && debug)
811 {
812 Info<< "map distribute..." << endl;
813 }
814
815 // Do all swapping
816 map.distribute(compactCoarseSf);
817 map.distribute(compactCoarseCf);
818 map.distribute(compactFineCf);
819 map.distribute(compactFineSf);
820 map.distribute(compactPoints);
821 map.distribute(compactPatchId);
822
823 // Plot all rays between visible faces.
824 if (dumpRays)
825 {
826 writeRays
827 (
828 runTime.path()/"allVisibleFaces.obj",
829 compactCoarseCf,
830 remoteCoarseCf[Pstream::myProcNo()],
831 visibleFaceFaces
832 );
833 }
834
835
836 // Fill local view factor matrix
837 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
839 (
841 (
842 "F",
843 mesh.facesInstance(),
844 mesh,
848 ),
849 nCoarseFaces
850 );
851
852 const label totalPatches =
853 returnReduce(coarsePatches.size(), maxOp<label>());
854
855 // Matrix sum in j(Fij) for each i (if enclosure sum = 1)
856 scalarSquareMatrix sumViewFactorPatch(totalPatches, Zero);
857
858 scalarList patchArea(totalPatches, Zero);
859
860 if (Pstream::master())
861 {
862 Info<< "\nCalculating view factors..." << endl;
863 }
864
865 FixedList<scalarList, 5> GaussPoints;
866 GaussPoints[0].setSize(1);
867 GaussPoints[0] = 0;
868
869 GaussPoints[1].setSize(2);
870 GaussPoints[1][0] = 1/std::sqrt(3);
871 GaussPoints[1][1] = -1/std::sqrt(3);
872
873 GaussPoints[2].setSize(3);
874 GaussPoints[2][0] = 0;
875 GaussPoints[2][1] = std::sqrt(3.0/5.0);
876 GaussPoints[2][2] = -std::sqrt(3.0/5.0);
877
878 GaussPoints[3].setSize(4);
879 GaussPoints[3][0] = std::sqrt(3.0/7.0 - (2.0/7.0)*std::sqrt(6.0/5.0));
880 GaussPoints[3][1] = -GaussPoints[3][0];
881 GaussPoints[3][2] = std::sqrt(3.0/7.0 + (2.0/7.0)*std::sqrt(6.0/5.0));
882 GaussPoints[3][3] = -GaussPoints[3][2];
883
884 GaussPoints[4].setSize(5);
885 GaussPoints[4][0] = 0;
886 GaussPoints[4][1] = (1.0/3.0)*std::sqrt(5.0 - 2.0*std::sqrt(10.0/7.0));
887 GaussPoints[4][2] = -GaussPoints[4][1];
888 GaussPoints[4][3] = (1.0/3.0)*std::sqrt(5.0 + 2.0*std::sqrt(10.0/7.0));
889 GaussPoints[4][4] = -GaussPoints[4][3];
890
891
892 FixedList<scalarList, 5> GaussWeights;
893 GaussWeights[0].setSize(1);
894 GaussWeights[0] = 2;
895
896 GaussWeights[1].setSize(2);
897 GaussWeights[1][0] = 1;
898 GaussWeights[1][1] = 1;
899
900 GaussWeights[2].setSize(3);
901 GaussWeights[2][0] = 8.0/9.0;
902 GaussWeights[2][1] = 5.0/9.0;
903 GaussWeights[2][2] = 5.0/9.0;
904
905 GaussWeights[3].setSize(4);
906 GaussWeights[3][0] = (18.0 + std::sqrt(30))/36.0;
907 GaussWeights[3][1] = (18.0 + std::sqrt(30))/36.0;
908 GaussWeights[3][2] = (18.0 - std::sqrt(30))/36.0;
909 GaussWeights[3][3] = (18.0 - std::sqrt(30))/36.0;
910
911 GaussWeights[4].setSize(5);
912 GaussWeights[4][0] = 128.0/225.0;
913 GaussWeights[4][1] = (322.0 + 13.0*std::sqrt(70))/900.0;
914 GaussWeights[4][2] = (322.0 + 13.0*std::sqrt(70))/900.0;
915 GaussWeights[4][3] = (322.0 - 13.0*std::sqrt(70))/900.0;
916 GaussWeights[4][4] = (322.0 - 13.0*std::sqrt(70))/900.0;
917
918 const label maxQuadOrder = 5;
919
920 if (mesh.nSolutionD() == 3)
921 {
922 forAll(localCoarseSf, coarseFaceI)
923 {
924 const List<point>& localFineSf = compactFineSf[coarseFaceI];
925 const vector Ai = sum(localFineSf);
926 const List<point>& localFineCf = compactFineCf[coarseFaceI];
927 const label fromPatchId = compactPatchId[coarseFaceI];
928
929 const List<point>& lPoints = compactPoints[coarseFaceI];
930
931 patchArea[fromPatchId] += mag(Ai);
932
933 const labelList& visCoarseFaces = visibleFaceFaces[coarseFaceI];
934
935 forAll(visCoarseFaces, visCoarseFaceI)
936 {
937 //F2AI[coarseFaceI].setSize(visCoarseFaces.size());
938 F2LI[coarseFaceI].setSize(visCoarseFaces.size());
939 label compactJ = visCoarseFaces[visCoarseFaceI];
940 const List<point>& remoteFineSj = compactFineSf[compactJ];
941 const List<point>& remoteFineCj = compactFineCf[compactJ];
942
943 const List<point>& rPointsCj = compactPoints[compactJ];
944
945 const label toPatchId = compactPatchId[compactJ];
946
947 bool far(false);
948 // Relative distance
949 forAll(localFineSf, i)
950 {
951 const scalar dAi =
953 (
954 mag(localFineSf[i])/constant::mathematical::pi
955 );
956 const vector& dCi = localFineCf[i];
957
958 forAll(remoteFineSj, j)
959 {
960 const scalar dAj =
962 (
963 mag(remoteFineSj[j])/constant::mathematical::pi
964 );
965 const vector& dCj = remoteFineCj[j];
966
967 const scalar dist = mag(dCi - dCj)/((dAi + dAj)/2);
968
969 if (dist > distTol)
970 {
971 far = true;
972 }
973 }
974 }
975
976 if (far)
977 {
978 // 2AI method
979 scalar F2AIij = 0;
980
981 forAll(localFineSf, i)
982 {
983 const vector& dAi = localFineSf[i];
984 const vector& dCi = localFineCf[i];
985
986 forAll(remoteFineSj, j)
987 {
988 const vector& dAj = remoteFineSj[j];
989 const vector& dCj = remoteFineCj[j];
990
991 scalar dIntFij = calculateViewFactorFij2AI
992 (
993 dCi,
994 dCj,
995 dAi,
996 dAj
997 );
998
999 F2AIij += dIntFij;
1000 }
1001 }
1002 F2LI[coarseFaceI][visCoarseFaceI] = F2AIij/mag(Ai);
1003 }
1004 else
1005 {
1006 // 2LI method
1007 label nLocal = lPoints.size();
1008 label nRemote = rPointsCj.size();
1009
1010 // Using sub-divisions (quadrature)
1011 scalar oldEToeInt = 0;
1012 for (label gi=0; gi < maxQuadOrder; gi++)
1013 {
1014 scalar F2LIij = 0;
1015 for(label i=0; i<nLocal; i++)
1016 {
1017 vector si;
1018 vector ci;
1019
1020 vector sj;
1021 vector cj;
1022
1023 if (i == 0)
1024 {
1025 si = lPoints[i] - lPoints[nLocal-1];
1026 ci = (lPoints[i] + lPoints[nLocal-1])/2;
1027 }
1028 else
1029 {
1030 si = lPoints[i] - lPoints[i-1];
1031 ci = (lPoints[i] + lPoints[i-1])/2;
1032 }
1033
1034 for(label j=0; j<nRemote; j++)
1035 {
1036 if (j == 0)
1037 {
1038 sj = rPointsCj[j]-rPointsCj[nRemote-1];
1039 cj = (rPointsCj[j]+rPointsCj[nRemote-1])/2;
1040 }
1041 else
1042 {
1043 sj = rPointsCj[j] - rPointsCj[j-1];
1044 cj = (rPointsCj[j] + rPointsCj[j-1])/2;
1045 }
1046
1047
1048 scalar magSi = mag(si);
1049 scalar magSj = mag(sj);
1050 scalar cosij = (si & sj)/(magSi * magSj);
1051
1052 vector di = si/magSi;
1053 vector dj = sj/magSj;
1054
1055 label quadOrder = gi;
1056 const vector r(ci - cj);
1057 // Common edges use n = 0
1058 if (mag(r) < SMALL)
1059 {
1060 quadOrder = 0;
1061 }
1062
1063 scalar dIntFij =
1064 GaussQuad
1065 (
1066 GaussWeights[gi],
1067 GaussPoints[gi],
1068 magSi,
1069 magSj,
1070 di,
1071 dj,
1072 ci,
1073 cj,
1074 cosij,
1075 alpha,
1076 quadOrder
1077 );
1078
1079 F2LIij += dIntFij;
1080 }
1081 }
1082
1083 const scalar err
1084 (
1085 mag(F2LIij) > ROOTVSMALL
1086 ? (F2LIij-oldEToeInt)/F2LIij
1087 : 0
1088 );
1089
1090 if
1091 (
1092 (mag(err) < GaussQuadTol && gi > 0)
1093 || gi == maxQuadOrder-1
1094 )
1095 {
1096 F2LI[coarseFaceI][visCoarseFaceI] =
1097 F2LIij/mag(Ai)/4/constant::mathematical::pi;
1098 break;
1099 }
1100 else
1101 {
1102 oldEToeInt = F2LIij;
1103 }
1104 }
1105 }
1106
1107 sumViewFactorPatch[fromPatchId][toPatchId] +=
1108 F2LI[coarseFaceI][visCoarseFaceI]*mag(Ai);
1109 }
1110 }
1111 }
1112 else if (mesh.nSolutionD() == 2)
1113 {
1114 const boundBox& box = mesh.bounds();
1115 const Vector<label>& dirs = mesh.geometricD();
1116 vector emptyDir = Zero;
1117 forAll(dirs, i)
1118 {
1119 if (dirs[i] == -1)
1120 {
1121 emptyDir[i] = 1.0;
1122 }
1123 }
1124
1125 scalar wideBy2 = (box.span() & emptyDir)*2.0;
1126
1127 forAll(localCoarseSf, coarseFaceI)
1128 {
1129 const vector& Ai = localCoarseSf[coarseFaceI];
1130 const vector& Ci = localCoarseCf[coarseFaceI];
1131 vector Ain = Ai/mag(Ai);
1132 vector R1i = Ci + (mag(Ai)/wideBy2)*(Ain ^ emptyDir);
1133 vector R2i = Ci - (mag(Ai)/wideBy2)*(Ain ^ emptyDir) ;
1134
1135 const label fromPatchId = compactPatchId[coarseFaceI];
1136 patchArea[fromPatchId] += mag(Ai);
1137
1138 const labelList& visCoarseFaces = visibleFaceFaces[coarseFaceI];
1139 forAll(visCoarseFaces, visCoarseFaceI)
1140 {
1141 F2LI[coarseFaceI].setSize(visCoarseFaces.size());
1142 label compactJ = visCoarseFaces[visCoarseFaceI];
1143 const vector& Aj = compactCoarseSf[compactJ];
1144 const vector& Cj = compactCoarseCf[compactJ];
1145
1146 const label toPatchId = compactPatchId[compactJ];
1147
1148 vector Ajn = Aj/mag(Aj);
1149 vector R1j = Cj + (mag(Aj)/wideBy2)*(Ajn ^ emptyDir);
1150 vector R2j = Cj - (mag(Aj)/wideBy2)*(Ajn ^ emptyDir);
1151
1152 scalar d1 = mag(R1i - R2j);
1153 scalar d2 = mag(R2i - R1j);
1154 scalar s1 = mag(R1i - R1j);
1155 scalar s2 = mag(R2i - R2j);
1156
1157 scalar Fij = mag((d1 + d2) - (s1 + s2))/(4.0*mag(Ai)/wideBy2);
1158
1159 F2LI[coarseFaceI][visCoarseFaceI] = Fij;
1160 sumViewFactorPatch[fromPatchId][toPatchId] += Fij*mag(Ai);
1161 }
1162 }
1163 }
1164 else
1165 {
1167 << " View factors are not available in 1D "
1168 << exit(FatalError);
1169 }
1170
1171 // Write view factors matrix in listlist form
1172 F2LI.write();
1173
1174 reduce(sumViewFactorPatch, sumOp<scalarSquareMatrix>());
1175 reduce(patchArea, sumOp<scalarList>());
1176
1177 if (Pstream::master() && debug)
1178 {
1179 forAll(viewFactorsPatches, i)
1180 {
1181 label patchI = viewFactorsPatches[i];
1182 for (label j=i; j<viewFactorsPatches.size(); j++)
1183 {
1184 label patchJ = viewFactorsPatches[j];
1185
1186 Info << "F" << patchI << patchJ << ": "
1187 << sumViewFactorPatch[patchI][patchJ]/patchArea[patchI]
1188 << endl;
1189 }
1190 }
1191 }
1192
1193 if (writeViewFactors)
1194 {
1195 if (Pstream::master())
1196 {
1197 Info << "Writing view factor matrix..." << endl;
1198 }
1199
1200 volScalarField viewFactorField
1201 (
1202 IOobject
1203 (
1204 "viewFactorField",
1205 mesh.time().timeName(),
1206 mesh,
1209 ),
1210 mesh,
1212 );
1213
1214 label compactI = 0;
1215
1216 volScalarField::Boundary& vfbf = viewFactorField.boundaryFieldRef();
1217 forAll(viewFactorsPatches, i)
1218 {
1219 label patchID = viewFactorsPatches[i];
1220 const labelList& agglom = finalAgglom[patchID];
1221 if (agglom.size() > 0)
1222 {
1223 label nAgglom = max(agglom)+1;
1224 labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
1225 const labelList& coarsePatchFace =
1226 coarseMesh.patchFaceMap()[patchID];
1227
1228 forAll(coarseToFine, coarseI)
1229 {
1230 const scalar FiSum = sum(F2LI[compactI]);
1231
1232 const label coarseFaceID = coarsePatchFace[coarseI];
1233 const labelList& fineFaces = coarseToFine[coarseFaceID];
1234 forAll(fineFaces, fineId)
1235 {
1236 const label faceID = fineFaces[fineId];
1237 vfbf[patchID][faceID] = FiSum;
1238 }
1239 compactI++;
1240 }
1241 }
1242 }
1243 viewFactorField.write();
1244 }
1245
1246
1247 // Invert compactMap (from processor+localface to compact) to go
1248 // from compact to processor+localface (expressed as a globalIndex)
1249 // globalIndex globalCoarFaceNum(coarseMesh.nFaces());
1250 labelList compactToGlobal(map.constructSize());
1251
1252 // Local indices first (note: are not in compactMap)
1253 for (label i = 0; i < globalNumbering.localSize(); i++)
1254 {
1255 compactToGlobal[i] = globalNumbering.toGlobal(i);
1256 }
1257
1258
1259 forAll(compactMap, procI)
1260 {
1261 const Map<label>& localToCompactMap = compactMap[procI];
1262
1263 forAllConstIters(localToCompactMap, iter)
1264 {
1265 compactToGlobal[*iter] = globalNumbering.toGlobal
1266 (
1267 procI,
1268 iter.key()
1269 );
1270 }
1271 }
1272
1273
1274 labelListList globalFaceFaces(visibleFaceFaces.size());
1275
1276 // Create globalFaceFaces needed to insert view factors
1277 // in F to the global matrix Fmatrix
1278 forAll(globalFaceFaces, faceI)
1279 {
1280 globalFaceFaces[faceI] = renumber
1281 (
1282 compactToGlobal,
1283 visibleFaceFaces[faceI]
1284 );
1285 }
1286
1287 labelListIOList IOglobalFaceFaces
1288 (
1289 IOobject
1290 (
1291 "globalFaceFaces",
1292 mesh.facesInstance(),
1293 mesh,
1297 ),
1298 std::move(globalFaceFaces)
1299 );
1300 IOglobalFaceFaces.write();
1301
1302
1303 IOmapDistribute IOmapDist
1304 (
1305 IOobject
1306 (
1307 "mapDist",
1308 mesh.facesInstance(),
1309 mesh,
1312 ),
1313 std::move(map)
1314 );
1315
1316 IOmapDist.write();
1317
1318 Info<< "End\n" << endl;
1319 return 0;
1320}
1321
1322
1323// ************************************************************************* //
CGAL::Segment_3< K > Segment
CGAL::Point_3< K > Point
CGAL::Triangle_3< K > Triangle
CGAL::Exact_predicates_exact_constructions_kernel K
constexpr scalar pi(M_PI)
Required Classes.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
void setSize(const label n)
Dummy function, to make FixedList consistent with List.
Definition FixedList.H:437
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
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
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...
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
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual communication.
static void broadcast(Type &value, const int communicator=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all communicator ranks. Does nothing in non-paral...
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 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
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
vector span() const
The bounding box span (from minimum to maximum).
Definition boundBoxI.H:192
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
A class for handling file names.
Definition fileName.H:75
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
label toGlobal(const label proci, const label i) const
From local to global on proci.
label localSize(const label proci) const
Size of proci data.
A triFace with additional (region) index.
Definition labelledTri.H:56
Class containing processor-to-processor mapping information.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
label start() const noexcept
The start label of boundary faces in the polyMesh face list.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
fvMesh as subset of other mesh. Consists of one cell and all original boundary faces....
Triangulated surface description with patch information.
Definition triSurface.H:74
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:66
volScalarField & p
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
engineTime & runTime
Required Classes.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
return returnReduce(nRefine-oldNRefine, sumOp< label >())
Mathematical constants.
constexpr scalar pi(M_PI)
Different types of constants.
Namespace for handling debugging switches.
Definition debug.C:45
const wordList surface
Standard surface field types (scalar, vector, tensor, etc).
const std::string patch
OpenFOAM patch number as a std::string.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition meshTools.C:196
Namespace for OpenFOAM.
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
IOList< scalarList > scalarListIOList
IO for a List of scalarList.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimless
Dimensionless.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
IntListType renumber(const labelUList &oldToNew, const IntListType &input)
Renumber the values within a list.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
dimensionedScalar log(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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).
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...
Field< vector > vectorField
Specialisation of Field<T> for vector.
SquareMatrix< scalar > scalarSquareMatrix
vector point
Point is a vector.
Definition point.H:37
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
PrimitivePatch< List< face >, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points).
Definition bMesh.H:39
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
vectorField pointField
pointField is a vectorField.
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.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
volScalarField & alpha
labelList triSurfaceToAgglom(5 *nFineFaces)
std::vector< Triangle > triangles
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Dispatch tag: Construct 'one-sided' from local sizes, using gather but no broadcast.
Foam::surfaceFields.