predicates.cpp 6.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2019 Qingnan Zhou <[email protected]>
  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. #include "./predicates.h"
  9. // This is a different file also called predicates.h
  10. #include <predicates.h>
  11. #include <type_traits>
  12. namespace igl {
  13. namespace predicates {
  14. using REAL = IGL_PREDICATES_REAL;
  15. #ifdef LIBIGL_PREDICATES_USE_FLOAT
  16. #define IGL_PREDICATES_ASSERT_SCALAR(Vector) \
  17. static_assert( \
  18. std::is_same<typename Vector::Scalar, float>::value, \
  19. "Shewchuk's exact predicates only support float")
  20. #else
  21. #define IGL_PREDICATES_ASSERT_SCALAR(Vector) \
  22. static_assert( \
  23. std::is_same<typename Vector::Scalar, double>::value || \
  24. std::is_same<typename Vector::Scalar, float>::value, \
  25. "Shewchuk's exact predicates only support float and double")
  26. #endif
  27. IGL_INLINE void exactinit() {
  28. // Thread-safe initialization using Meyers' singleton
  29. class MySingleton {
  30. public:
  31. static MySingleton& instance() {
  32. static MySingleton instance;
  33. return instance;
  34. }
  35. private:
  36. MySingleton() { ::exactinit(); }
  37. };
  38. MySingleton::instance();
  39. }
  40. template<typename Vector2D>
  41. IGL_INLINE Orientation orient2d(
  42. const Eigen::MatrixBase<Vector2D>& pa,
  43. const Eigen::MatrixBase<Vector2D>& pb,
  44. const Eigen::MatrixBase<Vector2D>& pc)
  45. {
  46. EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Vector2D, 2);
  47. IGL_PREDICATES_ASSERT_SCALAR(Vector2D);
  48. using Point = Eigen::Matrix<REAL, 2, 1>;
  49. Point a{pa[0], pa[1]};
  50. Point b{pb[0], pb[1]};
  51. Point c{pc[0], pc[1]};
  52. const auto r = ::orient2d(a.data(), b.data(), c.data());
  53. if (r > 0) return Orientation::POSITIVE;
  54. else if (r < 0) return Orientation::NEGATIVE;
  55. else return Orientation::COLLINEAR;
  56. }
  57. template<typename Vector3D>
  58. IGL_INLINE Orientation orient3d(
  59. const Eigen::MatrixBase<Vector3D>& pa,
  60. const Eigen::MatrixBase<Vector3D>& pb,
  61. const Eigen::MatrixBase<Vector3D>& pc,
  62. const Eigen::MatrixBase<Vector3D>& pd)
  63. {
  64. EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Vector3D, 3);
  65. IGL_PREDICATES_ASSERT_SCALAR(Vector3D);
  66. using Point = Eigen::Matrix<REAL, 3, 1>;
  67. Point a{pa[0], pa[1], pa[2]};
  68. Point b{pb[0], pb[1], pb[2]};
  69. Point c{pc[0], pc[1], pc[2]};
  70. Point d{pd[0], pd[1], pd[2]};
  71. const auto r = ::orient3d(a.data(), b.data(), c.data(), d.data());
  72. if (r > 0) return Orientation::POSITIVE;
  73. else if (r < 0) return Orientation::NEGATIVE;
  74. else return Orientation::COPLANAR;
  75. }
  76. template<typename Vector2D>
  77. IGL_INLINE Orientation incircle(
  78. const Eigen::MatrixBase<Vector2D>& pa,
  79. const Eigen::MatrixBase<Vector2D>& pb,
  80. const Eigen::MatrixBase<Vector2D>& pc,
  81. const Eigen::MatrixBase<Vector2D>& pd)
  82. {
  83. EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Vector2D, 2);
  84. IGL_PREDICATES_ASSERT_SCALAR(Vector2D);
  85. using Point = Eigen::Matrix<REAL, 2, 1>;
  86. Point a{pa[0], pa[1]};
  87. Point b{pb[0], pb[1]};
  88. Point c{pc[0], pc[1]};
  89. Point d{pd[0], pd[1]};
  90. const auto r = ::incircle(a.data(), b.data(), c.data(), d.data());
  91. if (r > 0) return Orientation::INSIDE;
  92. else if (r < 0) return Orientation::OUTSIDE;
  93. else return Orientation::COCIRCULAR;
  94. }
  95. template<typename Vector3D>
  96. IGL_INLINE Orientation insphere(
  97. const Eigen::MatrixBase<Vector3D>& pa,
  98. const Eigen::MatrixBase<Vector3D>& pb,
  99. const Eigen::MatrixBase<Vector3D>& pc,
  100. const Eigen::MatrixBase<Vector3D>& pd,
  101. const Eigen::MatrixBase<Vector3D>& pe)
  102. {
  103. EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Vector3D, 3);
  104. IGL_PREDICATES_ASSERT_SCALAR(Vector3D);
  105. using Point = Eigen::Matrix<REAL, 3, 1>;
  106. Point a{pa[0], pa[1], pa[2]};
  107. Point b{pb[0], pb[1], pb[2]};
  108. Point c{pc[0], pc[1], pc[2]};
  109. Point d{pd[0], pd[1], pd[2]};
  110. Point e{pe[0], pe[1], pe[2]};
  111. const auto r = ::insphere(a.data(), b.data(), c.data(), d.data(), e.data());
  112. if (r > 0) return Orientation::INSIDE;
  113. else if (r < 0) return Orientation::OUTSIDE;
  114. else return Orientation::COSPHERICAL;
  115. }
  116. }
  117. }
  118. #undef IGL_PREDICATES_ASSERT_SCALAR
  119. #ifdef IGL_STATIC_LIBRARY
  120. // Explicit template instantiation
  121. #define IGL_ORIENT2D(Vector) template igl::predicates::Orientation igl::predicates::orient2d<Vector>(const Eigen::MatrixBase<Vector>&, const Eigen::MatrixBase<Vector>&, const Eigen::MatrixBase<Vector>&)
  122. #define IGL_INCIRCLE(Vector) template igl::predicates::Orientation igl::predicates::incircle<Vector>(const Eigen::MatrixBase<Vector>&, const Eigen::MatrixBase<Vector>&, const Eigen::MatrixBase<Vector>&, const Eigen::MatrixBase<Vector>&)
  123. #define IGL_ORIENT3D(Vector) template igl::predicates::Orientation igl::predicates::orient3d<Vector>(const Eigen::MatrixBase<Vector>&, const Eigen::MatrixBase<Vector>&, const Eigen::MatrixBase<Vector>&, const Eigen::MatrixBase<Vector>&)
  124. #define IGL_INSPHERE(Vector) template igl::predicates::Orientation igl::predicates::insphere<Vector>(const Eigen::MatrixBase<Vector>&, const Eigen::MatrixBase<Vector>&, const Eigen::MatrixBase<Vector>&, const Eigen::MatrixBase<Vector>&, const Eigen::MatrixBase<Vector>&)
  125. #define IGL_MATRIX(T, R, C) Eigen::Matrix<T, R, C>
  126. IGL_ORIENT2D(IGL_MATRIX(float, 1, 2));
  127. IGL_INCIRCLE(IGL_MATRIX(float, 1, 2));
  128. IGL_ORIENT2D(IGL_MATRIX(float, 2, 1));
  129. IGL_INCIRCLE(IGL_MATRIX(float, 2, 1));
  130. IGL_ORIENT3D(IGL_MATRIX(float, 1, 3));
  131. IGL_INSPHERE(IGL_MATRIX(float, 1, 3));
  132. IGL_ORIENT3D(IGL_MATRIX(float, 3, 1));
  133. IGL_INSPHERE(IGL_MATRIX(float, 3, 1));
  134. #ifndef LIBIGL_PREDICATES_USE_FLOAT
  135. IGL_ORIENT2D(IGL_MATRIX(double, 1, 2));
  136. IGL_INCIRCLE(IGL_MATRIX(double, 1, 2));
  137. IGL_ORIENT2D(IGL_MATRIX(double, 2, 1));
  138. IGL_INCIRCLE(IGL_MATRIX(double, 2, 1));
  139. IGL_ORIENT3D(IGL_MATRIX(double, 1, 3));
  140. IGL_INSPHERE(IGL_MATRIX(double, 1, 3));
  141. IGL_ORIENT3D(IGL_MATRIX(double, 3, 1));
  142. IGL_INSPHERE(IGL_MATRIX(double, 3, 1));
  143. #endif
  144. #undef IGL_MATRIX
  145. #undef IGL_ORIENT2D
  146. #undef IGL_ORIENT3D
  147. #undef IGL_INCIRCLE
  148. #undef IGL_INSPHERE
  149. #endif