bijective_composite_harmonic_mapping.cpp 4.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2017 Alec Jacobson <alecjacobson@gmail.com>
  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 "bijective_composite_harmonic_mapping.h"
  9. #include "slice.h"
  10. #include "doublearea.h"
  11. #include "harmonic.h"
  12. //#include "matlab/MatlabWorkspace.h"
  13. #include <iostream>
  14. template <
  15. typename DerivedV,
  16. typename DerivedF,
  17. typename Derivedb,
  18. typename Derivedbc,
  19. typename DerivedU>
  20. IGL_INLINE bool igl::bijective_composite_harmonic_mapping(
  21. const Eigen::MatrixBase<DerivedV> & V,
  22. const Eigen::MatrixBase<DerivedF> & F,
  23. const Eigen::MatrixBase<Derivedb> & b,
  24. const Eigen::MatrixBase<Derivedbc> & bc,
  25. Eigen::PlainObjectBase<DerivedU> & U)
  26. {
  27. return bijective_composite_harmonic_mapping(V,F,b,bc,1,200,20,true,U);
  28. }
  29. template <
  30. typename DerivedV,
  31. typename DerivedF,
  32. typename Derivedb,
  33. typename Derivedbc,
  34. typename DerivedU>
  35. IGL_INLINE bool igl::bijective_composite_harmonic_mapping(
  36. const Eigen::MatrixBase<DerivedV> & V,
  37. const Eigen::MatrixBase<DerivedF> & F,
  38. const Eigen::MatrixBase<Derivedb> & b,
  39. const Eigen::MatrixBase<Derivedbc> & bc,
  40. const int min_steps,
  41. const int max_steps,
  42. const int num_inner_iters,
  43. const bool test_for_flips,
  44. Eigen::PlainObjectBase<DerivedU> & U)
  45. {
  46. typedef typename Derivedbc::Scalar Scalar;
  47. assert(V.cols() == 2 && bc.cols() == 2 && "Input should be 2D");
  48. assert(F.cols() == 3 && "F should contain triangles");
  49. int tries = 0;
  50. int nsteps = min_steps;
  51. Derivedbc bc0;
  52. slice(V,b,1,bc0);
  53. // It's difficult to check for flips "robustly" in the sense that the input
  54. // mesh might not have positive/consistent sign to begin with.
  55. while(nsteps<=max_steps)
  56. {
  57. U = V;
  58. int flipped = 0;
  59. int nans = 0;
  60. int step = 0;
  61. for(;step<=nsteps;step++)
  62. {
  63. const Scalar t = ((Scalar)step)/((Scalar)nsteps);
  64. // linearly interpolate boundary conditions
  65. // TODO: replace this with something that guarantees a homotopic "morph"
  66. // of the boundary conditions. Something like "Homotopic Morphing of
  67. // Planar Curves" [Dym et al. 2015] but also handling multiple connected
  68. // components.
  69. Derivedbc bct = bc0 + t*(bc - bc0);
  70. // Compute dsicrete harmonic map using metric of previous step
  71. for(int iter = 0;iter<num_inner_iters;iter++)
  72. {
  73. //std::cout<<nsteps<<" t: "<<t<<" iter: "<<iter;
  74. //igl::matlab::MatlabWorkspace mw;
  75. //mw.save(U,"U");
  76. //mw.save_index(F,"F");
  77. //mw.save_index(b,"b");
  78. //mw.save(bct,"bct");
  79. //mw.write("numerical.mat");
  80. harmonic(DerivedU(U),F,b,bct,1,U);
  81. igl::slice(U,b,1,bct);
  82. nans = (U.array() != U.array()).count();
  83. if(test_for_flips)
  84. {
  85. Eigen::Matrix<Scalar,Eigen::Dynamic,1> A;
  86. doublearea(U,F,A);
  87. flipped = (A.array() < 0 ).count();
  88. //std::cout<<" "<<flipped<<" nan? "<<(U.array() != U.array()).any()<<std::endl;
  89. if(flipped == 0 && nans == 0) break;
  90. }
  91. }
  92. if(flipped > 0 || nans>0) break;
  93. }
  94. if(flipped == 0 && nans == 0)
  95. {
  96. return step == nsteps+1;
  97. }
  98. nsteps *= 2;
  99. }
  100. //std::cout<<"failed to finish in "<<nsteps<<"..."<<std::endl;
  101. return false;
  102. }
  103. #ifdef IGL_STATIC_LIBRARY
  104. // Explicit template instantiation
  105. // generated by autoexplicit.sh
  106. template bool igl::bijective_composite_harmonic_mapping<Eigen::Matrix<double, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -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, 1, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -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&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> >&);
  107. // generated by autoexplicit.sh
  108. template bool igl::bijective_composite_harmonic_mapping<Eigen::Matrix<double, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -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, 1, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -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, int, int, bool, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> >&);
  109. #endif