ExactPredicate.h 5.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2016 Qingan Zhou <qnzhou@gmail.com>
  4. //
  5. // This Source Code Form is subject to the terms of the Mozilla Public License
  6. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  7. // obtain one at http://mozilla.org/MPL/2.0/.
  8. #ifndef IGL_COPYLEFT_CGAL_EXACT_PREDICATE
  9. #define IGL_COPYLEFT_CGAL_EXACT_PREDICATE
  10. #include <CGAL/Exact_predicates_exact_constructions_kernel.h>
  11. #include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
  12. namespace igl
  13. {
  14. namespace copyleft
  15. {
  16. namespace cgal
  17. {
  18. // This class extracts the exact predicates routines from CGAL kernels.
  19. //
  20. // Example:
  21. // Eigen::PlainObjectBase<DerivedV> V = ...;
  22. // typedef typename DerivedV::Scalar Scalar;
  23. // typedef igl::copyleft::cgal::ExactPredicate<Scalar> Predicate;
  24. // auto result = Predicate::orientation(
  25. // {V(0, 0), V(0, 1)},
  26. // {V(1, 0), V(1, 1)},
  27. // {V(2, 0), V(2, 1)});
  28. template<typename Scalar>
  29. class ExactPredicate {
  30. public:
  31. typedef CGAL::Exact_predicates_exact_constructions_kernel Epeck;
  32. typedef CGAL::Exact_predicates_inexact_constructions_kernel Epick;
  33. typedef typename std::conditional<std::is_same<Scalar, Epeck::FT>::value,
  34. Epeck, Epick>::type Kernel;
  35. // Inputs:
  36. // pa,pb,pc 2D points.
  37. // Output:
  38. // 1 if pa,pb,pc are counterclockwise oriented.
  39. // 0 if pa,pb,pc are counterclockwise oriented.
  40. // -1 if pa,pb,pc are clockwise oriented.
  41. static short orientation(const Scalar pa[2], const Scalar pb[2], const Scalar pc[2])
  42. {
  43. switch(CGAL::orientation(
  44. typename Kernel::Point_2(pa[0], pa[1]),
  45. typename Kernel::Point_2(pb[0], pb[1]),
  46. typename Kernel::Point_2(pc[0], pc[1]))) {
  47. case CGAL::LEFT_TURN:
  48. return 1;
  49. case CGAL::RIGHT_TURN:
  50. return -1;
  51. case CGAL::COLLINEAR:
  52. return 0;
  53. default:
  54. throw "Invalid orientation";
  55. }
  56. }
  57. // Inputs:
  58. // pa,pb,pc,pd 3D points.
  59. // Output:
  60. // 1 if pa,pb,pc,pd forms a tet of positive volume.
  61. // 0 if pa,pb,pc,pd are coplanar.
  62. // -1 if pa,pb,pc,pd forms a tet of negative volume.
  63. static short orientation(const Scalar pa[3], const Scalar pb[3], const Scalar pc[3], const Scalar pd[3])
  64. {
  65. switch(CGAL::orientation(
  66. typename Kernel::Point_3(pa[0], pa[1], pa[2]),
  67. typename Kernel::Point_3(pb[0], pb[1], pb[2]),
  68. typename Kernel::Point_3(pc[0], pc[1], pc[2]),
  69. typename Kernel::Point_3(pd[0], pd[1], pd[2]))) {
  70. case CGAL::POSITIVE:
  71. return 1;
  72. case CGAL::NEGATIVE:
  73. return -1;
  74. case CGAL::COPLANAR:
  75. return 0;
  76. default:
  77. throw "Invalid orientation";
  78. }
  79. }
  80. // Inputs:
  81. // pa,pb,pc,pd 2D points.
  82. // Output:
  83. // 1 if pd is inside of the oriented circle formed by pa,pb,pc.
  84. // 0 if pd is co-circular with pa,pb,pc.
  85. // -1 if pd is outside of the oriented circle formed by pa,pb,pc.
  86. static short incircle(const Scalar pa[2], const Scalar pb[2], const Scalar pc[2], const Scalar pd[2])
  87. {
  88. switch(CGAL::side_of_oriented_circle(
  89. typename Kernel::Point_2(pa[0], pa[1]),
  90. typename Kernel::Point_2(pb[0], pb[1]),
  91. typename Kernel::Point_2(pc[0], pc[1]),
  92. typename Kernel::Point_2(pd[0], pd[1]))) {
  93. case CGAL::ON_POSITIVE_SIDE:
  94. return 1;
  95. case CGAL::ON_NEGATIVE_SIDE:
  96. return -1;
  97. case CGAL::ON_ORIENTED_BOUNDARY:
  98. return 0;
  99. default:
  100. throw "Invalid incircle result";
  101. }
  102. }
  103. // Inputs:
  104. // pa,pb,pc,pd,pe 3D points.
  105. // Output:
  106. // 1 if pe is inside of the oriented sphere formed by pa,pb,pc,pd.
  107. // 0 if pe is co-spherical with pa,pb,pc,pd.
  108. // -1 if pe is outside of the oriented sphere formed by pa,pb,pc,pd.
  109. static short insphere(const Scalar pa[3], const Scalar pb[3], const Scalar pc[3], const Scalar pd[3],
  110. const Scalar pe[3])
  111. {
  112. switch(CGAL::side_of_oriented_sphere(
  113. typename Kernel::Point_3(pa[0], pa[1], pa[2]),
  114. typename Kernel::Point_3(pb[0], pb[1], pb[2]),
  115. typename Kernel::Point_3(pc[0], pc[1], pc[2]),
  116. typename Kernel::Point_3(pd[0], pd[1], pd[2]),
  117. typename Kernel::Point_3(pe[0], pe[1], pe[2]))) {
  118. case CGAL::ON_POSITIVE_SIDE:
  119. return 1;
  120. case CGAL::ON_NEGATIVE_SIDE:
  121. return -1;
  122. case CGAL::ON_ORIENTED_BOUNDARY:
  123. return 0;
  124. default:
  125. throw "Invalid incircle result";
  126. }
  127. }
  128. };
  129. }
  130. }
  131. }
  132. #endif