min_quad_with_fixed.cpp 16 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530
  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, typename Derivedknown>
  20. IGL_INLINE bool igl::min_quad_with_fixed_precompute(
  21. const Eigen::SparseMatrix<T>& A,
  22. const Eigen::PlainObjectBase<Derivedknown> & known,
  23. const Eigen::SparseMatrix<T>& Aeq,
  24. const bool pd,
  25. min_quad_with_fixed_data<T> & data
  26. )
  27. {
  28. #define MIN_QUAD_WITH_FIXED_CPP_DEBUG
  29. using namespace Eigen;
  30. using namespace std;
  31. using namespace igl;
  32. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  33. cout<<" pre"<<endl;
  34. #endif
  35. // number of rows
  36. int n = A.rows();
  37. // cache problem size
  38. data.n = n;
  39. int neq = Aeq.rows();
  40. // default is to have 0 linear equality constraints
  41. if(Aeq.size() != 0)
  42. {
  43. assert(n == Aeq.cols());
  44. }
  45. assert(A.rows() == n);
  46. assert(A.cols() == n);
  47. // number of known rows
  48. int kr = known.size();
  49. assert(kr == 0 || known.minCoeff() >= 0);
  50. assert(kr == 0 || known.maxCoeff() < n);
  51. assert(neq <= n);
  52. // cache known
  53. data.known = known;
  54. // get list of unknown indices
  55. data.unknown.resize(n-kr);
  56. std::vector<bool> unknown_mask;
  57. unknown_mask.resize(n,true);
  58. for(int i = 0;i<kr;i++)
  59. {
  60. unknown_mask[known(i)] = false;
  61. }
  62. int u = 0;
  63. for(int i = 0;i<n;i++)
  64. {
  65. if(unknown_mask[i])
  66. {
  67. data.unknown(u) = i;
  68. u++;
  69. }
  70. }
  71. // get list of lagrange multiplier indices
  72. data.lagrange.resize(neq);
  73. for(int i = 0;i<neq;i++)
  74. {
  75. data.lagrange(i) = n + i;
  76. }
  77. // cache unknown followed by lagrange indices
  78. data.unknown_lagrange.resize(data.unknown.size()+data.lagrange.size());
  79. data.unknown_lagrange << data.unknown, data.lagrange;
  80. SparseMatrix<T> Auu;
  81. slice(A,data.unknown,data.unknown,Auu);
  82. // Positive definiteness is *not* determined, rather it is given as a
  83. // parameter
  84. data.Auu_pd = pd;
  85. if(data.Auu_pd)
  86. {
  87. // PD implies symmetric
  88. data.Auu_sym = true;
  89. // This is an annoying assertion unless EPS can be chosen in a nicer way.
  90. //assert(is_symmetric(Auu,EPS<double>()));
  91. assert(is_symmetric(Auu,1.0));
  92. }else
  93. {
  94. // determine if A(unknown,unknown) is symmetric and/or positive definite
  95. data.Auu_sym = is_symmetric(Auu,EPS<double>());
  96. }
  97. // Determine number of linearly independent constraints
  98. int nc = 0;
  99. if(neq>0)
  100. {
  101. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  102. cout<<" qr"<<endl;
  103. #endif
  104. // QR decomposition to determine row rank in Aequ
  105. slice(Aeq,data.unknown,2,data.Aequ);
  106. assert(data.Aequ.rows() == neq);
  107. assert(data.Aequ.cols() == data.unknown.size());
  108. data.AeqTQR.compute(data.Aequ.transpose().eval());
  109. cout<<endl<<matlab_format(SparseMatrix<T>(data.Aequ.transpose().eval()),"AeqT")<<endl<<endl;
  110. switch(data.AeqTQR.info())
  111. {
  112. case Eigen::Success:
  113. break;
  114. case Eigen::NumericalIssue:
  115. cerr<<"Error: Numerical issue."<<endl;
  116. return false;
  117. case Eigen::InvalidInput:
  118. cerr<<"Error: Invalid input."<<endl;
  119. return false;
  120. default:
  121. cerr<<"Error: Other."<<endl;
  122. return false;
  123. }
  124. nc = data.AeqTQR.rank();
  125. assert(nc<=neq);
  126. data.Aeq_li = nc == neq;
  127. }else
  128. {
  129. data.Aeq_li = true;
  130. }
  131. if(data.Aeq_li)
  132. {
  133. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  134. cout<<" Aeq_li=true"<<endl;
  135. #endif
  136. // Append lagrange multiplier quadratic terms
  137. SparseMatrix<T> new_A;
  138. SparseMatrix<T> AeqT = Aeq.transpose();
  139. SparseMatrix<T> Z(neq,neq);
  140. // This is a bit slower. But why isn't cat fast?
  141. new_A = cat(1, cat(2, A, AeqT ),
  142. cat(2, Aeq, Z ));
  143. // precompute RHS builders
  144. if(kr > 0)
  145. {
  146. SparseMatrix<T> Aulk,Akul;
  147. // Slow
  148. slice(new_A,data.unknown_lagrange,data.known,Aulk);
  149. //// This doesn't work!!!
  150. //data.preY = Aulk + Akul.transpose();
  151. // Slow
  152. if(data.Auu_sym)
  153. {
  154. data.preY = Aulk*2;
  155. }else
  156. {
  157. slice(new_A,data.known,data.unknown_lagrange,Akul);
  158. SparseMatrix<T> AkulT = Akul.transpose();
  159. data.preY = Aulk + AkulT;
  160. }
  161. }else
  162. {
  163. data.preY.resize(data.unknown_lagrange.size(),0);
  164. }
  165. // Positive definite and no equality constraints (Postive definiteness
  166. // implies symmetric)
  167. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  168. cout<<" factorize"<<endl;
  169. #endif
  170. if(data.Auu_pd && neq == 0)
  171. {
  172. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  173. cout<<" llt"<<endl;
  174. #endif
  175. data.llt.compute(Auu);
  176. switch(data.llt.info())
  177. {
  178. case Eigen::Success:
  179. break;
  180. case Eigen::NumericalIssue:
  181. cerr<<"Error: Numerical issue."<<endl;
  182. return false;
  183. default:
  184. cerr<<"Error: Other."<<endl;
  185. return false;
  186. }
  187. data.solver_type = min_quad_with_fixed_data<T>::LLT;
  188. }else
  189. {
  190. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  191. cout<<" ldlt"<<endl;
  192. #endif
  193. // Either not PD or there are equality constraints
  194. SparseMatrix<T> NA;
  195. slice(new_A,data.unknown_lagrange,data.unknown_lagrange,NA);
  196. data.NA = NA;
  197. // Ideally we'd use LDLT but Eigen doesn't support positive semi-definite
  198. // matrices:
  199. // http://forum.kde.org/viewtopic.php?f=74&t=106962&p=291990#p291990
  200. if(data.Auu_sym && false)
  201. {
  202. data.ldlt.compute(NA);
  203. switch(data.ldlt.info())
  204. {
  205. case Eigen::Success:
  206. break;
  207. case Eigen::NumericalIssue:
  208. cerr<<"Error: Numerical issue."<<endl;
  209. return false;
  210. default:
  211. cerr<<"Error: Other."<<endl;
  212. return false;
  213. }
  214. data.solver_type = min_quad_with_fixed_data<T>::LDLT;
  215. }else
  216. {
  217. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  218. cout<<" lu"<<endl;
  219. #endif
  220. // Resort to LU
  221. // Bottleneck >1/2
  222. data.lu.compute(NA);
  223. //std::cout<<"NA=["<<std::endl<<NA<<std::endl<<"];"<<std::endl;
  224. switch(data.lu.info())
  225. {
  226. case Eigen::Success:
  227. break;
  228. case Eigen::NumericalIssue:
  229. cerr<<"Error: Numerical issue."<<endl;
  230. return false;
  231. case Eigen::InvalidInput:
  232. cerr<<"Error: Invalid Input."<<endl;
  233. return false;
  234. default:
  235. cerr<<"Error: Other."<<endl;
  236. return false;
  237. }
  238. data.solver_type = min_quad_with_fixed_data<T>::LU;
  239. }
  240. }
  241. }else
  242. {
  243. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  244. cout<<" Aeq_li=false"<<endl;
  245. #endif
  246. data.neq = neq;
  247. const int nu = data.unknown.size();
  248. //cout<<"nu: "<<nu<<endl;
  249. //cout<<"neq: "<<neq<<endl;
  250. //cout<<"nc: "<<nc<<endl;
  251. //cout<<" matrixR"<<endl;
  252. SparseMatrix<T> AeqTR,AeqTQ;
  253. AeqTR = data.AeqTQR.matrixR();
  254. // This shouldn't be necessary
  255. AeqTR.prune(0.0);
  256. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  257. cout<<" matrixQ"<<endl;
  258. #endif
  259. // THIS IS ESSENTIALLY DENSE AND THIS IS BY FAR THE BOTTLENECK
  260. // http://forum.kde.org/viewtopic.php?f=74&t=117500
  261. AeqTQ = data.AeqTQR.matrixQ();
  262. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  263. cout<<" prune"<<endl;
  264. cout<<" nnz: "<<AeqTQ.nonZeros()<<endl;
  265. #endif
  266. // This shouldn't be necessary
  267. AeqTQ.prune(0.0);
  268. //cout<<"AeqTQ: "<<AeqTQ.rows()<<" "<<AeqTQ.cols()<<endl;
  269. //cout<<matlab_format(AeqTQ,"AeqTQ")<<endl;
  270. //cout<<" perms"<<endl;
  271. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  272. cout<<" nnz: "<<AeqTQ.nonZeros()<<endl;
  273. cout<<" perm"<<endl;
  274. #endif
  275. SparseMatrix<double> I(neq,neq);
  276. I.setIdentity();
  277. data.AeqTE = data.AeqTQR.colsPermutation() * I;
  278. data.AeqTET = data.AeqTQR.colsPermutation().transpose() * I;
  279. assert(AeqTR.rows() == neq);
  280. assert(AeqTQ.rows() == nu);
  281. assert(AeqTQ.cols() == nu);
  282. assert(AeqTR.cols() == neq);
  283. //cout<<" slice"<<endl;
  284. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  285. cout<<" slice"<<endl;
  286. #endif
  287. data.AeqTQ1 = AeqTQ.topLeftCorner(nu,nc);
  288. data.AeqTQ1T = data.AeqTQ1.transpose().eval();
  289. // ALREADY TRIM (Not 100% sure about this)
  290. data.AeqTR1 = AeqTR.topLeftCorner(nc,nc);
  291. data.AeqTR1T = data.AeqTR1.transpose().eval();
  292. //cout<<"AeqTR1T.size() "<<data.AeqTR1T.rows()<<" "<<data.AeqTR1T.cols()<<endl;
  293. // Null space
  294. data.AeqTQ2 = AeqTQ.bottomRightCorner(nu,nu-nc);
  295. data.AeqTQ2T = data.AeqTQ2.transpose().eval();
  296. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  297. cout<<" proj"<<endl;
  298. #endif
  299. // Projected hessian
  300. SparseMatrix<T> QRAuu = data.AeqTQ2T * Auu * data.AeqTQ2;
  301. {
  302. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  303. cout<<" factorize"<<endl;
  304. #endif
  305. // QRAuu should always be PD
  306. data.llt.compute(QRAuu);
  307. switch(data.llt.info())
  308. {
  309. case Eigen::Success:
  310. break;
  311. case Eigen::NumericalIssue:
  312. cerr<<"Error: Numerical issue."<<endl;
  313. return false;
  314. default:
  315. cerr<<"Error: Other."<<endl;
  316. return false;
  317. }
  318. data.solver_type = min_quad_with_fixed_data<T>::QR_LLT;
  319. }
  320. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  321. cout<<" smash"<<endl;
  322. #endif
  323. // Known value multiplier
  324. SparseMatrix<T> Auk;
  325. slice(A,data.unknown,data.known,Auk);
  326. SparseMatrix<T> Aku;
  327. slice(A,data.known,data.unknown,Aku);
  328. SparseMatrix<T> AkuT = Aku.transpose();
  329. data.preY = Auk + AkuT;
  330. // Needed during solve
  331. data.Auu = Auu;
  332. slice(Aeq,data.known,2,data.Aeqk);
  333. assert(data.Aeqk.rows() == neq);
  334. assert(data.Aeqk.cols() == data.known.size());
  335. }
  336. return true;
  337. }
  338. template <
  339. typename T,
  340. typename DerivedB,
  341. typename DerivedY,
  342. typename DerivedBeq,
  343. typename DerivedZ,
  344. typename Derivedsol>
  345. IGL_INLINE bool igl::min_quad_with_fixed_solve(
  346. const min_quad_with_fixed_data<T> & data,
  347. const Eigen::PlainObjectBase<DerivedB> & B,
  348. const Eigen::PlainObjectBase<DerivedY> & Y,
  349. const Eigen::PlainObjectBase<DerivedBeq> & Beq,
  350. Eigen::PlainObjectBase<DerivedZ> & Z,
  351. Eigen::PlainObjectBase<Derivedsol> & sol)
  352. {
  353. using namespace std;
  354. using namespace Eigen;
  355. using namespace igl;
  356. // number of known rows
  357. int kr = data.known.size();
  358. if(kr!=0)
  359. {
  360. assert(kr == Y.rows());
  361. }
  362. // number of columns to solve
  363. int cols = Y.cols();
  364. assert(B.cols() == 1);
  365. assert(Beq.size() == 0 || Beq.cols() == 1);
  366. // resize output
  367. Z.resize(data.n,cols);
  368. // Set known values
  369. for(int i = 0;i < kr;i++)
  370. {
  371. for(int j = 0;j < cols;j++)
  372. {
  373. Z(data.known(i),j) = Y(i,j);
  374. }
  375. }
  376. if(data.Aeq_li)
  377. {
  378. // number of lagrange multipliers aka linear equality constraints
  379. int neq = data.lagrange.size();
  380. // append lagrange multiplier rhs's
  381. Eigen::Matrix<T,Eigen::Dynamic,1> BBeq(B.size() + Beq.size());
  382. BBeq << B, (Beq*-2.0);
  383. // Build right hand side
  384. Eigen::Matrix<T,Eigen::Dynamic,1> BBequl;
  385. igl::slice(BBeq,data.unknown_lagrange,BBequl);
  386. Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> BBequlcols;
  387. igl::repmat(BBequl,1,cols,BBequlcols);
  388. Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> NB;
  389. if(kr == 0)
  390. {
  391. NB = BBequlcols;
  392. }else
  393. {
  394. NB = data.preY * Y + BBequlcols;
  395. }
  396. //std::cout<<"NB=["<<std::endl<<NB<<std::endl<<"];"<<std::endl;
  397. //cout<<matlab_format(NB,"NB")<<endl;
  398. switch(data.solver_type)
  399. {
  400. case igl::min_quad_with_fixed_data<T>::LLT:
  401. sol = data.llt.solve(NB);
  402. break;
  403. case igl::min_quad_with_fixed_data<T>::LDLT:
  404. sol = data.ldlt.solve(NB);
  405. break;
  406. case igl::min_quad_with_fixed_data<T>::LU:
  407. // Not a bottleneck
  408. sol = data.lu.solve(NB);
  409. break;
  410. default:
  411. cerr<<"Error: invalid solver type"<<endl;
  412. return false;
  413. }
  414. //std::cout<<"sol=["<<std::endl<<sol<<std::endl<<"];"<<std::endl;
  415. // Now sol contains sol/-0.5
  416. sol *= -0.5;
  417. // Now sol contains solution
  418. // Place solution in Z
  419. for(int i = 0;i<(sol.rows()-neq);i++)
  420. {
  421. for(int j = 0;j<sol.cols();j++)
  422. {
  423. Z(data.unknown_lagrange(i),j) = sol(i,j);
  424. }
  425. }
  426. }else
  427. {
  428. assert(data.solver_type == min_quad_with_fixed_data<T>::QR_LLT);
  429. PlainObjectBase<DerivedBeq> eff_Beq;
  430. // Adjust Aeq rhs to include known parts
  431. eff_Beq =
  432. //data.AeqTQR.colsPermutation().transpose() * (-data.Aeqk * Y + Beq);
  433. data.AeqTET * (-data.Aeqk * Y + Beq);
  434. // Where did this -0.5 come from? Probably the same place as above.
  435. PlainObjectBase<DerivedB> Bu;
  436. slice(B,data.unknown,Bu);
  437. PlainObjectBase<DerivedB> NB;
  438. NB = -0.5*(Bu + data.preY * Y);
  439. // Trim eff_Beq
  440. const int nc = data.AeqTQR.rank();
  441. const int neq = Beq.rows();
  442. eff_Beq = eff_Beq.topLeftCorner(nc,1).eval();
  443. data.AeqTR1T.template triangularView<Lower>().solveInPlace(eff_Beq);
  444. // Now eff_Beq = (data.AeqTR1T \ (data.AeqTET * (-data.Aeqk * Y + Beq)))
  445. PlainObjectBase<DerivedB> lambda_0;
  446. lambda_0 = data.AeqTQ1 * eff_Beq;
  447. //cout<<matlab_format(lambda_0,"lambda_0")<<endl;
  448. PlainObjectBase<DerivedB> QRB;
  449. QRB = -data.AeqTQ2T * (data.Auu * lambda_0) + data.AeqTQ2T * NB;
  450. PlainObjectBase<Derivedsol> lambda;
  451. lambda = data.llt.solve(QRB);
  452. // prepare output
  453. PlainObjectBase<Derivedsol> solu;
  454. solu = data.AeqTQ2 * lambda + lambda_0;
  455. // http://www.math.uh.edu/~rohop/fall_06/Chapter3.pdf
  456. PlainObjectBase<Derivedsol> solLambda;
  457. {
  458. PlainObjectBase<Derivedsol> temp1,temp2;
  459. temp1 = (data.AeqTQ1T * NB - data.AeqTQ1T * data.Auu * solu);
  460. data.AeqTR1.template triangularView<Upper>().solveInPlace(temp1);
  461. //cout<<matlab_format(temp1,"temp1")<<endl;
  462. temp2 = PlainObjectBase<Derivedsol>::Zero(neq,1);
  463. temp2.topLeftCorner(nc,1) = temp1;
  464. //solLambda = data.AeqTQR.colsPermutation() * temp2;
  465. solLambda = data.AeqTE * temp2;
  466. }
  467. // sol is [Z(unknown);Lambda]
  468. assert(data.unknown.size() == solu.rows());
  469. assert(cols == solu.cols());
  470. assert(data.neq == neq);
  471. assert(data.neq == solLambda.rows());
  472. assert(cols == solLambda.cols());
  473. sol.resize(data.unknown.size()+data.neq,cols);
  474. sol.block(0,0,solu.rows(),solu.cols()) = solu;
  475. sol.block(solu.rows(),0,solLambda.rows(),solLambda.cols()) = solLambda;
  476. for(int u = 0;u<data.unknown.size();u++)
  477. {
  478. for(int j = 0;j<Z.cols();j++)
  479. {
  480. Z(data.unknown(u),j) = solu(u,j);
  481. }
  482. }
  483. }
  484. return true;
  485. }
  486. template <
  487. typename T,
  488. typename DerivedB,
  489. typename DerivedY,
  490. typename DerivedBeq,
  491. typename DerivedZ>
  492. IGL_INLINE bool igl::min_quad_with_fixed_solve(
  493. const min_quad_with_fixed_data<T> & data,
  494. const Eigen::PlainObjectBase<DerivedB> & B,
  495. const Eigen::PlainObjectBase<DerivedY> & Y,
  496. const Eigen::PlainObjectBase<DerivedBeq> & Beq,
  497. Eigen::PlainObjectBase<DerivedZ> & Z)
  498. {
  499. Eigen::PlainObjectBase<DerivedZ> sol;
  500. return min_quad_with_fixed_solve(data,B,Y,Beq,Z,sol);
  501. }
  502. #ifndef IGL_HEADER_ONLY
  503. // Explicit template specialization
  504. 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> >&);
  505. 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>&);
  506. 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> >&);
  507. #endif