main.cpp 4.3 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/edge_midpoints.h>
  10. #include <igl/edge_vectors.h>
  11. #include <igl/average_from_edges_onto_vertices.h>
  12. #include <igl/PI.h>
  13. #include <Eigen/Core>
  14. #include <Eigen/SparseCholesky>
  15. #include <Eigen/Geometry>
  16. #include <iostream>
  17. #include <set>
  18. #include <limits>
  19. #include <stdlib.h>
  20. #include "tutorial_shared_path.h"
  21. int main(int argc, char * argv[])
  22. {
  23. typedef Eigen::SparseMatrix<double> SparseMat;
  24. //Constants used for smoothing
  25. const double howMuchToSmoothBy = 1e-1;
  26. const int howManySmoothingInterations = 50;
  27. //Read our mesh
  28. Eigen::MatrixXd V;
  29. Eigen::MatrixXi F;
  30. if(!igl::read_triangle_mesh
  31. (argc>1?argv[1]: TUTORIAL_SHARED_PATH "/elephant.obj",V,F)) {
  32. std::cout << "Failed to load mesh." << std::endl;
  33. }
  34. //Orient edges for plotting
  35. Eigen::MatrixXi E, oE;
  36. igl::orient_halfedges(F, E, oE);
  37. //Compute edge midpoints & edge vectors
  38. Eigen::MatrixXd edgeMps, parVec, perpVec;
  39. igl::edge_midpoints(V, F, E, oE, edgeMps);
  40. igl::edge_vectors(V, F, E, oE, parVec, perpVec);
  41. //Constructing a function to add noise to
  42. const auto zraw_function = [] (const Eigen::Vector3d& x) {
  43. return Eigen::Vector3d(0.2*x(1) + cos(2*x(1)+0.2),
  44. 0.5*x(0) + 0.15,
  45. 0.3*cos(0.2+igl::PI*x(2)));
  46. };
  47. Eigen::VectorXd zraw(2*edgeMps.rows());
  48. for(int i=0; i<edgeMps.rows(); ++i) {
  49. const Eigen::Vector3d f = zraw_function(edgeMps.row(i));
  50. zraw(i) = f.dot(parVec.row(i));
  51. zraw(i+edgeMps.rows()) = f.dot(perpVec.row(i));
  52. }
  53. //Add noise
  54. srand(71);
  55. const double l = 15;
  56. Eigen::VectorXd znoisy = zraw + l*Eigen::VectorXd::Random(zraw.size());
  57. //Denoise function using the vector Dirichlet energy
  58. Eigen::VectorXd zsmoothed = znoisy;
  59. for(int i=0; i<howManySmoothingInterations; ++i) {
  60. //Compute Laplacian and mass matrix
  61. SparseMat L, M;
  62. igl::cr_vector_mass(V, F, E, oE, M);
  63. igl::cr_vector_laplacian(V, F, E, oE, L);
  64. //Implicit step
  65. Eigen::SimplicialLDLT<SparseMat> rhsSolver(M + howMuchToSmoothBy*L);
  66. zsmoothed = rhsSolver.solve(M*zsmoothed);
  67. }
  68. //Convert vector fields for plotting
  69. const auto cr_result_to_vecs_and_colors = [&]
  70. (const Eigen::VectorXd& z, Eigen::MatrixXd& vecs, Eigen::MatrixXd& colors) {
  71. vecs.resize(edgeMps.rows(), 3);
  72. for(int i=0; i<edgeMps.rows(); ++i) {
  73. vecs.row(i) = z(i)*parVec.row(i)
  74. + z(i+edgeMps.rows())*perpVec.row(i);
  75. }
  76. igl::average_from_edges_onto_vertices
  77. (F, E, oE, vecs.rowwise().norm(), colors);
  78. };
  79. Eigen::MatrixXd noisyvecs, noisycolors, smoothedvecs, smoothedcolors,
  80. rawvecs, rawcolors;
  81. cr_result_to_vecs_and_colors(znoisy, noisyvecs, noisycolors);
  82. cr_result_to_vecs_and_colors(zsmoothed, smoothedvecs, smoothedcolors);
  83. cr_result_to_vecs_and_colors(zraw, rawvecs, rawcolors);
  84. //Viewer that shows noisy and denoised functions
  85. igl::opengl::glfw::Viewer viewer;
  86. viewer.data().set_mesh(V,F);
  87. viewer.data().show_lines = false;
  88. viewer.callback_key_down =
  89. [&](igl::opengl::glfw::Viewer & viewer, unsigned char key, int mod)->bool
  90. {
  91. const Eigen::MatrixXd *vecs, *colors;
  92. switch(key) {
  93. case '1':
  94. vecs = &rawvecs;
  95. colors = &rawcolors;
  96. break;
  97. case '2':
  98. vecs = &noisyvecs;
  99. colors = &noisycolors;
  100. break;
  101. case '3':
  102. vecs = &smoothedvecs;
  103. colors = &smoothedcolors;
  104. break;
  105. default:
  106. return false;
  107. }
  108. viewer.data().set_data(*colors);
  109. viewer.data().clear_edges();
  110. const double s = 0.08; //How much to scale vectors during plotting
  111. viewer.data().add_edges(edgeMps, edgeMps + s*(*vecs),
  112. Eigen::RowVector3d(0.1, 0.1, 0.1));
  113. return true;
  114. };
  115. std::cout << R"(Usage:
  116. 1 Show raw function
  117. 2 Show noisy function
  118. 3 Show smoothed function
  119. )";
  120. Eigen::MatrixXd CM;
  121. igl::parula(Eigen::VectorXd::LinSpaced
  122. (500,znoisy.minCoeff(), znoisy.maxCoeff()).eval(), true, CM);
  123. viewer.data().set_colormap(CM);
  124. viewer.callback_key_down(viewer, '1', 0);
  125. viewer.launch();
  126. return 0;
  127. }