meshRefinementProblemCells.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-2015 OpenFOAM Foundation
9  Copyright (C) 2015-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "meshRefinement.H"
30 #include "fvMesh.H"
31 #include "syncTools.H"
32 #include "Time.H"
33 #include "refinementSurfaces.H"
34 #include "pointSet.H"
35 #include "faceSet.H"
36 #include "indirectPrimitivePatch.H"
37 #include "cellSet.H"
38 #include "searchableSurfaces.H"
39 #include "polyMeshGeometry.H"
40 #include "IOmanip.H"
41 #include "unitConversion.H"
42 #include "snappySnapDriver.H"
43 
44 #include "snapParameters.H"
45 #include "motionSmoother.H"
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::meshRefinement::markBoundaryFace
50 (
51  const label facei,
52  boolList& isBoundaryFace,
53  boolList& isBoundaryEdge,
54  boolList& isBoundaryPoint
55 ) const
56 {
57  isBoundaryFace[facei] = true;
58 
59  const labelList& fEdges = mesh_.faceEdges(facei);
60 
61  for (const label edgei : fEdges)
62  {
63  isBoundaryEdge[edgei] = true;
64  }
65 
66  const face& f = mesh_.faces()[facei];
67 
68  for (const label pointi : f)
69  {
70  isBoundaryPoint[pointi] = true;
71  }
72 }
73 
74 
75 void Foam::meshRefinement::findNearest
76 (
77  const labelList& meshFaces,
78  List<pointIndexHit>& nearestInfo,
79  labelList& nearestSurface,
80  labelList& nearestRegion,
81  vectorField& nearestNormal
82 ) const
83 {
84  pointField fc(meshFaces.size());
85  forAll(meshFaces, i)
86  {
87  fc[i] = mesh_.faceCentres()[meshFaces[i]];
88  }
89 
90  const labelList allSurfaces(identity(surfaces_.surfaces().size()));
91 
92  surfaces_.findNearest
93  (
94  allSurfaces,
95  fc,
96  scalarField(fc.size(), sqr(GREAT)), // sqr of attraction
97  nearestSurface,
98  nearestInfo
99  );
100 
101  // Do normal testing per surface.
102  nearestNormal.setSize(nearestInfo.size());
103  nearestRegion.setSize(nearestInfo.size());
104 
105  forAll(allSurfaces, surfI)
106  {
107  DynamicList<pointIndexHit> localHits;
108 
109  forAll(nearestSurface, i)
110  {
111  if (nearestSurface[i] == surfI)
112  {
113  localHits.append(nearestInfo[i]);
114  }
115  }
116 
117  label geomI = surfaces_.surfaces()[surfI];
118 
119  pointField localNormals;
120  surfaces_.geometry()[geomI].getNormal(localHits, localNormals);
121 
122  labelList localRegion;
123  surfaces_.geometry()[geomI].getRegion(localHits, localRegion);
124 
125  label localI = 0;
126  forAll(nearestSurface, i)
127  {
128  if (nearestSurface[i] == surfI)
129  {
130  nearestNormal[i] = localNormals[localI];
131  nearestRegion[i] = localRegion[localI];
132  localI++;
133  }
134  }
135  }
136 }
137 
138 
139 Foam::Map<Foam::label> Foam::meshRefinement::findEdgeConnectedProblemCells
140 (
141  const scalarField& perpendicularAngle,
142  const labelList& globalToPatch
143 ) const
144 {
145  // Construct addressing engine from all patches added for meshing.
146  autoPtr<indirectPrimitivePatch> ppPtr
147  (
149  (
150  mesh_,
151  meshedPatches()
152  )
153  );
154  const indirectPrimitivePatch& pp = ppPtr();
155 
156 
157  // 1. Collect faces to test
158  // ~~~~~~~~~~~~~~~~~~~~~~~~
159 
160  DynamicList<label> candidateFaces(pp.size()/20);
161 
162  const labelListList& edgeFaces = pp.edgeFaces();
163 
164  const labelList& cellLevel = meshCutter_.cellLevel();
165 
166  forAll(edgeFaces, edgeI)
167  {
168  const labelList& eFaces = edgeFaces[edgeI];
169 
170  if (eFaces.size() == 2)
171  {
172  label face0 = pp.addressing()[eFaces[0]];
173  label face1 = pp.addressing()[eFaces[1]];
174 
175  label cell0 = mesh_.faceOwner()[face0];
176  label cell1 = mesh_.faceOwner()[face1];
177 
178  if (cellLevel[cell0] > cellLevel[cell1])
179  {
180  // cell0 smaller.
181  const vector& n0 = pp.faceNormals()[eFaces[0]];
182  const vector& n1 = pp.faceNormals()[eFaces[1]];
183 
184  if (mag(n0 & n1) < 0.1)
185  {
186  candidateFaces.append(face0);
187  }
188  }
189  else if (cellLevel[cell1] > cellLevel[cell0])
190  {
191  // cell1 smaller.
192  const vector& n0 = pp.faceNormals()[eFaces[0]];
193  const vector& n1 = pp.faceNormals()[eFaces[1]];
194 
195  if (mag(n0 & n1) < 0.1)
196  {
197  candidateFaces.append(face1);
198  }
199  }
200  }
201  }
202  candidateFaces.shrink();
203 
204  Info<< "Testing " << returnReduce(candidateFaces.size(), sumOp<label>())
205  << " faces on edge-connected cells of differing level."
206  << endl;
207 
209  {
210  faceSet fSet(mesh_, "edgeConnectedFaces", candidateFaces);
211  fSet.instance() = timeName();
212  Pout<< "Writing " << fSet.size()
213  << " with problematic topology to faceSet "
214  << fSet.objectPath() << endl;
215  fSet.write();
216  }
217 
218 
219  // 2. Find nearest surface on candidate faces
220  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
221 
222  List<pointIndexHit> nearestInfo;
223  labelList nearestSurface;
224  labelList nearestRegion;
225  vectorField nearestNormal;
226  findNearest
227  (
228  candidateFaces,
229  nearestInfo,
230  nearestSurface,
231  nearestRegion,
232  nearestNormal
233  );
234 
235 
236  // 3. Test angle to surface
237  // ~~~~~~~~~~~~~~~~~~~~~~~~
238 
239  Map<label> candidateCells(candidateFaces.size());
240 
241  faceSet perpFaces(mesh_, "perpendicularFaces", pp.size()/100);
242 
243  forAll(candidateFaces, i)
244  {
245  label facei = candidateFaces[i];
246 
247  const vector n = normalised(mesh_.faceAreas()[facei]);
248 
249  label region = surfaces_.globalRegion
250  (
251  nearestSurface[i],
252  nearestRegion[i]
253  );
254 
255  scalar angle = degToRad(perpendicularAngle[region]);
256 
257  if (angle >= 0)
258  {
259  if (mag(n & nearestNormal[i]) < Foam::sin(angle))
260  {
261  perpFaces.insert(facei);
262  candidateCells.insert
263  (
264  mesh_.faceOwner()[facei],
265  globalToPatch[region]
266  );
267  }
268  }
269  }
270 
272  {
273  perpFaces.instance() = timeName();
274  Pout<< "Writing " << perpFaces.size()
275  << " faces that are perpendicular to the surface to set "
276  << perpFaces.objectPath() << endl;
277  perpFaces.write();
278  }
279  return candidateCells;
280 }
281 
282 
283 // Check if moving face to new points causes a 'collapsed' face.
284 // Uses new point position only for the face, not the neighbouring
285 // cell centres
286 bool Foam::meshRefinement::isCollapsedFace
287 (
288  const pointField& points,
289  const pointField& neiCc,
290  const scalar minFaceArea,
291  const scalar maxNonOrtho,
292  const label facei
293 ) const
294 {
295  // Severe nonorthogonality threshold
296  const scalar severeNonorthogonalityThreshold =
297  ::cos(degToRad(maxNonOrtho));
298 
299  vector s = mesh_.faces()[facei].areaNormal(points);
300  scalar magS = mag(s);
301 
302  // Check face area
303  if (magS < minFaceArea)
304  {
305  return true;
306  }
307 
308  // Check orthogonality
309  const point& ownCc = mesh_.cellCentres()[mesh_.faceOwner()[facei]];
310 
311  if (mesh_.isInternalFace(facei))
312  {
313  label nei = mesh_.faceNeighbour()[facei];
314  vector d = mesh_.cellCentres()[nei] - ownCc;
315 
316  scalar dDotS = (d & s)/(mag(d)*magS + VSMALL);
317 
318  if (dDotS < severeNonorthogonalityThreshold)
319  {
320  return true;
321  }
322  else
323  {
324  return false;
325  }
326  }
327  else
328  {
329  label patchi = mesh_.boundaryMesh().whichPatch(facei);
330 
331  if (mesh_.boundaryMesh()[patchi].coupled())
332  {
333  vector d = neiCc[facei-mesh_.nInternalFaces()] - ownCc;
334 
335  scalar dDotS = (d & s)/(mag(d)*magS + VSMALL);
336 
337  if (dDotS < severeNonorthogonalityThreshold)
338  {
339  return true;
340  }
341  else
342  {
343  return false;
344  }
345  }
346  else
347  {
348  // Collapsing normal boundary face does not cause problems with
349  // non-orthogonality
350  return false;
351  }
352  }
353 }
354 
355 
356 // Check if moving cell to new points causes it to collapse.
357 bool Foam::meshRefinement::isCollapsedCell
358 (
359  const pointField& points,
360  const scalar volFraction,
361  const label celli
362 ) const
363 {
364  scalar vol = mesh_.cells()[celli].mag(points, mesh_.faces());
365 
366  if (vol/mesh_.cellVolumes()[celli] < volFraction)
367  {
368  return true;
369  }
370  else
371  {
372  return false;
373  }
374 }
375 
376 
377 // Returns list with for every internal face -1 or the patch they should
378 // be baffled into. Gets run after createBaffles so all the unzoned surface
379 // intersections have already been turned into baffles. (Note: zoned surfaces
380 // are not baffled at this stage)
381 // Used to remove cells by baffling all their faces and have the
382 // splitMeshRegions chuck away non used regions.
383 void Foam::meshRefinement::markFacesOnProblemCells
384 (
385  const dictionary& motionDict,
386  const bool removeEdgeConnectedCells,
387  const scalarField& perpendicularAngle,
388  const labelList& globalToPatch,
389 
390  labelList& facePatch,
391  labelList& faceZone
392 ) const
393 {
394  const labelList& cellLevel = meshCutter_.cellLevel();
395  const labelList& pointLevel = meshCutter_.pointLevel();
396  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
397 
398 
399  // Mark all points and edges on baffle patches (so not on any inlets,
400  // outlets etc.)
401  boolList isBoundaryPoint(mesh_.nPoints(), false);
402  boolList isBoundaryEdge(mesh_.nEdges(), false);
403  boolList isBoundaryFace(mesh_.nFaces(), false);
404 
405  // Fill boundary data. All elements on meshed patches get marked.
406  // Get the labels of added patches.
407  const labelList adaptPatchIDs(meshedPatches());
408 
409  forAll(adaptPatchIDs, i)
410  {
411  const polyPatch& pp = patches[adaptPatchIDs[i]];
412 
413  label facei = pp.start();
414 
415  forAll(pp, j)
416  {
417  markBoundaryFace
418  (
419  facei,
420  isBoundaryFace,
421  isBoundaryEdge,
422  isBoundaryPoint
423  );
424 
425  facei++;
426  }
427  }
428 
429 
430  // Per mesh face the nearest adaptPatch and its faceZone (if any)
431  labelList nearestAdaptPatch;
432  labelList nearestAdaptZone;
433  nearestPatch(adaptPatchIDs, nearestAdaptPatch, nearestAdaptZone);
434 
435 
436  // Per internal face (boundary faces not used) the patch that the
437  // baffle should get (or -1)
438  facePatch.setSize(mesh_.nFaces());
439  facePatch = -1;
440  faceZone.setSize(mesh_.nFaces());
441  faceZone = -1;
442 
443  // Swap neighbouring cell centres and cell level
444  labelList neiLevel(mesh_.nBoundaryFaces());
445  pointField neiCc(mesh_.nBoundaryFaces());
446  calcNeighbourData(neiLevel, neiCc);
447 
448 
449  // Count of faces marked for baffling
450  label nBaffleFaces = 0;
451  bitSet isMasterFace(syncTools::getMasterFaces(mesh_));
452 
453  // Count of faces not baffled since would not cause a collapse
454  label nPrevented = 0;
455 
456  if (removeEdgeConnectedCells && max(perpendicularAngle) >= 0)
457  {
458  Info<< "markFacesOnProblemCells :"
459  << " Checking for edge-connected cells of highly differing sizes."
460  << endl;
461 
462  // Pick up the cells that need to be removed and (a guess for)
463  // the patch they should be patched with.
464  Map<label> problemCells
465  (
466  findEdgeConnectedProblemCells
467  (
468  perpendicularAngle,
469  globalToPatch
470  )
471  );
472 
473  // Baffle all faces of cells that need to be removed
474  forAllConstIters(problemCells, iter)
475  {
476  const cell& cFaces = mesh_.cells()[iter.key()];
477 
478  forAll(cFaces, i)
479  {
480  label facei = cFaces[i];
481 
482  if (facePatch[facei] == -1 && mesh_.isInternalFace(facei))
483  {
484  facePatch[facei] = nearestAdaptPatch[facei];
485  faceZone[facei] = nearestAdaptZone[facei];
486  nBaffleFaces++;
487 
488  // Mark face as a 'boundary'
489  markBoundaryFace
490  (
491  facei,
492  isBoundaryFace,
493  isBoundaryEdge,
494  isBoundaryPoint
495  );
496  }
497  }
498  }
499  Info<< "markFacesOnProblemCells : Marked "
500  << returnReduce(nBaffleFaces, sumOp<label>())
501  << " additional internal faces to be converted into baffles"
502  << " due to "
503  << returnReduce(problemCells.size(), sumOp<label>())
504  << " cells edge-connected to lower level cells." << endl;
505 
507  {
508  cellSet problemCellSet(mesh_, "problemCells", problemCells.toc());
509  problemCellSet.instance() = timeName();
510  Pout<< "Writing " << problemCellSet.size()
511  << " cells that are edge connected to coarser cell to set "
512  << problemCellSet.objectPath() << endl;
513  problemCellSet.write();
514  }
515  }
516 
518  (
519  mesh_,
520  isBoundaryPoint,
521  orEqOp<bool>(),
522  false // null value
523  );
524 
526  (
527  mesh_,
528  isBoundaryEdge,
529  orEqOp<bool>(),
530  false // null value
531  );
532 
534  (
535  mesh_,
536  isBoundaryFace,
537  orEqOp<bool>()
538  );
539 
540 
541  // See if checking for collapse
542  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
543 
544  // Collapse checking parameters
545  const scalar volFraction =
546  motionDict.getOrDefault<scalar>("minVolCollapseRatio", -1);
547 
548  const bool checkCollapse = (volFraction > 0);
549  scalar minArea = -1;
550  scalar maxNonOrtho = -1;
551 
552 
553  // Find nearest (non-baffle) surface
554  pointField newPoints;
555 
556  if (checkCollapse)
557  {
558  minArea = get<scalar>(motionDict, "minArea", dryRun_);
559  maxNonOrtho = get<scalar>(motionDict, "maxNonOrtho", dryRun_);
560 
561  Info<< "markFacesOnProblemCells :"
562  << " Deleting all-anchor surface cells only if"
563  << " snapping them violates mesh quality constraints:" << nl
564  << " snapped/original cell volume < " << volFraction << nl
565  << " face area < " << minArea << nl
566  << " non-orthogonality > " << maxNonOrtho << nl
567  << endl;
568 
569  // Construct addressing engine.
570  autoPtr<indirectPrimitivePatch> ppPtr
571  (
573  (
574  mesh_,
575  adaptPatchIDs
576  )
577  );
578  const indirectPrimitivePatch& pp = ppPtr();
579  const pointField& localPoints = pp.localPoints();
580  const labelList& meshPoints = pp.meshPoints();
581 
582  List<pointIndexHit> hitInfo;
583  labelList hitSurface;
584  surfaces_.findNearest
585  (
586  surfaceZonesInfo::getUnnamedSurfaces(surfaces_.surfZones()),
587  localPoints,
588  scalarField(localPoints.size(), sqr(GREAT)), // sqr of attraction
589  hitSurface,
590  hitInfo
591  );
592 
593  // Start off from current points
594  newPoints = mesh_.points();
595 
596  forAll(hitInfo, i)
597  {
598  if (hitInfo[i].hit())
599  {
600  newPoints[meshPoints[i]] = hitInfo[i].point();
601  }
602  }
603 
605  {
606  const_cast<Time&>(mesh_.time())++;
607  pointField oldPoints(mesh_.points());
608  mesh_.movePoints(newPoints);
609 
610  // Unset any moving state
611  mesh_.moving(false);
612 
613  Pout<< "Writing newPoints mesh to time " << timeName()
614  << endl;
615  write
616  (
617  debugType(debug),
618  writeType(writeLevel() | WRITEMESH),
619  mesh_.time().path()/"newPoints"
620  );
621  mesh_.movePoints(oldPoints);
622 
623  // Unset any moving state
624  mesh_.moving(false);
625  }
626  }
627 
628 
629 
630  // For each cell count the number of anchor points that are on
631  // the boundary:
632  // 8 : check the number of (baffle) boundary faces. If 3 or more block
633  // off the cell since the cell would get squeezed down to a diamond
634  // (probably; if the 3 or more faces are unrefined (only use the
635  // anchor points))
636  // 7 : store. Used to check later on whether there are points with
637  // 3 or more of these cells. (note that on a flat surface a boundary
638  // point will only have 4 cells connected to it)
639 
640 
641  // Does cell have exactly 7 of its 8 anchor points on the boundary?
642  bitSet hasSevenBoundaryAnchorPoints(mesh_.nCells());
643  // If so what is the remaining non-boundary anchor point?
644  labelHashSet nonBoundaryAnchors(mesh_.nCells()/10000);
645 
646  // On-the-fly addressing storage.
647  DynamicList<label> dynFEdges;
648  DynamicList<label> dynCPoints;
649  labelHashSet pSet;
650 
651  forAll(cellLevel, celli)
652  {
653  const labelList& cPoints = mesh_.cellPoints(celli, pSet, dynCPoints);
654 
655  // Get number of anchor points (pointLevel <= cellLevel)
656 
657  label nBoundaryAnchors = 0;
658  // label nNonAnchorBoundary = 0;
659  label nonBoundaryAnchor = -1;
660 
661  for (const label pointi : cPoints)
662  {
663  if (pointLevel[pointi] <= cellLevel[celli])
664  {
665  // Anchor point
666  if (isBoundaryPoint[pointi])
667  {
668  nBoundaryAnchors++;
669  }
670  else
671  {
672  // Anchor point which is not on the surface
673  nonBoundaryAnchor = pointi;
674  }
675  }
676  // else if (isBoundaryPoint[pointi])
677  // {
678  // ++nNonAnchorBoundary;
679  // }
680  }
681 
682  if (nBoundaryAnchors == 8)
683  {
684  const cell& cFaces = mesh_.cells()[celli];
685 
686  // Count boundary faces.
687  // label nBfaces = 0;
688  //
689  // forAll(cFaces, cFacei)
690  // {
691  // if (isBoundaryFace[cFaces[cFacei]])
692  // {
693  // ++nBfaces;
694  // }
695  // }
696 
697  // If nBfaces > 1 make all non-boundary non-baffle faces baffles.
698  // We assume that this situation is where there is a single
699  // cell sticking out which would get flattened.
700 
701  // Eugene: delete cell no matter what.
702  //if (nBfaces > 1)
703  {
704  if
705  (
706  checkCollapse
707  && !isCollapsedCell(newPoints, volFraction, celli)
708  )
709  {
710  nPrevented++;
711  //Pout<< "Preventing baffling/removal of 8 anchor point"
712  // << " cell "
713  // << celli << " at " << mesh_.cellCentres()[celli]
714  // << " since new volume "
715  // << mesh_.cells()[celli].mag(newPoints, mesh_.faces())
716  // << " old volume " << mesh_.cellVolumes()[celli]
717  // << endl;
718  }
719  else
720  {
721  // Block all faces of cell
722  forAll(cFaces, cf)
723  {
724  label facei = cFaces[cf];
725 
726  if
727  (
728  facePatch[facei] == -1
729  && mesh_.isInternalFace(facei)
730  )
731  {
732  facePatch[facei] = nearestAdaptPatch[facei];
733  faceZone[facei] = nearestAdaptZone[facei];
734  nBaffleFaces++;
735 
736  // Mark face as a 'boundary'
737  markBoundaryFace
738  (
739  facei,
740  isBoundaryFace,
741  isBoundaryEdge,
742  isBoundaryPoint
743  );
744  }
745  }
746  }
747  }
748  }
749  else if (nBoundaryAnchors == 7 && nonBoundaryAnchor != -1)
750  {
751  // Mark the cell. Store the (single!) non-boundary anchor point.
752  hasSevenBoundaryAnchorPoints.set(celli);
753  nonBoundaryAnchors.insert(nonBoundaryAnchor);
754  }
755  }
756 
757 
758  // Loop over all points. If a point is connected to 4 or more cells
759  // with 7 anchor points on the boundary set those cell's non-boundary faces
760  // to baffles
761 
762  DynamicList<label> dynPCells;
763 
764  for (const label pointi : nonBoundaryAnchors)
765  {
766  const labelList& pCells = mesh_.pointCells(pointi, dynPCells);
767 
768  // Count number of 'hasSevenBoundaryAnchorPoints' cells.
769  label n = 0;
770 
771  for (const label celli : pCells)
772  {
773  if (hasSevenBoundaryAnchorPoints.test(celli))
774  {
775  ++n;
776  }
777  }
778 
779  if (n > 3)
780  {
781  // Point in danger of being what? Remove all 7-cells.
782  for (const label celli : pCells)
783  {
784  if (hasSevenBoundaryAnchorPoints.test(celli))
785  {
786  if
787  (
788  checkCollapse
789  && !isCollapsedCell(newPoints, volFraction, celli)
790  )
791  {
792  nPrevented++;
793  //Pout<< "Preventing baffling of 7 anchor cell "
794  // << celli
795  // << " at " << mesh_.cellCentres()[celli]
796  // << " since new volume "
797  // << mesh_.cells()[celli].mag
798  // (newPoints, mesh_.faces())
799  // << " old volume " << mesh_.cellVolumes()[celli]
800  // << endl;
801  }
802  else
803  {
804  const cell& cFaces = mesh_.cells()[celli];
805 
806  for (const label facei : cFaces)
807  {
808  if
809  (
810  facePatch[facei] == -1
811  && mesh_.isInternalFace(facei)
812  )
813  {
814  facePatch[facei] = nearestAdaptPatch[facei];
815  faceZone[facei] = nearestAdaptZone[facei];
816  nBaffleFaces++;
817 
818  // Mark face as a 'boundary'
819  markBoundaryFace
820  (
821  facei,
822  isBoundaryFace,
823  isBoundaryEdge,
824  isBoundaryPoint
825  );
826  }
827  }
828  }
829  }
830  }
831  }
832  }
833 
834 
835  // Sync all. (note that pointdata and facedata not used anymore but sync
836  // anyway)
837 
839  (
840  mesh_,
841  isBoundaryPoint,
842  orEqOp<bool>(),
843  false // null value
844  );
845 
847  (
848  mesh_,
849  isBoundaryEdge,
850  orEqOp<bool>(),
851  false // null value
852  );
853 
855  (
856  mesh_,
857  isBoundaryFace,
858  orEqOp<bool>()
859  );
860 
861 
862  // Find faces with all edges on the boundary and make them baffles
863  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
864  {
865  if (facePatch[facei] == -1)
866  {
867  const labelList& fEdges = mesh_.faceEdges(facei, dynFEdges);
868  label nFaceBoundaryEdges = 0;
869 
870  for (const label edgei : fEdges)
871  {
872  if (isBoundaryEdge[edgei])
873  {
874  ++nFaceBoundaryEdges;
875  }
876  }
877 
878  if (nFaceBoundaryEdges == fEdges.size())
879  {
880  if
881  (
882  checkCollapse
883  && !isCollapsedFace
884  (
885  newPoints,
886  neiCc,
887  minArea,
888  maxNonOrtho,
889  facei
890  )
891  )
892  {
893  nPrevented++;
894  //Pout<< "Preventing baffling (to avoid collapse) of face "
895  // << facei
896  // << " with all boundary edges "
897  // << " at " << mesh_.faceCentres()[facei]
898  // << endl;
899  }
900  else
901  {
902  facePatch[facei] = nearestAdaptPatch[facei];
903  faceZone[facei] = nearestAdaptZone[facei];
904  nBaffleFaces++;
905 
906  // Do NOT update boundary data since this would grow blocked
907  // faces across gaps.
908  }
909  }
910  }
911  }
912 
913  for (const polyPatch& pp : patches)
914  {
915  if (pp.coupled())
916  {
917  label facei = pp.start();
918 
919  forAll(pp, i)
920  {
921  if (facePatch[facei] == -1)
922  {
923  const labelList& fEdges = mesh_.faceEdges(facei, dynFEdges);
924  label nFaceBoundaryEdges = 0;
925 
926  for (const label edgei : fEdges)
927  {
928  if (isBoundaryEdge[edgei])
929  {
930  ++nFaceBoundaryEdges;
931  }
932  }
933 
934  if (nFaceBoundaryEdges == fEdges.size())
935  {
936  if
937  (
938  checkCollapse
939  && !isCollapsedFace
940  (
941  newPoints,
942  neiCc,
943  minArea,
944  maxNonOrtho,
945  facei
946  )
947  )
948  {
949  nPrevented++;
950  //Pout<< "Preventing baffling of coupled face "
951  // << facei
952  // << " with all boundary edges "
953  // << " at " << mesh_.faceCentres()[facei]
954  // << endl;
955  }
956  else
957  {
958  facePatch[facei] = nearestAdaptPatch[facei];
959  faceZone[facei] = nearestAdaptZone[facei];
960  if (isMasterFace.test(facei))
961  {
962  ++nBaffleFaces;
963  }
964 
965  // Do NOT update boundary data since this would grow
966  // blocked faces across gaps.
967  }
968  }
969  }
970 
971  facei++;
972  }
973  }
974  }
975 
976 
977  // Because of isCollapsedFace one side can decide not to baffle whereas
978  // the other side does so sync. Baffling is preferred over not baffling.
979  if (checkCollapse) // Or always?
980  {
982  (
983  mesh_,
984  facePatch,
985  maxEqOp<label>()
986  );
988  (
989  mesh_,
990  faceZone,
991  maxEqOp<label>()
992  );
993  }
994 
995  Info<< "markFacesOnProblemCells : marked "
996  << returnReduce(nBaffleFaces, sumOp<label>())
997  << " additional internal faces to be converted into baffles."
998  << endl;
999 
1000  if (checkCollapse)
1001  {
1002  Info<< "markFacesOnProblemCells : prevented "
1003  << returnReduce(nPrevented, sumOp<label>())
1004  << " internal faces from getting converted into baffles."
1005  << endl;
1006  }
1007 }
1008 
1009 
1010 // Mark faces to be baffled to prevent snapping problems. Does
1011 // test to find nearest surface and checks which faces would get squashed.
1012 void Foam::meshRefinement::markFacesOnProblemCellsGeometric
1013 (
1014  const snapParameters& snapParams,
1015  const dictionary& motionDict,
1016  const labelList& globalToMasterPatch,
1017  const labelList& globalToSlavePatch,
1018 
1019  labelList& facePatch,
1020  labelList& faceZone
1021 ) const
1022 {
1023  pointField oldPoints(mesh_.points());
1024 
1025  // Repeat (most of) snappySnapDriver::doSnap
1026  {
1027  const labelList adaptPatchIDs(meshedPatches());
1028 
1029  // Construct addressing engine.
1030  autoPtr<indirectPrimitivePatch> ppPtr
1031  (
1033  (
1034  mesh_,
1035  adaptPatchIDs
1036  )
1037  );
1038  indirectPrimitivePatch& pp = ppPtr();
1039 
1040  // Distance to attract to nearest feature on surface
1041  const scalarField snapDist
1042  (
1043  snappySnapDriver::calcSnapDistance(mesh_, snapParams, pp)
1044  );
1045 
1046 
1047  // Construct iterative mesh mover.
1048  Info<< "Constructing mesh displacer ..." << endl;
1049  Info<< "Using mesh parameters " << motionDict << nl << endl;
1050 
1051  const pointMesh& pMesh = pointMesh::New(mesh_);
1052 
1053  motionSmoother meshMover
1054  (
1055  mesh_,
1056  pp,
1057  adaptPatchIDs,
1058  meshRefinement::makeDisplacementField(pMesh, adaptPatchIDs)(),
1059  motionDict
1060  );
1061 
1062 
1063  // Check initial mesh
1064  Info<< "Checking initial mesh ..." << endl;
1065  labelHashSet wrongFaces(mesh_.nFaces()/100);
1067  (
1068  false,
1069  mesh_,
1070  motionDict,
1071  wrongFaces,
1072  dryRun_
1073  );
1074  const label nInitErrors = returnReduce
1075  (
1076  wrongFaces.size(),
1077  sumOp<label>()
1078  );
1079 
1080  Info<< "Detected " << nInitErrors << " illegal faces"
1081  << " (concave, zero area or negative cell pyramid volume)"
1082  << endl;
1083 
1084 
1085  Info<< "Checked initial mesh in = "
1086  << mesh_.time().cpuTimeIncrement() << " s\n" << nl << endl;
1087 
1088  // Pre-smooth patch vertices (so before determining nearest)
1090  (
1091  *this,
1092  snapParams,
1093  nInitErrors,
1094  List<labelPair>(0), //baffles
1095  meshMover
1096  );
1097 
1098  pointField nearestPoint;
1099  vectorField nearestNormal;
1100  const vectorField disp
1101  (
1102  snappySnapDriver::calcNearestSurface
1103  (
1104  snapParams.strictRegionSnap(),
1105  *this,
1106  globalToMasterPatch,
1107  globalToSlavePatch,
1108  snapDist, // attraction
1109  pp,
1110  nearestPoint,
1111  nearestNormal
1112  )
1113  );
1114 
1115  const labelList& meshPoints = pp.meshPoints();
1116 
1117  pointField newPoints(mesh_.points());
1118  forAll(meshPoints, i)
1119  {
1120  newPoints[meshPoints[i]] += disp[i];
1121  }
1122 
1124  (
1125  mesh_,
1126  newPoints,
1127  minMagSqrEqOp<point>(), // combine op
1128  vector(GREAT, GREAT, GREAT) // null value (note: cannot use VGREAT)
1129  );
1130 
1131  mesh_.movePoints(newPoints);
1132 
1133  // Unset any moving state
1134  mesh_.moving(false);
1135  }
1136 
1137 
1138  // Per face the nearest adaptPatch
1139  //const labelList nearestAdaptPatch(nearestPatch(meshedPatches()));
1140  labelList nearestAdaptPatch;
1141  labelList nearestAdaptZone;
1142  nearestPatch(meshedPatches(), nearestAdaptPatch, nearestAdaptZone);
1143 
1144  // Per face (internal or coupled!) the patch that the
1145  // baffle should get (or -1).
1146  facePatch.setSize(mesh_.nFaces());
1147  facePatch = -1;
1148  faceZone.setSize(mesh_.nFaces());
1149  faceZone = -1;
1150  // Count of baffled faces
1151  label nBaffleFaces = 0;
1152 
1153  {
1154  faceSet wrongFaces(mesh_, "wrongFaces", 100);
1155  {
1156  //motionSmoother::checkMesh(false, mesh_, motionDict, wrongFaces);
1157 
1158  // Just check the errors from squashing
1159  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1160 
1161  const labelList allFaces(identity(mesh_.nFaces()));
1162  label nWrongFaces = 0;
1163 
1164  //const scalar minV
1165  //(motionDict.get<scalar>("minVol", keyType::REGEX_RECURSIVE));
1166  //if (minV > -GREAT)
1167  //{
1168  // polyMeshGeometry::checkFacePyramids
1169  // (
1170  // false,
1171  // minV,
1172  // mesh_,
1173  // mesh_.cellCentres(),
1174  // mesh_.points(),
1175  // allFaces,
1176  // List<labelPair>(0),
1177  // &wrongFaces
1178  // );
1179  //
1180  // label nNewWrongFaces = returnReduce
1181  // (
1182  // wrongFaces.size(),
1183  // sumOp<label>()
1184  // );
1185  //
1186  // Info<< " faces with pyramid volume < "
1187  // << setw(5) << minV
1188  // << " m^3 : "
1189  // << nNewWrongFaces-nWrongFaces << endl;
1190  //
1191  // nWrongFaces = nNewWrongFaces;
1192  //}
1193 
1194  scalar minArea(get<scalar>(motionDict, "minArea", dryRun_));
1195  if (minArea > -SMALL)
1196  {
1198  (
1199  false,
1200  minArea,
1201  mesh_,
1202  mesh_.faceAreas(),
1203  allFaces,
1204  &wrongFaces
1205  );
1206 
1207  label nNewWrongFaces = returnReduce
1208  (
1209  wrongFaces.size(),
1210  sumOp<label>()
1211  );
1212 
1213  Info<< " faces with area < "
1214  << setw(5) << minArea
1215  << " m^2 : "
1216  << nNewWrongFaces-nWrongFaces << endl;
1217 
1218  nWrongFaces = nNewWrongFaces;
1219  }
1220 
1221  scalar minDet(get<scalar>(motionDict, "minDeterminant", dryRun_));
1222  if (minDet > -1)
1223  {
1225  (
1226  false,
1227  minDet,
1228  mesh_,
1229  mesh_.faceAreas(),
1230  allFaces,
1231  polyMeshGeometry::affectedCells(mesh_, allFaces),
1232  &wrongFaces
1233  );
1234 
1235  label nNewWrongFaces = returnReduce
1236  (
1237  wrongFaces.size(),
1238  sumOp<label>()
1239  );
1240 
1241  Info<< " faces on cells with determinant < "
1242  << setw(5) << minDet << " : "
1243  << nNewWrongFaces-nWrongFaces << endl;
1244 
1245  nWrongFaces = nNewWrongFaces;
1246  }
1247  }
1248 
1249 
1250  for (const label facei : wrongFaces)
1251  {
1252  const label patchi = mesh_.boundaryMesh().whichPatch(facei);
1253 
1254  if (patchi == -1 || mesh_.boundaryMesh()[patchi].coupled())
1255  {
1256  facePatch[facei] = nearestAdaptPatch[facei];
1257  faceZone[facei] = nearestAdaptZone[facei];
1258  nBaffleFaces++;
1259 
1260  //Pout<< " " << facei
1261  // //<< " on patch " << mesh_.boundaryMesh()[patchi].name()
1262  // << " is destined for patch " << facePatch[facei]
1263  // << endl;
1264  }
1265  }
1266  }
1267 
1268 
1269  // Restore points.
1270  mesh_.movePoints(oldPoints);
1271 
1272  // Unset any moving state
1273  mesh_.moving(false);
1274 
1275 
1276  Info<< "markFacesOnProblemCellsGeometric : marked "
1277  << returnReduce(nBaffleFaces, sumOp<label>())
1278  << " additional internal and coupled faces"
1279  << " to be converted into baffles." << endl;
1280 
1282  (
1283  mesh_,
1284  facePatch,
1285  maxEqOp<label>()
1286  );
1288  (
1289  mesh_,
1290  faceZone,
1291  maxEqOp<label>()
1292  );
1293 }
1294 
1295 
1296 // ************************************************************************* //
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.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const Field< point_type > & localPoints() const
Return pointField of points in patch.
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
static void syncPointPositions(const polyMesh &mesh, List< point > &positions, const CombineOp &cop, const point &nullValue)
Synchronize locations on all mesh points.
Definition: syncTools.H:219
const labelListList & faceEdges() const
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new 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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:432
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field...
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
const labelList & meshPoints() const
Return labelList of mesh points in patch.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
word timeName
Definition: getTimeIndex.H:3
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dimensionedScalar cos(const dimensionedScalar &ds)
const pointField & points
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
const labelListList & edgeFaces() const
Return edge-face addressing.
static void preSmoothPatch(const meshRefinement &meshRefiner, const snapParameters &snapParams, const label nInitErrors, const List< labelPair > &baffles, motionSmoother &)
Smooth the mesh (patch and internal) to increase visibility.
Vector< scalar > vector
Definition: vector.H:57
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
static labelList getUnnamedSurfaces(const PtrList< surfaceZonesInfo > &surfList)
Get indices of unnamed surfaces (surfaces without faceZoneName)
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
static scalarField calcSnapDistance(const fvMesh &mesh, const snapParameters &snapParams, const indirectPrimitivePatch &)
Calculate edge length per patch point.
dimensionedScalar sin(const dimensionedScalar &ds)
Istream and Ostream manipulators taking arguments.
static bool checkFaceArea(const bool report, const scalar minArea, const polyMesh &, const vectorField &faceAreas, const labelList &checkFaces, labelHashSet *setPtr)
Small faces.
int debug
Static debugging option.
labelList f(nPoints)
static autoPtr< indirectPrimitivePatch > makePatch(const polyMesh &, const labelList &)
Create patch from set of patches.
vector point
Point is a vector.
Definition: point.H:37
const polyBoundaryMesh & patches
static labelList affectedCells(const polyMesh &, const labelList &changedFaces)
Helper function: get affected cells from faces.
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
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.
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.
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))
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
static bool checkCellDeterminant(const bool report, const scalar minDet, const polyMesh &, const vectorField &faceAreas, const labelList &checkFaces, const labelList &affectedCells, labelHashSet *setPtr)
Area of internal faces v.s. boundary faces.
List< bool > boolList
A List of bools.
Definition: List.H:60
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
A HashTable to objects of type <T> with a label key.