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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 "polyMeshGeometry.H"
37 #include "IOmanip.H"
38 #include "unitConversion.H"
39 #include "autoSnapDriver.H"
40 
41 #include "snapParameters.H"
42 #include "motionSmoother.H"
43 #include "topoDistanceData.H"
44 #include "FaceCellWave.H"
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
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 
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  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 
139 (
140  const scalarField& perpendicularAngle,
141  const labelList& globalToPatch
142 ) const
143 {
144  // Construct addressing engine from all patches added for meshing.
146  (
147  meshRefinement::makePatch
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  label face0 = pp.addressing()[eFaces[0]];
172  label face1 = pp.addressing()[eFaces[1]];
173 
174  label cell0 = mesh_.faceOwner()[face0];
175  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  label faceI = candidateFaces[i];
245 
246  vector n = mesh_.faceAreas()[faceI];
247  n /= mag(n);
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 
271  if (debug&meshRefinement::MESH)
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
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].normal(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.
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 
378 (
379  const labelList& adaptPatchIDs
380 ) const
381 {
382  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
383 
384  labelList nearestAdaptPatch;
385 
386  if (adaptPatchIDs.size())
387  {
388  nearestAdaptPatch.setSize(mesh_.nFaces(), adaptPatchIDs[0]);
389 
390 
391  // Count number of faces in adaptPatchIDs
392  label nFaces = 0;
393  forAll(adaptPatchIDs, i)
394  {
395  const polyPatch& pp = patches[adaptPatchIDs[i]];
396  nFaces += pp.size();
397  }
398 
399  // Field on cells and faces.
400  List<topoDistanceData> cellData(mesh_.nCells());
401  List<topoDistanceData> faceData(mesh_.nFaces());
402 
403  // Start of changes
404  labelList patchFaces(nFaces);
405  List<topoDistanceData> patchData(nFaces);
406  nFaces = 0;
407  forAll(adaptPatchIDs, i)
408  {
409  label patchI = adaptPatchIDs[i];
410  const polyPatch& pp = patches[patchI];
411 
412  forAll(pp, i)
413  {
414  patchFaces[nFaces] = pp.start()+i;
415  patchData[nFaces] = topoDistanceData(patchI, 0);
416  nFaces++;
417  }
418  }
419 
420  // Propagate information inwards
422  (
423  mesh_,
424  patchFaces,
425  patchData,
426  faceData,
427  cellData,
428  mesh_.globalData().nTotalCells()+1
429  );
430 
431  // And extract
432 
433  bool haveWarned = false;
434  forAll(faceData, faceI)
435  {
436  if (!faceData[faceI].valid(deltaCalc.data()))
437  {
438  if (!haveWarned)
439  {
441  << "Did not visit some faces, e.g. face " << faceI
442  << " at " << mesh_.faceCentres()[faceI] << endl
443  << "Assigning these cells to patch "
444  << adaptPatchIDs[0]
445  << endl;
446  haveWarned = true;
447  }
448  }
449  else
450  {
451  nearestAdaptPatch[faceI] = faceData[faceI].data();
452  }
453  }
454  }
455  else
456  {
457  // Use patch 0
458  nearestAdaptPatch.setSize(mesh_.nFaces(), 0);
459  }
460 
461  return nearestAdaptPatch;
462 }
463 
464 
465 // Returns list with for every internal face -1 or the patch they should
466 // be baffled into. Gets run after createBaffles so all the unzoned surface
467 // intersections have already been turned into baffles. (Note: zoned surfaces
468 // are not baffled at this stage)
469 // Used to remove cells by baffling all their faces and have the
470 // splitMeshRegions chuck away non used regions.
472 (
473  const dictionary& motionDict,
474  const bool removeEdgeConnectedCells,
475  const scalarField& perpendicularAngle,
476  const labelList& globalToPatch
477 ) const
478 {
479  const labelList& cellLevel = meshCutter_.cellLevel();
480  const labelList& pointLevel = meshCutter_.pointLevel();
481  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
482 
483 
484  // Mark all points and edges on baffle patches (so not on any inlets,
485  // outlets etc.)
486  boolList isBoundaryPoint(mesh_.nPoints(), false);
487  boolList isBoundaryEdge(mesh_.nEdges(), false);
488  boolList isBoundaryFace(mesh_.nFaces(), false);
489 
490  // Fill boundary data. All elements on meshed patches get marked.
491  // Get the labels of added patches.
492  labelList adaptPatchIDs(meshedPatches());
493 
494  forAll(adaptPatchIDs, i)
495  {
496  const polyPatch& pp = patches[adaptPatchIDs[i]];
497 
498  label faceI = pp.start();
499 
500  forAll(pp, j)
501  {
502  markBoundaryFace
503  (
504  faceI,
505  isBoundaryFace,
506  isBoundaryEdge,
507  isBoundaryPoint
508  );
509 
510  faceI++;
511  }
512  }
513 
514 
515  // Per face the nearest adaptPatch
516  const labelList nearestAdaptPatch(nearestPatch(adaptPatchIDs));
517 
518 
519  // Per internal face (boundary faces not used) the patch that the
520  // baffle should get (or -1)
521  labelList facePatch(mesh_.nFaces(), -1);
522 
523  // Swap neighbouring cell centres and cell level
524  labelList neiLevel(mesh_.nFaces()-mesh_.nInternalFaces());
525  pointField neiCc(mesh_.nFaces()-mesh_.nInternalFaces());
526  calcNeighbourData(neiLevel, neiCc);
527 
528 
529  // Count of faces marked for baffling
530  label nBaffleFaces = 0;
531  PackedBoolList isMasterFace(syncTools::getMasterFaces(mesh_));
532 
533  // Count of faces not baffled since would not cause a collapse
534  label nPrevented = 0;
535 
536  if (removeEdgeConnectedCells && max(perpendicularAngle) >= 0)
537  {
538  Info<< "markFacesOnProblemCells :"
539  << " Checking for edge-connected cells of highly differing sizes."
540  << endl;
541 
542  // Pick up the cells that need to be removed and (a guess for)
543  // the patch they should be patched with.
544  Map<label> problemCells
545  (
546  findEdgeConnectedProblemCells
547  (
548  perpendicularAngle,
549  globalToPatch
550  )
551  );
552 
553  // Baffle all faces of cells that need to be removed
554  forAllConstIter(Map<label>, problemCells, iter)
555  {
556  const cell& cFaces = mesh_.cells()[iter.key()];
557 
558  forAll(cFaces, i)
559  {
560  label faceI = cFaces[i];
561 
562  if (facePatch[faceI] == -1 && mesh_.isInternalFace(faceI))
563  {
564  facePatch[faceI] = nearestAdaptPatch[faceI];
565  nBaffleFaces++;
566 
567  // Mark face as a 'boundary'
568  markBoundaryFace
569  (
570  faceI,
571  isBoundaryFace,
572  isBoundaryEdge,
573  isBoundaryPoint
574  );
575  }
576  }
577  }
578  Info<< "markFacesOnProblemCells : Marked "
579  << returnReduce(nBaffleFaces, sumOp<label>())
580  << " additional internal faces to be converted into baffles"
581  << " due to "
582  << returnReduce(problemCells.size(), sumOp<label>())
583  << " cells edge-connected to lower level cells." << endl;
584 
585  if (debug&meshRefinement::MESH)
586  {
587  cellSet problemCellSet(mesh_, "problemCells", problemCells.toc());
588  problemCellSet.instance() = timeName();
589  Pout<< "Writing " << problemCellSet.size()
590  << " cells that are edge connected to coarser cell to set "
591  << problemCellSet.objectPath() << endl;
592  problemCellSet.write();
593  }
594  }
595 
596  syncTools::syncPointList
597  (
598  mesh_,
599  isBoundaryPoint,
600  orEqOp<bool>(),
601  false // null value
602  );
603 
604  syncTools::syncEdgeList
605  (
606  mesh_,
607  isBoundaryEdge,
608  orEqOp<bool>(),
609  false // null value
610  );
611 
612  syncTools::syncFaceList
613  (
614  mesh_,
615  isBoundaryFace,
616  orEqOp<bool>()
617  );
618 
619 
620  // See if checking for collapse
621  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
622 
623  // Collapse checking parameters
624  const scalar volFraction =
625  motionDict.lookupOrDefault<scalar>("minVolCollapseRatio", -1);
626 
627  const bool checkCollapse = (volFraction > 0);
628  scalar minArea = -1;
629  scalar maxNonOrtho = -1;
630 
631 
632  // Find nearest (non-baffle) surface
633  pointField newPoints;
634 
635  if (checkCollapse)
636  {
637  minArea = readScalar(motionDict.lookup("minArea"));
638  maxNonOrtho = readScalar(motionDict.lookup("maxNonOrtho"));
639 
640  Info<< "markFacesOnProblemCells :"
641  << " Deleting all-anchor surface cells only if"
642  << " snapping them violates mesh quality constraints:" << nl
643  << " snapped/original cell volume < " << volFraction << nl
644  << " face area < " << minArea << nl
645  << " non-orthogonality > " << maxNonOrtho << nl
646  << endl;
647 
648  // Construct addressing engine.
650  (
651  meshRefinement::makePatch
652  (
653  mesh_,
654  adaptPatchIDs
655  )
656  );
657  const indirectPrimitivePatch& pp = ppPtr();
658  const pointField& localPoints = pp.localPoints();
659  const labelList& meshPoints = pp.meshPoints();
660 
661  List<pointIndexHit> hitInfo;
662  labelList hitSurface;
663  surfaces_.findNearest
664  (
665  surfaceZonesInfo::getUnnamedSurfaces(surfaces_.surfZones()),
666  localPoints,
667  scalarField(localPoints.size(), sqr(GREAT)), // sqr of attraction
668  hitSurface,
669  hitInfo
670  );
671 
672  // Start off from current points
673  newPoints = mesh_.points();
674 
675  forAll(hitInfo, i)
676  {
677  if (hitInfo[i].hit())
678  {
679  newPoints[meshPoints[i]] = hitInfo[i].hitPoint();
680  }
681  }
682 
683  if (debug&meshRefinement::MESH)
684  {
685  const_cast<Time&>(mesh_.time())++;
686  pointField oldPoints(mesh_.points());
687  mesh_.movePoints(newPoints);
688  Pout<< "Writing newPoints mesh to time " << timeName()
689  << endl;
690  write
691  (
692  debugType(debug),
693  writeType(writeLevel() | WRITEMESH),
694  mesh_.time().path()/"newPoints"
695  );
696  mesh_.movePoints(oldPoints);
697  }
698  }
699 
700 
701 
702  // For each cell count the number of anchor points that are on
703  // the boundary:
704  // 8 : check the number of (baffle) boundary faces. If 3 or more block
705  // off the cell since the cell would get squeezed down to a diamond
706  // (probably; if the 3 or more faces are unrefined (only use the
707  // anchor points))
708  // 7 : store. Used to check later on whether there are points with
709  // 3 or more of these cells. (note that on a flat surface a boundary
710  // point will only have 4 cells connected to it)
711 
712 
713  // Does cell have exactly 7 of its 8 anchor points on the boundary?
714  PackedBoolList hasSevenBoundaryAnchorPoints(mesh_.nCells());
715  // If so what is the remaining non-boundary anchor point?
716  labelHashSet nonBoundaryAnchors(mesh_.nCells()/10000);
717 
718  // On-the-fly addressing storage.
719  DynamicList<label> dynFEdges;
720  DynamicList<label> dynCPoints;
721 
722  forAll(cellLevel, cellI)
723  {
724  const labelList& cPoints = mesh_.cellPoints(cellI, dynCPoints);
725 
726  // Get number of anchor points (pointLevel <= cellLevel)
727 
728  label nBoundaryAnchors = 0;
729  label nNonAnchorBoundary = 0;
730  label nonBoundaryAnchor = -1;
731 
732  forAll(cPoints, i)
733  {
734  label pointI = cPoints[i];
735 
736  if (pointLevel[pointI] <= cellLevel[cellI])
737  {
738  // Anchor point
739  if (isBoundaryPoint[pointI])
740  {
741  nBoundaryAnchors++;
742  }
743  else
744  {
745  // Anchor point which is not on the surface
746  nonBoundaryAnchor = pointI;
747  }
748  }
749  else if (isBoundaryPoint[pointI])
750  {
751  nNonAnchorBoundary++;
752  }
753  }
754 
755  if (nBoundaryAnchors == 8)
756  {
757  const cell& cFaces = mesh_.cells()[cellI];
758 
759  // Count boundary faces.
760  label nBfaces = 0;
761 
762  forAll(cFaces, cFaceI)
763  {
764  if (isBoundaryFace[cFaces[cFaceI]])
765  {
766  nBfaces++;
767  }
768  }
769 
770  // If nBfaces > 1 make all non-boundary non-baffle faces baffles.
771  // We assume that this situation is where there is a single
772  // cell sticking out which would get flattened.
773 
774  // Eugene: delete cell no matter what.
775  //if (nBfaces > 1)
776  {
777  if
778  (
779  checkCollapse
780  && !isCollapsedCell(newPoints, volFraction, cellI)
781  )
782  {
783  nPrevented++;
784  //Pout<< "Preventing baffling/removal of 8 anchor point"
785  // << " cell "
786  // << cellI << " at " << mesh_.cellCentres()[cellI]
787  // << " since new volume "
788  // << mesh_.cells()[cellI].mag(newPoints, mesh_.faces())
789  // << " old volume " << mesh_.cellVolumes()[cellI]
790  // << endl;
791  }
792  else
793  {
794  // Block all faces of cell
795  forAll(cFaces, cf)
796  {
797  label faceI = cFaces[cf];
798 
799  if
800  (
801  facePatch[faceI] == -1
802  && mesh_.isInternalFace(faceI)
803  )
804  {
805  facePatch[faceI] = nearestAdaptPatch[faceI];
806  nBaffleFaces++;
807 
808  // Mark face as a 'boundary'
809  markBoundaryFace
810  (
811  faceI,
812  isBoundaryFace,
813  isBoundaryEdge,
814  isBoundaryPoint
815  );
816  }
817  }
818  }
819  }
820  }
821  else if (nBoundaryAnchors == 7)
822  {
823  // Mark the cell. Store the (single!) non-boundary anchor point.
824  hasSevenBoundaryAnchorPoints.set(cellI, 1u);
825  nonBoundaryAnchors.insert(nonBoundaryAnchor);
826  }
827  }
828 
829 
830  // Loop over all points. If a point is connected to 4 or more cells
831  // with 7 anchor points on the boundary set those cell's non-boundary faces
832  // to baffles
833 
834  DynamicList<label> dynPCells;
835 
836  forAllConstIter(labelHashSet, nonBoundaryAnchors, iter)
837  {
838  label pointI = iter.key();
839 
840  const labelList& pCells = mesh_.pointCells(pointI, dynPCells);
841 
842  // Count number of 'hasSevenBoundaryAnchorPoints' cells.
843  label n = 0;
844 
845  forAll(pCells, i)
846  {
847  if (hasSevenBoundaryAnchorPoints.get(pCells[i]) == 1u)
848  {
849  n++;
850  }
851  }
852 
853  if (n > 3)
854  {
855  // Point in danger of being what? Remove all 7-cells.
856  forAll(pCells, i)
857  {
858  label cellI = pCells[i];
859 
860  if (hasSevenBoundaryAnchorPoints.get(cellI) == 1u)
861  {
862  if
863  (
864  checkCollapse
865  && !isCollapsedCell(newPoints, volFraction, cellI)
866  )
867  {
868  nPrevented++;
869  //Pout<< "Preventing baffling of 7 anchor cell "
870  // << cellI
871  // << " at " << mesh_.cellCentres()[cellI]
872  // << " since new volume "
873  // << mesh_.cells()[cellI].mag
874  // (newPoints, mesh_.faces())
875  // << " old volume " << mesh_.cellVolumes()[cellI]
876  // << endl;
877  }
878  else
879  {
880  const cell& cFaces = mesh_.cells()[cellI];
881 
882  forAll(cFaces, cf)
883  {
884  label faceI = cFaces[cf];
885 
886  if
887  (
888  facePatch[faceI] == -1
889  && mesh_.isInternalFace(faceI)
890  )
891  {
892  facePatch[faceI] = nearestAdaptPatch[faceI];
893  nBaffleFaces++;
894 
895  // Mark face as a 'boundary'
896  markBoundaryFace
897  (
898  faceI,
899  isBoundaryFace,
900  isBoundaryEdge,
901  isBoundaryPoint
902  );
903  }
904  }
905  }
906  }
907  }
908  }
909  }
910 
911 
912  // Sync all. (note that pointdata and facedata not used anymore but sync
913  // anyway)
914 
915  syncTools::syncPointList
916  (
917  mesh_,
918  isBoundaryPoint,
919  orEqOp<bool>(),
920  false // null value
921  );
922 
923  syncTools::syncEdgeList
924  (
925  mesh_,
926  isBoundaryEdge,
927  orEqOp<bool>(),
928  false // null value
929  );
930 
931  syncTools::syncFaceList
932  (
933  mesh_,
934  isBoundaryFace,
935  orEqOp<bool>()
936  );
937 
938 
939  // Find faces with all edges on the boundary and make them baffles
940  for (label faceI = 0; faceI < mesh_.nInternalFaces(); faceI++)
941  {
942  if (facePatch[faceI] == -1)
943  {
944  const labelList& fEdges = mesh_.faceEdges(faceI, dynFEdges);
945  label nFaceBoundaryEdges = 0;
946 
947  forAll(fEdges, fe)
948  {
949  if (isBoundaryEdge[fEdges[fe]])
950  {
951  nFaceBoundaryEdges++;
952  }
953  }
954 
955  if (nFaceBoundaryEdges == fEdges.size())
956  {
957  if
958  (
959  checkCollapse
960  && !isCollapsedFace
961  (
962  newPoints,
963  neiCc,
964  minArea,
965  maxNonOrtho,
966  faceI
967  )
968  )
969  {
970  nPrevented++;
971  //Pout<< "Preventing baffling (to avoid collapse) of face "
972  // << faceI
973  // << " with all boundary edges "
974  // << " at " << mesh_.faceCentres()[faceI]
975  // << endl;
976  }
977  else
978  {
979  facePatch[faceI] = nearestAdaptPatch[faceI];
980  nBaffleFaces++;
981 
982  // Do NOT update boundary data since this would grow blocked
983  // faces across gaps.
984  }
985  }
986  }
987  }
988 
989  forAll(patches, patchI)
990  {
991  const polyPatch& pp = patches[patchI];
992 
993  if (pp.coupled())
994  {
995  label faceI = pp.start();
996 
997  forAll(pp, i)
998  {
999  if (facePatch[faceI] == -1)
1000  {
1001  const labelList& fEdges = mesh_.faceEdges(faceI, dynFEdges);
1002  label nFaceBoundaryEdges = 0;
1003 
1004  forAll(fEdges, fe)
1005  {
1006  if (isBoundaryEdge[fEdges[fe]])
1007  {
1008  nFaceBoundaryEdges++;
1009  }
1010  }
1011 
1012  if (nFaceBoundaryEdges == fEdges.size())
1013  {
1014  if
1015  (
1016  checkCollapse
1017  && !isCollapsedFace
1018  (
1019  newPoints,
1020  neiCc,
1021  minArea,
1022  maxNonOrtho,
1023  faceI
1024  )
1025  )
1026  {
1027  nPrevented++;
1028  //Pout<< "Preventing baffling of coupled face "
1029  // << faceI
1030  // << " with all boundary edges "
1031  // << " at " << mesh_.faceCentres()[faceI]
1032  // << endl;
1033  }
1034  else
1035  {
1036  facePatch[faceI] = nearestAdaptPatch[faceI];
1037  if (isMasterFace[faceI])
1038  {
1039  nBaffleFaces++;
1040  }
1041 
1042  // Do NOT update boundary data since this would grow
1043  // blocked faces across gaps.
1044  }
1045  }
1046  }
1047 
1048  faceI++;
1049  }
1050  }
1051  }
1052 
1053 
1054  // Because of isCollapsedFace one side can decide not to baffle whereas
1055  // the other side does so sync. Baffling is prefered over not baffling.
1056  if (checkCollapse) // Or always?
1057  {
1058  syncTools::syncFaceList
1059  (
1060  mesh_,
1061  facePatch,
1062  maxEqOp<label>()
1063  );
1064  }
1065 
1066  Info<< "markFacesOnProblemCells : marked "
1067  << returnReduce(nBaffleFaces, sumOp<label>())
1068  << " additional internal faces to be converted into baffles."
1069  << endl;
1070 
1071  if (checkCollapse)
1072  {
1073  Info<< "markFacesOnProblemCells : prevented "
1074  << returnReduce(nPrevented, sumOp<label>())
1075  << " internal faces from getting converted into baffles."
1076  << endl;
1077  }
1078 
1079  return facePatch;
1080 }
1081 
1082 
1083 // Mark faces to be baffled to prevent snapping problems. Does
1084 // test to find nearest surface and checks which faces would get squashed.
1087  const snapParameters& snapParams,
1088  const dictionary& motionDict,
1089  const labelList& globalToMasterPatch,
1090  const labelList& globalToSlavePatch
1091 ) const
1092 {
1093  pointField oldPoints(mesh_.points());
1094 
1095  // Repeat (most of) autoSnapDriver::doSnap
1096  {
1097  labelList adaptPatchIDs(meshedPatches());
1098 
1099  // Construct addressing engine.
1101  (
1102  meshRefinement::makePatch
1103  (
1104  mesh_,
1105  adaptPatchIDs
1106  )
1107  );
1108  indirectPrimitivePatch& pp = ppPtr();
1109 
1110  // Distance to attract to nearest feature on surface
1111  const scalarField snapDist
1112  (
1113  autoSnapDriver::calcSnapDistance(mesh_, snapParams, pp)
1114  );
1115 
1116 
1117  // Construct iterative mesh mover.
1118  Info<< "Constructing mesh displacer ..." << endl;
1119  Info<< "Using mesh parameters " << motionDict << nl << endl;
1120 
1121  const pointMesh& pMesh = pointMesh::New(mesh_);
1122 
1123  motionSmoother meshMover
1124  (
1125  mesh_,
1126  pp,
1127  adaptPatchIDs,
1128  meshRefinement::makeDisplacementField(pMesh, adaptPatchIDs)(),
1129  motionDict
1130  );
1131 
1132 
1133  // Check initial mesh
1134  Info<< "Checking initial mesh ..." << endl;
1135  labelHashSet wrongFaces(mesh_.nFaces()/100);
1136  motionSmoother::checkMesh(false, mesh_, motionDict, wrongFaces);
1137  const label nInitErrors = returnReduce
1138  (
1139  wrongFaces.size(),
1140  sumOp<label>()
1141  );
1142 
1143  Info<< "Detected " << nInitErrors << " illegal faces"
1144  << " (concave, zero area or negative cell pyramid volume)"
1145  << endl;
1146 
1147 
1148  Info<< "Checked initial mesh in = "
1149  << mesh_.time().cpuTimeIncrement() << " s\n" << nl << endl;
1150 
1151  // Pre-smooth patch vertices (so before determining nearest)
1152  autoSnapDriver::preSmoothPatch
1153  (
1154  *this,
1155  snapParams,
1156  nInitErrors,
1157  List<labelPair>(0), //baffles
1158  meshMover
1159  );
1160 
1161  pointField nearestPoint;
1162  vectorField nearestNormal;
1163  const vectorField disp
1164  (
1165  autoSnapDriver::calcNearestSurface
1166  (
1167  snapParams.strictRegionSnap(),
1168  *this,
1169  globalToMasterPatch,
1170  globalToSlavePatch,
1171  snapDist, // attraction
1172  pp,
1173  nearestPoint,
1174  nearestNormal
1175  )
1176  );
1177 
1178  const labelList& meshPoints = pp.meshPoints();
1179 
1180  pointField newPoints(mesh_.points());
1181  forAll(meshPoints, i)
1182  {
1183  newPoints[meshPoints[i]] += disp[i];
1184  }
1185 
1186  syncTools::syncPointList
1187  (
1188  mesh_,
1189  newPoints,
1190  minMagSqrEqOp<point>(), // combine op
1191  vector(GREAT, GREAT, GREAT) // null value (note: cannot use VGREAT)
1192  );
1193 
1194  mesh_.movePoints(newPoints);
1195  }
1196 
1197 
1198  // Per face the nearest adaptPatch
1199  const labelList nearestAdaptPatch(nearestPatch(meshedPatches()));
1200 
1201  // Per face (internal or coupled!) the patch that the
1202  // baffle should get (or -1).
1203  labelList facePatch(mesh_.nFaces(), -1);
1204  // Count of baffled faces
1205  label nBaffleFaces = 0;
1206 
1207  {
1208  faceSet wrongFaces(mesh_, "wrongFaces", 100);
1209  {
1210  //motionSmoother::checkMesh(false, mesh_, motionDict, wrongFaces);
1211 
1212  // Just check the errors from squashing
1213  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1214 
1215  const labelList allFaces(identity(mesh_.nFaces()));
1216  label nWrongFaces = 0;
1217 
1218  //const scalar minV(readScalar(motionDict.lookup("minVol", true)));
1219  //if (minV > -GREAT)
1220  //{
1221  // polyMeshGeometry::checkFacePyramids
1222  // (
1223  // false,
1224  // minV,
1225  // mesh_,
1226  // mesh_.cellCentres(),
1227  // mesh_.points(),
1228  // allFaces,
1229  // List<labelPair>(0),
1230  // &wrongFaces
1231  // );
1232  //
1233  // label nNewWrongFaces = returnReduce
1234  // (
1235  // wrongFaces.size(),
1236  // sumOp<label>()
1237  // );
1238  //
1239  // Info<< " faces with pyramid volume < "
1240  // << setw(5) << minV
1241  // << " m^3 : "
1242  // << nNewWrongFaces-nWrongFaces << endl;
1243  //
1244  // nWrongFaces = nNewWrongFaces;
1245  //}
1246 
1247  scalar minArea(readScalar(motionDict.lookup("minArea")));
1248  if (minArea > -SMALL)
1249  {
1250  polyMeshGeometry::checkFaceArea
1251  (
1252  false,
1253  minArea,
1254  mesh_,
1255  mesh_.faceAreas(),
1256  allFaces,
1257  &wrongFaces
1258  );
1259 
1260  label nNewWrongFaces = returnReduce
1261  (
1262  wrongFaces.size(),
1263  sumOp<label>()
1264  );
1265 
1266  Info<< " faces with area < "
1267  << setw(5) << minArea
1268  << " m^2 : "
1269  << nNewWrongFaces-nWrongFaces << endl;
1270 
1271  nWrongFaces = nNewWrongFaces;
1272  }
1273 
1274  scalar minDet(readScalar(motionDict.lookup("minDeterminant")));
1275  if (minDet > -1)
1276  {
1277  polyMeshGeometry::checkCellDeterminant
1278  (
1279  false,
1280  minDet,
1281  mesh_,
1282  mesh_.faceAreas(),
1283  allFaces,
1284  polyMeshGeometry::affectedCells(mesh_, allFaces),
1285  &wrongFaces
1286  );
1287 
1288  label nNewWrongFaces = returnReduce
1289  (
1290  wrongFaces.size(),
1291  sumOp<label>()
1292  );
1293 
1294  Info<< " faces on cells with determinant < "
1295  << setw(5) << minDet << " : "
1296  << nNewWrongFaces-nWrongFaces << endl;
1297 
1298  nWrongFaces = nNewWrongFaces;
1299  }
1300  }
1301 
1302 
1303  forAllConstIter(faceSet, wrongFaces, iter)
1304  {
1305  label patchI = mesh_.boundaryMesh().whichPatch(iter.key());
1306 
1307  if (patchI == -1 || mesh_.boundaryMesh()[patchI].coupled())
1308  {
1309  facePatch[iter.key()] = nearestAdaptPatch[iter.key()];
1310  nBaffleFaces++;
1311 
1312  //Pout<< " " << iter.key()
1313  // //<< " on patch " << mesh_.boundaryMesh()[patchI].name()
1314  // << " is destined for patch " << facePatch[iter.key()]
1315  // << endl;
1316  }
1317  }
1318  }
1319 
1320 
1321  // Restore points.
1322  mesh_.movePoints(oldPoints);
1323 
1324 
1325  Info<< "markFacesOnProblemCellsGeometric : marked "
1326  << returnReduce(nBaffleFaces, sumOp<label>())
1327  << " additional internal and coupled faces"
1328  << " to be converted into baffles." << endl;
1329 
1330  syncTools::syncFaceList
1331  (
1332  mesh_,
1333  facePatch,
1334  maxEqOp<label>()
1335  );
1336 
1337  return facePatch;
1338 }
1339 
1340 
1341 // ************************************************************************* //
topoDistanceData.H
Foam::PackedBoolList
A bit-packed bool list.
Definition: PackedBoolList.H:63
Foam::topoDistanceData
For use with FaceCellWave. Determines topological distance to starting faces.
Definition: topoDistanceData.H:53
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Foam::meshRefinement::isCollapsedCell
bool isCollapsedCell(const pointField &, const scalar volFraction, const label cellI) const
Definition: meshRefinementProblemCells.C:358
Foam::PrimitivePatch::edgeFaces
const labelListList & edgeFaces() const
Return edge-face addressing.
Definition: PrimitivePatchTemplate.C:292
Foam::meshRefinement::isCollapsedFace
bool isCollapsedFace(const pointField &, const pointField &neiCc, const scalar minFaceArea, const scalar maxNonOrtho, const label faceI) const
Definition: meshRefinementProblemCells.C:287
Foam::polyMeshGenChecks::checkMesh
bool checkMesh(const polyMeshGen &mesh, const bool report)
Check mesh for correctness. Returns false for no error.
Definition: polyMeshGenChecks.C:104
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
polyMeshGeometry.H
Foam::PrimitivePatch::localPoints
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatchTemplate.C:432
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:255
Foam::PackedBoolList::set
void set(const PackedList< 1 > &)
Set specified bits.
Definition: PackedBoolList.C:169
Foam::meshRefinement::markBoundaryFace
void markBoundaryFace(const label faceI, boolList &isBoundaryFace, boolList &isBoundaryEdge, boolList &isBoundaryPoint) const
Helper function to mark face as being on 'boundary'. Used by.
Definition: meshRefinementProblemCells.C:49
Foam::polyPatch::coupled
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:322
motionSmoother.H
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::maxEqOp
Definition: ops.H:77
Foam::FaceCellWave::data
const TrackingData & data() const
Additional data to be passed into container.
Definition: FaceCellWave.H:331
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Definition: dictionaryTemplates.C:33
Foam::Map
A HashTable to objects of type <T> with a label key.
Definition: PrimitivePatchTemplate.H:68
unitConversion.H
Unit conversion functions.
Foam::regIOobject::write
virtual bool write() const
Write using setting from DB.
Definition: regIOobjectWrite.C:126
Foam::faceSet
A list of face labels.
Definition: faceSet.H:48
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::meshRefinement::markFacesOnProblemCellsGeometric
labelList markFacesOnProblemCellsGeometric(const snapParameters &snapParams, const dictionary &motionDict, const labelList &globalToMasterPatch, const labelList &globalToSlavePatch) const
Returns list with for every internal face -1 or the patch.
Definition: meshRefinementProblemCells.C:1086
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::IOobject::instance
const fileName & instance() const
Definition: IOobject.H:350
Foam::HashSet< label, Hash< label > >
patchFaces
labelList patchFaces(const polyBoundaryMesh &patches, const wordList &names)
Definition: extrudeMesh.C:148
syncTools.H
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::meshRefinement::nearestPatch
labelList nearestPatch(const labelList &adaptPatchIDs) const
Returns list with for every face the label of the nearest.
Definition: meshRefinementProblemCells.C:378
Foam::IOobject::objectPath
fileName objectPath() const
Return complete path + object name.
Definition: IOobject.H:376
searchableSurfaces.H
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::snapParameters
Simple container to keep together snap specific information.
Definition: snapParameters.H:50
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
refinementSurfaces.H
MESH
@ MESH
Definition: extrudeMesh.C:60
faceSet.H
Foam::snapParameters::strictRegionSnap
Switch strictRegionSnap() const
Attract point to corresponding surface region only.
Definition: snapParameters.H:165
nBfaces
label nBfaces
Definition: readKivaGrid.H:45
Foam::meshRefinement::findEdgeConnectedProblemCells
Map< label > findEdgeConnectedProblemCells(const scalarField &perpendicularAngle, const labelList &) const
Definition: meshRefinementProblemCells.C:139
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::motionSmoother
Given a displacement moves the mesh by scaling the displacement back until there are no more mesh err...
Definition: motionSmoother.H:89
Foam::DynamicList::shrink
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:258
Foam::orEqOp
Definition: ops.H:82
Foam::identity
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::HashTable::size
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
meshRefinement.H
fvMesh.H
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::cellSet
A collection of cell labels.
Definition: cellSet.H:48
Foam::PackedList::get
unsigned int get(const label) const
Get value at index I.
Definition: PackedListI.H:964
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
s
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))
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
readScalar
#define readScalar
Definition: doubleScalar.C:38
autoSnapDriver.H
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::FaceCellWave
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:75
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
nearestPatch
labelList nearestPatch(const polyMesh &mesh, const labelList &patchIDs)
Definition: subsetMesh.C:48
Foam::sumOp
Definition: ops.H:162
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::PrimitivePatch::faceNormals
const Field< PointType > & faceNormals() const
Return face normals for patch.
Definition: PrimitivePatchTemplate.C:520
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::meshRefinement::writeType
writeType
Definition: meshRefinement.H:130
Foam::List
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
scalarField
volScalarField scalarField(fieldObject, mesh)
Foam::meshRefinement::findNearest
void findNearest(const labelList &meshFaces, List< pointIndexHit > &nearestInfo, labelList &nearestSurface, labelList &nearestRegion, vectorField &nearestNormal) const
Definition: meshRefinementProblemCells.C:75
FaceCellWave.H
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
timeName
word timeName
Definition: getTimeIndex.H:3
snapParameters.H
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::minMagSqrEqOp
Definition: ops.H:79
cellSet.H
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatchTemplate.C:392
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
write
Tcoeff write()
Foam::meshRefinement::markFacesOnProblemCells
labelList markFacesOnProblemCells(const dictionary &motionDict, const bool removeEdgeConnectedCells, const scalarField &perpendicularAngle, const labelList &globalToMasterPatch) const
Returns list with for every internal face -1 or the patch.
Definition: meshRefinementProblemCells.C:472
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88
Foam::degToRad
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
Definition: unitConversion.H:45
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:256
pointSet.H
Foam::meshRefinement::debugType
debugType
Definition: meshRefinement.H:100