libigl v2.5.0
Loading...
Searching...
No Matches
min_quad_with_fixed.impl.h
Go to the documentation of this file.
1// This file is part of libigl, a simple c++ geometry processing library.
2//
3// Copyright (C) 2016 Alec Jacobson <alecjacobson@gmail.com>
4//
5// This Source Code Form is subject to the terms of the Mozilla Public License
6// v. 2.0. If a copy of the MPL was not distributed with this file, You can
7// obtain one at http://mozilla.org/MPL/2.0/.
8#pragma once
9
10#include "min_quad_with_fixed.h"
11
12#include "slice.h"
13#include "is_symmetric.h"
14#include "find.h"
15#include "sparse.h"
16#include "repmat.h"
17#include "EPS.h"
18#include "cat.h"
19
20//#include <Eigen/SparseExtra>
21// Bug in unsupported/Eigen/SparseExtra needs iostream first
22#include <iostream>
23#include <unsupported/Eigen/SparseExtra>
24#include <cassert>
25#include <cstdio>
26#include "matlab_format.h"
27#include <type_traits>
28
29template <typename T, typename Derivedknown>
31 const Eigen::SparseMatrix<T>& A2,
32 const Eigen::MatrixBase<Derivedknown> & known,
33 const Eigen::SparseMatrix<T>& Aeq,
34 const bool pd,
36 )
37{
38//#define MIN_QUAD_WITH_FIXED_CPP_DEBUG
39 using namespace Eigen;
40 using namespace std;
41 const Eigen::SparseMatrix<T> A = 0.5*A2;
42#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
43 cout<<" pre"<<endl;
44#endif
45 // number of rows
46 int n = A.rows();
47 // cache problem size
48 data.n = n;
49
50 int neq = Aeq.rows();
51 // default is to have 0 linear equality constraints
52 if(Aeq.size() != 0)
53 {
54 assert(n == Aeq.cols() && "#Aeq.cols() should match A.rows()");
55 }
56
57 assert(known.cols() == 1 && "known should be a vector");
58 assert(A.rows() == n && "A should be square");
59 assert(A.cols() == n && "A should be square");
60
61 // number of known rows
62 int kr = known.size();
63
64 assert((kr == 0 || known.minCoeff() >= 0)&& "known indices should be in [0,n)");
65 assert((kr == 0 || known.maxCoeff() < n) && "known indices should be in [0,n)");
66 assert(neq <= n && "Number of equality constraints should be less than DOFs");
67
68
69 // cache known
70 // FIXME: This is *NOT* generic and introduces a copy.
71 data.known = known.template cast<int>();
72
73 // get list of unknown indices
74 data.unknown.resize(n-kr);
75 std::vector<bool> unknown_mask;
76 unknown_mask.resize(n,true);
77 for(int i = 0;i<kr;i++)
78 {
79 unknown_mask[known(i, 0)] = false;
80 }
81 int u = 0;
82 for(int i = 0;i<n;i++)
83 {
84 if(unknown_mask[i])
85 {
86 data.unknown(u) = i;
87 u++;
88 }
89 }
90 // get list of lagrange multiplier indices
91 data.lagrange.resize(neq);
92 for(int i = 0;i<neq;i++)
93 {
94 data.lagrange(i) = n + i;
95 }
96 // cache unknown followed by lagrange indices
97 data.unknown_lagrange.resize(data.unknown.size()+data.lagrange.size());
98 // Would like to do:
99 //data.unknown_lagrange << data.unknown, data.lagrange;
100 // but Eigen can't handle empty vectors in comma initialization
101 // https://forum.kde.org/viewtopic.php?f=74&t=107974&p=364947#p364947
102 if(data.unknown.size() > 0)
103 {
104 data.unknown_lagrange.head(data.unknown.size()) = data.unknown;
105 }
106 if(data.lagrange.size() > 0)
107 {
108 data.unknown_lagrange.tail(data.lagrange.size()) = data.lagrange;
109 }
110
111 SparseMatrix<T> Auu;
112 slice(A,data.unknown,data.unknown,Auu);
113 assert(Auu.size() != 0 && Auu.rows() > 0 && "There should be at least one unknown.");
114
115 // Positive definiteness is *not* determined, rather it is given as a
116 // parameter
117 data.Auu_pd = pd;
118 if(data.Auu_pd)
119 {
120 // PD implies symmetric
121 data.Auu_sym = true;
122 // This is an annoying assertion unless EPS can be chosen in a nicer way.
123 //assert(is_symmetric(Auu,EPS<T>()));
124 assert(is_symmetric(Auu,1.0) &&
125 "Auu should be symmetric if positive definite");
126 }else
127 {
128 // determine if A(unknown,unknown) is symmetric and/or positive definite
129 VectorXi AuuI,AuuJ;
130 Matrix<T,Eigen::Dynamic,Eigen::Dynamic> AuuV;
131 find(Auu,AuuI,AuuJ,AuuV);
132 data.Auu_sym = is_symmetric(Auu,EPS<T>()*AuuV.maxCoeff());
133 }
134
135 // Determine number of linearly independent constraints
136 int nc = 0;
137 if(neq>0)
138 {
139#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
140 cout<<" qr"<<endl;
141#endif
142 // QR decomposition to determine row rank in Aequ
143 slice(Aeq,data.unknown,2,data.Aequ);
144 assert(data.Aequ.rows() == neq &&
145 "#Rows in Aequ should match #constraints");
146 assert(data.Aequ.cols() == data.unknown.size() &&
147 "#cols in Aequ should match #unknowns");
148 data.AeqTQR.compute(data.Aequ.transpose().eval());
149#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
150 //cout<<endl<<matlab_format(SparseMatrix<T>(data.Aequ.transpose().eval()),"AeqT")<<endl<<endl;
151#endif
152 switch(data.AeqTQR.info())
153 {
154 case Eigen::Success:
155 break;
156 case Eigen::NumericalIssue:
157 cerr<<"Error: Numerical issue."<<endl;
158 return false;
159 case Eigen::InvalidInput:
160 cerr<<"Error: Invalid input."<<endl;
161 return false;
162 default:
163 cerr<<"Error: Other."<<endl;
164 return false;
165 }
166 nc = data.AeqTQR.rank();
167 assert(nc<=neq &&
168 "Rank of reduced constraints should be <= #original constraints");
169 data.Aeq_li = nc == neq;
170 //cout<<"data.Aeq_li: "<<data.Aeq_li<<endl;
171 }else
172 {
173 data.Aeq_li = true;
174 }
175
176 if(data.Aeq_li)
177 {
178#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
179 cout<<" Aeq_li=true"<<endl;
180#endif
181 // Append lagrange multiplier quadratic terms
182 SparseMatrix<T> new_A;
183 SparseMatrix<T> AeqT = Aeq.transpose();
184 SparseMatrix<T> Z(neq,neq);
185 // This is a bit slower. But why isn't cat fast?
186 new_A = cat(1, cat(2, A, AeqT ),
187 cat(2, Aeq, Z ));
188
189 // precompute RHS builders
190 if(kr > 0)
191 {
192 SparseMatrix<T> Aulk,Akul;
193 // Slow
194 slice(new_A,data.unknown_lagrange,data.known,Aulk);
196 //data.preY = Aulk + Akul.transpose();
197 // Slow
198 if(data.Auu_sym)
199 {
200 data.preY = Aulk*2;
201 }else
202 {
203 slice(new_A,data.known,data.unknown_lagrange,Akul);
204 SparseMatrix<T> AkulT = Akul.transpose();
205 data.preY = Aulk + AkulT;
206 }
207 }else
208 {
209 data.preY.resize(data.unknown_lagrange.size(),0);
210 }
211
212 // Positive definite and no equality constraints (Positive definiteness
213 // implies symmetric)
214#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
215 cout<<" factorize"<<endl;
216#endif
217 if(data.Auu_pd && neq == 0)
218 {
219#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
220 cout<<" llt"<<endl;
221#endif
222 data.llt.compute(Auu);
223 switch(data.llt.info())
224 {
225 case Eigen::Success:
226 break;
227 case Eigen::NumericalIssue:
228 cerr<<"Error: Numerical issue."<<endl;
229 return false;
230 default:
231 cerr<<"Error: Other."<<endl;
232 return false;
233 }
235 }else
236 {
237#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
238 cout<<" ldlt/lu"<<endl;
239#endif
240 // Either not PD or there are equality constraints
241 SparseMatrix<T> NA;
242 slice(new_A,data.unknown_lagrange,data.unknown_lagrange,NA);
243 data.NA = NA;
244 if(data.Auu_pd)
245 {
246#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
247 cout<<" ldlt"<<endl;
248#endif
249 data.ldlt.compute(NA);
250 switch(data.ldlt.info())
251 {
252 case Eigen::Success:
253 break;
254 case Eigen::NumericalIssue:
255 cerr<<"Error: Numerical issue."<<endl;
256 return false;
257 default:
258 cerr<<"Error: Other."<<endl;
259 return false;
260 }
262 }else
263 {
264#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
265 cout<<" lu"<<endl;
266#endif
267 // Resort to LU
268 // Bottleneck >1/2
269 data.lu.compute(NA);
270 //std::cout<<"NA=["<<std::endl<<NA<<std::endl<<"];"<<std::endl;
271 switch(data.lu.info())
272 {
273 case Eigen::Success:
274 break;
275 case Eigen::NumericalIssue:
276 cerr<<"Error: Numerical issue."<<endl;
277 return false;
278 case Eigen::InvalidInput:
279 cerr<<"Error: Invalid Input."<<endl;
280 return false;
281 default:
282 cerr<<"Error: Other."<<endl;
283 return false;
284 }
286 }
287 }
288 }else
289 {
290#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
291 cout<<" Aeq_li=false"<<endl;
292#endif
293 data.neq = neq;
294 const int nu = data.unknown.size();
295 //cout<<"nu: "<<nu<<endl;
296 //cout<<"neq: "<<neq<<endl;
297 //cout<<"nc: "<<nc<<endl;
298 //cout<<" matrixR"<<endl;
299 SparseMatrix<T> AeqTR,AeqTQ;
300 AeqTR = data.AeqTQR.matrixR();
301 // This shouldn't be necessary
302 AeqTR.prune(static_cast<T>(0.0));
303#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
304 cout<<" matrixQ"<<endl;
305#endif
306 // THIS IS ESSENTIALLY DENSE AND THIS IS BY FAR THE BOTTLENECK
307 // http://forum.kde.org/viewtopic.php?f=74&t=117500
308 AeqTQ = data.AeqTQR.matrixQ();
309#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
310 cout<<" prune"<<endl;
311 cout<<" nnz: "<<AeqTQ.nonZeros()<<endl;
312#endif
313 // This shouldn't be necessary
314 AeqTQ.prune(static_cast<T>(0.0));
315 //cout<<"AeqTQ: "<<AeqTQ.rows()<<" "<<AeqTQ.cols()<<endl;
316 //cout<<matlab_format(AeqTQ,"AeqTQ")<<endl;
317 //cout<<" perms"<<endl;
318#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
319 cout<<" nnz: "<<AeqTQ.nonZeros()<<endl;
320 cout<<" perm"<<endl;
321#endif
322 SparseMatrix<T> I(neq,neq);
323 I.setIdentity();
324 data.AeqTE = data.AeqTQR.colsPermutation() * I;
325 data.AeqTET = data.AeqTQR.colsPermutation().transpose() * I;
326 assert(AeqTR.rows() == nu && "#rows in AeqTR should match #unknowns");
327 assert(AeqTR.cols() == neq && "#cols in AeqTR should match #constraints");
328 assert(AeqTQ.rows() == nu && "#rows in AeqTQ should match #unknowns");
329 assert(AeqTQ.cols() == nu && "#cols in AeqTQ should match #unknowns");
330 //cout<<" slice"<<endl;
331#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
332 cout<<" slice"<<endl;
333#endif
334 data.AeqTQ1 = AeqTQ.topLeftCorner(nu,nc);
335 data.AeqTQ1T = data.AeqTQ1.transpose().eval();
336 // ALREADY TRIM (Not 100% sure about this)
337 data.AeqTR1 = AeqTR.topLeftCorner(nc,nc);
338 data.AeqTR1T = data.AeqTR1.transpose().eval();
339 //cout<<"AeqTR1T.size() "<<data.AeqTR1T.rows()<<" "<<data.AeqTR1T.cols()<<endl;
340 // Null space
341 data.AeqTQ2 = AeqTQ.bottomRightCorner(nu,nu-nc);
342 data.AeqTQ2T = data.AeqTQ2.transpose().eval();
343#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
344 cout<<" proj"<<endl;
345#endif
346 // Projected hessian
347 SparseMatrix<T> QRAuu = data.AeqTQ2T * Auu * data.AeqTQ2;
348 {
349#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
350 cout<<" factorize"<<endl;
351#endif
352 // QRAuu should always be PD
353 data.llt.compute(QRAuu);
354 switch(data.llt.info())
355 {
356 case Eigen::Success:
357 break;
358 case Eigen::NumericalIssue:
359 cerr<<"Error: Numerical issue."<<endl;
360 return false;
361 default:
362 cerr<<"Error: Other."<<endl;
363 return false;
364 }
366 }
367#ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
368 cout<<" smash"<<endl;
369#endif
370 // Known value multiplier
371 SparseMatrix<T> Auk;
372 slice(A,data.unknown,data.known,Auk);
373 SparseMatrix<T> Aku;
374 slice(A,data.known,data.unknown,Aku);
375 SparseMatrix<T> AkuT = Aku.transpose();
376 data.preY = Auk + AkuT;
377 // Needed during solve
378 data.Auu = Auu;
379 slice(Aeq,data.known,2,data.Aeqk);
380 assert(data.Aeqk.rows() == neq);
381 assert(data.Aeqk.cols() == data.known.size());
382 }
383 return true;
384}
385
386
387template <
388 typename T,
389 typename DerivedB,
390 typename DerivedY,
391 typename DerivedBeq,
392 typename DerivedZ,
393 typename Derivedsol>
395 const min_quad_with_fixed_data<T> & data,
396 const Eigen::MatrixBase<DerivedB> & B,
397 const Eigen::MatrixBase<DerivedY> & Y,
398 const Eigen::MatrixBase<DerivedBeq> & Beq,
399 Eigen::PlainObjectBase<DerivedZ> & Z,
400 Eigen::PlainObjectBase<Derivedsol> & sol)
401{
402 using namespace std;
403 using namespace Eigen;
404 typedef Matrix<T,Dynamic,Dynamic> MatrixXT;
405 // number of known rows
406 int kr = data.known.size();
407 if(kr!=0)
408 {
409 assert(kr == Y.rows());
410 }
411 // number of columns to solve
412 int cols = Y.cols();
413 assert(B.cols() == 1 || B.cols() == cols);
414 assert(Beq.size() == 0 || Beq.cols() == 1 || Beq.cols() == cols);
415
416 // resize output
417 Z.resize(data.n,cols);
418 // Set known values
419 for(int i = 0;i < kr;i++)
420 {
421 for(int j = 0;j < cols;j++)
422 {
423 Z(data.known(i),j) = Y(i,j);
424 }
425 }
426
427 if(data.Aeq_li)
428 {
429 // number of lagrange multipliers aka linear equality constraints
430 int neq = data.lagrange.size();
431 // append lagrange multiplier rhs's
432 MatrixXT BBeq(B.rows() + Beq.rows(),cols);
433 if(B.size() > 0)
434 {
435 BBeq.topLeftCorner(B.rows(),cols) = B.replicate(1,B.cols()==cols?1:cols);
436 }
437 if(Beq.size() > 0)
438 {
439 BBeq.bottomLeftCorner(Beq.rows(),cols) = -2.0*Beq.replicate(1,Beq.cols()==cols?1:cols);
440 }
441
442 // Build right hand side
443 MatrixXT BBequlcols = BBeq(data.unknown_lagrange,Eigen::all);
444 MatrixXT NB;
445 if(kr == 0)
446 {
447 NB = BBequlcols;
448 }else
449 {
450 NB = data.preY * Y + BBequlcols;
451 }
452
453 //std::cout<<"NB=["<<std::endl<<NB<<std::endl<<"];"<<std::endl;
454 //cout<<matlab_format(NB,"NB")<<endl;
455 switch(data.solver_type)
456 {
458 sol = data.llt.solve(NB);
459 break;
461 sol = data.ldlt.solve(NB);
462 break;
464 // Not a bottleneck
465 sol = data.lu.solve(NB);
466 break;
467 default:
468 cerr<<"Error: invalid solver type"<<endl;
469 return false;
470 }
471 //std::cout<<"sol=["<<std::endl<<sol<<std::endl<<"];"<<std::endl;
472 // Now sol contains sol/-0.5
473 sol *= -0.5;
474 // Now sol contains solution
475 // Place solution in Z
476 for(int i = 0;i<(sol.rows()-neq);i++)
477 {
478 for(int j = 0;j<sol.cols();j++)
479 {
480 Z(data.unknown_lagrange(i),j) = sol(i,j);
481 }
482 }
483 }else
484 {
486 MatrixXT eff_Beq;
487 // Adjust Aeq rhs to include known parts
488 eff_Beq =
489 //data.AeqTQR.colsPermutation().transpose() * (-data.Aeqk * Y + Beq);
490 data.AeqTET * (-data.Aeqk * Y + Beq.replicate(1,Beq.cols()==cols?1:cols));
491 // Where did this -0.5 come from? Probably the same place as above.
492 MatrixXT Bu = B(data.unknown,Eigen::all);
493 MatrixXT NB;
494 NB = -0.5*(Bu.replicate(1,B.cols()==cols?1:cols) + data.preY * Y);
495 // Trim eff_Beq
496 const int nc = data.AeqTQR.rank();
497 const int neq = Beq.rows();
498 eff_Beq = eff_Beq.topLeftCorner(nc,cols).eval();
499 data.AeqTR1T.template triangularView<Lower>().solveInPlace(eff_Beq);
500 // Now eff_Beq = (data.AeqTR1T \ (data.AeqTET * (-data.Aeqk * Y + Beq)))
501 MatrixXT lambda_0;
502 lambda_0 = data.AeqTQ1 * eff_Beq;
503 //cout<<matlab_format(lambda_0,"lambda_0")<<endl;
504 MatrixXT QRB;
505 QRB = -data.AeqTQ2T * (data.Auu * lambda_0) + data.AeqTQ2T * NB;
506 Derivedsol lambda;
507 lambda = data.llt.solve(QRB);
508 // prepare output
509 Derivedsol solu;
510 solu = data.AeqTQ2 * lambda + lambda_0;
511 // http://www.math.uh.edu/~rohop/fall_06/Chapter3.pdf
512 Derivedsol solLambda;
513 {
514 Derivedsol temp1,temp2;
515 temp1 = (data.AeqTQ1T * NB - data.AeqTQ1T * data.Auu * solu);
516 data.AeqTR1.template triangularView<Upper>().solveInPlace(temp1);
517 //cout<<matlab_format(temp1,"temp1")<<endl;
518 temp2 = Derivedsol::Zero(neq,cols);
519 temp2.topLeftCorner(nc,cols) = temp1;
520 //solLambda = data.AeqTQR.colsPermutation() * temp2;
521 solLambda = data.AeqTE * temp2;
522 }
523 // sol is [Z(unknown);Lambda]
524 assert(data.unknown.size() == solu.rows());
525 assert(cols == solu.cols());
526 assert(data.neq == neq);
527 assert(data.neq == solLambda.rows());
528 assert(cols == solLambda.cols());
529 sol.resize(data.unknown.size()+data.neq,cols);
530 sol.block(0,0,solu.rows(),solu.cols()) = solu;
531 sol.block(solu.rows(),0,solLambda.rows(),solLambda.cols()) = solLambda;
532 for(int u = 0;u<data.unknown.size();u++)
533 {
534 for(int j = 0;j<Z.cols();j++)
535 {
536 Z(data.unknown(u),j) = solu(u,j);
537 }
538 }
539 }
540 return true;
541}
542
543template <
544 typename T,
545 typename DerivedB,
546 typename DerivedY,
547 typename DerivedBeq,
548 typename DerivedZ>
550 const min_quad_with_fixed_data<T> & data,
551 const Eigen::MatrixBase<DerivedB> & B,
552 const Eigen::MatrixBase<DerivedY> & Y,
553 const Eigen::MatrixBase<DerivedBeq> & Beq,
554 Eigen::PlainObjectBase<DerivedZ> & Z)
555{
556 Eigen::Matrix<typename DerivedZ::Scalar, Eigen::Dynamic, Eigen::Dynamic> sol;
557 return min_quad_with_fixed_solve(data,B,Y,Beq,Z,sol);
558}
559
560template <
561 typename T,
562 typename Derivedknown,
563 typename DerivedB,
564 typename DerivedY,
565 typename DerivedBeq,
566 typename DerivedZ>
568 const Eigen::SparseMatrix<T>& A,
569 const Eigen::MatrixBase<DerivedB> & B,
570 const Eigen::MatrixBase<Derivedknown> & known,
571 const Eigen::MatrixBase<DerivedY> & Y,
572 const Eigen::SparseMatrix<T>& Aeq,
573 const Eigen::MatrixBase<DerivedBeq> & Beq,
574 const bool pd,
575 Eigen::PlainObjectBase<DerivedZ> & Z)
576{
578 if(!min_quad_with_fixed_precompute(A,known,Aeq,pd,data))
579 {
580 return false;
581 }
582 return min_quad_with_fixed_solve(data,B,Y,Beq,Z);
583}
584
585
586template <typename Scalar, int n, int m, bool Hpd>
587IGL_INLINE Eigen::Matrix<Scalar,n,1> igl::min_quad_with_fixed(
588 const Eigen::Matrix<Scalar,n,n> & H,
589 const Eigen::Matrix<Scalar,n,1> & f,
590 const Eigen::Array<bool,n,1> & k,
591 const Eigen::Matrix<Scalar,n,1> & bc,
592 const Eigen::Matrix<Scalar,m,n> & A,
593 const Eigen::Matrix<Scalar,m,1> & b)
594{
595 const auto dyn_n = n == Eigen::Dynamic ? H.rows() : n;
596 const auto dyn_m = m == Eigen::Dynamic ? A.rows() : m;
597 constexpr const int nn = n == Eigen::Dynamic ? Eigen::Dynamic : n+m;
598 const auto dyn_nn = nn == Eigen::Dynamic ? dyn_n+dyn_m : nn;
599 if(dyn_m == 0)
600 {
601 return igl::min_quad_with_fixed<Scalar,n,Hpd>(H,f,k,bc);
602 }
603 // min_x ½ xᵀ H x + xᵀ f subject to A x = b and x(k) = bc(k)
604 // let zᵀ = [xᵀ λᵀ]
605 // min_z ½ zᵀ [H Aᵀ;A 0] z + zᵀ [f;-b] z(k) = bc(k)
606 const auto make_HH = [&]()
607 {
608 // Windows can't remember that nn is const.
609 constexpr const int nn = n == Eigen::Dynamic ? Eigen::Dynamic : n+m;
610 Eigen::Matrix<Scalar,nn,nn> HH =
611 Eigen::Matrix<Scalar,nn,nn>::Zero(dyn_nn,dyn_nn);
612 HH.topLeftCorner(dyn_n,dyn_n) = H;
613 HH.bottomLeftCorner(dyn_m,dyn_n) = A;
614 HH.topRightCorner(dyn_n,dyn_m) = A.transpose();
615 return HH;
616 };
617 const Eigen::Matrix<Scalar,nn,nn> HH = make_HH();
618 const auto make_ff = [&]()
619 {
620 // Windows can't remember that nn is const.
621 constexpr const int nn = n == Eigen::Dynamic ? Eigen::Dynamic : n+m;
622 Eigen::Matrix<Scalar,nn,1> ff(dyn_nn);
623 ff.head(dyn_n) = f;
624 ff.tail(dyn_m) = -b;
625 return ff;
626 };
627 const Eigen::Matrix<Scalar,nn,1> ff = make_ff();
628 const auto make_kk = [&]()
629 {
630 // Windows can't remember that nn is const.
631 constexpr const int nn = n == Eigen::Dynamic ? Eigen::Dynamic : n+m;
632 Eigen::Array<bool,nn,1> kk =
633 Eigen::Array<bool,nn,1>::Constant(dyn_nn,1,false);
634 kk.head(dyn_n) = k;
635 return kk;
636 };
637 const Eigen::Array<bool,nn,1> kk = make_kk();
638 const auto make_bcbc= [&]()
639 {
640 // Windows can't remember that nn is const.
641 constexpr const int nn = n == Eigen::Dynamic ? Eigen::Dynamic : n+m;
642 Eigen::Matrix<Scalar,nn,1> bcbc(dyn_nn);
643 bcbc.head(dyn_n) = bc;
644 return bcbc;
645 };
646 const Eigen::Matrix<Scalar,nn,1> bcbc = make_bcbc();
647 const Eigen::Matrix<Scalar,nn,1> xx =
648 min_quad_with_fixed<Scalar,nn,false>(HH,ff,kk,bcbc);
649 return xx.head(dyn_n);
650}
651
652template <typename Scalar, int n, bool Hpd>
653IGL_INLINE Eigen::Matrix<Scalar,n,1> igl::min_quad_with_fixed(
654 const Eigen::Matrix<Scalar,n,n> & H,
655 const Eigen::Matrix<Scalar,n,1> & f,
656 const Eigen::Array<bool,n,1> & k,
657 const Eigen::Matrix<Scalar,n,1> & bc)
658{
659 assert(H.isApprox(H.transpose(),1e-7));
660 assert(H.rows() == H.cols());
661 assert(H.rows() == f.size());
662 assert(H.rows() == k.size());
663 assert(H.rows() == bc.size());
664 const auto kcount = k.count();
665 // Everything fixed
666 if(kcount == (Eigen::Dynamic?H.rows():n))
667 {
668 return bc;
669 }
670 // Nothing fixed
671 if(kcount == 0)
672 {
673 // avoid function call
674 typedef Eigen::Matrix<Scalar,n,n> MatrixSn;
675 typedef typename
676 std::conditional<Hpd,Eigen::LLT<MatrixSn>,Eigen::CompleteOrthogonalDecomposition<MatrixSn>>::type
677 Solver;
678 return Solver(H).solve(-f);
679 }
680 // All-but-one fixed
681 if( (Eigen::Dynamic?H.rows():n)-kcount == 1)
682 {
683 // which one is not fixed?
684 int u = -1;
685 for(int i=0;i<k.size();i++){ if(!k(i)){ u=i; break; } }
686 assert(u>=0);
687 // min ½ x(u) Huu x(u) + x(u)(fu + H(u,k)bc(k))
688 // Huu x(u) = -(fu + H(u,k) bc(k))
689 // x(u) = (-fu + ∑ -Huj bcj)/Huu
690 Eigen::Matrix<Scalar,n,1> x = bc;
691 x(u) = -f(u);
692 for(int i=0;i<k.size();i++){ if(i!=u){ x(u)-=bc(i)*H(i,u); } }
693 x(u) /= H(u,u);
694 return x;
695 }
696 // Alec: Is there a smart template way to do this?
697 // jdumas: I guess you could do a templated for-loop starting from 16, and
698 // dispatching to the appropriate templated function when the argument matches
699 // (with a fallback to the dynamic version). Cf this example:
700 // https://gist.github.com/disconnect3d/13c2d035bb31b244df14
701 switch(kcount)
702 {
703 case 0: assert(false && "Handled above."); return Eigen::Matrix<Scalar,n,1>();
704 // % Matlibberish for generating these case statements:
705 // maxi=16;for i=1:maxi;fprintf(' case %d:\n {\n const bool D = (n-%d<=0)||(%d>=n)||(n>%d);\n return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:%d,Hpd>(H,f,k,bc);\n }\n',[i i i maxi i]);end
706 case 1:
707 {
708 const bool D = (n-1<=0)||(1>=n)||(n>16);
709 return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:1,Hpd>(H,f,k,bc);
710 }
711 case 2:
712 {
713 const bool D = (n-2<=0)||(2>=n)||(n>16);
714 return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:2,Hpd>(H,f,k,bc);
715 }
716 case 3:
717 {
718 const bool D = (n-3<=0)||(3>=n)||(n>16);
719 return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:3,Hpd>(H,f,k,bc);
720 }
721 case 4:
722 {
723 const bool D = (n-4<=0)||(4>=n)||(n>16);
724 return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:4,Hpd>(H,f,k,bc);
725 }
726 case 5:
727 {
728 const bool D = (n-5<=0)||(5>=n)||(n>16);
729 return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:5,Hpd>(H,f,k,bc);
730 }
731 case 6:
732 {
733 const bool D = (n-6<=0)||(6>=n)||(n>16);
734 return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:6,Hpd>(H,f,k,bc);
735 }
736 case 7:
737 {
738 const bool D = (n-7<=0)||(7>=n)||(n>16);
739 return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:7,Hpd>(H,f,k,bc);
740 }
741 case 8:
742 {
743 const bool D = (n-8<=0)||(8>=n)||(n>16);
744 return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:8,Hpd>(H,f,k,bc);
745 }
746 case 9:
747 {
748 const bool D = (n-9<=0)||(9>=n)||(n>16);
749 return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:9,Hpd>(H,f,k,bc);
750 }
751 case 10:
752 {
753 const bool D = (n-10<=0)||(10>=n)||(n>16);
754 return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:10,Hpd>(H,f,k,bc);
755 }
756 case 11:
757 {
758 const bool D = (n-11<=0)||(11>=n)||(n>16);
759 return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:11,Hpd>(H,f,k,bc);
760 }
761 case 12:
762 {
763 const bool D = (n-12<=0)||(12>=n)||(n>16);
764 return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:12,Hpd>(H,f,k,bc);
765 }
766 case 13:
767 {
768 const bool D = (n-13<=0)||(13>=n)||(n>16);
769 return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:13,Hpd>(H,f,k,bc);
770 }
771 case 14:
772 {
773 const bool D = (n-14<=0)||(14>=n)||(n>16);
774 return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:14,Hpd>(H,f,k,bc);
775 }
776 case 15:
777 {
778 const bool D = (n-15<=0)||(15>=n)||(n>16);
779 return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:15,Hpd>(H,f,k,bc);
780 }
781 case 16:
782 {
783 const bool D = (n-16<=0)||(16>=n)||(n>16);
784 return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:16,Hpd>(H,f,k,bc);
785 }
786 default:
787 return min_quad_with_fixed<Scalar,Eigen::Dynamic,Eigen::Dynamic,Hpd>(H,f,k,bc);
788 }
789}
790
791template <typename Scalar, int n, int kcount, bool Hpd>
792IGL_INLINE Eigen::Matrix<Scalar,n,1> igl::min_quad_with_fixed(
793 const Eigen::Matrix<Scalar,n,n> & H,
794 const Eigen::Matrix<Scalar,n,1> & f,
795 const Eigen::Array<bool,n,1> & k,
796 const Eigen::Matrix<Scalar,n,1> & bc)
797{
798 // 0 and n should be handle outside this function
799 static_assert(kcount==Eigen::Dynamic || kcount>0 ,"");
800 static_assert(kcount==Eigen::Dynamic || kcount<n ,"");
801 const int ucount = n==Eigen::Dynamic ? Eigen::Dynamic : n-kcount;
802 static_assert(kcount==Eigen::Dynamic || ucount+kcount == n ,"");
803 static_assert((n==Eigen::Dynamic) == (ucount==Eigen::Dynamic),"");
804 static_assert((kcount==Eigen::Dynamic) == (ucount==Eigen::Dynamic),"");
805 assert((n==Eigen::Dynamic) || n == H.rows());
806 assert((kcount==Eigen::Dynamic) || kcount == k.count());
807 typedef Eigen::Matrix<Scalar,ucount,ucount> MatrixSuu;
808 typedef Eigen::Matrix<Scalar,ucount,kcount> MatrixSuk;
809 typedef Eigen::Matrix<Scalar,n,1> VectorSn;
810 typedef Eigen::Matrix<Scalar,ucount,1> VectorSu;
811 typedef Eigen::Matrix<Scalar,kcount,1> VectorSk;
812 const auto dyn_n = n==Eigen::Dynamic ? H.rows() : n;
813 const auto dyn_kcount = kcount==Eigen::Dynamic ? k.count() : kcount;
814 const auto dyn_ucount = ucount==Eigen::Dynamic ? dyn_n- dyn_kcount : ucount;
815 // For ucount==2 or kcount==2 this calls the coefficient initiliazer rather
816 // than the size initilizer, but I guess that's ok.
817 MatrixSuu Huu(dyn_ucount,dyn_ucount);
818 MatrixSuk Huk(dyn_ucount,dyn_kcount);
819 VectorSu mrhs(dyn_ucount);
820 VectorSk bck(dyn_kcount);
821 {
822 int ui = 0;
823 int ki = 0;
824 for(int i = 0;i<dyn_n;i++)
825 {
826 if(k(i))
827 {
828 bck(ki) = bc(i);
829 ki++;
830 }else
831 {
832 mrhs(ui) = f(i);
833 int uj = 0;
834 int kj = 0;
835 for(int j = 0;j<dyn_n;j++)
836 {
837 if(k(j))
838 {
839 Huk(ui,kj) = H(i,j);
840 kj++;
841 }else
842 {
843 Huu(ui,uj) = H(i,j);
844 uj++;
845 }
846 }
847 ui++;
848 }
849 }
850 }
851 mrhs += Huk * bck;
852 typedef typename
853 std::conditional<Hpd,
854 Eigen::LLT<MatrixSuu>,
855 // LDLT should be faster for indefinite problems but already found some
856 // cases where it was too inaccurate when called via quadprog_primal.
857 // Ideally this function takes LLT,LDLT, or
858 // CompleteOrthogonalDecomposition as a template parameter. "template
859 // template" parameters did work because LLT,LDLT have different number of
860 // template parameters from CompleteOrthogonalDecomposition. Perhaps
861 // there's a way to take advantage of LLT and LDLT's default template
862 // parameters (I couldn't figure out how).
863 Eigen::CompleteOrthogonalDecomposition<MatrixSuu>>::type
864 Solver;
865 VectorSu xu = Solver(Huu).solve(-mrhs);
866 VectorSn x(dyn_n);
867 {
868 int ui = 0;
869 int ki = 0;
870 for(int i = 0;i<dyn_n;i++)
871 {
872 if(k(i))
873 {
874 x(i) = bck(ki);
875 ki++;
876 }else
877 {
878 x(i) = xu(ui);
879 ui++;
880 }
881 }
882 }
883 return x;
884}
#define IGL_INLINE
Definition igl_inline.h:15
bool is_symmetric(const Eigen::SparseMatrix< AT > &A)
Returns true if the given matrix is symmetric.
void slice(const Eigen::SparseMatrix< TX > &X, const Eigen::DenseBase< DerivedR > &R, const Eigen::DenseBase< DerivedC > &C, Eigen::SparseMatrix< TY > &Y)
Act like the matlab X(row_indices,col_indices) operator, where row_indices, col_indices are non-negat...
bool min_quad_with_fixed(const Eigen::SparseMatrix< T > &A, const Eigen::MatrixBase< DerivedB > &B, const Eigen::MatrixBase< Derivedknown > &known, const Eigen::MatrixBase< DerivedY > &Y, const Eigen::SparseMatrix< T > &Aeq, const Eigen::MatrixBase< DerivedBeq > &Beq, const bool pd, Eigen::PlainObjectBase< DerivedZ > &Z)
This is an overloaded member function, provided for convenience. It differs from the above function o...
Definition min_quad_with_fixed.impl.h:567
void find(const Eigen::SparseMatrix< T > &X, Eigen::DenseBase< DerivedI > &I, Eigen::DenseBase< DerivedJ > &J, Eigen::DenseBase< DerivedV > &V)
Find the non-zero entries and there respective indices in a sparse matrix.
void cat(const int dim, const Eigen::SparseMatrix< Scalar > &A, const Eigen::SparseMatrix< Scalar > &B, Eigen::SparseMatrix< Scalar > &C)
Perform concatenation of a two sparse matrices along a single dimension If dim == 1,...
bool min_quad_with_fixed_solve(const min_quad_with_fixed_data< T > &data, const Eigen::MatrixBase< DerivedB > &B, const Eigen::MatrixBase< DerivedY > &Y, const Eigen::MatrixBase< DerivedBeq > &Beq, Eigen::PlainObjectBase< DerivedZ > &Z, Eigen::PlainObjectBase< Derivedsol > &sol)
Solves a system previously factored using min_quad_with_fixed_precompute.
Definition min_quad_with_fixed.impl.h:394
bool min_quad_with_fixed_precompute(const Eigen::SparseMatrix< T > &A, const Eigen::MatrixBase< Derivedknown > &known, const Eigen::SparseMatrix< T > &Aeq, const bool pd, min_quad_with_fixed_data< T > &data)
Minimize a convex quadratic energy subject to fixed value and linear equality constraints.
Definition min_quad_with_fixed.impl.h:30
Parameters and precomputed values for min_quad_with_fixed.
Definition min_quad_with_fixed.h:166
Eigen::VectorXi known
Indices of known variables.
Definition min_quad_with_fixed.h:174
Eigen::SparseQR< Eigen::SparseMatrix< T >, Eigen::COLAMDOrdering< int > > AeqTQR
Definition min_quad_with_fixed.h:201
Eigen::SimplicialLDLT< Eigen::SparseMatrix< T > > ldlt
Definition min_quad_with_fixed.h:194
Eigen::SparseMatrix< T > preY
Matrix multiplied against Y when constructing right hand side.
Definition min_quad_with_fixed.h:182
Eigen::VectorXi unknown
Indices of unknown variables.
Definition min_quad_with_fixed.h:176
Eigen::SimplicialLLT< Eigen::SparseMatrix< T > > llt
Solver data (factorization)
Definition min_quad_with_fixed.h:193
Eigen::SparseMatrix< T > AeqTQ2
Definition min_quad_with_fixed.h:207
Eigen::SparseMatrix< T > AeqTQ1
Definition min_quad_with_fixed.h:205
Eigen::SparseLU< Eigen::SparseMatrix< T, Eigen::ColMajor >, Eigen::COLAMDOrdering< int > > lu
Definition min_quad_with_fixed.h:195
int n
Size of original system: number of unknowns + number of knowns.
Definition min_quad_with_fixed.h:168
bool Aeq_li
QR factorization Are rows of Aeq linearly independent?
Definition min_quad_with_fixed.h:198
Eigen::VectorXi lagrange
Indices of lagrange variables.
Definition min_quad_with_fixed.h:178
bool Auu_sym
Whether A(unknown,unknown) is symmetric.
Definition min_quad_with_fixed.h:172
int neq
Columns of Aeq corresponding to unknowns.
Definition min_quad_with_fixed.h:200
Eigen::SparseMatrix< T > AeqTR1
Definition min_quad_with_fixed.h:209
Eigen::VectorXi unknown_lagrange
Indices of unknown variable followed by Indices of lagrange variables.
Definition min_quad_with_fixed.h:180
Eigen::SparseMatrix< T > AeqTET
Definition min_quad_with_fixed.h:212
bool Auu_pd
Whether A(unknown,unknown) is positive definite.
Definition min_quad_with_fixed.h:170
Eigen::SparseMatrix< T > AeqTR1T
Definition min_quad_with_fixed.h:210
Eigen::SparseMatrix< T > AeqTQ2T
Definition min_quad_with_fixed.h:208
Eigen::SparseMatrix< T > Aeqk
Definition min_quad_with_fixed.h:202
Eigen::SparseMatrix< T > Auu
Definition min_quad_with_fixed.h:204
enum igl::min_quad_with_fixed_data::SolverType solver_type
Eigen::SparseMatrix< T > Aequ
Definition min_quad_with_fixed.h:203
Eigen::SparseMatrix< T > AeqTE
Definition min_quad_with_fixed.h:211
Eigen::SparseMatrix< T > AeqTQ1T
Definition min_quad_with_fixed.h:206