surfaceToCell.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-2021 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 \*---------------------------------------------------------------------------*/
25 
26 #include "surfaceToCell.H"
27 #include "polyMesh.H"
28 #include "meshSearch.H"
29 #include "triSurface.H"
30 #include "triSurfaceSearch.H"
31 #include "cellClassification.H"
32 #include "cpuTime.H"
33 #include "demandDrivenData.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(surfaceToCell, 0);
41  addToRunTimeSelectionTable(topoSetSource, surfaceToCell, word);
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 Foam::label Foam::surfaceToCell::getNearest
48 (
49  const triSurfaceSearch& querySurf,
50  const label pointi,
51  const point& pt,
52  const vector& span,
53  Map<label>& cache
54 )
55 {
56  Map<label>::const_iterator iter = cache.find(pointi);
57 
58  if (iter != cache.end())
59  {
60  // Found cached answer
61  return iter();
62  }
63  else
64  {
65  pointIndexHit inter = querySurf.nearest(pt, span);
66 
67  // Triangle label (can be -1)
68  label triI = inter.index();
69 
70  // Store triangle on point
71  cache.insert(pointi, triI);
72 
73  return triI;
74  }
75 }
76 
77 
78 bool Foam::surfaceToCell::differingPointNormals
79 (
80  const triSurfaceSearch& querySurf,
81 
82  const vector& span, // Current search span
83  const label celli,
84  const label cellTriI, // Nearest (to cell centre) surface triangle
85 
86  Map<label>& pointToNearest // Cache for nearest triangle to point
87 ) const
88 {
89  const triSurface& surf = querySurf.surface();
90  const vectorField& normals = surf.faceNormals();
91 
92  const faceList& faces = mesh().faces();
93  const pointField& points = mesh().points();
94 
95  const labelList& cFaces = mesh().cells()[celli];
96 
97  forAll(cFaces, cFacei)
98  {
99  const face& f = faces[cFaces[cFacei]];
100 
101  forAll(f, fp)
102  {
103  label pointi = f[fp];
104 
105  label pointTriI =
106  getNearest
107  (
108  querySurf,
109  pointi,
110  points[pointi],
111  span,
112  pointToNearest
113  );
114 
115  if (pointTriI != -1 && pointTriI != cellTriI)
116  {
117  scalar cosAngle = normals[pointTriI] & normals[cellTriI];
118 
119  if (cosAngle < 0.9)
120  {
121  return true;
122  }
123  }
124  }
125  }
126  return false;
127 }
128 
129 
130 void Foam::surfaceToCell::combine(topoSet& set, const bool add) const
131 {
132  cpuTime timer;
133 
134 
135  if (useSurfaceOrientation_ && (includeInside_ || includeOutside_))
136  {
137  const meshSearch queryMesh(mesh_);
138 
139  //- Calculate for each searchPoint inside/outside status.
140  boolList isInside(querySurf().calcInside(mesh_.cellCentres()));
141 
142  Info<< " Marked inside/outside using surface orientation in = "
143  << timer.cpuTimeIncrement() << " s" << endl << endl;
144 
145  forAll(isInside, celli)
146  {
147  if (isInside[celli] && includeInside_)
148  {
149  addOrDelete(set, celli, add);
150  }
151  else if (!isInside[celli] && includeOutside_)
152  {
153  addOrDelete(set, celli, add);
154  }
155  }
156  }
157  else if (includeCut_ || includeInside_ || includeOutside_)
158  {
159  //
160  // Cut cells with surface and classify cells
161  //
162 
163 
164  // Construct search engine on mesh
165 
166  const meshSearch queryMesh(mesh_);
167 
168 
169  // Check all 'outside' points
170  forAll(outsidePoints_, outsideI)
171  {
172  const point& outsidePoint = outsidePoints_[outsideI];
173 
174  // Find cell point is in. Linear search.
175  label celli = queryMesh.findCell(outsidePoint, -1, false);
176  if (returnReduce(celli, maxOp<label>()) == -1)
177  {
179  << "outsidePoint " << outsidePoint
180  << " is not inside any cell"
181  << exit(FatalError);
182  }
183  }
184 
185  // Cut faces with surface and classify cells
186 
187  cellClassification cellType
188  (
189  mesh_,
190  queryMesh,
191  querySurf(),
192  outsidePoints_
193  );
194 
195 
196  Info<< " Marked inside/outside using surface intersection in = "
197  << timer.cpuTimeIncrement() << " s" << endl << endl;
198 
199  //- Add/remove cells using set
200  forAll(cellType, celli)
201  {
202  label cType = cellType[celli];
203 
204  if
205  (
206  (
207  includeCut_
208  && (cType == cellClassification::CUT)
209  )
210  || (
211  includeInside_
212  && (cType == cellClassification::INSIDE)
213  )
214  || (
215  includeOutside_
216  && (cType == cellClassification::OUTSIDE)
217  )
218  )
219  {
220  addOrDelete(set, celli, add);
221  }
222  }
223  }
224 
225 
226  if (nearDist_ > 0)
227  {
228  //
229  // Determine distance to surface
230  //
231 
232  const pointField& ctrs = mesh_.cellCentres();
233 
234  // Box dimensions to search in octree.
235  const vector span(nearDist_, nearDist_, nearDist_);
236 
237 
238  if (curvature_ < -1)
239  {
240  Info<< " Selecting cells with cellCentre closer than "
241  << nearDist_ << " to surface" << endl;
242 
243  // No need to test curvature. Insert near cells into set.
244 
245  forAll(ctrs, celli)
246  {
247  const point& c = ctrs[celli];
248 
249  pointIndexHit inter = querySurf().nearest(c, span);
250 
251  if (inter.hit() && (mag(inter.hitPoint() - c) < nearDist_))
252  {
253  addOrDelete(set, celli, add);
254  }
255  }
256 
257  Info<< " Determined nearest surface point in = "
258  << timer.cpuTimeIncrement() << " s" << endl << endl;
259 
260  }
261  else
262  {
263  // Test near cells for curvature
264 
265  Info<< " Selecting cells with cellCentre closer than "
266  << nearDist_ << " to surface and curvature factor"
267  << " less than " << curvature_ << endl;
268 
269  // Cache for nearest surface triangle for a point
270  Map<label> pointToNearest(mesh_.nCells()/10);
271 
272  forAll(ctrs, celli)
273  {
274  const point& c = ctrs[celli];
275 
276  pointIndexHit inter = querySurf().nearest(c, span);
277 
278  if (inter.hit() && (mag(inter.hitPoint() - c) < nearDist_))
279  {
280  if
281  (
282  differingPointNormals
283  (
284  querySurf(),
285  span,
286  celli,
287  inter.index(), // nearest surface triangle
288  pointToNearest
289  )
290  )
291  {
292  addOrDelete(set, celli, add);
293  }
294  }
295  }
296 
297  Info<< " Determined nearest surface point in = "
298  << timer.cpuTimeIncrement() << " s" << endl << endl;
299  }
300  }
301 }
302 
303 
304 void Foam::surfaceToCell::checkSettings() const
305 {
306  if
307  (
308  (nearDist_ < 0)
309  && (curvature_ < -1)
310  && (
311  (includeCut_ && includeInside_ && includeOutside_)
312  || (!includeCut_ && !includeInside_ && !includeOutside_)
313  )
314  )
315  {
317  << "Illegal include cell specification."
318  << " Result would be either all or no cells." << endl
319  << "Please set one of includeCut, includeInside, includeOutside"
320  << " to true, set nearDistance to a value > 0"
321  << " or set curvature to a value -1 .. 1."
322  << exit(FatalError);
323  }
324 
325  if (useSurfaceOrientation_ && includeCut_)
326  {
328  << "Illegal include cell specification."
329  << " You cannot specify both 'useSurfaceOrientation'"
330  << " and 'includeCut'"
331  << " since 'includeCut' specifies a topological split"
332  << exit(FatalError);
333  }
334 }
335 
336 
337 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
338 
340 (
341  const polyMesh& mesh,
342  const fileName& surfName,
343  const pointField& outsidePoints,
344  const bool includeCut,
345  const bool includeInside,
346  const bool includeOutside,
347  const bool useSurfaceOrientation,
348  const scalar nearDist,
349  const scalar curvature
350 )
351 :
352  topoSetSource(mesh),
353  surfName_(surfName),
354  outsidePoints_(outsidePoints),
355  includeCut_(includeCut),
356  includeInside_(includeInside),
357  includeOutside_(includeOutside),
358  useSurfaceOrientation_(useSurfaceOrientation),
359  nearDist_(nearDist),
360  curvature_(curvature),
361  surfPtr_(new triSurface(surfName_)),
362  querySurfPtr_(new triSurfaceSearch(*surfPtr_)),
363  IOwnPtrs_(true)
364 {
365  checkSettings();
366 }
367 
368 
370 (
371  const polyMesh& mesh,
372  const fileName& surfName,
373  const triSurface& surf,
374  const triSurfaceSearch& querySurf,
375  const pointField& outsidePoints,
376  const bool includeCut,
377  const bool includeInside,
378  const bool includeOutside,
379  const bool useSurfaceOrientation,
380  const scalar nearDist,
381  const scalar curvature
382 )
383 :
384  topoSetSource(mesh),
385  surfName_(surfName),
386  outsidePoints_(outsidePoints),
387  includeCut_(includeCut),
388  includeInside_(includeInside),
389  includeOutside_(includeOutside),
390  useSurfaceOrientation_(useSurfaceOrientation),
391  nearDist_(nearDist),
392  curvature_(curvature),
393  surfPtr_(&surf),
394  querySurfPtr_(&querySurf),
395  IOwnPtrs_(false)
396 {
397  checkSettings();
398 }
399 
400 
402 (
403  const polyMesh& mesh,
404  const dictionary& dict
405 )
406 :
407  topoSetSource(mesh),
408  surfName_(fileName(dict.lookup("file")).expand()),
409  outsidePoints_(dict.lookup("outsidePoints")),
410  includeCut_(readBool(dict.lookup("includeCut"))),
411  includeInside_(readBool(dict.lookup("includeInside"))),
412  includeOutside_(readBool(dict.lookup("includeOutside"))),
413  useSurfaceOrientation_
414  (
415  dict.lookupOrDefault<bool>("useSurfaceOrientation", false)
416  ),
417  nearDist_(dict.lookup<scalar>("nearDistance")),
418  curvature_(dict.lookup<scalar>("curvature")),
419  surfPtr_(new triSurface(surfName_)),
420  querySurfPtr_(new triSurfaceSearch(*surfPtr_)),
421  IOwnPtrs_(true)
422 {
423  checkSettings();
424 }
425 
426 
427 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
428 
430 {
431  if (IOwnPtrs_)
432  {
433  deleteDemandDrivenData(surfPtr_);
434  deleteDemandDrivenData(querySurfPtr_);
435  }
436 }
437 
438 
439 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
440 
442 (
443  const topoSetSource::setAction action,
444  topoSet& set
445 ) const
446 {
447  if ((action == topoSetSource::NEW) || (action == topoSetSource::ADD))
448  {
449  Info<< " Adding cells in relation to surface " << surfName_
450  << " ..." << endl;
451 
452  combine(set, true);
453  }
454  else if (action == topoSetSource::DELETE)
455  {
456  Info<< " Removing cells in relation to surface " << surfName_
457  << " ..." << endl;
458 
459  combine(set, false);
460  }
461 }
462 
463 
464 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
string expand(const string &, const HashTable< string, word, string::hash > &mapping, const char sigil='$')
Expand occurrences of variables according to the mapping.
Definition: stringOps.C:69
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A class for handling file names.
Definition: fileName.H:79
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
virtual void applyToSet(const topoSetSource::setAction action, topoSet &) const
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
HashTable< label, label, Hash< label > >::const_iterator const_iterator
Definition: Map.H:59
const cellList & cells() const
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
bool readBool(Istream &)
Definition: boolIO.C:60
fvMesh & mesh
const dimensionedScalar c
Speed of light in a vacuum.
virtual ~surfaceToCell()
Destructor.
Macros for easy insertion into run-time selection tables.
Base class of a source for a topoSet.
Definition: topoSetSource.H:63
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
Helper class to search on triSurface.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
setAction
Enumeration defining the valid actions.
Definition: topoSetSource.H:82
List< label > labelList
A List of labels.
Definition: labelList.H:56
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1224
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
General set of labels of mesh quantity (points, cells, faces).
Definition: topoSet.H:61
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.
Template functions to aid in the implementation of demand driven data.
vector point
Point is a vector.
Definition: point.H:41
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
Triangulated surface description with patch information.
Definition: triSurface.H:66
void deleteDemandDrivenData(DataPtr &dataPtr)
Namespace for OpenFOAM.
treeBoundBox combine(const treeBoundBox &a, const treeBoundBox &b)
Definition: patchToPatch.C:77
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864