selectCells.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-2022 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 Application
25  selectCells
26 
27 Description
28  Select cells in relation to surface.
29 
30  Divides cells into three sets:
31  - cutCells : cells cut by surface or close to surface.
32  - outside : cells not in cutCells and reachable from set of
33  user-defined points (outsidePoints)
34  - inside : same but not reachable.
35 
36  Finally the wanted sets are combined into a cellSet 'selected'. Apart
37  from straightforward adding the contents there are a few extra rules to
38  make sure that the surface of the 'outside' of the mesh is singly
39  connected.
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #include "argList.H"
44 #include "Time.H"
45 #include "polyMesh.H"
46 #include "IOdictionary.H"
47 #include "twoDPointCorrector.H"
48 #include "OFstream.H"
49 #include "meshTools.H"
50 
51 #include "triSurface.H"
52 #include "triSurfaceSearch.H"
53 #include "meshSearch.H"
54 #include "cellClassification.H"
55 #include "cellSet.H"
56 #include "cellInfo.H"
57 #include "FaceCellWave.H"
58 #include "edgeStats.H"
59 #include "treeDataTriSurface.H"
60 #include "indexedOctree.H"
61 #include "globalMeshData.H"
62 
63 using namespace Foam;
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 // cellType for cells included/not included in mesh.
68 static const label MESH = cellClassification::INSIDE;
69 static const label NONMESH = cellClassification::OUTSIDE;
70 
71 
72 void writeSet(const cellSet& cells, const string& msg)
73 {
74  Info<< "Writing " << msg << " (" << cells.size() << ") to cellSet "
75  << cells.instance()/cells.local()/cells.name()
76  << endl << endl;
77  cells.write();
78 }
79 
80 
81 void getType(const labelList& elems, const label type, labelHashSet& set)
82 {
83  forAll(elems, i)
84  {
85  if (elems[i] == type)
86  {
87  set.insert(i);
88  }
89  }
90 }
91 
92 
93 void cutBySurface
94 (
95  const polyMesh& mesh,
96  const meshSearch& queryMesh,
97  const triSurfaceSearch& querySurf,
98 
99  const pointField& outsidePts,
100  const bool selectCut,
101  const bool selectInside,
102  const bool selectOutside,
103  const scalar nearDist,
104 
105  cellClassification& cellType
106 )
107 {
108  // Cut with surface and classify as inside/outside/cut
109  cellType =
111  (
112  mesh,
113  queryMesh,
114  querySurf,
115  outsidePts
116  );
117 
118  // Get inside/outside/cutCells cellSets.
119  cellSet inside(mesh, "inside", mesh.nCells()/10);
120  getType(cellType, cellClassification::INSIDE, inside);
121  writeSet(inside, "inside cells");
122 
123  cellSet outside(mesh, "outside", mesh.nCells()/10);
124  getType(cellType, cellClassification::OUTSIDE, outside);
125  writeSet(outside, "outside cells");
126 
127  cellSet cutCells(mesh, "cutCells", mesh.nCells()/10);
128  getType(cellType, cellClassification::CUT, cutCells);
129  writeSet(cutCells, "cells cut by surface");
130 
131 
132  // Change cellType to reflect selected part of mesh. Use
133  // MESH to denote selected part, NONMESH for all
134  // other cells.
135  // Is a bit of a hack but allows us to reuse all the functionality
136  // in cellClassification.
137 
138  forAll(cellType, celli)
139  {
140  label cType = cellType[celli];
141 
142  if (cType == cellClassification::CUT)
143  {
144  if (selectCut)
145  {
146  cellType[celli] = MESH;
147  }
148  else
149  {
150  cellType[celli] = NONMESH;
151  }
152  }
153  else if (cType == cellClassification::INSIDE)
154  {
155  if (selectInside)
156  {
157  cellType[celli] = MESH;
158  }
159  else
160  {
161  cellType[celli] = NONMESH;
162  }
163  }
164  else if (cType == cellClassification::OUTSIDE)
165  {
166  if (selectOutside)
167  {
168  cellType[celli] = MESH;
169  }
170  else
171  {
172  cellType[celli] = NONMESH;
173  }
174  }
175  else
176  {
178  << "Multiple mesh regions in original mesh" << endl
179  << "Please use splitMeshRegions to separate these"
180  << exit(FatalError);
181  }
182  }
183 
184 
185  if (nearDist > 0)
186  {
187  Info<< "Removing cells with points closer than " << nearDist
188  << " to the surface ..." << nl << endl;
189 
190  const pointField& pts = mesh.points();
191  const indexedOctree<treeDataTriSurface>& tree = querySurf.tree();
192 
193  label nRemoved = 0;
194 
195  forAll(pts, pointi)
196  {
197  const point& pt = pts[pointi];
198 
199  pointIndexHit hitInfo = tree.findNearest(pt, sqr(nearDist));
200 
201  if (hitInfo.hit())
202  {
203  const labelList& pCells = mesh.pointCells()[pointi];
204 
205  forAll(pCells, i)
206  {
207  if (cellType[pCells[i]] != NONMESH)
208  {
209  cellType[pCells[i]] = NONMESH;
210  nRemoved++;
211  }
212  }
213  }
214  }
215 
216 // tmp<pointField> tnearest = querySurf.calcNearest(pts);
217 // const pointField& nearest = tnearest();
218 //
219 // label nRemoved = 0;
220 //
221 // forAll(nearest, pointi)
222 // {
223 // if (mag(nearest[pointi] - pts[pointi]) < nearDist)
224 // {
225 // const labelList& pCells = mesh.pointCells()[pointi];
226 //
227 // forAll(pCells, i)
228 // {
229 // if (cellType[pCells[i]] != NONMESH)
230 // {
231 // cellType[pCells[i]] = NONMESH;
232 // nRemoved++;
233 // }
234 // }
235 // }
236 // }
237 
238  Info<< "Removed " << nRemoved << " cells since too close to surface"
239  << nl << endl;
240  }
241 }
242 
243 
244 
245 // We're meshing the outside. Subset the currently selected mesh cells with the
246 // ones reachable from the outsidepoints.
247 label selectOutsideCells
248 (
249  const polyMesh& mesh,
250  const meshSearch& queryMesh,
251  const pointField& outsidePts,
252  cellClassification& cellType
253 )
254 {
255  //
256  // Check all outsidePts and for all of them inside a mesh cell
257  // collect the faces to start walking from
258  //
259 
260  // Outside faces
261  labelHashSet outsideFacesMap(outsidePts.size() * 6 * 2);
262  DynamicList<label> outsideFaces(outsideFacesMap.size());
263  // CellInfo on outside faces
264  DynamicList<cellInfo> outsideFacesInfo(outsideFacesMap.size());
265 
266  // cellInfo for mesh cell
267  const cellInfo meshInfo(MESH);
268 
269  forAll(outsidePts, outsidePtI)
270  {
271  // Find cell containing point. Linear search.
272  label celli = queryMesh.findCell(outsidePts[outsidePtI], -1, false);
273 
274  if (celli != -1 && cellType[celli] == MESH)
275  {
276  Info<< "Marking cell " << celli << " containing outside point "
277  << outsidePts[outsidePtI] << " with type " << cellType[celli]
278  << " ..." << endl;
279 
280  //
281  // Mark this cell and its faces to start walking from
282  //
283 
284  // Mark faces of celli
285  const labelList& cFaces = mesh.cells()[celli];
286  forAll(cFaces, i)
287  {
288  label facei = cFaces[i];
289 
290  if (outsideFacesMap.insert(facei))
291  {
292  outsideFaces.append(facei);
293  outsideFacesInfo.append(meshInfo);
294  }
295  }
296  }
297  }
298 
299  // Floodfill starting from outsideFaces (of type meshInfo)
300  List<cellInfo> faceInfoList(mesh.nFaces()), cellInfoList(mesh.nCells());
301  FaceCellWave<cellInfo> regionCalc
302  (
303  mesh,
304  outsideFaces.shrink(),
305  outsideFacesInfo.shrink(),
306  faceInfoList,
307  cellInfoList,
308  mesh.globalData().nTotalCells() + 1 // max iterations
309  );
310 
311  // Now cellInfoList should hold info on cells that are reachable from
312  // changedFaces. Use these to subset cellType
313  label nChanged = 0;
314  forAll(cellInfoList, celli)
315  {
316  if (cellType[celli] == MESH)
317  {
318  // Original cell was selected for meshing. Check if cell was
319  // reached from outsidePoints
320  if (cellInfoList[celli].type() != MESH)
321  {
322  cellType[celli] = NONMESH;
323  nChanged++;
324  }
325  }
326  }
327 
328  return nChanged;
329 }
330 
331 
332 
333 int main(int argc, char *argv[])
334 {
336 
337  #include "setRootCase.H"
338  #include "createTime.H"
339  #include "createPolyMesh.H"
340 
341  // Mesh edge statistics calculator
342  edgeStats edgeCalc(mesh);
343 
344 
345  IOdictionary refineDict
346  (
347  IOobject
348  (
349  "selectCellsDict",
350  runTime.system(),
351  mesh,
354  )
355  );
356 
357  fileName surfName(refineDict.lookup("surface"));
358  pointField outsidePts(refineDict.lookup("outsidePoints"));
359  bool useSurface(readBool(refineDict.lookup("useSurface")));
360  bool selectCut(readBool(refineDict.lookup("selectCut")));
361  bool selectInside(readBool(refineDict.lookup("selectInside")));
362  bool selectOutside(readBool(refineDict.lookup("selectOutside")));
363  scalar nearDist(refineDict.lookup<scalar>("nearDistance"));
364 
365 
366  if (useSurface)
367  {
368  Info<< "Cells to be used for meshing (0=false, 1=true):" << nl
369  << " cells cut by surface : " << selectCut << nl
370  << " cells inside of surface : " << selectInside << nl
371  << " cells outside of surface : " << selectOutside << nl
372  << " cells with points further than : " << nearDist << nl
373  << endl;
374  }
375  else
376  {
377  Info<< "Cells to be used for meshing (0=false, 1=true):" << nl
378  << " cells reachable from outsidePoints:" << selectOutside << nl
379  << endl;
380  }
381 
382  // Print edge stats on original mesh.
383  (void)edgeCalc.minLen(Info);
384 
385  // Search engine on mesh. Face decomposition since faces might be warped.
386  meshSearch queryMesh(mesh);
387 
388  // Check all 'outside' points
389  forAll(outsidePts, outsideI)
390  {
391  const point& outsidePoint = outsidePts[outsideI];
392 
393  label celli = queryMesh.findCell(outsidePoint, -1, false);
394  if (returnReduce(celli, maxOp<label>()) == -1)
395  {
397  << "outsidePoint " << outsidePoint
398  << " is not inside any cell"
399  << exit(FatalError);
400  }
401  }
402 
403  // Cell status (compared to surface if provided): inside/outside/cut.
404  // Start off from everything selected and cut later.
405  cellClassification cellType
406  (
407  mesh,
408  labelList
409  (
410  mesh.nCells(),
412  )
413  );
414 
415 
416  // Surface
417  autoPtr<triSurface> surf(nullptr);
418  // Search engine on surface.
419  autoPtr<triSurfaceSearch> querySurf(nullptr);
420 
421  if (useSurface)
422  {
423  surf.reset(new triSurface(surfName));
424 
425  // Dump some stats
426  surf().writeStats(Info);
427 
428  // Search engine on surface.
429  querySurf.reset(new triSurfaceSearch(surf));
430 
431  // Set cellType[celli] according to relation to surface
432  cutBySurface
433  (
434  mesh,
435  queryMesh,
436  querySurf,
437 
438  outsidePts,
439  selectCut,
440  selectInside,
441  selectOutside,
442  nearDist,
443 
444  cellType
445  );
446  }
447 
448 
449  // Now 'trim' all the corners from the mesh so meshing/surface extraction
450  // becomes easier.
451 
452  label nHanging, nRegionEdges, nRegionPoints, nOutside;
453 
454  do
455  {
456  Info<< "Removing cells which after subsetting would have all points"
457  << " on outside ..." << nl << endl;
458 
459  nHanging = cellType.fillHangingCells
460  (
461  MESH, // meshType
462  NONMESH, // fill type
463  mesh.nCells()
464  );
465 
466 
467  Info<< "Removing edges connecting cells unconnected by faces ..."
468  << nl << endl;
469 
470  nRegionEdges = cellType.fillRegionEdges
471  (
472  MESH, // meshType
473  NONMESH, // fill type
474  mesh.nCells()
475  );
476 
477 
478  Info<< "Removing points connecting cells unconnected by faces ..."
479  << nl << endl;
480 
481  nRegionPoints = cellType.fillRegionPoints
482  (
483  MESH, // meshType
484  NONMESH, // fill type
485  mesh.nCells()
486  );
487 
488  nOutside = 0;
489  if (selectOutside)
490  {
491  // Since we're selecting the cells reachable from outsidePoints
492  // and the set might have changed, redo the outsideCells
493  // calculation
494  nOutside = selectOutsideCells
495  (
496  mesh,
497  queryMesh,
498  outsidePts,
499  cellType
500  );
501  }
502  } while
503  (
504  nHanging != 0
505  || nRegionEdges != 0
506  || nRegionPoints != 0
507  || nOutside != 0
508  );
509 
510  cellSet selectedCells(mesh, "selected", mesh.nCells()/10);
511  getType(cellType, MESH, selectedCells);
512 
513  writeSet(selectedCells, "cells selected for meshing");
514 
515 
516  Info<< "End\n" << endl;
517 
518  return 0;
519 }
520 
521 
522 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:79
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:119
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:54
bool hit() const
Is there a hit.
static void noParallel()
Remove the parallel options.
Definition: argList.C:175
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
'Cuts' a mesh with a surface.
label fillRegionPoints(const label meshType, const label fillType, const label maxIter)
Find regionPoints and fill all neighbours. Iterate until nothing.
label fillHangingCells(const label meshType, const label fillType, const label maxIter)
Find hanging cells (cells with all points on outside) and set their.
label fillRegionEdges(const label meshType, const label fillType, const label maxIter)
Find regionEdges and fill one neighbour. Iterate until nothing.
Holds information regarding type of cell. Used in inside/outside determination in cellClassification.
Definition: cellInfo.H:67
A collection of cell labels.
Definition: cellSet.H:51
Helper class to calculate minimum edge length on mesh.
Definition: edgeStats.H:55
A class for handling file names.
Definition: fileName.H:82
label nTotalCells() const
Return total number of cells in decomposed mesh.
Non-pointer based hierarchical recursive searching.
Definition: indexedOctree.H:72
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:58
label findCell(const point &location, const label seedCelli=-1, const bool useTreeSearch=true) const
Find cell containing location.
Definition: meshSearch.C:750
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1562
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1360
label nCells() const
const labelListList & pointCells() const
label nFaces() const
const cellList & cells() const
Helper class to search on triSurface.
const indexedOctree< treeDataTriSurface > & tree() const
Demand driven construction of the octree.
Triangulated surface description with patch information.
Definition: triSurface.H:69
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
int main(int argc, char *argv[])
Definition: financialFoam.C:44
const cellShapeList & cells
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
bool readBool(Istream &)
Definition: boolIO.C:66
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:251
messageStream Info
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
error FatalError
static const char nl
Definition: Ostream.H:260
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488