meshRefinement.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-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2024 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 \*---------------------------------------------------------------------------*/
28 
29 #include "meshRefinement.H"
30 #include "volMesh.H"
31 #include "volFields.H"
32 #include "surfaceMesh.H"
33 #include "syncTools.H"
34 #include "Time.H"
35 #include "refinementSurfaces.H"
36 #include "refinementFeatures.H"
37 #include "decompositionMethod.H"
38 #include "regionSplit.H"
39 #include "fvMeshDistribute.H"
40 #include "indirectPrimitivePatch.H"
41 #include "polyTopoChange.H"
42 #include "removeCells.H"
43 #include "mapDistributePolyMesh.H"
44 #include "localPointRegion.H"
45 #include "pointMesh.H"
46 #include "pointFields.H"
47 #include "slipPointPatchFields.H"
51 #include "processorPointPatch.H"
52 #include "globalIndex.H"
53 #include "meshTools.H"
54 #include "OFstream.H"
55 #include "Random.H"
56 #include "searchableSurfaces.H"
57 #include "treeBoundBox.H"
59 #include "fvMeshTools.H"
60 #include "motionSmoother.H"
61 #include "faceSet.H"
62 #include "topoDistanceData.H"
63 #include "FaceCellWave.H"
64 #include "PackedBoolList.H"
65 
66 // Leak path
67 #include "shortestPathSet.H"
68 #include "meshSearch.H"
69 
70 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
71 
72 namespace Foam
73 {
74  defineTypeNameAndDebug(meshRefinement, 0);
75 }
76 
77 
78 const Foam::Enum
79 <
81 >
83 ({
84  { debugType::MESH, "mesh" },
85  { debugType::OBJINTERSECTIONS, "intersections" },
86  { debugType::FEATURESEEDS, "featureSeeds" },
87  { debugType::ATTRACTION, "attraction" },
88  { debugType::LAYERINFO, "layerInfo" },
89 });
90 
91 
92 //const Foam::Enum
93 //<
94 // Foam::meshRefinement::outputType
95 //>
96 //Foam::meshRefinement::outputTypeNames
97 //({
98 // { outputType::OUTPUTLAYERINFO, "layerInfo" }
99 //});
100 
101 
102 const Foam::Enum
103 <
105 >
107 ({
108  { writeType::WRITEMESH, "mesh" },
109  { writeType::NOWRITEREFINEMENT, "noRefinement" },
110  { writeType::WRITELEVELS, "scalarLevels" },
111  { writeType::WRITELAYERSETS, "layerSets" },
112  { writeType::WRITELAYERFIELDS, "layerFields" },
113 });
114 
115 
116 Foam::meshRefinement::writeType Foam::meshRefinement::writeLevel_;
117 
118 //Foam::meshRefinement::outputType Foam::meshRefinement::outputLevel_;
119 
120 // Inside/outside test for polyMesh:.findCell()
121 // 2.4.x : default = polyMesh::FACE_DIAG_TRIS
122 // 1712 : default = polyMesh::CELL_TETS
123 //
124 // - CELL_TETS is better with concave cells, but much slower.
125 // - use faster method (FACE_DIAG_TRIS) here
126 
129 
130 
131 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
132 
133 Foam::label Foam::meshRefinement::globalFaceCount(const labelList& elems) const
134 {
135  //- Check for duplicates
136  const bitSet isElem(mesh_.nFaces(), elems);
137  if (label(isElem.count()) != elems.size())
138  {
139  FatalErrorInFunction << "Problem Duplicates:"
140  << " isElem:" << isElem.count()
141  << " elems:" << elems.size()
142  << exit(FatalError);
143  }
144 
145  //- Check for same entries on coupled faces
146  {
147  bitSet isElem2(isElem);
148  syncTools::swapFaceList(mesh_, isElem2);
149 
150  for
151  (
152  label facei = mesh_.nInternalFaces();
153  facei < mesh_.nFaces();
154  facei++
155  )
156  {
157  if (isElem2[facei] != isElem[facei])
158  {
160  << "at face:" << facei
161  << " at:" << mesh_.faceCentres()[facei]
162  << " patch:" << mesh_.boundaryMesh().whichPatch(facei)
163  << " isElem:" << isElem[facei]
164  << " isElem2:" << isElem2[facei]
165  << exit(FatalError);
166  }
167  }
168  }
169 
170  //- Count number of master elements
171  const bitSet isMaster(syncTools::getMasterFaces(mesh_));
172  label count = 0;
173  for (const label i : isElem)
174  {
175  if (isMaster[i])
176  {
177  count++;
178  }
179  }
180  return returnReduce(count, sumOp<label>());
181 }
182 
183 
184 void Foam::meshRefinement::calcNeighbourData
185 (
186  labelList& neiLevel,
187  pointField& neiCc
188 ) const
189 {
190  const labelList& cellLevel = meshCutter_.cellLevel();
191  const pointField& cellCentres = mesh_.cellCentres();
192 
193  const label nBoundaryFaces = mesh_.nBoundaryFaces();
194 
195  if (neiLevel.size() != nBoundaryFaces || neiCc.size() != nBoundaryFaces)
196  {
198  << nBoundaryFaces << " neiLevel:" << neiLevel.size()
199  << abort(FatalError);
200  }
201 
202  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
203 
204  labelHashSet addedPatchIDSet(meshedPatches());
205 
206  forAll(patches, patchi)
207  {
208  const polyPatch& pp = patches[patchi];
209 
210  const labelUList& faceCells = pp.faceCells();
211  const vectorField::subField faceCentres = pp.faceCentres();
212  const vectorField::subField faceAreas = pp.faceAreas();
213 
214  label bFacei = pp.start()-mesh_.nInternalFaces();
215 
216  if (pp.coupled())
217  {
218  forAll(faceCells, i)
219  {
220  neiLevel[bFacei] = cellLevel[faceCells[i]];
221  neiCc[bFacei] = cellCentres[faceCells[i]];
222  bFacei++;
223  }
224  }
225  else if (addedPatchIDSet.found(patchi))
226  {
227  // Face was introduced from cell-cell intersection. Try to
228  // reconstruct other side cell(centre). Three possibilities:
229  // - cells same size.
230  // - preserved cell smaller. Not handled.
231  // - preserved cell larger.
232  forAll(faceCells, i)
233  {
234  // Extrapolate the face centre.
235  const vector fn = normalised(faceAreas[i]);
236 
237  label own = faceCells[i];
238  label ownLevel = cellLevel[own];
239  label faceLevel = meshCutter_.faceLevel(pp.start()+i);
240  if (faceLevel < 0)
241  {
242  // Due to e.g. face merging no longer a consistent
243  // refinementlevel of face. Assume same as cell
244  faceLevel = ownLevel;
245  }
246 
247  // Normal distance from face centre to cell centre
248  scalar d = ((faceCentres[i] - cellCentres[own]) & fn);
249  if (faceLevel > ownLevel)
250  {
251  // Other cell more refined. Adjust normal distance
252  d *= 0.5;
253  }
254  neiLevel[bFacei] = faceLevel;
255  // Calculate other cell centre by extrapolation
256  neiCc[bFacei] = faceCentres[i] + d*fn;
257  bFacei++;
258  }
259  }
260  else
261  {
262  forAll(faceCells, i)
263  {
264  neiLevel[bFacei] = cellLevel[faceCells[i]];
265  neiCc[bFacei] = faceCentres[i];
266  bFacei++;
267  }
268  }
269  }
270 
271  // Swap coupled boundaries. Apply separation to cc since is coordinate.
273  syncTools::swapBoundaryFaceList(mesh_, neiLevel);
274 }
275 
276 
277 void Foam::meshRefinement::calcCellCellRays
278 (
279  const pointField& neiCc,
280  const labelList& neiLevel,
281  const labelList& testFaces,
282  pointField& start,
283  pointField& end,
284  labelList& minLevel
285 ) const
286 {
287  const labelList& cellLevel = meshCutter_.cellLevel();
288  const pointField& cellCentres = mesh_.cellCentres();
289 
290 
291  // Mark all non-coupled or coupled+master faces. Leaves only slave of
292  // coupled unset.
293  bitSet isMaster(mesh_.nBoundaryFaces(), true);
294  {
295  for (const polyPatch& pp : mesh_.boundaryMesh())
296  {
297  if (pp.coupled() && !refCast<const coupledPolyPatch>(pp).owner())
298  {
299  isMaster.unset(labelRange(pp.offset(), pp.size()));
300  }
301  }
302  }
303 
304 
305  start.setSize(testFaces.size());
306  end.setSize(testFaces.size());
307  minLevel.setSize(testFaces.size());
308 
309  forAll(testFaces, i)
310  {
311  const label facei = testFaces[i];
312  const label own = mesh_.faceOwner()[facei];
313 
314  if (mesh_.isInternalFace(facei))
315  {
316  const label nei = mesh_.faceNeighbour()[facei];
317 
318  start[i] = cellCentres[own];
319  end[i] = cellCentres[nei];
320  minLevel[i] = min(cellLevel[own], cellLevel[nei]);
321  }
322  else
323  {
324  const label bFacei = facei - mesh_.nInternalFaces();
325 
326  if (isMaster[bFacei])
327  {
328  start[i] = cellCentres[own];
329  end[i] = neiCc[bFacei];
330  }
331  else
332  {
333  // Slave face
334  start[i] = neiCc[bFacei];
335  end[i] = cellCentres[own];
336  }
337  minLevel[i] = min(cellLevel[own], neiLevel[bFacei]);
338  }
339  }
340 
341  // Extend segments a bit
342  {
343  const vectorField smallVec(ROOTSMALL*(end-start));
344  start -= smallVec;
345  end += smallVec;
346  }
347 }
348 
351 {
352  // Stats on edges to test. Count proc faces only once.
353  bitSet isMasterFace(syncTools::getMasterFaces(mesh_));
354 
355  {
356  const label nMasterFaces =
357  returnReduce(isMasterFace.count(), sumOp<label>());
358 
359  label nChangedFaces = 0;
360  forAll(changedFaces, i)
361  {
362  if (isMasterFace.test(changedFaces[i]))
363  {
364  ++nChangedFaces;
365  }
366  }
367  reduce(nChangedFaces, sumOp<label>());
368 
369  if (!dryRun_)
370  {
371  Info<< "Edge intersection testing:" << nl
372  << " Number of edges : " << nMasterFaces << nl
373  << " Number of edges to retest : " << nChangedFaces
374  << endl;
375  }
376  }
377 
378 
379  // Get boundary face centre and level. Coupled aware.
380  labelList neiLevel(mesh_.nBoundaryFaces());
381  pointField neiCc(mesh_.nBoundaryFaces());
382  calcNeighbourData(neiLevel, neiCc);
383 
384  // Collect segments we want to test for
385  pointField start(changedFaces.size());
386  pointField end(changedFaces.size());
387  {
388  labelList minLevel;
389  calcCellCellRays
390  (
391  neiCc,
392  neiLevel,
393  changedFaces,
394  start,
395  end,
396  minLevel
397  );
398  }
399 
400 
401  // Do tests in one go
402  labelList surfaceHit;
403  {
404  labelList surfaceLevel;
405  surfaces_.findHigherIntersection
406  (
407  shells_,
408  start,
409  end,
410  labelList(start.size(), -1), // accept any intersection
411  surfaceHit,
412  surfaceLevel
413  );
414  }
415 
416  // Keep just surface hit
417  forAll(surfaceHit, i)
418  {
419  surfaceIndex_[changedFaces[i]] = surfaceHit[i];
420  }
421 
422  // Make sure both sides have same information. This should be
423  // case in general since same vectors but just to make sure.
424  syncTools::syncFaceList(mesh_, surfaceIndex_, maxEqOp<label>());
425 
426  label nTotHits = returnReduce(countHits(), sumOp<label>());
427 
428  if (!dryRun_)
429  {
430  Info<< " Number of intersected edges : " << nTotHits << endl;
431  }
432 
433  // Set files to same time as mesh
434  setInstance(mesh_.facesInstance());
435 }
436 
437 
438 void Foam::meshRefinement::nearestFace
439 (
440  const labelUList& startFaces,
441  const bitSet& isBlockedFace,
442 
443  autoPtr<mapDistribute>& mapPtr,
444  labelList& faceToStart,
445  const label nIter
446 ) const
447 {
448  // From startFaces walk out (but not through isBlockedFace). Returns
449  // faceToStart which is the index into startFaces (or rather distributed
450  // version of it). E.g.
451  // pointField startFc(mesh.faceCentres(), startFaces);
452  // mapPtr().distribute(startFc);
453  // forAll(faceToStart, facei)
454  // {
455  // const label sloti = faceToStart[facei];
456  // if (sloti != -1)
457  // {
458  // Pout<< "face:" << mesh.faceCentres()[facei]
459  // << " nearest:" << startFc[sloti]
460  // << endl;
461  // }
462  // }
463 
464 
465  // Consecutive global numbering for start elements
466  const globalIndex globalStart(startFaces.size());
467 
468  // Field on cells and faces.
469  List<topoDistanceData<label>> cellData(mesh_.nCells());
470  List<topoDistanceData<label>> faceData(mesh_.nFaces());
471 
472  // Mark blocked faces to there not visited
473  for (const label facei : isBlockedFace)
474  {
475  faceData[facei] = topoDistanceData<label>(0, -1);
476  }
477 
478  List<topoDistanceData<label>> startData(startFaces.size());
479  forAll(startFaces, i)
480  {
481  const label facei = startFaces[i];
482  if (isBlockedFace[facei])
483  {
484  FatalErrorInFunction << "Start face:" << facei
485  << " at:" << mesh_.faceCentres()[facei]
486  << " is also blocked" << exit(FatalError);
487  }
488  startData[i] = topoDistanceData<label>(0, globalStart.toGlobal(i));
489  }
490 
491 
492  // Propagate information inwards
493  FaceCellWave<topoDistanceData<label>> deltaCalc
494  (
495  mesh_,
496  startFaces,
497  startData,
498  faceData,
499  cellData,
500  0
501  );
502  deltaCalc.iterate(nIter);
503 
504  // And extract
505 
506  faceToStart.setSize(mesh_.nFaces());
507  faceToStart = -1;
508  bool haveWarned = false;
509  forAll(faceData, facei)
510  {
511  if (!faceData[facei].valid(deltaCalc.data()))
512  {
513  if (!haveWarned)
514  {
516  << "Did not visit some faces, e.g. face " << facei
517  << " at " << mesh_.faceCentres()[facei]
518  << endl;
519  haveWarned = true;
520  }
521  }
522  else
523  {
524  faceToStart[facei] = faceData[facei].data();
525  }
526  }
527 
528  // Compact
529  List<Map<label>> compactMap;
530  mapPtr.reset(new mapDistribute(globalStart, faceToStart, compactMap));
531 }
532 
533 
534 void Foam::meshRefinement::nearestPatch
535 (
536  const labelList& adaptPatchIDs,
537  labelList& nearestPatch,
538  labelList& nearestZone
539 ) const
540 {
541  // Determine nearest patch for all mesh faces. Used when removing cells
542  // to give some reasonable patch to exposed faces.
543 
544  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
545 
546  nearestZone.setSize(mesh_.nFaces(), -1);
547 
548  if (adaptPatchIDs.size())
549  {
550  nearestPatch.setSize(mesh_.nFaces(), adaptPatchIDs[0]);
551 
552 
553  // Get per-face the zone or -1
554  labelList faceToZone(mesh_.nFaces(), -1);
555  {
556  for (const faceZone& fz : mesh_.faceZones())
557  {
558  UIndirectList<label>(faceToZone, fz) = fz.index();
559  }
560  }
561 
562 
563  // Count number of faces in adaptPatchIDs
564  label nFaces = 0;
565  forAll(adaptPatchIDs, i)
566  {
567  const polyPatch& pp = patches[adaptPatchIDs[i]];
568  nFaces += pp.size();
569  }
570 
571  // Field on cells and faces.
572  List<topoDistanceData<labelPair>> cellData(mesh_.nCells());
573  List<topoDistanceData<labelPair>> faceData(mesh_.nFaces());
574 
575  // Start of changes
576  labelList patchFaces(nFaces);
577  List<topoDistanceData<labelPair>> patchData(nFaces);
578  nFaces = 0;
579  forAll(adaptPatchIDs, i)
580  {
581  label patchi = adaptPatchIDs[i];
582  const polyPatch& pp = patches[patchi];
583 
584  forAll(pp, i)
585  {
586  patchFaces[nFaces] = pp.start()+i;
587  patchData[nFaces] = topoDistanceData<labelPair>
588  (
589  0,
590  labelPair
591  (
592  patchi,
593  faceToZone[pp.start()+i]
594  )
595  );
596  nFaces++;
597  }
598  }
599 
600  // Propagate information inwards
601  FaceCellWave<topoDistanceData<labelPair>> deltaCalc
602  (
603  mesh_,
604  patchFaces,
605  patchData,
606  faceData,
607  cellData,
608  mesh_.globalData().nTotalCells()+1
609  );
610 
611  // And extract
612 
613  bool haveWarned = false;
614  forAll(faceData, facei)
615  {
616  if (!faceData[facei].valid(deltaCalc.data()))
617  {
618  if (!haveWarned)
619  {
621  << "Did not visit some faces, e.g. face " << facei
622  << " at " << mesh_.faceCentres()[facei] << endl
623  << "Assigning these faces to patch "
624  << adaptPatchIDs[0]
625  << endl;
626  haveWarned = true;
627  }
628  }
629  else
630  {
631  const labelPair& data = faceData[facei].data();
632  nearestPatch[facei] = data.first();
633  nearestZone[facei] = data.second();
634  }
635  }
636  }
637  else
638  {
639  // Use patch 0
640  nearestPatch.setSize(mesh_.nFaces(), 0);
641  }
642 }
643 
644 
645 Foam::labelList Foam::meshRefinement::nearestPatch
646 (
647  const labelList& adaptPatchIDs
648 ) const
649 {
650  labelList nearestAdaptPatch;
651  labelList nearestAdaptZone;
652  nearestPatch(adaptPatchIDs, nearestAdaptPatch, nearestAdaptZone);
653  return nearestAdaptPatch;
654 }
655 
656 
657 void Foam::meshRefinement::nearestIntersection
658 (
659  const labelList& surfacesToTest,
660  const labelList& testFaces,
661 
662  labelList& surface1,
663  List<pointIndexHit>& hit1,
664  labelList& region1,
665  labelList& surface2,
666  List<pointIndexHit>& hit2,
667  labelList& region2
668 ) const
669 {
670  // Swap neighbouring cell centres and cell level
671  labelList neiLevel(mesh_.nBoundaryFaces());
672  pointField neiCc(mesh_.nBoundaryFaces());
673  calcNeighbourData(neiLevel, neiCc);
674 
675 
676  // Collect segments
677 
678  pointField start(testFaces.size());
679  pointField end(testFaces.size());
680  labelList minLevel(testFaces.size());
681 
682  calcCellCellRays
683  (
684  neiCc,
685  neiLevel,
686  testFaces,
687  start,
688  end,
689  minLevel
690  );
691 
692  // Do tests in one go
693 
694  surfaces_.findNearestIntersection
695  (
696  surfacesToTest,
697  start,
698  end,
699 
700  surface1,
701  hit1,
702  region1,
703  surface2,
704  hit2,
705  region2
706  );
707 }
708 
709 
710 Foam::labelList Foam::meshRefinement::nearestIntersection
711 (
712  const labelList& surfacesToTest,
713  const label defaultRegion
714 ) const
715 {
716  // Determine nearest intersection for all mesh faces. Used when removing
717  // cells to give some reasonable patch to exposed faces. Use this
718  // function instead of nearestPatch if you don't have patches yet.
719 
720 
721  // All faces with any intersection
722  const labelList testFaces(intersectedFaces());
723 
724  // Find intersection (if any)
725  labelList surface1;
726  List<pointIndexHit> hit1;
727  labelList region1;
728  labelList surface2;
729  List<pointIndexHit> hit2;
730  labelList region2;
731  nearestIntersection
732  (
733  surfacesToTest,
734  testFaces,
735 
736  surface1,
737  hit1,
738  region1,
739  surface2,
740  hit2,
741  region2
742  );
743 
744 
745  labelList nearestRegion(mesh_.nFaces(), defaultRegion);
746 
747  // Field on cells and faces.
748  List<topoDistanceData<label>> cellData(mesh_.nCells());
749  List<topoDistanceData<label>> faceData(mesh_.nFaces());
750 
751  // Start walking from all intersected faces
752  DynamicList<label> patchFaces(testFaces.size());
753  DynamicList<topoDistanceData<label>> patchData(testFaces.size());
754  forAll(testFaces, i)
755  {
756  label facei = testFaces[i];
757  if (surface1[i] != -1)
758  {
759  patchFaces.append(facei);
760  label regioni = surfaces_.globalRegion(surface1[i], region1[i]);
761  patchData.append(topoDistanceData<label>(0, regioni));
762  }
763  else if (surface2[i] != -1)
764  {
765  patchFaces.append(facei);
766  label regioni = surfaces_.globalRegion(surface2[i], region2[i]);
767  patchData.append(topoDistanceData<label>(0, regioni));
768  }
769  }
770 
771  // Propagate information inwards
772  FaceCellWave<topoDistanceData<label>> deltaCalc
773  (
774  mesh_,
775  patchFaces,
776  patchData,
777  faceData,
778  cellData,
779  mesh_.globalData().nTotalCells()+1
780  );
781 
782  // And extract
783 
784  bool haveWarned = false;
785  forAll(faceData, facei)
786  {
787  if (!faceData[facei].valid(deltaCalc.data()))
788  {
789  if (!haveWarned)
790  {
792  << "Did not visit some faces, e.g. face " << facei
793  << " at " << mesh_.faceCentres()[facei] << endl
794  << "Assigning these faces to global region "
795  << defaultRegion << endl;
796  haveWarned = true;
797  }
798  }
799  else
800  {
801  nearestRegion[facei] = faceData[facei].data();
802  }
803  }
804 
805  return nearestRegion;
806 }
807 
808 
810 (
811  const string& msg,
812  const polyMesh& mesh,
813  const List<scalar>& fld
814 )
815 {
816  if (fld.size() != mesh.nPoints())
817  {
819  << msg << endl
820  << "fld size:" << fld.size() << " mesh points:" << mesh.nPoints()
821  << abort(FatalError);
822  }
823 
824  Pout<< "Checking field " << msg << endl;
825  scalarField minFld(fld);
827  (
828  mesh,
829  minFld,
830  minEqOp<scalar>(),
831  GREAT
832  );
833  scalarField maxFld(fld);
835  (
836  mesh,
837  maxFld,
838  maxEqOp<scalar>(),
839  -GREAT
840  );
841  forAll(minFld, pointi)
842  {
843  const scalar& minVal = minFld[pointi];
844  const scalar& maxVal = maxFld[pointi];
845  if (mag(minVal-maxVal) > SMALL)
846  {
847  Pout<< msg << " at:" << mesh.points()[pointi] << nl
848  << " minFld:" << minVal << nl
849  << " maxFld:" << maxVal << nl
850  << endl;
851  }
852  }
853 }
854 
855 
857 (
858  const string& msg,
859  const polyMesh& mesh,
860  const List<point>& fld
861 )
862 {
863  if (fld.size() != mesh.nPoints())
864  {
866  << msg << endl
867  << "fld size:" << fld.size() << " mesh points:" << mesh.nPoints()
868  << abort(FatalError);
869  }
870 
871  Pout<< "Checking field " << msg << endl;
872  pointField minFld(fld);
874  (
875  mesh,
876  minFld,
878  point(GREAT, GREAT, GREAT)
879  );
880  pointField maxFld(fld);
882  (
883  mesh,
884  maxFld,
886  vector::zero
887  );
888  forAll(minFld, pointi)
889  {
890  const point& minVal = minFld[pointi];
891  const point& maxVal = maxFld[pointi];
892  if (mag(minVal-maxVal) > SMALL)
893  {
894  Pout<< msg << " at:" << mesh.points()[pointi] << nl
895  << " minFld:" << minVal << nl
896  << " maxFld:" << maxVal << nl
897  << endl;
898  }
899  }
900 }
901 
904 {
905  Pout<< "meshRefinement::checkData() : Checking refinement structure."
906  << endl;
907  meshCutter_.checkMesh();
908 
909  Pout<< "meshRefinement::checkData() : Checking refinement levels."
910  << endl;
911  meshCutter_.checkRefinementLevels(1, labelList(0));
912 
913 
914  const label nBnd = mesh_.nBoundaryFaces();
915 
916  Pout<< "meshRefinement::checkData() : Checking synchronization."
917  << endl;
918 
919  // Check face centres
920  {
921  // Boundary face centres
922  pointField::subList boundaryFc
923  (
924  mesh_.faceCentres(),
925  nBnd,
926  mesh_.nInternalFaces()
927  );
928 
929  // Get neighbouring face centres
930  pointField neiBoundaryFc(boundaryFc);
932  (
933  mesh_,
934  neiBoundaryFc,
935  eqOp<point>()
936  );
937 
938  // Compare
939  testSyncBoundaryFaceList
940  (
941  mergeDistance_,
942  "testing faceCentres : ",
943  boundaryFc,
944  neiBoundaryFc
945  );
946  }
947 
948  // Check meshRefinement
949  const labelList& surfIndex = surfaceIndex();
950 
951 
952  {
953  // Get boundary face centre and level. Coupled aware.
954  labelList neiLevel(nBnd);
955  pointField neiCc(nBnd);
956  calcNeighbourData(neiLevel, neiCc);
957 
958  // Collect segments we want to test for
959  pointField start(mesh_.nFaces());
960  pointField end(mesh_.nFaces());
961 
962  forAll(start, facei)
963  {
964  start[facei] = mesh_.cellCentres()[mesh_.faceOwner()[facei]];
965 
966  if (mesh_.isInternalFace(facei))
967  {
968  end[facei] = mesh_.cellCentres()[mesh_.faceNeighbour()[facei]];
969  }
970  else
971  {
972  end[facei] = neiCc[facei-mesh_.nInternalFaces()];
973  }
974  }
975 
976  // Extend segments a bit
977  {
978  const vectorField smallVec(ROOTSMALL*(end-start));
979  start -= smallVec;
980  end += smallVec;
981  }
982 
983 
984  // Do tests in one go
985  labelList surfaceHit;
986  {
987  labelList surfaceLevel;
988  surfaces_.findHigherIntersection
989  (
990  shells_,
991  start,
992  end,
993  labelList(start.size(), -1), // accept any intersection
994  surfaceHit,
995  surfaceLevel
996  );
997  }
998  // Get the coupled hit
999  labelList neiHit
1000  (
1002  (
1003  surfaceHit,
1004  nBnd,
1005  mesh_.nInternalFaces()
1006  )
1007  );
1008  syncTools::swapBoundaryFaceList(mesh_, neiHit);
1009 
1010  // Check
1011  forAll(surfaceHit, facei)
1012  {
1013  if (surfIndex[facei] != surfaceHit[facei])
1014  {
1015  if (mesh_.isInternalFace(facei))
1016  {
1018  << "Internal face:" << facei
1019  << " fc:" << mesh_.faceCentres()[facei]
1020  << " cached surfaceIndex_:" << surfIndex[facei]
1021  << " current:" << surfaceHit[facei]
1022  << " ownCc:"
1023  << mesh_.cellCentres()[mesh_.faceOwner()[facei]]
1024  << " neiCc:"
1025  << mesh_.cellCentres()[mesh_.faceNeighbour()[facei]]
1026  << endl;
1027  }
1028  else if
1029  (
1030  surfIndex[facei]
1031  != neiHit[facei-mesh_.nInternalFaces()]
1032  )
1033  {
1035  << "Boundary face:" << facei
1036  << " fc:" << mesh_.faceCentres()[facei]
1037  << " cached surfaceIndex_:" << surfIndex[facei]
1038  << " current:" << surfaceHit[facei]
1039  << " ownCc:"
1040  << mesh_.cellCentres()[mesh_.faceOwner()[facei]]
1041  << " neiCc:"
1042  << neiCc[facei-mesh_.nInternalFaces()]
1043  << " end:" << end[facei]
1044  << " ownLevel:"
1045  << meshCutter_.cellLevel()[mesh_.faceOwner()[facei]]
1046  << " faceLevel:"
1047  << meshCutter_.faceLevel(facei)
1048  << endl;
1049  }
1050  }
1051  }
1052  }
1053  {
1054  labelList::subList boundarySurface
1055  (
1056  surfaceIndex_,
1057  mesh_.nBoundaryFaces(),
1058  mesh_.nInternalFaces()
1059  );
1060 
1061  labelList neiBoundarySurface(boundarySurface);
1063  (
1064  mesh_,
1065  neiBoundarySurface
1066  );
1067 
1068  // Compare
1069  testSyncBoundaryFaceList
1070  (
1071  0, // tolerance
1072  "testing surfaceIndex() : ",
1073  boundarySurface,
1074  neiBoundarySurface
1075  );
1076  }
1077 
1078 
1079  // Find duplicate faces
1080  Pout<< "meshRefinement::checkData() : Counting duplicate faces."
1081  << endl;
1082 
1083  labelList duplicateFace
1084  (
1086  (
1087  mesh_,
1088  identity(mesh_.nBoundaryFaces(), mesh_.nInternalFaces())
1089  )
1090  );
1091 
1092  // Count
1093  {
1094  label nDup = 0;
1095 
1096  forAll(duplicateFace, i)
1097  {
1098  if (duplicateFace[i] != -1)
1099  {
1100  nDup++;
1101  }
1102  }
1103  nDup /= 2; // will have counted both faces of duplicate
1104  Pout<< "meshRefinement::checkData() : Found " << nDup
1105  << " duplicate pairs of faces." << endl;
1106  }
1107 }
1108 
1111 {
1112  meshCutter_.setInstance(inst);
1113  surfaceIndex_.instance() = inst;
1114 }
1115 
1116 
1118 (
1119  const labelList& cellsToRemove,
1120  const labelList& exposedFaces,
1121  const labelList& exposedPatchIDs,
1122  removeCells& cellRemover
1123 )
1124 {
1125  polyTopoChange meshMod(mesh_);
1126 
1127  // Arbitrary: put exposed faces into last patch.
1128  cellRemover.setRefinement
1129  (
1130  cellsToRemove,
1131  exposedFaces,
1132  exposedPatchIDs,
1133  meshMod
1134  );
1135 
1136  // Remove any unnecessary fields
1137  mesh_.clearOut();
1138  mesh_.moving(false);
1139 
1140  // Change the mesh (no inflation)
1141  autoPtr<mapPolyMesh> mapPtr = meshMod.changeMesh(mesh_, false, true);
1142  mapPolyMesh& map = *mapPtr;
1143 
1144  // Update fields
1145  mesh_.updateMesh(map);
1146 
1147  // Move mesh (since morphing might not do this)
1148  if (map.hasMotionPoints())
1149  {
1150  mesh_.movePoints(map.preMotionPoints());
1151  }
1152  else
1153  {
1154  // Delete mesh volumes. No other way to do this?
1155  mesh_.clearOut();
1156  }
1157 
1158  // Reset the instance for if in overwrite mode
1159  mesh_.setInstance(timeName());
1160  setInstance(mesh_.facesInstance());
1161 
1162  // Update local mesh data
1163  cellRemover.updateMesh(map);
1164 
1165  // Update intersections. Recalculate intersections for exposed faces.
1166  labelList newExposedFaces = renumber
1167  (
1168  map.reverseFaceMap(),
1169  exposedFaces
1170  );
1171 
1172  //Pout<< "removeCells : updating intersections for "
1173  // << newExposedFaces.size() << " newly exposed faces." << endl;
1174 
1175  updateMesh(map, newExposedFaces);
1176 
1177  return mapPtr;
1178 }
1179 
1180 
1182 (
1183  const labelList& splitFaces,
1184  const labelPairList& splits,
1185  //const List<Pair<point>>& splitPoints,
1186  polyTopoChange& meshMod
1187 ) const
1188 {
1189  forAll(splitFaces, i)
1190  {
1191  label facei = splitFaces[i];
1192  const face& f = mesh_.faces()[facei];
1193 
1194  // Split as start and end index in face
1195  const labelPair& split = splits[i];
1196 
1197  label nVerts = split[1]-split[0];
1198  if (nVerts < 0)
1199  {
1200  nVerts += f.size();
1201  }
1202  nVerts += 1;
1203 
1204 
1205  // Split into f0, f1
1206  face f0(nVerts);
1207 
1208  label fp = split[0];
1209  forAll(f0, i)
1210  {
1211  f0[i] = f[fp];
1212  fp = f.fcIndex(fp);
1213  }
1214 
1215  face f1(f.size()-f0.size()+2);
1216  fp = split[1];
1217  forAll(f1, i)
1218  {
1219  f1[i] = f[fp];
1220  fp = f.fcIndex(fp);
1221  }
1222 
1223 
1224  // Determine face properties
1225  label own = mesh_.faceOwner()[facei];
1226  label nei = -1;
1227  label patchi = -1;
1228  if (facei >= mesh_.nInternalFaces())
1229  {
1230  patchi = mesh_.boundaryMesh().whichPatch(facei);
1231  }
1232  else
1233  {
1234  nei = mesh_.faceNeighbour()[facei];
1235  }
1236 
1237  label zonei = mesh_.faceZones().whichZone(facei);
1238  bool zoneFlip = false;
1239  if (zonei != -1)
1240  {
1241  const faceZone& fz = mesh_.faceZones()[zonei];
1242  zoneFlip = fz.flipMap()[fz.whichFace(facei)];
1243  }
1244 
1245 
1246  if (debug)
1247  {
1248  Pout<< "face:" << facei << " verts:" << f
1249  << " split into f0:" << f0
1250  << " f1:" << f1 << endl;
1251  }
1252 
1253  // Change/add faces
1254  meshMod.modifyFace
1255  (
1256  f0, // modified face
1257  facei, // label of face
1258  own, // owner
1259  nei, // neighbour
1260  false, // face flip
1261  patchi, // patch for face
1262  zonei, // zone for face
1263  zoneFlip // face flip in zone
1264  );
1265 
1266  meshMod.addFace
1267  (
1268  f1, // modified face
1269  own, // owner
1270  nei, // neighbour
1271  -1, // master point
1272  -1, // master edge
1273  facei, // master face
1274  false, // face flip
1275  patchi, // patch for face
1276  zonei, // zone for face
1277  zoneFlip // face flip in zone
1278  );
1279 
1280 
1282  //meshMod.modifyPoint
1283  //(
1284  // f[split[0]],
1285  // splitPoints[i][0],
1286  // -1,
1287  // true
1288  //);
1289  //meshMod.modifyPoint
1290  //(
1291  // f[split[1]],
1292  // splitPoints[i][1],
1293  // -1,
1294  // true
1295  //);
1296  }
1297 }
1298 
1299 
1301 (
1302  const labelList& splitFaces,
1303  const labelPairList& splits,
1304  const dictionary& motionDict,
1305 
1306  labelList& duplicateFace,
1307  List<labelPair>& baffles
1308 )
1309 {
1310  label nSplit = returnReduce(splitFaces.size(), sumOp<label>());
1311  Info<< nl
1312  << "Splitting " << nSplit
1313  << " faces across diagonals" << "." << nl << endl;
1314 
1315  // To undo: store original faces
1316  faceList originalFaces(UIndirectList<face>(mesh_.faces(), splitFaces));
1317  labelPairList facePairs(splitFaces.size(), labelPair(-1, -1));
1318 
1319 
1320  {
1321  polyTopoChange meshMod(mesh_);
1322  meshMod.setCapacity
1323  (
1324  meshMod.points().size(),
1325  meshMod.faces().size()+splitFaces.size(),
1326  mesh_.nCells()
1327  );
1328 
1329  // Insert the mesh changes
1330  doSplitFaces(splitFaces, splits, meshMod);
1331 
1332  // Remove any unnecessary fields
1333  mesh_.clearOut();
1334  mesh_.moving(false);
1335 
1336  // Change the mesh (no inflation)
1337  autoPtr<mapPolyMesh> mapPtr = meshMod.changeMesh(mesh_, false, true);
1338  mapPolyMesh& map = *mapPtr;
1339 
1340  // Update fields
1341  mesh_.updateMesh(map);
1342 
1343  // Move mesh (since morphing might not do this)
1344  if (map.hasMotionPoints())
1345  {
1346  mesh_.movePoints(map.preMotionPoints());
1347  }
1348  else
1349  {
1350  mesh_.clearOut();
1351  }
1352 
1353  // Reset the instance for if in overwrite mode
1354  mesh_.setInstance(timeName());
1355  setInstance(mesh_.facesInstance());
1356 
1357 
1358  // Update local mesh data
1359  // ~~~~~~~~~~~~~~~~~~~~~~
1360 
1361  forAll(originalFaces, i)
1362  {
1363  inplaceRenumber(map.reversePointMap(), originalFaces[i]);
1364  }
1365 
1366  {
1367  Map<label> splitFaceToIndex(invertToMap(splitFaces));
1368 
1369  forAll(map.faceMap(), facei)
1370  {
1371  label oldFacei = map.faceMap()[facei];
1372 
1373  const auto oldFaceFnd = splitFaceToIndex.cfind(oldFacei);
1374 
1375  if (oldFaceFnd.good())
1376  {
1377  labelPair& twoFaces = facePairs[oldFaceFnd.val()];
1378  if (twoFaces[0] == -1)
1379  {
1380  twoFaces[0] = facei;
1381  }
1382  else if (twoFaces[1] == -1)
1383  {
1384  twoFaces[1] = facei;
1385  }
1386  else
1387  {
1389  << "problem: twoFaces:" << twoFaces
1390  << exit(FatalError);
1391  }
1392  }
1393  }
1394  }
1395 
1396 
1397  // Update baffle data
1398  // ~~~~~~~~~~~~~~~~~~
1399 
1400  if (duplicateFace.size())
1401  {
1403  (
1404  map.faceMap(),
1405  label(-1),
1406  duplicateFace
1407  );
1408  }
1409 
1410  const labelList& oldToNewFaces = map.reverseFaceMap();
1411  forAll(baffles, i)
1412  {
1413  labelPair& baffle = baffles[i];
1414  baffle.first() = oldToNewFaces[baffle.first()];
1415  baffle.second() = oldToNewFaces[baffle.second()];
1416 
1417  if (baffle.first() == -1 || baffle.second() == -1)
1418  {
1420  << "Removed baffle : faces:" << baffle
1421  << exit(FatalError);
1422  }
1423  }
1424 
1425 
1426  // Update intersections
1427  // ~~~~~~~~~~~~~~~~~~~~
1428 
1429  {
1430  DynamicList<label> changedFaces(facePairs.size());
1431  forAll(facePairs, i)
1432  {
1433  changedFaces.append(facePairs[i].first());
1434  changedFaces.append(facePairs[i].second());
1435  }
1436 
1437  // Update intersections on changed faces
1438  updateMesh(map, growFaceCellFace(changedFaces));
1439  }
1440  }
1441 
1442 
1443 
1444  // Undo loop
1445  // ~~~~~~~~~
1446  // Maintains two bits of information:
1447  // facePairs : two faces originating from the same face
1448  // originalFaces : original face in current vertices
1449 
1450 
1451  for (label iteration = 0; iteration < 100; iteration++)
1452  {
1453  Info<< nl
1454  << "Undo iteration " << iteration << nl
1455  << "----------------" << endl;
1456 
1457 
1458  // Check mesh for errors
1459  // ~~~~~~~~~~~~~~~~~~~~~
1460 
1461  faceSet errorFaces(mesh_, "errorFaces", mesh_.nBoundaryFaces());
1462  bool hasErrors = motionSmoother::checkMesh
1463  (
1464  false, // report
1465  mesh_,
1466  motionDict,
1467  errorFaces,
1468  dryRun_
1469  );
1470  if (!hasErrors)
1471  {
1472  break;
1473  }
1474 
1475  // Extend faces
1476  {
1477  const labelList grownFaces(growFaceCellFace(errorFaces));
1478  errorFaces.clear();
1479  errorFaces.insert(grownFaces);
1480  }
1481 
1482 
1483  // Merge face pairs
1484  // ~~~~~~~~~~~~~~~~
1485  // (if one of the faces is in the errorFaces set)
1486 
1487  polyTopoChange meshMod(mesh_);
1488 
1489  // Indices (in facePairs) of merged faces
1490  labelHashSet mergedIndices(facePairs.size());
1491  forAll(facePairs, index)
1492  {
1493  const labelPair& twoFaces = facePairs[index];
1494 
1495  if
1496  (
1497  errorFaces.found(twoFaces.first())
1498  || errorFaces.found(twoFaces.second())
1499  )
1500  {
1501  const face& originalFace = originalFaces[index];
1502 
1503 
1504  // Determine face properties
1505  label own = mesh_.faceOwner()[twoFaces[0]];
1506  label nei = -1;
1507  label patchi = -1;
1508  if (twoFaces[0] >= mesh_.nInternalFaces())
1509  {
1510  patchi = mesh_.boundaryMesh().whichPatch(twoFaces[0]);
1511  }
1512  else
1513  {
1514  nei = mesh_.faceNeighbour()[twoFaces[0]];
1515  }
1516 
1517  label zonei = mesh_.faceZones().whichZone(twoFaces[0]);
1518  bool zoneFlip = false;
1519  if (zonei != -1)
1520  {
1521  const faceZone& fz = mesh_.faceZones()[zonei];
1522  zoneFlip = fz.flipMap()[fz.whichFace(twoFaces[0])];
1523  }
1524 
1525  // Modify first face
1526  meshMod.modifyFace
1527  (
1528  originalFace, // modified face
1529  twoFaces[0], // label of face
1530  own, // owner
1531  nei, // neighbour
1532  false, // face flip
1533  patchi, // patch for face
1534  zonei, // zone for face
1535  zoneFlip // face flip in zone
1536  );
1537  // Merge second face into first
1538  meshMod.removeFace(twoFaces[1], twoFaces[0]);
1539 
1540  mergedIndices.insert(index);
1541  }
1542  }
1543 
1544  label n = returnReduce(mergedIndices.size(), sumOp<label>());
1545 
1546  Info<< "Detected " << n
1547  << " split faces that will be restored to their original faces."
1548  << nl << endl;
1549 
1550  if (n == 0)
1551  {
1552  // Nothing to be restored
1553  break;
1554  }
1555 
1556  nSplit -= n;
1557 
1558 
1559  // Remove any unnecessary fields
1560  mesh_.clearOut();
1561  mesh_.moving(false);
1562 
1563  // Change the mesh (no inflation)
1564  autoPtr<mapPolyMesh> mapPtr = meshMod.changeMesh(mesh_, false, true);
1565  mapPolyMesh& map = *mapPtr;
1566 
1567  // Update fields
1568  mesh_.updateMesh(map);
1569 
1570  // Move mesh (since morphing might not do this)
1571  if (map.hasMotionPoints())
1572  {
1573  mesh_.movePoints(map.preMotionPoints());
1574  }
1575  else
1576  {
1577  mesh_.clearOut();
1578  }
1579 
1580  // Reset the instance for if in overwrite mode
1581  mesh_.setInstance(timeName());
1582  setInstance(mesh_.facesInstance());
1583 
1584 
1585  // Update local mesh data
1586  // ~~~~~~~~~~~~~~~~~~~~~~
1587 
1588  {
1589  const labelList& oldToNewFaces = map.reverseFaceMap();
1590  const labelList& oldToNewPoints = map.reversePointMap();
1591 
1592  // Compact out merged faces
1593  DynamicList<label> changedFaces(mergedIndices.size());
1594 
1595  label newIndex = 0;
1596  forAll(facePairs, index)
1597  {
1598  const labelPair& oldSplit = facePairs[index];
1599  label new0 = oldToNewFaces[oldSplit[0]];
1600  label new1 = oldToNewFaces[oldSplit[1]];
1601 
1602  if (!mergedIndices.found(index))
1603  {
1604  // Faces still split
1605  if (new0 < 0 || new1 < 0)
1606  {
1608  << "Problem: oldFaces:" << oldSplit
1609  << " newFaces:" << labelPair(new0, new1)
1610  << exit(FatalError);
1611  }
1612 
1613  facePairs[newIndex] = labelPair(new0, new1);
1614  originalFaces[newIndex] = renumber
1615  (
1616  oldToNewPoints,
1617  originalFaces[index]
1618  );
1619  newIndex++;
1620  }
1621  else
1622  {
1623  // Merged face. Only new0 kept.
1624  if (new0 < 0 || new1 == -1)
1625  {
1627  << "Problem: oldFaces:" << oldSplit
1628  << " newFace:" << labelPair(new0, new1)
1629  << exit(FatalError);
1630  }
1631  changedFaces.append(new0);
1632  }
1633  }
1634 
1635  facePairs.setSize(newIndex);
1636  originalFaces.setSize(newIndex);
1637 
1638 
1639  // Update intersections
1640  updateMesh(map, growFaceCellFace(changedFaces));
1641  }
1642 
1643  // Update baffle data
1644  // ~~~~~~~~~~~~~~~~~~
1645  {
1646  if (duplicateFace.size())
1647  {
1649  (
1650  map.faceMap(),
1651  label(-1),
1652  duplicateFace
1653  );
1654  }
1655 
1656  const labelList& reverseFaceMap = map.reverseFaceMap();
1657  forAll(baffles, i)
1658  {
1659  labelPair& baffle = baffles[i];
1660  baffle.first() = reverseFaceMap[baffle.first()];
1661  baffle.second() = reverseFaceMap[baffle.second()];
1662 
1663  if (baffle.first() == -1 || baffle.second() == -1)
1664  {
1666  << "Removed baffle : faces:" << baffle
1667  << exit(FatalError);
1668  }
1669  }
1670  }
1671 
1672  }
1673 
1674  return nSplit;
1675 }
1676 
1677 
1678 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1679 
1680 Foam::meshRefinement::meshRefinement
1681 (
1682  fvMesh& mesh,
1683  const scalar mergeDistance,
1684  const bool overwrite,
1685  const refinementSurfaces& surfaces,
1686  const refinementFeatures& features,
1687  const shellSurfaces& shells,
1688  const shellSurfaces& limitShells,
1689  const labelUList& checkFaces,
1690  const bool dryRun
1691 )
1692 :
1693  mesh_(mesh),
1694  mergeDistance_(mergeDistance),
1695  overwrite_(overwrite),
1696  oldInstance_(mesh.pointsInstance()),
1697  surfaces_(surfaces),
1698  features_(features),
1699  shells_(shells),
1700  limitShells_(limitShells),
1701  dryRun_(dryRun),
1702  meshCutter_
1703  (
1704  mesh,
1705  false // do not try to read history.
1706  ),
1707  surfaceIndex_
1708  (
1709  IOobject
1710  (
1711  "surfaceIndex",
1712  mesh_.facesInstance(),
1713  polyMesh::meshSubDir,
1714  mesh_,
1715  IOobject::NO_READ,
1716  IOobject::NO_WRITE,
1717  IOobject::NO_REGISTER
1718  ),
1719  labelList(mesh_.nFaces(), -1)
1720  ),
1721  userFaceData_(0)
1722 {
1723  // recalculate intersections for all faces
1724  updateIntersections(checkFaces);
1725 }
1726 
1727 
1728 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1731 {
1732  if (surfaceIndex_.size() != mesh_.nFaces())
1733  {
1734  const_cast<meshRefinement&>(*this).updateIntersections
1735  (
1736  identity(mesh_.nFaces())
1737  );
1738  }
1739  return surfaceIndex_;
1740 }
1741 
1744 {
1745  if (surfaceIndex_.size() != mesh_.nFaces())
1746  {
1747  updateIntersections(identity(mesh_.nFaces()));
1748  }
1749  return surfaceIndex_;
1750 }
1751 
1753 Foam::label Foam::meshRefinement::countHits() const
1754 {
1755  // Stats on edges to test. Count proc faces only once.
1756  bitSet isMasterFace(syncTools::getMasterFaces(mesh_));
1757 
1758  label nHits = 0;
1759 
1760  const labelList& surfIndex = surfaceIndex();
1761 
1762  forAll(surfIndex, facei)
1763  {
1764  if (surfIndex[facei] >= 0 && isMasterFace.test(facei))
1765  {
1766  ++nHits;
1767  }
1768  }
1769  return nHits;
1770 }
1771 
1772 
1774 (
1775  const bool keepZoneFaces,
1776  const bool keepBaffles,
1777  const scalarField& cellWeights,
1778  decompositionMethod& decomposer,
1779  fvMeshDistribute& distributor
1780 )
1781 {
1783 
1784  if (Pstream::parRun())
1785  {
1786  // Wanted distribution
1788 
1789 
1790  // Faces where owner and neighbour are not 'connected' so can
1791  // go to different processors.
1792  boolList blockedFace;
1793  label nUnblocked = 0;
1794 
1795  // Faces that move as block onto single processor
1796  PtrList<labelList> specifiedProcessorFaces;
1797  labelList specifiedProcessor;
1798 
1799  // Pairs of baffles
1800  List<labelPair> couples;
1801 
1802  // Constraints from decomposeParDict
1803  decomposer.setConstraints
1804  (
1805  mesh_,
1806  blockedFace,
1807  specifiedProcessorFaces,
1808  specifiedProcessor,
1809  couples
1810  );
1811 
1812 
1813  if (keepZoneFaces || keepBaffles)
1814  {
1815  if (keepZoneFaces)
1816  {
1817  // Determine decomposition to keep/move surface zones
1818  // on one processor. The reason is that snapping will make these
1819  // into baffles, move and convert them back so if they were
1820  // proc boundaries after baffling&moving the points might be no
1821  // longer synchronised so recoupling will fail. To prevent this
1822  // keep owner&neighbour of such a surface zone on the same
1823  // processor.
1824 
1825  const faceZoneMesh& fZones = mesh_.faceZones();
1826  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
1827 
1828  // Get faces whose owner and neighbour should stay together,
1829  // i.e. they are not 'blocked'.
1830 
1831  forAll(fZones, zonei)
1832  {
1833  const faceZone& fZone = fZones[zonei];
1834 
1835  forAll(fZone, i)
1836  {
1837  label facei = fZone[i];
1838  if (blockedFace[facei])
1839  {
1840  if
1841  (
1842  mesh_.isInternalFace(facei)
1843  || pbm[pbm.whichPatch(facei)].coupled()
1844  )
1845  {
1846  blockedFace[facei] = false;
1847  nUnblocked++;
1848  }
1849  }
1850  }
1851  }
1852 
1853 
1854  // If the faceZones are not synchronised the blockedFace
1855  // might not be synchronised. If you are sure the faceZones
1856  // are synchronised remove below check.
1858  (
1859  mesh_,
1860  blockedFace,
1861  andEqOp<bool>() // combine operator
1862  );
1863  }
1864  reduce(nUnblocked, sumOp<label>());
1865 
1866  if (keepZoneFaces)
1867  {
1868  Info<< "Found " << nUnblocked
1869  << " zoned faces to keep together." << endl;
1870  }
1871 
1872 
1873  // Extend un-blockedFaces with any cyclics
1874  {
1875  boolList separatedCoupledFace(mesh_.nFaces(), false);
1876  selectSeparatedCoupledFaces(separatedCoupledFace);
1877 
1878  label nSeparated = 0;
1879  forAll(separatedCoupledFace, facei)
1880  {
1881  if (separatedCoupledFace[facei])
1882  {
1883  if (blockedFace[facei])
1884  {
1885  blockedFace[facei] = false;
1886  nSeparated++;
1887  }
1888  }
1889  }
1890  reduce(nSeparated, sumOp<label>());
1891  Info<< "Found " << nSeparated
1892  << " separated coupled faces to keep together." << endl;
1893 
1894  nUnblocked += nSeparated;
1895  }
1896 
1897 
1898  if (keepBaffles)
1899  {
1900  const label nBnd = mesh_.nBoundaryFaces();
1901 
1902  labelList coupledFace(mesh_.nFaces(), -1);
1903  {
1904  // Get boundary baffles that need to stay together
1905  List<labelPair> allCouples
1906  (
1908  );
1909 
1910  // Merge with any couples from
1911  // decompositionMethod::setConstraints
1912  forAll(couples, i)
1913  {
1914  const labelPair& baffle = couples[i];
1915  coupledFace[baffle.first()] = baffle.second();
1916  coupledFace[baffle.second()] = baffle.first();
1917  }
1918  forAll(allCouples, i)
1919  {
1920  const labelPair& baffle = allCouples[i];
1921  coupledFace[baffle.first()] = baffle.second();
1922  coupledFace[baffle.second()] = baffle.first();
1923  }
1924  }
1925 
1926  couples.setSize(nBnd);
1927  label nCpl = 0;
1928  forAll(coupledFace, facei)
1929  {
1930  if (coupledFace[facei] != -1 && facei < coupledFace[facei])
1931  {
1932  couples[nCpl++] = labelPair(facei, coupledFace[facei]);
1933  }
1934  }
1935  couples.setSize(nCpl);
1936  }
1937  label nCouples = returnReduce(couples.size(), sumOp<label>());
1938 
1939  if (keepBaffles)
1940  {
1941  Info<< "Found " << nCouples << " baffles to keep together."
1942  << endl;
1943  }
1944  }
1945 
1946 
1947  // Make sure blockedFace not set on couples
1948  forAll(couples, i)
1949  {
1950  const labelPair& baffle = couples[i];
1951  blockedFace[baffle.first()] = false;
1952  blockedFace[baffle.second()] = false;
1953  }
1954 
1955  distribution = decomposer.decompose
1956  (
1957  mesh_,
1958  cellWeights,
1959  blockedFace,
1960  specifiedProcessorFaces,
1961  specifiedProcessor,
1962  couples // explicit connections
1963  );
1964 
1965  if (debug)
1966  {
1967  labelList nProcCells(distributor.countCells(distribution));
1968  Pout<< "Wanted distribution:" << nProcCells << endl;
1969 
1971 
1972  Pout<< "Wanted resulting decomposition:" << endl;
1973  forAll(nProcCells, proci)
1974  {
1975  Pout<< " " << proci << '\t' << nProcCells[proci] << endl;
1976  }
1977  Pout<< endl;
1978  }
1979 
1980  // Do actual sending/receiving of mesh
1981  map = distributor.distribute(distribution);
1982 
1983  // Update numbering of meshRefiner
1984  distribute(map());
1985 
1986  // Set correct instance (for if overwrite)
1987  mesh_.setInstance(timeName());
1988  setInstance(mesh_.facesInstance());
1989 
1990 
1991  if (debug && keepZoneFaces)
1992  {
1993  const faceZoneMesh& fZones = mesh_.faceZones();
1994  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
1995 
1996  // Check that faceZone faces are indeed internal
1997  forAll(fZones, zonei)
1998  {
1999  const faceZone& fZone = fZones[zonei];
2000 
2001  forAll(fZone, i)
2002  {
2003  label facei = fZone[i];
2004  label patchi = pbm.whichPatch(facei);
2005 
2006  if (patchi >= 0 && pbm[patchi].coupled())
2007  {
2009  << "Face at " << mesh_.faceCentres()[facei]
2010  << " on zone " << fZone.name()
2011  << " is on coupled patch " << pbm[patchi].name()
2012  << endl;
2013  }
2014  }
2015  }
2016  }
2017  }
2018 
2019  return map;
2020 }
2021 
2024 {
2025  label nBoundaryFaces = 0;
2026 
2027  const labelList& surfIndex = surfaceIndex();
2028 
2029  forAll(surfIndex, facei)
2030  {
2031  if (surfIndex[facei] != -1)
2032  {
2033  nBoundaryFaces++;
2034  }
2035  }
2036 
2037  labelList surfaceFaces(nBoundaryFaces);
2038  nBoundaryFaces = 0;
2039 
2040  forAll(surfIndex, facei)
2041  {
2042  if (surfIndex[facei] != -1)
2043  {
2044  surfaceFaces[nBoundaryFaces++] = facei;
2045  }
2046  }
2047  return surfaceFaces;
2048 }
2049 
2052 {
2053  const faceList& faces = mesh_.faces();
2054 
2055  // Mark all points on faces that will become baffles
2056  bitSet isBoundaryPoint(mesh_.nPoints());
2057  label nBoundaryPoints = 0;
2058 
2059  const labelList& surfIndex = surfaceIndex();
2060 
2061  forAll(surfIndex, facei)
2062  {
2063  if (surfIndex[facei] != -1)
2064  {
2065  const face& f = faces[facei];
2066 
2067  forAll(f, fp)
2068  {
2069  if (isBoundaryPoint.set(f[fp]))
2070  {
2071  nBoundaryPoints++;
2072  }
2073  }
2074  }
2075  }
2076 
2078  //labelList adaptPatchIDs(meshedPatches());
2079  //forAll(adaptPatchIDs, i)
2080  //{
2081  // label patchi = adaptPatchIDs[i];
2082  //
2083  // if (patchi != -1)
2084  // {
2085  // const polyPatch& pp = mesh_.boundaryMesh()[patchi];
2086  //
2087  // label facei = pp.start();
2088  //
2089  // forAll(pp, i)
2090  // {
2091  // const face& f = faces[facei];
2092  //
2093  // forAll(f, fp)
2094  // {
2095  // if (isBoundaryPoint.set(f[fp]))
2096  // nBoundaryPoints++;
2097  // }
2098  // }
2099  // facei++;
2100  // }
2101  // }
2102  //}
2103 
2104 
2105  // Pack
2106  labelList boundaryPoints(nBoundaryPoints);
2107  nBoundaryPoints = 0;
2108  forAll(isBoundaryPoint, pointi)
2109  {
2110  if (isBoundaryPoint.test(pointi))
2111  {
2112  boundaryPoints[nBoundaryPoints++] = pointi;
2113  }
2114  }
2115 
2116  return boundaryPoints;
2117 }
2118 
2119 
2120 //- Create patch from set of patches
2122 (
2123  const polyMesh& mesh,
2124  const labelList& patchIDs
2125 )
2126 {
2128 
2129  // Count faces.
2130  label nFaces = 0;
2131 
2132  forAll(patchIDs, i)
2133  {
2134  const polyPatch& pp = patches[patchIDs[i]];
2135 
2136  nFaces += pp.size();
2137  }
2138 
2139  // Collect faces.
2140  labelList addressing(nFaces);
2141  nFaces = 0;
2142 
2143  forAll(patchIDs, i)
2144  {
2145  const polyPatch& pp = patches[patchIDs[i]];
2146 
2147  label meshFacei = pp.start();
2148 
2149  forAll(pp, i)
2150  {
2151  addressing[nFaces++] = meshFacei++;
2152  }
2153  }
2154 
2156  (
2157  IndirectList<face>(mesh.faces(), addressing),
2158  mesh.points()
2159  );
2160 }
2161 
2162 
2164 (
2165  const pointMesh& pMesh,
2166  const labelList& adaptPatchIDs
2167 )
2168 {
2169  const polyMesh& mesh = pMesh();
2170 
2171  // Construct displacement field.
2172  const pointBoundaryMesh& pointPatches = pMesh.boundary();
2173 
2174  wordList patchFieldTypes
2175  (
2176  pointPatches.size(),
2177  slipPointPatchVectorField::typeName
2178  );
2179 
2180  forAll(adaptPatchIDs, i)
2181  {
2182  patchFieldTypes[adaptPatchIDs[i]] =
2183  fixedValuePointPatchVectorField::typeName;
2184  }
2185 
2186  forAll(pointPatches, patchi)
2187  {
2188  if (isA<processorPointPatch>(pointPatches[patchi]))
2189  {
2190  patchFieldTypes[patchi] = calculatedPointPatchVectorField::typeName;
2191  }
2192  else if (isA<cyclicPointPatch>(pointPatches[patchi]))
2193  {
2194  patchFieldTypes[patchi] = cyclicSlipPointPatchVectorField::typeName;
2195  }
2196  }
2197 
2198  // Note: time().timeName() instead of meshRefinement::timeName() since
2199  // postprocessable field.
2201  (
2202  IOobject
2203  (
2204  "pointDisplacement",
2205  mesh.time().timeName(),
2206  mesh,
2209  ),
2210  pMesh,
2212  patchFieldTypes
2213  );
2214 }
2215 
2216 
2219  const faceZoneMesh& fZones = mesh.faceZones();
2220 
2221  // Check any zones are present anywhere and in same order
2222 
2223  {
2224  List<wordList> zoneNames(Pstream::nProcs());
2225  zoneNames[Pstream::myProcNo()] = fZones.names();
2226  Pstream::allGatherList(zoneNames);
2227 
2228  // All have same data now. Check.
2229  forAll(zoneNames, proci)
2230  {
2231  if
2232  (
2233  proci != Pstream::myProcNo()
2234  && (zoneNames[proci] != zoneNames[Pstream::myProcNo()])
2235  )
2236  {
2238  << "faceZones are not synchronised on processors." << nl
2239  << "Processor " << proci << " has faceZones "
2240  << zoneNames[proci] << nl
2241  << "Processor " << Pstream::myProcNo()
2242  << " has faceZones "
2243  << zoneNames[Pstream::myProcNo()] << nl
2244  << exit(FatalError);
2245  }
2246  }
2247  }
2248 
2249  // Check that coupled faces are present on both sides.
2250 
2251  labelList faceToZone(mesh.nBoundaryFaces(), -1);
2252 
2253  forAll(fZones, zonei)
2254  {
2255  const faceZone& fZone = fZones[zonei];
2256 
2257  forAll(fZone, i)
2258  {
2259  label bFacei = fZone[i]-mesh.nInternalFaces();
2260 
2261  if (bFacei >= 0)
2262  {
2263  if (faceToZone[bFacei] == -1)
2264  {
2265  faceToZone[bFacei] = zonei;
2266  }
2267  else if (faceToZone[bFacei] == zonei)
2268  {
2270  << "Face " << fZone[i] << " in zone "
2271  << fZone.name()
2272  << " is twice in zone!"
2273  << abort(FatalError);
2274  }
2275  else
2276  {
2278  << "Face " << fZone[i] << " in zone "
2279  << fZone.name()
2280  << " is also in zone "
2281  << fZones[faceToZone[bFacei]].name()
2282  << abort(FatalError);
2283  }
2284  }
2285  }
2286  }
2287 
2288  labelList neiFaceToZone(faceToZone);
2289  syncTools::swapBoundaryFaceList(mesh, neiFaceToZone);
2290 
2291  forAll(faceToZone, i)
2292  {
2293  if (faceToZone[i] != neiFaceToZone[i])
2294  {
2296  << "Face " << mesh.nInternalFaces()+i
2297  << " is in zone " << faceToZone[i]
2298  << ", its coupled face is in zone " << neiFaceToZone[i]
2299  << abort(FatalError);
2300  }
2301  }
2302 }
2303 
2304 
2306 (
2307  const polyMesh& mesh,
2308  const bitSet& isMasterEdge,
2309  const labelList& meshPoints,
2310  const edgeList& edges,
2311  scalarField& edgeWeights,
2312  scalarField& invSumWeight
2313 )
2314 {
2315  const pointField& pts = mesh.points();
2316 
2317  // Calculate edgeWeights and inverse sum of edge weights
2318  edgeWeights.setSize(isMasterEdge.size());
2319  invSumWeight.setSize(meshPoints.size());
2320 
2321  forAll(edges, edgei)
2322  {
2323  const edge& e = edges[edgei];
2324  scalar eMag = max
2325  (
2326  SMALL,
2327  mag
2328  (
2329  pts[meshPoints[e[1]]]
2330  - pts[meshPoints[e[0]]]
2331  )
2332  );
2333  edgeWeights[edgei] = 1.0/eMag;
2334  }
2335 
2336  // Sum per point all edge weights
2337  weightedSum
2338  (
2339  mesh,
2340  isMasterEdge,
2341  meshPoints,
2342  edges,
2343  edgeWeights,
2344  scalarField(meshPoints.size(), 1.0), // data
2345  invSumWeight
2346  );
2347 
2348  // Inplace invert
2349  forAll(invSumWeight, pointi)
2350  {
2351  scalar w = invSumWeight[pointi];
2352 
2353  if (w > 0.0)
2354  {
2355  invSumWeight[pointi] = 1.0/w;
2356  }
2357  }
2358 }
2359 
2360 
2362 (
2364  const label insertPatchi,
2365  const word& patchName,
2366  const dictionary& patchDict
2367 )
2368 {
2369  // Clear local fields and e.g. polyMesh parallelInfo.
2370  mesh.clearOut();
2371 
2372  polyBoundaryMesh& polyPatches =
2373  const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
2374  fvBoundaryMesh& fvPatches = const_cast<fvBoundaryMesh&>(mesh.boundary());
2375 
2376  // The new location (at the end)
2377  const label patchi = polyPatches.size();
2378 
2379  // Add polyPatch at the end
2380  polyPatches.push_back
2381  (
2383  (
2384  patchName,
2385  patchDict,
2386  insertPatchi,
2387  polyPatches
2388  )
2389  );
2390 
2391  fvPatches.push_back
2392  (
2393  fvPatch::New
2394  (
2395  polyPatches[patchi], // point to newly added polyPatch
2396  mesh.boundary()
2397  )
2398  );
2399 
2400  addPatchFields<volScalarField>
2401  (
2402  mesh,
2404  );
2405  addPatchFields<volVectorField>
2406  (
2407  mesh,
2409  );
2410  addPatchFields<volSphericalTensorField>
2411  (
2412  mesh,
2414  );
2415  addPatchFields<volSymmTensorField>
2416  (
2417  mesh,
2419  );
2420  addPatchFields<volTensorField>
2421  (
2422  mesh,
2424  );
2425 
2426  // Surface fields
2427 
2428  addPatchFields<surfaceScalarField>
2429  (
2430  mesh,
2432  );
2433  addPatchFields<surfaceVectorField>
2434  (
2435  mesh,
2437  );
2438  addPatchFields<surfaceSphericalTensorField>
2439  (
2440  mesh,
2442  );
2443  addPatchFields<surfaceSymmTensorField>
2444  (
2445  mesh,
2447  );
2448  addPatchFields<surfaceTensorField>
2449  (
2450  mesh,
2452  );
2453  return patchi;
2454 }
2455 
2456 
2458 (
2460  const word& patchName,
2461  const dictionary& patchInfo
2462 )
2463 {
2464  polyBoundaryMesh& polyPatches =
2465  const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
2466  fvBoundaryMesh& fvPatches = const_cast<fvBoundaryMesh&>(mesh.boundary());
2467 
2468  const label patchi = polyPatches.findPatchID(patchName);
2469  if (patchi != -1)
2470  {
2471  // Already there
2472  return patchi;
2473  }
2474 
2475 
2476  label insertPatchi = polyPatches.size();
2477  label startFacei = mesh.nFaces();
2478 
2479  forAll(polyPatches, patchi)
2480  {
2481  const polyPatch& pp = polyPatches[patchi];
2482 
2483  if (isA<processorPolyPatch>(pp))
2484  {
2485  insertPatchi = patchi;
2486  startFacei = pp.start();
2487  break;
2488  }
2489  }
2490 
2491  dictionary patchDict(patchInfo);
2492  patchDict.set("nFaces", 0);
2493  patchDict.set("startFace", startFacei);
2494 
2495  // Below is all quite a hack. Feel free to change once there is a better
2496  // mechanism to insert and reorder patches.
2497 
2498  label addedPatchi = appendPatch(mesh, insertPatchi, patchName, patchDict);
2499 
2500 
2501  // Create reordering list
2502  // patches before insert position stay as is
2503  labelList oldToNew(addedPatchi+1);
2504  for (label i = 0; i < insertPatchi; i++)
2505  {
2506  oldToNew[i] = i;
2507  }
2508  // patches after insert position move one up
2509  for (label i = insertPatchi; i < addedPatchi; i++)
2510  {
2511  oldToNew[i] = i+1;
2512  }
2513  // appended patch gets moved to insert position
2514  oldToNew[addedPatchi] = insertPatchi;
2515 
2516  // Shuffle into place
2517  polyPatches.reorder(oldToNew, true);
2518  fvPatches.reorder(oldToNew);
2519 
2520  reorderPatchFields<volScalarField>(mesh, oldToNew);
2521  reorderPatchFields<volVectorField>(mesh, oldToNew);
2522  reorderPatchFields<volSphericalTensorField>(mesh, oldToNew);
2523  reorderPatchFields<volSymmTensorField>(mesh, oldToNew);
2524  reorderPatchFields<volTensorField>(mesh, oldToNew);
2525  reorderPatchFields<surfaceScalarField>(mesh, oldToNew);
2526  reorderPatchFields<surfaceVectorField>(mesh, oldToNew);
2527  reorderPatchFields<surfaceSphericalTensorField>(mesh, oldToNew);
2528  reorderPatchFields<surfaceSymmTensorField>(mesh, oldToNew);
2529  reorderPatchFields<surfaceTensorField>(mesh, oldToNew);
2530 
2531  return insertPatchi;
2532 }
2533 
2534 
2536 (
2537  const word& name,
2538  const dictionary& patchInfo
2539 )
2540 {
2541  label meshedi = meshedPatches_.find(name);
2542 
2543  if (meshedi != -1)
2544  {
2545  // Already there. Get corresponding polypatch
2546  return mesh_.boundaryMesh().findPatchID(name);
2547  }
2548  else
2549  {
2550  // Add patch
2551  label patchi = addPatch(mesh_, name, patchInfo);
2552 
2553 // dictionary patchDict(patchInfo);
2554 // patchDict.set("nFaces", 0);
2555 // patchDict.set("startFace", 0);
2556 // autoPtr<polyPatch> ppPtr
2557 // (
2558 // polyPatch::New
2559 // (
2560 // name,
2561 // patchDict,
2562 // 0,
2563 // mesh_.boundaryMesh()
2564 // )
2565 // );
2566 // label patchi = fvMeshTools::addPatch
2567 // (
2568 // mesh_,
2569 // ppPtr(),
2570 // dictionary(), // optional field values
2571 // fvPatchFieldBase::calculatedType(),
2572 // true
2573 // );
2574 
2575  // Store
2576  meshedPatches_.append(name);
2577 
2578  // Clear patch based addressing
2579  faceToCoupledPatch_.clear();
2580 
2581  return patchi;
2582  }
2583 }
2584 
2585 
2588  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
2589 
2590  DynamicList<label> patchIDs(meshedPatches_.size());
2591  forAll(meshedPatches_, i)
2592  {
2593  label patchi = patches.findPatchID(meshedPatches_[i]);
2594 
2595  if (patchi == -1)
2596  {
2598  << "Problem : did not find patch " << meshedPatches_[i]
2599  << endl << "Valid patches are " << patches.names()
2600  << endl; //abort(FatalError);
2601  }
2602  else if (!polyPatch::constraintType(patches[patchi].type()))
2603  {
2604  patchIDs.append(patchi);
2605  }
2606  }
2607 
2608  return patchIDs;
2609 }
2610 
2611 
2613 (
2614  const word& fzName,
2615  const word& masterPatch,
2616  const word& slavePatch,
2617  const surfaceZonesInfo::faceZoneType& fzType
2618 )
2619 {
2620  label zonei = surfaceZonesInfo::addFaceZone
2621  (
2622  fzName, //name
2623  labelList(0), //addressing
2624  boolList(0), //flipmap
2625  mesh_
2626  );
2627 
2628  faceZoneToMasterPatch_.insert(fzName, masterPatch);
2629  faceZoneToSlavePatch_.insert(fzName, slavePatch);
2630  faceZoneToType_.insert(fzName, fzType);
2631 
2632  return zonei;
2633 }
2634 
2635 
2637 (
2638  const word& fzName,
2639  label& masterPatchID,
2640  label& slavePatchID,
2642 ) const
2643 {
2644  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
2645 
2646  if (!faceZoneToMasterPatch_.found(fzName))
2647  {
2648  return false;
2649  }
2650  else
2651  {
2652  const word& masterName = faceZoneToMasterPatch_[fzName];
2653  masterPatchID = pbm.findPatchID(masterName);
2654 
2655  const word& slaveName = faceZoneToSlavePatch_[fzName];
2656  slavePatchID = pbm.findPatchID(slaveName);
2657 
2658  fzType = faceZoneToType_[fzName];
2659 
2660  return true;
2661  }
2662 }
2663 
2664 
2665 Foam::label Foam::meshRefinement::addPointZone(const word& name)
2667  pointZoneMesh& pointZones = mesh_.pointZones();
2668 
2669  label zoneI = pointZones.findZoneID(name);
2670 
2671  if (zoneI == -1)
2672  {
2673  zoneI = pointZones.size();
2674  pointZones.clearAddressing();
2675 
2676  pointZones.emplace_back
2677  (
2678  name, // name
2679  zoneI, // index
2680  pointZones // pointZoneMesh
2681  );
2682  }
2683  return zoneI;
2684 }
2685 
2686 
2689  for (const polyPatch& pp : mesh_.boundaryMesh())
2690  {
2691  // Check all coupled. Avoid using .coupled() so we also pick up AMI.
2692  const auto* cpp = isA<coupledPolyPatch>(pp);
2693 
2694  if (cpp && (cpp->separated() || !cpp->parallel()))
2695  {
2696  SubList<bool>(selected, pp.size(), pp.start()) = true;
2697  }
2698  }
2699 }
2700 
2701 
2703 (
2705 ) const
2706 {
2707  // Count number of faces per edge. Parallel consistent.
2708 
2709  const labelListList& edgeFaces = pp.edgeFaces();
2710  labelList nEdgeFaces(edgeFaces.size());
2711  forAll(edgeFaces, edgei)
2712  {
2713  nEdgeFaces[edgei] = edgeFaces[edgei].size();
2714  }
2715 
2716  // Sync across processor patches
2717  if (Pstream::parRun())
2718  {
2719  const globalMeshData& globalData = mesh_.globalData();
2720  const mapDistribute& map = globalData.globalEdgeSlavesMap();
2721  const indirectPrimitivePatch& cpp = globalData.coupledPatch();
2722 
2723  // Match pp edges to coupled edges
2724  labelList patchEdges;
2725  labelList coupledEdges;
2726  PackedBoolList sameEdgeOrientation;
2728  (
2729  pp,
2730  cpp,
2731  patchEdges,
2732  coupledEdges,
2733  sameEdgeOrientation
2734  );
2735 
2736  // Convert patch-edge data into cpp-edge data
2737  labelList coupledNEdgeFaces(map.constructSize(), Zero);
2738  UIndirectList<label>(coupledNEdgeFaces, coupledEdges) =
2739  UIndirectList<label>(nEdgeFaces, patchEdges);
2740 
2741  // Synchronise
2742  globalData.syncData
2743  (
2744  coupledNEdgeFaces,
2745  globalData.globalEdgeSlaves(),
2746  globalData.globalEdgeTransformedSlaves(),
2747  map,
2748  plusEqOp<label>()
2749  );
2750 
2751  // Convert back from cpp-edge to patch-edge
2752  UIndirectList<label>(nEdgeFaces, patchEdges) =
2753  UIndirectList<label>(coupledNEdgeFaces, coupledEdges);
2754  }
2755  return nEdgeFaces;
2756 }
2757 
2758 
2760 (
2761  const polyMesh& mesh,
2762  const vector& perturbVec,
2763  const point& p
2764 )
2765 {
2766  // Force calculation of base points (needs to be synchronised)
2767  (void)mesh.tetBasePtIs();
2768 
2769  label celli = mesh.findCell(p, findCellMode);
2770 
2771  if (returnReduceAnd(celli < 0))
2772  {
2773  // See if we can perturb a bit
2774  celli = mesh.findCell(p+perturbVec, findCellMode);
2775  }
2776  return celli;
2777 }
2778 
2779 
2781 (
2782  const polyMesh& mesh,
2783  const labelList& cellToRegion,
2784  const vector& perturbVec,
2785  const point& p
2786 )
2787 {
2788  label regioni = -1;
2789 
2790  // Force calculation of base points (needs to be synchronised)
2791  (void)mesh.tetBasePtIs();
2792 
2793  label celli = mesh.findCell(p, findCellMode);
2794  if (celli != -1)
2795  {
2796  regioni = cellToRegion[celli];
2797  }
2798  reduce(regioni, maxOp<label>());
2799 
2800  if (regioni == -1)
2801  {
2802  // See if we can perturb a bit
2803  celli = mesh.findCell(p+perturbVec, findCellMode);
2804  if (celli != -1)
2805  {
2806  regioni = cellToRegion[celli];
2807  }
2808  reduce(regioni, maxOp<label>());
2809  }
2810  return regioni;
2811 }
2812 
2813 
2814 Foam::fileName Foam::meshRefinement::writeLeakPath
2815 (
2816  const polyMesh& mesh,
2817  const pointField& locationsInMesh,
2818  const pointField& locationsOutsideMesh,
2819  const boolList& blockedFace,
2821 )
2822 {
2824 
2825  fileName outputDir
2826  (
2827  mesh.time().globalPath()
2829  / mesh.pointsInstance()
2830  );
2831  outputDir.clean(); // Remove unneeded ".."
2832 
2833  // Write the leak path
2834 
2835  meshSearch searchEngine(mesh);
2836  shortestPathSet leakPath
2837  (
2838  "leakPath",
2839  mesh,
2840  searchEngine,
2842  false, //true,
2843  50, // tbd. Number of iterations
2844  pbm.groupPatchIDs()["wall"],
2845  locationsInMesh,
2846  locationsOutsideMesh,
2847  blockedFace
2848  );
2849 
2850  // Split leak path according to segment. Note: segment index
2851  // is global (= index in locationsInsideMesh)
2852  List<pointList> segmentPoints;
2853  List<scalarList> segmentDist;
2854  {
2855  label nSegments = 0;
2856  if (leakPath.segments().size())
2857  {
2858  nSegments = max(leakPath.segments())+1;
2859  }
2860  reduce(nSegments, maxOp<label>());
2861 
2862  labelList nElemsPerSegment(nSegments, Zero);
2863  for (label segmenti : leakPath.segments())
2864  {
2865  nElemsPerSegment[segmenti]++;
2866  }
2867  segmentPoints.setSize(nElemsPerSegment.size());
2868  segmentDist.setSize(nElemsPerSegment.size());
2869  forAll(nElemsPerSegment, i)
2870  {
2871  segmentPoints[i].setSize(nElemsPerSegment[i]);
2872  segmentDist[i].setSize(nElemsPerSegment[i]);
2873  }
2874  nElemsPerSegment = 0;
2875 
2876  forAll(leakPath, elemi)
2877  {
2878  label segmenti = leakPath.segments()[elemi];
2879  pointList& points = segmentPoints[segmenti];
2880  scalarList& dist = segmentDist[segmenti];
2881  label& n = nElemsPerSegment[segmenti];
2882 
2883  points[n] = leakPath[elemi];
2884  dist[n] = leakPath.distance()[elemi];
2885  n++;
2886  }
2887  }
2888 
2889  PtrList<coordSet> allLeakPaths(segmentPoints.size());
2890  forAll(allLeakPaths, segmenti)
2891  {
2892  // Collect data from all processors
2893  List<pointList> gatheredPts(Pstream::nProcs());
2894  gatheredPts[Pstream::myProcNo()] = std::move(segmentPoints[segmenti]);
2895  Pstream::gatherList(gatheredPts);
2896 
2897  List<scalarList> gatheredDist(Pstream::nProcs());
2898  gatheredDist[Pstream::myProcNo()] = std::move(segmentDist[segmenti]);
2899  Pstream::gatherList(gatheredDist);
2900 
2901  // Combine processor lists into one big list.
2902  pointList allPts
2903  (
2904  ListListOps::combine<pointList>
2905  (
2906  gatheredPts, accessOp<pointList>()
2907  )
2908  );
2909  scalarList allDist
2910  (
2911  ListListOps::combine<scalarList>
2912  (
2913  gatheredDist, accessOp<scalarList>()
2914  )
2915  );
2916 
2917  // Sort according to distance
2918  labelList indexSet(Foam::sortedOrder(allDist));
2919 
2920  allLeakPaths.set
2921  (
2922  segmenti,
2923  new coordSet
2924  (
2925  leakPath.name(),
2926  leakPath.axis(),
2927  pointList(allPts, indexSet),
2928  //scalarList(allDist, indexSet)
2929  scalarList(allPts.size(), scalar(segmenti))
2930  )
2931  );
2932  }
2933 
2934  fileName fName;
2935  if (Pstream::master())
2936  {
2937  List<scalarField> allLeakData(allLeakPaths.size());
2938  forAll(allLeakPaths, segmenti)
2939  {
2940  allLeakData[segmenti] = allLeakPaths[segmenti].distance();
2941  }
2942 
2943  writer.nFields(1);
2944 
2945  writer.open
2946  (
2947  allLeakPaths,
2948  (outputDir / allLeakPaths[0].name())
2949  );
2950 
2951  fName = writer.write("leakPath", allLeakData);
2952 
2953  // Force writing to finish before FatalError exit
2954  writer.close(true);
2955  }
2956 
2957  // Probably do not need to broadcast name (only written on master anyhow)
2958  Pstream::broadcast(fName);
2959 
2960  return fName;
2961 }
2962 
2963 
2964 // Modify cellRegion to be consistent with locationsInMesh.
2965 // - all regions not in locationsInMesh are set to -1
2966 // - check that all regions inside locationsOutsideMesh are set to -1
2968 (
2969  const polyMesh& mesh,
2970  const vector& perturbVec,
2971  const pointField& locationsInMesh,
2972  const pointField& locationsOutsideMesh,
2973  const label nRegions,
2974  labelList& cellRegion,
2975  const boolList& blockedFace,
2976  // Leak-path
2977  const bool exitIfLeakPath,
2978  const refPtr<coordSetWriter>& leakPathFormatter
2979 )
2980 {
2981  bitSet insideCell(mesh.nCells());
2982 
2983  // Mark all cells reachable from locationsInMesh
2984  labelList insideRegions(locationsInMesh.size());
2985  forAll(insideRegions, i)
2986  {
2987  // Find the region containing the point
2988  label regioni = findRegion
2989  (
2990  mesh,
2991  cellRegion,
2992  perturbVec,
2993  locationsInMesh[i]
2994  );
2995 
2996  insideRegions[i] = regioni;
2997 
2998  // Mark all cells in the region as being inside
2999  forAll(cellRegion, celli)
3000  {
3001  if (cellRegion[celli] == regioni)
3002  {
3003  insideCell.set(celli);
3004  }
3005  }
3006  }
3007 
3008 
3009 
3010  // Check that all the locations outside the
3011  // mesh do not conflict with those inside
3012  forAll(locationsOutsideMesh, i)
3013  {
3014  // Find the region containing the point,
3015  // and the corresponding inside region index
3016 
3017  label indexi;
3018  label regioni = findRegion
3019  (
3020  mesh,
3021  cellRegion,
3022  perturbVec,
3023  locationsOutsideMesh[i]
3024  );
3025 
3026  if (regioni == -1 && (indexi = insideRegions.find(regioni)) != -1)
3027  {
3028  if (leakPathFormatter)
3029  {
3030  const fileName fName
3031  (
3032  writeLeakPath
3033  (
3034  mesh,
3035  locationsInMesh,
3036  locationsOutsideMesh,
3037  blockedFace,
3038  leakPathFormatter.constCast()
3039  )
3040  );
3041  Info<< "Dumped leak path to " << fName << endl;
3042  }
3043 
3044  auto& err =
3045  (
3046  exitIfLeakPath
3049  );
3050 
3051  err << "Location in mesh " << locationsInMesh[indexi]
3052  << " is inside same mesh region " << regioni
3053  << " as one of the locations outside mesh "
3054  << locationsOutsideMesh << endl;
3055 
3056  if (exitIfLeakPath)
3057  {
3059  }
3060  }
3061  }
3062 
3063 
3064  label nRemove = 0;
3065 
3066  // Now update cellRegion to -1 for unreachable cells
3067  forAll(insideCell, celli)
3068  {
3069  if (!insideCell.test(celli))
3070  {
3071  cellRegion[celli] = -1;
3072  ++nRemove;
3073  }
3074  else if (cellRegion[celli] == -1)
3075  {
3076  ++nRemove;
3077  }
3078  }
3079 
3080  return nRemove;
3081 }
3082 
3083 
3085 (
3086  const labelList& globalToMasterPatch,
3087  const labelList& globalToSlavePatch,
3088  const pointField& locationsInMesh,
3089  const pointField& locationsOutsideMesh,
3090  const bool exitIfLeakPath,
3091  const refPtr<coordSetWriter>& leakPathFormatter
3092 )
3093 {
3094  // Force calculation of face decomposition (used in findCell)
3095  (void)mesh_.tetBasePtIs();
3096 
3097  // Determine connected regions. regionSplit is the labelList with the
3098  // region per cell.
3099 
3100  boolList blockedFace(mesh_.nFaces(), false);
3101  selectSeparatedCoupledFaces(blockedFace);
3102 
3103  regionSplit cellRegion(mesh_, blockedFace);
3104 
3105  label nRemove = findRegions
3106  (
3107  mesh_,
3108  vector::uniform(mergeDistance_), // perturbVec
3109  locationsInMesh,
3110  locationsOutsideMesh,
3111  cellRegion.nRegions(),
3112  cellRegion,
3113  blockedFace,
3114  // Leak-path
3115  exitIfLeakPath,
3116  leakPathFormatter
3117  );
3118 
3119  // Subset
3120  // ~~~~~~
3121 
3122  // Get cells to remove
3123  DynamicList<label> cellsToRemove(nRemove);
3124  forAll(cellRegion, celli)
3125  {
3126  if (cellRegion[celli] == -1)
3127  {
3128  cellsToRemove.append(celli);
3129  }
3130  }
3131  cellsToRemove.shrink();
3132 
3133  label nTotCellsToRemove = returnReduce
3134  (
3135  cellsToRemove.size(),
3136  sumOp<label>()
3137  );
3138 
3139 
3140  autoPtr<mapPolyMesh> mapPtr;
3141  if (nTotCellsToRemove > 0)
3142  {
3143  label nCellsToKeep =
3144  mesh_.globalData().nTotalCells()
3145  - nTotCellsToRemove;
3146 
3147  Info<< "Keeping all cells containing points " << locationsInMesh << endl
3148  << "Selected for keeping : "
3149  << nCellsToKeep
3150  << " cells." << endl;
3151 
3152 
3153  // Remove cells
3154  removeCells cellRemover(mesh_);
3155 
3156  labelList exposedFaces(cellRemover.getExposedFaces(cellsToRemove));
3157  labelList exposedPatch;
3158 
3159  label nExposedFaces = returnReduce(exposedFaces.size(), sumOp<label>());
3160  if (nExposedFaces)
3161  {
3162  // FatalErrorInFunction
3163  // << "Removing non-reachable cells should only expose"
3164  // << " boundary faces" << nl
3165  // << "ExposedFaces:" << exposedFaces << abort(FatalError);
3166 
3167  // Patch for exposed faces for lack of anything sensible.
3168  label defaultPatch = 0;
3169  if (globalToMasterPatch.size())
3170  {
3171  defaultPatch = globalToMasterPatch[0];
3172  }
3173 
3175  << "Removing non-reachable cells exposes "
3176  << nExposedFaces << " internal or coupled faces." << endl
3177  << " These get put into patch " << defaultPatch << endl;
3178  exposedPatch.setSize(exposedFaces.size(), defaultPatch);
3179  }
3180 
3181  mapPtr = doRemoveCells
3182  (
3183  cellsToRemove,
3184  exposedFaces,
3185  exposedPatch,
3186  cellRemover
3187  );
3188  }
3189  return mapPtr;
3190 }
3191 
3192 
3195  // mesh_ already distributed; distribute my member data
3196 
3197  // surfaceQueries_ ok.
3198 
3199  // refinement
3200  meshCutter_.distribute(map);
3201 
3202  // surfaceIndex is face data.
3203  map.distributeFaceData(surfaceIndex_);
3204 
3205  // faceToPatch (baffles that were on coupled faces) is not maintained
3206  // (since baffling also disconnects points)
3207  faceToCoupledPatch_.clear();
3208 
3209  // maintainedFaces are indices of faces.
3210  forAll(userFaceData_, i)
3211  {
3212  map.distributeFaceData(userFaceData_[i].second());
3213  }
3214 
3215  // Redistribute surface and any fields on it.
3216  {
3217  Random rndGen(653213);
3218 
3219  // Get local mesh bounding box. Single box for now.
3221  (
3222  1,
3223  treeBoundBox(mesh_.points()).extend(rndGen, 1e-4)
3224  );
3225 
3226  // Distribute all geometry (so refinementSurfaces and shellSurfaces)
3227  searchableSurfaces& geometry =
3228  const_cast<searchableSurfaces&>(surfaces_.geometry());
3229 
3230  forAll(geometry, i)
3231  {
3233  autoPtr<mapDistribute> pointMap;
3234  geometry[i].distribute
3235  (
3236  meshBb,
3237  false, // do not keep outside triangles
3238  faceMap,
3239  pointMap
3240  );
3241 
3242  if (faceMap)
3243  {
3244  // (ab)use the instance() to signal current modification time
3245  geometry[i].instance() = geometry[i].time().timeName();
3246  }
3247 
3248  faceMap.clear();
3249  pointMap.clear();
3250  }
3251  }
3252 }
3253 
3254 
3256 (
3257  const mapPolyMesh& map,
3258  const labelList& changedFaces
3259 )
3260 {
3261  Map<label> dummyMap(0);
3262 
3263  updateMesh(map, changedFaces, dummyMap, dummyMap, dummyMap);
3264 }
3265 
3266 
3268 (
3269  const labelList& pointsToStore,
3270  const labelList& facesToStore,
3271  const labelList& cellsToStore
3272 )
3273 {
3274  // For now only meshCutter has storable/retrievable data.
3275  meshCutter_.storeData
3276  (
3277  pointsToStore,
3278  facesToStore,
3279  cellsToStore
3280  );
3281 }
3282 
3283 
3285 (
3286  const mapPolyMesh& map,
3287  const labelList& changedFaces,
3288  const Map<label>& pointsToRestore,
3289  const Map<label>& facesToRestore,
3290  const Map<label>& cellsToRestore
3291 )
3292 {
3293  // For now only meshCutter has storable/retrievable data.
3294 
3295  // Update numbering of cells/vertices.
3296  meshCutter_.updateMesh
3297  (
3298  map,
3299  pointsToRestore,
3300  facesToRestore,
3301  cellsToRestore
3302  );
3303 
3304  // Update surfaceIndex
3305  updateList(map.faceMap(), label(-1), surfaceIndex_);
3306 
3307  // Update faceToCoupledPatch_
3308  {
3309  Map<label> newFaceToPatch(faceToCoupledPatch_.size());
3310  forAllConstIters(faceToCoupledPatch_, iter)
3311  {
3312  const label newFacei = map.reverseFaceMap()[iter.key()];
3313 
3314  if (newFacei >= 0)
3315  {
3316  newFaceToPatch.insert(newFacei, iter.val());
3317  }
3318  }
3319  faceToCoupledPatch_.transfer(newFaceToPatch);
3320  }
3321 
3322 
3323  // Update cached intersection information
3324  updateIntersections(changedFaces);
3325 
3326  // Update maintained faces
3327  forAll(userFaceData_, i)
3328  {
3329  labelList& data = userFaceData_[i].second();
3330 
3331  if (userFaceData_[i].first() == KEEPALL)
3332  {
3333  // extend list with face-from-face data
3334  updateList(map.faceMap(), label(-1), data);
3335  }
3336  else if (userFaceData_[i].first() == MASTERONLY)
3337  {
3338  // keep master only
3339  labelList newFaceData(map.faceMap().size(), -1);
3340 
3341  forAll(newFaceData, facei)
3342  {
3343  label oldFacei = map.faceMap()[facei];
3344 
3345  if (oldFacei >= 0 && map.reverseFaceMap()[oldFacei] == facei)
3346  {
3347  newFaceData[facei] = data[oldFacei];
3348  }
3349  }
3350  data.transfer(newFaceData);
3351  }
3352  else
3353  {
3354  // remove any face that has been refined i.e. referenced more than
3355  // once.
3356 
3357  // 1. Determine all old faces that get referenced more than once.
3358  // These get marked with -1 in reverseFaceMap
3359  labelList reverseFaceMap(map.reverseFaceMap());
3360 
3361  forAll(map.faceMap(), facei)
3362  {
3363  label oldFacei = map.faceMap()[facei];
3364 
3365  if (oldFacei >= 0)
3366  {
3367  if (reverseFaceMap[oldFacei] != facei)
3368  {
3369  // facei is slave face. Mark old face.
3370  reverseFaceMap[oldFacei] = -1;
3371  }
3372  }
3373  }
3374 
3375  // 2. Map only faces with intact reverseFaceMap
3376  labelList newFaceData(map.faceMap().size(), -1);
3377  forAll(newFaceData, facei)
3378  {
3379  label oldFacei = map.faceMap()[facei];
3380 
3381  if (oldFacei >= 0)
3382  {
3383  if (reverseFaceMap[oldFacei] == facei)
3384  {
3385  newFaceData[facei] = data[oldFacei];
3386  }
3387  }
3388  }
3389  data.transfer(newFaceData);
3390  }
3391  }
3392 }
3393 
3394 
3395 bool Foam::meshRefinement::write() const
3397  bool writeOk = mesh_.write();
3398 
3399  // Make sure that any distributed surfaces (so ones which probably have
3400  // been changed) get written as well.
3401  // Note: should ideally have some 'modified' flag to say whether it
3402  // has been changed or not.
3403  searchableSurfaces& geometry =
3404  const_cast<searchableSurfaces&>(surfaces_.geometry());
3405 
3406  forAll(geometry, i)
3407  {
3408  searchableSurface& s = geometry[i];
3409 
3410  // Check if instance() of surface is not constant or system.
3411  // Is good hint that surface is distributed.
3412  if
3413  (
3414  s.instance() != s.time().system()
3415  && s.instance() != s.time().caseSystem()
3416  && s.instance() != s.time().constant()
3417  && s.instance() != s.time().caseConstant()
3418  )
3419  {
3420  // Make sure it gets written to current time, not constant.
3421  s.instance() = s.time().timeName();
3422  writeOk = writeOk && s.write();
3423  }
3424  }
3425 
3426  return writeOk;
3427 }
3428 
3429 
3431 (
3432  const polyMesh& mesh,
3433  const labelList& meshPoints
3434 )
3435 {
3436  const label myProci = UPstream::myProcNo();
3437  const globalIndex globalPoints(meshPoints.size());
3438 
3439  labelList myPoints
3440  (
3441  Foam::identity(globalPoints.range(myProci))
3442  );
3443 
3445  (
3446  mesh,
3447  meshPoints,
3448  myPoints,
3449  minEqOp<label>(),
3450  labelMax
3451  );
3452 
3453 
3454  bitSet isPatchMasterPoint(meshPoints.size());
3455  forAll(meshPoints, pointi)
3456  {
3457  if (myPoints[pointi] == globalPoints.toGlobal(myProci, pointi))
3458  {
3459  isPatchMasterPoint.set(pointi);
3460  }
3461  }
3462 
3463  return isPatchMasterPoint;
3464 }
3465 
3466 
3468 (
3469  const polyMesh& mesh,
3470  const labelList& meshEdges
3471 )
3472 {
3473  const label myProci = UPstream::myProcNo();
3474  const globalIndex globalEdges(meshEdges.size());
3475 
3476  labelList myEdges
3477  (
3478  Foam::identity(globalEdges.range(myProci))
3479  );
3480 
3482  (
3483  mesh,
3484  meshEdges,
3485  myEdges,
3486  minEqOp<label>(),
3487  labelMax
3488  );
3489 
3490 
3491  bitSet isMasterEdge(meshEdges.size());
3492  forAll(meshEdges, edgei)
3493  {
3494  if (myEdges[edgei] == globalEdges.toGlobal(myProci, edgei))
3495  {
3496  isMasterEdge.set(edgei);
3497  }
3498  }
3499 
3500  return isMasterEdge;
3501 }
3502 
3503 
3505 (
3506  const bool debug,
3507  const string& msg,
3508  const bool printCellLevel
3509 )
3510 const
3511 {
3512  const globalMeshData& pData = mesh_.globalData();
3513 
3514  if (debug)
3515  {
3516  Pout<< msg.c_str()
3517  << " : cells(local):" << mesh_.nCells()
3518  << " faces(local):" << mesh_.nFaces()
3519  << " points(local):" << mesh_.nPoints()
3520  << endl;
3521  }
3522 
3523  {
3524  bitSet isMasterFace(syncTools::getMasterFaces(mesh_));
3525  label nMasterFaces = isMasterFace.count();
3526 
3527  bitSet isMeshMasterPoint(syncTools::getMasterPoints(mesh_));
3528  label nMasterPoints = isMeshMasterPoint.count();
3529 
3530  Info<< msg.c_str()
3531  << " : cells:" << pData.nTotalCells()
3532  << " faces:" << returnReduce(nMasterFaces, sumOp<label>())
3533  << " points:" << returnReduce(nMasterPoints, sumOp<label>());
3534 
3535  if (UPstream::parRun())
3536  {
3537  const scalar nIdealCells =
3538  scalar(pData.nTotalCells())/Pstream::nProcs();
3539  const scalar unbalance = returnReduce
3540  (
3541  mag(1.0-mesh_.nCells()/nIdealCells),
3542  maxOp<scalar>()
3543  );
3544  Info<< " unbalance:" << unbalance;
3545  }
3546  Info<< endl;
3547  }
3548 
3549  if (printCellLevel)
3550  {
3551  const labelList& cellLevel = meshCutter_.cellLevel();
3552 
3553  labelList nCells(gMax(cellLevel)+1, Zero);
3554 
3555  forAll(cellLevel, celli)
3556  {
3557  nCells[cellLevel[celli]]++;
3558  }
3559 
3561 
3563  if (Pstream::master())
3564  {
3565  Info<< "Cells per refinement level:" << endl;
3566  forAll(nCells, leveli)
3567  {
3568  Info<< " " << leveli << '\t' << nCells[leveli]
3569  << endl;
3570  }
3571  }
3572  }
3573 }
3574 
3575 
3578  if (overwrite_ && mesh_.time().timeIndex() == 0)
3579  {
3580  return oldInstance_;
3581  }
3582 
3583  return mesh_.time().timeName();
3584 }
3585 
3586 
3589  // Note: use time().timeName(), not meshRefinement::timeName()
3590  // so as to dump the fields to 0, not to constant.
3591  {
3592  volScalarField volRefLevel
3593  (
3594  IOobject
3595  (
3596  "cellLevel",
3597  mesh_.time().timeName(),
3598  mesh_,
3602  ),
3603  mesh_,
3606  );
3607 
3608  const labelList& cellLevel = meshCutter_.cellLevel();
3609 
3610  forAll(volRefLevel, celli)
3611  {
3612  volRefLevel[celli] = cellLevel[celli];
3613  }
3614 
3615  volRefLevel.write();
3616  }
3617 
3618  // Dump pointLevel
3619  {
3620  const pointMesh& pMesh = pointMesh::New(mesh_);
3621 
3622  pointScalarField pointRefLevel
3623  (
3624  IOobject
3625  (
3626  "pointLevel",
3627  mesh_.time().timeName(),
3628  mesh_,
3632  ),
3633  pMesh,
3635  );
3636 
3637  const labelList& pointLevel = meshCutter_.pointLevel();
3638 
3639  forAll(pointRefLevel, pointi)
3640  {
3641  pointRefLevel[pointi] = pointLevel[pointi];
3642  }
3643 
3644  pointRefLevel.write();
3645  }
3646 }
3647 
3648 
3649 void Foam::meshRefinement::dumpIntersections(const fileName& prefix) const
3651  {
3652  OFstream str(prefix + "_edges.obj");
3653  label verti = 0;
3654  Pout<< "meshRefinement::dumpIntersections :"
3655  << " Writing cellcentre-cellcentre intersections to file "
3656  << str.name() << endl;
3657 
3658 
3659  // Redo all intersections
3660  // ~~~~~~~~~~~~~~~~~~~~~~
3661 
3662  // Get boundary face centre and level. Coupled aware.
3663  labelList neiLevel(mesh_.nBoundaryFaces());
3664  pointField neiCc(mesh_.nBoundaryFaces());
3665  calcNeighbourData(neiLevel, neiCc);
3666 
3667  labelList intersectionFaces(intersectedFaces());
3668 
3669  // Collect segments we want to test for
3670  pointField start(intersectionFaces.size());
3671  pointField end(intersectionFaces.size());
3672  {
3673  labelList minLevel;
3674  calcCellCellRays
3675  (
3676  neiCc,
3677  labelList(neiCc.size(), -1),
3678  intersectionFaces,
3679  start,
3680  end,
3681  minLevel
3682  );
3683  }
3684 
3685 
3686  // Do tests in one go
3687  labelList surfaceHit;
3688  List<pointIndexHit> surfaceHitInfo;
3689  surfaces_.findAnyIntersection
3690  (
3691  start,
3692  end,
3693  surfaceHit,
3694  surfaceHitInfo
3695  );
3696 
3697  forAll(intersectionFaces, i)
3698  {
3699  if (surfaceHit[i] != -1)
3700  {
3701  meshTools::writeOBJ(str, start[i]);
3702  verti++;
3703  meshTools::writeOBJ(str, surfaceHitInfo[i].point());
3704  verti++;
3705  meshTools::writeOBJ(str, end[i]);
3706  verti++;
3707  str << "l " << verti-2 << ' ' << verti-1 << nl
3708  << "l " << verti-1 << ' ' << verti << nl;
3709  }
3710  }
3711  }
3712 
3713  Pout<< endl;
3714 }
3715 
3716 
3718 (
3719  const debugType debugFlags,
3720  const writeType writeFlags,
3721  const fileName& prefix
3722 ) const
3723 {
3724  if (writeFlags & WRITEMESH)
3725  {
3726  write();
3727  }
3728 
3729  if (writeFlags && !(writeFlags & NOWRITEREFINEMENT))
3730  {
3731  meshCutter_.write();
3732 
3733  // Force calculation before writing
3734  (void)surfaceIndex();
3735  surfaceIndex_.write();
3736  }
3737 
3738  if (writeFlags & WRITELEVELS)
3739  {
3740  dumpRefinementLevel();
3741  }
3742 
3743  if ((debugFlags & OBJINTERSECTIONS) && prefix.size())
3744  {
3745  dumpIntersections(prefix);
3746  }
3747 }
3748 
3749 
3752  IOobject io
3753  (
3754  "dummy",
3755  mesh.facesInstance(),
3757  mesh
3758  );
3759  fileName setsDir(io.path());
3760 
3761  if (topoSet::debug) DebugVar(setsDir);
3762 
3763  // Remove local files
3764  if (exists(setsDir/"surfaceIndex"))
3765  {
3766  rm(setsDir/"surfaceIndex");
3767  }
3768 
3769  // Remove other files
3771 }
3772 
3773 
3776  return writeLevel_;
3777 }
3778 
3779 
3780 void Foam::meshRefinement::writeLevel(const writeType flags)
3782  writeLevel_ = flags;
3783 }
3784 
3785 
3786 //Foam::meshRefinement::outputType Foam::meshRefinement::outputLevel()
3787 //{
3788 // return outputLevel_;
3789 //}
3790 //
3791 //
3792 //void Foam::meshRefinement::outputLevel(const outputType flags)
3793 //{
3794 // outputLevel_ = flags;
3795 //}
3796 
3797 
3799 (
3801  const word& keyword,
3802  const bool noExit,
3803  enum keyType::option matchOpt
3804 )
3805 {
3806  const dictionary* dictptr = dict.findDict(keyword, matchOpt);
3807 
3808  if (!dictptr)
3809  {
3811  << "Entry '" << keyword
3812  << "' not found (or not a dictionary) in dictionary "
3813  << dict.relativeName() << nl;
3814 
3815  if (noExit)
3816  {
3817  // Dummy return
3818  return dictionary::null;
3819  }
3820  else
3821  {
3823  }
3824  }
3825 
3826  return *dictptr;
3827 }
3828 
3829 
3831 (
3833  const word& keyword,
3834  const bool noExit,
3835  enum keyType::option matchOpt
3836 )
3837 {
3838  const entry* eptr = dict.findEntry(keyword, matchOpt);
3839 
3840  if (!eptr)
3841  {
3843  << "Entry '" << keyword << "' not found in dictionary "
3844  << dict.relativeName() << nl;
3845 
3846  if (noExit)
3847  {
3848  // Dummy return
3849  return ITstream::empty_stream();
3850  }
3851  else
3852  {
3854  }
3855  }
3856 
3857  return eptr->stream();
3858 }
3859 
3860 
3861 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
static void listCombineGather(UList< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Combines List elements.
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:136
static bool checkMesh(const bool report, const polyMesh &mesh, const dictionary &dict, labelHashSet &wrongFaces, const bool dryRun=false)
Check mesh with mesh settings in dict. Collects incorrect faces.
const labelList patchIDs(pbm.indices(polyPatchNames, true))
const Field< point_type > & faceAreas() const
Return face area vectors for patch.
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search...
Definition: meshSearch.H:56
const polyBoundaryMesh & pbm
static void swapFaceList(const polyMesh &mesh, UList< T > &faceValues, const bool parRun=UPstream::parRun())
Swap coupled face values. Uses eqOp.
Definition: syncTools.H:567
const T & first() const noexcept
Access the first element.
Definition: Pair.H:137
label addFaceZone(const word &fzName, const word &masterPatch, const word &slavePatch, const surfaceZonesInfo::faceZoneType &fzType)
Add/lookup faceZone and update information. Return index of.
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const bool parRun=UPstream::parRun())
Synchronize values on all mesh faces.
Definition: syncTools.H:465
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
labelList getExposedFaces(const bitSet &removedCell) const
Get labels of faces exposed after cells removal.
Definition: removeCells.C:77
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type unset(const label i)
Unset the bool entry at specified position, always false for out-of-range access. ...
Definition: UList.H:805
unsigned int count(const bool on=true) const
Count number of bits set.
Definition: bitSetI.H:420
static const dictionary & subDict(const dictionary &dict, const word &keyword, const bool noExit, enum keyType::option matchOpt=keyType::REGEX)
Wrapper around dictionary::subDict which does not exit.
void clearAddressing()
Clear addressing.
Definition: ZoneMesh.C:905
Definition: ops.H:67
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:498
A class for handling file names.
Definition: fileName.H:72
A list of face labels.
Definition: faceSet.H:47
const DynamicList< face > & faces() const
label addPointZone(const word &name)
Add pointZone if does not exist. Return index of zone.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Foam::pointBoundaryMesh.
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:859
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/"finiteArea-edgesCentres"))
static bitSet getMasterFaces(const polyMesh &mesh)
Get per face whether it is uncoupled or a master of a coupled set of faces.
Definition: syncTools.C:119
label countHits() const
Count number of intersections (local)
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
label nPoints() const noexcept
Number of mesh points.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
static void removeFiles(const polyMesh &)
Helper: remove all relevant files from mesh instance.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:521
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:899
const bitSet isBlockedFace(intersectedFaces())
bool getFaceZoneInfo(const word &fzName, label &masterPatchID, label &slavePatchID, surfaceZonesInfo::faceZoneType &fzType) const
Lookup faceZone information. Return false if no information.
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition: OSstream.H:134
static label findRegions(const polyMesh &, const vector &perturbVec, const pointField &locationsInMesh, const pointField &locationsOutsideMesh, const label nRegions, labelList &cellRegion, const boolList &blockedFace, const bool exitIfLeakPath, const refPtr< coordSetWriter > &leakPathFormatter)
Find regions points are in.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
static void removeFiles(const polyMesh &)
Helper: remove all relevant files from mesh instance.
Definition: hexRef8.C:5810
const mapDistribute & globalEdgeSlavesMap() const
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create MeshObject registered with typeName.
Definition: MeshObject.C:53
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
static label findRegion(const polyMesh &, const labelList &cellRegion, const vector &perturbVec, const point &p)
Find region point is in. Uses optional perturbation to re-test.
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
IntListType renumber(const labelUList &oldToNew, const IntListType &input)
Renumber the values within a list.
Output to file stream as an OSstream, normally using std::ofstream for the actual output...
Definition: OFstream.H:71
void printMeshInfo(const bool debug, const string &msg, const bool printCellLevel) const
Print some mesh stats.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
autoPtr< mapPolyMesh > doRemoveCells(const labelList &cellsToRemove, const labelList &exposedFaces, const labelList &exposedPatchIDs, removeCells &cellRemover)
Remove cells. Put exposedFaces into exposedPatchIDs.
static ITstream & lookup(const dictionary &dict, const word &keyword, const bool noExit, enum keyType::option matchOpt=keyType::REGEX)
Wrapper around dictionary::lookup which does not exit.
label splitFacesUndo(const labelList &splitFaces, const labelPairList &splits, const dictionary &motionDict, labelList &duplicateFace, List< labelPair > &baffles)
Split faces along diagonal. Maintain mesh quality. Return.
static void syncBoundaryFacePositions(const polyMesh &mesh, UList< point > &positions, const CombineOp &cop)
Synchronize locations on boundary faces only.
Definition: syncTools.H:445
static void matchEdges(const PrimitivePatch< FaceList1, PointField1 > &p1, const PrimitivePatch< FaceList2, PointField2 > &p2, labelList &p1EdgeLabels, labelList &p2EdgeLabels, bitSet &sameOrientation)
Find corresponding edges on patches sharing the same points.
T * data() noexcept
Return pointer to the underlying array serving as data storage.
Definition: UListI.H:265
Given list of cells to remove, insert all the topology changes.
Definition: removeCells.H:59
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const bool parRun=UPstream::parRun())
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:524
Random rndGen
Definition: createFields.H:23
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:411
void reorder(const labelUList &oldToNew, const bool validBoundary)
Reorders patches. Ordering does not have to be done in.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static void gatherList(const UList< commsStruct > &comms, UList< T > &values, const int tag, const label comm)
Gather data, but keep individual values separate. Uses the specified communication schedule...
void setConstraints(const polyMesh &mesh, boolList &blockedFace, PtrList< labelList > &specifiedProcessorFaces, labelList &specifiedProcessor, List< labelPair > &explicitConnections) const
Helper: extract constraints:
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1061
void reorder(const labelUList &oldToNew, const bool check=false)
Reorder elements. Reordering must be unique (ie, shuffle).
Definition: UPtrList.C:79
void clear() noexcept
Same as reset(nullptr)
Definition: autoPtr.H:255
label constructSize() const noexcept
Constructed data size.
labelList intersectedFaces() const
Get faces with intersection.
static writeType writeLevel()
Get/set write level.
void clearOut(const bool isMeshUpdate=false)
Clear all geometry and addressing.
Definition: fvMesh.C:227
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
static const word & calculatedType() noexcept
The type name for calculated patch fields.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:196
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:397
static bitSet getMasterEdges(const polyMesh &mesh, const labelList &meshEdges)
Determine master edge for subset of edges. If coupled.
Base class of (analytical or triangulated) surface. Encapsulates all the search routines. WIP.
Ignore writing from objectRegistry::writeObject()
const labelListList & globalEdgeSlaves() const
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
const dimensionSet dimless
Dimensionless.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1086
static void calculateEdgeWeights(const polyMesh &mesh, const bitSet &isMasterEdge, const labelList &meshPoints, const edgeList &edges, scalarField &edgeWeights, scalarField &invSumWeight)
Helper: calculate edge weights (1/length)
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
label nFaces() const noexcept
Number of mesh faces.
labelList intersectedPoints() const
Get points on surfaces with intersection and boundary faces.
autoPtr< mapPolyMesh > changeMesh(polyMesh &mesh, const labelUList &patchMap, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
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.
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:104
Finds shortest path (in terms of cell centres) to walk on mesh from any point in insidePoints to any ...
void write(const word &fieldName, const UList< Type > &field)
Write primitive field of CellData (Poly or Line) or PointData values.
A class for managing references or pointers (no reference counting)
Definition: HashPtrTable.H:49
T & emplace_back(Args &&... args)
Construct and append an element to the end of the list, return reference to the new list element...
Definition: PtrListI.H:105
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
static void broadcast(Type &value, const label comm=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all communicator ranks. Does nothing in non-paral...
static labelPairList findDuplicateFacePairs(const polyMesh &)
Helper routine to find all baffles (two boundary faces.
fileName path() const
The complete path for the object (with instance, local,...).
Definition: IOobject.C:480
Container for data on surfaces used for surface-driven refinement. Contains all the data about the le...
static label addPatch(fvMesh &, const word &name, const dictionary &)
Helper:add patch to mesh. Update all registered fields.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:45
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static void syncData(List< Type > &elems, const labelListList &slaves, const labelListList &transformedSlaves, const mapDistribute &slavesMap, const globalIndexAndTransform &, const CombineOp &cop, const TransformOp &top)
Helper: synchronise data with transforms.
const labelList & reverseFaceMap() const noexcept
Reverse face map.
Definition: mapPolyMesh.H:620
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Definition: UListI.H:90
T & constCast() const
Return non-const reference to the object or to the contents of a (non-null) managed pointer...
Definition: refPtr.H:278
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:37
static const Enum< writeType > writeTypeNames
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Encapsulates queries for features.
void distributeFaceData(List< T > &values) const
Distribute list of face data.
static label appendPatch(fvMesh &, const label insertPatchi, const word &, const dictionary &)
Helper:append patch to end of mesh.
void doSplitFaces(const labelList &splitFaces, const labelPairList &splits, polyTopoChange &meshMod) const
Split faces into two.
void distribute(const mapDistributePolyMesh &)
Update local numbering for mesh redistribution.
static void updateList(const labelList &newToOld, const T &nullValue, List< T > &elems)
Helper: reorder list according to map.
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
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.
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of &#39;true&#39; entries.
Definition: BitOps.H:73
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:61
A List obtained as a section of another List.
Definition: SubList.H:50
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1077
Sends/receives parts of mesh+fvfields to neighbouring processors. Used in load balancing.
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:853
void setSize(const label n)
Alias for resize()
Definition: List.H:320
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
label addFace(const face &f, const label own, const label nei, const label masterPointID, const label masterEdgeID, const label masterFaceID, const bool flipFaceFlux, const label patchID, const label zoneID, const bool zoneFlip)
Add face to cells. Return new face label.
const pointField & points
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:718
const_iterator cfind(const Key &key) const
Find and return an const_iterator set at the hashed entry.
Definition: HashTableI.H:113
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
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
fileName globalPath() const
Return global path for the case = rootPath/globalCaseName. Same as TimePaths::globalPath() ...
Definition: Time.H:512
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static bitSet getMasterPoints(const polyMesh &mesh)
Get per point whether it is uncoupled or a master of a coupled set of points.
Definition: syncTools.C:61
const DynamicList< point > & points() const
Points. Shrunk after constructing mesh (or calling of compact())
wordList names() const
Return a list of patch names.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
const Field< point_type > & faceCentres() const
Return face centres for patch.
Base class for writing coordSet(s) and tracks with fields.
static void checkCoupledFaceZones(const polyMesh &)
Helper function: check that face zones are synced.
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:486
Container for searchableSurfaces. The collection is specified as a dictionary. For example...
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
bool exists(const fileName &name, const bool checkGzip=true, const bool followLink=true)
Does the name exist (as DIRECTORY or FILE) in the file system?
Definition: POSIX.C:835
const labelList & surfaceIndex() const
Per start-end edge the index of the surface hit.
static void swapBoundaryFacePositions(const polyMesh &mesh, UList< point > &positions, const bool parRun=UPstream::parRun())
Swap coupled positions. Uses eqOp.
Definition: syncTools.H:545
void modifyFace(const face &f, const label facei, const label own, const label nei, const bool flipFaceFlux, const label patchID, const label zoneID, const bool zoneFlip, const bool multiZone=false)
Modify vertices or cell of face.
const labelListList & edgeFaces() const
Return edge-face addressing.
Abstract base class for domain decomposition.
label nInternalFaces() const noexcept
Number of internal faces.
Encapsulates queries for volume refinement (&#39;refine all cells within shell&#39;).
Definition: shellSurfaces.H:53
const labelList & reversePointMap() const noexcept
Reverse point map.
Definition: mapPolyMesh.H:582
autoPtr< mapDistributePolyMesh > balance(const bool keepZoneFaces, const bool keepBaffles, const scalarField &cellWeights, decompositionMethod &decomposer, fvMeshDistribute &distributor)
Redecompose according to cell count.
Vector< scalar > vector
Definition: vector.H:57
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
label whichPatch(const label meshFacei) const
Return patch index for a given mesh face index. Uses binary search.
Accumulating histogram of values. Specified bin resolution automatic generation of bins...
Definition: distribution.H:57
const HashTable< labelList > & groupPatchIDs() const
The patch indices per patch group.
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition: removeCells.H:149
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
bool test(const label pos) const
Test for True value at specified position, never auto-vivify entries.
Definition: bitSet.H:329
Random number generator.
Definition: Random.H:55
void setRefinement(const bitSet &removedCell, const labelUList &facesToExpose, const labelUList &patchIDs, polyTopoChange &) const
Play commands into polyTopoChange to remove cells.
Definition: removeCells.C:181
word timeName() const
Replacement for Time::timeName() that returns oldInstance (if overwrite_)
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
void setCapacity(const label nPoints, const label nFaces, const label nCells)
Explicitly pre-size the dynamic storage for expected mesh size for if construct-without-mesh.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
void checkData()
Debugging: check that all faces still obey start()>end()
int debug
Static debugging option.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:51
Type gMax(const FieldField< Field, Type > &f)
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:671
virtual bool open(const fileName &file, bool parallel=UPstream::parRun())
Open file for writing (creates parent directory).
defineTypeNameAndDebug(combustionModel, 0)
static bool clean(std::string &str)
Cleanup filename string, possibly applies other transformations such as changing the path separator e...
Definition: fileName.C:192
label addMeshedPatch(const word &name, const dictionary &)
Add patch originating from meshing. Update meshedPatches_.
labelList meshedPatches() const
Get patchIDs for patches added in addMeshedPatch.
autoPtr< mapDistributePolyMesh > distribute(const labelList &dist)
Send cells to neighbours according to distribution.
labelList f(nPoints)
void setSize(const label newLen)
Same as resize()
Definition: PtrList.H:337
label nTotalCells() const noexcept
Total global number of mesh cells.
static autoPtr< indirectPrimitivePatch > makePatch(const polyMesh &, const labelList &)
Create patch from set of patches.
fileName relativeName(const bool caseTag=false) const
The dictionary name relative to the case.
Definition: dictionary.C:179
void updateMesh(const mapPolyMesh &, const labelList &changedFaces)
Update for external change to mesh. changedFaces are in new mesh.
void updateIntersections(const labelList &changedFaces)
Find any intersection of surface. Store in surfaceIndex_.
static const Foam::polyMesh::cellDecomposition findCellMode(Foam::polyMesh::FACE_DIAG_TRIS)
bool hasMotionPoints() const noexcept
Has valid preMotionPoints?
Definition: mapPolyMesh.H:767
autoPtr< mapPolyMesh > splitMeshRegions(const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const pointField &locationsInMesh, const pointField &locationsOutsideMesh, const bool exitIfLeakPath, const refPtr< coordSetWriter > &leakPathFormatter)
Split mesh. Keep part containing point. Return empty map if.
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;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Helper class which maintains intersections of (changing) mesh with (static) surfaces.
Use additional distance field for (scalar) axis.
static bool split(const std::string &line, std::string &key, std::string &val)
Definition: cpuInfo.C:32
const indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual linear/tree communicat...
void setInstance(const fileName &)
Set instance of all local IOobjects.
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:98
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values within a list.
Class containing processor-to-processor mapping information.
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition: polyPatch.C:255
void storeData(const labelList &pointsToStore, const labelList &facesToStore, const labelList &cellsToStore)
Signal points/face/cells for which to store data.
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.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
static void testSyncPointList(const string &msg, const polyMesh &mesh, const List< scalar > &fld)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:637
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label nCells() const noexcept
Number of mesh cells.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Foam::fvBoundaryMesh.
static bitSet getMasterPoints(const polyMesh &mesh, const labelList &meshPoints)
Determine master point for subset of points. If coupled.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:52
Direct mesh changes based on v1.3 polyTopoChange syntax.
void close()
End the file contents and close the file after writing.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1094
const polyBoundaryMesh & patches
const word & name() const noexcept
The zone name.
Nothing to be read.
option
Enumeration for the data type and search/match modes (bitmask)
Definition: keyType.H:79
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: fvPatchField.H:204
Automatically write from objectRegistry::writeObject()
static word outputPrefix
Directory prefix.
void push_back(T *ptr)
Append an element to the end of the list.
Definition: PtrListI.H:114
For use with FaceCellWave. Determines topological distance to starting faces. Templated on passive tr...
label nRegions() const
Return total number of regions.
Definition: regionSplit.H:337
const entry * findEntry(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:84
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1511
static const Enum< debugType > debugTypeNames
const boolList & flipMap() const noexcept
Return face flip map.
Definition: faceZone.H:367
static labelList findDuplicateFaces(const primitiveMesh &, const labelList &)
Helper routine to find baffles (two boundary faces using the.
void selectSeparatedCoupledFaces(boolList &) const
Select coupled faces that are not collocated.
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:451
messageStream Info
Information stream (stdout output on master, null elsewhere)
constexpr label labelMax
Definition: label.H:55
static ITstream & empty_stream()
Return reference to an empty ITstream, for functions needing to return an ITstream reference but whic...
Definition: ITstream.C:68
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
SubField< vector > subField
Declare type of subField.
Definition: Field.H:128
writeType
Enumeration for what to write. Used as a bit-pattern.
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return pointer to a new patch created on freestore from components.
Definition: polyPatchNew.C:28
faceZoneType
What to do with faceZone faces.
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
#define DebugVar(var)
Report a variable name and value.
T * data() noexcept
Return pointer to the underlying array serving as data storage.
Definition: FixedListI.H:121
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:765
static autoPtr< fvPatch > New(const polyPatch &, const fvBoundaryMesh &)
Return a pointer to a new patch created on freestore from polyPatch.
Definition: fvPatchNew.C:28
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
const labelList & faceMap() const noexcept
Old face map.
Definition: mapPolyMesh.H:503
static void syncEdgeList(const polyMesh &mesh, List< T > &edgeValues, const CombineOp &cop, const T &nullValue, const TransformOp &top, const FlipOp &fop)
Synchronize values on all mesh edges.
bool write() const
Write mesh and all data.
List< point > pointList
List of point.
Definition: pointList.H:32
List< label > labelList
A List of labels.
Definition: List.H:62
static tmp< pointVectorField > makeDisplacementField(const pointMesh &pMesh, const labelList &adaptPatchIDs)
Helper function to make a pointVectorField with correct.
volScalarField & p
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
bool coupled
const T & second() const noexcept
Access the second element.
Definition: Pair.H:147
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))
const pointBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: pointMesh.H:140
void dumpIntersections(const fileName &prefix) const
Debug: Write intersection information to OBJ format.
virtual ITstream & stream() const =0
Return token stream, if entry is a primitive entry.
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
static labelList countCells(const labelList &)
Helper function: count cells per processor in wanted distribution.
void removeFace(const label facei, const label mergeFacei)
Remove/merge face.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Map< label > invertToMap(const labelUList &values)
Create inverse mapping, which is a lookup table into the given list.
Definition: ListOps.C:107
List< bool > boolList
A List of bools.
Definition: List.H:60
A List with indirect addressing.
Definition: IndirectList.H:60
void dumpRefinementLevel() const
Write refinement level as volScalarFields for postprocessing.
const pointField & preMotionPoints() const noexcept
Pre-motion point positions.
Definition: mapPolyMesh.H:759
labelList countEdgeFaces(const uindirectPrimitivePatch &pp) const
Count number of faces per patch edge. Parallel consistent.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
label size() const noexcept
Number of entries.
Definition: PackedList.H:393
Do not request registration (bool: false)
void reduce(T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) using linear/tree communication schedule.
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
An input stream of tokens.
Definition: ITstream.H:52
static label addFaceZone(const word &name, const labelList &addressing, const boolList &flipMap, polyMesh &mesh)
static label findCell(const polyMesh &, const vector &perturbVec, const point &p)
Find cell point is in. Uses optional perturbation to re-test.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
A keyword and a list of tokens is an &#39;entry&#39;.
Definition: entry.H:63
virtual labelList decompose(const pointField &points, const scalarField &pointWeights=scalarField::null()) const
Return the wanted processor number for every coordinate, using uniform or specified point weights...
static void listCombineReduce(List< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Combines List elements. After completion all processors have the same data.
static const Enum< coordFormat > coordFormatNames
String representation of coordFormat enum.
Definition: coordSet.H:74
const labelListList & globalEdgeTransformedSlaves() const
bool rm(const fileName &file)
Remove a file (or its gz equivalent), returning true if successful.
Definition: POSIX.C:1404
const pointField & pts
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and it is a dictionary) otherwise return nullptr...
Definition: dictionaryI.H:124
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127