meshRefinementBlock.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) 2018-2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "meshRefinement.H"
29 #include "fvMesh.H"
30 #include "Time.H"
31 #include "refinementSurfaces.H"
32 #include "removeCells.H"
33 #include "unitConversion.H"
34 #include "bitSet.H"
35 #include "volFields.H"
36 
37 // Leak path
38 #include "shortestPathSet.H"
39 #include "meshSearch.H"
40 #include "topoDistanceData.H"
41 #include "FaceCellWave.H"
42 #include "removeCells.H"
43 #include "regionSplit.H"
44 
45 #include "volFields.H"
46 #include "wallPoints.H"
47 #include "searchableSurfaces.H"
49 
50 #include "holeToFace.H"
51 #include "refinementParameters.H"
52 #include "indirectPrimitivePatch.H"
53 #include "OBJstream.H"
54 #include "PatchTools.H"
55 
56 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
57 
58 //Foam::label Foam::meshRefinement::markFakeGapRefinement
59 //(
60 // const scalar planarCos,
61 //
62 // const label nAllowRefine,
63 // const labelList& neiLevel,
64 // const pointField& neiCc,
65 //
66 // labelList& refineCell,
67 // label& nRefine
68 //) const
69 //{
70 // label oldNRefine = nRefine;
71 //
72 //
73 // // Collect candidate faces (i.e. intersecting any surface and
74 // // owner/neighbour not yet refined.
75 // const labelList testFaces(getRefineCandidateFaces(refineCell));
76 //
77 // // Collect segments
78 // pointField start(testFaces.size());
79 // pointField end(testFaces.size());
80 // labelList minLevel(testFaces.size());
81 //
82 // calcCellCellRays
83 // (
84 // neiCc,
85 // neiLevel,
86 // testFaces,
87 // start,
88 // end,
89 // minLevel
90 // );
91 //
92 //
93 // // Re-use the gap shooting methods. This needs:
94 // // - shell gapLevel : faked. Set to 0,labelMax
95 // // - surface gapLevel : faked by overwriting
96 //
97 //
98 // List<FixedList<label, 3>>& surfGapLevel = const_cast
99 // <
100 // List<FixedList<label, 3>>&
101 // >(surfaces_.extendedGapLevel());
102 //
103 // List<volumeType>& surfGapMode = const_cast
104 // <
105 // List<volumeType>&
106 // >(surfaces_.extendedGapMode());
107 //
108 // const List<FixedList<label, 3>> surfOldLevel(surfGapLevel);
109 // const List<volumeType> surfOldMode(surfGapMode);
110 //
111 // // Set the extended gap levels
112 // forAll(surfaces_.gapLevel(), regioni)
113 // {
114 // surfGapLevel[regioni] = FixedList<label, 3>
115 // ({
116 // 3,
117 // -1,
118 // surfaces_.gapLevel()[regioni]+1
119 // });
120 // }
121 // surfGapMode = volumeType::MIXED;
122 //
123 //Pout<< "gapLevel was:" << surfOldLevel << endl;
124 //Pout<< "gapLevel now:" << surfGapLevel << endl;
125 //Pout<< "gapMode was:" << surfOldMode << endl;
126 //Pout<< "gapMode now:" << surfGapMode << endl;
127 //Pout<< "nRefine was:" << oldNRefine << endl;
128 //
129 //
130 //
131 // List<List<FixedList<label, 3>>>& shellGapLevel = const_cast
132 // <
133 // List<List<FixedList<label, 3>>>&
134 // >(shells_.extendedGapLevel());
135 //
136 // List<List<volumeType>>& shellGapMode = const_cast
137 // <
138 // List<List<volumeType>>&
139 // >(shells_.extendedGapMode());
140 //
141 // const List<List<FixedList<label, 3>>> shellOldLevel(shellGapLevel);
142 // const List<List<volumeType>> shellOldMode(shellGapMode);
143 //
144 // // Set the extended gap levels
145 // forAll(shellGapLevel, shelli)
146 // {
147 // shellGapLevel[shelli] = FixedList<label, 3>({3, -1, labelMax});
148 // shellGapMode[shelli] = volumeType::MIXED;
149 // }
150 //Pout<< "shellLevel was:" << shellOldLevel << endl;
151 //Pout<< "shellLevel now:" << shellGapLevel << endl;
152 //
153 // const label nAdditionalRefined = markSurfaceGapRefinement
154 // (
155 // planarCos,
156 //
157 // nAllowRefine,
158 // neiLevel,
159 // neiCc,
160 //
161 // refineCell,
162 // nRefine
163 // );
164 //
165 //Pout<< "nRefine now:" << nRefine << endl;
166 //
167 // // Restore
168 // surfGapLevel = surfOldLevel;
169 // surfGapMode = surfOldMode;
170 // shellGapLevel = shellOldLevel;
171 // shellGapMode = shellOldMode;
172 //
173 // return nAdditionalRefined;
174 //}
175 
176 
178 (
179  const labelList& cellLevel,
180  const labelList& neiLevel,
181  const labelList& refineCell,
182  bitSet& isOutsideFace
183 ) const
184 {
185  // Get faces:
186  // - on outside of cell set
187  // - inbetween same cell level (i.e. quads)
188 
189  isOutsideFace.setSize(mesh_.nFaces());
190  isOutsideFace = Zero;
191 
192  forAll(mesh_.faceNeighbour(), facei)
193  {
194  label own = mesh_.faceOwner()[facei];
195  label nei = mesh_.faceNeighbour()[facei];
196  if
197  (
198  (cellLevel[own] == cellLevel[nei])
199  && (
200  (refineCell[own] != -1)
201  != (refineCell[nei] != -1)
202  )
203  )
204  {
205  isOutsideFace.set(facei);
206  }
207  }
208  {
209 
210  const label nBnd = mesh_.nBoundaryFaces();
211 
212  labelList neiRefineCell(nBnd);
213  syncTools::swapBoundaryCellList(mesh_, refineCell, neiRefineCell);
214  for (label bFacei = 0; bFacei < nBnd; ++bFacei)
215  {
216  label facei = mesh_.nInternalFaces()+bFacei;
217  label own = mesh_.faceOwner()[facei];
218 
219  if
220  (
221  (cellLevel[own] == neiLevel[bFacei])
222  && (
223  (refineCell[own] != -1)
224  != (neiRefineCell[bFacei] != -1)
225  )
226  )
227  {
228  isOutsideFace.set(facei);
229  }
230  }
231  }
232 }
233 
234 
236 (
237  const bitSet& isOutsideFace,
238  const label celli
239 ) const
240 {
241  const cell& cFaces = mesh_.cells()[celli];
242  const vectorField& faceAreas = mesh_.faceAreas();
243 
244  Vector<bool> haveDirs(vector::uniform(false));
245 
246  forAll(cFaces, cFacei)
247  {
248  const label facei = cFaces[cFacei];
249 
250  if (isOutsideFace[facei])
251  {
252  const vector& n = faceAreas[facei];
253  scalar magSqrN = magSqr(n);
254 
255  if (magSqrN > ROOTVSMALL)
256  {
257  for
258  (
259  direction dir = 0;
261  dir++
262  )
263  {
264  if (Foam::sqr(n[dir]) > 0.99*magSqrN)
265  {
266  haveDirs[dir] = true;
267  break;
268  }
269  }
270  }
271  }
272  }
273 
274  label nDirs = 0;
275  forAll(haveDirs, dir)
276  {
277  if (haveDirs[dir])
278  {
279  nDirs++;
280  }
281  }
282  return nDirs;
283 }
284 
285 
287 (
288  const labelList& neiLevel,
289  const bitSet& isOutsideFace,
290  labelList& refineCell,
291  label& nRefine
292 ) const
293 {
294  // Get cells with three or more outside faces
295  const cellList& cells = mesh_.cells();
296  forAll(cells, celli)
297  {
298  if (refineCell[celli] == -1)
299  {
300  if (countFaceDirs(isOutsideFace, celli) == 3)
301  {
302  // Mark cell with any value
303  refineCell[celli] = 0;
304  nRefine++;
305  }
306  }
307  }
308 }
309 
310 
311 //void Foam::meshRefinement::markMultiRegionCell
312 //(
313 // const label celli,
314 // const FixedList<label, 3>& surface,
315 //
316 // Map<FixedList<label, 3>>& cellToRegions,
317 // bitSet& isMultiRegion
318 //) const
319 //{
320 // if (!isMultiRegion[celli])
321 // {
322 // Map<FixedList<label, 3>>::iterator fnd = cellToRegions.find(celli);
323 // if (!fnd.good())
324 // {
325 // cellToRegions.insert(celli, surface);
326 // }
327 // else if (fnd() != surface)
328 // {
329 // // Found multiple intersections on cell
330 // isMultiRegion.set(celli);
331 // }
332 // }
333 //}
334 
335 
336 //void Foam::meshRefinement::detectMultiRegionCells
337 //(
338 // const labelListList& faceZones,
339 // const labelList& testFaces,
340 //
341 // const labelList& surface1,
342 // const List<pointIndexHit>& hit1,
343 // const labelList& region1,
344 //
345 // const labelList& surface2,
346 // const List<pointIndexHit>& hit2,
347 // const labelList& region2,
348 //
349 // bitSet& isMultiRegion
350 //) const
351 //{
352 // isMultiRegion.clear();
353 // isMultiRegion.setSize(mesh_.nCells());
354 //
355 // Map<FixedList<label, 3>> cellToRegions(testFaces.size());
356 //
357 // forAll(testFaces, i)
358 // {
359 // const pointIndexHit& info1 = hit1[i];
360 // if (info1.hit())
361 // {
362 // const label facei = testFaces[i];
363 // const labelList& fz1 = faceZones[surface1[i]];
364 // const FixedList<label, 3> surfaceInfo1
365 // ({
366 // surface1[i],
367 // region1[i],
368 // (fz1.size() ? fz1[info1.index()] : region1[i])
369 // });
370 //
371 // markMultiRegionCell
372 // (
373 // mesh_.faceOwner()[facei],
374 // surfaceInfo1,
375 // cellToRegions,
376 // isMultiRegion
377 // );
378 // if (mesh_.isInternalFace(facei))
379 // {
380 // markMultiRegionCell
381 // (
382 // mesh_.faceNeighbour()[facei],
383 // surfaceInfo1,
384 // cellToRegions,
385 // isMultiRegion
386 // );
387 // }
388 //
389 // const pointIndexHit& info2 = hit2[i];
390 //
391 // if (info2.hit() && info1 != info2)
392 // {
393 // const labelList& fz2 = faceZones[surface2[i]];
394 // const FixedList<label, 3> surfaceInfo2
395 // ({
396 // surface2[i],
397 // region2[i],
398 // (fz2.size() ? fz2[info2.index()] : region2[i])
399 // });
400 //
401 // markMultiRegionCell
402 // (
403 // mesh_.faceOwner()[facei],
404 // surfaceInfo2,
405 // cellToRegions,
406 // isMultiRegion
407 // );
408 // if (mesh_.isInternalFace(facei))
409 // {
410 // markMultiRegionCell
411 // (
412 // mesh_.faceNeighbour()[facei],
413 // surfaceInfo2,
414 // cellToRegions,
415 // isMultiRegion
416 // );
417 // }
418 // }
419 // }
420 // }
421 //
422 //
423 // if (debug&meshRefinement::MESH)
424 // {
425 // volScalarField multiCell
426 // (
427 // IOobject
428 // (
429 // "multiCell",
430 // mesh_.time().timeName(),
431 // mesh_,
432 // IOobject::NO_READ,
433 // IOobject::NO_WRITE,
434 // IOobject::NO_REGISTER
435 // ),
436 // mesh_,
437 // dimensionedScalar
438 // (
439 // "zero",
440 // dimensionSet(0, 1, 0, 0, 0),
441 // 0.0
442 // )
443 // );
444 // forAll(isMultiRegion, celli)
445 // {
446 // if (isMultiRegion[celli])
447 // {
448 // multiCell[celli] = 1.0;
449 // }
450 // }
451 //
452 // Info<< "Writing all multi cells to " << multiCell.name() << endl;
453 // multiCell.write();
454 // }
455 //}
456 
457 
458 Foam::label Foam::meshRefinement::markProximityRefinementWave
459 (
460  const scalar planarCos,
461  const labelList& blockedSurfaces,
462  const label nAllowRefine,
463  const labelList& neiLevel,
464  const pointField& neiCc,
465 
466  labelList& refineCell,
467  label& nRefine
468 ) const
469 {
470  labelListList faceZones(surfaces_.surfaces().size());
471  {
472  // If triSurface do additional zoning based on connectivity
473  for (const label surfi : blockedSurfaces)
474  {
475  const label geomi = surfaces_.surfaces()[surfi];
476  const searchableSurface& s = surfaces_.geometry()[geomi];
477  if (isA<triSurfaceMesh>(s) && !isA<distributedTriSurfaceMesh>(s))
478  {
479  const triSurfaceMesh& surf = refCast<const triSurfaceMesh>(s);
480  const labelListList& edFaces = surf.edgeFaces();
481  boolList isOpenEdge(edFaces.size(), false);
482  forAll(edFaces, i)
483  {
484  if (edFaces[i].size() == 1)
485  {
486  isOpenEdge[i] = true;
487  }
488  }
489 
490  labelList faceZone;
491  const label nZones = surf.markZones(isOpenEdge, faceZone);
492  if (nZones > 1)
493  {
494  faceZones[surfi].transfer(faceZone);
495  }
496  }
497  }
498  }
499 
500 
501  // Re-work the blockLevel of the blockedSurfaces into a length scale
502  // and a number of cells to cross
503  List<scalarList> regionToBlockSize(surfaces_.surfaces().size());
504  for (const label surfi : blockedSurfaces)
505  {
506  const label geomi = surfaces_.surfaces()[surfi];
507  const searchableSurface& s = surfaces_.geometry()[geomi];
508  const label nRegions = s.regions().size();
509  regionToBlockSize[surfi].setSize(nRegions);
510  for (label regioni = 0; regioni < nRegions; regioni++)
511  {
512  const label globalRegioni = surfaces_.globalRegion(surfi, regioni);
513  const label bLevel = surfaces_.blockLevel()[globalRegioni];
514  regionToBlockSize[surfi][regioni] =
515  meshCutter_.level0EdgeLength()/pow(2.0, bLevel);
516 
517  //const label mLevel = surfaces_.maxLevel()[globalRegioni];
521  //if (isA<triSurfaceMesh>(s) && !isA<distributedTriSurfaceMesh>(s))
522  //{
523  // const triSurfaceMesh& surf = refCast<const triSurfaceMesh>(s);
524  //}
525 
526  //nIters = max(nIters, (2<<(mLevel-bLevel)));
527  }
528  }
529 
530  // Clever limiting of the number of iterations (= max cells in the channel)
531  // since it has too many problematic issues, e.g. with volume refinement
532  // and the real check uses the proper distance anyway just disable.
533  const label nIters = mesh_.globalData().nTotalCells();
534 
535 
536  // Collect candidate faces (i.e. intersecting any surface and
537  // owner/neighbour not yet refined)
538  const labelList testFaces(getRefineCandidateFaces(refineCell));
539 
540  // Collect segments
541  pointField start(testFaces.size());
542  pointField end(testFaces.size());
543  labelList minLevel(testFaces.size());
544 
545  calcCellCellRays
546  (
547  neiCc,
548  neiLevel,
549  testFaces,
550  start,
551  end,
552  minLevel
553  );
554  // TBD. correct nIters for actual cellLevel (since e.g. volume refinement
555  // might add to cell level). Unfortunately we only have minLevel,
556  // not maxLevel. Problem: what if volume refinement only in middle of
557  // channel? Say channel 1m wide with a 0.1m sphere of refinement
558  // Workaround: have dummy surface with e.g. maxLevel 100 to
559  // force nIters to be high enough.
560 
561 
562  // Test for all intersections (with surfaces of higher gap level than
563  // minLevel) and cache per cell the max surface level and the local normal
564  // on that surface.
565 
566  labelList surface1;
567  List<pointIndexHit> hit1;
568  labelList region1;
569  vectorField normal1;
570 
571  labelList surface2;
572  List<pointIndexHit> hit2;
573  labelList region2;
574  vectorField normal2;
575 
576  surfaces_.findNearestIntersection
577  (
578  blockedSurfaces,
579  start,
580  end,
581 
582  surface1,
583  hit1,
584  region1, // local region
585  normal1,
586 
587  surface2,
588  hit2,
589  region2, // local region
590  normal2
591  );
592 
593 
594  // Detect cells that are using multiple surface regions
595  //bitSet isMultiRegion;
596  //detectMultiRegionCells
597  //(
598  // faceZones,
599  // testFaces,
600  //
601  // surface1,
602  // hit1,
603  // region1,
604  //
605  // surface2,
606  // hit2,
607  // region2,
608  //
609  // isMultiRegion
610  //);
611 
612 
613  label n = 0;
614  forAll(testFaces, i)
615  {
616  if (hit1[i].hit())
617  {
618  n++;
619  }
620  }
621 
622  List<wallPoints> faceDist(n);
623  labelList changedFaces(n);
624  n = 0;
625 
626  DynamicList<point> originLocation(2);
627  DynamicList<scalar> originDistSqr(2);
628  DynamicList<FixedList<label, 3>> originSurface(2);
629  //DynamicList<point> originNormal(2);
630 
631 
632  //- To avoid walking through surfaces we mark all faces that have been
633  // intersected. We can either mark only those faces intersecting
634  // blockedSurfaces (i.e. with a 'blockLevel') or mark faces intersecting
635  // any (refinement) surface (this includes e.g. faceZones). This is
636  // much easier since that information is already cached
637  // (meshRefinement::intersectedFaces())
638 
639  //bitSet isBlockedFace(mesh_.nFaces());
640  forAll(testFaces, i)
641  {
642  if (hit1[i].hit())
643  {
644  const label facei = testFaces[i];
645  //isBlockedFace.set(facei);
646  const point& fc = mesh_.faceCentres()[facei];
647  const labelList& fz1 = faceZones[surface1[i]];
648 
649  originLocation.clear();
650  originDistSqr.clear();
651  //originNormal.clear();
652  originSurface.clear();
653 
654  originLocation.append(hit1[i].point());
655  originDistSqr.append(hit1[i].point().distSqr(fc));
656  //originNormal.append(normal1[i]);
657  originSurface.append
658  (
659  FixedList<label, 3>
660  ({
661  surface1[i],
662  region1[i],
663  (fz1.size() ? fz1[hit1[i].index()] : region1[i])
664  })
665  );
666 
667  if (hit2[i].hit() && hit1[i] != hit2[i])
668  {
669  const labelList& fz2 = faceZones[surface2[i]];
670  originLocation.append(hit2[i].point());
671  originDistSqr.append(hit2[i].point().distSqr(fc));
672  //originNormal.append(normal2[i]);
673  originSurface.append
674  (
675  FixedList<label, 3>
676  ({
677  surface2[i],
678  region2[i],
679  (fz2.size() ? fz2[hit2[i].index()] : region2[i])
680  })
681  );
682  }
683 
684  // Collect all seed data. Currently walking does not look at
685  // surface direction - if so pass in surface normal as well
686  faceDist[n] = wallPoints
687  (
688  originLocation, // origin
689  originDistSqr, // distance to origin
690  originSurface // surface+region+zone
691  //originNormal // normal at origin
692  );
693  changedFaces[n] = facei;
694  n++;
695  }
696  }
697 
698 
699  // Clear intersection info
700  surface1.clear();
701  hit1.clear();
702  region1.clear();
703  normal1.clear();
704  surface2.clear();
705  hit2.clear();
706  region2.clear();
707  normal2.clear();
708 
709 
710  List<wallPoints> allFaceInfo(mesh_.nFaces());
711  List<wallPoints> allCellInfo(mesh_.nCells());
712 
713  // Any refinement surface (even a faceZone) should stop the gap walking.
714  // This is exactly the information which is cached in the surfaceIndex_
715  // field.
716  const bitSet isBlockedFace(intersectedFaces());
717 
718  wallPoints::trackData td(isBlockedFace, regionToBlockSize);
719  FaceCellWave<wallPoints, wallPoints::trackData> wallDistCalc
720  (
721  mesh_,
722  changedFaces,
723  faceDist,
724  allFaceInfo,
725  allCellInfo,
726  0, // max iterations
727  td
728  );
729  wallDistCalc.iterate(nIters);
730 
731 
733  {
734  // Dump current nearest opposite surfaces
736  (
737  IOobject
738  (
739  "gapSize",
740  mesh_.time().timeName(),
741  mesh_,
745  ),
746  mesh_,
748  (
749  "zero",
750  dimLength, //dimensionSet(0, 1, 0, 0, 0),
751  -1
752  )
753  );
754 
755  forAll(allCellInfo, celli)
756  {
757  if (allCellInfo[celli].valid(wallDistCalc.data()))
758  {
759  const point& cc = mesh_.cellCentres()[celli];
760  // Nearest surface points
761  const List<point>& origin = allCellInfo[celli].origin();
762 
763  // Find 'opposite' pair with minimum distance
764  for (label i = 0; i < origin.size(); i++)
765  {
766  for (label j = i + 1; j < origin.size(); j++)
767  {
768  if (((cc-origin[i]) & (cc-origin[j])) < 0)
769  {
770  const scalar d(mag(origin[i]-origin[j]));
771  if (distance[celli] < 0)
772  {
773  distance[celli] = d;
774  }
775  else
776  {
777  distance[celli] = min(distance[celli], d);
778  }
779  }
780  }
781  }
782  }
783  }
784  distance.correctBoundaryConditions();
785 
786  Info<< "Writing measured gap distance to "
787  << distance.name() << endl;
788  distance.write();
789  }
790 
791 
792 
793  // Detect tight gaps:
794  // - cell is inbetween the two surfaces
795  // - two surfaces are planarish
796  // - two surfaces are not too far apart
797  // (number of walking iterations is a too-coarse measure)
798 
799  scalarField smallGapDistance(mesh_.nCells(), 0.0);
800  label nMulti = 0;
801  label nSmallGap = 0;
802 
803  //OBJstream str(mesh_.time().timePath()/"multiRegion.obj");
804 
805 
806  forAll(allCellInfo, celli)
807  {
808  if (allCellInfo[celli].valid(wallDistCalc.data()))
809  {
810  const point& cc = mesh_.cellCentres()[celli];
811 
812  const List<point>& origin = allCellInfo[celli].origin();
813  const List<FixedList<label, 3>>& surface =
814  allCellInfo[celli].surface();
815 
816  // Find pair with minimum distance
817  for (label i = 0; i < origin.size(); i++)
818  {
819  for (label j = i + 1; j < origin.size(); j++)
820  {
821  //if (isMultiRegion[celli])
822  //{
823  // // The intersection locations are too inaccurate
824  // // (since not proper nearest, just a cell-cell ray
825  // // intersection) so include always
826  //
827  // smallGapDistance[celli] =
828  // max(smallGapDistance[celli], maxDist);
829  //
830  //
831  // str.writeLine(cc, origin[i]);
832  // str.writeLine(cc, origin[j]);
833  //
834  // nMulti++;
835  //}
836  //else
837  if (((cc-origin[i]) & (cc-origin[j])) < 0)
838  {
839  const label surfi = surface[i][0];
840  const label regioni = surface[i][1];
841 
842  const label surfj = surface[j][0];
843  const label regionj = surface[j][1];
844 
845  const scalar maxSize = max
846  (
847  regionToBlockSize[surfi][regioni],
848  regionToBlockSize[surfj][regionj]
849  );
850 
851  if
852  (
853  magSqr(origin[i]-origin[j])
854  < Foam::sqr(2*maxSize)
855  )
856  {
857  const scalar maxDist
858  (
859  max
860  (
861  mag(cc-origin[i]),
862  mag(cc-origin[j])
863  )
864  );
865 
866  smallGapDistance[celli] =
867  max(smallGapDistance[celli], maxDist);
868  nSmallGap++;
869  }
870  }
871  }
872  }
873  }
874  }
875 
876 
877  if (debug)
878  {
879  Info<< "Marked for blocking due to intersecting multiple surfaces : "
880  << returnReduce(nMulti, sumOp<label>()) << " cells." << endl;
881  Info<< "Marked for blocking due to close opposite surfaces : "
882  << returnReduce(nSmallGap, sumOp<label>()) << " cells." << endl;
883  }
884 
886  {
888  (
889  IOobject
890  (
891  "smallGapDistance",
892  mesh_.time().timeName(),
893  mesh_,
897  ),
898  mesh_,
900  );
901  distance.field() = smallGapDistance;
902  distance.correctBoundaryConditions();
903 
904  Info<< "Writing all small-gap cells to "
905  << distance.name() << endl;
906  distance.write();
907  }
908 
909 
910  // Mark refinement
911  const label oldNRefine = nRefine;
912  forAll(smallGapDistance, celli)
913  {
914  if (smallGapDistance[celli] > SMALL)
915  {
916  if
917  (
918  !markForRefine
919  (
920  0, // mark level
921  nAllowRefine,
922  refineCell[celli],
923  nRefine
924  )
925  )
926  {
927  if (debug)
928  {
929  Pout<< "Stopped refining since reaching my cell"
930  << " limit of " << mesh_.nCells()+7*nRefine
931  << endl;
932  }
933  break;
934  }
935  }
936  }
937 
938  if
939  (
940  returnReduce(nRefine, sumOp<label>())
941  > returnReduce(nAllowRefine, sumOp<label>())
942  )
943  {
944  Info<< "Reached refinement limit." << endl;
945  }
946 
947  return returnReduce(nRefine-oldNRefine, sumOp<label>());
948 }
949 
950 
952 (
953  const scalar planarAngle,
954  const labelList& minSurfaceLevel,
955  const labelList& globalToMasterPatch,
956  const label growIter
957 )
958 {
959  // Swap neighbouring cell centres and cell level
960  labelList neiLevel(mesh_.nBoundaryFaces());
961  pointField neiCc(mesh_.nBoundaryFaces());
962  calcNeighbourData(neiLevel, neiCc);
963 
964  labelList refineCell(mesh_.nCells(), -1);
965  label nRefine = 0;
966  //markProximityRefinement
967  //(
968  // Foam::cos(degToRad(planarAngle)),
969  //
970  // minSurfaceLevel, // surface min level
971  // labelList(minSurfaceLevel.size(), labelMax), // surfaceGapLevel
972  //
973  // labelMax/Pstream::nProcs(), //nAllowRefine,
974  // neiLevel,
975  // neiCc,
976  //
977  // refineCell,
978  // nRefine
979  //);
980 
981 
982  // Determine minimum blockLevel per surface
983  Map<label> surfToBlockLevel;
984 
985  forAll(surfaces_.surfaces(), surfi)
986  {
987  const label geomi = surfaces_.surfaces()[surfi];
988  const searchableSurface& s = surfaces_.geometry()[geomi];
989  const label nRegions = s.regions().size();
990 
991  label minBlockLevel = labelMax;
992  for (label regioni = 0; regioni < nRegions; regioni++)
993  {
994  const label globalRegioni = surfaces_.globalRegion(surfi, regioni);
995  minBlockLevel = min
996  (
997  minBlockLevel,
998  surfaces_.blockLevel()[globalRegioni]
999  );
1000  }
1001 
1002  if (minBlockLevel < labelMax)
1003  {
1004  surfToBlockLevel.insert(surfi, minBlockLevel);
1005  }
1006  }
1007 
1008 
1009  markProximityRefinementWave
1010  (
1011  Foam::cos(degToRad(planarAngle)),
1012  surfToBlockLevel.sortedToc(),
1013 
1014  labelMax/Pstream::nProcs(), //nAllowRefine,
1015  neiLevel,
1016  neiCc,
1017 
1018  refineCell,
1019  nRefine
1020  );
1021 
1022 
1024  //markFakeGapRefinement
1025  //(
1026  // Foam::cos(degToRad(planarAngle)),
1027  //
1028  // labelMax/Pstream::nProcs(), //nAllowRefine,
1029  // neiLevel,
1030  // neiCc,
1031  //
1032  // refineCell,
1033  // nRefine
1034  //);
1035 
1036 
1037  Info<< "Marked for blocking due to close opposite surfaces : "
1038  << returnReduce(nRefine, sumOp<label>()) << " cells." << endl;
1039 
1040  // Remove outliers, i.e. cells with all points exposed
1041  if (growIter)
1042  {
1043  labelList oldRefineCell(refineCell);
1044 
1045  // Pass1: extend the set to fill in gaps
1046  bitSet isOutsideFace;
1047  for (label iter = 0; iter < growIter; iter++)
1048  {
1049  // Get outside faces
1050  markOutsideFaces
1051  (
1052  meshCutter_.cellLevel(),
1053  neiLevel,
1054  refineCell,
1055  isOutsideFace
1056  );
1057  // Extend with cells with three outside faces
1058  growSet(neiLevel, isOutsideFace, refineCell, nRefine);
1059  }
1060 
1061 
1062  // Pass2: erode back to original set if pass1 didn't help
1063  for (label iter = 0; iter < growIter; iter++)
1064  {
1065  // Get outside faces. Ignore cell level.
1066  markOutsideFaces
1067  (
1068  labelList(mesh_.nCells(), 0),
1069  labelList(neiLevel.size(), 0),
1070  refineCell,
1071  isOutsideFace
1072  );
1073 
1074  // Unmark cells with three or more outside faces
1075  for (label celli = 0; celli < mesh_.nCells(); celli++)
1076  {
1077  if (refineCell[celli] != -1 && oldRefineCell[celli] == -1)
1078  {
1079  if (countFaceDirs(isOutsideFace, celli) >= 3)
1080  {
1081  refineCell[celli] = -1;
1082  --nRefine;
1083  }
1084  }
1085  }
1086  }
1087 
1088  Info<< "Marked for blocking after filtering : "
1089  << returnReduce(nRefine, sumOp<label>()) << " cells." << endl;
1090  }
1091 
1092 
1093  // Determine patch for every mesh face
1094  const PtrList<surfaceZonesInfo>& surfZones = surfaces_.surfZones();
1095  labelList unnamedSurfaces(surfaceZonesInfo::getUnnamedSurfaces(surfZones));
1096  const label defaultRegion(surfaces_.globalRegion(unnamedSurfaces[0], 0));
1097 
1098  const labelList nearestRegion
1099  (
1100  nearestIntersection
1101  (
1102  unnamedSurfaces,
1103  defaultRegion
1104  )
1105  );
1106 
1107  // Pack
1108  labelList cellsToRemove(nRefine);
1109  nRefine = 0;
1110 
1111  forAll(refineCell, cellI)
1112  {
1113  if (refineCell[cellI] != -1)
1114  {
1115  cellsToRemove[nRefine++] = cellI;
1116  }
1117  }
1118 
1119  // Remove cells
1120  removeCells cellRemover(mesh_);
1121  labelList exposedFaces(cellRemover.getExposedFaces(cellsToRemove));
1122 
1123  labelList exposedPatches(exposedFaces.size());
1124  forAll(exposedFaces, i)
1125  {
1126  label facei = exposedFaces[i];
1127  exposedPatches[i] = globalToMasterPatch[nearestRegion[facei]];
1128  }
1129 
1130  return doRemoveCells
1131  (
1132  cellsToRemove,
1133  exposedFaces,
1134  exposedPatches,
1135  cellRemover
1136  );
1138 
1139 
1141 (
1142  const labelList& selectedSurfaces,
1143  boolList& isBlockedFace
1144 ) const
1145 {
1146  // Like meshRefinement::selectSeparatedCoupledFaces. tbd: convert to bitSet
1147 
1148  // Check that no connection between inside and outside points
1149  isBlockedFace.setSize(mesh_.nFaces(), false);
1150 
1151  // Block off separated couples.
1152  selectSeparatedCoupledFaces(isBlockedFace);
1153 
1154  // Block off intersections with selected surfaces
1155 
1156  // Mark per face (for efficiency)
1157  boolList isSelectedSurf(surfaces_.surfaces().size(), false);
1158  UIndirectList<bool>(isSelectedSurf, selectedSurfaces) = true;
1159 
1160  forAll(surfaceIndex_, facei)
1161  {
1162  const label surfi = surfaceIndex_[facei];
1163  if (surfi != -1 && isSelectedSurf[surfi])
1164  {
1165  isBlockedFace[facei] = true;
1166  }
1167  }
1168 }
1169 
1170 
1171 //Foam::labelList Foam::meshRefinement::detectLeakCells
1172 //(
1173 // const boolList& isBlockedFace,
1174 // const labelList& leakFaces,
1175 // const labelList& seedCells
1176 //) const
1177 //{
1178 // int dummyTrackData = 0;
1179 // List<topoDistanceData<label>> allFaceInfo(mesh_.nFaces());
1180 // List<topoDistanceData<label>> allCellInfo(mesh_.nCells());
1181 //
1182 // // Block faces
1183 // forAll(isBlockedFace, facei)
1184 // {
1185 // if (isBlockedFace[facei])
1186 // {
1187 // allFaceInfo[facei] = topoDistanceData<label>(labelMax, 123);
1188 // }
1189 // }
1190 // for (const label facei : leakFaces)
1191 // {
1192 // allFaceInfo[facei] = topoDistanceData<label>(labelMax, 123);
1193 // }
1194 //
1195 // // Walk from inside cell
1196 // DynamicList<topoDistanceData<label>> faceDist;
1197 // DynamicList<label> cFaces1;
1198 // for (const label celli : seedCells)
1199 // {
1200 // if (celli != -1)
1201 // {
1202 // const labelList& cFaces = mesh_.cells()[celli];
1203 // faceDist.reserve(cFaces.size());
1204 // cFaces1.reserve(cFaces.size());
1205 //
1206 // for (const label facei : cFaces)
1207 // {
1208 // if (!allFaceInfo[facei].valid(dummyTrackData))
1209 // {
1210 // cFaces1.append(facei);
1211 // faceDist.append(topoDistanceData<label>(0, 123));
1212 // }
1213 // }
1214 // }
1215 // }
1216 //
1217 // // Walk through face-cell wave till all cells are reached
1218 // FaceCellWave<topoDistanceData<label>> wallDistCalc
1219 // (
1220 // mesh_,
1221 // cFaces1,
1222 // faceDist,
1223 // allFaceInfo,
1224 // allCellInfo,
1225 // mesh_.globalData().nTotalCells()+1 // max iterations
1226 // );
1227 //
1228 // label nRemove = 0;
1229 // for (const label facei : leakFaces)
1230 // {
1231 // const label own = mesh_.faceOwner()[facei];
1232 // if (!allCellInfo[own].valid(dummyTrackData))
1233 // {
1234 // nRemove++;
1235 // }
1236 // if (mesh_.isInternalFace(facei))
1237 // {
1238 // const label nei = mesh_.faceNeighbour()[facei];
1239 // if (!allCellInfo[nei].valid(dummyTrackData))
1240 // {
1241 // nRemove++;
1242 // }
1243 // }
1244 // }
1245 //
1246 // labelList cellsToRemove(nRemove);
1247 // nRemove = 0;
1248 // for (const label facei : leakFaces)
1249 // {
1250 // const label own = mesh_.faceOwner()[facei];
1251 // if (!allCellInfo[own].valid(dummyTrackData))
1252 // {
1253 // cellsToRemove[nRemove++] = own;
1254 // }
1255 // if (mesh_.isInternalFace(facei))
1256 // {
1257 // const label nei = mesh_.faceNeighbour()[facei];
1258 // if (!allCellInfo[nei].valid(dummyTrackData))
1259 // {
1260 // cellsToRemove[nRemove++] = nei;
1261 // }
1262 // }
1263 // }
1264 //
1265 // if (debug)
1266 // {
1267 // volScalarField fld
1268 // (
1269 // IOobject
1270 // (
1271 // "cellsToKeep",
1272 // mesh_.time().timeName(),
1273 // mesh_,
1274 // IOobject::NO_READ,
1275 // IOobject::NO_WRITE
1276 // ),
1277 // mesh_,
1278 // dimensionedScalar(dimless, Zero)
1279 // );
1280 // forAll(allCellInfo, celli)
1281 // {
1282 // if (allCellInfo[celli].valid(dummyTrackData))
1283 // {
1284 // fld[celli] = allCellInfo[celli].distance();
1285 // }
1286 // }
1287 // forAll(fld.boundaryField(), patchi)
1288 // {
1289 // const polyPatch& pp = mesh_.boundaryMesh()[patchi];
1290 // SubList<topoDistanceData<label>> p(pp.patchSlice(allFaceInfo));
1291 // scalarField pfld
1292 // (
1293 // fld.boundaryField()[patchi].size(),
1294 // Zero
1295 // );
1296 // forAll(pfld, i)
1297 // {
1298 // if (p[i].valid(dummyTrackData))
1299 // {
1300 // pfld[i] = 1.0*p[i].distance();
1301 // }
1302 // }
1303 // fld.boundaryFieldRef()[patchi] == pfld;
1304 // }
1305 // //Note: do not swap cell values so do not do
1306 // //fld.correctBoundaryConditions();
1307 // Pout<< "Writing distance field for initial cells "
1308 // << seedCells << " to " << fld.objectPath() << endl;
1309 // fld.write();
1310 // }
1311 //
1312 // return cellsToRemove;
1313 //}
1314 //
1315 //
1316 //Foam::autoPtr<Foam::mapPolyMesh> Foam::meshRefinement::removeLeakCells
1317 //(
1318 // const labelList& globalToMasterPatch,
1319 // const labelList& globalToSlavePatch,
1320 // const pointField& locationsInMesh,
1321 // const wordList& zonesInMesh,
1322 // const pointField& locationsOutsideMesh,
1323 // const labelList& selectedSurfaces
1324 //)
1325 //{
1326 // boolList isBlockedFace;
1327 // selectIntersectedFaces(selectedSurfaces, isBlockedFace);
1328 //
1329 // // Determine cell regions
1330 // const regionSplit cellRegion(mesh_, isBlockedFace);
1331 //
1332 // // Detect locationsInMesh regions
1333 // labelList insideCells(locationsInMesh.size(), -1);
1334 // labelList insideRegions(locationsInMesh.size(), -1);
1335 // forAll(locationsInMesh, i)
1336 // {
1337 // insideCells[i] = findCell
1338 // (
1339 // mesh_,
1340 // mergeDistance_*vector::one, //perturbVec,
1341 // locationsInMesh[i]
1342 // );
1343 // if (insideCells[i] != -1)
1344 // {
1345 // insideRegions[i] = cellRegion[insideCells[i]];
1346 // }
1347 // reduce(insideRegions[i], maxOp<label>());
1348 //
1349 // if (insideRegions[i] == -1)
1350 // {
1351 // // See if we can perturb a bit
1352 // insideCells[i] = findCell
1353 // (
1354 // mesh_,
1355 // mergeDistance_*vector::one, //perturbVec,
1356 // locationsInMesh[i]+mergeDistance_*vector::one
1357 // );
1358 // if (insideCells[i] != -1)
1359 // {
1360 // insideRegions[i] = cellRegion[insideCells[i]];
1361 // }
1362 // reduce(insideRegions[i], maxOp<label>());
1363 //
1364 // if (insideRegions[i] == -1)
1365 // {
1366 // FatalErrorInFunction
1367 // << "Cannot find locationInMesh " << locationsInMesh[i]
1368 // << " on any processor" << exit(FatalError);
1369 // }
1370 // }
1371 // }
1372 //
1373 //
1374 // // Check that all the locations outside the
1375 // // mesh do not conflict with those inside
1376 //
1377 // bool haveLeak = false;
1378 // forAll(locationsOutsideMesh, i)
1379 // {
1380 // // Find the region containing the point
1381 // label regioni = findRegion
1382 // (
1383 // mesh_,
1384 // cellRegion,
1385 // mergeDistance_*vector::one, //perturbVec,
1386 // locationsOutsideMesh[i]
1387 // );
1388 //
1389 // if (regioni != -1)
1390 // {
1391 // // Check for locationsOutsideMesh overlapping with inside ones
1392 // if (insideRegions.find(regioni) != -1)
1393 // {
1394 // haveLeak = true;
1395 // WarningInFunction
1396 // << "Outside location " << locationsOutsideMesh[i]
1397 // << " in region " << regioni
1398 // << " is connected to one of the inside points "
1399 // << locationsInMesh << endl;
1400 // }
1401 // }
1402 // }
1403 //
1404 //
1405 // autoPtr<mapPolyMesh> mapPtr;
1406 // if (returnReduceOr(haveLeak))
1407 // {
1408 // // Use shortestPathSet to provide a minimum set of faces needed
1409 // // to close hole. Tbd: maybe directly use wave?
1410 // meshSearch searchEngine(mesh_);
1411 // shortestPathSet leakPath
1412 // (
1413 // "leakPath",
1414 // mesh_,
1415 // searchEngine,
1416 // coordSet::coordFormatNames[coordSet::coordFormat::DISTANCE],
1417 // true,
1418 // 50, // tbd. Number of iterations. This is the maximum
1419 // // number of faces in the leak hole
1420 //
1421 // //pbm.groupPatchIDs()["wall"], // patch to grow from
1422 // meshedPatches(), // patch to grow from
1423 //
1424 // locationsInMesh,
1425 // locationsOutsideMesh,
1426 // isBlockedFace
1427 // );
1428 //
1429 //
1430 // // Use leak path to find minimum set of cells to delete
1431 // const labelList cellsToRemove
1432 // (
1433 // detectLeakCells
1434 // (
1435 // isBlockedFace,
1436 // leakPath.leakFaces(),
1437 // insideCells
1438 // )
1439 // );
1440 //
1441 // // Re-do intersections to find nearest unnamed surface
1442 // const label defaultRegion
1443 // (
1444 // surfaces().globalRegion
1445 // (
1446 // selectedSurfaces[0],
1447 // 0
1448 // )
1449 // );
1450 //
1451 // const labelList nearestRegion
1452 // (
1453 // nearestIntersection
1454 // (
1455 // selectedSurfaces,
1456 // defaultRegion
1457 // )
1458 // );
1459 //
1460 //
1461 // // Remove cells
1462 // removeCells cellRemover(mesh_);
1463 // const labelList exposedFaces
1464 // (
1465 // cellRemover.getExposedFaces(cellsToRemove)
1466 // );
1467 //
1468 // labelList exposedPatches(exposedFaces.size());
1469 // forAll(exposedFaces, i)
1470 // {
1471 // label facei = exposedFaces[i];
1472 // exposedPatches[i] = globalToMasterPatch[nearestRegion[facei]];
1473 // }
1474 //
1475 // mapPtr = doRemoveCells
1476 // (
1477 // cellsToRemove,
1478 // exposedFaces,
1479 // exposedPatches,
1480 // cellRemover
1481 // );
1482 //
1483 //
1484 // // Put the exposed faces into a special faceZone
1485 // {
1486 // // Add to "frozenFaces" zone
1487 // faceZoneMesh& faceZones = mesh_.faceZones();
1488 //
1489 // // Get current frozen faces (if any)
1490 // bitSet isFrozenFace(mesh_.nFaces());
1491 // label zonei = faceZones.findZoneID("frozenFaces");
1492 // if (zonei != -1)
1493 // {
1494 // const bitSet oldSet(mesh_.nFaces(), faceZones[zonei]);
1495 // isFrozenFace.set(oldSet);
1496 // }
1497 //
1498 // // Add newly exposed faces (if not yet in any faceZone!)
1499 // const labelList exposed
1500 // (
1501 // renumber
1502 // (
1503 // mapPtr().reverseFaceMap(),
1504 // exposedFaces
1505 // )
1506 // );
1507 // bitSet isZonedFace(mesh_.nFaces(), faceZones.zoneMap().toc());
1508 // for (const label facei : exposed)
1509 // {
1510 // if (!isZonedFace[facei])
1511 // {
1512 // isFrozenFace.set(facei);
1513 // }
1514 // }
1515 //
1516 // syncTools::syncFaceList
1517 // (
1518 // mesh_,
1519 // isFrozenFace,
1520 // orEqOp<unsigned int>(),
1521 // 0u
1522 // );
1523 //
1524 // // Add faceZone if non existing
1525 // faceZones.clearAddressing();
1526 // if (zonei == -1)
1527 // {
1528 // zonei = faceZones.size();
1529 //
1530 // faceZones.emplace_back
1531 // (
1532 // "frozenFaces", // name
1533 // zonei, // index
1534 // faceZones
1535 // );
1536 // }
1537 //
1538 // // Update faceZone with new contents
1539 // labelList frozenFaces(isFrozenFace.toc());
1540 // boolList frozenFlip(frozenFaces.size(), false);
1541 //
1542 // faceZones[zonei].resetAddressing
1543 // (
1544 // std::move(frozenFaces),
1545 // std::move(frozenFlip)
1546 // );
1547 // }
1548 //
1549 //
1550 // //// Put the exposed points into a special pointZone
1551 // //if (false)
1552 // //{
1553 // // const labelList meshFaceIDs
1554 // // (
1555 // // renumber
1556 // // (
1557 // // mapPtr().reverseFaceMap(),
1558 // // exposedFaces
1559 // // )
1560 // // );
1561 // // const uindirectPrimitivePatch pp
1562 // // (
1563 // // UIndirectList<face>(mesh_.faces(), meshFaceIDs),
1564 // // mesh_.points()
1565 // // );
1566 // //
1567 // // // Count number of faces per edge
1568 // // const labelListList& edgeFaces = pp.edgeFaces();
1569 // // labelList nEdgeFaces(edgeFaces.size());
1570 // // forAll(edgeFaces, edgei)
1571 // // {
1572 // // nEdgeFaces[edgei] = edgeFaces[edgei].size();
1573 // // }
1574 // //
1575 // // // Sync across processor patches
1576 // // if (Pstream::parRun())
1577 // // {
1578 // // const globalMeshData& globalData = mesh_.globalData();
1579 // // const mapDistribute& map = globalData.globalEdgeSlavesMap();
1580 // // const indirectPrimitivePatch& cpp =
1581 // // globalData.coupledPatch();
1582 // //
1583 // // // Match pp edges to coupled edges
1584 // // labelList patchEdges;
1585 // // labelList coupledEdges;
1586 // // PackedBoolList sameEdgeOrientation;
1587 // // PatchTools::matchEdges
1588 // // (
1589 // // pp,
1590 // // cpp,
1591 // // patchEdges,
1592 // // coupledEdges,
1593 // // sameEdgeOrientation
1594 // // );
1595 // //
1596 // //
1597 // // // Convert patch-edge data into cpp-edge data
1598 // // labelList coupledNEdgeFaces(map.constructSize(), Zero);
1599 // // UIndirectList<label>(coupledNEdgeFaces, coupledEdges) =
1600 // // UIndirectList<label>(nEdgeFaces, patchEdges);
1601 // //
1602 // // // Synchronise
1603 // // globalData.syncData
1604 // // (
1605 // // coupledNEdgeFaces,
1606 // // globalData.globalEdgeSlaves(),
1607 // // globalData.globalEdgeTransformedSlaves(),
1608 // // map,
1609 // // plusEqOp<label>()
1610 // // );
1611 // //
1612 // // // Convert back from cpp-edge to patch-edge
1613 // // UIndirectList<label>(nEdgeFaces, patchEdges) =
1614 // // UIndirectList<label>(coupledNEdgeFaces, coupledEdges);
1615 // // }
1616 // //
1617 // // // Freeze all internal points
1618 // // bitSet isFrozenPoint(mesh_.nPoints());
1619 // // forAll(nEdgeFaces, edgei)
1620 // // {
1621 // // if (nEdgeFaces[edgei] != 1)
1622 // // {
1623 // // const edge& e = pp.edges()[edgei];
1624 // // isFrozenPoint.set(pp.meshPoints()[e[0]]);
1625 // // isFrozenPoint.set(pp.meshPoints()[e[1]]);
1626 // // }
1627 // // }
1628 // // // Add to "frozenPoints" zone
1629 // // pointZoneMesh& pointZones = mesh_.pointZones();
1630 // // pointZones.clearAddressing();
1631 // // label zonei = pointZones.findZoneID("frozenPoints");
1632 // // if (zonei != -1)
1633 // // {
1634 // // const bitSet oldSet(mesh_.nPoints(), pointZones[zonei]);
1635 // // // Add to isFrozenPoint
1636 // // isFrozenPoint.set(oldSet);
1637 // // }
1638 // //
1639 // // syncTools::syncPointList
1640 // // (
1641 // // mesh_,
1642 // // isFrozenPoint,
1643 // // orEqOp<unsigned int>(),
1644 // // 0u
1645 // // );
1646 // //
1647 // // if (zonei == -1)
1648 // // {
1649 // // zonei = pointZones.size();
1650 // //
1651 // // pointZones.emplace_back
1652 // // (
1653 // // "frozenPoints", // name
1654 // // isFrozenPoint.sortedToc(), // addressing
1655 // // zonei, // index
1656 // // pointZones // pointZoneMesh
1657 // // );
1658 // // }
1659 // //}
1660 //
1661 //
1662 // if (debug&meshRefinement::MESH)
1663 // {
1664 // const_cast<Time&>(mesh_.time())++;
1665 //
1666 // Pout<< "Writing current mesh to time "
1667 // << timeName() << endl;
1668 // write
1669 // (
1670 // meshRefinement::debugType(debug),
1671 // meshRefinement::writeType
1672 // (
1673 // meshRefinement::writeLevel()
1674 // | meshRefinement::WRITEMESH
1675 // ),
1676 // mesh_.time().path()/timeName()
1677 // );
1678 // Pout<< "Dumped mesh in = "
1679 // << mesh_.time().cpuTimeIncrement() << " s\n" << nl << endl;
1680 // }
1681 // }
1682 // return mapPtr;
1683 //}
1684 
1685 
1686 //Foam::autoPtr<Foam::mapDistribute> Foam::meshRefinement::nearestFace
1687 //(
1688 // const globalIndex& globalSeedFaces,
1689 // const labelList& seedFaces,
1690 // const labelList& closureFaces
1691 //) const
1692 //{
1693 // // Takes a set of faces for which we have information (seedFaces,
1694 // // globalSeedFaces - these are e.g. intersected faces) and walks out
1695 // // across nonSeedFace. Returns for
1696 // // every nonSeedFace the nearest seed face (in global indexing).
1697 // // Used e.g. in hole closing. Assumes that seedFaces, closureFaces
1698 // // are a small subset of the master
1699 // // so are not large - it uses edge addressing on the closureFaces
1700 //
1701 //
1702 // if (seedFaces.size() != globalSeedFaces.localSize())
1703 // {
1704 // FatalErrorInFunction << "problem : seedFaces:" << seedFaces.size()
1705 // << " globalSeedFaces:" << globalSeedFaces.localSize()
1706 // << exit(FatalError);
1707 // }
1708 //
1709 // //// Mark mesh points that are used by any closureFaces. This is for
1710 // //// initial filtering
1711 // //bitSet isNonSeedPoint(mesh.nPoints());
1712 // //for (const label facei : closureFaces)
1713 // //{
1714 // // isNonSeedPoint.set(mesh_.faces()[facei]);
1715 // //}
1716 // //syncTools::syncPointList
1717 // //(
1718 // // mesh_,
1719 // // isNonSeedPoint,
1720 // // orEqOp<unsigned int>(),
1721 // // 0u
1722 // //);
1723 //
1724 // // Make patch
1725 // const uindirectPrimitivePatch pp
1726 // (
1727 // IndirectList<face>(mesh_.faces(), closureFaces),
1728 // mesh_.points()
1729 // );
1730 // const edgeList& edges = pp.edges();
1731 // const labelList& mp = pp.meshPoints();
1732 // const label nBndEdges = pp.nEdges() - pp.nInternalEdges();
1733 //
1734 // // For all faces in seedFaces mark the edge with a face. No special
1735 // // handling for multiple faces sharing the edge - first one wins
1736 // EdgeMap<label> edgeMap(pp.nEdges());
1737 // for (const label facei : seedFaces)
1738 // {
1739 // const label globalFacei = globalSeedFaces.toGlobal(facei);
1740 // const face& f = mesh_.faces()[facei];
1741 // forAll(f, fp)
1742 // {
1743 // label nextFp = f.fcIndex(fp);
1744 // edgeMap.insert(edge(f[fp], f[nextFp]), globalFacei);
1745 // }
1746 // }
1747 // syncTools::syncEdgeMap(mesh_, edgeMap, maxEqOp<label>());
1748 //
1749 //
1750 //
1751 // // Seed
1752 // DynamicList<label> initialEdges(2*nBndEdges);
1753 // DynamicList<edgeTopoDistanceData<label, uindirectPrimitivePatch>>
1754 // initialEdgesInfo(2*nBndEdges);
1755 // forAll(edges, edgei)
1756 // {
1757 // const edge& e = edges[edgei];
1758 // const edge meshE = edge(mp[e[0]], mp[e[1]]);
1759 //
1760 // const auto iter = edgeMap.cfind(meshE);
1761 // if (iter.good())
1762 // {
1763 // initialEdges.append(edgei);
1764 // initialEdgesInfo.append
1765 // (
1766 // edgeTopoDistanceData<label, uindirectPrimitivePatch>
1767 // (
1768 // 0, // distance
1769 // iter.val() // globalFacei
1770 // )
1771 // );
1772 // }
1773 // }
1774 //
1775 // // Data on all edges and faces
1776 // List<edgeTopoDistanceData<label, uindirectPrimitivePatch>> allEdgeInfo
1777 // (
1778 // pp.nEdges()
1779 // );
1780 // List<edgeTopoDistanceData<label, uindirectPrimitivePatch>> allFaceInfo
1781 // (
1782 // pp.size()
1783 // );
1784 //
1785 // // Walk
1786 // PatchEdgeFaceWave
1787 // <
1788 // uindirectPrimitivePatch,
1789 // edgeTopoDistanceData<label, uindirectPrimitivePatch>
1790 // > calc
1791 // (
1792 // mesh_,
1793 // pp,
1794 // initialEdges,
1795 // initialEdgesInfo,
1796 // allEdgeInfo,
1797 // allFaceInfo,
1798 // returnReduce(pp.nEdges(), sumOp<label>())
1799 // );
1800 //
1801 //
1802 // // Per non-seed face the seed face
1803 // labelList closureToSeed(pp.size(), -1);
1804 // forAll(allFaceInfo, facei)
1805 // {
1806 // if (allFaceInfo[facei].valid(calc.data()))
1807 // {
1808 // closureToSeed[facei] = allFaceInfo[facei].data();
1809 // }
1810 // }
1811 //
1812 // List<Map<label>> compactMap;
1813 // return autoPtr<mapDistribute>::New
1814 // (
1815 // globalSeedFaces,
1816 // closureToSeed,
1817 // compactMap
1818 // );
1819 //}
1820 
1821 
1823 (
1824  const labelList& globalToMasterPatch,
1825  const labelList& globalToSlavePatch,
1826  const pointField& locationsInMesh,
1827  const wordList& zonesInMesh,
1828  const pointField& locationsOutsideMesh,
1829  const labelList& selectedSurfaces
1830 )
1831 {
1832  // Problem: this is based on cached intersection information. This might
1833  // not include the surface we actually want to use. In which case the
1834  // surface would not be seen as intersected!
1835  boolList isBlockedFace;
1836  selectIntersectedFaces(selectedSurfaces, isBlockedFace);
1837 
1838  // Determine cell regions
1839  const regionSplit cellRegion(mesh_, isBlockedFace);
1840 
1841  // Detect locationsInMesh regions
1842  labelList insideCells(locationsInMesh.size(), -1);
1843  labelList insideRegions(locationsInMesh.size(), -1);
1844  forAll(locationsInMesh, i)
1845  {
1846  insideCells[i] = findCell
1847  (
1848  mesh_,
1849  mergeDistance_*vector::one, //perturbVec,
1850  locationsInMesh[i]
1851  );
1852  if (insideCells[i] != -1)
1853  {
1854  insideRegions[i] = cellRegion[insideCells[i]];
1855  }
1856  reduce(insideRegions[i], maxOp<label>());
1857 
1858  if (insideRegions[i] == -1)
1859  {
1860  // See if we can perturb a bit
1861  insideCells[i] = findCell
1862  (
1863  mesh_,
1864  mergeDistance_*vector::one, //perturbVec,
1865  locationsInMesh[i]+mergeDistance_*vector::one
1866  );
1867  if (insideCells[i] != -1)
1868  {
1869  insideRegions[i] = cellRegion[insideCells[i]];
1870  }
1871  reduce(insideRegions[i], maxOp<label>());
1872 
1873  if (insideRegions[i] == -1)
1874  {
1876  << "Cannot find locationInMesh " << locationsInMesh[i]
1877  << " on any processor" << exit(FatalError);
1878  }
1879  }
1880  }
1881 
1882 
1883  // Check that all the locations outside the
1884  // mesh do not conflict with those inside
1885 
1886  bool haveLeak = false;
1887  forAll(locationsOutsideMesh, i)
1888  {
1889  // Find the region containing the point
1890  label regioni = findRegion
1891  (
1892  mesh_,
1893  cellRegion,
1894  mergeDistance_*vector::one, //perturbVec,
1895  locationsOutsideMesh[i]
1896  );
1897 
1898  if (regioni != -1)
1899  {
1900  // Check for locationsOutsideMesh overlapping with inside ones
1901  if (insideRegions.find(regioni) != -1)
1902  {
1903  haveLeak = true;
1905  << "Outside location " << locationsOutsideMesh[i]
1906  << " in region " << regioni
1907  << " is connected to one of the inside points "
1908  << locationsInMesh << endl;
1909  }
1910  }
1911  }
1912 
1913 
1914  autoPtr<mapPolyMesh> mapPtr;
1915  if (returnReduceOr(haveLeak))
1916  {
1917  // Use holeToFace to provide a minimum set of faces needed
1918  // to close hole.
1919 
1920  const List<pointField> allLocations
1921  (
1923  (
1924  locationsInMesh,
1925  zonesInMesh,
1926  locationsOutsideMesh
1927  )
1928  );
1929 
1930  const labelList blockingFaces(findIndices(isBlockedFace, true));
1931 
1932  labelList closureFaces;
1933  labelList closureToBlocked;
1934  autoPtr<mapDistribute> closureMapPtr;
1935  {
1936  const globalIndex globalBlockingFaces(blockingFaces.size());
1937 
1938  closureMapPtr = holeToFace::calcClosure
1939  (
1940  mesh_,
1941  allLocations,
1942  blockingFaces,
1943  globalBlockingFaces,
1944  true, // allow erosion
1945 
1946  closureFaces,
1947  closureToBlocked
1948  );
1949 
1950  if (debug)
1951  {
1952  Pout<< "meshRefinement::blockLeakFaces :"
1953  << " found closure faces:" << closureFaces.size()
1954  << " map:" << bool(closureMapPtr) << endl;
1955  }
1956 
1957  if (!closureMapPtr)
1958  {
1960  << "have leak but did not find any closure faces"
1961  << exit(FatalError);
1962  }
1963  }
1964 
1965  // Baffle faces
1966  labelList ownPatch(mesh_.nFaces(), -1);
1967  labelList neiPatch(mesh_.nFaces(), -1);
1968 
1969  // Keep (global) boundary faces in their patch
1970  {
1971  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
1972  for (label patchi = 0; patchi < patches.nNonProcessor(); ++patchi)
1973  {
1974  const polyPatch& pp = patches[patchi];
1975 
1976  forAll(pp, i)
1977  {
1978  ownPatch[pp.start()+i] = patchi;
1979  neiPatch[pp.start()+i] = patchi;
1980  }
1981  }
1982  }
1983 
1984  const faceZoneMesh& fzs = mesh_.faceZones();
1985 
1986  // Mark zone per face
1987  labelList faceToZone(mesh_.nFaces(), -1);
1988  boolList faceToFlip(mesh_.nFaces(), false);
1989  forAll(fzs, zonei)
1990  {
1991  const labelList& addressing = fzs[zonei];
1992  const boolList& flipMap = fzs[zonei].flipMap();
1993 
1994  forAll(addressing, i)
1995  {
1996  faceToZone[addressing[i]] = zonei;
1997  faceToFlip[addressing[i]] = flipMap[i];
1998  }
1999  }
2000 
2001 
2002  // Fetch patch and zone info from blockingFaces
2003  labelList packedOwnPatch(labelUIndList(ownPatch, blockingFaces));
2004  closureMapPtr->distribute(packedOwnPatch);
2005  labelList packedNeiPatch(labelUIndList(neiPatch, blockingFaces));
2006  closureMapPtr->distribute(packedNeiPatch);
2007  labelList packedZone(labelUIndList(faceToZone, blockingFaces));
2008  closureMapPtr->distribute(packedZone);
2009  boolList packedFlip(UIndirectList<bool>(faceToFlip, blockingFaces));
2010  closureMapPtr->distribute(packedFlip);
2011  forAll(closureFaces, i)
2012  {
2013  const label facei = closureFaces[i];
2014  const label sloti = closureToBlocked[i];
2015 
2016  if (sloti != -1)
2017  {
2018  // TBD. how to orient own/nei patch?
2019  ownPatch[facei] = packedOwnPatch[sloti];
2020  neiPatch[facei] = packedNeiPatch[sloti];
2021  faceToZone[facei] = packedZone[sloti];
2022  faceToFlip[facei] = packedFlip[sloti];
2023  }
2024  }
2025 
2026 
2027  // Add faces to faceZone. For now do this outside of createBaffles
2028  // below
2029  {
2030  List<DynamicList<label>> zoneToFaces(fzs.size());
2031  List<DynamicList<bool>> zoneToFlip(fzs.size());
2032 
2033  // Add any to-be-patched face
2034  forAll(faceToZone, facei)
2035  {
2036  const label zonei = faceToZone[facei];
2037  if (zonei != -1)
2038  {
2039  zoneToFaces[zonei].append(facei);
2040  zoneToFlip[zonei].append(faceToFlip[facei]);
2041  }
2042  }
2043 
2044  forAll(zoneToFaces, zonei)
2045  {
2047  (
2048  fzs[zonei].name(),
2049  zoneToFaces[zonei],
2050  zoneToFlip[zonei],
2051  mesh_
2052  );
2053  }
2054  }
2055 
2056  // Put the points of closureFaces into a special pointZone
2057  {
2059  (
2060  UIndirectList<face>(mesh_.faces(), closureFaces),
2061  mesh_.points()
2062  );
2063 
2064  // Count number of faces per edge
2065  const labelList nEdgeFaces(countEdgeFaces(pp));
2066 
2067  // Freeze all internal points
2068  bitSet isFrozenPoint(mesh_.nPoints());
2069  forAll(nEdgeFaces, edgei)
2070  {
2071  if (nEdgeFaces[edgei] != 1)
2072  {
2073  const edge& e = pp.edges()[edgei];
2074  isFrozenPoint.set(pp.meshPoints()[e[0]]);
2075  isFrozenPoint.set(pp.meshPoints()[e[1]]);
2076  }
2077  }
2078 
2079  // Lookup/add pointZone and include its points
2080  pointZoneMesh& pointZones =
2081  const_cast<pointZoneMesh&>(mesh_.pointZones());
2082  const label zonei = addPointZone("frozenPoints");
2083  const bitSet oldSet(mesh_.nPoints(), pointZones[zonei]);
2084  isFrozenPoint.set(oldSet);
2085 
2087  (
2088  mesh_,
2089  isFrozenPoint,
2091  0u
2092  );
2093 
2094  // Override addressing
2095  pointZones.clearAddressing();
2096  pointZones[zonei] = isFrozenPoint.sortedToc();
2097  }
2098 
2099 
2100 
2101  // Create baffles for faces
2102  mapPtr = createBaffles(ownPatch, neiPatch);
2103 
2105  //{
2106  // // Add newly exposed faces (if not yet in any faceZone!)
2107  // const labelList exposed
2108  // (
2109  // renumber
2110  // (
2111  // mapPtr().reverseFaceMap(),
2112  // blockingFaces
2113  // )
2114  // );
2115  //
2116  // surfaceZonesInfo::addFaceZone
2117  // (
2118  // "frozenFaces",
2119  // exposed,
2120  // boolList(exposed.size(), false),
2121  // mesh_
2122  // );
2123  //}
2124 
2125 
2127  {
2128  const_cast<Time&>(mesh_.time())++;
2129 
2130  Pout<< "Writing current mesh to time "
2131  << timeName() << endl;
2132  write
2133  (
2136  (
2139  ),
2140  mesh_.time().path()/timeName()
2141  );
2142  Pout<< "Dumped mesh in = "
2143  << mesh_.time().cpuTimeIncrement() << " s\n" << nl << endl;
2144  }
2145  }
2146  return mapPtr;
2147 }
2148 
2149 
2150 // ************************************************************************* //
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:136
void markOutsideFaces(const labelList &cellLevel, const labelList &neiLevel, const labelList &refineCell, bitSet &isOutsideFace) const
Mark faces on interface between set and rest.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
void clearAddressing()
Clear addressing.
Definition: ZoneMesh.C:815
uint8_t direction
Definition: direction.H:46
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:504
List< cell > cellList
List of cell.
Definition: cellListFwd.H:39
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
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:598
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:517
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute List data using default commsType, default flip/negate operator.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:489
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
autoPtr< mapPolyMesh > removeGapCells(const scalar planarAngle, const labelList &minSurfaceLevel, const labelList &globalToMasterPatch, const label growIter)
Detect gapRefinement cells and remove them.
void setSize(const label n, unsigned int val=0u)
Alias for resize()
Definition: PackedList.H:785
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static writeType writeLevel()
Get/set write level.
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
label countFaceDirs(const bitSet &isOutsideFace, const label celli) const
Count number of faces on cell that are in set.
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
Base class of (analytical or triangulated) surface. Encapsulates all the search routines. WIP.
Ignore writing from objectRegistry::writeObject()
label nFaces() const noexcept
Number of mesh faces.
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.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
const labelList & meshPoints() const
Return labelList of mesh points in patch.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:152
word timeName
Definition: getTimeIndex.H:3
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
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
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:1065
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dimensionedScalar cos(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition: edge.H:59
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:137
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
autoPtr< mapPolyMesh > blockLeakFaces(const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const pointField &locationsInMesh, const wordList &zonesInMesh, const pointField &locationsOutsideMesh, const labelList &selectedSurfaces)
Baffle faces to break any leak from inside to outside.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
void growSet(const labelList &neiLevel, const bitSet &isOutsideFace, labelList &refineCell, label &nRefine) const
Add one layer of cells to set.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nInternalFaces() const noexcept
Number of internal faces.
Vector< scalar > vector
Definition: vector.H:57
Container with cells to refine. Refinement given as single direction.
Definition: refineCell.H:52
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
static labelList getUnnamedSurfaces(const PtrList< surfaceZonesInfo > &surfList)
Get indices of unnamed surfaces (surfaces without faceZoneName)
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
const wordList surface
Standard surface field types (scalar, vector, tensor, etc)
int debug
Static debugging option.
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
vector point
Point is a vector.
Definition: point.H:37
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:139
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:56
const polyBoundaryMesh & patches
Nothing to be read.
void selectIntersectedFaces(const labelList &surfaces, boolList &isBlockedFace) const
Faces currently on boundary or intersected by surface.
static List< pointField > zonePoints(const pointField &locationsInMesh, const wordList &zonesInMesh, const pointField &locationsOutsideMesh)
Helper: per zone (entry in zonesInMesh) the locations with.
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)
constexpr label labelMax
Definition: label.H:55
writeType
Enumeration for what to write. Used as a bit-pattern.
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
static autoPtr< mapDistribute > calcClosure(const polyMesh &mesh, const List< pointField > &zonePoints, const labelList &blockedFaces, const globalIndex &globalBlockedFaces, const bool erode, labelList &closureFaces, labelList &closureToBlocked)
Optional direct use to generate the set of faces and the method to.
Definition: holeToFace.C:1190
debugType
Enumeration for what to debug. Used as a bit-pattern.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
List< label > labelList
A List of labels.
Definition: List.H:62
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
label nNonProcessor() const
The number of patches before the first processor patch.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
List< bool > boolList
A List of bools.
Definition: List.H:60
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Do not request registration (bool: false)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
static label addFaceZone(const word &name, const labelList &addressing, const boolList &flipMap, polyMesh &mesh)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127