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