volume.cpp 5.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 Alec Jacobson <alecjacobson@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. #include "volume.h"
  9. #include "cross.h"
  10. #include <Eigen/Geometry>
  11. template <
  12. typename DerivedV,
  13. typename DerivedT,
  14. typename Derivedvol>
  15. IGL_INLINE void igl::volume(
  16. const Eigen::MatrixBase<DerivedV>& V,
  17. const Eigen::MatrixBase<DerivedT>& T,
  18. Eigen::PlainObjectBase<Derivedvol>& vol)
  19. {
  20. using namespace Eigen;
  21. const int m = T.rows();
  22. vol.resize(m,1);
  23. for(int t = 0;t<m;t++)
  24. {
  25. typedef Eigen::Matrix<typename DerivedV::Scalar,1,3> RowVector3S;
  26. const RowVector3S & a = V.row(T(t,0));
  27. const RowVector3S & b = V.row(T(t,1));
  28. const RowVector3S & c = V.row(T(t,2));
  29. const RowVector3S & d = V.row(T(t,3));
  30. vol(t) = -(a-d).dot((b-d).cross(c-d))/6.;
  31. }
  32. }
  33. template <
  34. typename DerivedA,
  35. typename DerivedB,
  36. typename DerivedC,
  37. typename DerivedD,
  38. typename Derivedvol>
  39. IGL_INLINE void igl::volume(
  40. const Eigen::MatrixBase<DerivedA> & A,
  41. const Eigen::MatrixBase<DerivedB> & B,
  42. const Eigen::MatrixBase<DerivedC> & C,
  43. const Eigen::MatrixBase<DerivedD> & D,
  44. Eigen::PlainObjectBase<Derivedvol> & vol)
  45. {
  46. const auto & AmD = A-D;
  47. const auto & BmD = B-D;
  48. const auto & CmD = C-D;
  49. DerivedA BmDxCmD;
  50. cross(BmD.eval(),CmD.eval(),BmDxCmD);
  51. const auto & AmDdx = (AmD.array() * BmDxCmD.array()).rowwise().sum();
  52. vol = -AmDdx/6.;
  53. }
  54. template <
  55. typename VecA,
  56. typename VecB,
  57. typename VecC,
  58. typename VecD>
  59. IGL_INLINE typename VecA::Scalar igl::volume_single(
  60. const VecA & a,
  61. const VecB & b,
  62. const VecC & c,
  63. const VecD & d)
  64. {
  65. return -(a-d).dot((b-d).cross(c-d))/6.;
  66. }
  67. template <
  68. typename DerivedL,
  69. typename Derivedvol>
  70. IGL_INLINE void igl::volume(
  71. const Eigen::MatrixBase<DerivedL>& L,
  72. Eigen::PlainObjectBase<Derivedvol>& vol)
  73. {
  74. using namespace Eigen;
  75. const int m = L.rows();
  76. typedef typename Derivedvol::Scalar ScalarS;
  77. vol.resize(m,1);
  78. for(int t = 0;t<m;t++)
  79. {
  80. const ScalarS u = L(t,0);
  81. const ScalarS v = L(t,1);
  82. const ScalarS w = L(t,2);
  83. const ScalarS U = L(t,3);
  84. const ScalarS V = L(t,4);
  85. const ScalarS W = L(t,5);
  86. const ScalarS X = (w - U + v)*(U + v + w);
  87. const ScalarS x = (U - v + w)*(v - w + U);
  88. const ScalarS Y = (u - V + w)*(V + w + u);
  89. const ScalarS y = (V - w + u)*(w - u + V);
  90. const ScalarS Z = (v - W + u)*(W + u + v);
  91. const ScalarS z = (W - u + v)*(u - v + W);
  92. const ScalarS a = sqrt(x*Y*Z);
  93. const ScalarS b = sqrt(y*Z*X);
  94. const ScalarS c = sqrt(z*X*Y);
  95. const ScalarS d = sqrt(x*y*z);
  96. vol(t) = sqrt(
  97. (-a + b + c + d)*
  98. ( a - b + c + d)*
  99. ( a + b - c + d)*
  100. ( a + b + c - d))/
  101. (192.*u*v*w);
  102. }
  103. }
  104. #ifdef IGL_STATIC_LIBRARY
  105. // Explicit template instantiation
  106. // Nonsense template
  107. namespace igl{ template<> void volume<Eigen::Matrix<double, -1, 2, 0, -1, 2>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, 2, 0, -1, 2> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&){} }
  108. // generated by autoexplicit.sh
  109. template void igl::volume<Eigen::Matrix<float, -1, 3, 1, -1, 3>, Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<float, -1, 3, 1, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  110. // generated by autoexplicit.sh
  111. template void igl::volume<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::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<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  112. // generated by autoexplicit.sh
  113. template void igl::volume<Eigen::Matrix<double, -1, 6, 0, -1, 6>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, 6, 0, -1, 6> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  114. // generated by autoexplicit.sh
  115. template void igl::volume<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> >&);
  116. template Eigen::Matrix<double, 1, 3, 1, 1, 3>::Scalar igl::volume_single<Eigen::Matrix<double, 1, 3, 1, 1, 3>, Eigen::Matrix<double, 1, 3, 1, 1, 3>, Eigen::Matrix<double, 1, 3, 1, 1, 3>, Eigen::Matrix<double, 1, 3, 1, 1, 3> >(Eigen::Matrix<double, 1, 3, 1, 1, 3> const&, Eigen::Matrix<double, 1, 3, 1, 1, 3> const&, Eigen::Matrix<double, 1, 3, 1, 1, 3> const&, Eigen::Matrix<double, 1, 3, 1, 1, 3> const&);
  117. template void igl::volume<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::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> > &);
  118. #endif