bijective_composite_harmonic_mapping.cpp 4.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2017 Alec Jacobson <[email protected]>
  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 nsteps = min_steps;
  50. Eigen::Matrix<typename Derivedbc::Scalar, Eigen::Dynamic, Eigen::Dynamic> bc0;
  51. slice(V,b.col(0),1,bc0);
  52. // It's difficult to check for flips "robustly" in the sense that the input
  53. // mesh might not have positive/consistent sign to begin with.
  54. while(nsteps<=max_steps)
  55. {
  56. U = V;
  57. int flipped = 0;
  58. int nans = 0;
  59. int step = 0;
  60. for(;step<=nsteps;step++)
  61. {
  62. const Scalar t = ((Scalar)step)/((Scalar)nsteps);
  63. // linearly interpolate boundary conditions
  64. // TODO: replace this with something that guarantees a homotopic "morph"
  65. // of the boundary conditions. Something like "Homotopic Morphing of
  66. // Planar Curves" [Dym et al. 2015] but also handling multiple connected
  67. // components.
  68. Eigen::Matrix<typename Derivedbc::Scalar, Eigen::Dynamic, Eigen::Dynamic> bct = bc0 + t * (bc - bc0);
  69. // Compute dsicrete harmonic map using metric of previous step
  70. for(int iter = 0;iter<num_inner_iters;iter++)
  71. {
  72. //std::cout<<nsteps<<" t: "<<t<<" iter: "<<iter;
  73. //igl::matlab::MatlabWorkspace mw;
  74. //mw.save(U,"U");
  75. //mw.save_index(F,"F");
  76. //mw.save_index(b,"b");
  77. //mw.save(bct,"bct");
  78. //mw.write("numerical.mat");
  79. harmonic(Eigen::Matrix<typename DerivedU::Scalar, Eigen::Dynamic, Eigen::Dynamic>(U), F, b, bct, 1, U);
  80. igl::slice(U,b.col(0),1,bct);
  81. nans = (U.array() != U.array()).count();
  82. if(test_for_flips)
  83. {
  84. Eigen::Matrix<Scalar,Eigen::Dynamic,1> A;
  85. doublearea(U,F,A);
  86. flipped = (A.array() < 0 ).count();
  87. //std::cout<<" "<<flipped<<" nan? "<<(U.array() != U.array()).any()<<std::endl;
  88. if(flipped == 0 && nans == 0) break;
  89. }
  90. }
  91. if(flipped > 0 || nans>0) break;
  92. }
  93. if(flipped == 0 && nans == 0)
  94. {
  95. return step == nsteps+1;
  96. }
  97. nsteps *= 2;
  98. }
  99. //std::cout<<"failed to finish in "<<nsteps<<"..."<<std::endl;
  100. return false;
  101. }
  102. #ifdef IGL_STATIC_LIBRARY
  103. // Explicit template instantiation
  104. // generated by autoexplicit.sh
  105. 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> >&);
  106. // generated by autoexplicit.sh
  107. 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> >&);
  108. #endif