min_quad_with_fixed.cpp 8.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321
  1. #include "min_quad_with_fixed.h"
  2. #include "slice.h"
  3. #include "is_symmetric.h"
  4. #include "find.h"
  5. #include "sparse.h"
  6. #include "repmat.h"
  7. #include "lu_lagrange.h"
  8. #include "full.h"
  9. #include "matlab_format.h"
  10. #include "EPS.h"
  11. #include "cat.h"
  12. //#include <Eigen/SparseExtra>
  13. // Bug in unsupported/Eigen/SparseExtra needs iostream first
  14. #include <iostream>
  15. #include <unsupported/Eigen/SparseExtra>
  16. #include <cassert>
  17. #include <cstdio>
  18. #include <iostream>
  19. template <typename T>
  20. IGL_INLINE bool igl::min_quad_with_fixed_precompute(
  21. const Eigen::SparseMatrix<T>& A,
  22. const Eigen::Matrix<int,Eigen::Dynamic,1> & known,
  23. const Eigen::SparseMatrix<T>& Aeq,
  24. const bool pd,
  25. igl::min_quad_with_fixed_data<T> & data
  26. )
  27. {
  28. using namespace Eigen;
  29. using namespace std;
  30. using namespace igl;
  31. // number of rows
  32. int n = A.rows();
  33. // cache problem size
  34. data.n = n;
  35. int neq = Aeq.rows();
  36. // default is to have 0 linear equality constraints
  37. if(Aeq.size() != 0)
  38. {
  39. assert(n == Aeq.cols());
  40. }
  41. assert(A.rows() == n);
  42. assert(A.cols() == n);
  43. // number of known rows
  44. int kr = known.size();
  45. assert(kr == 0 || known.minCoeff() >= 0);
  46. assert(kr == 0 || known.maxCoeff() < n);
  47. assert(neq <= n);
  48. // cache known
  49. data.known = known;
  50. // get list of unknown indices
  51. data.unknown.resize(n-kr);
  52. std::vector<bool> unknown_mask;
  53. unknown_mask.resize(n,true);
  54. for(int i = 0;i<kr;i++)
  55. {
  56. unknown_mask[known(i)] = false;
  57. }
  58. int u = 0;
  59. for(int i = 0;i<n;i++)
  60. {
  61. if(unknown_mask[i])
  62. {
  63. data.unknown(u) = i;
  64. u++;
  65. }
  66. }
  67. // get list of lagrange multiplier indices
  68. data.lagrange.resize(neq);
  69. for(int i = 0;i<neq;i++)
  70. {
  71. data.lagrange(i) = n + i;
  72. }
  73. // cache unknown followed by lagrange indices
  74. data.unknown_lagrange.resize(data.unknown.size()+data.lagrange.size());
  75. data.unknown_lagrange << data.unknown, data.lagrange;
  76. SparseMatrix<T> Auu;
  77. igl::slice(A,data.unknown,data.unknown,Auu);
  78. // Positive definiteness is *not* determined, rather it is given as a
  79. // parameter
  80. data.Auu_pd = pd;
  81. if(data.Auu_pd)
  82. {
  83. // PD implies symmetric
  84. data.Auu_sym = true;
  85. assert(igl::is_symmetric(Auu,igl::EPS<double>()));
  86. }else
  87. {
  88. // determine if A(unknown,unknown) is symmetric and/or positive definite
  89. data.Auu_sym = igl::is_symmetric(Auu,igl::EPS<double>());
  90. }
  91. // Append lagrange multiplier quadratic terms
  92. SparseMatrix<T> new_A;
  93. SparseMatrix<T> AeqT = Aeq.transpose();
  94. SparseMatrix<T> Z(neq,neq);
  95. // This is a bit slower. But why isn't cat fast?
  96. new_A =
  97. cat(1,
  98. cat(2, A, AeqT ),
  99. cat(2, Aeq, Z ));
  100. //cout<<matlab_format(new_A,"new_A")<<endl;
  101. //Matrix<int,Dynamic,1> AI;
  102. //Matrix<int,Dynamic,1> AJ;
  103. //Matrix<T,Dynamic,1> AV;
  104. //// Slow
  105. //igl::find(A,AI,AJ,AV);
  106. //Matrix<int,Dynamic,1> AeqI(0,1);
  107. //Matrix<int,Dynamic,1> AeqJ(0,1);
  108. //Matrix<T,Dynamic,1> AeqV(0,1);
  109. //if(neq > 0)
  110. //{
  111. // // Slow
  112. // igl::find(Aeq,AeqI,AeqJ,AeqV);
  113. //}
  114. //Matrix<int,Dynamic,1> new_AI(AV.size()+AeqV.size()*2);
  115. //Matrix<int,Dynamic,1> new_AJ(AV.size()+AeqV.size()*2);
  116. //Matrix<T,Dynamic,1> new_AV(AV.size()+AeqV.size()*2);
  117. //new_AI << AI, (AeqI.array()+n).matrix(), AeqJ;
  118. //new_AJ << AJ, AeqJ, (AeqI.array()+n).matrix();
  119. //new_AV << AV, AeqV, AeqV;
  120. //// Slow
  121. //igl::sparse(new_AI,new_AJ,new_AV,n+neq,n+neq,new_A);
  122. //cout<<matlab_format(new_A,"new_A")<<endl;
  123. // precompute RHS builders
  124. if(kr > 0)
  125. {
  126. SparseMatrix<T> Aulk,Akul;
  127. // Slow
  128. igl::slice(new_A,data.unknown_lagrange,data.known,Aulk);
  129. //// This doesn't work!!!
  130. //data.preY = Aulk + Akul.transpose();
  131. // Slow
  132. if(data.Auu_sym)
  133. {
  134. data.preY = Aulk*2;
  135. }else
  136. {
  137. igl::slice(new_A,data.known,data.unknown_lagrange,Akul);
  138. SparseMatrix<T> AkulT = Akul.transpose();
  139. data.preY = Aulk + AkulT;
  140. }
  141. }else
  142. {
  143. data.preY.resize(data.unknown_lagrange.size(),0);
  144. }
  145. // Positive definite and no equality constraints (Postive definiteness
  146. // implies symmetric)
  147. if(data.Auu_pd && neq == 0)
  148. {
  149. data.llt.compute(Auu);
  150. switch(data.ldlt.info())
  151. {
  152. case Eigen::Success:
  153. break;
  154. case Eigen::NumericalIssue:
  155. cerr<<"Error: Numerical issue."<<endl;
  156. return false;
  157. default:
  158. cerr<<"Error: Other."<<endl;
  159. return false;
  160. }
  161. data.solver_type = igl::min_quad_with_fixed_data<T>::LLT;
  162. }else
  163. {
  164. // Either not PD or there are equality constraints
  165. SparseMatrix<T> NA;
  166. igl::slice(new_A,data.unknown_lagrange,data.unknown_lagrange,NA);
  167. data.NA = NA;
  168. // Ideally we'd use LDLT but Eigen doesn't support positive semi-definite
  169. // matrices:
  170. // http://forum.kde.org/viewtopic.php?f=74&t=106962&p=291990#p291990
  171. if(data.Auu_sym && false)
  172. {
  173. data.ldlt.compute(NA);
  174. switch(data.ldlt.info())
  175. {
  176. case Eigen::Success:
  177. break;
  178. case Eigen::NumericalIssue:
  179. cerr<<"Error: Numerical issue."<<endl;
  180. return false;
  181. default:
  182. cerr<<"Error: Other."<<endl;
  183. return false;
  184. }
  185. data.solver_type = igl::min_quad_with_fixed_data<T>::LDLT;
  186. }else
  187. {
  188. // Resort to LU
  189. // Bottleneck >1/2
  190. data.lu.compute(NA);
  191. MatrixXd NAf;
  192. full(NA,NAf);
  193. //std::cout<<"NA=["<<std::endl<<NA<<std::endl<<"];"<<std::endl;
  194. switch(data.lu.info())
  195. {
  196. case Eigen::Success:
  197. break;
  198. case Eigen::NumericalIssue:
  199. cerr<<"Error: Numerical issue."<<endl;
  200. return false;
  201. case Eigen::InvalidInput:
  202. cerr<<"Error: Invalid Input."<<endl;
  203. return false;
  204. default:
  205. cerr<<"Error: Other."<<endl;
  206. return false;
  207. }
  208. data.solver_type = igl::min_quad_with_fixed_data<T>::LU;
  209. }
  210. }
  211. return true;
  212. }
  213. template <typename T, typename DerivedY, typename DerivedZ>
  214. IGL_INLINE bool igl::min_quad_with_fixed_solve(
  215. const igl::min_quad_with_fixed_data<T> & data,
  216. const Eigen::Matrix<T,Eigen::Dynamic,1> & B,
  217. const Eigen::PlainObjectBase<DerivedY> & Y,
  218. const Eigen::Matrix<T,Eigen::Dynamic,1> & Beq,
  219. Eigen::PlainObjectBase<DerivedZ> & Z)
  220. {
  221. using namespace std;
  222. // number of known rows
  223. int kr = data.known.size();
  224. if(kr!=0)
  225. {
  226. assert(kr == Y.rows());
  227. }
  228. // number of columns to solve
  229. int cols = Y.cols();
  230. // number of lagrange multipliers aka linear equality constraints
  231. int neq = data.lagrange.size();
  232. // append lagrange multiplier rhs's
  233. Eigen::Matrix<T,Eigen::Dynamic,1> BBeq(B.size() + Beq.size());
  234. BBeq << B, (Beq*-2.0);
  235. // Build right hand side
  236. Eigen::Matrix<T,Eigen::Dynamic,1> BBequl;
  237. igl::slice(BBeq,data.unknown_lagrange,BBequl);
  238. Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> BBequlcols;
  239. igl::repmat(BBequl,1,cols,BBequlcols);
  240. Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> NB;
  241. if(kr == 0)
  242. {
  243. NB = BBequlcols;
  244. }else
  245. {
  246. NB = data.preY * Y + BBequlcols;
  247. }
  248. // resize output
  249. Z.resize(data.n,cols);
  250. // Set known values
  251. for(int i = 0;i < kr;i++)
  252. {
  253. for(int j = 0;j < cols;j++)
  254. {
  255. Z(data.known(i),j) = Y(i,j);
  256. }
  257. }
  258. //std::cout<<"NB=["<<std::endl<<NB<<std::endl<<"];"<<std::endl;
  259. Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> sol;
  260. //cout<<matlab_format(NB,"NB")<<endl;
  261. switch(data.solver_type)
  262. {
  263. case igl::min_quad_with_fixed_data<T>::LLT:
  264. sol = data.llt.solve(NB);
  265. break;
  266. case igl::min_quad_with_fixed_data<T>::LDLT:
  267. sol = data.ldlt.solve(NB);
  268. break;
  269. case igl::min_quad_with_fixed_data<T>::LU:
  270. // Not a bottleneck
  271. sol = data.lu.solve(NB);
  272. break;
  273. default:
  274. cerr<<"Error: invalid solver type"<<endl;
  275. return false;
  276. }
  277. //std::cout<<"sol=["<<std::endl<<sol<<std::endl<<"];"<<std::endl;
  278. // Now sol contains sol/-0.5
  279. sol *= -0.5;
  280. // Now sol contains solution
  281. // Place solution in Z
  282. for(int i = 0;i<(sol.rows()-neq);i++)
  283. {
  284. for(int j = 0;j<sol.cols();j++)
  285. {
  286. Z(data.unknown_lagrange(i),j) = sol(i,j);
  287. }
  288. }
  289. return true;
  290. }
  291. #ifndef IGL_HEADER_ONLY
  292. // Explicit template specialization
  293. 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> >(igl::min_quad_with_fixed_data<double> const&, Eigen::Matrix<double, -1, 1, 0, -1, 1> const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::Matrix<double, -1, 1, 0, -1, 1> const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  294. template bool igl::min_quad_with_fixed_precompute<double>(Eigen::SparseMatrix<double, 0, int> const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, Eigen::SparseMatrix<double, 0, int> const&, bool, igl::min_quad_with_fixed_data<double>&);
  295. #endif