main.cpp 3.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111
  1. #include <igl/read_triangle_mesh.h>
  2. #include <igl/hessian_energy.h>
  3. #include <igl/massmatrix.h>
  4. #include <igl/cotmatrix.h>
  5. #include <igl/jet.h>
  6. #include <igl/edges.h>
  7. #include <igl/components.h>
  8. #include <igl/remove_unreferenced.h>
  9. #include <igl/opengl/glfw/Viewer.h>
  10. #include <Eigen/Core>
  11. #include <Eigen/SparseCholesky>
  12. #include <iostream>
  13. #include <set>
  14. #include <limits>
  15. #include <stdlib.h>
  16. #include "tutorial_shared_path.h"
  17. #include <igl/isolines.h>
  18. int main(int argc, char * argv[])
  19. {
  20. typedef Eigen::SparseMatrix<double> SparseMat;
  21. //Read our mesh
  22. Eigen::MatrixXd V;
  23. Eigen::MatrixXi F, E;
  24. if(!igl::read_triangle_mesh(
  25. argc>1?argv[1]: TUTORIAL_SHARED_PATH "/beetle.off",V,F)) {
  26. std::cout << "Failed to load mesh." << std::endl;
  27. }
  28. igl::edges(F,E);
  29. //Constructing an exact function to smooth
  30. Eigen::VectorXd zexact = V.block(0,2,V.rows(),1).array()
  31. + 0.5*V.block(0,1,V.rows(),1).array()
  32. + V.block(0,1,V.rows(),1).array().pow(2)
  33. + V.block(0,2,V.rows(),1).array().pow(3);
  34. //Make the exact function noisy
  35. srand(5);
  36. const double s = 0.2*(zexact.maxCoeff() - zexact.minCoeff());
  37. Eigen::VectorXd znoisy = zexact + s*Eigen::VectorXd::Random(zexact.size());
  38. //Constructing the squared Laplacian and squared Hessian energy
  39. SparseMat L, M;
  40. igl::cotmatrix(V, F, L);
  41. igl::massmatrix(V, F, igl::MASSMATRIX_TYPE_BARYCENTRIC, M);
  42. Eigen::SimplicialLDLT<SparseMat> solver(M);
  43. SparseMat MinvL = solver.solve(L);
  44. SparseMat QL = L.transpose()*MinvL;
  45. SparseMat QH;
  46. igl::hessian_energy(V, F, QH);
  47. //Solve to find Laplacian-smoothed and Hessian-smoothed solutions
  48. const double al = 8e-4;
  49. Eigen::SimplicialLDLT<SparseMat> lapSolver(al*QL + (1.-al)*M);
  50. Eigen::VectorXd zl = lapSolver.solve(al*M*znoisy);
  51. const double ah = 5e-6;
  52. Eigen::SimplicialLDLT<SparseMat> hessSolver(ah*QH + (1.-ah)*M);
  53. Eigen::VectorXd zh = hessSolver.solve(ah*M*znoisy);
  54. //Viewer that shows all functions: zexact, znoisy, zl, zh
  55. igl::opengl::glfw::Viewer viewer;
  56. viewer.data().set_mesh(V,F);
  57. viewer.data().show_lines = false;
  58. viewer.callback_key_down =
  59. [&](igl::opengl::glfw::Viewer & viewer, unsigned char key, int mod)->bool
  60. {
  61. //Graduate result to show isolines, then compute color matrix
  62. const Eigen::VectorXd* z;
  63. switch(key) {
  64. case '1':
  65. z = &zexact;
  66. break;
  67. case '2':
  68. z = &znoisy;
  69. break;
  70. case '3':
  71. z = &zl;
  72. break;
  73. case '4':
  74. z = &zh;
  75. break;
  76. default:
  77. return false;
  78. }
  79. Eigen::MatrixXd isoV;
  80. Eigen::MatrixXi isoE;
  81. if(key!='2')
  82. igl::isolines(V, F, *z, 30, isoV, isoE);
  83. viewer.data().set_edges(isoV,isoE,Eigen::RowVector3d(0,0,0));
  84. Eigen::MatrixXd colors;
  85. igl::jet(*z, true, colors);
  86. viewer.data().set_colors(colors);
  87. return true;
  88. };
  89. std::cout << R"(Usage:
  90. 1 Show original function
  91. 2 Show noisy function
  92. 3 Biharmonic smoothing (zero Neumann boundary)
  93. 4 Biharmonic smoothing (natural Hessian boundary)
  94. )";
  95. viewer.launch();
  96. return 0;
  97. }