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