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