triIntersect.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) 2021-2025 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 "triIntersect.H"
27 #include "boundBox.H"
28 #include "cubicEqn.H"
29 #include "unitConversion.H"
30 #include "OFstream.H"
31 #include "tensor2D.H"
32 #include "vtkWritePolyData.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace triIntersect
39 {
40 
41 
42 //- The maximum dot product between a source point normal and a target plane
43 // considered to be a valid, forward projection
44 const scalar maxDot = - cos(degToRad(80));
45 
46 
47 //- Print 3x3 FixedListLists on one line
48 template <class Type>
50 {
51  os << token::BEGIN_LIST;
52  forAll(l, i)
53  {
54  if (i) os << token::SPACE;
55  os << l[i];
56  }
57  os << token::END_LIST;
58  return os;
59 }
60 
61 
62 //- Clip the given vector between values of 0 and 1, and also clip one minus
63 // its component sum. Clipping is applied to groups of components. It is done
64 // by moving the value linearly towards the value where all components in the
65 // group, and one minus their sum, share the same value.
67 {
68  vector y(x);
69 
70  for (label group = 0; group < groups[findMax(groups)] + 1; ++ group)
71  {
72  label n = 1;
73  forAll(x, i)
74  {
75  if (groups[i] == group)
76  {
77  n ++;
78  }
79  }
80 
81  if (n == 1)
82  {
83  continue;
84  }
85  else if (n == 2)
86  {
87  forAll(x, i)
88  {
89  if (groups[i] == group)
90  {
91  y[i] = min(max(x[i], 0), 1);
92  }
93  }
94  }
95  else
96  {
97  scalar xn = 1;
98  forAll(x, i)
99  {
100  if (groups[i] == group)
101  {
102  xn -= x[i];
103  }
104  }
105 
106  scalar phi = 0;
107  forAll(x, i)
108  {
109  if (groups[i] == group)
110  {
111  if (x[i] < 0)
112  {
113  phi = max(phi, n*x[i]/(n*x[i] - 1));
114  }
115  }
116  }
117  if (xn < 0)
118  {
119  phi = max(phi, n*xn/(n*xn - 1));
120  }
121 
122  forAll(x, i)
123  {
124  if (groups[i] == group)
125  {
126  y[i] = min(max((1 - phi)*x[i] + phi/n, 0), 1);
127  }
128  }
129  }
130  }
131 
132  return y;
133 }
134 
135 
136 //- Solve a projection equation given a value of the t variable
138 (
139  const vector& C,
140  const vector& Ct,
141  const vector& Cu,
142  const vector& Cv,
143  const vector& Ctu,
144  const vector& Ctv,
145  const FixedList<label, 3> groups,
146  const scalar t
147 )
148 {
149  // Solve a least squares problem for u and v
150  const vector CCtt = C + Ct*t;
151  const vector CuCtut = Cu + Ctu*t;
152  const vector CvCtvt = Cv + Ctv*t;
153 
154  const tensor2D A
155  (
156  CuCtut & CuCtut, CuCtut & CvCtvt,
157  CvCtvt & CuCtut, CvCtvt & CvCtvt
158  );
159  const vector2D B
160  (
161  - CuCtut & CCtt,
162  - CvCtvt & CCtt
163  );
164 
165  const scalar detA = det(A);
166  const vector2D detAuv = cof(A) & B;
167 
168  const vector tuv
169  (
170  t,
171  protectedDivide(detAuv.x(), detA),
172  protectedDivide(detAuv.y(), detA)
173  );
174 
175  // Apply group clipping
176  return clipped01(tuv, groups);
177 }
178 
179 
180 //- Solve a projection equation
182 (
183  const vector& C,
184  const vector& Ct,
185  const vector& Cu,
186  const vector& Cv,
187  const vector& Ctu,
188  const vector& Ctv,
189  const FixedList<label, 3> groups
190 )
191 {
192  // Solve the cubic projection equation for t
193  const Roots<3> tRoots =
194  cubicEqn
195  (
196  (Ct ^ Ctu) & Ctv,
197  ((C ^ Ctu) & Ctv) + ((Ct ^ Cu) & Ctv) + ((Ct ^ Ctu) & Cv),
198  ((C ^ Cu) & Ctv) + ((C ^ Ctu) & Cv) + ((Ct ^ Cu) & Cv),
199  (C ^ Cu) & Cv
200  ).roots();
201 
202  // Solve the remaining problem for u and v
203  label nTuvs = 0;
205  forAll(tRoots, tRooti)
206  {
207  if (tRoots.type(tRooti) != rootType::real) continue;
208 
209  if (mag(tRoots[tRooti]) > great) continue;
210 
211  const vector tuv =
213  (
214  C,
215  Ct,
216  Cu,
217  Cv,
218  Ctu,
219  Ctv,
220  {-1, -1, -1},
221  tRoots[tRooti]
222  );
223 
224  if (cmptMax(cmptMag(tuv)) > rootVGreat) continue;
225 
226  tuvs[nTuvs ++] = tuv;
227  }
228 
229  // Apply clipping
230  FixedList<scalar, 3> tuvClippage(NaN);
231  for (label i = 0; i < nTuvs; ++ i)
232  {
233  const vector tuvOld = tuvs[i];
234  tuvs[i] = clipped01(tuvs[i], groups);
235  tuvClippage[i] = cmptSum(cmptMag(tuvs[i] - tuvOld));
236  }
237 
238  // Sort so the least clipped roots come first
239  for (label i = 0; i < nTuvs - 1; ++ i)
240  {
241  for (label j = 0; j < nTuvs - 1; ++ j)
242  {
243  if (tuvClippage[j] > tuvClippage[j + 1])
244  {
245  Swap(tuvs[j], tuvs[j + 1]);
246  Swap(tuvClippage[j], tuvClippage[j + 1]);
247  }
248  }
249  }
250 
251  // Analyse each clipped solution value versus estimated error. If the
252  // value is small relative to the error, then return success.
253  for (label i = 0; i < nTuvs; ++ i)
254  {
255  const scalar t = tuvs[i].x(), u = tuvs[i].y(), v = tuvs[i].z();
256 
257  const scalar magSqrF =
258  magSqr(C + Ct*t + Cu*u + Cv*v + Ctu*t*u + Ctv*t*v);
259  const scalar magSqrErrF =
260  magSqr
261  (
262  1*cmptMag(C)
263  + 2*cmptMag(Ct)*mag(t)
264  + 2*cmptMag(Cu)*mag(u)
265  + 2*cmptMag(Cv)*mag(v)
266  + 3*cmptMag(Ctu)*mag(t)*mag(u)
267  + 3*cmptMag(Ctv)*mag(t)*mag(v)
268  )*small;
269 
270  if (magSqrF < magSqrErrF)
271  {
272  return Tuple2<bool, vector>(true, tuvs[i]);
273  }
274  }
275 
276  // No suitable roots were found. Return failure.
277  return Tuple2<bool, vector>(false, vector::uniform(NaN));
278 }
279 
280 
281 //- Calculate the non-dimensional offsets of the source points from the target
282 // edges. These values are considered indicative only. The calculation is not
283 // as reliable as that for the target point offsets from the source edges (see
284 // below). The target point offsets should take precedence where possible.
286 (
287  const FixedList<point, 3>& srcPs,
288  const FixedList<vector, 3>& srcNs,
289  const FixedList<bool, 3>& srcOwns,
290  const FixedList<point, 3>& tgtPs,
291  const FixedList<bool, 3>& tgtOwns
292 )
293 {
295 
296  const triPointRef tgtTri(tgtPs[0], tgtPs[1], tgtPs[2]);
297  const vector tgtN = tgtTri.normal();
298 
299  // Check whether the intersections occur in a forward direction
300  FixedList<bool, 3> srcFwd;
301  forAll(srcPs, srcPi)
302  {
303  srcFwd[srcPi] = (srcNs[srcPi] & tgtN) < maxDot;
304  }
305 
306  // For all forward projecting source points, determine the offset from the
307  // target edges
308  forAll(srcPs, srcPi)
309  {
310  if (srcFwd[srcPi])
311  {
312  forAll(tgtPs, tgtEi)
313  {
314  const label tgtPi0 = tgtEi, tgtPi1 = (tgtEi + 1) % 3;
315  result[srcPi][tgtEi] =
317  (
318  srcPs[srcPi],
319  srcNs[srcPi],
320  {tgtPs[tgtPi0], tgtPs[tgtPi1]},
321  tgtOwns[tgtEi]
322  );
323  }
324  }
325  }
326 
327  // For all backward projecting source points, initialise to outside
328  // everything
329  forAll(srcPs, srcPi)
330  {
331  if (!srcFwd[srcPi])
332  {
333  result[srcPi] = {-vGreat, -vGreat, -vGreat};
334  }
335  }
336 
337  // For source edges with one forward projecting and one backward
338  // projecting point, compute the point and normal where the projection
339  // direction changes and use this to determine the offset of the backward
340  // projecting point
341  forAll(srcPs, srcEi)
342  {
343  const label srcPi0 = srcEi, srcPi1 = (srcEi + 1) % 3;
344 
345  if (srcFwd[srcPi0] != srcFwd[srcPi1])
346  {
347  // Get the source edge parameter and normal at the asymptote where
348  // the normal switches sign relative to the target
349  const scalar srcT =
351  (
352  maxDot - (srcNs[srcPi0] & tgtN),
353  (srcNs[srcPi1] - srcNs[srcPi0]) & tgtN
354  );
355  const point srcP = (1 - srcT)*srcPs[srcPi0] + srcT*srcPs[srcPi1];
356  const vector srcN = (1 - srcT)*srcNs[srcPi0] + srcT*srcNs[srcPi1];
357 
358  forAll(tgtPs, tgtEi)
359  {
360  const label tgtPi0 = tgtEi, tgtPi1 = (tgtEi + 1) % 3;
361  result[srcFwd[srcPi0] ? srcPi1 : srcPi0][tgtEi] =
363  (
364  srcP,
365  srcN,
366  {tgtPs[tgtPi0], tgtPs[tgtPi1]},
367  tgtOwns[tgtEi]
368  );
369  }
370  }
371  }
372 
373  return result;
374 }
375 
376 
377 //- Calculate the non-dimensional offsets the target points from the source
378 // edges. These values are considered definitive, and should take precedence
379 // over the source point target edge offsets.
381 (
382  const FixedList<point, 3>& srcPs,
383  const FixedList<vector, 3>& srcNs,
384  const FixedList<bool, 3>& srcOwns,
385  const FixedList<point, 3>& tgtPs,
386  const FixedList<bool, 3>& tgtOwns
387 )
388 {
390 
391  // For all target points, determine the offset from each source edge
392  forAll(tgtPs, tgtPi)
393  {
394  forAll(srcPs, srcEi)
395  {
396  const label srcPi0 = srcEi, srcPi1 = (srcEi + 1) % 3;
397  result[tgtPi][srcEi] =
399  (
400  {srcPs[srcPi0], srcPs[srcPi1]},
401  {srcNs[srcPi0], srcNs[srcPi1]},
402  tgtPs[tgtPi],
403  srcOwns[srcEi]
404  );
405  }
406  }
407 
408  return result;
409 }
410 
411 
412 //- Construct point-inside/outside-edge topology from a set of point-edge
413 // offsets. Uses the sign of the offsets.
415 (
416  const FixedList<FixedList<scalar, 3>, 3>& thisOtherEdgeOffset
417 )
418 {
420 
421  // Determine the edge association from the sign of the offset
422  forAll(thisOtherEdgeOffset, thisPi)
423  {
424  forAll(thisOtherEdgeOffset[thisPi], otherEi)
425  {
426  result[thisPi][otherEi] =
427  thisOtherEdgeOffset[thisPi][otherEi] > 0 ? +1 : -1;
428  }
429  }
430 
431  return result;
432 }
433 
434 
435 //- Construct point-inside/outside-triangle topology from a set of
436 // point-inside/outside-edge topology
438 (
440 )
441 {
442  FixedList<label, 3> result;
443 
444  // Combine edge associations to get triangle associations
445  forAll(thisInOtherEdge, thisPi)
446  {
447  result[thisPi] = count(thisInOtherEdge[thisPi], 1) == 3 ? +1 : -1;
448  }
449 
450  return result;
451 }
452 
453 
454 //- Construct target-point-inside/outside-source-triangle topology from a set
455 // of target-point-inside/outside-source-edge topology, and some additional
456 // geometric information to handle cases where the source normal direction
457 // is reversed relative to the target triangle
459 (
460  const FixedList<point, 3>& srcPs,
461  const FixedList<vector, 3>& srcNs,
462  const FixedList<point, 3>& tgtPs,
463  const FixedList<FixedList<label, 3>, 3>& tgtInSrcEdge
464 )
465 {
466  const triPointRef tgtTri(tgtPs[0], tgtPs[1], tgtPs[2]);
467  const vector tgtN = tgtTri.normal();
468 
469  // Combine edge associations to get triangle associations
470  FixedList<label, 3> result = thisInOtherTri(tgtInSrcEdge);
471 
472  // Filter to only include forward intersections
473  forAll(tgtInSrcEdge, tgtPi)
474  {
475  if (result[tgtPi] == 1)
476  {
477  const barycentric2D srcTs =
478  srcTriTgtPointIntersection(srcPs, srcNs, tgtPs[tgtPi]);
479  const vector srcN = srcTriInterpolate(srcTs, srcNs);
480  const bool tgtFwd = (srcN & tgtN) < maxDot;
481  result[tgtPi] = tgtFwd ? +1 : -1;
482  }
483  }
484 
485  return result;
486 }
487 
488 
489 //- Override results of the srcInTgt/tgtInSrc calculations with explicit
490 // connections between points on either side
492 (
493  const FixedList<label, 3>& thisOtherPis,
496 )
497 {
498  forAll(thisOtherPis, thisPi)
499  {
500  const label otherPi = thisOtherPis[thisPi];
501 
502  if (otherPi != -1)
503  {
504  const label otherEi0 = (otherPi + 2) % 3, otherEi1 = otherPi;
505 
506  thisInOtherTri[thisPi] = 0;
507 
508  thisInOtherEdge[thisPi][otherEi0] = 0;
509  thisInOtherEdge[thisPi][otherEi1] = 0;
510  }
511  }
512 }
513 
514 
515 //- Calculate whether the points of the given source triangle project inside or
516 // outside the opposing target triangle and its edges
518 (
519  const FixedList<point, 3>& srcPs,
520  const FixedList<vector, 3>& srcNs,
521  const FixedList<bool, 3>& srcOwns,
522  const FixedList<label, 3>& srcTgtPis,
523  const FixedList<point, 3>& tgtPs,
524  const FixedList<bool, 3>& tgtOwns,
525  FixedList<FixedList<label, 3>, 3>& srcInTgtEdge,
526  FixedList<label, 3>& srcInTgtTri
527 )
528 {
530  triIntersect::srcTgtEdgeOffset(srcPs, srcNs, srcOwns, tgtPs, tgtOwns);
531 
532  srcInTgtEdge = thisInOtherEdge(srcTgtEdgeOffset);
533 
534  srcInTgtTri = thisInOtherTri(srcInTgtEdge);
535 
536  thisIsOther(srcTgtPis, srcInTgtEdge, srcInTgtTri);
537 }
538 
539 
540 //- Calculate whether the points of the given target triangle project inside or
541 // outside the opposing source triangle and its edges
543 (
544  const FixedList<point, 3>& srcPs,
545  const FixedList<vector, 3>& srcNs,
546  const FixedList<bool, 3>& srcOwns,
547  const FixedList<point, 3>& tgtPs,
548  const FixedList<bool, 3>& tgtOwns,
549  const FixedList<label, 3>& tgtSrcPis,
550  FixedList<FixedList<label, 3>, 3>& tgtInSrcEdge,
552 )
553 {
555  triIntersect::tgtSrcEdgeOffset(srcPs, srcNs, srcOwns, tgtPs, tgtOwns);
556 
557  tgtInSrcEdge = thisInOtherEdge(tgtSrcEdgeOffset);
558 
559  tgtInSrcTri = triIntersect::tgtInSrcTri(srcPs, srcNs, tgtPs, tgtInSrcEdge);
560 
561  thisIsOther(tgtSrcPis, tgtInSrcEdge, tgtInSrcTri);
562 }
563 
564 
565 //- Order intersection locations into a polygon
567 (
568  const UList<location>& locations,
569  bool isSrcEdge,
570  const label i0,
571  label& nVisited,
572  boolList& visited,
574 )
575 {
576  // Mark this location as visited
577  order[nVisited ++] = i0;
578  visited[i0] = true;
579 
580  // Get the index of the edge attached to this point
581  const location& l0 = locations[i0];
582  const label ei0 =
583  isSrcEdge
584  ? (l0.isSrcPoint() ? l0.srcPointi() : l0.srcEdgei())
585  : (l0.isTgtPoint() ? (l0.tgtPointi() + 2) % 3 : l0.tgtEdgei());
586 
587  // Terminate if connected back to the first location
588  {
589  const label i1 = order.first();
590 
591  const location& l1 = locations[i1];
592 
593  if
594  (
595  i0 != order.first()
596  && !(l1.isSrcNotTgtPoint() && !isSrcEdge)
597  && !(l1.isTgtNotSrcPoint() && isSrcEdge)
598  )
599  {
600  const location& l1 = locations[i1];
601  const label ei1 =
602  isSrcEdge
603  ? (l1.isSrcPoint() ? (l1.srcPointi() + 2) % 3 : l1.srcEdgei())
604  : (l1.isTgtPoint() ? l1.tgtPointi() : l1.tgtEdgei());
605 
606  if (ei0 == ei1)
607  {
608  return true;
609  }
610  }
611  }
612 
613  // Search for the next connected location and recurse if found
614  forAll(locations, i1)
615  {
616  if (!visited[i1])
617  {
618  const location& l1 = locations[i1];
619 
620  if
621  (
622  !(l1.isSrcNotTgtPoint() && !isSrcEdge)
623  && !(l1.isTgtNotSrcPoint() && isSrcEdge)
624  )
625  {
626  const label ei1 =
627  isSrcEdge
628  ? (l1.isSrcPoint() ? (l1.srcPointi() + 2) % 3 : l1.srcEdgei())
629  : (l1.isTgtPoint() ? l1.tgtPointi() : l1.tgtEdgei());
630 
631  if (ei0 == ei1)
632  {
633  auto branch = [&](const bool isSrcEdge)
634  {
635  return orderLocations
636  (
637  locations,
638  isSrcEdge,
639  i1,
640  nVisited,
641  visited,
642  order
643  );
644  };
645 
646  if
647  (
648  (
649  !l1.isSrcAndTgtPoint()
650  && branch(l1.isIntersection() != isSrcEdge)
651  )
652  || (
653  l1.isSrcAndTgtPoint()
654  && (branch(true) || branch(false))
655  )
656  )
657  {
658  return true;
659  }
660  }
661  }
662  }
663  }
664 
665  // This branch failed to find a connected location. Un-visit this location.
666  order[-- nVisited] = -1;
667  visited[i0] = false;
668 
669  return false;
670 }
671 
672 
673 //- Construct the intersection topology
675 (
676  const FixedList<label, 3>& tgtSrcPis,
677  const FixedList<FixedList<label, 3>, 3>& srcInTgtEdge,
678  const FixedList<FixedList<label, 3>, 3>& tgtInSrcEdge,
679  const FixedList<label, 3>& srcInTgtTri,
681  DynamicList<location>& pointLocations
682 )
683 {
684  // Step 1: Process trivial rejection cases
685 
686  // If the entire target triangle is outside or on the same source edge
687  // then there can be no intersection.
688  forAll(srcInTgtEdge, srcEi)
689  {
690  bool outside = true;
691  forAll(tgtInSrcEdge, tgtPi)
692  {
693  if (tgtInSrcEdge[tgtPi][srcEi] == 1)
694  {
695  outside = false;
696  break;
697  }
698  }
699  if (outside)
700  {
701  return true;
702  }
703  }
704 
705  // If all source points are outside all target edges this indicates
706  // that the triangles are oppositely oriented, in which case there can
707  // also be no intersection.
708  if (count(srcInTgtEdge, {-1, -1, -1}) == 3)
709  {
710  return true;
711  }
712 
713  // Step 2: Define point addition/checking functions
714 
715  // Add crossing point locations, inserting source points as necessary
716  auto addPointLocations = [&pointLocations]
717  (
718  const location l1,
719  const location l2 = location(),
720  const bool add = true
721  )
722  {
723  if (!pointLocations.empty())
724  {
725  const location l0 = pointLocations.last();
726 
727  if (l0.isIntersection() || l0.isSrcAndTgtPoint())
728  {
729  const label srcEi0 =
730  l0.isIntersection()
731  ? l0.srcEdgei()
732  : (l0.srcPointi() + 2) % 3;
733  const label tgtEi0 =
734  l0.isIntersection()
735  ? l0.tgtEdgei()
736  : l0.tgtPointi();
737 
738  const label srcEi1 =
739  l1.isIntersection()
740  ? l1.srcEdgei()
741  : l1.srcPointi();
742  const label tgtEi1 =
743  l1.isIntersection()
744  ? l1.tgtEdgei()
745  : (l1.tgtPointi() + 2) % 3;
746 
747  if
748  (
749  (l0.isIntersection() && l1.isIntersection())
750  || tgtEi0 != tgtEi1
751  )
752  {
753  for
754  (
755  label srcEj = srcEi0;
756  srcEj != srcEi1;
757  srcEj = (srcEj + 2) % 3
758  )
759  {
760  pointLocations.append(location::srcPoint(srcEj));
761  }
762  }
763  }
764  }
765 
766  if (add)
767  {
768  pointLocations.append(l1);
769 
770  if (!l2.isNull())
771  {
772  pointLocations.append(l2);
773  }
774  }
775  };
776 
777  // One target point is within the source triangle and one is not
778  auto inTriToOut = [&addPointLocations,&srcInTgtEdge]
779  (
780  const label tgtEi,
781  const label tgtOutSrcEi1,
782  const label tgtOutSrcPi1,
783  const bool reverse
784  )
785  {
786  const label srcEi =
787  tgtOutSrcEi1 != -1
788  ? tgtOutSrcEi1
789  : srcInTgtEdge[tgtOutSrcPi1][tgtEi] == 1
790  ? (tgtOutSrcPi1 + 2*reverse) % 3
791  : (tgtOutSrcPi1 + 2*!reverse) % 3;
792 
793  addPointLocations(location::intersection(srcEi, tgtEi));
794  };
795 
796  // One target point is a source point and the other is outside a source edge
797  auto isPointToOutEdge = [&addPointLocations,&srcInTgtEdge]
798  (
799  const label tgtEi,
800  const label tgtIsSrcPi0,
801  const label tgtOutSrcEi1,
802  const bool reverse
803  )
804  {
805  const label srcEi0Next = (tgtIsSrcPi0 + 2*reverse) % 3;
806  const label srcEi0Opp = (tgtIsSrcPi0 + 1) % 3;
807 
808  if
809  (
810  srcInTgtEdge[(tgtIsSrcPi0 + 1) % 3][tgtEi] == -1
811  && srcInTgtEdge[(tgtIsSrcPi0 + 2) % 3][tgtEi] == -1
812  )
813  {
814  return false;
815  }
816 
817  if (tgtOutSrcEi1 == srcEi0Next)
818  {
819  return false;
820  }
821 
822  if (tgtOutSrcEi1 == srcEi0Opp)
823  {
824  addPointLocations(location::intersection(srcEi0Opp, tgtEi));
825  }
826 
827  return true;
828  };
829 
830  // One target point is a source point and the other is outside a source
831  // corner
832  auto isPointToOutCorner = []
833  (
834  const label tgtEi,
835  const label tgtIsSrcPi0,
836  const label tgtOutSrcPi1,
837  const bool reverse
838  )
839  {
840  return tgtOutSrcPi1 != (tgtIsSrcPi0 + 1 + reverse) % 3;
841  };
842 
843  // Both target points are outside source edges
844  auto outEdgeToOutEdge = [&addPointLocations,&srcInTgtEdge]
845  (
846  const label tgtEi,
847  const label tgtOutSrcEi0,
848  const label tgtOutSrcEi1
849  )
850  {
851  const label srcPi = (5 - tgtOutSrcEi0 - tgtOutSrcEi1) % 3;
852 
853  if
854  (
855  (tgtOutSrcEi0 != (tgtOutSrcEi1 + 1) % 3)
856  && (srcInTgtEdge[srcPi][tgtEi] != 1)
857  )
858  {
859  return false;
860  }
861 
862  if
863  (
864  (tgtOutSrcEi0 == (tgtOutSrcEi1 + 1) % 3)
865  != (srcInTgtEdge[srcPi][tgtEi] == 1)
866  )
867  {
868  addPointLocations
869  (
870  location::intersection(tgtOutSrcEi0, tgtEi),
871  location::intersection(tgtOutSrcEi1, tgtEi)
872  );
873  }
874 
875  return true;
876  };
877 
878  // One target point is outside a source edge and the other is outside a
879  // source corner
880  auto outEdgeToOutCorner = [&addPointLocations,&srcInTgtEdge]
881  (
882  const label tgtEi,
883  const label tgtOutSrcEi0,
884  const label tgtOutSrcPi1,
885  const bool reverse
886  )
887  {
888  if (tgtOutSrcEi0 == tgtOutSrcPi1)
889  {
890  return !reverse;
891  }
892 
893  if ((tgtOutSrcEi0 + 1) % 3 == tgtOutSrcPi1)
894  {
895  return reverse;
896  }
897 
898  const label srcPi1Prev = (tgtOutSrcPi1 + 1 + !reverse) % 3;
899 
900  if (srcInTgtEdge[srcPi1Prev][tgtEi] == -1)
901  {
902  return false;
903  }
904 
905  const label srcPi1Next = (tgtOutSrcPi1 + 1 + reverse) % 3;
906 
907  if (srcInTgtEdge[srcPi1Next][tgtEi] == 1)
908  {
909  return true;
910  }
911 
912  location l1 = location::intersection(tgtOutSrcEi0, tgtEi);
913 
914  const label srcEi =
915  srcInTgtEdge[tgtOutSrcPi1][tgtEi] == 1
916  ? (tgtOutSrcPi1 + 2*reverse) % 3
917  : (tgtOutSrcPi1 + 2*!reverse) % 3;
918 
919  location l2 = location::intersection(srcEi, tgtEi);
920 
921  if (reverse)
922  {
923  Swap(l1, l2);
924  }
925 
926  addPointLocations(l1, l2);
927 
928  return true;
929  };
930 
931  // Both target points are outside source corners
932  auto outCornerToOutCorner = []
933  (
934  const label tgtEi,
935  const label tgtOutSrcPi0,
936  const label tgtOutSrcPi1
937  )
938  {
939  return tgtOutSrcPi0 != (tgtOutSrcPi1 + 2) % 3;
940  };
941 
942  // Step 3: Walk around the target edges to form the intersection polygon
943  for (label tgtEi = 0; tgtEi < 3; tgtEi ++)
944  {
945  const label tgtPi0 = tgtEi, tgtPi1 = (tgtEi + 1) % 3;
946 
947  const bool tgtInSrcTri0 = tgtInSrcTri[tgtPi0] == 1;
948  const bool tgtInSrcTri1 = tgtInSrcTri[tgtPi1] == 1;
949 
950  const label tgtIsSrcPi0 =
951  tgtInSrcTri[tgtPi0] == 0 ? tgtSrcPis[tgtPi0] : -1;
952  const label tgtIsSrcPi1 =
953  tgtInSrcTri[tgtPi1] == 0 ? tgtSrcPis[tgtPi1] : -1;
954 
955  const label tgtOutSrcEi0 =
956  count(tgtInSrcEdge[tgtPi0], -1) == 1
957  ? findIndex(tgtInSrcEdge[tgtPi0], -1)
958  : -1;
959  const label tgtOutSrcEi1 =
960  count(tgtInSrcEdge[tgtPi1], -1) == 1
961  ? findIndex(tgtInSrcEdge[tgtPi1], -1)
962  : -1;
963 
964  const label tgtOutSrcPi0 =
965  count(tgtInSrcEdge[tgtPi0], -1) == 2
966  ? (findIndex(tgtInSrcEdge[tgtPi0], 1) + 2) % 3
967  : -1;
968  const label tgtOutSrcPi1 =
969  count(tgtInSrcEdge[tgtPi1], -1) == 2
970  ? (findIndex(tgtInSrcEdge[tgtPi1], 1) + 2) % 3
971  : -1;
972 
973  // Add the first point if it within or part of the source triangle
974  if (tgtInSrcTri0)
975  {
976  pointLocations.append(location::tgtPoint(tgtPi0));
977  }
978  if (tgtIsSrcPi0 != -1)
979  {
980  addPointLocations(location::srcTgtPoint(tgtIsSrcPi0, tgtPi0));
981  }
982 
983  // Add crossings
984  if
985  (
986  (tgtInSrcTri0 && tgtInSrcTri1)
987  || (tgtOutSrcEi0 != -1 && tgtOutSrcEi0 == tgtOutSrcEi1)
988  || (tgtOutSrcPi0 != -1 && tgtOutSrcPi0 == tgtOutSrcPi1)
989  )
990  {
991  // Both target points are in the same source quadrant. There is
992  // nothing to check or to add.
993  }
994  else if
995  (
996  (tgtInSrcTri0 && tgtIsSrcPi1 != -1)
997  || (tgtIsSrcPi0 != -1 && tgtInSrcTri1)
998  )
999  {
1000  // One target point is within the source triangle and one is a
1001  // source point. There is nothing to check or to add.
1002  }
1003  else if (tgtInSrcTri0 && (tgtOutSrcEi1 != -1 || tgtOutSrcPi1 != -1))
1004  {
1005  // The first target point is within the source triangle and the
1006  // second is outside
1007  inTriToOut(tgtEi, tgtOutSrcEi1, tgtOutSrcPi1, 0);
1008  }
1009  else if ((tgtOutSrcEi0 != -1 || tgtOutSrcPi0 != -1) && tgtInSrcTri1)
1010  {
1011  // (reverse of previous clause)
1012  inTriToOut(tgtEi, tgtOutSrcEi0, tgtOutSrcPi0, 1);
1013  }
1014  else if (tgtIsSrcPi0 != -1 && tgtIsSrcPi1 != -1)
1015  {
1016  // Both target points are source points. Check the ordering is
1017  // compatible with an intersection.
1018  if (tgtIsSrcPi0 != (tgtIsSrcPi1 + 1) % 3)
1019  {
1020  pointLocations.clear();
1021  return false;
1022  }
1023  }
1024  else if (tgtIsSrcPi0 != -1 && tgtOutSrcEi1 != -1)
1025  {
1026  // The first target point is a source point and the second is
1027  // outside a source edge
1028  if (!isPointToOutEdge(tgtEi, tgtIsSrcPi0, tgtOutSrcEi1, 0))
1029  {
1030  pointLocations.clear();
1031  return false;
1032  }
1033  }
1034  else if (tgtOutSrcEi0 != -1 && tgtIsSrcPi1 != -1)
1035  {
1036  // (reverse of previous clause)
1037  if (!isPointToOutEdge(tgtEi, tgtIsSrcPi1, tgtOutSrcEi0, 1))
1038  {
1039  pointLocations.clear();
1040  return false;
1041  }
1042  }
1043  else if (tgtIsSrcPi0 != -1 && tgtOutSrcPi1 != -1)
1044  {
1045  // The first target point is a source point and the second is
1046  // outside a source corner
1047  if (!isPointToOutCorner(tgtEi, tgtIsSrcPi0, tgtOutSrcPi1, 0))
1048  {
1049  pointLocations.clear();
1050  return false;
1051  }
1052  }
1053  else if (tgtOutSrcPi0 != -1 && tgtIsSrcPi1 != -1)
1054  {
1055  // (reverse of previous clause)
1056  if (!isPointToOutCorner(tgtEi, tgtIsSrcPi1, tgtOutSrcPi0, 1))
1057  {
1058  pointLocations.clear();
1059  return false;
1060  }
1061  }
1062  else if (tgtOutSrcEi0 != -1 && tgtOutSrcEi1 != -1)
1063  {
1064  // Both target points are outside source edges
1065  if (!outEdgeToOutEdge(tgtEi, tgtOutSrcEi0, tgtOutSrcEi1))
1066  {
1067  pointLocations.clear();
1068  return false;
1069  }
1070  }
1071  else if (tgtOutSrcEi0 != -1 && tgtOutSrcPi1 != -1)
1072  {
1073  // The first target point is outside a source edge and the
1074  // second is outside a source corner
1075  if (!outEdgeToOutCorner(tgtEi, tgtOutSrcEi0, tgtOutSrcPi1, 0))
1076  {
1077  pointLocations.clear();
1078  return false;
1079  }
1080  }
1081  else if (tgtOutSrcPi0 != -1 && tgtOutSrcEi1 != -1)
1082  {
1083  // (reverse of previous clause)
1084  if (!outEdgeToOutCorner(tgtEi, tgtOutSrcEi1, tgtOutSrcPi0, 1))
1085  {
1086  pointLocations.clear();
1087  return false;
1088  }
1089  }
1090  else if (tgtOutSrcPi0 != -1 && tgtOutSrcPi1 != -1)
1091  {
1092  // Both target points are outside source corners
1093  if (!outCornerToOutCorner(tgtEi, tgtOutSrcPi0, tgtOutSrcPi1))
1094  {
1095  pointLocations.clear();
1096  return false;
1097  }
1098  }
1099  else
1100  {
1101  // A target point is outside all source edges. The projection
1102  // has collapsed.
1103  pointLocations.clear();
1104  return false;
1105  }
1106  }
1107 
1108  // Step 4: Complete the polygon by adding any remaining source points that
1109  // were not traversed during the walk of the target edges
1110  if (!pointLocations.empty())
1111  {
1112  const location& l = pointLocations.first();
1113 
1114  if (l.isIntersection() || l.isSrcAndTgtPoint())
1115  {
1116  addPointLocations(l, location(), false);
1117  }
1118  }
1119  else
1120  {
1121  forAllReverse(srcInTgtEdge, srcPi)
1122  {
1123  pointLocations.append(location::srcPoint(srcPi));
1124  }
1125  }
1126 
1127  // Step 5: The above walk was done around the target triangle, but the
1128  // result should be ordered in the direction of the source triangle, so the
1129  // list of locations must be reversed
1130  inplaceReverseList(pointLocations);
1131 
1132  return true;
1133 }
1134 
1135 
1136 //- Construct the intersection geometry
1138 (
1139  const FixedList<point, 3>& srcPs,
1140  const FixedList<vector, 3>& srcNs,
1141  const FixedList<point, 3>& tgtPs,
1142  DynamicList<point>& srcPoints,
1143  DynamicList<vector>& srcPointNormals,
1144  DynamicList<point>& tgtPoints,
1145  const DynamicList<location>& pointLocations
1146 )
1147 {
1148  srcPoints.resize(pointLocations.size());
1149  srcPointNormals.resize(pointLocations.size());
1150  tgtPoints.resize(pointLocations.size());
1151 
1152  forAll(pointLocations, pointi)
1153  {
1154  const location& l = pointLocations[pointi];
1155 
1156  if (l.isSrcAndTgtPoint())
1157  {
1158  const point& srcP = srcPs[l.srcPointi()];
1159  const vector& srcN = srcNs[l.srcPointi()];
1160  const point& tgtP = tgtPs[l.tgtPointi()];
1161 
1162  srcPoints[pointi] = srcP;
1163  srcPointNormals[pointi] = srcN;
1164  tgtPoints[pointi] = tgtP;
1165  }
1166  else if (l.isSrcPoint())
1167  {
1168  const point& srcP = srcPs[l.srcPointi()];
1169  const vector& srcN = srcNs[l.srcPointi()];
1170  barycentric2D tgtTs =
1171  srcPointTgtTriIntersection(srcP, srcN, tgtPs);
1172 
1173  // Force inside the target triangle
1174  if (cmptMin(tgtTs) < 0)
1175  {
1176  const direction iMin = findMin(tgtTs);
1177  const direction iMax = findMax(tgtTs);
1178  const direction iMid = 3 - iMin - iMax;
1179 
1180  if (tgtTs[iMid] < 0)
1181  {
1182  tgtTs[iMin] = 0;
1183  tgtTs[iMax] = 1;
1184  tgtTs[iMid] = 0;
1185  }
1186  else
1187  {
1188  const scalar t = tgtTs[iMax] + tgtTs[iMid];
1189  tgtTs[iMin] = 0;
1190  tgtTs[iMax] /= t;
1191  tgtTs[iMid] /= t;
1192  }
1193  }
1194 
1195  srcPoints[pointi] = srcP;
1196  srcPointNormals[pointi] = srcN;
1197  tgtPoints[pointi] = tgtTriInterpolate(tgtTs, tgtPs);
1198  }
1199  else if (l.isTgtPoint())
1200  {
1201  const point& tgtP = tgtPs[l.tgtPointi()];
1202  const barycentric2D srcTs =
1203  srcTriTgtPointIntersection(srcPs, srcNs, tgtP);
1204 
1205  srcPoints[pointi] = srcTriInterpolate(srcTs, srcPs);
1206  srcPointNormals[pointi] = srcTriInterpolate(srcTs, srcNs);
1207  tgtPoints[pointi] = tgtP;
1208  }
1209  else // if (l.isIntersection())
1210  {
1211  const label srcPi0 = l.srcEdgei(), srcPi1 = (srcPi0 + 1) % 3;
1212  const label tgtPi0 = l.tgtEdgei(), tgtPi1 = (tgtPi0 + 1) % 3;
1213 
1214  const Pair<scalar> ts =
1216  (
1217  {srcPs[srcPi0], srcPs[srcPi1]},
1218  {srcNs[srcPi0], srcNs[srcPi1]},
1219  {tgtPs[tgtPi0], tgtPs[tgtPi1]}
1220  );
1221  const scalar srcT = ts.first(), tgtT = ts.second();
1222 
1223  srcPoints[pointi] =
1224  (1 - srcT)*srcPs[srcPi0] + srcT*srcPs[srcPi1];
1225  srcPointNormals[pointi] =
1226  (1 - srcT)*srcNs[srcPi0] + srcT*srcNs[srcPi1];
1227  tgtPoints[pointi] =
1228  (1 - tgtT)*tgtPs[tgtPi0] + tgtT*tgtPs[tgtPi1];
1229  }
1230  }
1231 }
1232 
1233 
1234 } // End namespace triIntersect
1235 } // End namespace Foam
1236 
1237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1238 
1240 (
1241  const word& name,
1242  const FixedList<point, 3>& srcPs,
1243  const FixedList<vector, 3>& srcNs,
1244  const label nEdge,
1245  const label nNormal,
1246  const scalar lNormal
1247 )
1248 {
1249  scalar lengthScale = 0;
1250  for (label i = 0; i < 3; ++ i)
1251  {
1252  lengthScale = max(lengthScale, mag(srcPs[i] - srcPs[(i + 1) % 3]));
1253  }
1254 
1255  const label nu = nEdge, nv = nNormal;
1256  const scalar u0 = 0, u1 = 1;
1257  const scalar v0 = -lNormal/2*lengthScale, v1 = lNormal/2*lengthScale;
1258 
1259  pointField ps(3*(nu + 1)*(nv + 1));
1260  for (label i = 0; i < 3; ++ i)
1261  {
1262  const point& p0 = srcPs[i], & p1 = srcPs[(i + 1) % 3];
1263  const vector& n0 = srcNs[i], & n1 = srcNs[(i + 1) % 3];
1264 
1265  for (label iu = 0; iu <= nu; ++ iu)
1266  {
1267  const scalar u = u0 + (u1 - u0)*scalar(iu)/nu;
1268  for (label iv = 0; iv <= nv; ++ iv)
1269  {
1270  const scalar v = v0 + (v1 - v0)*scalar(iv)/nv;
1271  const vector x = p0 + (p1 - p0)*u + (n0 + (n1 - n0)*u)*v;
1272  ps[i*(nu + 1)*(nv + 1) + iu*(nv + 1) + iv] = x;
1273  }
1274  }
1275  }
1276 
1277  faceList fs(3*nu*nv);
1278  for (label i = 0; i < 3; ++ i)
1279  {
1280  for (label iu = 0; iu < nu; ++ iu)
1281  {
1282  for (label iv = 0; iv < nv; ++ iv)
1283  {
1284  fs[i*nu*nv + iu*nv + iv] =
1285  face
1286  ({
1287  i*(nu + 1)*(nv + 1) + (nv + 1)*iu + iv,
1288  i*(nu + 1)*(nv + 1) + (nv + 1)*iu + iv + 1,
1289  i*(nu + 1)*(nv + 1) + (nv + 1)*(iu + 1) + iv + 1,
1290  i*(nu + 1)*(nv + 1) + (nv + 1)*(iu + 1) + iv
1291  });
1292  }
1293  }
1294  }
1295 
1296  Info<< indent << "Writing face to " << name + ".vtk" << endl;
1298  (
1299  name + ".vtk",
1300  name,
1301  false,
1302  ps,
1303  labelList(),
1304  labelListList(),
1305  fs
1306  );
1307 }
1308 
1309 
1311 (
1312  const Pair<point>& srcPs,
1313  const Pair<vector>& srcNs,
1314  const point& tgtP
1315 )
1316 {
1317  const tensor A(srcPs[1] - srcPs[0], srcNs[0], srcNs[1]);
1318  const scalar detA = det(A);
1319 
1320  const tensor T(A.y()^A.z(), A.z()^A.x(), A.x()^A.y());
1321 
1322  const vector detAY = T & (tgtP - srcPs[0]);
1323  const scalar Yx = protectedDivideAndClip01(detAY.x(), detA);
1324 
1325  const scalar offset = Yx*detAY.y() - (1 - Yx)*detAY.z();
1326 
1327  return offset == 0 ? - vSmall : offset;
1328 }
1329 
1330 
1332 (
1333  const Pair<point>& srcPs,
1334  const Pair<vector>& srcNs,
1335  const point& tgtP,
1336  const bool srcDirection
1337 )
1338 {
1339  if (srcDirection)
1340  {
1341  return srcEdgeTgtPointOffset(srcPs, srcNs, tgtP);
1342  }
1343  else
1344  {
1345  return - srcEdgeTgtPointOffset(reverse(srcPs), reverse(srcNs), tgtP);
1346  }
1347 }
1348 
1349 
1351 (
1352  const point& srcP,
1353  const vector& srcN,
1354  const Pair<point>& tgtPs
1355 )
1356 {
1357  const tensor A(srcN, tgtPs[1] - tgtPs[0], srcN^(tgtPs[1] - tgtPs[0]));
1358  const scalar detA = det(A);
1359 
1360  const tensor T(A.y()^A.z(), A.z()^A.x(), A.x()^A.y());
1361 
1362  const vector detAY = T & (tgtPs[0] - srcP);
1363 
1364  const scalar offset = protectedDivide(detAY.z(), detA);
1365 
1366  return offset == 0 ? - vSmall : offset;
1367 }
1368 
1369 
1371 (
1372  const point& srcP,
1373  const vector& srcN,
1374  const Pair<point>& tgtPs,
1375  const bool tgtDirection
1376 )
1377 {
1378  if (tgtDirection)
1379  {
1380  return srcPointTgtEdgeOffset(srcP, srcN, tgtPs);
1381  }
1382  else
1383  {
1384  return - srcPointTgtEdgeOffset(srcP, srcN, reverse(tgtPs));
1385  }
1386 }
1387 
1388 
1390 (
1391  const Pair<point>& srcPs,
1392  const Pair<vector>& srcNs,
1393  const Pair<point>& tgtPs
1394 )
1395 {
1396  scalar srcT, tgtT;
1397 
1398  #ifndef BISECT
1401  (
1402  srcPs[0] - tgtPs[0],
1403  srcPs[1] - srcPs[0],
1404  tgtPs[0] - tgtPs[1],
1405  srcNs[0],
1406  Zero,
1407  srcNs[1] - srcNs[0],
1408  {0, 1, -1}
1409  );
1410  if (solution.first())
1411  {
1412  // If the analytical solution succeeds, then use the result
1413  srcT = solution.second().x();
1414  tgtT = solution.second().y();
1415  }
1416  else
1417  #endif
1418  {
1419  // If the analytical solution fails, then solve by bisection
1420 
1421  // !!! This method, whilst elegant, isn't sufficiently robust. The
1422  // srcPointTgtEdgeOffset calculation is not as reliable as
1423  // srcEdgeTgtPointOffset. So, we can't bisect the source edge to get
1424  // srcT and then reuse the solveProjectionGivenT stuff. We need to
1425  // bisect the target edge to get tgtT, and then have a specific process
1426  // for calculating srcT from tgtT.
1427  /*
1428  scalar srcT0 = 0, srcT1 = 1;
1429 
1430  const scalar o0 = srcPointTgtEdgeOffset(srcPs[0], srcNs[0], tgtPs);
1431  const scalar o1 = srcPointTgtEdgeOffset(srcPs[1], srcNs[1], tgtPs);
1432 
1433  const scalar s = o0 > o1 ? +1 : -1;
1434 
1435  for (label i = 0; i < ceil(std::log2(1/small)); ++ i)
1436  {
1437  const scalar srcT = (srcT0 + srcT1)/2;
1438  const vector srcP = (1 - srcT)*srcPs[0] + srcT*srcPs[1];
1439  const vector srcN = (1 - srcT)*srcNs[0] + srcT*srcNs[1];
1440 
1441  const scalar o = s*srcPointTgtEdgeOffset(srcP, srcN, tgtPs);
1442 
1443  if (o > 0)
1444  {
1445  srcT0 = srcT;
1446  }
1447  else
1448  {
1449  srcT1 = srcT;
1450  }
1451  }
1452 
1453  srcT = (srcT0 + srcT1)/2;
1454  tgtT =
1455  solveProjectionGivenT
1456  (
1457  srcPs[0] - tgtPs[0],
1458  srcPs[1] - srcPs[0],
1459  tgtPs[0] - tgtPs[1],
1460  srcNs[0],
1461  Zero,
1462  srcNs[1] - srcNs[0],
1463  {0, 1, -1},
1464  srcT
1465  ).y();
1466  */
1467 
1468  // !!! This method appears robust
1469  scalar tgtT0 = 0, tgtT1 = 1;
1470 
1471  const scalar o0 = srcEdgeTgtPointOffset(srcPs, srcNs, tgtPs[0]);
1472  const scalar o1 = srcEdgeTgtPointOffset(srcPs, srcNs, tgtPs[1]);
1473 
1474  const scalar s = o0 > o1 ? +1 : -1;
1475 
1476  for (label i = 0; i < ceil(std::log2(1/small)); ++ i)
1477  {
1478  const scalar tgtT = (tgtT0 + tgtT1)/2;
1479  const vector tgtP = (1 - tgtT)*tgtPs[0] + tgtT*tgtPs[1];
1480 
1481  const scalar o = s*srcEdgeTgtPointOffset(srcPs, srcNs, tgtP);
1482 
1483  if (o > 0)
1484  {
1485  tgtT0 = tgtT;
1486  }
1487  else
1488  {
1489  tgtT1 = tgtT;
1490  }
1491  }
1492 
1493  tgtT = (tgtT0 + tgtT1)/2;
1494 
1495  // Solve for the corresponding source edge coordinate
1496  const vector srcDP = srcPs[1] - srcPs[0];
1497  const vector tgtP = (1 - tgtT)*tgtPs[0] + tgtT*tgtPs[1];
1498 
1499  const tensor A(srcDP, srcNs[0], srcNs[1]);
1500  const scalar detA = det(A);
1501 
1502  const vector Tx = A.y()^A.z();
1503 
1504  const scalar magDetAYx = sign(detA)*(Tx & (tgtP - srcPs[0]));
1505 
1506  const scalar srcTStar = protectedDivideAndClip01(magDetAYx, mag(detA));
1507 
1508  const vector srcN = (1 - srcTStar)*srcNs[0] + srcTStar*srcNs[1];
1509 
1510  const vector srcDPPerpN = srcDP - (srcDP & srcN)*srcN;
1511 
1512  srcT =
1514  (
1515  (tgtP - srcPs[0]) & srcDPPerpN,
1516  srcDP & srcDPPerpN
1517  );
1518  }
1519 
1520  /*
1521  // Check that the points match
1522  {
1523  const point srcP = (1 - srcT)*srcPs[0] + srcT*srcPs[1];
1524  const point srcN = (1 - srcT)*srcNs[0] + srcT*srcNs[1];
1525  const point tgtP = (1 - tgtT)*tgtPs[0] + tgtT*tgtPs[1];
1526  const scalar srcU = ((tgtP - srcP) & srcN)/magSqr(srcN);
1527  const point srcQ = srcP + srcU*srcN;
1528  Info<< "srcT=" << srcT << ", tgtT=" << tgtT
1529  << ", err=" << magSqr(srcQ - tgtP) << endl;
1530  }
1531  */
1532 
1533  return Pair<scalar>(srcT, tgtT);
1534 }
1535 
1536 
1538 (
1539  const FixedList<point, 3>& srcPs,
1540  const FixedList<vector, 3>& srcNs,
1541  const point& tgtP
1542 )
1543 {
1544  auto srcPN = []
1545  (
1546  const FixedList<vector, 3>& srcPNs,
1547  const vector2D& srcTs
1548  )
1549  {
1550  const scalar srcT0 = 1 - srcTs.x() - srcTs.y();
1551  return srcT0*srcPNs[0] + srcTs.x()*srcPNs[1] + srcTs.y()*srcPNs[2];
1552  };
1553 
1554  auto srcEdgePNs = [srcPN]
1555  (
1556  const FixedList<vector, 3>& srcPNs,
1557  const vector2D& srcTs0,
1558  const vector2D& srcTs1
1559  )
1560  {
1561  return Pair<vector>(srcPN(srcPNs, srcTs0), srcPN(srcPNs, srcTs1));
1562  };
1563 
1564  auto offset = [&](const vector2D& srcTs0, const vector2D& srcTs1)
1565  {
1566  return srcEdgeTgtPointOffset
1567  (
1568  srcEdgePNs(srcPs, srcTs0, srcTs1),
1569  srcEdgePNs(srcNs, srcTs0, srcTs1),
1570  tgtP
1571  );
1572  };
1573 
1574  const scalar oA = offset(vector2D(0, 0), vector2D(1, 0));
1575  const scalar oB = offset(vector2D(1, 0), vector2D(0, 1));
1576  const scalar oC = offset(vector2D(0, 1), vector2D(0, 0));
1577  const FixedList<scalar, 3> offsets({oA, oB, oC});
1578 
1579  // If inside the triangle (or outside if the triangle is inverted) ...
1580  if (offsets[findMin(offsets)] >= 0 || offsets[findMax(offsets)] <= 0)
1581  {
1582  scalar srcT, srcU;
1583 
1584  #ifndef BISECT
1587  (
1588  srcPs[0] - tgtP,
1589  srcNs[0],
1590  srcPs[1] - srcPs[0],
1591  srcPs[2] - srcPs[0],
1592  srcNs[1] - srcNs[0],
1593  srcNs[2] - srcNs[0],
1594  {-1, 0, 0}
1595  );
1596  if (solution.first())
1597  {
1598  // If the analytical solution succeeds, then use the result
1599  srcT = solution.second().y();
1600  srcU = solution.second().z();
1601  }
1602  else
1603  #endif
1604  {
1605  // If the analytical solution fails, then solve by bisection
1606  const scalar sign = offsets[findMin(offsets)] > 0 ? +1 : -1;
1607 
1608  vector2D srcTsA(0, 0), srcTsB(1, 0), srcTsC(0, 1);
1609 
1610  for (label i = 0; i < ceil(std::log2(1/small)); ++ i)
1611  {
1612  const vector2D srcTsAB = (srcTsA + srcTsB)/2;
1613  const vector2D srcTsBC = (srcTsB + srcTsC)/2;
1614  const vector2D srcTsCA = (srcTsC + srcTsA)/2;
1615 
1616  const scalar oA = sign*offset(srcTsCA, srcTsAB);
1617  const scalar oB = sign*offset(srcTsAB, srcTsBC);
1618  const scalar oC = sign*offset(srcTsBC, srcTsCA);
1619  const FixedList<scalar, 3> offsets({oA, oB, oC});
1620 
1621  const label offsetMini = findMin(offsets);
1622  if (offsets[offsetMini] > 0)
1623  {
1624  srcTsA = srcTsAB;
1625  srcTsB = srcTsBC;
1626  srcTsC = srcTsCA;
1627  }
1628  else if (offsetMini == 0)
1629  {
1630  srcTsC = srcTsCA;
1631  srcTsB = srcTsAB;
1632  }
1633  else if (offsetMini == 1)
1634  {
1635  srcTsA = srcTsAB;
1636  srcTsC = srcTsBC;
1637  }
1638  else if (offsetMini == 2)
1639  {
1640  srcTsB = srcTsBC;
1641  srcTsA = srcTsCA;
1642  }
1643  }
1644 
1645  srcT = (srcTsA[0] + srcTsB[0] + srcTsC[0])/3;
1646  srcU = (srcTsA[1] + srcTsB[1] + srcTsC[1])/3;
1647  }
1648 
1649  return barycentric2D(1 - srcT - srcU, srcT, srcU);
1650  }
1651 
1652  // If outside an edge ...
1653  forAll(srcPs, srcEi)
1654  {
1655  const label srcEi0 = (srcEi + 2) % 3, srcEi1 = (srcEi + 1) % 3;
1656 
1657  if
1658  (
1659  offsets[srcEi] <= 0
1660  && offsets[srcEi0] >= 0
1661  && offsets[srcEi1] >= 0
1662  )
1663  {
1664  const label srcPi0 = srcEi, srcPi1 = (srcEi + 1) % 3;
1665  const label srcPiOpp = (srcEi + 2) % 3;
1666 
1667  scalar srcT, srcU;
1668 
1669  #ifndef BISECT
1672  (
1673  srcPs[srcPi0] - tgtP,
1674  srcPs[srcPi1] - srcPs[srcPi0],
1675  srcPs[srcPi0] - srcPs[srcPiOpp],
1676  srcNs[srcPi0],
1677  srcPs[srcPi1] - srcPs[srcPi0],
1678  srcNs[srcPi1] - srcNs[srcPi0],
1679  {0, -1, -1}
1680  );
1681  if (solution.first())
1682  {
1683  // If the analytical solution succeeds, then use the result
1684  srcT = solution.second().x();
1685  srcU = solution.second().y();
1686  }
1687  else
1688  #endif
1689  {
1690  // If the analytical solution fails, then solve by bisection
1691  const vector2D srcTsOpp(srcPiOpp == 1, srcPiOpp == 2);
1692  const vector2D srcTs0(srcPi0 == 1, srcPi0 == 2);
1693  const vector2D srcTs1(srcPi1 == 1, srcPi1 == 2);
1694 
1695  scalar srcT0 = 0, srcT1 = 1;
1696 
1697  for (label i = 0; i < ceil(std::log2(1/small)); ++ i)
1698  {
1699  const scalar srcT = (srcT0 + srcT1)/2;
1700  const vector2D srcTs01(srcTs0*(1 - srcT) + srcTs1*srcT);
1701 
1702  const scalar o = offset(srcTsOpp, srcTs01);
1703 
1704  if (o > 0)
1705  {
1706  srcT0 = srcT;
1707  }
1708  else
1709  {
1710  srcT1 = srcT;
1711  }
1712  }
1713 
1714  srcT = (srcT0 + srcT1)/2;
1715  srcU =
1717  (
1718  srcPs[srcPi0] - tgtP,
1719  srcPs[srcPi1] - srcPs[srcPi0],
1720  srcPs[srcPi0] - srcPs[srcPiOpp],
1721  srcNs[srcPi0],
1722  srcPs[srcPi1] - srcPs[srcPi0],
1723  srcNs[srcPi1] - srcNs[srcPi0],
1724  {0, -1, -1},
1725  srcT
1726  ).y();
1727  }
1728 
1729  // Convert to the triangle's coordinate system
1730  barycentric2D y;
1731  y[srcPiOpp] = - srcU;
1732  y[srcPi0] = (1 + srcU)*(1 - srcT);
1733  y[srcPi1] = (1 + srcU)*srcT;
1734  return y;
1735  }
1736  }
1737 
1738  // If outside a corner ...
1739  forAll(srcPs, srcPi)
1740  {
1741  const label srcEiOpp = (srcPi + 1) % 3;
1742  const label srcEi0 = (srcPi + 2) % 3, srcEi1 = srcPi;
1743 
1744  if
1745  (
1746  offsets[srcEiOpp] >= 0
1747  && offsets[srcEi0] <= 0
1748  && offsets[srcEi1] <= 0
1749  )
1750  {
1751  // Solve for the intersection coordinates directly
1752  const label srcPi0 = (srcPi + 2) % 3, srcPi1 = (srcPi + 1) % 3;
1753 
1754  const tensor A
1755  (
1756  srcPs[srcPi] - srcPs[srcPi0],
1757  srcPs[srcPi] - srcPs[srcPi1],
1758  srcNs[srcPi]
1759  );
1760  const vector T0(A.y()^A.z()), T1(A.z()^A.x());
1761  const scalar detA = A.x() & T0;
1762 
1763  const scalar srcT =
1764  protectedDivide(T0 & (tgtP - srcPs[srcPi]), detA);
1765  const scalar srcU =
1766  protectedDivide(T1 & (tgtP - srcPs[srcPi]), detA);
1767 
1768  // Convert to the triangle's coordinate system
1769  barycentric2D y;
1770  y[srcPi0] = - srcT;
1771  y[srcPi1] = - srcU;
1772  y[srcPi] = 1 + srcT + srcU;
1773  return y;
1774  }
1775  }
1776 
1777  // Above logic means we should never reach here
1779  << "Point " << tgtP << " could not be classified within triangle "
1780  << srcPs << " with projection normals " << srcNs << exit(FatalError);
1781  return barycentric2D::uniform(NaN);
1782 }
1783 
1784 
1786 (
1787  const point& srcP,
1788  const vector& srcN,
1789  const FixedList<point, 3>& tgtPs
1790 )
1791 {
1792  const tensor A(tgtPs[1] - tgtPs[0], tgtPs[2] - tgtPs[0], - srcN);
1793  const scalar detA = det(A);
1794 
1795  const vector T0(A.y()^A.z()), T1(A.z()^A.x());
1796  const tensor T(- T0 - T1, T0, T1);
1797 
1798  const vector detAY = (T & (srcP - tgtPs[0])) + vector(detA, 0, 0);
1799  const scalar maxMagDetAY = mag(detAY[findMax(cmptMag(detAY))]);
1800 
1801  // Usual case. The source normal is not parallel to the target triangle.
1802  // The intersection is a single unambiguous point.
1803  if (maxMagDetAY/vGreat < mag(detA))
1804  {
1805  const vector y = detAY/detA;
1806  return barycentric2D(y.x(), y.y(), y.z());
1807  }
1808 
1809  // Degenerate case. The source normal is parallel to, and the source point
1810  // is out of the plane of, the target triangle. Really, there is no
1811  // intersection, but for the purposes of this function we can say there is
1812  // an intersection and it is arbitrarily far away.
1813  if (maxMagDetAY > 0)
1814  {
1815  const vector y = detAY/maxMagDetAY*vGreat;
1816  return barycentric2D(y.x(), y.y(), y.z());
1817  }
1818 
1819  const tensor2D A2
1820  (
1821  A.x() & A.x(), A.x() & A.y(),
1822  A.y() & A.x(), A.y() & A.y()
1823  );
1824  const scalar detA2(det(A2));
1825 
1826  const tensor2D T2(cof(A2));
1827 
1828  const vector2D detAY2 =
1829  T2 & vector2D(A.x() & (srcP - tgtPs[0]), A.y() & (srcP - tgtPs[0]));
1830  const scalar maxMagDetAY2 = mag(detAY[findMax(cmptMag(detAY2))]);
1831 
1832  // Very degenerate case. The source normal is parallel to, and the source
1833  // point is on the plane of, the target triangle. The intersection is a
1834  // line. Choose the point on the line corresponding to the source point.
1835  if (maxMagDetAY2/vGreat < mag(detA2))
1836  {
1837  const vector2D y2 = detAY2/detA2;
1838  return barycentric2D(1 - cmptSum(y2), y2.x(), y2.y());
1839  }
1840 
1841  // Most degenerate case. The target triangle has collapsed to a line.
1842  // Choose an arbitrary point a long way away. It's possible there's more we
1843  // could do here, but a need has yet to present itself.
1844  return barycentric2D(-vGreat, vGreat, vGreat);
1845 }
1846 
1847 
1849 (
1850  const FixedList<point, 3>& srcPs,
1851  const FixedList<vector, 3>& srcNs,
1852  const FixedList<bool, 3>& srcOwns,
1853  const FixedList<label, 3>& srcTgtPis,
1854  const FixedList<point, 3>& tgtPs,
1855  const FixedList<bool, 3>& tgtOwns,
1856  const FixedList<label, 3>& tgtSrcPis,
1857  DynamicList<point>& srcPoints,
1858  DynamicList<vector>& srcPointNormals,
1859  DynamicList<point>& tgtPoints,
1860  DynamicList<location>& pointLocations,
1861  const bool debug,
1862  const word& writePrefix
1863 )
1864 {
1865  const bool write = writePrefix != word::null;
1866 
1867  if (debug || write)
1868  {
1869  Info<< indent << "Intersecting triangles" << incrIndent << endl;
1870  }
1871 
1872  if (write)
1873  {
1874  writePolygon(writePrefix + "_srcTri", srcPs);
1875  writePolygon(writePrefix + "_tgtTri", tgtPs);
1876  writeTriProjection(writePrefix + "_srcPrj", srcPs, srcNs);
1877  }
1878 
1879  // Determine what source points lie within target edges and vice-versa
1880  FixedList<FixedList<label, 3>, 3> srcInTgtEdge, tgtInSrcEdge;
1881  FixedList<label, 3> srcInTgtTri, tgtInSrcTri;
1882  srcInTgt
1883  (
1884  srcPs, srcNs, srcOwns, srcTgtPis,
1885  tgtPs, tgtOwns,
1886  srcInTgtEdge,
1887  srcInTgtTri
1888  );
1889  tgtInSrc
1890  (
1891  srcPs, srcNs, srcOwns,
1892  tgtPs, tgtOwns, tgtSrcPis,
1893  tgtInSrcEdge,
1894  tgtInSrcTri
1895  );
1896 
1897  if (debug)
1898  {
1899  if (count(srcTgtPis, -1) != 3)
1900  {
1901  Info<< indent << "srcTgtPis=" << srcTgtPis << endl;
1902  }
1903  Info<< indent << "srcInTgtTri=" << srcInTgtTri << endl
1904  << indent << "srcInTgtEdge=" << srcInTgtEdge << endl;
1905  if (count(tgtSrcPis, -1) != 3)
1906  {
1907  Info<< indent << "tgtSrcPis=" << tgtSrcPis << endl;
1908  }
1909  Info<< indent << "tgtInSrcTri=" << tgtInSrcTri << endl
1910  << indent << "tgtInSrcEdge=" << tgtInSrcEdge << endl;
1911  }
1912 
1913  // Generate the locations
1915  (
1916  tgtSrcPis,
1917  srcInTgtEdge,
1918  tgtInSrcEdge,
1919  srcInTgtTri,
1920  tgtInSrcTri,
1921  pointLocations
1922  );
1923 
1924  // Generate the geometry
1926  (
1927  srcPs,
1928  srcNs,
1929  tgtPs,
1930  srcPoints,
1931  srcPointNormals,
1932  tgtPoints,
1933  pointLocations
1934  );
1935 
1936  if (write)
1937  {
1938  writePolygon(writePrefix + "_srcIctFace", srcPoints);
1939  writePolygon(writePrefix + "_tgtIctFace", tgtPoints);
1940  }
1941 
1942  if (debug || write)
1943  {
1944  Info<< decrIndent;
1945  }
1946 }
1947 
1948 // ************************************************************************* //
scalar Cv(const scalar p, const scalar T) const
Definition: HtoEthermo.H:2
scalar y
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:445
Templated 2D Barycentric derived from VectorSpace. Has 3 components, one of which is redundant.
Definition: Barycentric2D.H:53
Graphite solid properties.
Definition: C.H:51
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
void resize(const label)
Alter the addressed list size.
Definition: DynamicListI.H:216
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:236
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:78
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
const Type & second() const
Return second.
Definition: PairI.H:121
const Type & first() const
Return first.
Definition: PairI.H:107
Templated storage for the roots of polynomial equations, plus flags to indicate the nature of the roo...
Definition: Roots.H:67
void type(const direction i, const rootType t)
Set the type of the i-th root.
Definition: RootsI.H:120
Templated 2D tensor derived from VectorSpace adding construction from 4 components,...
Definition: Tensor2D.H:59
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
T * first()
Return the first entry.
Definition: UILList.H:109
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:74
T & first()
Return the first element of the list.
Definition: UListI.H:114
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
T & last()
Return the last element of the list.
Definition: UListI.H:128
const Cmpt & y() const
Definition: Vector2DI.H:74
const Cmpt & x() const
Definition: Vector2DI.H:68
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
const Cmpt & z() const
Definition: VectorI.H:87
const Cmpt & y() const
Definition: VectorI.H:81
const Cmpt & x() const
Definition: VectorI.H:75
Cubic equation of the form a*x^3 + b*x^2 + c*x + d = 0.
Definition: cubicEqn.H:52
Roots< 3 > roots(const bool real=false) const
Get the roots.
Definition: cubicEqn.C:32
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:51
@ BEGIN_LIST
Definition: token.H:109
@ END_LIST
Definition: token.H:110
bool isTgtNotSrcPoint() const
Return whether the location is a target point and not a source.
static location intersection(const label srcEi, const label tgtEi)
Construct an intersection location between a source and a.
bool isTgtPoint() const
Return whether the location is a target point.
label tgtEdgei() const
Return the target edge index.
static location srcTgtPoint(const label srcPi, const label tgtPi)
Construct a source and target point location.
bool isSrcNotTgtPoint() const
Return whether the location is a source point and not a target.
static location tgtPoint(const label tgtPi)
Construct a target point location.
label srcPointi() const
Return the source point index.
label tgtPointi() const
Return the target point index.
bool isSrcAndTgtPoint() const
Return whether the location is a source point and a target point.
bool isIntersection() const
Return whether the location is an intersection.
bool isSrcPoint() const
Return whether the location is a source point.
label srcEdgei() const
Return the source edge index.
static location srcPoint(const label srcPi)
Construct a source point location.
A triangle primitive used to calculate face areas and swept volumes.
Definition: triangle.H:80
vector normal() const
Return unit normal.
Definition: triangleI.H:110
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const char *const group
Group name for atomic constants.
static const coefficient B("B", dimless, 18.678)
static const coefficient A("A", dimPressure, 611.21)
barycentric2D srcPointTgtTriIntersection(const point &srcP, const vector &srcN, const FixedList< point, 3 > &tgtPs)
Calculate the intersection of a projected source point with a target.
scalar srcPointTgtEdgeOffset(const point &srcP, const vector &srcN, const Pair< point > &tgtPs)
Calculate the signed offset of a projected source point in relation to a.
void tgtInSrc(const FixedList< point, 3 > &srcPs, const FixedList< vector, 3 > &srcNs, const FixedList< bool, 3 > &srcOwns, const FixedList< point, 3 > &tgtPs, const FixedList< bool, 3 > &tgtOwns, const FixedList< label, 3 > &tgtSrcPis, FixedList< FixedList< label, 3 >, 3 > &tgtInSrcEdge, FixedList< label, 3 > &tgtInSrcTri)
Calculate whether the points of the given target triangle project inside or.
Definition: triIntersect.C:543
void generateGeometry(const FixedList< point, 3 > &srcPs, const FixedList< vector, 3 > &srcNs, const FixedList< point, 3 > &tgtPs, DynamicList< point > &srcPoints, DynamicList< vector > &srcPointNormals, DynamicList< point > &tgtPoints, const DynamicList< location > &pointLocations)
Construct the intersection geometry.
barycentric2D srcTriTgtPointIntersection(const FixedList< point, 3 > &srcPs, const FixedList< vector, 3 > &srcNs, const point &tgtP)
Calculate the intersection of a target point with a source triangle's.
Tuple2< bool, vector > solveProjection(const vector &C, const vector &Ct, const vector &Cu, const vector &Cv, const vector &Ctu, const vector &Ctv, const FixedList< label, 3 > groups)
Solve a projection equation.
Definition: triIntersect.C:182
void writeTriProjection(const word &name, const FixedList< point, 3 > &srcTriPs, const FixedList< vector, 3 > &srcTriNs, const label nEdge=20, const label nNormal=20, const scalar lNormal=0.5)
Write a VTK file of a triangle projection.
FixedList< FixedList< scalar, 3 >, 3 > tgtSrcEdgeOffset(const FixedList< point, 3 > &srcPs, const FixedList< vector, 3 > &srcNs, const FixedList< bool, 3 > &srcOwns, const FixedList< point, 3 > &tgtPs, const FixedList< bool, 3 > &tgtOwns)
Calculate the non-dimensional offsets the target points from the source.
Definition: triIntersect.C:381
void srcInTgt(const FixedList< point, 3 > &srcPs, const FixedList< vector, 3 > &srcNs, const FixedList< bool, 3 > &srcOwns, const FixedList< label, 3 > &srcTgtPis, const FixedList< point, 3 > &tgtPs, const FixedList< bool, 3 > &tgtOwns, FixedList< FixedList< label, 3 >, 3 > &srcInTgtEdge, FixedList< label, 3 > &srcInTgtTri)
Calculate whether the points of the given source triangle project inside or.
Definition: triIntersect.C:518
void thisIsOther(const FixedList< label, 3 > &thisOtherPis, FixedList< FixedList< label, 3 >, 3 > &thisInOtherEdge, FixedList< label, 3 > &thisInOtherTri)
Override results of the srcInTgt/tgtInSrc calculations with explicit.
Definition: triIntersect.C:492
FixedList< label, 3 > tgtInSrcTri(const FixedList< point, 3 > &srcPs, const FixedList< vector, 3 > &srcNs, const FixedList< point, 3 > &tgtPs, const FixedList< FixedList< label, 3 >, 3 > &tgtInSrcEdge)
Construct target-point-inside/outside-source-triangle topology from a set.
Definition: triIntersect.C:459
FixedList< FixedList< label, 3 >, 3 > thisInOtherEdge(const FixedList< FixedList< scalar, 3 >, 3 > &thisOtherEdgeOffset)
Construct point-inside/outside-edge topology from a set of point-edge.
Definition: triIntersect.C:415
Ostream & operator<<(Ostream &os, const FixedList< FixedList< Type, 3 >, 3 > &l)
Print 3x3 FixedListLists on one line.
Definition: triIntersect.C:49
bool generateLocations(const FixedList< label, 3 > &tgtSrcPis, const FixedList< FixedList< label, 3 >, 3 > &srcInTgtEdge, const FixedList< FixedList< label, 3 >, 3 > &tgtInSrcEdge, const FixedList< label, 3 > &srcInTgtTri, const FixedList< label, 3 > &tgtInSrcTri, DynamicList< location > &pointLocations)
Construct the intersection topology.
Definition: triIntersect.C:675
vector clipped01(const vector x, const FixedList< label, 3 > groups)
Clip the given vector between values of 0 and 1, and also clip one minus.
Definition: triIntersect.C:66
const scalar maxDot
The maximum dot product between a source point normal and a target plane.
Definition: triIntersect.C:44
FixedList< label, 3 > thisInOtherTri(const FixedList< FixedList< label, 3 >, 3 > &thisInOtherEdge)
Construct point-inside/outside-triangle topology from a set of.
Definition: triIntersect.C:438
scalar srcEdgeTgtPointOffset(const Pair< point > &srcPs, const Pair< vector > &srcNs, const point &tgtP)
Calculate the signed offset of a target point in relation to a projected.
vector solveProjectionGivenT(const vector &C, const vector &Ct, const vector &Cu, const vector &Cv, const vector &Ctu, const vector &Ctv, const FixedList< label, 3 > groups, const scalar t)
Solve a projection equation given a value of the t variable.
Definition: triIntersect.C:138
Type tgtTriInterpolate(const barycentric2D &y, const FixedList< Type, 3 > &tgtPsis)
Use the coordinates obtained from srcPointTgtTriIntersection to interpolate.
void writePolygon(const word &name, const PointField &ps)
Write a VTK file of a polygon.
FixedList< FixedList< scalar, 3 >, 3 > srcTgtEdgeOffset(const FixedList< point, 3 > &srcPs, const FixedList< vector, 3 > &srcNs, const FixedList< bool, 3 > &srcOwns, const FixedList< point, 3 > &tgtPs, const FixedList< bool, 3 > &tgtOwns)
Calculate the non-dimensional offsets of the source points from the target.
Definition: triIntersect.C:286
bool orderLocations(const UList< location > &locations, bool isSrcEdge, const label i0, label &nVisited, boolList &visited, labelList &order)
Order intersection locations into a polygon.
Definition: triIntersect.C:567
Pair< scalar > srcEdgeTgtEdgeIntersection(const Pair< point > &srcPs, const Pair< vector > &srcNs, const Pair< point > &tgtPs)
Calculate the intersection of a target edge with a source edge's.
Type srcTriInterpolate(const barycentric2D &y, const FixedList< Type, 3 > &tgtPsis)
Use the coordinates obtained from srcTriTgtPointIntersection to interpolate.
void intersectTris(const FixedList< point, 3 > &srcPs, const FixedList< vector, 3 > &srcNs, const FixedList< bool, 3 > &srcOwns, const FixedList< label, 3 > &srcTgtPis, const FixedList< point, 3 > &tgtPs, const FixedList< bool, 3 > &tgtOwns, const FixedList< label, 3 > &tgtSrcPis, DynamicList< point > &srcPoints, DynamicList< vector > &srcPointNormals, DynamicList< point > &tgtPoints, DynamicList< location > &pointLocations, const bool debug, const word &writePrefix=word::null)
Construct the intersection of a source triangle's projected volume and a.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
Vector2D< scalar > vector2D
vector2D obtained from generic Vector2D
Definition: vector2D.H:49
Scalar protectedDivideAndClip01(const Scalar a, const Scalar b)
Divide two numbers, protect the result from overflowing, and clip the.
Definition: Scalar.H:389
Scalar protectedDivide(const Scalar a, const Scalar b)
Divide two numbers and protect the result from overflowing.
Definition: Scalar.H:381
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
scalar degToRad(const scalar deg)
Convert degrees to radians.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:242
List< label > labelList
A List of labels.
Definition: labelList.H:56
label log2(label i)
Return the log base 2 by successive bit-shifting of the given label.
Definition: label.H:79
dimensionedScalar sign(const dimensionedScalar &ds)
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
label findMax(const ListType &, const label start=0)
Find index of max element (and larger than given element).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
void inplaceReverseList(ListType &list)
Inplace reversal of a list using Swap.
int order(const scalar s)
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
messageStream Info
Barycentric2D< scalar > barycentric2D
A scalar version of the templated Barycentric2D.
Definition: barycentric2D.H:45
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:235
void det(LagrangianPatchField< scalar > &f, const LagrangianPatchField< tensor > &f1)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void reverse(UList< T > &, const label n)
Definition: UListI.H:334
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedSymmTensor cof(const dimensionedSymmTensor &dt)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
label findMin(const ListType &, const label start=0)
Find index of min element (and less than given element).
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
void Swap(T &a, T &b)
Definition: Swap.H:43
label count(const ListType &l, typename ListType::const_reference x)
Count the number of occurrences of a value in a list.
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
Cmpt cmptSum(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:507
void add(LagrangianPatchField< typename typeOfSum< Type1, Type2 >::type > &f, const LagrangianPatchField< Type1 > &f1, const LagrangianPatchField< Type2 > &f2)
void offset(label &lst, const label o)
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:228
void cmptMag(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
dimensionedScalar cos(const dimensionedScalar &ds)
uint8_t direction
Definition: direction.H:45
scalar u0
scalar T0
Definition: createFields.H:22
Functions with which to perform an intersection of a pair of triangles; the source and target....