circumradius.cpp 7.8 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2016 Alec Jacobson <[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 "circumradius.h"
  9. #include "edge_lengths.h"
  10. #include "doublearea.h"
  11. #include "placeholders.h"
  12. #include <Eigen/QR>
  13. template <
  14. typename DerivedV,
  15. typename DerivedF,
  16. typename DerivedR>
  17. IGL_INLINE void igl::circumradius(
  18. const Eigen::MatrixBase<DerivedV> & V,
  19. const Eigen::MatrixBase<DerivedF> & F,
  20. Eigen::PlainObjectBase<DerivedR> & R)
  21. {
  22. // cols at compiletime should be Dynamic or 3
  23. static_assert(DerivedF::ColsAtCompileTime == Eigen::Dynamic || DerivedF::ColsAtCompileTime == 3,"F should contain triangles");
  24. assert(F.cols() == 3 && "F should contain triangles");
  25. Eigen::Matrix<typename DerivedV::Scalar,Eigen::Dynamic,3> l;
  26. igl::edge_lengths(V,F,l);
  27. DerivedR A;
  28. igl::doublearea(l,0.,A);
  29. // use formula: R=abc/(4*area) to compute the circum radius
  30. R = l.col(0).array() * l.col(1).array() * l.col(2).array() / (2.0*A.array());
  31. }
  32. template <
  33. typename DerivedV,
  34. typename DerivedT,
  35. typename DerivedR,
  36. typename DerivedC,
  37. typename DerivedB>
  38. IGL_INLINE void igl::circumradius(
  39. const Eigen::MatrixBase<DerivedV> & V,
  40. const Eigen::MatrixBase<DerivedT> & T,
  41. Eigen::PlainObjectBase<DerivedR> & R,
  42. Eigen::PlainObjectBase<DerivedC> & C,
  43. Eigen::PlainObjectBase<DerivedB> & B)
  44. {
  45. using Scalar = typename DerivedV::Scalar;
  46. const int TCOLS = DerivedT::ColsAtCompileTime;
  47. const int VCOLS = DerivedV::ColsAtCompileTime;
  48. const int ACOLS = TCOLS==Eigen::Dynamic?Eigen::Dynamic:TCOLS+1;
  49. const int ss = T.cols();
  50. R.resize(T.rows(),1);
  51. C.resize(T.rows(),V.cols());
  52. B.resize(T.rows(),T.cols());
  53. for(int i = 0;i<T.rows();i++)
  54. {
  55. Eigen::Matrix<Scalar,ACOLS,ACOLS> A(T.cols()+1,T.cols()+1);
  56. // Not sure if this .eval() is a good idea
  57. const auto Vi = V(T.row(i),igl::placeholders::all).eval();
  58. A.topLeftCorner(T.cols(),T.cols()) = 2*(Vi*Vi.transpose());
  59. A.block(0,T.cols(),T.cols(),1).setConstant(1);
  60. A.block(T.cols(),0,1,T.cols()).setConstant(1);
  61. A(T.cols(),T.cols()) = 0;
  62. Eigen::Matrix<Scalar,ACOLS,1> b(T.cols()+1,1);
  63. b.head(T.cols()) = (Vi.array().square().rowwise().sum()).eval();
  64. b(T.cols()) = 1;
  65. // [B;λ] = A\b
  66. const auto x = (A.colPivHouseholderQr().solve(b)).eval();
  67. B.row(i) = x.head(T.cols());
  68. C.row(i) = B.row(i) * Vi;
  69. const Scalar lambda = x(T.cols());
  70. R(i) = std::sqrt(lambda + C.row(i).squaredNorm());
  71. }
  72. }
  73. #ifdef IGL_STATIC_LIBRARY
  74. // Explicit template instantiation
  75. // generated by autoexplicit.sh
  76. template void igl::circumradius<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  77. template void igl::circumradius<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  78. template void igl::circumradius<Eigen::Matrix<double, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 4, 0, -1, 4> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 4, 0, -1, 4> >&);
  79. template void igl::circumradius<Eigen::Matrix<double, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 3, 0, -1, 3> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&);
  80. template void igl::circumradius<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 4, 0, -1, 4> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 4, 0, -1, 4> >&);
  81. template void igl::circumradius<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 3, 0, -1, 3> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&);
  82. template void igl::circumradius<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 4, 0, -1, 4> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 4, 0, -1, 4> >&);
  83. template void igl::circumradius<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 4, 0, -1, 4>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 4, 0, -1, 4> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 4, 0, -1, 4> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 4, 0, -1, 4> >&);
  84. template void igl::circumradius<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 3, 0, -1, 3> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&);
  85. template void igl::circumradius<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  86. #endif