17#include <Eigen/Geometry>
27extern double orient2d(
double*,
double*,
double*);
36 const_cast<double*
>(b.
data()),
37 const_cast<double*
>(c.data()));
45 const_cast<double*
>(b.
data()),
46 const_cast<double*
>(c.data()));
72 double const orientation =
85bool parallel(Eigen::Vector3d v, Eigen::Vector3d w)
87 const double eps(std::numeric_limits<double>::epsilon());
88 double const eps_squared = eps * eps;
91 if (v.squaredNorm() < eps_squared)
96 if (w.squaredNorm() < eps_squared)
105 if (std::abs(v[0] - w[0]) > eps)
109 if (std::abs(v[1] - w[1]) > eps)
113 if (std::abs(v[2] - w[2]) > eps)
124 if (std::abs(v[0] - w[0]) > eps)
128 if (std::abs(v[1] - w[1]) > eps)
132 if (std::abs(v[2] - w[2]) > eps)
150 if (!isCoplanar(pa, pb, pc, pd))
160 Eigen::Vector3d
const v = b - a;
161 Eigen::Vector3d
const w = d - c;
162 Eigen::Vector3d
const qp = c - a;
163 Eigen::Vector3d
const pq = a - c;
165 double const eps = std::numeric_limits<double>::epsilon();
166 double const squared_eps = eps * eps;
168 if (qp.squaredNorm() < squared_eps || (d - a).squaredNorm() < squared_eps)
173 if ((c - b).squaredNorm() < squared_eps ||
174 (d - b).squaredNorm() < squared_eps)
180 auto isLineSegmentIntersectingAB =
181 [&v](Eigen::Vector3d
const& ap, std::size_t i)
184 return 0.0 <= ap[i] / v[i] && ap[i] / v[i] <= 1.0;
203 std::size_t i_max(std::abs(v[0]) <= std::abs(v[1]) ? 1 : 0);
204 i_max = std::abs(v[i_max]) <= std::abs(v[2]) ? 2 : i_max;
205 if (isLineSegmentIntersectingAB(qp, i_max))
210 Eigen::Vector3d
const ad = d - a;
211 if (isLineSegmentIntersectingAB(ad, i_max))
222 const double sqr_len_v(v.squaredNorm());
223 const double sqr_len_w(w.squaredNorm());
226 mat(0, 0) = sqr_len_v;
227 mat(0, 1) = -v.dot(w);
228 mat(1, 1) = sqr_len_w;
229 mat(1, 0) = mat(0, 1);
231 Eigen::Vector2d rhs{v.dot(qp), w.dot(pq)};
233 rhs = mat.partialPivLu().solve(rhs);
237 const double l(-1.0 * std::numeric_limits<float>::epsilon());
239 const double u(1.0 + std::numeric_limits<float>::epsilon());
240 if (rhs[0] < l || u < rhs[0] || rhs[1] < l || u < rhs[1])
246 GeoLib::Point const p0(a[0] + rhs[0] * v[0], a[1] + rhs[0] * v[1],
247 a[2] + rhs[0] * v[2]);
248 GeoLib::Point const p1(c[0] + rhs[1] * w[0], c[1] + rhs[1] * w[1],
249 c[2] + rhs[1] * w[2]);
252 double const min_seg_len(
253 std::min(std::sqrt(sqr_len_v), std::sqrt(sqr_len_w)));
254 if (min_dist < min_seg_len * 1e-6)
256 s[0] = 0.5 * (p0[0] + p1[0]);
257 s[1] = 0.5 * (p0[1] + p1[1]);
258 s[2] = 0.5 * (p0[2] + p1[2]);
273 for (seg_it0 = ply->
begin(); seg_it0 != ply->
end() - 2; ++seg_it0)
275 seg_it1 = seg_it0 + 2;
277 for (; seg_it1 != ply->
end(); ++seg_it1)
294 std::vector<GeoLib::Point*>& pnts)
301 Eigen::Matrix3d rot_mat = Eigen::Matrix3d::Zero();
303 if (n[0] == 0 && n[1] == 0)
316 rot_mat(0, 0) = -1.0;
317 rot_mat(2, 2) = -1.0;
324 double const h0(std::sqrt(n[0] * n[0] + n[2] * n[2]));
330 if (h0 < std::numeric_limits<double>::epsilon())
335 rot_mat(1, 2) = -1.0;
341 rot_mat(2, 1) = -1.0;
346 double const h1(1 / n.norm());
349 rot_mat(0, 0) = n[2] / h0;
351 rot_mat(0, 2) = -n[0] / h0;
352 rot_mat(1, 0) = -n[1] * n[0] / h0 * h1;
353 rot_mat(1, 1) = h0 * h1;
354 rot_mat(1, 2) = -n[1] * n[2] / h0 * h1;
355 rot_mat(2, 0) = n[0] * h1;
356 rot_mat(2, 1) = n[1] * h1;
357 rot_mat(2, 2) = n[2] * h1;
364 return rotatePointsToXY(pnts.begin(), pnts.end(), pnts.begin(), pnts.end());
375 Eigen::Vector3d
const pc = c.asEigenVector3d() - p.asEigenVector3d();
377 double u = pq.cross(pc).dot(pb);
382 double v = pq.cross(pa).dot(pc);
387 double w = pq.cross(pb).dot(pa);
393 const double denom(1.0 / (u + v + w));
397 return std::make_unique<GeoLib::Point>(u * a[0] + v * b[0] + w * c[0],
398 u * a[1] + v * b[1] + w * c[1],
399 u * a[2] + v * b[2] + w * c[2]);
403 std::vector<GeoLib::Polyline*>& plys)
405 auto computeSegmentIntersections =
408 for (
auto seg0_it(poly0.
begin()); seg0_it != poly0.
end(); ++seg0_it)
410 for (
auto seg1_it(poly1.begin()); seg1_it != poly1.end(); ++seg1_it)
415 std::size_t
const id(
417 poly0.
insertPoint(seg0_it.getSegmentNumber() + 1,
id);
418 poly1.insertPoint(seg1_it.getSegmentNumber() + 1,
id);
424 for (
auto it0(plys.begin()); it0 != plys.end(); ++it0)
428 for (; it1 != plys.end(); ++it1)
430 computeSegmentIntersections(*(*it0), *(*it1));
435std::tuple<std::vector<GeoLib::Point*>, Eigen::Vector3d>
439 std::vector<GeoLib::Point*> polygon_points;
448 Eigen::Matrix3d
const rot_mat =
453 std::for_each(polygon_points.begin(), polygon_points.end(),
456 return {polygon_points, plane_normal};
471 if ((orient_abc > 0 && orient_abd > 0) ||
472 (orient_abc < 0 && orient_abd < 0))
474 return std::vector<MathLib::Point3d>();
478 if (orient_abc == 0.0 && orient_abd == 0.0)
480 double const eps(std::numeric_limits<double>::epsilon());
492 auto isPointOnSegment = [](
double q,
double p0,
double p1)
494 double const t((q - p0) / (p1 - p0));
495 return 0 <= t && t <= 1;
499 if (isPointOnSegment(c[0], a[0], b[0]))
502 if (isPointOnSegment(a[0], c[0], d[0]))
512 if (isPointOnSegment(b[0], c[0], d[0]))
517 if (isPointOnSegment(d[0], a[0], b[0]))
521 std::stringstream err;
522 err.precision(std::numeric_limits<double>::max_digits10);
523 err << ab <<
" x " << cd;
525 "The case of parallel line segments ({:s}) is not handled yet. "
531 if (isPointOnSegment(d[0], a[0], b[0]))
534 if (isPointOnSegment(a[0], c[0], d[0]))
544 if (isPointOnSegment(b[0], c[0], d[0]))
549 if (isPointOnSegment(c[0], a[0], b[0]))
554 std::stringstream err;
555 err.precision(std::numeric_limits<double>::max_digits10);
556 err << ab <<
" x " << cd;
558 "The case of parallel line segments ({:s}) is not handled yet. "
562 return std::vector<MathLib::Point3d>();
571 if (b[0] - a[0] != 0)
573 double const t = (c[0] - a[0]) / (b[0] - a[0]);
574 return 0.0 <= t && t <= 1.0;
576 if (b[1] - a[1] != 0)
578 double const t = (c[1] - a[1]) / (b[1] - a[1]);
579 return 0.0 <= t && t <= 1.0;
581 if (b[2] - a[2] != 0)
583 double const t = (c[2] - a[2]) / (b[2] - a[2]);
584 return 0.0 <= t && t <= 1.0;
589 if (orient_abc == 0.0)
591 if (isCollinearPointOntoLineSegment(a, b, c))
595 return std::vector<MathLib::Point3d>();
598 if (orient_abd == 0.0)
600 if (isCollinearPointOntoLineSegment(a, b, d))
604 return std::vector<MathLib::Point3d>();
610 if ((orient_cda > 0 && orient_cdb > 0) ||
611 (orient_cda < 0 && orient_cdb < 0))
613 return std::vector<MathLib::Point3d>();
620 mat(0, 0) = b[0] - a[0];
621 mat(0, 1) = c[0] - d[0];
622 mat(1, 0) = b[1] - a[1];
623 mat(1, 1) = c[1] - d[1];
624 Eigen::Vector2d rhs{c[0] - a[0], c[1] - a[1]};
626 rhs = mat.partialPivLu().solve(rhs);
627 if (0 <= rhs[1] && rhs[1] <= 1.0)
630 {c[0] + rhs[1] * (d[0] - c[0]), c[1] + rhs[1] * (d[1] - c[1]),
631 c[2] + rhs[1] * (d[2] - c[2])}}}};
633 return std::vector<MathLib::Point3d>();
638 std::vector<GeoLib::LineSegment>& sub_segments)
640 double const eps(std::numeric_limits<double>::epsilon());
642 auto findNextSegment =
644 std::vector<GeoLib::LineSegment>& sub_segments,
645 std::vector<GeoLib::LineSegment>::iterator& sub_seg_it)
647 if (sub_seg_it == sub_segments.end())
652 auto act_beg_seg_it = std::find_if(
653 sub_seg_it, sub_segments.end(),
656 return MathLib::sqrDist(seg_beg_pnt, seg.getBeginPoint()) <
658 MathLib::sqrDist(seg_beg_pnt, seg.getEndPoint()) < eps;
660 if (act_beg_seg_it == sub_segments.end())
669 std::swap(act_beg_seg_it->getBeginPoint(),
670 act_beg_seg_it->getEndPoint());
672 assert(sub_seg_it != sub_segments.end());
674 if (sub_seg_it != act_beg_seg_it)
676 std::swap(*sub_seg_it, *act_beg_seg_it);
681 auto seg_it = sub_segments.begin();
682 findNextSegment(seg_beg_pnt, sub_segments, seg_it);
684 while (seg_it != sub_segments.end())
688 if (seg_it != sub_segments.end())
690 findNextSegment(new_seg_beg_pnt, sub_segments, seg_it);
697 Eigen::Matrix3d rot_mat = Eigen::Matrix3d::Zero();
698 const double cos_theta = v[0];
699 const double sin_theta = v[1];
700 rot_mat(0, 0) = rot_mat(1, 1) = cos_theta;
701 rot_mat(0, 1) = sin_theta;
702 rot_mat(1, 0) = -sin_theta;
710 Eigen::Vector3d yy = Eigen::Vector3d::Zero();
711 auto const eps = std::numeric_limits<double>::epsilon();
712 if (std::abs(v[0]) > 0.0 && std::abs(v[1]) + std::abs(v[2]) < eps)
716 else if (std::abs(v[1]) > 0.0 && std::abs(v[0]) + std::abs(v[2]) < eps)
720 else if (std::abs(v[2]) > 0.0 && std::abs(v[0]) + std::abs(v[1]) < eps)
726 for (
unsigned i = 0; i < 3; i++)
728 if (std::abs(v[i]) > 0.0)
736 Eigen::Vector3d
const zz = v.cross(yy).normalized();
738 yy = zz.cross(v).normalized();
740 Eigen::Matrix3d rot_mat;
double orient2dfast(double *, double *, double *)
double orient2d(double *, double *, double *)
Definition of analytical geometry functions.
Definition of the PointVec class.
Definition of the PolyLine class.
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...
std::size_t push_back(Point *pnt)
std::size_t getSegmentNumber() const
Class Polyline consists mainly of a reference to a point vector and a vector that stores the indices ...
std::size_t getNumberOfSegments() const
std::size_t getNumberOfPoints() const
const Point * getPoint(std::size_t i) const
returns the i-th point contained in the polyline
virtual bool insertPoint(std::size_t pos, std::size_t pnt_id)
SegmentIterator begin() const
SegmentIterator end() const
Eigen::Vector3d const & asEigenVector3d() const
const double * data() const
double getOrientation2d(MathLib::Point3d const &a, MathLib::Point3d const &b, MathLib::Point3d const &c)
double getOrientation2dFast(MathLib::Point3d const &a, MathLib::Point3d const &b, MathLib::Point3d const &c)
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)
double sqrDist(MathLib::Point3d const &p0, MathLib::Point3d const &p1)