OGS
AnalyticalGeometry.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
5
6#include <Eigen/Geometry>
7#include <algorithm>
8#include <cmath>
9#include <limits>
10
11#include "BaseLib/StringTools.h"
13#include "PointVec.h"
14#include "Polyline.h"
15#include "predicates.h"
16
17namespace
18{
20 MathLib::Point3d const& c)
21{
22 return orient2d(const_cast<double*>(a.data()),
23 const_cast<double*>(b.data()),
24 const_cast<double*>(c.data()));
25}
26
28 MathLib::Point3d const& b,
29 MathLib::Point3d const& c)
30{
31 return orient2dfast(const_cast<double*>(a.data()),
32 const_cast<double*>(b.data()),
33 const_cast<double*>(c.data()));
34}
35} // namespace
36
37namespace GeoLib
38{
40 MathLib::Point3d const& p1,
41 MathLib::Point3d const& p2)
42{
43 double const orientation = getOrientation2d(p0, p1, p2);
44 if (orientation > 0)
45 {
46 return CCW;
47 }
48 if (orientation < 0)
49 {
50 return CW;
51 }
52 return COLLINEAR;
53}
54
56 MathLib::Point3d const& p1,
57 MathLib::Point3d const& p2)
58{
59 double const orientation = getOrientation2dFast(p0, p1, p2);
60 if (orientation > 0)
61 {
62 return CCW;
63 }
64 if (orientation < 0)
65 {
66 return CW;
67 }
68 return COLLINEAR;
69}
70
71bool parallel(Eigen::Vector3d v, Eigen::Vector3d w)
72{
73 const double eps(std::numeric_limits<double>::epsilon());
74 double const eps_squared = eps * eps;
75
76 // check degenerated cases
77 if (v.squaredNorm() < eps_squared)
78 {
79 return false;
80 }
81
82 if (w.squaredNorm() < eps_squared)
83 {
84 return false;
85 }
86
87 v.normalize();
88 w.normalize();
89
90 bool parallel(true);
91 if (std::abs(v[0] - w[0]) > eps)
92 {
93 parallel = false;
94 }
95 if (std::abs(v[1] - w[1]) > eps)
96 {
97 parallel = false;
98 }
99 if (std::abs(v[2] - w[2]) > eps)
100 {
101 parallel = false;
102 }
103
104 if (!parallel)
105 {
106 parallel = true;
107 // change sense of direction of v_normalised
108 v *= -1.0;
109 // check again
110 if (std::abs(v[0] - w[0]) > eps)
111 {
112 parallel = false;
113 }
114 if (std::abs(v[1] - w[1]) > eps)
115 {
116 parallel = false;
117 }
118 if (std::abs(v[2] - w[2]) > eps)
119 {
120 parallel = false;
121 }
122 }
123
124 return parallel;
125}
126
128 GeoLib::LineSegment const& s1,
129 GeoLib::Point& s)
130{
131 GeoLib::Point const& pa{s0.getBeginPoint()};
132 GeoLib::Point const& pb{s0.getEndPoint()};
133 GeoLib::Point const& pc{s1.getBeginPoint()};
134 GeoLib::Point const& pd{s1.getEndPoint()};
135
136 if (!isCoplanar(pa, pb, pc, pd))
137 {
138 return false;
139 }
140
141 auto const& a = s0.getBeginPoint().asEigenVector3d();
142 auto const& b = s0.getEndPoint().asEigenVector3d();
143 auto const& c = s1.getBeginPoint().asEigenVector3d();
144 auto const& d = s1.getEndPoint().asEigenVector3d();
145
146 Eigen::Vector3d const v = b - a;
147 Eigen::Vector3d const w = d - c;
148 Eigen::Vector3d const qp = c - a;
149 Eigen::Vector3d const pq = a - c;
150
151 double const eps = std::numeric_limits<double>::epsilon();
152 double const squared_eps = eps * eps;
153 // handle special cases here to avoid computing intersection numerical
154 if (qp.squaredNorm() < squared_eps || (d - a).squaredNorm() < squared_eps)
155 {
156 s = pa;
157 return true;
158 }
159 if ((c - b).squaredNorm() < squared_eps ||
160 (d - b).squaredNorm() < squared_eps)
161 {
162 s = pb;
163 return true;
164 }
165
166 auto isLineSegmentIntersectingAB =
167 [&v](Eigen::Vector3d const& ap, std::size_t i)
168 {
169 // check if p is located at v=(a,b): (ap = t*v, t in [0,1])
170 return 0.0 <= ap[i] / v[i] && ap[i] / v[i] <= 1.0;
171 };
172
173 if (parallel(v, w))
174 { // original line segments (a,b) and (c,d) are parallel
175 if (parallel(pq, v))
176 { // line segment (a,b) and (a,c) are also parallel
177 // Here it is already checked that the line segments (a,b) and (c,d)
178 // are parallel. At this point it is also known that the line
179 // segment (a,c) is also parallel to (a,b). In that case it is
180 // possible to express c as c(t) = a + t * (b-a) (analog for the
181 // point d). Since the evaluation of all three coordinate equations
182 // (x,y,z) have to lead to the same solution for the parameter t it
183 // is sufficient to evaluate t only once.
184
185 // Search id of coordinate with largest absolute value which is will
186 // be used in the subsequent computations. This prevents division by
187 // zero in case the line segments are parallel to one of the
188 // coordinate axis.
189 std::size_t i_max(std::abs(v[0]) <= std::abs(v[1]) ? 1 : 0);
190 i_max = std::abs(v[i_max]) <= std::abs(v[2]) ? 2 : i_max;
191 if (isLineSegmentIntersectingAB(qp, i_max))
192 {
193 s = pc;
194 return true;
195 }
196 Eigen::Vector3d const ad = d - a;
197 if (isLineSegmentIntersectingAB(ad, i_max))
198 {
199 s = pd;
200 return true;
201 }
202 return false;
203 }
204 return false;
205 }
206
207 // general case
208 const double sqr_len_v(v.squaredNorm());
209 const double sqr_len_w(w.squaredNorm());
210
211 Eigen::Matrix2d mat;
212 mat(0, 0) = sqr_len_v;
213 mat(0, 1) = -v.dot(w);
214 mat(1, 1) = sqr_len_w;
215 mat(1, 0) = mat(0, 1);
216
217 Eigen::Vector2d rhs{v.dot(qp), w.dot(pq)};
218
219 rhs = mat.partialPivLu().solve(rhs);
220
221 // no theory for the following tolerances, determined by testing
222 // lower tolerance: little bit smaller than zero
223 const double l(-1.0 * std::numeric_limits<float>::epsilon());
224 // upper tolerance a little bit greater than one
225 const double u(1.0 + std::numeric_limits<float>::epsilon());
226 if (rhs[0] < l || u < rhs[0] || rhs[1] < l || u < rhs[1])
227 {
228 return false;
229 }
230
231 // compute points along line segments with minimal distance
232 GeoLib::Point const p0(a[0] + rhs[0] * v[0], a[1] + rhs[0] * v[1],
233 a[2] + rhs[0] * v[2]);
234 GeoLib::Point const p1(c[0] + rhs[1] * w[0], c[1] + rhs[1] * w[1],
235 c[2] + rhs[1] * w[2]);
236
237 double const min_dist(std::sqrt(MathLib::sqrDist(p0, p1)));
238 double const min_seg_len(
239 std::min(std::sqrt(sqr_len_v), std::sqrt(sqr_len_w)));
240 if (min_dist < min_seg_len * 1e-6)
241 {
242 s[0] = 0.5 * (p0[0] + p1[0]);
243 s[1] = 0.5 * (p0[1] + p1[1]);
244 s[2] = 0.5 * (p0[2] + p1[2]);
245 return true;
246 }
247
248 return false;
249}
250
254 GeoLib::Point& intersection_pnt)
255{
256 std::size_t const n_segs(ply->getNumberOfSegments());
257 // Neighbouring segments always intersects at a common vertex. The algorithm
258 // checks for intersections of non-neighbouring segments.
259 for (seg_it0 = ply->begin(); seg_it0 != ply->end() - 2; ++seg_it0)
260 {
261 seg_it1 = seg_it0 + 2;
262 std::size_t const seg_num_0 = seg_it0.getSegmentNumber();
263 for (; seg_it1 != ply->end(); ++seg_it1)
264 {
265 // Do not check first and last segment, because they are
266 // neighboured.
267 if (!(seg_num_0 == 0 && seg_it1.getSegmentNumber() == n_segs - 1))
268 {
269 if (lineSegmentIntersect(*seg_it0, *seg_it1, intersection_pnt))
270 {
271 return true;
272 }
273 }
274 }
275 }
276 return false;
277}
278
279void rotatePoints(Eigen::Matrix3d const& rot_mat,
280 std::vector<GeoLib::Point*>& pnts)
281{
282 rotatePoints(rot_mat, pnts.begin(), pnts.end());
283}
284
285Eigen::Matrix3d computeRotationMatrixToXY(Eigen::Vector3d const& n)
286{
287 Eigen::Matrix3d rot_mat = Eigen::Matrix3d::Zero();
288 // check if normal points already in the right direction
289 if (n[0] == 0 && n[1] == 0)
290 {
291 rot_mat(1, 1) = 1.0;
292
293 if (n[2] > 0)
294 {
295 // identity matrix
296 rot_mat(0, 0) = 1.0;
297 rot_mat(2, 2) = 1.0;
298 }
299 else
300 {
301 // rotate by pi about the y-axis
302 rot_mat(0, 0) = -1.0;
303 rot_mat(2, 2) = -1.0;
304 }
305
306 return rot_mat;
307 }
308
309 // sqrt (n_1^2 + n_3^2)
310 double const h0(std::sqrt(n[0] * n[0] + n[2] * n[2]));
311
312 // In case the x and z components of the normal are both zero the rotation
313 // to the x-z-plane is not required, i.e. only the rotation in the z-axis is
314 // required. The angle is either pi/2 or 3/2*pi. Thus the components of
315 // rot_mat are as follows.
316 if (h0 < std::numeric_limits<double>::epsilon())
317 {
318 rot_mat(0, 0) = 1.0;
319 if (n[1] > 0)
320 {
321 rot_mat(1, 2) = -1.0;
322 rot_mat(2, 1) = 1.0;
323 }
324 else
325 {
326 rot_mat(1, 2) = 1.0;
327 rot_mat(2, 1) = -1.0;
328 }
329 return rot_mat;
330 }
331
332 double const h1(1 / n.norm());
333
334 // general case: calculate entries of rotation matrix
335 rot_mat(0, 0) = n[2] / h0;
336 rot_mat(0, 1) = 0;
337 rot_mat(0, 2) = -n[0] / h0;
338 rot_mat(1, 0) = -n[1] * n[0] / h0 * h1;
339 rot_mat(1, 1) = h0 * h1;
340 rot_mat(1, 2) = -n[1] * n[2] / h0 * h1;
341 rot_mat(2, 0) = n[0] * h1;
342 rot_mat(2, 1) = n[1] * h1;
343 rot_mat(2, 2) = n[2] * h1;
344
345 return rot_mat;
346}
347
348Eigen::Matrix3d rotatePointsToXY(std::vector<GeoLib::Point*>& pnts)
349{
350 return rotatePointsToXY(pnts.begin(), pnts.end(), pnts.begin(), pnts.end());
351}
352
353std::unique_ptr<GeoLib::Point> triangleLineIntersection(
354 MathLib::Point3d const& a, MathLib::Point3d const& b,
355 MathLib::Point3d const& c, MathLib::Point3d const& p,
356 MathLib::Point3d const& q)
357{
358 Eigen::Vector3d const pq = q.asEigenVector3d() - p.asEigenVector3d();
359 Eigen::Vector3d const pa = a.asEigenVector3d() - p.asEigenVector3d();
360 Eigen::Vector3d const pb = b.asEigenVector3d() - p.asEigenVector3d();
361 Eigen::Vector3d const pc = c.asEigenVector3d() - p.asEigenVector3d();
362
363 double u = pq.cross(pc).dot(pb);
364 if (u < 0)
365 {
366 return nullptr;
367 }
368 double v = pq.cross(pa).dot(pc);
369 if (v < 0)
370 {
371 return nullptr;
372 }
373 double w = pq.cross(pb).dot(pa);
374 if (w < 0)
375 {
376 return nullptr;
377 }
378
379 const double denom(1.0 / (u + v + w));
380 u *= denom;
381 v *= denom;
382 w *= denom;
383 return std::make_unique<GeoLib::Point>(u * a[0] + v * b[0] + w * c[0],
384 u * a[1] + v * b[1] + w * c[1],
385 u * a[2] + v * b[2] + w * c[2]);
386}
387
389 std::vector<GeoLib::Polyline*>& plys)
390{
391 auto computeSegmentIntersections =
392 [&pnt_vec](GeoLib::Polyline& poly0, GeoLib::Polyline& poly1)
393 {
394 for (auto seg0_it(poly0.begin()); seg0_it != poly0.end(); ++seg0_it)
395 {
396 for (auto seg1_it(poly1.begin()); seg1_it != poly1.end(); ++seg1_it)
397 {
398 GeoLib::Point s(0.0, 0.0, 0.0, pnt_vec.size());
399 if (lineSegmentIntersect(*seg0_it, *seg1_it, s))
400 {
401 std::size_t const id(
402 pnt_vec.push_back(new GeoLib::Point(s)));
403 poly0.insertPoint(seg0_it.getSegmentNumber() + 1, id);
404 poly1.insertPoint(seg1_it.getSegmentNumber() + 1, id);
405 }
406 }
407 }
408 };
409
410 for (auto it0(plys.begin()); it0 != plys.end(); ++it0)
411 {
412 auto it1(it0);
413 ++it1;
414 for (; it1 != plys.end(); ++it1)
415 {
416 computeSegmentIntersections(*(*it0), *(*it1));
417 }
418 }
419}
420
421std::tuple<std::vector<GeoLib::Point*>, Eigen::Vector3d>
423{
424 // 1 copy all points
425 std::vector<GeoLib::Point*> polygon_points;
426 polygon_points.reserve(polygon_in.getNumberOfPoints());
427 for (std::size_t k(0); k < polygon_in.getNumberOfPoints(); k++)
428 {
429 polygon_points.push_back(new GeoLib::Point(*(polygon_in.getPoint(k))));
430 }
431
432 // 2 rotate points
433 auto [plane_normal, d_polygon] = GeoLib::getNewellPlane(polygon_points);
434 Eigen::Matrix3d const rot_mat =
436 GeoLib::rotatePoints(rot_mat, polygon_points);
437
438 // 3 set z coord to zero
439 std::for_each(polygon_points.begin(), polygon_points.end(),
440 [](GeoLib::Point* p) { (*p)[2] = 0.0; });
441
442 return {polygon_points, plane_normal};
443}
444
445std::vector<MathLib::Point3d> lineSegmentIntersect2d(
446 GeoLib::LineSegment const& ab, GeoLib::LineSegment const& cd)
447{
448 GeoLib::Point const& a{ab.getBeginPoint()};
449 GeoLib::Point const& b{ab.getEndPoint()};
450 GeoLib::Point const& c{cd.getBeginPoint()};
451 GeoLib::Point const& d{cd.getEndPoint()};
452
453 double const orient_abc(getOrientation(a, b, c));
454 double const orient_abd(getOrientation(a, b, d));
455
456 // check if the segment (cd) lies on the left or on the right of (ab)
457 if ((orient_abc > 0 && orient_abd > 0) ||
458 (orient_abc < 0 && orient_abd < 0))
459 {
460 return std::vector<MathLib::Point3d>();
461 }
462
463 // check: (cd) and (ab) are on the same line
464 if (orient_abc == 0.0 && orient_abd == 0.0)
465 {
466 double const eps(std::numeric_limits<double>::epsilon());
467 if (MathLib::sqrDist2d(a, c) < eps && MathLib::sqrDist2d(b, d) < eps)
468 {
469 return {{a, b}};
470 }
471 if (MathLib::sqrDist2d(a, d) < eps && MathLib::sqrDist2d(b, c) < eps)
472 {
473 return {{a, b}};
474 }
475
476 // Since orient_ab and orient_abd vanish, a, b, c, d are on the same
477 // line and for this reason it is enough to check the x-component.
478 auto isPointOnSegment = [](double q, double p0, double p1)
479 {
480 double const t((q - p0) / (p1 - p0));
481 return 0 <= t && t <= 1;
482 };
483
484 // check if c in (ab)
485 if (isPointOnSegment(c[0], a[0], b[0]))
486 {
487 // check if a in (cd)
488 if (isPointOnSegment(a[0], c[0], d[0]))
489 {
490 return {{a, c}};
491 }
492 // check b == c
493 if (MathLib::sqrDist2d(b, c) < eps)
494 {
495 return {{b}};
496 }
497 // check if b in (cd)
498 if (isPointOnSegment(b[0], c[0], d[0]))
499 {
500 return {{b, c}};
501 }
502 // check d in (ab)
503 if (isPointOnSegment(d[0], a[0], b[0]))
504 {
505 return {{c, d}};
506 }
507 std::stringstream err;
508 err.precision(std::numeric_limits<double>::max_digits10);
509 err << ab << " x " << cd;
510 OGS_FATAL(
511 "The case of parallel line segments ({:s}) is not handled yet. "
512 "Aborting.",
513 err.str());
514 }
515
516 // check if d in (ab)
517 if (isPointOnSegment(d[0], a[0], b[0]))
518 {
519 // check if a in (cd)
520 if (isPointOnSegment(a[0], c[0], d[0]))
521 {
522 return {{a, d}};
523 }
524 // check if b==d
525 if (MathLib::sqrDist2d(b, d) < eps)
526 {
527 return {{b}};
528 }
529 // check if b in (cd)
530 if (isPointOnSegment(b[0], c[0], d[0]))
531 {
532 return {{b, d}};
533 }
534 // d in (ab), b not in (cd): check c in (ab)
535 if (isPointOnSegment(c[0], a[0], b[0]))
536 {
537 return {{c, d}};
538 }
539
540 std::stringstream err;
541 err.precision(std::numeric_limits<double>::max_digits10);
542 err << ab << " x " << cd;
543 OGS_FATAL(
544 "The case of parallel line segments ({:s}) is not handled yet. "
545 "Aborting.",
546 err.str());
547 }
548 return std::vector<MathLib::Point3d>();
549 }
550
551 // precondition: points a, b, c are collinear
552 // the function checks if the point c is onto the line segment (a,b)
553 auto isCollinearPointOntoLineSegment = [](MathLib::Point3d const& a,
554 MathLib::Point3d const& b,
555 MathLib::Point3d const& c)
556 {
557 if (b[0] - a[0] != 0)
558 {
559 double const t = (c[0] - a[0]) / (b[0] - a[0]);
560 return 0.0 <= t && t <= 1.0;
561 }
562 if (b[1] - a[1] != 0)
563 {
564 double const t = (c[1] - a[1]) / (b[1] - a[1]);
565 return 0.0 <= t && t <= 1.0;
566 }
567 if (b[2] - a[2] != 0)
568 {
569 double const t = (c[2] - a[2]) / (b[2] - a[2]);
570 return 0.0 <= t && t <= 1.0;
571 }
572 return false;
573 };
574
575 if (orient_abc == 0.0)
576 {
577 if (isCollinearPointOntoLineSegment(a, b, c))
578 {
579 return {{c}};
580 }
581 return std::vector<MathLib::Point3d>();
582 }
583
584 if (orient_abd == 0.0)
585 {
586 if (isCollinearPointOntoLineSegment(a, b, d))
587 {
588 return {{d}};
589 }
590 return std::vector<MathLib::Point3d>();
591 }
592
593 // check if the segment (ab) lies on the left or on the right of (cd)
594 double const orient_cda(getOrientation(c, d, a));
595 double const orient_cdb(getOrientation(c, d, b));
596 if ((orient_cda > 0 && orient_cdb > 0) ||
597 (orient_cda < 0 && orient_cdb < 0))
598 {
599 return std::vector<MathLib::Point3d>();
600 }
601
602 // at this point it is sure that there is an intersection and the system of
603 // linear equations will be invertible
604 // solve the two linear equations (b-a, c-d) (t, s)^T = (c-a) simultaneously
605 Eigen::Matrix2d mat;
606 mat(0, 0) = b[0] - a[0];
607 mat(0, 1) = c[0] - d[0];
608 mat(1, 0) = b[1] - a[1];
609 mat(1, 1) = c[1] - d[1];
610 Eigen::Vector2d rhs{c[0] - a[0], c[1] - a[1]};
611
612 rhs = mat.partialPivLu().solve(rhs);
613 if (0 <= rhs[1] && rhs[1] <= 1.0)
614 {
615 return {MathLib::Point3d{std::array<double, 3>{
616 {c[0] + rhs[1] * (d[0] - c[0]), c[1] + rhs[1] * (d[1] - c[1]),
617 c[2] + rhs[1] * (d[2] - c[2])}}}};
618 }
619 return std::vector<MathLib::Point3d>(); // parameter s not in the valid
620 // range
621}
622
623void sortSegments(MathLib::Point3d const& seg_beg_pnt,
624 std::vector<GeoLib::LineSegment>& sub_segments)
625{
626 double const eps(std::numeric_limits<double>::epsilon());
627
628 auto findNextSegment =
629 [&eps](MathLib::Point3d const& seg_beg_pnt,
630 std::vector<GeoLib::LineSegment>& sub_segments,
631 std::vector<GeoLib::LineSegment>::iterator& sub_seg_it)
632 {
633 if (sub_seg_it == sub_segments.end())
634 {
635 return;
636 }
637 // find appropriate segment for the given segment begin point
638 auto act_beg_seg_it = std::find_if(
639 sub_seg_it, sub_segments.end(),
640 [&seg_beg_pnt, &eps](GeoLib::LineSegment const& seg)
641 {
642 return MathLib::sqrDist(seg_beg_pnt, seg.getBeginPoint()) <
643 eps ||
644 MathLib::sqrDist(seg_beg_pnt, seg.getEndPoint()) < eps;
645 });
646 if (act_beg_seg_it == sub_segments.end())
647 {
648 return;
649 }
650 // if necessary correct orientation of segment, i.e. swap beg and
651 // end
652 if (MathLib::sqrDist(seg_beg_pnt, act_beg_seg_it->getEndPoint()) <
653 MathLib::sqrDist(seg_beg_pnt, act_beg_seg_it->getBeginPoint()))
654 {
655 std::swap(act_beg_seg_it->getBeginPoint(),
656 act_beg_seg_it->getEndPoint());
657 }
658 assert(sub_seg_it != sub_segments.end());
659 // exchange segments within the container
660 if (sub_seg_it != act_beg_seg_it)
661 {
662 std::swap(*sub_seg_it, *act_beg_seg_it);
663 }
664 };
665
666 // find start segment
667 auto seg_it = sub_segments.begin();
668 findNextSegment(seg_beg_pnt, sub_segments, seg_it);
669
670 while (seg_it != sub_segments.end())
671 {
672 MathLib::Point3d& new_seg_beg_pnt(seg_it->getEndPoint());
673 seg_it++;
674 if (seg_it != sub_segments.end())
675 {
676 findNextSegment(new_seg_beg_pnt, sub_segments, seg_it);
677 }
678 }
679}
680
681Eigen::Matrix3d compute2DRotationMatrixToX(Eigen::Vector3d const& v)
682{
683 Eigen::Matrix3d rot_mat = Eigen::Matrix3d::Zero();
684 const double cos_theta = v[0];
685 const double sin_theta = v[1];
686 rot_mat(0, 0) = rot_mat(1, 1) = cos_theta;
687 rot_mat(0, 1) = sin_theta;
688 rot_mat(1, 0) = -sin_theta;
689 rot_mat(2, 2) = 1.0;
690 return rot_mat;
691}
692
693Eigen::Matrix3d compute3DRotationMatrixToX(Eigen::Vector3d const& v)
694{
695 // a vector on the plane
696 Eigen::Vector3d yy = Eigen::Vector3d::Zero();
697 auto const eps = std::numeric_limits<double>::epsilon();
698 if (std::abs(v[0]) > 0.0 && std::abs(v[1]) + std::abs(v[2]) < eps)
699 {
700 yy[2] = 1.0;
701 }
702 else if (std::abs(v[1]) > 0.0 && std::abs(v[0]) + std::abs(v[2]) < eps)
703 {
704 yy[0] = 1.0;
705 }
706 else if (std::abs(v[2]) > 0.0 && std::abs(v[0]) + std::abs(v[1]) < eps)
707 {
708 yy[1] = 1.0;
709 }
710 else
711 {
712 for (unsigned i = 0; i < 3; i++)
713 {
714 if (std::abs(v[i]) > 0.0)
715 {
716 yy[i] = -v[i];
717 break;
718 }
719 }
720 }
721 // z"_vec
722 Eigen::Vector3d const zz = v.cross(yy).normalized();
723 // y"_vec
724 yy = zz.cross(v).normalized();
725
726 Eigen::Matrix3d rot_mat;
727 rot_mat.row(0) = v;
728 rot_mat.row(1) = yy;
729 rot_mat.row(2) = zz;
730 return rot_mat;
731}
732
733} // end namespace GeoLib
#define OGS_FATAL(...)
Definition Error.h:19
GeoLib::Point const & getBeginPoint() const
GeoLib::Point const & getEndPoint() const
This class manages pointers to Points in a std::vector along with a name. It also handles the deletio...
Definition PointVec.h:25
std::size_t push_back(Point *pnt)
Definition PointVec.cpp:124
std::size_t getSegmentNumber() const
Definition Polyline.cpp:376
Class Polyline consists mainly of a reference to a point vector and a vector that stores the indices ...
Definition Polyline.h:29
std::size_t getNumberOfSegments() const
Definition Polyline.cpp:103
std::size_t getNumberOfPoints() const
Definition Polyline.cpp:98
const Point * getPoint(std::size_t i) const
returns the i-th point contained in the polyline
Definition Polyline.cpp:168
virtual bool insertPoint(std::size_t pos, std::size_t pnt_id)
Definition Polyline.cpp:44
SegmentIterator begin() const
Definition Polyline.h:177
SegmentIterator end() const
Definition Polyline.h:179
std::size_t size() const
Definition TemplateVec.h:88
Eigen::Vector3d const & asEigenVector3d() const
Definition Point3d.h:55
const double * data() const
Definition Point3d.h:51
Orientation getOrientationFast(MathLib::Point3d const &p0, MathLib::Point3d const &p1, MathLib::Point3d const &p2)
Eigen::Matrix3d compute3DRotationMatrixToX(Eigen::Vector3d const &v)
bool lineSegmentsIntersect(const GeoLib::Polyline *ply, GeoLib::Polyline::SegmentIterator &seg_it0, GeoLib::Polyline::SegmentIterator &seg_it1, GeoLib::Point &intersection_pnt)
void computeAndInsertAllIntersectionPoints(GeoLib::PointVec &pnt_vec, std::vector< GeoLib::Polyline * > &plys)
void rotatePoints(Eigen::Matrix3d const &rot_mat, InputIterator pnts_begin, InputIterator pnts_end)
bool parallel(Eigen::Vector3d v, Eigen::Vector3d w)
std::tuple< std::vector< GeoLib::Point * >, Eigen::Vector3d > rotatePolygonPointsToXY(GeoLib::Polygon const &polygon_in)
Eigen::Matrix3d computeRotationMatrixToXY(Eigen::Vector3d const &n)
Eigen::Matrix3d rotatePointsToXY(InputIterator1 p_pnts_begin, InputIterator1 p_pnts_end, InputIterator2 r_pnts_begin, InputIterator2 r_pnts_end)
void sortSegments(MathLib::Point3d const &seg_beg_pnt, std::vector< GeoLib::LineSegment > &sub_segments)
std::vector< MathLib::Point3d > lineSegmentIntersect2d(GeoLib::LineSegment const &ab, GeoLib::LineSegment const &cd)
std::unique_ptr< GeoLib::Point > triangleLineIntersection(MathLib::Point3d const &a, MathLib::Point3d const &b, MathLib::Point3d const &c, MathLib::Point3d const &p, MathLib::Point3d const &q)
std::pair< Eigen::Vector3d, double > getNewellPlane(InputIterator pnts_begin, InputIterator pnts_end)
Eigen::Matrix3d compute2DRotationMatrixToX(Eigen::Vector3d const &v)
bool lineSegmentIntersect(GeoLib::LineSegment const &s0, GeoLib::LineSegment const &s1, GeoLib::Point &s)
Orientation getOrientation(MathLib::Point3d const &p0, MathLib::Point3d const &p1, MathLib::Point3d const &p2)
double sqrDist2d(MathLib::Point3d const &p0, MathLib::Point3d const &p1)
Definition Point3d.h:114
double sqrDist(MathLib::Point3d const &p0, MathLib::Point3d const &p1)
Definition Point3d.cpp:19
double getOrientation2dFast(MathLib::Point3d const &a, MathLib::Point3d const &b, MathLib::Point3d const &c)
double getOrientation2d(MathLib::Point3d const &a, MathLib::Point3d const &b, MathLib::Point3d const &c)
double orient2dfast(double *, double *, double *)
double orient2d(double *, double *, double *)