// This file is part of libigl, a simple c++ geometry processing library. // // Copyright (C) 2013 Alec Jacobson // // This Source Code Form is subject to the terms of the Mozilla Public License // v. 2.0. If a copy of the MPL was not distributed with this file, You can // obtain one at http://mozilla.org/MPL/2.0/. #include "lu_lagrange.h" // Cholesky LLT decomposition for symmetric positive definite //#include // Bug in unsupported/Eigen/SparseExtra needs iostream first #include #include #include #include #include "find.h" #include "sparse.h" template IGL_INLINE bool igl::lu_lagrange( const Eigen::SparseMatrix & ATA, const Eigen::SparseMatrix & C, Eigen::SparseMatrix & L, Eigen::SparseMatrix & U) { #if EIGEN_VERSION_AT_LEAST(3,0,92) #if defined(_WIN32) #pragma message("lu_lagrange has not yet been implemented for your Eigen Version") #else #warning lu_lagrange has not yet been implemented for your Eigen Version #endif return false; #else // number of unknowns int n = ATA.rows(); // number of lagrange multipliers int m = C.cols(); assert(ATA.cols() == n); if(m != 0) { assert(C.rows() == n); if(C.nonZeros() == 0) { // See note above about empty columns in C fprintf(stderr,"Error: lu_lagrange() C has columns but no entries\n"); return false; } } // Check that each column of C has at least one entry std::vector has_entry; has_entry.resize(C.cols(),false); // Iterate over outside for(int k=0; k::InnerIterator it (C,k); it; ++it) { has_entry[it.col()] = true; } } for(int i=0;i<(int)has_entry.size();i++) { if(!has_entry[i]) { // See note above about empty columns in C fprintf(stderr,"Error: lu_lagrange() C(:,%d) has no entries\n",i); return false; } } // Cholesky factorization of ATA //// Eigen fails if you give a full view of the matrix like this: //Eigen::SparseLLT > ATA_LLT(ATA); Eigen::SparseMatrix ATA_LT = ATA.template triangularView(); Eigen::SparseLLT > ATA_LLT(ATA_LT); Eigen::SparseMatrix J = ATA_LLT.matrixL(); //if(!ATA_LLT.succeeded()) if(!((J*0).eval().nonZeros() == 0)) { fprintf(stderr,"Error: lu_lagrange() failed to factor ATA\n"); return false; } if(m == 0) { // If there are no constraints (C is empty) then LU decomposition is just L // and L' from cholesky decomposition L = J; U = J.transpose(); }else { // Construct helper matrix M Eigen::SparseMatrix M = C; J.template triangularView().solveInPlace(M); // Compute cholesky factorizaiton of M'*M Eigen::SparseMatrix MTM = M.transpose() * M; Eigen::SparseLLT > MTM_LLT(MTM.template triangularView()); Eigen::SparseMatrix K = MTM_LLT.matrixL(); //if(!MTM_LLT.succeeded()) if(!((K*0).eval().nonZeros() == 0)) { fprintf(stderr,"Error: lu_lagrange() failed to factor MTM\n"); return false; } // assemble LU decomposition of Q Eigen::Matrix MI; Eigen::Matrix MJ; Eigen::Matrix MV; igl::find(M,MI,MJ,MV); Eigen::Matrix KI; Eigen::Matrix KJ; Eigen::Matrix KV; igl::find(K,KI,KJ,KV); Eigen::Matrix JI; Eigen::Matrix JJ; Eigen::Matrix JV; igl::find(J,JI,JJ,JV); int nnz = JV.size() + MV.size() + KV.size(); Eigen::Matrix UI(nnz); Eigen::Matrix UJ(nnz); Eigen::Matrix UV(nnz); UI << JJ, MI, (KJ.array() + n).matrix(); UJ << JI, (MJ.array() + n).matrix(), (KI.array() + n).matrix(); UV << JV, MV, KV*-1; igl::sparse(UI,UJ,UV,U); Eigen::Matrix LI(nnz); Eigen::Matrix LJ(nnz); Eigen::Matrix LV(nnz); LI << JI, (MJ.array() + n).matrix(), (KI.array() + n).matrix(); LJ << JJ, MI, (KJ.array() + n).matrix(); LV << JV, MV, KV; igl::sparse(LI,LJ,LV,L); } return true; #endif } #ifdef IGL_STATIC_LIBRARY // Explicit template specialization template bool igl::lu_lagrange(Eigen::SparseMatrix const&, Eigen::SparseMatrix const&, Eigen::SparseMatrix&, Eigen::SparseMatrix&); #endif