isolines.cpp 4.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2017 Oded Stein <oded.stein@columbia.edu>
  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 "isolines.h"
  9. #include <vector>
  10. #include <array>
  11. #include <iostream>
  12. #include "remove_duplicate_vertices.h"
  13. template <typename DerivedV,
  14. typename DerivedF,
  15. typename DerivedZ,
  16. typename DerivedIsoV,
  17. typename DerivedIsoE>
  18. IGL_INLINE void igl::isolines(
  19. const Eigen::MatrixBase<DerivedV>& V,
  20. const Eigen::MatrixBase<DerivedF>& F,
  21. const Eigen::MatrixBase<DerivedZ>& z,
  22. const int n,
  23. Eigen::PlainObjectBase<DerivedIsoV>& isoV,
  24. Eigen::PlainObjectBase<DerivedIsoE>& isoE)
  25. {
  26. //Constants
  27. const int dim = V.cols();
  28. assert(dim==2 || dim==3);
  29. const int nVerts = V.rows();
  30. assert(z.rows() == nVerts &&
  31. "There must be as many function entries as vertices");
  32. const int nFaces = F.rows();
  33. const int np1 = n+1;
  34. const double min = z.minCoeff(), max = z.maxCoeff();
  35. //Following http://www.alecjacobson.com/weblog/?p=2529
  36. typedef typename DerivedZ::Scalar Scalar;
  37. typedef Eigen::Matrix<Scalar, Eigen::Dynamic, 1> Vec;
  38. Vec iso(np1);
  39. for(int i=0; i<np1; ++i)
  40. iso(i) = Scalar(i)/Scalar(n)*(max-min) + min;
  41. typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> Matrix;
  42. std::array<Matrix,3> t{{Matrix(nFaces, np1),
  43. Matrix(nFaces, np1), Matrix(nFaces, np1)}};
  44. for(int i=0; i<nFaces; ++i) {
  45. for(int k=0; k<3; ++k) {
  46. const Scalar z1=z(F(i,k)), z2=z(F(i,(k+1)%3));
  47. for(int j=0; j<np1; ++j) {
  48. t[k](i,j) = (iso(j)-z1) / (z2-z1);
  49. if(t[k](i,j)<0 || t[k](i,j)>1)
  50. t[k](i,j) = std::numeric_limits<Scalar>::quiet_NaN();
  51. }
  52. }
  53. }
  54. std::array<std::vector<int>,3> Fij, Iij;
  55. for(int i=0; i<nFaces; ++i) {
  56. for(int j=0; j<np1; ++j) {
  57. for(int k=0; k<3; ++k) {
  58. const int kp1=(k+1)%3, kp2=(k+2)%3;
  59. if(std::isfinite(t[kp1](i,j)) && std::isfinite(t[kp2](i,j))) {
  60. Fij[k].push_back(i);
  61. Iij[k].push_back(j);
  62. }
  63. }
  64. }
  65. }
  66. const int K = Fij[0].size()+Fij[1].size()+Fij[2].size();
  67. isoV.resize(2*K, dim);
  68. int b = 0;
  69. for(int k=0; k<3; ++k) {
  70. const int kp1=(k+1)%3, kp2=(k+2)%3;
  71. for(int i=0; i<Fij[k].size(); ++i) {
  72. isoV.row(b+i) = (1.-t[kp1](Fij[k][i],Iij[k][i]))*
  73. V.row(F(Fij[k][i],kp1)) +
  74. t[kp1](Fij[k][i],Iij[k][i])*V.row(F(Fij[k][i],kp2));
  75. isoV.row(K+b+i) = (1.-t[kp2](Fij[k][i],Iij[k][i]))*
  76. V.row(F(Fij[k][i],kp2)) +
  77. t[kp2](Fij[k][i],Iij[k][i])*V.row(F(Fij[k][i],k));
  78. }
  79. b += Fij[k].size();
  80. }
  81. isoE.resize(K,2);
  82. for(int i=0; i<K; ++i)
  83. isoE.row(i) << i, K+i;
  84. //Remove double entries
  85. typedef typename DerivedIsoV::Scalar LScalar;
  86. typedef typename DerivedIsoE::Scalar LInt;
  87. typedef Eigen::Matrix<LInt, Eigen::Dynamic, 1> LIVec;
  88. typedef Eigen::Matrix<LScalar, Eigen::Dynamic, Eigen::Dynamic> LMat;
  89. typedef Eigen::Matrix<LInt, Eigen::Dynamic, Eigen::Dynamic> LIMat;
  90. LIVec dummy1, dummy2;
  91. igl::remove_duplicate_vertices(LMat(isoV), LIMat(isoE),
  92. 2.2204e-15, isoV, dummy1, dummy2, isoE);
  93. }
  94. #ifdef IGL_STATIC_LIBRARY
  95. // Explicit template instantiation
  96. template void igl::isolines<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<int, -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::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, int const, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > &, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > &);
  97. #endif