ExactPredicate.h 3.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106
  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. template<typename Scalar>
  19. class ExactPredicate {
  20. public:
  21. typedef CGAL::Exact_predicates_exact_constructions_kernel Epeck;
  22. typedef CGAL::Exact_predicates_inexact_constructions_kernel Epick;
  23. typedef typename std::conditional<std::is_same<Scalar, Epeck::FT>::value,
  24. Epeck, Epick>::type Kernel;
  25. static short orientation(const Scalar pa[2], const Scalar pb[2], const Scalar pc[2])
  26. {
  27. switch(CGAL::orientation(
  28. typename Kernel::Point_2(pa[0], pa[1]),
  29. typename Kernel::Point_2(pb[0], pb[1]),
  30. typename Kernel::Point_2(pc[0], pc[1]))) {
  31. case CGAL::LEFT_TURN:
  32. return 1;
  33. case CGAL::RIGHT_TURN:
  34. return -1;
  35. case CGAL::COLLINEAR:
  36. return 0;
  37. default:
  38. throw "Invalid orientation";
  39. }
  40. }
  41. static short orientation(const Scalar pa[3], const Scalar pb[3], const Scalar pc[3], const Scalar pd[3])
  42. {
  43. switch(CGAL::orientation(
  44. typename Kernel::Point_3(pa[0], pa[1], pa[2]),
  45. typename Kernel::Point_3(pb[0], pb[1], pb[2]),
  46. typename Kernel::Point_3(pc[0], pc[1], pc[2]),
  47. typename Kernel::Point_3(pd[0], pd[1], pd[2]))) {
  48. case CGAL::POSITIVE:
  49. return 1;
  50. case CGAL::NEGATIVE:
  51. return -1;
  52. case CGAL::COPLANAR:
  53. return 0;
  54. default:
  55. throw "Invalid orientation";
  56. }
  57. }
  58. static short incircle(const Scalar pa[2], const Scalar pb[2], const Scalar pc[2], const Scalar pd[2])
  59. {
  60. switch(CGAL::side_of_oriented_circle(
  61. typename Kernel::Point_2(pa[0], pa[1]),
  62. typename Kernel::Point_2(pb[0], pb[1]),
  63. typename Kernel::Point_2(pc[0], pc[1]),
  64. typename Kernel::Point_2(pd[0], pd[1]))) {
  65. case CGAL::ON_POSITIVE_SIDE:
  66. return 1;
  67. case CGAL::ON_NEGATIVE_SIDE:
  68. return -1;
  69. case CGAL::ON_ORIENTED_BOUNDARY:
  70. return 0;
  71. default:
  72. throw "Invalid incircle result";
  73. }
  74. }
  75. static short insphere(const Scalar pa[3], const Scalar pb[3], const Scalar pc[3], const Scalar pd[3],
  76. const Scalar pe[3])
  77. {
  78. switch(CGAL::side_of_oriented_sphere(
  79. typename Kernel::Point_3(pa[0], pa[1], pa[2]),
  80. typename Kernel::Point_3(pb[0], pb[1], pb[2]),
  81. typename Kernel::Point_3(pc[0], pc[1], pc[2]),
  82. typename Kernel::Point_3(pd[0], pd[1], pd[2]),
  83. typename Kernel::Point_3(pe[0], pe[1], pe[2]))) {
  84. case CGAL::ON_POSITIVE_SIDE:
  85. return 1;
  86. case CGAL::ON_NEGATIVE_SIDE:
  87. return -1;
  88. case CGAL::ON_ORIENTED_BOUNDARY:
  89. return 0;
  90. default:
  91. throw "Invalid incircle result";
  92. }
  93. }
  94. };
  95. }
  96. }
  97. }
  98. #endif