123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530 |
- #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 "matlab_format.h"
- #include "EPS.h"
- #include "cat.h"
- //#include <Eigen/SparseExtra>
- // Bug in unsupported/Eigen/SparseExtra needs iostream first
- #include <iostream>
- #include <unsupported/Eigen/SparseExtra>
- #include <cassert>
- #include <cstdio>
- #include <iostream>
- template <typename T, typename Derivedknown>
- IGL_INLINE bool igl::min_quad_with_fixed_precompute(
- const Eigen::SparseMatrix<T>& A,
- const Eigen::PlainObjectBase<Derivedknown> & known,
- const Eigen::SparseMatrix<T>& Aeq,
- const bool pd,
- min_quad_with_fixed_data<T> & data
- )
- {
- #define MIN_QUAD_WITH_FIXED_CPP_DEBUG
- using namespace Eigen;
- using namespace std;
- using namespace igl;
- #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
- cout<<" pre"<<endl;
- #endif
- // 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<bool> unknown_mask;
- unknown_mask.resize(n,true);
- for(int i = 0;i<kr;i++)
- {
- unknown_mask[known(i)] = false;
- }
- int u = 0;
- for(int i = 0;i<n;i++)
- {
- if(unknown_mask[i])
- {
- data.unknown(u) = i;
- u++;
- }
- }
- // get list of lagrange multiplier indices
- data.lagrange.resize(neq);
- for(int i = 0;i<neq;i++)
- {
- data.lagrange(i) = n + i;
- }
- // cache unknown followed by lagrange indices
- data.unknown_lagrange.resize(data.unknown.size()+data.lagrange.size());
- data.unknown_lagrange << data.unknown, data.lagrange;
- SparseMatrix<T> Auu;
- slice(A,data.unknown,data.unknown,Auu);
- // Positive definiteness is *not* determined, rather it is given as a
- // parameter
- data.Auu_pd = pd;
- if(data.Auu_pd)
- {
- // PD implies symmetric
- data.Auu_sym = true;
- // This is an annoying assertion unless EPS can be chosen in a nicer way.
- //assert(is_symmetric(Auu,EPS<double>()));
- assert(is_symmetric(Auu,1.0));
- }else
- {
- // determine if A(unknown,unknown) is symmetric and/or positive definite
- data.Auu_sym = is_symmetric(Auu,EPS<double>());
- }
- // Determine number of linearly independent constraints
- int nc = 0;
- if(neq>0)
- {
- #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
- cout<<" qr"<<endl;
- #endif
- // QR decomposition to determine row rank in Aequ
- slice(Aeq,data.unknown,2,data.Aequ);
- assert(data.Aequ.rows() == neq);
- assert(data.Aequ.cols() == data.unknown.size());
- data.AeqTQR.compute(data.Aequ.transpose().eval());
- cout<<endl<<matlab_format(SparseMatrix<T>(data.Aequ.transpose().eval()),"AeqT")<<endl<<endl;
- switch(data.AeqTQR.info())
- {
- case Eigen::Success:
- break;
- case Eigen::NumericalIssue:
- cerr<<"Error: Numerical issue."<<endl;
- return false;
- case Eigen::InvalidInput:
- cerr<<"Error: Invalid input."<<endl;
- return false;
- default:
- cerr<<"Error: Other."<<endl;
- return false;
- }
- nc = data.AeqTQR.rank();
- assert(nc<=neq);
- data.Aeq_li = nc == neq;
- }else
- {
- data.Aeq_li = true;
- }
- if(data.Aeq_li)
- {
- #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
- cout<<" Aeq_li=true"<<endl;
- #endif
- // Append lagrange multiplier quadratic terms
- SparseMatrix<T> new_A;
- SparseMatrix<T> AeqT = Aeq.transpose();
- SparseMatrix<T> Z(neq,neq);
- // This is a bit slower. But why isn't cat fast?
- new_A = cat(1, cat(2, A, AeqT ),
- cat(2, Aeq, Z ));
- // precompute RHS builders
- if(kr > 0)
- {
- SparseMatrix<T> Aulk,Akul;
- // Slow
- slice(new_A,data.unknown_lagrange,data.known,Aulk);
- //// This doesn't work!!!
- //data.preY = Aulk + Akul.transpose();
- // Slow
- if(data.Auu_sym)
- {
- data.preY = Aulk*2;
- }else
- {
- slice(new_A,data.known,data.unknown_lagrange,Akul);
- SparseMatrix<T> AkulT = Akul.transpose();
- data.preY = Aulk + AkulT;
- }
- }else
- {
- data.preY.resize(data.unknown_lagrange.size(),0);
- }
- // Positive definite and no equality constraints (Postive definiteness
- // implies symmetric)
- #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
- cout<<" factorize"<<endl;
- #endif
- if(data.Auu_pd && neq == 0)
- {
- #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
- cout<<" llt"<<endl;
- #endif
- data.llt.compute(Auu);
- switch(data.llt.info())
- {
- case Eigen::Success:
- break;
- case Eigen::NumericalIssue:
- cerr<<"Error: Numerical issue."<<endl;
- return false;
- default:
- cerr<<"Error: Other."<<endl;
- return false;
- }
- data.solver_type = min_quad_with_fixed_data<T>::LLT;
- }else
- {
- #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
- cout<<" ldlt"<<endl;
- #endif
- // Either not PD or there are equality constraints
- SparseMatrix<T> NA;
- slice(new_A,data.unknown_lagrange,data.unknown_lagrange,NA);
- data.NA = NA;
- // Ideally we'd use LDLT but Eigen doesn't support positive semi-definite
- // matrices:
- // http://forum.kde.org/viewtopic.php?f=74&t=106962&p=291990#p291990
- if(data.Auu_sym && false)
- {
- data.ldlt.compute(NA);
- switch(data.ldlt.info())
- {
- case Eigen::Success:
- break;
- case Eigen::NumericalIssue:
- cerr<<"Error: Numerical issue."<<endl;
- return false;
- default:
- cerr<<"Error: Other."<<endl;
- return false;
- }
- data.solver_type = min_quad_with_fixed_data<T>::LDLT;
- }else
- {
- #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
- cout<<" lu"<<endl;
- #endif
- // Resort to LU
- // Bottleneck >1/2
- data.lu.compute(NA);
- //std::cout<<"NA=["<<std::endl<<NA<<std::endl<<"];"<<std::endl;
- switch(data.lu.info())
- {
- case Eigen::Success:
- break;
- case Eigen::NumericalIssue:
- cerr<<"Error: Numerical issue."<<endl;
- return false;
- case Eigen::InvalidInput:
- cerr<<"Error: Invalid Input."<<endl;
- return false;
- default:
- cerr<<"Error: Other."<<endl;
- return false;
- }
- data.solver_type = min_quad_with_fixed_data<T>::LU;
- }
- }
- }else
- {
- #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
- cout<<" Aeq_li=false"<<endl;
- #endif
- data.neq = neq;
- const int nu = data.unknown.size();
- //cout<<"nu: "<<nu<<endl;
- //cout<<"neq: "<<neq<<endl;
- //cout<<"nc: "<<nc<<endl;
- //cout<<" matrixR"<<endl;
- SparseMatrix<T> AeqTR,AeqTQ;
- AeqTR = data.AeqTQR.matrixR();
- // This shouldn't be necessary
- AeqTR.prune(0.0);
- #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
- cout<<" matrixQ"<<endl;
- #endif
- // THIS IS ESSENTIALLY DENSE AND THIS IS BY FAR THE BOTTLENECK
- // http://forum.kde.org/viewtopic.php?f=74&t=117500
- AeqTQ = data.AeqTQR.matrixQ();
- #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
- cout<<" prune"<<endl;
- cout<<" nnz: "<<AeqTQ.nonZeros()<<endl;
- #endif
- // This shouldn't be necessary
- AeqTQ.prune(0.0);
- //cout<<"AeqTQ: "<<AeqTQ.rows()<<" "<<AeqTQ.cols()<<endl;
- //cout<<matlab_format(AeqTQ,"AeqTQ")<<endl;
- //cout<<" perms"<<endl;
- #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
- cout<<" nnz: "<<AeqTQ.nonZeros()<<endl;
- cout<<" perm"<<endl;
- #endif
- SparseMatrix<double> I(neq,neq);
- I.setIdentity();
- data.AeqTE = data.AeqTQR.colsPermutation() * I;
- data.AeqTET = data.AeqTQR.colsPermutation().transpose() * I;
- assert(AeqTR.rows() == neq);
- assert(AeqTQ.rows() == nu);
- assert(AeqTQ.cols() == nu);
- assert(AeqTR.cols() == neq);
- //cout<<" slice"<<endl;
- #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
- cout<<" slice"<<endl;
- #endif
- data.AeqTQ1 = AeqTQ.topLeftCorner(nu,nc);
- data.AeqTQ1T = data.AeqTQ1.transpose().eval();
- // ALREADY TRIM (Not 100% sure about this)
- data.AeqTR1 = AeqTR.topLeftCorner(nc,nc);
- data.AeqTR1T = data.AeqTR1.transpose().eval();
- //cout<<"AeqTR1T.size() "<<data.AeqTR1T.rows()<<" "<<data.AeqTR1T.cols()<<endl;
- // Null space
- data.AeqTQ2 = AeqTQ.bottomRightCorner(nu,nu-nc);
- data.AeqTQ2T = data.AeqTQ2.transpose().eval();
- #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
- cout<<" proj"<<endl;
- #endif
- // Projected hessian
- SparseMatrix<T> QRAuu = data.AeqTQ2T * Auu * data.AeqTQ2;
- {
- #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
- cout<<" factorize"<<endl;
- #endif
- // QRAuu should always be PD
- data.llt.compute(QRAuu);
- switch(data.llt.info())
- {
- case Eigen::Success:
- break;
- case Eigen::NumericalIssue:
- cerr<<"Error: Numerical issue."<<endl;
- return false;
- default:
- cerr<<"Error: Other."<<endl;
- return false;
- }
- data.solver_type = min_quad_with_fixed_data<T>::QR_LLT;
- }
- #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
- cout<<" smash"<<endl;
- #endif
- // Known value multiplier
- SparseMatrix<T> Auk;
- slice(A,data.unknown,data.known,Auk);
- SparseMatrix<T> Aku;
- slice(A,data.known,data.unknown,Aku);
- SparseMatrix<T> AkuT = Aku.transpose();
- data.preY = Auk + AkuT;
- // Needed during solve
- data.Auu = Auu;
- slice(Aeq,data.known,2,data.Aeqk);
- assert(data.Aeqk.rows() == neq);
- assert(data.Aeqk.cols() == data.known.size());
- }
- return true;
- }
- template <
- typename T,
- typename DerivedB,
- typename DerivedY,
- typename DerivedBeq,
- typename DerivedZ,
- typename Derivedsol>
- IGL_INLINE bool igl::min_quad_with_fixed_solve(
- const min_quad_with_fixed_data<T> & data,
- const Eigen::PlainObjectBase<DerivedB> & B,
- const Eigen::PlainObjectBase<DerivedY> & Y,
- const Eigen::PlainObjectBase<DerivedBeq> & Beq,
- Eigen::PlainObjectBase<DerivedZ> & Z,
- Eigen::PlainObjectBase<Derivedsol> & sol)
- {
- using namespace std;
- using namespace Eigen;
- using namespace igl;
- // 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();
- assert(B.cols() == 1);
- assert(Beq.size() == 0 || Beq.cols() == 1);
- // 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);
- }
- }
- if(data.Aeq_li)
- {
- // number of lagrange multipliers aka linear equality constraints
- int neq = data.lagrange.size();
- // append lagrange multiplier rhs's
- Eigen::Matrix<T,Eigen::Dynamic,1> BBeq(B.size() + Beq.size());
- BBeq << B, (Beq*-2.0);
- // Build right hand side
- Eigen::Matrix<T,Eigen::Dynamic,1> BBequl;
- igl::slice(BBeq,data.unknown_lagrange,BBequl);
- Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> BBequlcols;
- igl::repmat(BBequl,1,cols,BBequlcols);
- Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> NB;
- if(kr == 0)
- {
- NB = BBequlcols;
- }else
- {
- NB = data.preY * Y + BBequlcols;
- }
- //std::cout<<"NB=["<<std::endl<<NB<<std::endl<<"];"<<std::endl;
- //cout<<matlab_format(NB,"NB")<<endl;
- switch(data.solver_type)
- {
- case igl::min_quad_with_fixed_data<T>::LLT:
- sol = data.llt.solve(NB);
- break;
- case igl::min_quad_with_fixed_data<T>::LDLT:
- sol = data.ldlt.solve(NB);
- break;
- case igl::min_quad_with_fixed_data<T>::LU:
- // Not a bottleneck
- sol = data.lu.solve(NB);
- break;
- default:
- cerr<<"Error: invalid solver type"<<endl;
- return false;
- }
- //std::cout<<"sol=["<<std::endl<<sol<<std::endl<<"];"<<std::endl;
- // Now sol contains sol/-0.5
- sol *= -0.5;
- // Now sol contains solution
- // Place solution in Z
- for(int i = 0;i<(sol.rows()-neq);i++)
- {
- for(int j = 0;j<sol.cols();j++)
- {
- Z(data.unknown_lagrange(i),j) = sol(i,j);
- }
- }
- }else
- {
- assert(data.solver_type == min_quad_with_fixed_data<T>::QR_LLT);
- PlainObjectBase<DerivedBeq> eff_Beq;
- // Adjust Aeq rhs to include known parts
- eff_Beq =
- //data.AeqTQR.colsPermutation().transpose() * (-data.Aeqk * Y + Beq);
- data.AeqTET * (-data.Aeqk * Y + Beq);
- // Where did this -0.5 come from? Probably the same place as above.
- PlainObjectBase<DerivedB> Bu;
- slice(B,data.unknown,Bu);
- PlainObjectBase<DerivedB> NB;
- NB = -0.5*(Bu + data.preY * Y);
- // Trim eff_Beq
- const int nc = data.AeqTQR.rank();
- const int neq = Beq.rows();
- eff_Beq = eff_Beq.topLeftCorner(nc,1).eval();
- data.AeqTR1T.template triangularView<Lower>().solveInPlace(eff_Beq);
- // Now eff_Beq = (data.AeqTR1T \ (data.AeqTET * (-data.Aeqk * Y + Beq)))
- PlainObjectBase<DerivedB> lambda_0;
- lambda_0 = data.AeqTQ1 * eff_Beq;
- //cout<<matlab_format(lambda_0,"lambda_0")<<endl;
- PlainObjectBase<DerivedB> QRB;
- QRB = -data.AeqTQ2T * (data.Auu * lambda_0) + data.AeqTQ2T * NB;
- PlainObjectBase<Derivedsol> lambda;
- lambda = data.llt.solve(QRB);
- // prepare output
- PlainObjectBase<Derivedsol> solu;
- solu = data.AeqTQ2 * lambda + lambda_0;
- // http://www.math.uh.edu/~rohop/fall_06/Chapter3.pdf
- PlainObjectBase<Derivedsol> solLambda;
- {
- PlainObjectBase<Derivedsol> temp1,temp2;
- temp1 = (data.AeqTQ1T * NB - data.AeqTQ1T * data.Auu * solu);
- data.AeqTR1.template triangularView<Upper>().solveInPlace(temp1);
- //cout<<matlab_format(temp1,"temp1")<<endl;
- temp2 = PlainObjectBase<Derivedsol>::Zero(neq,1);
- temp2.topLeftCorner(nc,1) = temp1;
- //solLambda = data.AeqTQR.colsPermutation() * temp2;
- solLambda = data.AeqTE * temp2;
- }
- // sol is [Z(unknown);Lambda]
- assert(data.unknown.size() == solu.rows());
- assert(cols == solu.cols());
- assert(data.neq == neq);
- assert(data.neq == solLambda.rows());
- assert(cols == solLambda.cols());
- sol.resize(data.unknown.size()+data.neq,cols);
- sol.block(0,0,solu.rows(),solu.cols()) = solu;
- sol.block(solu.rows(),0,solLambda.rows(),solLambda.cols()) = solLambda;
- for(int u = 0;u<data.unknown.size();u++)
- {
- for(int j = 0;j<Z.cols();j++)
- {
- Z(data.unknown(u),j) = solu(u,j);
- }
- }
- }
- return true;
- }
- template <
- typename T,
- typename DerivedB,
- typename DerivedY,
- typename DerivedBeq,
- typename DerivedZ>
- IGL_INLINE bool igl::min_quad_with_fixed_solve(
- const min_quad_with_fixed_data<T> & data,
- const Eigen::PlainObjectBase<DerivedB> & B,
- const Eigen::PlainObjectBase<DerivedY> & Y,
- const Eigen::PlainObjectBase<DerivedBeq> & Beq,
- Eigen::PlainObjectBase<DerivedZ> & Z)
- {
- Eigen::PlainObjectBase<DerivedZ> sol;
- return min_quad_with_fixed_solve(data,B,Y,Beq,Z,sol);
- }
- #ifndef IGL_HEADER_ONLY
- // Explicit template specialization
- template bool igl::min_quad_with_fixed_solve<double, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(igl::min_quad_with_fixed_data<double> const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
- template bool igl::min_quad_with_fixed_precompute<double, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::SparseMatrix<double, 0, int> const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::SparseMatrix<double, 0, int> const&, bool, igl::min_quad_with_fixed_data<double>&);
- template bool igl::min_quad_with_fixed_solve<double, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(igl::min_quad_with_fixed_data<double> const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
- #endif
|