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