#include "min_quad_with_fixed.h" #include "slice.h" #include "is_symmetric.h" #include "find.h" #include "sparse.h" #include "repmat.h" #include "lu_lagrange.h" #include "full.h" #include "EPS.h" //#include // Bug in unsupported/Eigen/SparseExtra needs iostream first #include #include #include #include #include template IGL_INLINE bool igl::min_quad_with_fixed_precompute( const Eigen::SparseMatrix& A, const Eigen::Matrix & known, const Eigen::SparseMatrix& Aeq, const bool pd, igl::min_quad_with_fixed_data & data ) { using namespace Eigen; using namespace std; // number of rows int n = A.rows(); // cache problem size data.n = n; int neq = Aeq.rows(); // default is to have 0 linear equality constraints if(Aeq.size() != 0) { assert(n == Aeq.cols()); } assert(A.rows() == n); assert(A.cols() == n); // number of known rows int kr = known.size(); assert(kr == 0 || known.minCoeff() >= 0); assert(kr == 0 || known.maxCoeff() < n); assert(neq <= n); // cache known data.known = known; // get list of unknown indices data.unknown.resize(n-kr); std::vector unknown_mask; unknown_mask.resize(n,true); for(int i = 0;i Auu; igl::slice(A,data.unknown,data.unknown,Auu); // determine if A(unknown,unknown) is symmetric and/or positive definite data.Auu_sym = igl::is_symmetric(Auu,igl::EPS()); // Positive definiteness is *not* determined, rather it is given as a // parameter data.Auu_pd = pd; // Append lagrange multiplier quadratic terms SparseMatrix new_A; Matrix AI; Matrix AJ; Matrix AV; igl::find(A,AI,AJ,AV); Matrix AeqI(0,1); Matrix AeqJ(0,1); Matrix AeqV(0,1); if(neq > 0) { igl::find(Aeq,AeqI,AeqJ,AeqV); } Matrix new_AI(AV.size()+AeqV.size()*2); Matrix new_AJ(AV.size()+AeqV.size()*2); Matrix new_AV(AV.size()+AeqV.size()*2); new_AI << AI, (AeqI.array()+n).matrix(), AeqJ; new_AJ << AJ, AeqJ, (AeqI.array()+n).matrix(); new_AV << AV, AeqV, AeqV; igl::sparse(new_AI,new_AJ,new_AV,n+neq,n+neq,new_A); // precompute RHS builders if(kr > 0) { SparseMatrix Aulk; igl::slice(new_A,data.unknown_lagrange,data.known,Aulk); SparseMatrix Akul; igl::slice(new_A,data.known,data.unknown_lagrange,Akul); //// This doesn't work!!! //data.preY = Aulk + Akul.transpose(); SparseMatrix AkulT = Akul.transpose(); data.preY = Aulk + AkulT; }else { data.preY.resize(data.unknown_lagrange.size(),0); } // Create factorization //if(data.Auu_sym && data.Auu_pd) //{ // data.sparse = true; // SparseMatrix Aequ(0,0); // if(neq>0) // { // Matrix Aeq_rows; // Aeq_rows.resize(neq); // for(int i = 0;i AequT = Aequ.transpose(); // // Compute LU decomposition // // TODO: This should be using SimplicialLDLT // // Q: Does SimplicialLDLT work on "Hermitian indefinite matrices" the way // // that matlabs ldl does? // // A: Maybe not. The eigen documentation writes, "This class provides a // // LDL^T Cholesky factorizations without square root of sparse matrices // // that are selfadjoint and positive definite" // bool lu_success = igl::lu_lagrange(Auu,AequT,data.L,data.U); // if(!lu_success) // { // return false; // } //}else //{ SparseMatrix NA; igl::slice(new_A,data.unknown_lagrange,data.unknown_lagrange,NA); if(data.Auu_sym) { SimplicialLDLT > ldlt_solver; ldlt_solver.compute(NA); if(ldlt_solver.info() != Eigen::Success) { return false; } // Implementation incomplete assert(false); cerr<<"ERROR min_quad_with_fixed_precompute implementation incomplete"< NAfull; igl::full(NA,NAfull); data.lu = FullPivLU >(NAfull); if(!data.lu.isInvertible()) { fprintf(stderr, "Error: min_quad_with_fixed_precompute() LU not invertible\n"); return false; } } //} return true; } template IGL_INLINE bool igl::min_quad_with_fixed_solve( const igl::min_quad_with_fixed_data & data, const Eigen::Matrix & B, const Eigen::PlainObjectBase & Y, const Eigen::Matrix & Beq, Eigen::PlainObjectBase & Z) { // number of known rows int kr = data.known.size(); if(kr!=0) { assert(kr == Y.rows()); } // number of columns to solve int cols = Y.cols(); // number of lagrange multipliers aka linear equality constraints int neq = data.lagrange.size(); // append lagrange multiplier rhs's Eigen::Matrix BBeq(B.size() + Beq.size()); BBeq << B, (Beq*-2.0); // Build right hand side Eigen::Matrix BBequl; igl::slice(BBeq,data.unknown_lagrange,BBequl); Eigen::Matrix BBequlcols; igl::repmat(BBequl,1,cols,BBequlcols); Eigen::Matrix NB; if(kr == 0) { NB = BBequlcols; }else { NB = data.preY * Y + BBequlcols; } // resize output Z.resize(data.n,cols); // Set known values for(int i = 0;i < kr;i++) { for(int j = 0;j < cols;j++) { Z(data.known(i),j) = Y(i,j); } } //std::cout<<"NB=["<().solveInPlace(NB); data.U.template triangularView().solveInPlace(NB); }else { NB = data.lu.solve(NB); } // Now NB contains sol/-0.5 NB *= -0.5; // Now NB contains solution // Place solution in Z for(int i = 0;i<(NB.rows()-neq);i++) { for(int j = 0;j, Eigen::Matrix >(igl::min_quad_with_fixed_data const&, Eigen::Matrix const&, Eigen::PlainObjectBase > const&, Eigen::Matrix const&, Eigen::PlainObjectBase >&); template bool igl::min_quad_with_fixed_precompute(Eigen::SparseMatrix const&, Eigen::Matrix const&, Eigen::SparseMatrix const&, bool, igl::min_quad_with_fixed_data&); #endif