lscm.cpp 4.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 Daniele Panozzo <[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 "lscm.h"
  9. #include "lscm_hessian.h"
  10. #include "massmatrix.h"
  11. #include "repdiag.h"
  12. #include "eigs.h"
  13. #include "min_quad_with_fixed.h"
  14. #include <iostream>
  15. template <
  16. typename DerivedV,
  17. typename DerivedF,
  18. typename Derivedb,
  19. typename Derivedbc,
  20. typename DerivedV_uv,
  21. typename QScalar>
  22. IGL_INLINE bool igl::lscm(
  23. const Eigen::MatrixBase<DerivedV> & V,
  24. const Eigen::MatrixBase<DerivedF> & F,
  25. const Eigen::MatrixBase<Derivedb> & b,
  26. const Eigen::MatrixBase<Derivedbc> & bc,
  27. Eigen::PlainObjectBase<DerivedV_uv> & V_uv,
  28. Eigen::SparseMatrix<QScalar> & Q)
  29. {
  30. igl::lscm_hessian(V,F,Q);
  31. Eigen::Matrix<typename Derivedb::Scalar,Eigen::Dynamic,1> b_flat(b.size()*bc.cols(),1);
  32. Eigen::Matrix<typename Derivedbc::Scalar,Eigen::Dynamic,1> bc_flat(bc.size(),1);
  33. for(int c = 0;c<bc.cols();c++)
  34. {
  35. b_flat.block(c*b.size(),0,b.rows(),1) = c*V.rows() + b.array();
  36. bc_flat.block(c*bc.rows(),0,bc.rows(),1) = bc.col(c);
  37. }
  38. const Eigen::VectorXd B_flat = Eigen::VectorXd::Zero(V.rows()*2);
  39. igl::min_quad_with_fixed_data<QScalar> data;
  40. if(!igl::min_quad_with_fixed_precompute(Q,b_flat,Eigen::SparseMatrix<QScalar>(),true,data))
  41. {
  42. return false;
  43. }
  44. Eigen::MatrixXd W_flat;
  45. if(!min_quad_with_fixed_solve(data,B_flat,bc_flat,Eigen::VectorXd(),W_flat))
  46. {
  47. return false;
  48. }
  49. assert(W_flat.rows() == V.rows()*2);
  50. V_uv.resize(V.rows(),2);
  51. for (unsigned i=0;i<V_uv.cols();++i)
  52. {
  53. V_uv.col(i) = W_flat.block(V_uv.rows()*i,0,V_uv.rows(),1);
  54. }
  55. return true;
  56. }
  57. template <
  58. typename DerivedV,
  59. typename DerivedF,
  60. typename Derivedb,
  61. typename Derivedbc,
  62. typename DerivedV_uv>
  63. IGL_INLINE bool igl::lscm(
  64. const Eigen::MatrixBase<DerivedV> & V,
  65. const Eigen::MatrixBase<DerivedF> & F,
  66. const Eigen::MatrixBase<Derivedb> & b,
  67. const Eigen::MatrixBase<Derivedbc> & bc,
  68. Eigen::PlainObjectBase<DerivedV_uv> & V_uv)
  69. {
  70. Eigen::SparseMatrix<typename DerivedV_uv::Scalar> Q;
  71. return lscm(V, F, b, bc, V_uv, Q);
  72. }
  73. template <
  74. typename DerivedV,
  75. typename DerivedF,
  76. typename DerivedV_uv>
  77. IGL_INLINE bool igl::lscm(
  78. const Eigen::MatrixBase<DerivedV> & V,
  79. const Eigen::MatrixBase<DerivedF> & F,
  80. Eigen::PlainObjectBase<DerivedV_uv> & V_uv)
  81. {
  82. using Scalar = typename DerivedV_uv::Scalar;
  83. Eigen::SparseMatrix<Scalar> Q;
  84. igl::lscm_hessian(V,F,Q);
  85. Eigen::SparseMatrix<Scalar> M;
  86. igl::massmatrix(V,F,igl::MASSMATRIX_TYPE_DEFAULT,M);
  87. Eigen::SparseMatrix<Scalar> M2;
  88. igl::repdiag(M,2,M2);
  89. Eigen::Matrix<Scalar,Eigen::Dynamic,Eigen::Dynamic> U;
  90. Eigen::Matrix<Scalar,Eigen::Dynamic,1> S;
  91. bool success = igl::eigs(Q,M2,3,igl::EIGS_TYPE_SM,U,S);
  92. if(!success)
  93. {
  94. return false;
  95. }
  96. V_uv.resize(V.rows(),2);
  97. V_uv<< U.col(0).head(V.rows()),U.col(0).tail(V.rows());
  98. return true;
  99. }
  100. #ifdef IGL_STATIC_LIBRARY
  101. // Explicit template instantiation
  102. // generated by autoexplicit.sh
  103. template bool igl::lscm<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  104. // generated by autoexplicit.sh
  105. template bool igl::lscm<Eigen::Matrix<double, -1, -1, 0, -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, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -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, 0, -1, -1> >&);
  106. #endif