Loading...
Searching...
No Matches
selectCells.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Application
28 selectCells
29
30Group
31 grpMeshAdvancedUtilities
32
33Description
34 Select cells in relation to surface.
35
36 Divides cells into three sets:
37 - cutCells : cells cut by surface or close to surface.
38 - outside : cells not in cutCells and reachable from set of
39 user-defined points (outsidePoints)
40 - inside : same but not reachable.
41
42 Finally the wanted sets are combined into a cellSet 'selected'. Apart
43 from straightforward adding the contents there are a few extra rules to
44 make sure that the surface of the 'outside' of the mesh is singly
45 connected.
46
47\*---------------------------------------------------------------------------*/
48
49#include "argList.H"
50#include "Time.H"
51#include "polyMesh.H"
52#include "IOdictionary.H"
53#include "twoDPointCorrector.H"
54#include "OFstream.H"
55#include "meshTools.H"
56
57#include "triSurface.H"
58#include "triSurfaceSearch.H"
59#include "meshSearch.H"
60#include "cellClassification.H"
61#include "cellSet.H"
62#include "cellInfo.H"
63#include "MeshWave.H"
64#include "edgeStats.H"
65#include "treeDataTriSurface.H"
66#include "indexedOctree.H"
67#include "globalMeshData.H"
68
69using namespace Foam;
70
71// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72
73// cellType for cells included/not included in mesh.
74static const label MESH = cellClassification::INSIDE;
75static const label NONMESH = cellClassification::OUTSIDE;
76
77
78void writeSet(const cellSet& cells, const string& msg)
79{
80 Info<< "Writing " << msg << " (" << cells.size() << ") to cellSet "
81 << cells.instance()/cells.local()/cells.name()
82 << endl << endl;
83 cells.write();
84}
85
86
87void getType(const labelList& elems, const label type, labelHashSet& set)
88{
89 forAll(elems, i)
90 {
91 if (elems[i] == type)
92 {
93 set.insert(i);
94 }
95 }
96}
97
98
99void cutBySurface
100(
101 const polyMesh& mesh,
102 const meshSearch& queryMesh,
103 const triSurfaceSearch& querySurf,
104
105 const pointField& outsidePts,
106 const bool selectCut,
107 const bool selectInside,
108 const bool selectOutside,
109 const scalar nearDist,
110
111 cellClassification& cellType
112)
113{
114 // Cut with surface and classify as inside/outside/cut
115 cellType =
117 (
118 mesh,
119 queryMesh,
120 querySurf,
121 outsidePts
122 );
123
124 // Get inside/outside/cutCells cellSets.
125 cellSet inside(mesh, "inside", mesh.nCells()/10);
126 getType(cellType, cellClassification::INSIDE, inside);
127 writeSet(inside, "inside cells");
128
129 cellSet outside(mesh, "outside", mesh.nCells()/10);
130 getType(cellType, cellClassification::OUTSIDE, outside);
131 writeSet(outside, "outside cells");
132
133 cellSet cutCells(mesh, "cutCells", mesh.nCells()/10);
134 getType(cellType, cellClassification::CUT, cutCells);
135 writeSet(cutCells, "cells cut by surface");
136
137
138 // Change cellType to reflect selected part of mesh. Use
139 // MESH to denote selected part, NONMESH for all
140 // other cells.
141 // Is a bit of a hack but allows us to reuse all the functionality
142 // in cellClassification.
143
144 forAll(cellType, celli)
145 {
146 label cType = cellType[celli];
147
148 if (cType == cellClassification::CUT)
149 {
150 if (selectCut)
151 {
152 cellType[celli] = MESH;
153 }
154 else
155 {
156 cellType[celli] = NONMESH;
157 }
158 }
159 else if (cType == cellClassification::INSIDE)
160 {
161 if (selectInside)
162 {
163 cellType[celli] = MESH;
164 }
165 else
166 {
167 cellType[celli] = NONMESH;
168 }
169 }
170 else if (cType == cellClassification::OUTSIDE)
171 {
172 if (selectOutside)
173 {
174 cellType[celli] = MESH;
175 }
176 else
177 {
178 cellType[celli] = NONMESH;
179 }
180 }
181 else
182 {
184 << "Multiple mesh regions in original mesh" << endl
185 << "Please use splitMeshRegions to separate these"
186 << exit(FatalError);
187 }
188 }
189
190
191 if (nearDist > 0)
192 {
193 Info<< "Removing cells with points closer than " << nearDist
194 << " to the surface ..." << nl << endl;
195
196 const pointField& pts = mesh.points();
197 const indexedOctree<treeDataTriSurface>& tree = querySurf.tree();
198
199 label nRemoved = 0;
200
201 forAll(pts, pointi)
202 {
203 const point& pt = pts[pointi];
204
205 pointIndexHit hitInfo = tree.findNearest(pt, sqr(nearDist));
206
207 if (hitInfo.hit())
208 {
209 const labelList& pCells = mesh.pointCells()[pointi];
210
211 forAll(pCells, i)
212 {
213 if (cellType[pCells[i]] != NONMESH)
214 {
215 cellType[pCells[i]] = NONMESH;
216 nRemoved++;
217 }
218 }
219 }
220 }
221
222// tmp<pointField> tnearest = querySurf.calcNearest(pts);
223// const pointField& nearest = tnearest();
224//
225// label nRemoved = 0;
226//
227// forAll(nearest, pointi)
228// {
229// if (mag(nearest[pointi] - pts[pointi]) < nearDist)
230// {
231// const labelList& pCells = mesh.pointCells()[pointi];
232//
233// forAll(pCells, i)
234// {
235// if (cellType[pCells[i]] != NONMESH)
236// {
237// cellType[pCells[i]] = NONMESH;
238// nRemoved++;
239// }
240// }
241// }
242// }
243
244 Info<< "Removed " << nRemoved << " cells since too close to surface"
245 << nl << endl;
246 }
247}
248
249
250
251// We're meshing the outside. Subset the currently selected mesh cells with the
252// ones reachable from the outsidepoints.
253label selectOutsideCells
254(
255 const polyMesh& mesh,
256 const meshSearch& queryMesh,
257 const pointField& outsidePts,
258 cellClassification& cellType
259)
260{
261 //
262 // Check all outsidePts and for all of them inside a mesh cell
263 // collect the faces to start walking from
264 //
265
266 // Outside faces
267 labelHashSet outsideFacesMap(outsidePts.size() * 6 * 2);
268 DynamicList<label> outsideFaces(outsideFacesMap.size());
269 // CellInfo on outside faces
270 DynamicList<cellInfo> outsideFacesInfo(outsideFacesMap.size());
271
272 // cellInfo for mesh cell
273 const cellInfo meshInfo(MESH);
274
275 forAll(outsidePts, outsidePtI)
276 {
277 // Find cell containing point. Linear search.
278 label celli = queryMesh.findCell(outsidePts[outsidePtI], -1, false);
279
280 if (celli != -1 && cellType[celli] == MESH)
281 {
282 Info<< "Marking cell " << celli << " containing outside point "
283 << outsidePts[outsidePtI] << " with type " << cellType[celli]
284 << " ..." << endl;
285
286 //
287 // Mark this cell and its faces to start walking from
288 //
289
290 // Mark faces of celli
291 const labelList& cFaces = mesh.cells()[celli];
292 forAll(cFaces, i)
293 {
294 label facei = cFaces[i];
295
296 if (outsideFacesMap.insert(facei))
297 {
298 outsideFaces.append(facei);
299 outsideFacesInfo.append(meshInfo);
300 }
301 }
302 }
303 }
304
305 // Floodfill starting from outsideFaces (of type meshInfo)
306 MeshWave<cellInfo> regionCalc
307 (
308 mesh,
309 outsideFaces.shrink(),
310 outsideFacesInfo.shrink(),
311 mesh.globalData().nTotalCells()+1 // max iterations
312 );
313
314 // Now regionCalc should hold info on cells that are reachable from
315 // changedFaces. Use these to subset cellType
316 const List<cellInfo>& allCellInfo = regionCalc.allCellInfo();
317
318 label nChanged = 0;
319
320 forAll(allCellInfo, celli)
321 {
322 if (cellType[celli] == MESH)
323 {
324 // Original cell was selected for meshing. Check if cell was
325 // reached from outsidePoints
326 if (allCellInfo[celli].type() != MESH)
327 {
328 cellType[celli] = NONMESH;
329 nChanged++;
330 }
331 }
332 }
333
334 return nChanged;
335}
336
337
338
339int main(int argc, char *argv[])
340{
342 (
343 "Select cells in relation to surface"
344 );
345
347
348 #include "setRootCase.H"
349 #include "createTime.H"
350 #include "createPolyMesh.H"
351
352 // Mesh edge statistics calculator
353 edgeStats edgeCalc(mesh);
354
355
356 IOdictionary refineDict
357 (
359 (
360 "selectCellsDict",
361 runTime.system(),
362 mesh,
365 )
366 );
367
368 fileName surfName(refineDict.get<fileName>("surface"));
369 pointField outsidePts(refineDict.lookup("outsidePoints"));
370 const bool useSurface(refineDict.get<bool>("useSurface"));
371 const bool selectCut(refineDict.get<bool>("selectCut"));
372 const bool selectInside(refineDict.get<bool>("selectInside"));
373 const bool selectOutside(refineDict.get<bool>("selectOutside"));
374 const scalar nearDist(refineDict.get<scalar>("nearDistance"));
375
376
377 if (useSurface)
378 {
379 Info<< "Cells to be used for meshing (0=false, 1=true):" << nl
380 << " cells cut by surface : " << selectCut << nl
381 << " cells inside of surface : " << selectInside << nl
382 << " cells outside of surface : " << selectOutside << nl
383 << " cells with points further than : " << nearDist << nl
384 << endl;
385 }
386 else
387 {
388 Info<< "Cells to be used for meshing (0=false, 1=true):" << nl
389 << " cells reachable from outsidePoints:" << selectOutside << nl
390 << endl;
391 }
392
393 // Print edge stats on original mesh.
394 (void)edgeCalc.minLen(Info);
395
396 // Search engine on mesh. Face decomposition since faces might be warped.
397 meshSearch queryMesh(mesh);
398
399 // Check all 'outside' points
400 for (const point& outsidePoint : outsidePts)
401 {
402 const label celli = queryMesh.findCell(outsidePoint, -1, false);
403
404 if (returnReduceAnd(celli < 0))
405 {
407 << "outsidePoint " << outsidePoint
408 << " is not inside any cell"
409 << exit(FatalError);
410 }
411 }
412
413 // Cell status (compared to surface if provided): inside/outside/cut.
414 // Start off from everything selected and cut later.
416 (
417 mesh,
419 (
420 mesh.nCells(),
422 )
423 );
424
425
426 if (useSurface)
427 {
428 triSurface surf(surfName);
429
430 // Dump some stats
431 surf.writeStats(Info);
432
433 triSurfaceSearch querySurf(surf);
434
435 // Set cellType[celli] according to relation to surface
436 cutBySurface
437 (
438 mesh,
439 queryMesh,
440 querySurf,
441
442 outsidePts,
443 selectCut,
444 selectInside,
445 selectOutside,
446 nearDist,
447
448 cellType
449 );
450 }
451
452
453 // Now 'trim' all the corners from the mesh so meshing/surface extraction
454 // becomes easier.
455
456 label nHanging, nRegionEdges, nRegionPoints, nOutside;
457
458 do
459 {
460 Info<< "Removing cells which after subsetting would have all points"
461 << " on outside ..." << nl << endl;
462
463 nHanging = cellType.fillHangingCells
464 (
465 MESH, // meshType
466 NONMESH, // fill type
467 mesh.nCells()
468 );
469
470
471 Info<< "Removing edges connecting cells unconnected by faces ..."
472 << nl << endl;
473
474 nRegionEdges = cellType.fillRegionEdges
475 (
476 MESH, // meshType
477 NONMESH, // fill type
478 mesh.nCells()
479 );
480
481
482 Info<< "Removing points connecting cells unconnected by faces ..."
483 << nl << endl;
484
485 nRegionPoints = cellType.fillRegionPoints
486 (
487 MESH, // meshType
488 NONMESH, // fill type
489 mesh.nCells()
490 );
491
492 nOutside = 0;
493 if (selectOutside)
494 {
495 // Since we're selecting the cells reachable from outsidePoints
496 // and the set might have changed, redo the outsideCells
497 // calculation
498 nOutside = selectOutsideCells
499 (
500 mesh,
501 queryMesh,
502 outsidePts,
503 cellType
504 );
505 }
506 } while
507 (
508 nHanging != 0
509 || nRegionEdges != 0
510 || nRegionPoints != 0
511 || nOutside != 0
512 );
513
514 cellSet selectedCells(mesh, "selected", mesh.nCells()/10);
515 getType(cellType, MESH, selectedCells);
516
517 writeSet(selectedCells, "cells selected for meshing");
518
519
520 Info<< "End\n" << endl;
521
522 return 0;
523}
524
525
526// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
FaceCellWave plus data.
Definition MeshWave.H:58
bool hit() const noexcept
Is there a hit?
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static void noParallel()
Remove the parallel options.
Definition argList.C:599
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
'Cuts' a mesh with a surface.
Holds information regarding type of cell. Used in inside/outside determination in cellClassification.
Definition cellInfo.H:61
A collection of cell labels.
Definition cellSet.H:50
Helper class to calculate minimum edge length on mesh.
Definition edgeStats.H:52
A class for handling file names.
Definition fileName.H:75
Non-pointer based hierarchical recursive searching.
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition meshSearch.H:57
label findCell(const point &location, const label seedCelli=-1, const bool useTreeSearch=true) const
Find cell containing location.
Definition meshSearch.C:759
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
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:74
dynamicFvMesh & mesh
engineTime & runTime
Required Variables.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const cellShapeList & cells
List< wallPoints > allCellInfo(mesh_.nCells())
void set(List< bool > &bools, const labelUList &locations)
Set the listed locations (assign 'true').
Definition BitOps.C:35
cellType
Equivalent to enumeration in "vtkCellType.h" (should be uint8_t).
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
bool returnReduceAnd(const bool value, const int communicator=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
vector point
Point is a vector.
Definition point.H:37
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
const pointField & pts
Tree tree(triangles.begin(), triangles.end())
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299