main.cpp 4.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147
  1. #include <igl/read_triangle_mesh.h>
  2. #include <igl/parula.h>
  3. #include <igl/remove_unreferenced.h>
  4. #include <igl/opengl/glfw/Viewer.h>
  5. #include <igl/per_face_normals.h>
  6. #include <igl/orient_halfedges.h>
  7. #include <igl/cr_vector_laplacian.h>
  8. #include <igl/cr_vector_mass.h>
  9. #include <igl/crouzeix_raviart_cotmatrix.h>
  10. #include <igl/crouzeix_raviart_massmatrix.h>
  11. #include <igl/edge_midpoints.h>
  12. #include <igl/edge_vectors.h>
  13. #include <igl/average_from_edges_onto_vertices.h>
  14. #include <igl/min_quad_with_fixed.h>
  15. #include <igl/heat_geodesics.h>
  16. #include <Eigen/Core>
  17. #include <Eigen/SparseCholesky>
  18. #include <Eigen/Geometry>
  19. #include <iostream>
  20. #include <set>
  21. #include <limits>
  22. #include <stdlib.h>
  23. int main(int argc, char * argv[])
  24. {
  25. typedef Eigen::SparseMatrix<double> SparseMat;
  26. typedef Eigen::Matrix<double, 1, 1> Vector1d;
  27. typedef Eigen::Matrix<int, 1, 1> Vector1i;
  28. //Constants used for smoothing
  29. const double howMuchToSmoothBy = 1e-1;
  30. const int howManySmoothingInterations = 50;
  31. //Read our mesh
  32. Eigen::MatrixXd V;
  33. Eigen::MatrixXi F;
  34. if(!igl::read_triangle_mesh
  35. (argc>1?argv[1]: TUTORIAL_SHARED_PATH "/cheburashka.off",V,F)) {
  36. std::cout << "Failed to load mesh." << std::endl;
  37. }
  38. //Compute vector Laplacian and mass matrix
  39. Eigen::MatrixXi E, oE;//Compute Laplacian and mass matrix
  40. SparseMat vecL, vecM;
  41. igl::cr_vector_mass(V, F, E, vecM);
  42. igl::cr_vector_laplacian(V, F, E, oE, vecL);
  43. const int m = vecL.rows()/2; //The number of edges in the mesh
  44. //Convert the E / oE matrix format to list of edges / EMAP format required
  45. // by the functions constructing scalar Crouzeix-Raviart functions
  46. Eigen::MatrixXi Elist(m,2), EMAP(3*F.rows(),1);
  47. for(int i=0; i<F.rows(); ++i) {
  48. for(int j=0; j<3; ++j) {
  49. const int e = E(i,j);
  50. EMAP(i+j*F.rows()) = e;
  51. if(oE(i,j)>0) {
  52. Elist.row(e) << F(i, (j+1)%3), F(i, (j+2)%3);
  53. }
  54. }
  55. }
  56. SparseMat scalarL, scalarM;
  57. igl::crouzeix_raviart_massmatrix(V, F, Elist, EMAP, scalarM);
  58. igl::crouzeix_raviart_cotmatrix(V, F, Elist, EMAP, scalarL);
  59. //Compute edge midpoints & edge vectors
  60. Eigen::MatrixXd edgeMps, parVec, perpVec;
  61. igl::edge_midpoints(V, F, E, oE, edgeMps);
  62. igl::edge_vectors(V, F, E, oE, parVec, perpVec);
  63. //Perform the vector heat method
  64. const int initialIndex = 14319;
  65. const double initialPara=0.95, initialPerp=0.08;
  66. const double t = 0.01;
  67. SparseMat Aeq;
  68. Eigen::VectorXd Beq;
  69. Eigen::VectorXi known = Eigen::Vector2i(initialIndex, initialIndex+m);
  70. Eigen::VectorXd knownVals = Eigen::Vector2d(initialPara, initialPerp);
  71. Eigen::VectorXd Y0 = Eigen::VectorXd::Zero(2*m), Yt;
  72. Y0(initialIndex) = initialPara; Y0(initialIndex+m) = initialPerp;
  73. igl::min_quad_with_fixed
  74. (SparseMat(vecM+t*vecL), Eigen::VectorXd(-vecM*Y0), known, knownVals,
  75. Aeq, Beq, false, Yt);
  76. Eigen::VectorXd u0 = Eigen::VectorXd::Zero(m), ut;
  77. u0(initialIndex) = sqrt(initialPara*initialPara + initialPerp*initialPerp);
  78. Eigen::VectorXi knownScal = Vector1i(initialIndex);
  79. Eigen::VectorXd knownScalVals = Vector1d(u0(initialIndex));
  80. igl::min_quad_with_fixed
  81. (SparseMat(scalarM+t*scalarL), Eigen::VectorXd(-scalarM*u0), knownScal,
  82. knownScalVals, Aeq, Beq, false, ut);
  83. Eigen::VectorXd phi0 = Eigen::VectorXd::Zero(m), phit;
  84. phi0(initialIndex) = 1;
  85. Eigen::VectorXd knownScalValsPhi = Vector1d(1);
  86. igl::min_quad_with_fixed
  87. (SparseMat(scalarM+t*scalarL), Eigen::VectorXd(-scalarM*phi0), knownScal,
  88. knownScalValsPhi, Aeq, Beq, false, phit);
  89. Eigen::ArrayXd Xtfactor = ut.array() /
  90. (phit.array() * (Yt.array().segment(0,m)*Yt.array().segment(0,m)
  91. + Yt.array().segment(m,m)*Yt.array().segment(m,m)).sqrt());
  92. Eigen::VectorXd Xt(2*m);
  93. Xt.segment(0,m) = Xtfactor * Yt.segment(0,m).array();
  94. Xt.segment(m,m) = Xtfactor * Yt.segment(m,m).array();
  95. //Compute scalar heat colors
  96. igl::HeatGeodesicsData<double> hgData;
  97. igl::heat_geodesics_precompute(V, F, hgData);
  98. Eigen::VectorXd heatColor;
  99. Eigen::VectorXi gamma = Elist.row(initialIndex);
  100. igl::heat_geodesics_solve(hgData, gamma, heatColor);
  101. //Convert vector field for plotting
  102. Eigen::MatrixXd vecs(m, 3);
  103. for(int i=0; i<edgeMps.rows(); ++i) {
  104. vecs.row(i) = Xt(i)*parVec.row(i) + Xt(i+edgeMps.rows())*perpVec.row(i);
  105. }
  106. //Viewer that shows parallel transported vector
  107. igl::opengl::glfw::Viewer viewer;
  108. viewer.data().set_mesh(V,F);
  109. viewer.data().show_lines = false;
  110. viewer.data().set_data(heatColor.maxCoeff()-heatColor.array(), //invert colormap
  111. igl::COLOR_MAP_TYPE_VIRIDIS);
  112. const double s = 0.012; //How much to scale vectors during plotting
  113. Eigen::MatrixXd vecColors(m, 3);
  114. for(int i=0; i<m; ++i) {
  115. vecColors.row(i) << 0.1, 0.1, 0.1;
  116. }
  117. vecColors.row(initialIndex) << 0.9, 0.1, 0.1;
  118. viewer.data().add_edges(edgeMps, edgeMps + s*vecs, vecColors);
  119. std::cout << R"(The red vector is parallel transported to every point on the surface.
  120. The surface is shaded by geodesic distance from the red vector.
  121. )"
  122. << std::endl;
  123. viewer.launch();
  124. return 0;
  125. }