lu_lagrange.cpp 3.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137
  1. #include "lu_lagrange.h"
  2. // Cholesky LLT decomposition for symmetric positive definite
  3. #include <Eigen/SparseExtra>
  4. #include <cassert>
  5. #include "find.h"
  6. #include "sparse.h"
  7. template <typename T>
  8. IGL_INLINE bool igl::lu_lagrange(
  9. const Eigen::SparseMatrix<T> & ATA,
  10. const Eigen::SparseMatrix<T> & C,
  11. Eigen::SparseMatrix<T> & L,
  12. Eigen::SparseMatrix<T> & U)
  13. {
  14. // number of unknowns
  15. int n = ATA.rows();
  16. // number of lagrange multipliers
  17. int m = C.cols();
  18. assert(ATA.cols() == n);
  19. if(m != 0)
  20. {
  21. assert(C.rows() == n);
  22. if(C.nonZeros() == 0)
  23. {
  24. // See note above about empty columns in C
  25. fprintf(stderr,"Error: lu_lagrange() C has columns but no entries\n");
  26. return false;
  27. }
  28. }
  29. // Check that each column of C has at least one entry
  30. std::vector<bool> has_entry; has_entry.resize(C.cols(),false);
  31. // Iterate over outside
  32. for(int k=0; k<C.outerSize(); ++k)
  33. {
  34. // Iterate over inside
  35. for(typename Eigen::SparseMatrix<T>::InnerIterator it (C,k); it; ++it)
  36. {
  37. has_entry[it.col()] = true;
  38. }
  39. }
  40. for(int i=0;i<(int)has_entry.size();i++)
  41. {
  42. if(!has_entry[i])
  43. {
  44. // See note above about empty columns in C
  45. fprintf(stderr,"Error: lu_lagrange() C(:,%d) has no entries\n",i);
  46. return false;
  47. }
  48. }
  49. // Cholesky factorization of ATA
  50. //// Eigen fails if you give a full view of the matrix like this:
  51. //Eigen::SparseLLT<SparseMatrix<T> > ATA_LLT(ATA);
  52. Eigen::SparseMatrix<T> ATA_LT = ATA.template triangularView<Eigen::Lower>();
  53. Eigen::SparseLLT<Eigen::SparseMatrix<T> > ATA_LLT(ATA_LT);
  54. Eigen::SparseMatrix<T> J = ATA_LLT.matrixL();
  55. //if(!ATA_LLT.succeeded())
  56. if(!((J*0).eval().nonZeros() == 0))
  57. {
  58. fprintf(stderr,"Error: lu_lagrange() failed to factor ATA\n");
  59. return false;
  60. }
  61. if(m == 0)
  62. {
  63. // If there are no constraints (C is empty) then LU decomposition is just L
  64. // and L' from cholesky decomposition
  65. L = J;
  66. U = J.transpose();
  67. }else
  68. {
  69. // Construct helper matrix M
  70. Eigen::SparseMatrix<T> M = C;
  71. J.template triangularView<Eigen::Lower>().solveInPlace(M);
  72. // Compute cholesky factorizaiton of M'*M
  73. Eigen::SparseMatrix<T> MTM = M.transpose() * M;
  74. Eigen::SparseLLT<Eigen::SparseMatrix<T> > MTM_LLT(MTM.template triangularView<Eigen::Lower>());
  75. Eigen::SparseMatrix<T> K = MTM_LLT.matrixL();
  76. //if(!MTM_LLT.succeeded())
  77. if(!((K*0).eval().nonZeros() == 0))
  78. {
  79. fprintf(stderr,"Error: lu_lagrange() failed to factor MTM\n");
  80. return false;
  81. }
  82. // assemble LU decomposition of Q
  83. Eigen::Matrix<int,Eigen::Dynamic,1> MI;
  84. Eigen::Matrix<int,Eigen::Dynamic,1> MJ;
  85. Eigen::Matrix<T,Eigen::Dynamic,1> MV;
  86. igl::find(M,MI,MJ,MV);
  87. Eigen::Matrix<int,Eigen::Dynamic,1> KI;
  88. Eigen::Matrix<int,Eigen::Dynamic,1> KJ;
  89. Eigen::Matrix<T,Eigen::Dynamic,1> KV;
  90. igl::find(K,KI,KJ,KV);
  91. Eigen::Matrix<int,Eigen::Dynamic,1> JI;
  92. Eigen::Matrix<int,Eigen::Dynamic,1> JJ;
  93. Eigen::Matrix<T,Eigen::Dynamic,1> JV;
  94. igl::find(J,JI,JJ,JV);
  95. int nnz = JV.size() + MV.size() + KV.size();
  96. Eigen::Matrix<int,Eigen::Dynamic,1> UI(nnz);
  97. Eigen::Matrix<int,Eigen::Dynamic,1> UJ(nnz);
  98. Eigen::Matrix<T,Eigen::Dynamic,1> UV(nnz);
  99. UI << JJ, MI, (KJ.array() + n).matrix();
  100. UJ << JI, (MJ.array() + n).matrix(), (KI.array() + n).matrix();
  101. UV << JV, MV, KV*-1;
  102. igl::sparse(UI,UJ,UV,U);
  103. Eigen::Matrix<int,Eigen::Dynamic,1> LI(nnz);
  104. Eigen::Matrix<int,Eigen::Dynamic,1> LJ(nnz);
  105. Eigen::Matrix<T,Eigen::Dynamic,1> LV(nnz);
  106. LI << JI, (MJ.array() + n).matrix(), (KI.array() + n).matrix();
  107. LJ << JJ, MI, (KJ.array() + n).matrix();
  108. LV << JV, MV, KV;
  109. igl::sparse(LI,LJ,LV,L);
  110. }
  111. return true;
  112. }
  113. #ifndef IGL_HEADER_ONLY
  114. // Explicit template specialization
  115. #endif