random_points_on_mesh.cpp 5.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 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 "random_points_on_mesh.h"
  9. #include "doublearea.h"
  10. #include "cumsum.h"
  11. #include "histc.h"
  12. #include <iostream>
  13. #include <cassert>
  14. template <typename DerivedV, typename DerivedF, typename DerivedB, typename DerivedFI>
  15. IGL_INLINE void igl::random_points_on_mesh(
  16. const int n,
  17. const Eigen::MatrixBase<DerivedV > & V,
  18. const Eigen::MatrixBase<DerivedF > & F,
  19. Eigen::PlainObjectBase<DerivedB > & B,
  20. Eigen::PlainObjectBase<DerivedFI > & FI)
  21. {
  22. using namespace Eigen;
  23. using namespace std;
  24. typedef typename DerivedV::Scalar Scalar;
  25. typedef Matrix<Scalar,Dynamic,1> VectorXs;
  26. VectorXs A;
  27. doublearea(V,F,A);
  28. // Should be traingle mesh. Although Turk's method 1 generalizes...
  29. assert(F.cols() == 3);
  30. VectorXs C;
  31. VectorXs A0(A.size()+1);
  32. A0(0) = 0;
  33. A0.bottomRightCorner(A.size(),1) = A;
  34. // Even faster would be to use the "Alias Table Method"
  35. cumsum(A0,1,C);
  36. const Scalar Cmax = C(C.size()-1);
  37. assert(Cmax > 0 && "Total surface area should be positive");
  38. // Why is this more accurate than `C /= C(C.size()-1)` ?
  39. for(int i = 0;i<C.size();i++) { C(i) = C(i)/Cmax; }
  40. const VectorXs R = (VectorXs::Random(n,1).array() + 1.)/2.;
  41. assert(R.minCoeff() >= 0);
  42. assert(R.maxCoeff() <= 1);
  43. histc(R,C,FI);
  44. const VectorXs S = (VectorXs::Random(n,1).array() + 1.)/2.;
  45. const VectorXs T = (VectorXs::Random(n,1).array() + 1.)/2.;
  46. B.resize(n,3);
  47. B.col(0) = 1.-T.array().sqrt();
  48. B.col(1) = (1.-S.array()) * T.array().sqrt();
  49. B.col(2) = S.array() * T.array().sqrt();
  50. }
  51. template <
  52. typename DerivedV,
  53. typename DerivedF,
  54. typename DerivedB,
  55. typename DerivedFI,
  56. typename DerivedX>
  57. IGL_INLINE void igl::random_points_on_mesh(
  58. const int n,
  59. const Eigen::MatrixBase<DerivedV > & V,
  60. const Eigen::MatrixBase<DerivedF > & F,
  61. Eigen::PlainObjectBase<DerivedB > & B,
  62. Eigen::PlainObjectBase<DerivedFI > & FI,
  63. Eigen::PlainObjectBase<DerivedX> & X)
  64. {
  65. random_points_on_mesh(n,V,F,B,FI);
  66. X = DerivedX::Zero(B.rows(),V.cols());
  67. for(int x = 0;x<B.rows();x++)
  68. {
  69. for(int b = 0;b<B.cols();b++)
  70. {
  71. X.row(x) += B(x,b)*V.row(F(FI(x),b));
  72. }
  73. }
  74. }
  75. template <typename DerivedV, typename DerivedF, typename ScalarB, typename DerivedFI>
  76. IGL_INLINE void igl::random_points_on_mesh(
  77. const int n,
  78. const Eigen::MatrixBase<DerivedV > & V,
  79. const Eigen::MatrixBase<DerivedF > & F,
  80. Eigen::SparseMatrix<ScalarB > & B,
  81. Eigen::PlainObjectBase<DerivedFI > & FI)
  82. {
  83. using namespace Eigen;
  84. using namespace std;
  85. Matrix<ScalarB,Dynamic,3> BC;
  86. random_points_on_mesh(n,V,F,BC,FI);
  87. vector<Triplet<ScalarB> > BIJV;
  88. BIJV.reserve(n*3);
  89. for(int s = 0;s<n;s++)
  90. {
  91. for(int c = 0;c<3;c++)
  92. {
  93. assert(FI(s) < F.rows());
  94. assert(FI(s) >= 0);
  95. const int v = F(FI(s),c);
  96. BIJV.push_back(Triplet<ScalarB>(s,v,BC(s,c)));
  97. }
  98. }
  99. B.resize(n,V.rows());
  100. B.reserve(n*3);
  101. B.setFromTriplets(BIJV.begin(),BIJV.end());
  102. }
  103. #ifdef IGL_STATIC_LIBRARY
  104. // Explicit template instantiation
  105. template void igl::random_points_on_mesh<Eigen::Matrix<float, -1, 3, 1, -1, 3>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<float, -1, 3, 1, -1, 3>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<float, -1, 3, 1, -1, 3> >(int, Eigen::MatrixBase<Eigen::Matrix<float, -1, 3, 1, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, 3, 1, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, 3, 1, -1, 3> >&);
  106. template void igl::random_points_on_mesh<Eigen::Matrix<float, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, float, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(int, Eigen::MatrixBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::SparseMatrix<float, 0, int>&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  107. template void igl::random_points_on_mesh<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<int, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(int, 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<int, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  108. template void igl::random_points_on_mesh<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, double, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(int, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::SparseMatrix<double, 0, int>&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  109. template void igl::random_points_on_mesh<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, double, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(int, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::SparseMatrix<double, 0, int>&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  110. #endif