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 void 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  scalar& area,
230  vector& centroid
231 ) const
232 {
233  // Work storage
234  FixedList<triPoints, 10> workTris1;
235  label nWorkTris1 = 0;
236 
237  FixedList<triPoints, 10> workTris2;
238  label nWorkTris2 = 0;
239 
240  // Cut source triangle with all inward pointing faces of target triangle
241  // - triangles in workTris1 are inside target triangle
242 
243  const scalar srcArea(triArea(src));
244  if (srcArea < ROOTVSMALL)
245  {
246  return;
247  }
248 
249  // Typical length scale
250  const scalar t = sqrt(srcArea);
251 
252  // Edge 0
253  {
254  // Cut triangle src with plane and put resulting sub-triangles in
255  // workTris1 list
256 
257  scalar s = mag(tgt1 - tgt0);
258  if (s < ROOTVSMALL)
259  {
260  return;
261  }
262 
263  // Note: outer product with n pre-scaled with edge length. This is
264  // purely to avoid numerical errors because of drastically
265  // different vector lengths.
266  const vector n0((tgt0 - tgt1)^(-s*n));
267  const scalar magSqrN0(magSqr(n0));
268  if (magSqrN0 < ROOTVSMALL)
269  {
270  // Triangle either zero edge length (s == 0) or
271  // perpendicular to face normal n. In either case zero
272  // overlap area
273  return;
274  }
275  else
276  {
277  plane pl0(tgt0, n0/Foam::sqrt(magSqrN0), false);
278  triSliceWithPlane(src, pl0, workTris1, nWorkTris1, t);
279  }
280  }
281 
282  if (nWorkTris1 == 0)
283  {
284  return;
285  }
286 
287  // Edge 1
288  {
289  // Cut workTris1 with plane and put resulting sub-triangles in
290  // workTris2 list (re-use tris storage)
291 
292  scalar s = mag(tgt2 - tgt1);
293  if (s < ROOTVSMALL)
294  {
295  return;
296  }
297  const vector n1((tgt1 - tgt2)^(-s*n));
298  const scalar magSqrN1(magSqr(n1));
299 
300  if (magSqrN1 < ROOTVSMALL)
301  {
302  // Triangle either zero edge length (s == 0) or
303  // perpendicular to face normal n. In either case zero
304  // overlap area
305  return;
306  }
307  else
308  {
309  plane pl1(tgt1, n1/Foam::sqrt(magSqrN1), false);
310 
311  nWorkTris2 = 0;
312 
313  for (label i = 0; i < nWorkTris1; ++i)
314  {
315  triSliceWithPlane(workTris1[i], pl1, workTris2, nWorkTris2, t);
316  }
317 
318  if (nWorkTris2 == 0)
319  {
320  return;
321  }
322  }
323  }
324 
325  // Edge 2
326  {
327  // Cut workTris2 with plane and put resulting sub-triangles in
328  // workTris1 list (re-use workTris1 storage)
329 
330  scalar s = mag(tgt2 - tgt0);
331  if (s < ROOTVSMALL)
332  {
333  return;
334  }
335  const vector n2((tgt2 - tgt0)^(-s*n));
336  const scalar magSqrN2(magSqr(n2));
337 
338  if (magSqrN2 < ROOTVSMALL)
339  {
340  // Triangle either zero edge length (s == 0) or
341  // perpendicular to face normal n. In either case zero
342  // overlap area
343  return;
344  }
345  else
346  {
347  plane pl2(tgt2, n2/Foam::sqrt(magSqrN2), false);
348 
349  nWorkTris1 = 0;
350 
351  for (label i = 0; i < nWorkTris2; ++i)
352  {
353  triSliceWithPlane(workTris2[i], pl2, workTris1, nWorkTris1, t);
354  }
355 
356  if (nWorkTris1 == 0)
357  {
358  return;
359  }
360  else
361  {
362  // Calculate area of sub-triangles
363  for (label i = 0; i < nWorkTris1; ++i)
364  {
365  // Area of intersection
366  const scalar currArea = triArea(workTris1[i]);
367  area += currArea;
368 
369  // Area-weighted centroid of intersection
370  centroid += currArea*triCentroid(workTris1[i]);
371 
372  if (cacheTriangulation_)
373  {
374  triangles_.append(workTris1[i]);
375  }
376  }
377  }
378  }
379  }
380 }
381 
382 
383 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
384 
386 (
387  const pointField& pointsA,
388  const pointField& pointsB,
389  const DynamicList<face>& trisA,
390  const DynamicList<face>& trisB,
391  const bool reverseB,
392  const bool cacheTriangulation
393 )
394 :
395  pointsA_(pointsA),
396  pointsB_(pointsB),
397  trisA_(trisA),
398  trisB_(trisB),
399  reverseB_(reverseB),
400  cacheTriangulation_(cacheTriangulation),
401  triangles_(cacheTriangulation ? 10 : 0)
402 {}
403 
404 
405 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
406 
408 (
409  const face& f,
410  const pointField& points,
411  const triangulationMode& triMode,
412  faceList& faceTris
413 )
414 {
415  faceTris.resize(f.nTriangles());
416 
417  switch (triMode)
418  {
419  case triangulationMode::tmFan:
420  {
421  for (label i = 0; i < f.nTriangles(); ++i)
422  {
423  faceTris[i] = face(3);
424  faceTris[i][0] = f[0];
425  faceTris[i][1] = f[i + 1];
426  faceTris[i][2] = f[i + 2];
427  }
428 
429  break;
430  }
431  case triangulationMode::tmMesh:
432  {
433  const label nFaceTris = f.nTriangles();
434 
435  label nFaceTris1 = 0;
436  const label nFaceTris2 = f.triangles(points, nFaceTris1, faceTris);
437 
438  if (nFaceTris != nFaceTris1 || nFaceTris != nFaceTris2)
439  {
441  << "The numbers of reported triangles in the face do not "
442  << "match that generated by the triangulation"
443  << exit(FatalError);
444  }
445  }
446  }
447 }
448 
449 
451 (
452  const face& faceA,
453  const face& faceB,
454  const vector& n,
455  scalar& area,
456  vector& centroid
457 ) const
458 {
459  if (cacheTriangulation_)
460  {
461  triangles_.clear();
462  }
463 
464  area = 0.0;
465  centroid = vector::zero;
466 
467  // Intersect triangles
468  for (const face& triA : trisA_)
469  {
470  triPoints tpA = getTriPoints(pointsA_, triA, false);
471 
472  for (const face& triB : trisB_)
473  {
474  if (reverseB_)
475  {
476  triangleIntersect
477  (
478  tpA,
479  pointsB_[triB[0]],
480  pointsB_[triB[1]],
481  pointsB_[triB[2]],
482  n,
483  area,
484  centroid
485  );
486  }
487  else
488  {
489  triangleIntersect
490  (
491  tpA,
492  pointsB_[triB[2]],
493  pointsB_[triB[1]],
494  pointsB_[triB[0]],
495  n,
496  area,
497  centroid
498  );
499  }
500  }
501  }
502 
503  // Area weighed centroid
504  if (area > 0)
505  {
506  centroid /= area;
507  }
508 }
509 
510 
512 (
513  const face& faceA,
514  const face& faceB,
515  const vector& n,
516  const scalar threshold
517 ) const
518 {
519  scalar area = 0.0;
520  vector centroid(Zero);
521 
522  // Intersect triangles
523  for (const face& triA : trisA_)
524  {
525  const triPoints tpA = getTriPoints(pointsA_, triA, false);
526 
527  for (const face& triB : trisB_)
528  {
529  if (reverseB_)
530  {
531  triangleIntersect
532  (
533  tpA,
534  pointsB_[triB[0]],
535  pointsB_[triB[1]],
536  pointsB_[triB[2]],
537  n,
538  area,
539  centroid
540  );
541  }
542  else
543  {
544  triangleIntersect
545  (
546  tpA,
547  pointsB_[triB[2]],
548  pointsB_[triB[1]],
549  pointsB_[triB[0]],
550  n,
551  area,
552  centroid
553  );
554  }
555 
556  if (area > threshold)
557  {
558  return true;
559  }
560  }
561  }
562 
563  return false;
564 }
565 
566 
567 // ************************************************************************* //
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
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::List::resize
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
Foam::faceAreaIntersect::triangulationMode
triangulationMode
Definition: faceAreaIntersect.H:63
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:55
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:408
Foam::faceAreaIntersect::calc
void calc(const face &faceA, const face &faceB, const vector &n, scalar &area, vector &centroid) const
Definition: faceAreaIntersect.C:451
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
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:69
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:386
Foam::Field< vector >
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:512
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:453
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::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::point
vector point
Point is a vector.
Definition: point.H:43