faceAreaIntersect.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-2017 OpenFOAM Foundation
9  Copyright (C) 2018 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 "faceAreaIntersect.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 const Foam::Enum
34 <
36 >
38 ({
39  { triangulationMode::tmFan, "fan" },
40  { triangulationMode::tmMesh, "mesh" },
41 });
42 
43 Foam::scalar Foam::faceAreaIntersect::tol = 1e-6;
44 
45 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
46 
47 void Foam::faceAreaIntersect::triSliceWithPlane
48 (
49  const triPoints& tri,
50  const plane& pln,
51  FixedList<triPoints, 10>& tris,
52  label& nTris,
53  const scalar len
54 ) const
55 {
56  // distance to cutting plane
57  FixedList<scalar, 3> d;
58 
59  // determine how many of the points are above the cutting plane
60  label nCoPlanar = 0;
61  label nPos = 0;
62  label posI = -1;
63  label negI = -1;
64  label copI = -1;
65  forAll(tri, i)
66  {
67  d[i] = pln.signedDistance(tri[i]);
68 
69  if (mag(d[i]) < tol*len)
70  {
71  nCoPlanar++;
72  copI = i;
73  d[i] = 0.0;
74  }
75  else
76  {
77  if (d[i] > 0)
78  {
79  nPos++;
80  posI = i;
81  }
82  else
83  {
84  negI = i;
85  }
86  }
87  }
88 
89 
90  // Determine triangle area contribution
91 
92  if
93  (
94  (nPos == 3)
95  || ((nPos == 2) && (nCoPlanar == 1))
96  || ((nPos == 1) && (nCoPlanar == 2))
97  )
98  {
99  /*
100  /\ _____
101  / \ \ / /\
102  /____\ \ / / \
103  __________ ____v____ __/____\__
104 
105  all points above cutting plane
106  - add complete triangle to list
107  */
108 
109  tris[nTris++] = tri;
110  }
111  else if ((nPos == 2) && (nCoPlanar == 0))
112  {
113  /*
114  i1________i2
115  \ /
116  --\----/--
117  \ /
118  \/
119  i0
120 
121  2 points above plane, 1 below
122  - resulting quad above plane split into 2 triangles
123  - forget triangle below plane
124  */
125 
126  // point under the plane
127  label i0 = negI;
128 
129  // indices of remaining points
130  label i1 = d.fcIndex(i0);
131  label i2 = d.fcIndex(i1);
132 
133  // determine the two intersection points
134  point p01 = planeIntersection(d, tri, i0, i1);
135  point p02 = planeIntersection(d, tri, i0, i2);
136 
137  // forget triangle below plane
138  // - decompose quad above plane into 2 triangles and add to list
139  setTriPoints(tri[i1], tri[i2], p02, nTris, tris);
140  setTriPoints(tri[i1], p02, p01, nTris, tris);
141  }
142  else if (nPos == 1)
143  {
144  // point above the plane
145  label i0 = posI;
146 
147  if (nCoPlanar == 0)
148  {
149  /*
150  i0
151  /\
152  / \
153  --/----\--
154  /______\
155  i2 i1
156 
157  1 point above plane, 2 below
158  - keep triangle above intersection plane
159  - forget quad below plane
160  */
161 
162  // indices of remaining points
163  label i1 = d.fcIndex(i0);
164  label i2 = d.fcIndex(i1);
165 
166  // determine the two intersection points
167  point p01 = planeIntersection(d, tri, i1, i0);
168  point p02 = planeIntersection(d, tri, i2, i0);
169 
170  // add triangle above plane to list
171  setTriPoints(tri[i0], p01, p02, nTris, tris);
172  }
173  else
174  {
175  /*
176  i0
177  |\
178  | \
179  __|__\_i2_
180  | /
181  | /
182  |/
183  i1
184 
185  1 point above plane, 1 on plane, 1 below
186  - keep triangle above intersection plane
187  */
188 
189  // point indices
190  label i1 = negI;
191  label i2 = copI;
192 
193  // determine the intersection point
194  point p01 = planeIntersection(d, tri, i1, i0);
195 
196  // add triangle above plane to list - clockwise points
197  if (d.fcIndex(i0) == i1)
198  {
199  setTriPoints(tri[i0], p01, tri[i2], nTris, tris);
200  }
201  else
202  {
203  setTriPoints(tri[i0], tri[i2], p01, nTris, tris);
204  }
205  }
206  }
207  else
208  {
209  /*
210  _________ __________ ___________
211  /\ \ /
212  /\ / \ \ /
213  / \ /____\ \/
214  /____\
215 
216  all points below cutting plane - forget
217  */
218  }
219 }
220 
221 
222 Foam::scalar Foam::faceAreaIntersect::triangleIntersect
223 (
224  const triPoints& src,
225  const point& tgt0,
226  const point& tgt1,
227  const point& tgt2,
228  const vector& n
229 ) const
230 {
231  // Work storage
232  FixedList<triPoints, 10> workTris1;
233  label nWorkTris1 = 0;
234 
235  FixedList<triPoints, 10> workTris2;
236  label nWorkTris2 = 0;
237 
238  // Cut source triangle with all inward pointing faces of target triangle
239  // - triangles in workTris1 are inside target triangle
240 
241  const scalar srcArea(triArea(src));
242  if (srcArea < ROOTVSMALL)
243  {
244  return 0.0;
245  }
246 
247  // Typical length scale
248  const scalar t = sqrt(srcArea);
249 
250  // Edge 0
251  {
252  // Cut triangle src with plane and put resulting sub-triangles in
253  // workTris1 list
254 
255  scalar s = mag(tgt1 - tgt0);
256  if (s < ROOTVSMALL)
257  {
258  return 0.0;
259  }
260 
261  // Note: outer product with n pre-scaled with edge length. This is
262  // purely to avoid numerical errors because of drastically
263  // different vector lengths.
264  const vector n0((tgt0 - tgt1)^(-s*n));
265  const scalar magSqrN0(magSqr(n0));
266  if (magSqrN0 < ROOTVSMALL)
267  {
268  // Triangle either zero edge length (s == 0) or
269  // perpendicular to face normal n. In either case zero
270  // overlap area
271  return 0.0;
272  }
273  else
274  {
275  plane pl0(tgt0, n0/Foam::sqrt(magSqrN0), false);
276  triSliceWithPlane(src, pl0, workTris1, nWorkTris1, t);
277  }
278  }
279 
280  if (nWorkTris1 == 0)
281  {
282  return 0.0;
283  }
284 
285  // Edge 1
286  {
287  // Cut workTris1 with plane and put resulting sub-triangles in
288  // workTris2 list (re-use tris storage)
289 
290  scalar s = mag(tgt2 - tgt1);
291  if (s < ROOTVSMALL)
292  {
293  return 0.0;
294  }
295  const vector n1((tgt1 - tgt2)^(-s*n));
296  const scalar magSqrN1(magSqr(n1));
297 
298  if (magSqrN1 < ROOTVSMALL)
299  {
300  // Triangle either zero edge length (s == 0) or
301  // perpendicular to face normal n. In either case zero
302  // overlap area
303  return 0.0;
304  }
305  else
306  {
307  plane pl1(tgt1, n1/Foam::sqrt(magSqrN1), false);
308 
309  nWorkTris2 = 0;
310 
311  for (label i = 0; i < nWorkTris1; ++i)
312  {
313  triSliceWithPlane(workTris1[i], pl1, workTris2, nWorkTris2, t);
314  }
315 
316  if (nWorkTris2 == 0)
317  {
318  return 0.0;
319  }
320  }
321  }
322 
323  // Edge 2
324  {
325  // Cut workTris2 with plane and put resulting sub-triangles in
326  // workTris1 list (re-use workTris1 storage)
327 
328  scalar s = mag(tgt2 - tgt0);
329  if (s < ROOTVSMALL)
330  {
331  return 0.0;
332  }
333  const vector n2((tgt2 - tgt0)^(-s*n));
334  const scalar magSqrN2(magSqr(n2));
335 
336  if (magSqrN2 < ROOTVSMALL)
337  {
338  // Triangle either zero edge length (s == 0) or
339  // perpendicular to face normal n. In either case zero
340  // overlap area
341  return 0.0;
342  }
343  else
344  {
345  plane pl2(tgt2, n2/Foam::sqrt(magSqrN2), false);
346 
347  nWorkTris1 = 0;
348 
349  for (label i = 0; i < nWorkTris2; ++i)
350  {
351  triSliceWithPlane(workTris2[i], pl2, workTris1, nWorkTris1, t);
352  }
353 
354  if (nWorkTris1 == 0)
355  {
356  return 0.0;
357  }
358  else
359  {
360  // Calculate area of sub-triangles
361  scalar area = 0.0;
362  for (label i = 0; i < nWorkTris1; ++i)
363  {
364  area += triArea(workTris1[i]);
365 
366  if (cacheTriangulation_)
367  {
368  triangles_.append(workTris1[i]);
369  }
370  }
371 
372  return area;
373  }
374  }
375  }
376 }
377 
378 
379 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
380 
382 (
383  const pointField& pointsA,
384  const pointField& pointsB,
385  const DynamicList<face>& trisA,
386  const DynamicList<face>& trisB,
387  const bool reverseB,
388  const bool cacheTriangulation
389 )
390 :
391  pointsA_(pointsA),
392  pointsB_(pointsB),
393  trisA_(trisA),
394  trisB_(trisB),
395  reverseB_(reverseB),
396  cacheTriangulation_(cacheTriangulation),
397  triangles_(cacheTriangulation ? 10 : 0)
398 {}
399 
400 
401 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
402 
404 (
405  const face& f,
406  const pointField& points,
407  const triangulationMode& triMode,
408  faceList& faceTris
409 )
410 {
411  faceTris.resize(f.nTriangles());
412 
413  switch (triMode)
414  {
415  case triangulationMode::tmFan:
416  {
417  for (label i = 0; i < f.nTriangles(); ++i)
418  {
419  faceTris[i] = face(3);
420  faceTris[i][0] = f[0];
421  faceTris[i][1] = f[i + 1];
422  faceTris[i][2] = f[i + 2];
423  }
424 
425  break;
426  }
427  case triangulationMode::tmMesh:
428  {
429  const label nFaceTris = f.nTriangles();
430 
431  label nFaceTris1 = 0;
432  const label nFaceTris2 = f.triangles(points, nFaceTris1, faceTris);
433 
434  if (nFaceTris != nFaceTris1 || nFaceTris != nFaceTris2)
435  {
437  << "The numbers of reported triangles in the face do not "
438  << "match that generated by the triangulation"
439  << exit(FatalError);
440  }
441  }
442  }
443 }
444 
445 
446 
448 (
449  const face& faceA,
450  const face& faceB,
451  const vector& n
452 ) const
453 {
454  if (cacheTriangulation_)
455  {
456  triangles_.clear();
457  }
458 
459  // Intersect triangles
460  scalar totalArea = 0.0;
461  for (const face& triA : trisA_)
462  {
463  triPoints tpA = getTriPoints(pointsA_, triA, false);
464 
465  for (const face& triB : trisB_)
466  {
467  if (reverseB_)
468  {
469  totalArea +=
470  triangleIntersect
471  (
472  tpA,
473  pointsB_[triB[0]],
474  pointsB_[triB[1]],
475  pointsB_[triB[2]],
476  n
477  );
478  }
479  else
480  {
481  totalArea +=
482  triangleIntersect
483  (
484  tpA,
485  pointsB_[triB[2]],
486  pointsB_[triB[1]],
487  pointsB_[triB[0]],
488  n
489  );
490  }
491  }
492  }
493 
494  return totalArea;
495 }
496 
497 
499 (
500  const face& faceA,
501  const face& faceB,
502  const vector& n,
503  const scalar threshold
504 ) const
505 {
506  // Intersect triangles
507  scalar totalArea = 0.0;
508  for (const face& triA : trisA_)
509  {
510  const triPoints tpA = getTriPoints(pointsA_, triA, false);
511 
512  for (const face& triB : trisB_)
513  {
514  if (reverseB_)
515  {
516  totalArea +=
517  triangleIntersect
518  (
519  tpA,
520  pointsB_[triB[0]],
521  pointsB_[triB[1]],
522  pointsB_[triB[2]],
523  n
524  );
525  }
526  else
527  {
528  totalArea +=
529  triangleIntersect
530  (
531  tpA,
532  pointsB_[triB[2]],
533  pointsB_[triB[1]],
534  pointsB_[triB[0]],
535  n
536  );
537  }
538 
539  if (totalArea > threshold)
540  {
541  return true;
542  }
543  }
544  }
545 
546  return false;
547 }
548 
549 
550 // ************************************************************************* //
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:51
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
Foam::faceAreaIntersect::triangulationMode
triangulationMode
Definition: faceAreaIntersect.H:62
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:57
Foam::faceAreaIntersect::triangulate
static void triangulate(const face &f, const pointField &points, const triangulationMode &triMode, faceList &faceTris)
Triangulate a face using the given triangulation mode.
Definition: faceAreaIntersect.C:404
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
faceAreaIntersect.H
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::faceAreaIntersect::triangulationModeNames_
static const Enum< triangulationMode > triangulationModeNames_
Definition: faceAreaIntersect.H:68
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::triPoints
Triangle storage. Null constructable (unfortunately triangle<point, point> is not)
Definition: triPoints.H:52
Foam::faceAreaIntersect::faceAreaIntersect
faceAreaIntersect(const pointField &pointsA, const pointField &pointsB, const DynamicList< face > &trisA, const DynamicList< face > &trisB, const bool reverseB=false, const bool cacheTriangulation=false)
Construct from components.
Definition: faceAreaIntersect.C:382
Foam::label
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:62
Foam::Field< vector >
Foam::List::resize
void resize(const label newSize)
Adjust allocated size of list.
Definition: ListI.H:139
Foam::FatalError
error FatalError
Foam::faceAreaIntersect::overlaps
bool overlaps(const face &faceA, const face &faceB, const vector &n, const scalar threshold) const
Return area of intersection of faceA with faceB.
Definition: faceAreaIntersect.C:499
Foam::faceAreaIntersect::calc
scalar calc(const face &faceA, const face &faceB, const vector &n) const
Return area of intersection of faceA with faceB.
Definition: faceAreaIntersect.C:448
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< face >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::fieldTypes::area
const wordList area
Standard area field types (scalar, vector, tensor, etc)
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:74
Foam::point
vector point
Point is a vector.
Definition: point.H:43