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