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