Loading...
Searching...
No Matches
surfaceToCell.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) 2017-2024 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
27\*---------------------------------------------------------------------------*/
28
29#include "surfaceToCell.H"
30#include "polyMesh.H"
31#include "meshSearch.H"
32#include "triSurface.H"
33#include "triSurfaceSearch.H"
34#include "cellClassification.H"
35#include "cpuTime.H"
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40namespace Foam
41{
47}
48
49
50Foam::topoSetSource::addToUsageTable Foam::surfaceToCell::usage_
51(
52 surfaceToCell::typeName,
53 "\n Usage: surfaceToCell"
54 "<surface> <outsidePoints> <cut> <inside> <outside> <near> <curvature>\n\n"
55 " <surface> name of triSurface\n"
56 " <outsidePoints> list of points that define outside\n"
57 " <cut> boolean whether to include cells cut by surface\n"
58 " <inside> ,, ,, inside surface\n"
59 " <outside> ,, ,, outside surface\n"
60 " <near> scalar; include cells with centre <= near to surface\n"
61 " <curvature> scalar; include cells close to strong curvature"
62 " on surface\n"
63 " (curvature defined as difference in surface normal at nearest"
64 " point on surface for each vertex of cell)\n\n"
65);
66
67
68// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
69
70Foam::label Foam::surfaceToCell::getNearest
71(
72 const triSurfaceSearch& querySurf,
73 const label pointi,
74 const point& pt,
75 const vector& span,
76 Map<label>& cache
77)
78{
79 const auto iter = cache.cfind(pointi);
80
81 if (iter.good())
82 {
83 return iter.val(); // Return cached value
84 }
85
86 pointIndexHit inter = querySurf.nearest(pt, span);
87
88 // Triangle label (can be -1)
89 const label trii = inter.index();
90
91 // Store triangle on point
92 cache.insert(pointi, trii);
93
94 return trii;
95}
96
97
98bool Foam::surfaceToCell::differingPointNormals
99(
100 const triSurfaceSearch& querySurf,
101
102 const vector& span, // Current search span
103 const label celli,
104 const label cellTriI, // Nearest (to cell centre) surface triangle
105
106 Map<label>& pointToNearest // Cache for nearest triangle to point
107) const
108{
109 const triSurface& surf = querySurf.surface();
110 const vectorField& normals = surf.faceNormals();
111
112 const faceList& faces = mesh().faces();
113 const pointField& points = mesh().points();
114
115 const labelList& cFaces = mesh().cells()[celli];
116
117 for (const label facei : cFaces)
118 {
119 const face& f = faces[facei];
120
121 for (const label pointi : f)
122 {
123 label pointTriI =
124 getNearest
125 (
126 querySurf,
127 pointi,
128 points[pointi],
129 span,
130 pointToNearest
131 );
132
133 if (pointTriI != -1 && pointTriI != cellTriI)
134 {
135 scalar cosAngle = normals[pointTriI] & normals[cellTriI];
136
137 if (cosAngle < 0.9)
138 {
139 return true;
140 }
141 }
142 }
143 }
144 return false;
145}
146
147
148void Foam::surfaceToCell::combine(topoSet& set, const bool add) const
149{
151
152 if (useSurfaceOrientation_ && (includeInside_ || includeOutside_))
153 {
154 const meshSearch queryMesh(mesh_);
155
156 //- Calculate for each searchPoint inside/outside status.
157 boolList isInside(querySurf().calcInside(mesh_.cellCentres()));
158
159 if (verbose_)
160 {
161 Info<< " Marked inside/outside using surface orientation in = "
162 << timer.cpuTimeIncrement() << " s" << nl << endl;
163 }
164
165 forAll(isInside, celli)
166 {
167 if (isInside[celli] ? includeInside_ : includeOutside_)
168 {
169 addOrDelete(set, celli, add);
170 }
171 }
172 }
173 else if (includeCut_ || includeInside_ || includeOutside_)
174 {
175 //
176 // Cut cells with surface and classify cells
177 //
178
179
180 // Construct search engine on mesh
181
182 const meshSearch queryMesh(mesh_);
183
184
185 // Check all 'outside' points
186 for (const point& outsidePoint : outsidePoints_)
187 {
188 // Find cell point is in. Linear search.
189 label celli = queryMesh.findCell(outsidePoint, -1, false);
190
191 if (returnReduceAnd(celli < 0))
192 {
194 << "outsidePoint " << outsidePoint
195 << " is not inside any cell" << nl
196 << exit(FatalError);
197 }
198 }
199
200 // Cut faces with surface and classify cells
201
203 (
204 mesh_,
205 queryMesh,
206 querySurf(),
207 outsidePoints_
208 );
209
210
211 if (verbose_)
212 {
213 Info<< " Marked inside/outside using surface intersection in = "
214 << timer.cpuTimeIncrement() << " s" << nl << endl;
215 }
216
217 //- Add/remove cells using set
218 forAll(cellType, celli)
219 {
220 label cType = cellType[celli];
221
222 if
223 (
224 (
225 includeCut_
226 && (cType == cellClassification::CUT)
227 )
228 || (
229 includeInside_
230 && (cType == cellClassification::INSIDE)
231 )
232 || (
233 includeOutside_
234 && (cType == cellClassification::OUTSIDE)
235 )
236 )
237 {
238 addOrDelete(set, celli, add);
239 }
240 }
241 }
242
243
244 if (nearDist_ > 0)
245 {
246 //
247 // Determine distance to surface
248 //
249
250 const pointField& ctrs = mesh_.cellCentres();
251
252 // Box dimensions to search in octree.
253 const vector span(nearDist_, nearDist_, nearDist_);
254
255
256 if (curvature_ < -1)
257 {
258 if (verbose_)
259 {
260 Info<< " Selecting cells with cellCentre closer than "
261 << nearDist_ << " to surface" << endl;
262 }
263
264 // No need to test curvature. Insert near cells into set.
265
266 forAll(ctrs, celli)
267 {
268 const point& c = ctrs[celli];
269
270 pointIndexHit inter = querySurf().nearest(c, span);
271
272 if (inter.hit() && inter.point().dist(c) < nearDist_)
273 {
274 addOrDelete(set, celli, add);
275 }
276 }
277
278 if (verbose_)
279 {
280 Info<< " Determined nearest surface point in = "
281 << timer.cpuTimeIncrement() << " s" << nl << endl;
282 }
283 }
284 else
285 {
286 // Test near cells for curvature
287
288 if (verbose_)
289 {
290 Info<< " Selecting cells with cellCentre closer than "
291 << nearDist_ << " to surface and curvature factor"
292 << " less than " << curvature_ << endl;
293 }
294
295 // Cache for nearest surface triangle for a point
296 Map<label> pointToNearest(mesh_.nCells()/10);
297
298 forAll(ctrs, celli)
299 {
300 const point& c = ctrs[celli];
301
302 pointIndexHit inter = querySurf().nearest(c, span);
303
304 if (inter.hit() && inter.point().dist(c) < nearDist_)
305 {
306 if
307 (
308 differingPointNormals
309 (
310 querySurf(),
311 span,
312 celli,
313 inter.index(), // nearest surface triangle
314 pointToNearest
315 )
316 )
317 {
318 addOrDelete(set, celli, add);
319 }
320 }
321 }
322
323 if (verbose_)
324 {
325 Info<< " Determined nearest surface point in = "
326 << timer.cpuTimeIncrement() << " s" << nl << endl;
327 }
328 }
329 }
330}
331
332
333void Foam::surfaceToCell::checkSettings() const
334{
335 if
336 (
337 (nearDist_ < 0)
338 && (curvature_ < -1)
339 && (
340 (includeCut_ && includeInside_ && includeOutside_)
341 || (!includeCut_ && !includeInside_ && !includeOutside_)
342 )
343 )
344 {
346 << "Illegal include cell specification."
347 << " Result would be either all or no cells." << endl
348 << "Please set one of includeCut, includeInside, includeOutside"
349 << " to true, set nearDistance to a value > 0"
350 << " or set curvature to a value -1 .. 1."
351 << exit(FatalError);
352 }
353
354 if (useSurfaceOrientation_ && includeCut_)
355 {
357 << "Illegal include cell specification."
358 << " You cannot specify both 'useSurfaceOrientation'"
359 << " and 'includeCut'"
360 << " since 'includeCut' specifies a topological split"
361 << exit(FatalError);
362 }
363}
364
366// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
367
369(
370 const polyMesh& mesh,
371 const fileName& surfName,
372 const pointField& outsidePoints,
373 const bool includeCut,
374 const bool includeInside,
375 const bool includeOutside,
376 const bool useSurfaceOrientation,
377 const scalar nearDist,
378 const scalar curvature
379)
380:
382 surfName_(surfName),
383 outsidePoints_(outsidePoints),
384 includeCut_(includeCut),
385 includeInside_(includeInside),
386 includeOutside_(includeOutside),
387 useSurfaceOrientation_(useSurfaceOrientation),
388 nearDist_(nearDist),
389 curvature_(curvature),
390 surfPtr_(new triSurface(surfName_)),
391 querySurfPtr_(new triSurfaceSearch(*surfPtr_))
392{
393 checkSettings();
395
396
398(
399 const polyMesh& mesh,
400 const fileName& surfName,
401 const triSurface& surf,
402 const triSurfaceSearch& querySurf,
403 const pointField& outsidePoints,
404 const bool includeCut,
405 const bool includeInside,
406 const bool includeOutside,
407 const bool useSurfaceOrientation,
408 const scalar nearDist,
409 const scalar curvature
410)
411:
413 surfName_(surfName),
414 outsidePoints_(outsidePoints),
415 includeCut_(includeCut),
416 includeInside_(includeInside),
417 includeOutside_(includeOutside),
418 useSurfaceOrientation_(useSurfaceOrientation),
419 nearDist_(nearDist),
420 curvature_(curvature),
421 surfPtr_(surf),
422 querySurfPtr_(querySurf)
423{
424 checkSettings();
426
427
429(
430 const polyMesh& mesh,
431 const dictionary& dict
432)
433:
435 surfName_(dict.get<fileName>("file").expand()),
436 outsidePoints_(dict.get<pointField>("outsidePoints")),
437 includeCut_(dict.get<bool>("includeCut")),
438 includeInside_(dict.get<bool>("includeInside")),
439 includeOutside_(dict.get<bool>("includeOutside")),
440 useSurfaceOrientation_
441 (
442 dict.getOrDefault("useSurfaceOrientation", false)
443 ),
444 nearDist_(dict.get<scalar>("nearDistance")),
445 curvature_(dict.get<scalar>("curvature")),
446 surfPtr_
447 (
448 new triSurface
449 (
450 surfName_,
451 dict.getOrDefault<word>("fileType", word::null),
452 dict.getOrDefault<scalar>("scale", -1)
453 )
454 ),
455 querySurfPtr_(new triSurfaceSearch(*surfPtr_))
456{
457 checkSettings();
459
460
462(
463 const polyMesh& mesh,
464 Istream& is
465)
466:
468 surfName_(checkIs(is)),
469 outsidePoints_(checkIs(is)),
470 includeCut_(readBool(checkIs(is))),
471 includeInside_(readBool(checkIs(is))),
472 includeOutside_(readBool(checkIs(is))),
473 useSurfaceOrientation_(false),
474 nearDist_(readScalar(checkIs(is))),
475 curvature_(readScalar(checkIs(is))),
476 surfPtr_(new triSurface(surfName_)),
477 querySurfPtr_(new triSurfaceSearch(*surfPtr_))
478{
479 checkSettings();
480}
481
483// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
484
486{} // Define here (incomplete type in header)
487
489// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
490
492(
493 const topoSetSource::setAction action,
494 topoSet& set
495) const
496{
497 if (action == topoSetSource::ADD || action == topoSetSource::NEW)
498 {
499 if (verbose_)
500 {
501 Info<< " Adding cells in relation to surface " << surfName_
502 << " ..." << endl;
503 }
504
505 combine(set, true);
506 }
507 else if (action == topoSetSource::SUBTRACT)
508 {
509 if (verbose_)
510 {
511 Info<< " Removing cells in relation to surface " << surfName_
512 << " ..." << endl;
513 }
514
515 combine(set, false);
516 }
517}
518
519
520// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
label index() const noexcept
Return the hit index.
'Cuts' a mesh with a surface.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
A class for handling file names.
Definition fileName.H:75
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition meshSearch.H:57
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
const cellList & cells() const
A topoSetCellSource to select cells based on relation to a surface given by an external file.
virtual void applyToSet(const topoSetSource::setAction action, topoSet &) const
Apply specified action to the topoSet.
surfaceToCell(const polyMesh &mesh, const fileName &surfName, const pointField &outsidePoints, const bool includeCut, const bool includeInside, const bool includeOutside, const bool useSurfaceOrientation, const scalar nearDist, const scalar curvature)
Construct from components.
virtual ~surfaceToCell()
Destructor.
Implements a timeout mechanism via sigalarm.
Definition timer.H:83
The topoSetCellSource is a intermediate class for handling topoSet sources for selecting cells.
topoSetCellSource(const polyMesh &mesh)
Construct from mesh.
Class with constructor to add usage string to table.
Base class of a source for a topoSet.
setAction
Enumeration defining various actions.
@ SUBTRACT
Subtract elements from current set.
@ ADD
Add elements to current set.
@ NEW
Create a new set and ADD elements to it.
bool verbose_
Output verbosity (default: true).
const polyMesh & mesh() const noexcept
Reference to the mesh.
static Istream & checkIs(Istream &is)
Check state of stream.
General set of labels of mesh quantity (points, cells, faces).
Definition topoSet.H:63
Helper class to search on triSurface.
Triangulated surface description with patch information.
Definition triSurface.H:74
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
const dimensionedScalar c
Speed of light in a vacuum.
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
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void add(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
bool returnReduceAnd(const bool value, const int communicator=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
cpuTimePosix cpuTime
Selection of preferred clock mechanism for the elapsed cpu time.
Definition cpuTimeFwd.H:38
bool readBool(Istream &is)
Read bool from stream using Foam::Switch(Istream&).
Definition bool.C:72
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299