Updating the changelog in the VERSION file, and version_sync.
[shapes.git] / source / solve_qp_dist_gen.cc
blobcc9576fd9ac79a3f908405daf4193615c735d325
1 /* This file is part of Shapes.
3 * Shapes is free software: you can redistribute it and/or modify
4 * it under the terms of the GNU General Public License as published by
5 * the Free Software Foundation, either version 3 of the License, or
6 * any later version.
8 * Shapes is distributed in the hope that it will be useful,
9 * but WITHOUT ANY WARRANTY; without even the implied warranty of
10 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11 * GNU General Public License for more details.
13 * You should have received a copy of the GNU General Public License
14 * along with Shapes. If not, see <http://www.gnu.org/licenses/>.
16 * Copyright 2009 Henrik Tidefelt
19 #include "quadratic_programs.h"
21 #include <gsl/gsl_linalg.h>
22 #include <gsl/gsl_blas.h>
23 #include <gsl/gsl_cblas.h>
25 using namespace Shapes;
27 #include <vector>
28 #include <cstdlib>
29 #include <cmath>
31 //#define DISPLAY
32 //#define TRACEFIGS // generate output on stdout that gives a graphical picture of the operation
34 #ifdef DISPLAY
35 #include "mathematicarepresentation.h"
36 size_t polytope_distance_generator_form_CallCount = 0;
37 #endif
38 #ifdef TRACEFIGS
39 #include <iomanip>
40 #endif
43 namespace Shapes
45 namespace Computation
47 /* Since an active constraint implies that the corresponding variable is zero, the
48 * subproblems are not concerned with these variables at all.
50 void polytope_distance_generator_form_sub( const size_t dim,
51 const double ** g1Begin, const double ** g1End,
52 const double ** g2Begin, const double ** g2End,
53 const double * yDiff,
54 double sigmaTol,
55 double * d1, double * d2,
56 double * work ); /* allocated memory of size max( dim * ( n1 + n2 - 2 ) + ( n1 + n2 - 2 )^2 + 2 * dim + n1 + n2 )*/
58 void polytope_distance_generator_form_sub_Gilbert( const size_t dim,
59 const double ** g1Begin, const double ** g1End,
60 const double ** g2Begin, const double ** g2End,
61 const double * yDiff,
62 double * d1, double * d2 );
64 void polytope_distance_generator_form_diff( const size_t dim,
65 const double ** g1Begin, const double ** g1End,
66 const double ** g2Begin, const double ** g2End,
67 const double * x1, const double * x2,
68 double * yDiff );
70 void polytope_distance_generator_form_costs( const size_t dim,
71 const double ** g1Begin, const double ** g1End,
72 const double ** g2Begin, const double ** g2End,
73 const double * yDiff,
74 double * primalCost,
75 double * dualCostAccum, /* This value is updated as an accumulated maximum. */
76 double * work ); /* Just needs to be of size (dim) */
78 /* This function will compute orthogonal vectors spanning the nullspaces of vectors of ones.
79 * Memorey will be reused, but never returned.
80 * Given argument n, it returns a pointer to an array of ( n - 1 ) pointers, each pointing to an
81 * array of doubles of length n.
83 const double ** polytope_distance_generator_form_kernel( size_t n );
87 template<class T>
88 class DeleteOnExit
90 T * mem_;
91 public:
92 DeleteOnExit( T * mem ) : mem_( mem ) { }
93 ~DeleteOnExit( ){ delete mem_; }
98 void
99 Computation::polytope_distance_generator_form( const size_t dim,
100 const size_t n1, const double * g1,
101 const size_t n2, const double * g2,
102 double * costLower, double * costUpper,
103 double stopLower, double stopUpper,
104 double gapTolAbs, double gapTolRel,
105 double lmTol,
106 double xTol,
107 double sigmaTol,
108 double * y1, double * y2,
109 QPSolverStatus * status )
111 /* In case we for some reason would return without setting status properly... */
112 status->code = QP_ERROR;
113 status->error = QP_ERROR_UNDEFINED;
115 if( stopLower < 0 )
117 stopLower = std::numeric_limits< double >::max( );
119 /* stopUpper < 0 gives a test that will never be true. */
121 if( lmTol <= 0 )
123 lmTol = 1e-14;
125 if( xTol <= 0 )
127 xTol = 1e-8;
129 if( sigmaTol <= 0 )
131 sigmaTol = 1e-3;
134 if( costLower != 0 )
136 /* Since the dual cost is not strictly increasing, we must compute the accumulative max of all dual costs.
137 * Since the cost cannot be negative, it is safe to initialize with 0.
139 *costLower = 0;
142 bool checkCost =
143 ( costLower != 0 && stopLower < std::numeric_limits< double >::max( ) )
144 || ( costUpper != 0 && stopUpper > 0 )
145 || ( costLower != 0 && costUpper != 0 && ( gapTolAbs >= 0 || gapTolRel >= 0 ) );
147 /*-------------------------------------------------------------------------*/
148 /* Allocate memory: this section should be removed later and free memory
149 should be passed as a pointer */
150 /*-------------------------------------------------------------------------*/
151 const double ** gMem = new const double * [ n1 + n2 ];
152 DeleteOnExit< const double * > gMemDeleter( gMem );
153 const double ** g1InactiveBegin = gMem;
154 const double ** g1InactiveEnd = g1InactiveBegin + n1;
155 const double ** g1ActiveEnd = g1InactiveEnd;
156 const double ** g2InactiveBegin = gMem + n1;
157 const double ** g2InactiveEnd = g2InactiveBegin + n2;
158 const double ** g2ActiveEnd = g2InactiveEnd;
159 double * xMem = new double[ 5 * ( n1 + n2 ) ];
160 DeleteOnExit< double > xMemDeleter( xMem );
161 /* The x1 and x2 variables correspond to the generators in the order they appear in the sets defined above.
162 * The elements that corresponds to generators with active constraints are zero by definition, and
163 * do not have to be zero in memory.
165 double * x1 = xMem;
166 double * x2 = x1 + n1;
167 /* Similarly, only the last part of the Lagrange multiplier vectors are actually used.
169 double * lm1 = x2 + n2;
170 double * lm2 = lm1 + n1;
171 /* ... and only the first part of the subproblem step vectors.
173 double * d1 = lm2 + n2;
174 double * d2 = d1 + n1;
175 double * x1short = d2 + n2;
176 double * x2short = x1short + n1;
177 double * x1long = x2short + n2;
178 double * x2long = x1long + n1;
180 const double xInit = 1. / n1;
181 double * xDst = x1;
182 const double * src = g1;
183 for( const double ** dst = g1InactiveBegin; dst != g1InactiveEnd; ++dst, src += dim, ++xDst )
185 *dst = src;
186 *xDst = xInit;
190 const double xInit = 1. / n2;
191 double * xDst = x2;
192 const double * src = g2;
193 for( const double ** dst = g2InactiveBegin; dst != g2InactiveEnd; ++dst, src += dim, ++xDst )
195 *dst = src;
196 *xDst = xInit;
199 double * yDiff = new double[ dim ];
200 DeleteOnExit< double > xDiffDeleter( yDiff );
201 double * yDiffEnd = yDiff + dim;
203 double * sub_work = new double[ ( dim + n1 + n2 - 2 ) * ( n1 + n2 - 2 ) + 2 * dim + n1 + n2 ];
204 DeleteOnExit< double > sub_workDeleter( sub_work );
206 #ifdef DISPLAY
207 ++polytope_distance_generator_form_CallCount;
208 std::cerr << "Polytope distance call count: " << polytope_distance_generator_form_CallCount << std::endl ;
209 if( polytope_distance_generator_form_CallCount <= 10 )
211 std::cerr << "=== Problem number " << polytope_distance_generator_form_CallCount << " ===" << std::endl ;
212 Computation::polytope_distance_generator_form_write_data( std::cerr,
213 dim,
214 n1, g1,
215 n2, g2 );
216 std::cerr << "========================" << std::endl ;
218 #endif
219 #ifdef TRACEFIGS
220 if( dim != 2 )
222 std::cout << "##echo Invalid program, dimension of metric space is not two!" << std::endl ;
224 bool tracefigPageOpen = false;
225 std::cout << "#catalog.[setbboxgroup 'a]" << std::endl ;
227 std::cout << std::fixed << "g1: ";
228 const double * begin = g1;
229 const double * end = begin + dim * n1;
230 for( const double * gSrc = begin; gSrc != end; gSrc += dim )
232 if( gSrc != begin )
234 std::cout << "--" ;
236 std::cout << "(" ;
237 for( const double * src = gSrc; src != gSrc + dim; ++src )
239 if( src != gSrc )
241 std::cout << "," ;
243 std::cout << "(" << *src << "bp)" ;
245 std::cout << ")" ;
247 std::cout << std::endl ;
250 std::cout << "g2: ";
251 const double * begin = g2;
252 const double * end = begin + dim * n2;
253 for( const double * gSrc = begin; gSrc != end; gSrc += dim )
255 if( gSrc != begin )
257 std::cout << "--" ;
259 std::cout << "(" ;
260 for( const double * src = gSrc; src != gSrc + dim; ++src )
262 if( src != gSrc )
264 std::cout << "," ;
266 std::cout << "(" << *src << "bp)" ;
268 std::cout << ")" ;
270 std::cout << std::endl ;
272 #endif
274 bool checkNextCost = true;
275 /* The number of iterations is theoretically bounded by the number of possible combinations of activated constraints. Here, we know that
276 * one can never activate all constraints in a polytope, leading to the following bound.
278 const size_t maxIter = ( (2<<(n1-1)) - 1 ) * ( (2<<(n2-1)) - 1 );
279 for( size_t iterCount = 0; ; ++iterCount )
281 if( iterCount > maxIter )
283 status->code = QP_FAIL;
284 status->failure = QP_FAIL_ITER;
285 #ifdef DISPLAY
286 std::cerr << "Early return: Maximum number of iterations exceeded." << std::endl ;
287 #endif
288 #ifdef TRACEFIGS
289 if( tracefigPageOpen )
291 std::cout << ")" << std::endl ;
292 tracefigPageOpen = false;
294 #endif
295 return;
297 #ifdef TRACEFIGS
298 if( tracefigPageOpen )
300 std::cout << ")" << std::endl ;
301 tracefigPageOpen = false;
303 std::cout << "#catalog << @text_size:5bp | ( newGroup" << std::endl ;
304 tracefigPageOpen = true;
305 std::cout << " << @stroking:RGB_RED&@width:2bp | [[range '0 [duration g1]].foldl \\ p e → p & [spot [g1 1*e].p] null]" << std::endl
306 << " << @stroking:RGB_GREEN&@width:2bp | [[range '0 [duration g2]].foldl \\ p e → p & [spot [g2 1*e].p] null]" << std::endl
307 << " << @dash:[dashpattern 2bp 2bp]&@width:0.3bp | [stroke [controlling_hull g1] & [controlling_hull g2]]" << std::endl ;
308 for( const double ** gSrc = g1InactiveEnd; gSrc != g2ActiveEnd; ++gSrc )
310 if( gSrc == g1ActiveEnd )
312 gSrc = g2InactiveEnd - 1;
313 continue;
315 std::cout << " << [stroke [circle 2bp] >> [shift (" ;
316 for( const double * src = *gSrc; src != *gSrc + dim; ++src )
318 if( src != *gSrc )
320 std::cout << "," ;
322 std::cout << "(" << *src << "bp)" ;
324 std::cout << ")]]" << std::endl ;
326 std::cout << " << @stroking:RGB_BLUE | [stroke (" ;
327 polytope_distance_generator_form_diff( dim,
328 g1InactiveBegin, g1InactiveEnd,
329 g2InactiveBegin, g2InactiveBegin,
330 x1, x2,
331 yDiff );
332 for( const double * src = yDiff; src != yDiff + dim; ++src )
334 if( src != yDiff )
336 std::cout << "," ;
338 std::cout << "(" << -*src << "bp)" ;
340 std::cout << ")--(" ;
341 polytope_distance_generator_form_diff( dim,
342 g1InactiveBegin, g1InactiveBegin,
343 g2InactiveBegin, g2InactiveEnd,
344 x1, x2,
345 yDiff );
346 for( const double * src = yDiff; src != yDiff + dim; ++src )
348 if( src != yDiff )
350 std::cout << "," ;
352 std::cout << "(" << *src << "bp)" ;
354 std::cout << ")]" << std::endl ;
355 #endif
356 size_t nInactive1 = g1InactiveEnd - g1InactiveBegin;
357 size_t nInactive2 = g2InactiveEnd - g2InactiveBegin;
358 polytope_distance_generator_form_diff( dim,
359 g1InactiveBegin, g1InactiveEnd,
360 g2InactiveBegin, g2InactiveEnd,
361 x1, x2,
362 yDiff );
363 if( checkCost && checkNextCost )
365 /* All generators are used here, not only those of inactive constraints. */
366 polytope_distance_generator_form_costs( dim,
367 g1InactiveBegin, g1ActiveEnd,
368 g2InactiveBegin, g2ActiveEnd,
369 yDiff,
370 costUpper,
371 costLower,
372 sub_work );
373 #ifdef DISPLAY
374 std::cerr << "Interation: " << iterCount + 1 << std::endl ;
375 std::cerr << "Cost: " << *costLower << " -- " << *costUpper << std::endl ;
376 #endif
378 if( costUpper != 0 && costLower != 0 )
380 if( *costUpper - *costLower < gapTolAbs )
382 status->code = QP_OK;
383 status->termination = QP_OK_GAP_ABS;
384 status->iterations = iterCount;
385 break;
387 if( *costUpper - *costLower < gapTolRel * *costUpper )
389 status->code = QP_OK;
390 status->termination = QP_OK_GAP_REL;
391 status->iterations = iterCount;
392 break;
395 if( costUpper != 0 && *costUpper < stopUpper )
397 status->code = QP_OK;
398 status->termination = QP_OK_UPPER;
399 status->iterations = iterCount;
400 break;
402 if( costLower != 0 && *costLower > stopLower )
404 status->code = QP_OK;
405 status->termination = QP_OK_LOWER;
406 status->iterations = iterCount;
407 break;
410 checkNextCost = true;
411 polytope_distance_generator_form_sub( dim,
412 g1InactiveBegin, g1InactiveEnd,
413 g2InactiveBegin, g2InactiveEnd,
414 yDiff,
415 sigmaTol,
416 d1, d2,
417 sub_work );
418 #ifdef DISPLAY
419 std::cerr << "Delta: " << Helpers::mathematicaFormat( nInactive1, 1, d1 )
420 << " " << Helpers::mathematicaFormat( nInactive2, 1, d2 ) << std::endl ;
421 #endif
423 /* Compute squared norm of subproblem step. */
424 double stepNorm2 = 0;
426 double * src = d1;
427 double * end = d1 + nInactive1;
428 for( ; src != end; ++src )
430 stepNorm2 += (*src) * (*src);
434 double * src = d2;
435 double * end = d2 + nInactive2;
436 for( ; src != end; ++src )
438 stepNorm2 += (*src) * (*src);
441 if( stepNorm2 < xTol )
443 /* No step taken in sub-problem. Are there any constraints to drop? */
444 /* Compute lagrange multipliers.
445 * Note that the lagrange multipliers of the inactive constraints are only used for temporary storage.
446 * Normalizing yDiff will only scale the Langrange multipliers, and by normalizing yDiff we give
447 * the Lagrange multiplier tolerance (lmTol) a geometric interpretation; lmTol becomes a gap tolerance which
448 * is linear in the distance (not distance squared).
450 polytope_distance_generator_form_diff( dim,
451 g1InactiveBegin, g1InactiveEnd,
452 g2InactiveBegin, g2InactiveEnd,
453 x1short, x2short,
454 yDiff );
456 /* Normalize yDiff. */
457 double lengthSquared = 0;
458 for( const double * src = yDiff; src != yDiffEnd; ++src )
460 lengthSquared += (*src) * (*src);
462 const double yDiffScaling = 1. / sqrt( lengthSquared );
463 for( double * dst = yDiff; dst != yDiffEnd; ++dst )
465 *dst *= yDiffScaling;
469 /* Initialize with negative gradient. */
470 double * lmDst = lm1;
471 for( const double ** gSrc = g1InactiveBegin; gSrc != g1ActiveEnd; ++gSrc, ++lmDst )
473 double tmp = 0;
474 const double * ySrc = yDiff;
475 const double * src = *gSrc;
476 for( ; ySrc != yDiffEnd; ++ySrc, ++src )
478 tmp -= (*ySrc) * (*src);
480 *lmDst = tmp;
482 /* Let nu be the dual variable of the given equality constraint (all variables sum to 1). */
483 double nuMax = -1;
484 double nu = 0; /* Initialize to avoid compiler warnings. */
485 double * lmInactiveEnd = lm1 + nInactive1;
486 double * lmPtr = lm1;
487 for( ; lmPtr != lmInactiveEnd; ++lmPtr )
489 double lmAbs = fabs( *lmPtr );
490 if( lmAbs > nuMax )
492 nuMax = lmAbs;
493 nu = *lmPtr;
496 /* Then the remaining Lagrange multipliers follow. */
497 double * lmActiveEnd = lm1 + n1;
498 for( ; lmPtr != lmActiveEnd; ++lmPtr )
500 *lmPtr -= nu;
505 /* Compute Lagrange multipliers for the second polytope in the same way... up to sign changes. */
506 double * lmDst = lm2;
507 for( const double ** gSrc = g2InactiveBegin; gSrc != g2ActiveEnd; ++gSrc, ++lmDst )
509 double tmp = 0;
510 const double * ySrc = yDiff;
511 const double * src = *gSrc;
512 for( ; ySrc != yDiffEnd; ++ySrc, ++src )
514 tmp += (*ySrc) * (*src); /* Note sign! */
516 *lmDst = tmp;
518 double nuMax = -1;
519 double nu = 0; /* Initialize to avoid compiler warnings. */
520 double * lmInactiveEnd = lm2 + nInactive2;
521 double * lmPtr = lm2;
522 for( ; lmPtr != lmInactiveEnd; ++lmPtr )
524 double lmAbs = fabs( *lmPtr );
525 if( lmAbs > nuMax )
527 nuMax = lmAbs;
528 nu = *lmPtr;
531 double * lmActiveEnd = lm2 + n2;
532 for( ; lmPtr != lmActiveEnd; ++lmPtr )
534 *lmPtr -= nu;
538 #ifdef DISPLAY
539 std::cerr << "lambda: " << Helpers::mathematicaFormat( n1 - nInactive1, 1, lm1 + nInactive1 )
540 << " " << Helpers::mathematicaFormat( n2 - nInactive2, 1, lm2 + nInactive2 ) << std::endl ;
541 #endif
542 #ifdef TRACEFIGS
544 const double * lambdaSrc = lm1 + nInactive1;
545 for( const double ** gSrc = g1InactiveEnd; gSrc != g2ActiveEnd; ++gSrc, ++lambdaSrc )
547 if( gSrc == g1ActiveEnd )
549 gSrc = g2InactiveEnd - 1;
550 lambdaSrc = lm2 + nInactive2 - 1;
551 continue;
553 std::cout << " << @nonstroking:[gray 0.5] | (newText << `" << std::scientific << std::setprecision(1) << *lambdaSrc << std::fixed << "´) >> [shift (0,3bp)+(" ;
554 for( const double * src = *gSrc; src != *gSrc + dim; ++src )
556 if( src != *gSrc )
558 std::cout << "," ;
560 std::cout << "(" << *src << "bp)" ;
562 std::cout << ")]" << std::endl ;
565 #endif
567 bool dropped = false;
569 const double * lmSrc = lm1 + nInactive1;
570 const double * lmEnd = lm1 + n1;
571 const double ** gSrc = g1InactiveEnd;
572 double * xDst = x1 + nInactive1;
573 for( ; lmSrc != lmEnd; ++lmSrc, ++gSrc )
575 if( *lmSrc < -lmTol )
577 /* Drop constraint.
578 * Swap the generator out of the active set, so that it appears at the end
579 * of the inactive set.
581 const double * gTmp = *g1InactiveEnd;
582 *g1InactiveEnd = *gSrc;
583 *gSrc = gTmp;
584 ++g1InactiveEnd;
586 /* Since the constraint has been active, the corresponding variable has been zero.
588 *xDst = 0;
589 ++xDst;
591 dropped = true;
596 const double * lmSrc = lm2 + nInactive2;
597 const double * lmEnd = lm2 + n2;
598 const double ** gSrc = g2InactiveEnd;
599 double * xDst = x2 + nInactive2;
600 for( ; lmSrc != lmEnd; ++lmSrc, ++gSrc )
602 if( *lmSrc < -lmTol )
604 /* Drop constraint.
605 * Swap the generator out of the active set, so that it appears at the end
606 * of the inactive set.
608 const double * gTmp = *g2InactiveEnd;
609 *g2InactiveEnd = *gSrc;
610 *gSrc = gTmp;
611 ++g2InactiveEnd;
613 /* Since the constraint has been active, the corresponding variable has been zero.
615 *xDst = 0;
616 ++xDst;
618 dropped = true;
622 if( ! dropped )
624 /* We reached optimum!
626 status->code = QP_OK;
627 status->termination = QP_OK_LM;
628 status->iterations = iterCount;
629 break;
631 /* If we dropped some constraint, we just return to the top of the loop.
632 * No need to re-normalize variables since added variables are exactly zero.
634 checkNextCost = false; /* No need to check the cost in the next iteration since no step was made. */
636 else
638 /* Find inactive constraints that limit the step length in each polytope. */
639 double s1 = 1;
640 const double ** gStop1 = 0;
642 const double ** gSrc = g1InactiveBegin;
643 const double ** gEnd = g1InactiveEnd;
644 const double * dSrc = d1;
645 const double * xSrc = x1;
646 for( ; gSrc != gEnd; ++gSrc, ++dSrc, ++xSrc )
648 double denom = -*dSrc;
649 if( denom > 0 )
651 double sTmp = *xSrc / denom;
652 if( sTmp < s1 )
654 s1 = sTmp;
655 gStop1 = gSrc;
660 double s2 = 1;
661 const double ** gStop2 = 0;
663 const double ** gSrc = g2InactiveBegin;
664 const double ** gEnd = g2InactiveEnd;
665 const double * dSrc = d2;
666 const double * xSrc = x2;
667 for( ; gSrc != gEnd; ++gSrc, ++dSrc, ++xSrc )
669 double denom = -*dSrc;
670 if( denom > 0 )
672 double sTmp = *xSrc / denom;
673 if( sTmp < s2 )
675 s2 = sTmp;
676 gStop2 = gSrc;
682 /* In a basic active set algorithm, one uses only the smaller of s1 and s2, but here
683 * we actually use both, if that gives further decrease in the objective function.
685 double sShort = ( s1 < s2 ) ? s1 : s2;
687 const double * xSrc = x1;
688 const double * xEnd = x1 + nInactive1;
689 const double * dSrc = d1;
690 double * dst = x1short;
691 for( ; xSrc != xEnd; ++xSrc, ++dSrc, ++dst )
693 *dst = *xSrc + sShort * (*dSrc);
697 const double * xSrc = x2;
698 const double * xEnd = x2 + nInactive2;
699 const double * dSrc = d2;
700 double * dst = x2short;
701 for( ; xSrc != xEnd; ++xSrc, ++dSrc, ++dst )
703 *dst = *xSrc + sShort * (*dSrc);
706 double shortCost = 0;
707 polytope_distance_generator_form_diff( dim,
708 g1InactiveBegin, g1InactiveEnd,
709 g2InactiveBegin, g2InactiveEnd,
710 x1short, x2short,
711 yDiff );
712 for( const double * src = yDiff; src != yDiffEnd; ++src )
714 shortCost += (*src) * (*src);
717 if( s1 > s2 )
719 const double * xSrc = x1;
720 const double * xEnd = x1 + nInactive1;
721 const double * dSrc = d1;
722 double * dst = x1long;
723 for( ; xSrc != xEnd; ++xSrc, ++dSrc, ++dst )
725 *dst = *xSrc + s1 * (*dSrc);
727 double longCost = 0;
728 polytope_distance_generator_form_diff( dim,
729 g1InactiveBegin, g1InactiveEnd,
730 g2InactiveBegin, g2InactiveEnd,
731 x1long, x2short,
732 yDiff );
733 for( const double * src = yDiff; src != yDiffEnd; ++src )
735 longCost += (*src) * (*src);
738 if( longCost > shortCost )
740 /* Dont make long step. */
741 gStop1 = 0;
742 s1 = s2;
745 else if( s2 > s1 )
747 const double * xSrc = x2;
748 const double * xEnd = x2 + nInactive2;
749 const double * dSrc = d2;
750 double * dst = x2long;
751 for( ; xSrc != xEnd; ++xSrc, ++dSrc, ++dst )
753 *dst = *xSrc + s2 * (*dSrc);
755 double longCost = 0;
756 polytope_distance_generator_form_diff( dim,
757 g1InactiveBegin, g1InactiveEnd,
758 g2InactiveBegin, g2InactiveEnd,
759 x1short, x2long,
760 yDiff );
761 for( const double * src = yDiff; src != yDiffEnd; ++src )
763 longCost += (*src) * (*src);
766 if( longCost > shortCost )
768 /* Dont make long step. */
769 gStop2 = 0;
770 s2 = s1;
774 /* Update variables and active constraints.
775 * Important: variables first, or we will get a mess!
778 const double * xSrc = x1;
779 const double * xEnd = x1 + nInactive1;
780 const double * dSrc = d1;
781 double * dst = x1;
782 for( ; xSrc != xEnd; ++xSrc, ++dSrc, ++dst )
784 *dst = *xSrc + s1 * (*dSrc);
787 if( gStop1 != 0 )
789 /* Swap constraints, and swap variables in the same way. New variable is zero by construction. */
790 --g1InactiveEnd;
791 const double * gTmp = *gStop1;
792 *gStop1 = *g1InactiveEnd;
793 *g1InactiveEnd = gTmp;
794 *( x1 + ( gStop1 - g1InactiveBegin ) ) = *( x1 + ( g1InactiveEnd - g1InactiveBegin ) );
797 const double * xSrc = x2;
798 const double * xEnd = x2 + nInactive2;
799 const double * dSrc = d2;
800 double * dst = x2;
801 for( ; xSrc != xEnd; ++xSrc, ++dSrc, ++dst )
803 *dst = *xSrc + s2 * (*dSrc);
806 if( gStop2 != 0 )
808 /* Swap constraints, and swap variables in the same way. New variable is zero by construction. */
809 --g2InactiveEnd;
810 const double * gTmp = *gStop2;
811 *gStop2 = *g2InactiveEnd;
812 *g2InactiveEnd = gTmp;
813 *( x2 + ( gStop2 - g2InactiveBegin ) ) = *( x2 + ( g2InactiveEnd - g2InactiveBegin ) );
818 #ifdef TRACEFIGS
819 if( tracefigPageOpen )
821 std::cout << ")" << std::endl ;
822 tracefigPageOpen = false;
824 #endif
826 /* Compute final cost. */
828 polytope_distance_generator_form_diff( dim,
829 g1InactiveBegin, g1InactiveEnd,
830 g2InactiveBegin, g2InactiveEnd,
831 x1, x2,
832 yDiff );
833 polytope_distance_generator_form_costs( dim,
834 g1InactiveBegin, g1ActiveEnd,
835 g2InactiveBegin, g2ActiveEnd,
836 yDiff,
837 costUpper,
838 costLower,
839 sub_work );
840 #ifdef DISPLAY
841 std::cerr << "Terminating with final cost: " << *costLower << " -- " << *costUpper << std::endl ;
842 #endif
845 /* Compute final point in each polytope. */
846 if( y1 != 0 )
848 double * y = y1;
849 double * yEnd = y + dim;
850 for( double * dst = y; dst != yEnd; ++dst )
852 *dst = 0;
854 const double * xSrc = x1;
855 for( const double ** gSrc = g1InactiveBegin; gSrc != g1InactiveEnd; ++gSrc, ++xSrc )
857 double xVal = *xSrc;
858 const double * src = *gSrc;
859 for( double * dst = y; dst != yEnd; ++dst, ++src )
861 *dst += xVal * (*src);
865 if( y2 != 0 )
867 double * y = y2;
868 double * yEnd = y + dim;
869 for( double * dst = y; dst != yEnd; ++dst )
871 *dst = 0;
873 const double * xSrc = x2;
874 for( const double ** gSrc = g2InactiveBegin; gSrc != g2InactiveEnd; ++gSrc, ++xSrc )
876 double xVal = *xSrc;
877 const double * src = *gSrc;
878 for( double * dst = y; dst != yEnd; ++dst, ++src )
880 *dst += xVal * (*src);
886 void
887 Computation::polytope_distance_generator_form_diff( const size_t dim,
888 const double ** g1Begin, const double ** g1End,
889 const double ** g2Begin, const double ** g2End,
890 const double * x1, const double * x2,
891 double * yDiff )
893 double * yDiffEnd = yDiff + dim;
894 for( double * dst = yDiff; dst != yDiffEnd; ++dst )
896 *dst = 0;
899 const double * xSrc = x2;
900 for( const double ** gSrc = g2Begin; gSrc != g2End; ++gSrc, ++xSrc )
902 double xVal = *xSrc;
903 const double * src = *gSrc;
904 for( double * dst = yDiff; dst != yDiffEnd; ++dst, ++src )
906 *dst += xVal * (*src);
911 const double * xSrc = x1;
912 for( const double ** gSrc = g1Begin; gSrc != g1End; ++gSrc, ++xSrc )
914 double xVal = *xSrc;
915 const double * src = *gSrc;
916 for( double * dst = yDiff; dst != yDiffEnd; ++dst, ++src )
918 *dst -= xVal * (*src);
924 void
925 Computation::polytope_distance_generator_form_costs( const size_t dim,
926 const double ** g1Begin, const double ** g1End,
927 const double ** g2Begin, const double ** g2End,
928 const double * yDiff,
929 double * primalCost,
930 double * dualCostAccum,
931 double * work )
933 const double * yDiffEnd = yDiff + dim;
935 double pCost = 0;
936 for( const double * src = yDiff; src != yDiffEnd; ++src )
938 pCost += (*src) * (*src);
941 if( primalCost != 0 )
943 *primalCost = pCost;
946 if( dualCostAccum == 0 )
948 return;
951 double * n = work; /* Normalaized normal to planes separating the two polytopes. */
952 double * nEnd = work + dim;
954 /* Normalize yDiff and store in n. */
955 double yDiffScale = 1. / sqrt( pCost );
956 double * dst = n;
957 for( const double * src = yDiff; src != yDiffEnd; ++src, ++dst )
959 *dst = yDiffScale * (*src);
963 /* Find dual cost by projecting generators onto n. */
964 double min2 = std::numeric_limits< double >::max( );
966 for( const double ** gSrc = g2Begin; gSrc != g2End; ++gSrc )
968 double tmp = 0;
969 const double * src = *gSrc;
970 for( const double * nSrc = n; nSrc != nEnd; ++nSrc, ++src )
972 tmp += (*nSrc) * (*src);
974 min2 = std::min( min2, tmp );
978 double max1 = -std::numeric_limits< double >::max( );
980 for( const double ** gSrc = g1Begin; gSrc != g1End; ++gSrc )
982 double tmp = 0;
983 const double * src = *gSrc;
984 for( const double * nSrc = n; nSrc != nEnd; ++nSrc, ++src )
986 tmp += (*nSrc) * (*src);
988 max1 = std::max( max1, tmp );
992 double diff = min2 - max1;
993 if( diff > 0 )
995 *dualCostAccum = std::max( *dualCostAccum, diff * diff );
1000 void
1001 Computation::polytope_distance_generator_form_sub( const size_t dim,
1002 const double ** g1Begin, const double ** g1End,
1003 const double ** g2Begin, const double ** g2End,
1004 const double * yDiff,
1005 double sigmaTol,
1006 double * d1, double * d2,
1007 double * work )
1009 size_t n1 = g1End - g1Begin;
1010 size_t n2 = g2End - g2Begin;
1011 double * d1End = d1 + n1;
1012 double * d2End = d2 + n2;
1013 for( double * dst = d1; dst != d1End; ++dst )
1015 *dst = 0;
1017 for( double * dst = d2; dst != d2End; ++dst )
1019 *dst = 0;
1021 if( n1 == 1 && n2 == 1 ) /* They must never be less than 1! */
1023 /* If there is only one inactive constraint in a polytope, the corresponding variable is in fact locked to 1.
1024 * Henece, we cannot move in that polytope. If we can't move in either polytope, there is nothing to optimize.
1025 * We have already set the result to 0.
1027 return;
1029 const double ** N1 = Computation::polytope_distance_generator_form_kernel( n1 );
1030 const double ** N2 = Computation::polytope_distance_generator_form_kernel( n2 );
1032 gsl_vector y_gsl;
1033 y_gsl.data = const_cast< double * >( yDiff );
1034 y_gsl.size = dim;
1035 y_gsl.stride = 1;
1037 gsl_vector d_gsl;
1038 d_gsl.data = work;
1039 work += ( d_gsl.size = n1 + n2 - 2 );
1040 d_gsl.stride = 1;
1042 if( dim >= n1 + n2 - 2 )
1044 /* Over-determined case. SVD of A can be computed as usual. */
1045 gsl_matrix A;
1046 A.data = work;
1047 work += ( A.size1 = dim ) * ( A.size2 = A.tda = n1 + n2 - 2 );
1048 gsl_matrix V;
1049 V.data = work;
1050 work += ( V.size1 = n1 + n2 - 2 ) * ( V.size2 = V.tda = n1 + n2 - 2 );
1051 gsl_vector S;
1052 S.data = work;
1053 work += ( S.size = n1 + n2 - 2 );
1054 S.stride = 1;
1055 gsl_vector SV_work;
1056 SV_work.data = work;
1057 work += ( SV_work.size = n1 + n2 - 2 );
1058 SV_work.stride = 1;
1060 /* Populate A with data.
1061 * The data access here can be improved. */
1063 /* Polytope 1 */
1064 double * ARowStart = A.data;
1065 for( size_t row = 0; row < dim; ++row, ARowStart += A.tda )
1067 const double ** NcolStart = N1;
1068 const double ** NcolEnd = N1 + n1 - 1;
1069 double * Adst = ARowStart;
1070 for( ; NcolStart != NcolEnd; ++NcolStart, ++Adst )
1072 double tmp = 0;
1073 const double ** gcolStart = g1Begin;
1074 const double ** gcolEnd = g1End;
1075 const double * Nsrc = *NcolStart;
1076 for( ; gcolStart != gcolEnd; ++gcolStart, ++Nsrc )
1078 tmp += (*((*gcolStart) + row)) * (*Nsrc);
1080 *Adst = tmp;
1085 /* Polytope 2 */
1086 double * ARowStart = A.data + ( n1 - 1 );
1087 for( size_t row = 0; row < dim; ++row, ARowStart += A.tda )
1089 const double ** NcolStart = N2;
1090 const double ** NcolEnd = N2 + n2 - 1;
1091 double * Adst = ARowStart;
1092 for( ; NcolStart != NcolEnd; ++NcolStart, ++Adst )
1094 double tmp = 0;
1095 const double ** gcolStart = g2Begin;
1096 const double ** gcolEnd = g2End;
1097 const double * Nsrc = *NcolStart;
1098 for( ; gcolStart != gcolEnd; ++gcolStart, ++Nsrc )
1100 tmp -= (*((*gcolStart) + row)) * (*Nsrc); /* Note the sign! */
1102 *Adst = tmp;
1107 gsl_linalg_SV_decomp( & A, & V, & S, & SV_work );
1109 /* Truncate small singular values. */
1110 double * dst = S.data;
1111 double * end = S.data + S.size;
1112 const double sigmaCut = sigmaTol * (*dst);
1113 ++dst;
1114 for( ; dst != end; ++dst )
1116 if( *dst < sigmaCut )
1118 break;
1121 for( ; dst != end; ++dst )
1123 *dst = 0;
1126 gsl_linalg_SV_solve( & A, & V, & S, & y_gsl, & d_gsl );
1128 else
1130 /* Under-determined case. Need to transpose A to compute SVD. */
1131 gsl_matrix A;
1132 A.data = work;
1133 work += ( A.size1 = n1 + n2 - 2 ) * ( A.size2 = A.tda = dim );
1134 gsl_matrix V;
1135 V.data = work;
1136 work += ( V.size1 = dim ) * ( V.size2 = V.tda = dim );
1137 gsl_vector S;
1138 S.data = work;
1139 work += ( S.size = dim );
1140 S.stride = 1;
1141 gsl_vector SV_work;
1142 SV_work.data = work;
1143 work += ( SV_work.size = dim );
1144 SV_work.stride = 1;
1146 /* Populate A with data.
1147 * The data access here can be improved. */
1149 /* Polytope 1 */
1150 double * AColStart = A.data;
1151 for( size_t row = 0; row < dim; ++row, ++AColStart )
1153 const double ** NcolStart = N1;
1154 const double ** NcolEnd = N1 + n1 - 1;
1155 double * Adst = AColStart;
1156 for( ; NcolStart != NcolEnd; ++NcolStart, Adst += A.tda )
1158 double tmp = 0;
1159 const double ** gcolStart = g1Begin;
1160 const double ** gcolEnd = g1End;
1161 const double * Nsrc = *NcolStart;
1162 for( ; gcolStart != gcolEnd; ++gcolStart, ++Nsrc )
1164 tmp += (*((*gcolStart) + row)) * (*Nsrc);
1166 *Adst = tmp;
1171 /* Polytope 2 */
1172 double * AColStart = A.data + ( n1 - 1 ) * A.tda;
1173 for( size_t row = 0; row < dim; ++row, ++AColStart )
1175 const double ** NcolStart = N2;
1176 const double ** NcolEnd = N2 + n2 - 1;
1177 double * Adst = AColStart;
1178 for( ; NcolStart != NcolEnd; ++NcolStart, Adst += A.tda )
1180 double tmp = 0;
1181 const double ** gcolStart = g2Begin;
1182 const double ** gcolEnd = g2End;
1183 const double * Nsrc = *NcolStart;
1184 for( ; gcolStart != gcolEnd; ++gcolStart, ++Nsrc )
1186 tmp -= (*((*gcolStart) + row)) * (*Nsrc); /* Note the sign! */
1188 *Adst = tmp;
1193 gsl_linalg_SV_decomp( & A, & V, & S, & SV_work );
1195 /* Truncate small singular values. */
1196 double * dst = S.data;
1197 double * end = S.data + S.size;
1198 const double sigmaCut = sigmaTol * (*dst);
1199 ++dst;
1200 for( ; dst != end; ++dst )
1202 if( *dst < sigmaCut )
1204 break;
1207 for( ; dst != end; ++dst )
1209 *dst = 0;
1212 /* Now, instead of using
1213 * gsl_linalg_SV_solve( & A, & V, & S, & y_gsl, & d_gsl );
1214 * we must solve the equation ourselves since the decomposition is transposed.
1216 gsl_blas_dgemv( CblasTrans, 1, & V, & y_gsl, 0, & SV_work );
1218 const double * src = S.data;
1219 const double * srcEnd = S.data + S.size;
1220 double * dst = SV_work.data;
1221 for( ; src != srcEnd; ++src, ++dst )
1223 if( *src > 0 )
1225 *dst /= *src;
1227 else
1229 *dst = 0;
1233 gsl_blas_dgemv( CblasNoTrans, 1, & A, & SV_work, 0, & d_gsl );
1237 /* Compute output variables. */
1238 /* First clear. */
1240 double * dst = d1;
1241 double * dstEnd = d1 + n1;
1242 for( ; dst != dstEnd; ++dst )
1244 *dst = 0;
1248 double * dst = d2;
1249 double * dstEnd = d2 + n2;
1250 for( ; dst != dstEnd; ++dst )
1252 *dst = 0;
1256 /* Then add linear combination of basis vectors. */
1257 const double * dSrc = d_gsl.data;
1259 const double ** NSrc = N1;
1260 const double ** NSrcEnd = N1 + ( n1 - 1 );
1261 for( ; NSrc != NSrcEnd; ++NSrc, ++dSrc )
1263 double w = *dSrc;
1264 double * dst = d1;
1265 double * dstEnd = d1 + n1;
1266 const double * src = *NSrc;
1267 for( ; dst != dstEnd; ++dst, ++src )
1269 *dst += w * *src;
1274 const double ** NSrc = N2;
1275 const double ** NSrcEnd = N2 + ( n2 - 1 );
1276 for( ; NSrc != NSrcEnd; ++NSrc, ++dSrc )
1278 double w = *dSrc;
1279 double * dst = d2;
1280 double * dstEnd = d2 + n2;
1281 const double * src = *NSrc;
1282 for( ; dst != dstEnd; ++dst, ++src )
1284 *dst += w * *src;
1291 void
1292 Computation::polytope_distance_generator_form_sub_Gilbert( const size_t dim,
1293 const double ** g1Begin, const double ** g1End,
1294 const double ** g2Begin, const double ** g2End,
1295 const double * yDiff,
1296 double * d1, double * d2 )
1301 const double **
1302 Computation::polytope_distance_generator_form_kernel( size_t n )
1304 static std::vector< const double ** > mem( 0 );
1305 if( n < mem.size( ) )
1307 return mem[ n ];
1309 mem.reserve( n + 1 );
1310 while( mem.size( ) <= n )
1312 size_t m = mem.size( );
1313 if( m == 0 )
1315 /* Invalid dimension, just skip. */
1316 mem.push_back( 0 );
1317 continue;
1319 if( m == 1 )
1321 /* Trivial case; kernel has dimension 0 */
1322 mem.push_back( 0 );
1323 continue;
1325 const double ** res = new const double * [ m - 1 ];
1326 gsl_matrix * qr = gsl_matrix_alloc( m, m );
1327 gsl_vector * tau = gsl_vector_alloc( m );
1328 gsl_matrix_set_zero( qr );
1329 for( size_t i = 0; i < m; ++i )
1331 gsl_matrix_set( qr, i, 0, 1. );
1333 gsl_matrix * q = gsl_matrix_alloc( m, m );
1334 gsl_matrix * r = gsl_matrix_alloc( m, m );
1336 gsl_linalg_QR_decomp( qr, tau );
1337 gsl_linalg_QR_unpack( qr, tau, q, r );
1339 const double ** dst = res;
1340 const double ** dstEnd = res + ( m - 1 );
1341 double * srcQStart = q->data + 1; /* Skip first column. */
1342 for( ; dst != dstEnd; ++dst, ++srcQStart )
1344 double * v = new double[ m ];
1345 *dst = v;
1346 double * vEnd = v + m;
1347 double * srcQ = srcQStart;
1348 for( ; v != vEnd; ++v, srcQ += q->tda )
1350 *v = *srcQ;
1355 mem.push_back( res );
1357 gsl_matrix_free( qr );
1358 gsl_vector_free( tau );
1359 gsl_matrix_free( q );
1360 gsl_matrix_free( r );
1362 return mem[ n ];