cuttingPlane.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 | Copyright (C) 2011-2016 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 "cuttingPlane.H"
27 #include "primitiveMesh.H"
28 #include "linePointRef.H"
29 #include "meshTools.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 // Set values for what is close to zero and what is considered to
34 // be positive (and not just rounding noise)
36 const Foam::scalar zeroish = Foam::SMALL;
37 const Foam::scalar positive = Foam::SMALL * 1E3;
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 void Foam::cuttingPlane::calcCutCells
43 (
44  const primitiveMesh& mesh,
45  const scalarField& dotProducts,
46  const labelUList& cellIdLabels
47 )
48 {
49  const labelListList& cellEdges = mesh.cellEdges();
50  const edgeList& edges = mesh.edges();
51 
52  label listSize = cellEdges.size();
53  if (notNull(cellIdLabels))
54  {
55  listSize = cellIdLabels.size();
56  }
57 
58  cutCells_.setSize(listSize);
59  label cutcelli(0);
60 
61  // Find the cut cells by detecting any cell that uses points with
62  // opposing dotProducts.
63  for (label listI = 0; listI < listSize; ++listI)
64  {
65  label celli = listI;
66 
67  if (notNull(cellIdLabels))
68  {
69  celli = cellIdLabels[listI];
70  }
71 
72  const labelList& cEdges = cellEdges[celli];
73 
74  label nCutEdges = 0;
75 
76  forAll(cEdges, i)
77  {
78  const edge& e = edges[cEdges[i]];
79 
80  if
81  (
82  (dotProducts[e[0]] < zeroish && dotProducts[e[1]] > positive)
83  || (dotProducts[e[1]] < zeroish && dotProducts[e[0]] > positive)
84  )
85  {
86  nCutEdges++;
87 
88  if (nCutEdges > 2)
89  {
90  cutCells_[cutcelli++] = celli;
91 
92  break;
93  }
94  }
95  }
96  }
97 
98  // Set correct list size
99  cutCells_.setSize(cutcelli);
100 }
101 
102 
103 void Foam::cuttingPlane::intersectEdges
104 (
105  const primitiveMesh& mesh,
106  const scalarField& dotProducts,
107  List<label>& edgePoint
108 )
109 {
110  // Use the dotProducts to find out the cut edges.
111  const edgeList& edges = mesh.edges();
112  const pointField& points = mesh.points();
113 
114  // Per edge -1 or the label of the intersection point
115  edgePoint.setSize(edges.size());
116 
117  DynamicList<point> dynCuttingPoints(4*cutCells_.size());
118 
119  forAll(edges, edgeI)
120  {
121  const edge& e = edges[edgeI];
122 
123  if
124  (
125  (dotProducts[e[0]] < zeroish && dotProducts[e[1]] > positive)
126  || (dotProducts[e[1]] < zeroish && dotProducts[e[0]] > positive)
127  )
128  {
129  // Edge is cut
130  edgePoint[edgeI] = dynCuttingPoints.size();
131 
132  const point& p0 = points[e[0]];
133  const point& p1 = points[e[1]];
134 
135  scalar alpha = lineIntersect(linePointRef(p0, p1));
136 
137  if (alpha < zeroish)
138  {
139  dynCuttingPoints.append(p0);
140  }
141  else if (alpha >= 1.0)
142  {
143  dynCuttingPoints.append(p1);
144  }
145  else
146  {
147  dynCuttingPoints.append((1-alpha)*p0 + alpha*p1);
148  }
149  }
150  else
151  {
152  edgePoint[edgeI] = -1;
153  }
154  }
155 
156  this->storedPoints().transfer(dynCuttingPoints);
157 }
158 
159 
160 bool Foam::cuttingPlane::walkCell
161 (
162  const primitiveMesh& mesh,
163  const labelUList& edgePoint,
164  const label celli,
165  const label startEdgeI,
166  DynamicList<label>& faceVerts
167 )
168 {
169  label facei = -1;
170  label edgeI = startEdgeI;
171 
172  label nIter = 0;
173 
174  faceVerts.clear();
175  do
176  {
177  faceVerts.append(edgePoint[edgeI]);
178 
179  // Cross edge to other face
180  facei = meshTools::otherFace(mesh, celli, facei, edgeI);
181 
182  // Find next cut edge on face.
183  const labelList& fEdges = mesh.faceEdges()[facei];
184 
185  label nextEdgeI = -1;
186 
187  //Note: here is where we should check for whether there are more
188  // than 2 intersections with the face (warped/non-convex face).
189  // If so should e.g. decompose the cells on both faces and redo
190  // the calculation.
191 
192  forAll(fEdges, i)
193  {
194  label edge2I = fEdges[i];
195 
196  if (edge2I != edgeI && edgePoint[edge2I] != -1)
197  {
198  nextEdgeI = edge2I;
199  break;
200  }
201  }
202 
203  if (nextEdgeI == -1)
204  {
205  // Did not find another cut edge on facei. Do what?
207  << "Did not find closed walk along surface of cell " << celli
208  << " starting from edge " << startEdgeI
209  << " in " << nIter << " iterations." << nl
210  << "Collected cutPoints so far:" << faceVerts
211  << endl;
212 
213  return false;
214  }
215 
216  edgeI = nextEdgeI;
217 
218  nIter++;
219 
220  if (nIter > 1000)
221  {
223  << "Did not find closed walk along surface of cell " << celli
224  << " starting from edge " << startEdgeI
225  << " in " << nIter << " iterations." << nl
226  << "Collected cutPoints so far:" << faceVerts
227  << endl;
228  return false;
229  }
230 
231  } while (edgeI != startEdgeI);
232 
233 
234  if (faceVerts.size() >= 3)
235  {
236  return true;
237  }
238  else
239  {
241  << "Did not find closed walk along surface of cell " << celli
242  << " starting from edge " << startEdgeI << nl
243  << "Collected cutPoints so far:" << faceVerts
244  << endl;
245 
246  return false;
247  }
248 }
249 
250 
251 void Foam::cuttingPlane::walkCellCuts
252 (
253  const primitiveMesh& mesh,
254  const bool triangulate,
255  const labelUList& edgePoint
256 )
257 {
258  const pointField& cutPoints = this->points();
259 
260  // use dynamic lists to handle triangulation and/or missed cuts
261  DynamicList<face> dynCutFaces(cutCells_.size());
262  DynamicList<label> dynCutCells(cutCells_.size());
263 
264  // scratch space for calculating the face vertices
265  DynamicList<label> faceVerts(10);
266 
267  forAll(cutCells_, i)
268  {
269  label celli = cutCells_[i];
270 
271  // Find the starting edge to walk from.
272  const labelList& cEdges = mesh.cellEdges()[celli];
273 
274  label startEdgeI = -1;
275 
276  forAll(cEdges, cEdgeI)
277  {
278  label edgeI = cEdges[cEdgeI];
279 
280  if (edgePoint[edgeI] != -1)
281  {
282  startEdgeI = edgeI;
283  break;
284  }
285  }
286 
287  // Check for the unexpected ...
288  if (startEdgeI == -1)
289  {
291  << "Cannot find cut edge for cut cell " << celli
292  << abort(FatalError);
293  }
294 
295  // Walk from starting edge around the circumference of the cell.
296  bool okCut = walkCell
297  (
298  mesh,
299  edgePoint,
300  celli,
301  startEdgeI,
302  faceVerts
303  );
304 
305  if (okCut)
306  {
307  face f(faceVerts);
308 
309  // Orient face to point in the same direction as the plane normal
310  if ((f.normal(cutPoints) & normal()) < 0)
311  {
312  f.flip();
313  }
314 
315  // the cut faces are usually quite ugly, so optionally triangulate
316  if (triangulate)
317  {
318  label nTri = f.triangles(cutPoints, dynCutFaces);
319  while (nTri--)
320  {
321  dynCutCells.append(celli);
322  }
323  }
324  else
325  {
326  dynCutFaces.append(f);
327  dynCutCells.append(celli);
328  }
329  }
330  }
331 
332  this->storedFaces().transfer(dynCutFaces);
333  cutCells_.transfer(dynCutCells);
334 }
335 
336 
337 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
338 
339 // Construct without cutting
341 :
342  plane(pln)
343 {}
344 
345 
346 // Construct from plane and mesh reference, restricted to a list of cells
348 (
349  const plane& pln,
350  const primitiveMesh& mesh,
351  const bool triangulate,
352  const labelUList& cellIdLabels
353 )
354 :
355  plane(pln)
356 {
357  reCut(mesh, triangulate, cellIdLabels);
358 }
359 
360 
361 
362 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
363 
365 (
366  const primitiveMesh& mesh,
367  const bool triangulate,
368  const labelUList& cellIdLabels
369 )
370 {
372  cutCells_.clear();
373 
374  const scalarField dotProducts((mesh.points() - refPoint()) & normal());
375 
376  // Determine cells that are (probably) cut.
377  calcCutCells(mesh, dotProducts, cellIdLabels);
378 
379  // Determine cutPoints and return list of edge cuts.
380  // per edge -1 or the label of the intersection point
381  labelList edgePoint;
382  intersectEdges(mesh, dotProducts, edgePoint);
383 
384  // Do topological walk around cell to find closed loop.
385  walkCellCuts(mesh, triangulate, edgePoint);
386 }
387 
388 
390 (
391  const labelUList& faceMap
392 )
393 {
394  // recalculate the cells cut
395  if (notNull(faceMap) && faceMap.size())
396  {
397  MeshStorage::remapFaces(faceMap);
398 
399  List<label> newCutCells(faceMap.size());
400  forAll(faceMap, facei)
401  {
402  newCutCells[facei] = cutCells_[faceMap[facei]];
403  }
404  cutCells_.transfer(newCutCells);
405  }
406 }
407 
408 
409 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
410 
412 {
413  // Check for assignment to self
414  if (this == &rhs)
415  {
417  << "Attempted assignment to self"
418  << abort(FatalError);
419  }
420 
421  static_cast<MeshStorage&>(*this) = rhs;
422  static_cast<plane&>(*this) = rhs;
423  cutCells_ = rhs.cutCells();
424 }
425 
426 
427 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
void reCut(const primitiveMesh &, const bool triangulate, const labelUList &cellIdLabels=labelUList::null())
Recut mesh with existing planeDesc, restricted to a list of cells.
Definition: cuttingPlane.C:365
error FatalError
const labelList & cutCells() const
Return List of cells cut by the plane.
Definition: cuttingPlane.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const point & refPoint() const
Return or return plane base point.
Definition: plane.C:246
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:74
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
virtual void clear()
Clear all storage.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const Field< PointType > & points() const
Return reference to global points.
virtual const pointField & points() const =0
Return mesh points.
void operator=(const cuttingPlane &)
Definition: cuttingPlane.C:411
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
const vector & normal() const
Return plane normal.
Definition: plane.C:240
UList< label > labelUList
Definition: UList.H:64
plane(const vector &normalVector)
Construct from normal vector through the origin.
Definition: plane.C:112
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Constructs plane through mesh.
Definition: cuttingPlane.H:60
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
List< edge > edgeList
Definition: edgeList.H:38
List< face > & storedFaces()
Non-const access to the faces.
void clear()
Clear the list, i.e. set size to zero.
Definition: List.C:356
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
pointField & storedPoints()
Non-const access to global points.
virtual void remapFaces(const labelUList &faceMap)
Remap action on triangulation or cleanup.
Definition: cuttingPlane.C:390
List< label > labelList
A List of labels.
Definition: labelList.H:56
scalar lineIntersect(const line< Point, PointRef > &l) const
Return the cutting point between the plane and.
Definition: plane.H:180
errorManip< error > abort(error &err)
Definition: errorManip.H:131
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
volScalarField scalarField(fieldObject, mesh)
static const char nl
Definition: Ostream.H:262
labelList f(nPoints)
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:46
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
void setSize(const label)
Reset size of List.
Definition: List.C:295
virtual void remapFaces(const labelUList &faceMap)
Set new zones from faceMap.
vector point
Point is a vector.
Definition: point.H:41
#define WarningInFunction
Report a warning using Foam::Warning.
cuttingPlane(const plane &)
Construct plane description without cutting.
Definition: cuttingPlane.C:340
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:365
label otherFace(const primitiveMesh &, const label celli, const label facei, const label edgeI)
Return face on cell using edgeI but not facei. Throws error.
Definition: meshTools.C:536