condensing_based_cp_solver.cpp
Go to the documentation of this file.
1 /*
2  * This file is part of ACADO Toolkit.
3  *
4  * ACADO Toolkit -- A Toolkit for Automatic Control and Dynamic Optimization.
5  * Copyright (C) 2008-2014 by Boris Houska, Hans Joachim Ferreau,
6  * Milan Vukov, Rien Quirynen, KU Leuven.
7  * Developed within the Optimization in Engineering Center (OPTEC)
8  * under supervision of Moritz Diehl. All rights reserved.
9  *
10  * ACADO Toolkit is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public
12  * License as published by the Free Software Foundation; either
13  * version 3 of the License, or (at your option) any later version.
14  *
15  * ACADO Toolkit is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with ACADO Toolkit; if not, write to the Free Software
22  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
23  *
24  */
25 
26 
35 
36 using namespace Eigen;
37 using namespace std;
38 
40 
41 
42 //
43 // PUBLIC MEMBER FUNCTIONS:
44 //
45 
47 {
48  nConstraints = 0;
49  blockDims = 0;
50 
52 
53  cpSolver = 0;
54  cpSolverRelaxed = 0;
55 }
56 
57 
59  uint nConstraints_,
60  const DVector& blockDims_
61  ) : BandedCPsolver( _userInteraction )
62 {
63  nConstraints = nConstraints_;
64  blockDims = blockDims_;
65 
67 
68  cpSolver = new QPsolver_qpOASES( _userInteraction );
69  cpSolverRelaxed = new QPsolver_qpOASES( _userInteraction );
70 }
71 
72 
74  :BandedCPsolver( rhs )
75 {
77  blockDims = rhs.blockDims;
78 
80 
81  if( rhs.cpSolver != 0 ) cpSolver = rhs.cpSolver->clone();
82  else cpSolver = 0 ;
83 
85  else cpSolverRelaxed = 0;
86 
87  denseCP = rhs.denseCP;
88 
89  deltaX = rhs.deltaX;
90  deltaP = rhs.deltaP;
91 }
92 
93 
95 
96  if( cpSolver != 0 ) delete cpSolver ;
97  if( cpSolverRelaxed != 0 ) delete cpSolverRelaxed;
98 }
99 
100 
102 
103  if ( this != &rhs ){
104 
105  if( cpSolver != 0 ) delete cpSolver ;
106  if( cpSolverRelaxed != 0 ) delete cpSolverRelaxed;
107 
109 
111  blockDims = rhs.blockDims;
112 
114 
115  if( rhs.cpSolver != 0 ) cpSolver = rhs.cpSolver->clone();
116  else cpSolver = 0 ;
117 
119  else cpSolverRelaxed = 0;
120 
121  denseCP = rhs.denseCP;
122 
123  deltaX = rhs.deltaX;
124  deltaP = rhs.deltaP;
125  }
126  return *this;
127 }
128 
129 
131 {
132  return new CondensingBasedCPsolver(*this);
133 }
134 
135 
136 
138  )
139 {
140  iter = iter_;
141 
142  // INITIALIZE THE CONDENSING MATRICES:
143  // -----------------------------------
146 
147 
148  // INITIALIZE THE CP SOLVER:
149  // -------------------------
150  int infeasibleQPhandling;
151  get( INFEASIBLE_QP_HANDLING,infeasibleQPhandling );
152 
153  if ( initializeCPsolver( (InfeasibleQPhandling)infeasibleQPhandling ) != SUCCESSFUL_RETURN )
155 
156 
157  return SUCCESSFUL_RETURN;
158 }
159 
160 
161 
163  )
164 {
165  RealClock clock;
166 
167  // CONDENSE THE KKT-SYSTEM:
168  // ------------------------
169 
170  int printLevel;
171  get( PRINTLEVEL,printLevel );
172 
173  if ( (PrintLevel)printLevel >= HIGH )
174  cout << "--> Condesing banded QP ...\n";
175 
176  clock.reset( );
177  clock.start( );
178 
179  returnValue returnvalue = condense( cp );
180  if( returnvalue != SUCCESSFUL_RETURN ) return ACADOERROR(returnvalue);
181 
182  clock.stop( );
184 
185  if ( (PrintLevel)printLevel >= HIGH )
186  cout << "<-- Condesing banded QP done.\n";
187 
188  return SUCCESSFUL_RETURN;
189 }
190 
191 
193  )
194 {
196  prepareSolve( cp );
197 
198 
199  // ADD THE FEEDBACK DATA TO THE DENSE QP (IF SPECIFIED):
200  // -----------------------------------------------------
201 
202  if ( deltaX.isEmpty( ) == BT_FALSE )
203  {
204  for( uint run1 = 0; run1 < getNX(); run1++ )
205  {
206  denseCP.lb(run1) = deltaX(run1);
207  denseCP.ub(run1) = deltaX(run1);
208  }
209  }
210 
211  if ( deltaP.isEmpty( ) == BT_FALSE )
212  {
213  for( uint run1 = 0; run1 < getNP(); run1++ )
214  {
215  denseCP.lb(getNX()+getNumPoints()*getNXA()+run1) = deltaP(run1);
216  denseCP.ub(getNX()+getNumPoints()*getNXA()+run1) = deltaP(run1);
217  }
218  }
219 
220 // cp.deltaX.print();
221 
222  // Solve QP subproblem
223  // ------------------------------------
224  int printLevel;
225  get( PRINTLEVEL,printLevel );
226 
227  if ( (PrintLevel)printLevel >= HIGH )
228  cout << "--> Solving condesed QP ...\n";
229 
230 
231 
232  returnValue returnvalue = solveCPsubproblem( );
233  if( returnvalue != SUCCESSFUL_RETURN ) return ACADOERROR( RET_BANDED_CP_SOLUTION_FAILED );
234 
235  if ( (PrintLevel)printLevel >= HIGH )
236  cout << "<-- Solving condesed QP done.\n";
237 
238  // Expand the KKT-System if neccessary:
239  // ------------------------------------
240 
242  return finalizeSolve( cp );
243  else
244  return SUCCESSFUL_RETURN;
245 }
246 
247 
249  )
250 {
251  RealClock clock;
252 
253  int printLevel;
254  get( PRINTLEVEL,printLevel );
255 
256  if ( (PrintLevel)printLevel >= HIGH )
257  cout << "--> Expanding condensed QP solution ...\n";
258 
259  // Expand the KKT-System if neccessary:
260  // ------------------------------------
261  clock.reset( );
262  clock.start( );
263 
264  returnValue returnvalue = expand( cp );
265  if( returnvalue != SUCCESSFUL_RETURN ) return ACADOERROR(returnvalue);
266 
267  clock.stop( );
268  setLast( LOG_TIME_EXPAND,clock.getTime() );
269 
270  if ( (PrintLevel)printLevel >= HIGH )
271  cout << "<-- Expanding condensed QP solution done.\n";
272 
273  return returnvalue;
274 }
275 
276 
277 
279 {
280  if ( p_.getDim( ) != getNP( ) )
282 
283  uint startIdx = getNX() + getNumPoints()*getNXA();
284 
285  for( uint i=0; i<getNP(); ++i )
286  p_( i ) = (*denseCP.x)( startIdx+i );
287 
288  return SUCCESSFUL_RETURN;
289 }
290 
291 
293 {
294  if ( u0_.getDim( ) != getNU( ) )
296 
297  uint startIdx = getNX() + getNumPoints()*getNXA() + getNP();
298 
299  for( uint i=0; i<getNU(); ++i )
300  u0_( i ) = (*denseCP.x)( startIdx+i );
301 
302  return SUCCESSFUL_RETURN;
303 }
304 
305 
306 
308 {
309  if ( cpSolver == 0 )
311 
312  return cpSolver->getVarianceCovariance( denseCP.H,var );
313 }
314 
315 
316 
318  const DVector& DeltaP
319  )
320 {
321  deltaX = DeltaX;
322  deltaP = DeltaP;
323 
324  return SUCCESSFUL_RETURN;
325 }
326 
327 
328 
330 {
333 
335 
336  return SUCCESSFUL_RETURN;
337 }
338 
339 
341 {
342  if ( condensingStatus == COS_FROZEN )
344 
345  return SUCCESSFUL_RETURN;
346 }
347 
348 
349 
350 //
351 // PROTECTED MEMBER FUNCTIONS:
352 //
353 
355 
356  if( dampingFactor < 0.0 ) return SUCCESSFUL_RETURN;
357 
358  // COMPUTE THE EIGENVALUES OF THE HESSIAN:
359  // ---------------------------------------
360 
362  MatrixXd V = es.eigenvectors();
363  VectorXd D = es.eigenvalues();
364 
365  // OVER-PROJECT THE EIGENVALUES BASED ON THE DAMPING TECHNIQUE:
366  // ------------------------------------------------------------
367 
368  for (unsigned el = 0; el < D.size(); el++)
369  if (D( el ) <= 0.1 * dampingFactor)
370  {
371  if (fabs(D( el )) >= dampingFactor)
372  D( el ) = fabs(D( el ));
373  else
374  D( el ) = dampingFactor;
375  }
376 
377  // RECONSTRUCT THE PROJECTED HESSIAN MATRIX:
378  // -----------------------------------------
379 
380  H_ = V * D.asDiagonal() * V.inverse();
381 
382  return SUCCESSFUL_RETURN;
383 }
384 
385 
386 
388 {
389  if( denseCP.isQP() == BT_FALSE )
391 
392  returnValue returnvalue;
393 
394  // ensure that Hessian matrix is symmetric
395  if ( denseCP.H.isSymmetric( ) == BT_FALSE )
396  {
398  denseCP.H.symmetrize();
399  }
400 
401 
402  // PROJECT HESSIAN TO POSITIVE DEFINITE CONE IF NECESSARY:
403  // -------------------------------------------------------
404  int hessianMode;
405  get( HESSIAN_APPROXIMATION,hessianMode );
406 
407  if ( (HessianApproximationMode)hessianMode == EXACT_HESSIAN )
408  {
409  double hessianProjectionFactor;
410  get( HESSIAN_PROJECTION_FACTOR, hessianProjectionFactor );
411  projectHessian( denseCP.H, hessianProjectionFactor );
412  }
413 
414  // APPLY LEVENBERG-MARQUARD REGULARISATION IF DESIRED:
415  // -------------------------------------------------------
416  double levenbergMarquard;
417  get(LEVENBERG_MARQUARDT, levenbergMarquard );
418 
419  if( levenbergMarquard > EPS )
420  denseCP.H += eye<double>( denseCP.H.rows() ) * levenbergMarquard;
421 
423  {
424  // Check condition number of the condensed Hessian.
425  double denseHConditionNumber = denseCP.H.getConditionNumber();
426  if (denseHConditionNumber > 1.0e16)
427  {
428  LOG( LVL_WARNING )
429  << "Condition number of the condensed Hessian is quite high: log_10(kappa( H )) = "
430  << log10( denseHConditionNumber ) << endl;
432  }
433  // Check for max and min entry in the condensed Hessian:
434  if (denseCP.H.getMin() < -1.0e16 || denseCP.H.getMax() > 1.0e16)
435  {
436  LOG( LVL_WARNING ) << "Ill formed condensed Hessian: min(.) < -1e16 or max(.) > 1e16" << endl;
438  }
439  }
440 
441  // SOLVE QP ALLOWING THE GIVEN NUMBER OF ITERATIONS:
442  // -------------------------------------------------------
443  int maxQPiter;
444  get( MAX_NUM_QP_ITERATIONS, maxQPiter );
445 
446 
447  RealClock clock;
448  clock.start( );
449 
450 // denseCP.print( "",PS_MATLAB );
451  returnvalue = solveQP( maxQPiter );
452 
453  clock.stop( );
454  setLast( LOG_TIME_QP,clock.getTime() );
455 
456  switch( returnvalue )
457  {
458  case SUCCESSFUL_RETURN:
461  break;
462 
467  break;
468 
469  case RET_QP_UNBOUNDED:
472  return ACADOERROR( RET_QP_UNBOUNDED );
473  break;
474 
475  default: //case: RET_QP_INFEASIBLE:
476  int infeasibleQPhandling;
477  get( INFEASIBLE_QP_HANDLING,infeasibleQPhandling );
478 
479  switch( (InfeasibleQPhandling) infeasibleQPhandling )
480  {
481 
482  case IQH_STOP:
483  return ACADOERROR( RET_QP_INFEASIBLE );
484 
485  case IQH_IGNORE:
488  break;
489 
490  default:
491  clock.reset( );
492  clock.start( );
493 
495 
496  if ( solveQP( maxQPiter - cpSolver->getNumberOfIterations( ),
497  (InfeasibleQPhandling) infeasibleQPhandling
498  ) != SUCCESSFUL_RETURN )
500 
501  clock.stop( );
504  }
505  break;
506  }
507 
508  return SUCCESSFUL_RETURN;
509 }
510 
511 
512 
514  )
515 {
518 
519  uint run1, run2, run3;
520 
521  uint nF = getNF();
522  uint nA = getNA();
523  uint N = getNumPoints();
524 
525 
526  if( getNX() + getNXA() > 0 )
527  {
530 
531  int rowOffset = 0;
532  uint rowOffset1 = 0;
533 
534  if ( condensingStatus != COS_FROZEN )
535  {
536  // generate H
537  hT = cp.hessian*T;
538  HDense = T^hT;
539 
540  if( getNX() != 0 ) generateHessianBlockLine( getNX(), rowOffset, rowOffset1 );
541  rowOffset++;
542 
543  for( run3 = 0; run3 < N; run3++ ){
544  if( getNXA() != 0 ) generateHessianBlockLine( getNXA(), rowOffset, rowOffset1 );
545  rowOffset++;
546  }
547 
548  if( getNP() != 0 ) generateHessianBlockLine( getNP(), rowOffset, rowOffset1 );
549  rowOffset++;
550 
551  for( run3 = 0; run3 < N-1; run3++ ){
552  if( getNU() != 0 ) generateHessianBlockLine( getNU(), rowOffset, rowOffset1 );
553  rowOffset++;
554  }
555 
556  for( run3 = 0; run3 < N-1; run3++ ){
557  if( getNW() != 0 ) generateHessianBlockLine( getNW(), rowOffset, rowOffset1 );
558  rowOffset++;
559  }
560 
561 
562  // generate A
564 
565  denseCP.A.setZero();
566 
567  rowOffset1 = 0;
568  for( run3 = 0; run3 < ADense.getNumRows(); run3++ ){
569  ASSERT( getNC() != 0 );
570  generateConstraintBlockLine( (int) blockDims(run3), run3, rowOffset1 );
571  rowOffset++;
572  }
573 
574  for( run3 = 1; run3 < N; run3++ ){
575  generateStateBoundBlockLine( getNX(), run3, rowOffset1 );
576  }
577  }
578 
579 
580  // generate g
581  gDense = (cp.objectiveGradient*T) + (d^hT);
582 
584 
585 
586  // generate lb, ub
587  BlockMatrix dCut(4*N+1,1);
588  DMatrix tmp;
589  for( run1 = 0; run1 < N; run1++ ){
590  d.getSubBlock(run1,0,tmp);
591  if( tmp.getDim() != 0 )
592  dCut.setDense(run1,0,tmp);
593  }
594 
595  lbDense = cp.lowerBoundResiduum - dCut;
596  ubDense = cp.upperBoundResiduum - dCut;
597 
599 
600 
601  // generate lbA, ubA
604 
605  rowOffset1 = 0;
606  for( run3 = 0; run3 < ADense.getNumRows(); run3++ ){
607  ASSERT( getNC() != 0 );
608  generateConstraintVectors( (int) blockDims(run3), run3, rowOffset1 );
609  rowOffset++;
610  }
611 
612  for( run3 = 1; run3 < N; run3++ ){
613  generateStateBoundVectors( getNX(), run3, rowOffset1 );
614  }
615 
616 
617 // acadoPrintf("denseCP.H = \n");
618 // denseCP.H.print();
619 //
620 // acadoPrintf("denseCP.A = \n");
621 // denseCP.A.print();
622 //
623 // acadoPrintf("denseCP.g = \n");
624 // denseCP.g.print();
625 //
626 // acadoPrintf("denseCP.lbA = \n");
627 // denseCP.lbA.print();
628 //
629 // acadoPrintf("denseCP.ubA = \n");
630 // denseCP.ubA.print();
631 //
632 // acadoPrintf("denseCP.lb = \n");
633 // denseCP.lb.print();
634 //
635 // acadoPrintf("denseCP.ub = \n");
636 // denseCP.ub.print();
637  }
638  else{
639 
640  ASSERT( N == 1 );
641  ASSERT( getNU() == 0 );
642  ASSERT( getNW() == 0 );
643 
644  denseCP.H.init( nF, nF );
645  denseCP.g.init( nF );
646 
647  denseCP.A.init( nA, nF );
648  denseCP.lbA.init( nA );
649  denseCP.ubA.init( nA );
650 
651  denseCP.lb.init( nF );
652  denseCP.ub.init( nF );
653 
654  DMatrix tmp;
655 
656  cp.hessian .getSubBlock( 2, 2, denseCP.H , getNP(), getNP() );
657  cp.objectiveGradient .getSubBlock( 0, 2, tmp, 1 , getNP() );
658  for( run1 = 0; run1 < getNP(); run1++ )
659  denseCP.g(run1) = tmp(0,run1);
660 
661  uint offset = 0;
662  for( run2 = 0; run2 < cp.constraintGradient.getNumRows(); run2++ ){
663 
664  cp.constraintGradient.getSubBlock( run2, 2, tmp );
665 
666  for( run1 = 0; run1 < tmp.getNumRows(); run1++ )
667  for( run3 = 0; run3 < tmp.getNumCols(); run3++ )
668  denseCP.A(offset + run1,run3) = tmp(run1,run3);
669 
670  cp.lowerConstraintResiduum.getSubBlock( run2, 0, tmp );
671  for( run1 = 0; run1 < tmp.getNumRows(); run1++ )
672  denseCP.lbA(offset + run1) = tmp(run1,0);
673 
674  cp.upperConstraintResiduum.getSubBlock( run2, 0, tmp );
675  for( run1 = 0; run1 < tmp.getNumRows(); run1++ )
676  denseCP.ubA(offset + run1) = tmp(run1,0);
677 
678  offset += tmp.getNumRows();
679  }
680 
681  cp.lowerBoundResiduum.getSubBlock( 2, 0, tmp, getNP(), 1 );
682  for( run1 = 0; run1 < getNP(); run1++ )
683  denseCP.lb(run1) = tmp(run1,0);
684  cp.upperBoundResiduum.getSubBlock( 2, 0, tmp, getNP(), 1 );
685  for( run1 = 0; run1 < getNP(); run1++ )
686  denseCP.ub(run1) = tmp(run1,0);
687  }
688 
689  if ( condensingStatus != COS_FROZEN )
691 
692  return SUCCESSFUL_RETURN;
693 }
694 
695 
696 
698 
699  uint run1, run2, run3;
700 
701  uint colOffset = 0;
702  uint colOffset1 = 0;
703 
704  uint N = getNumPoints();
705 
706  DMatrix tmp;
707 
708  if( getNX() != 0 ){
709  HDense.getSubBlock( rowOffset, 0, tmp, nn, getNX() );
710 
711  for( run1 = 0; run1 < nn; run1++ )
712  for( run2 = 0; run2 < getNX(); run2++ )
713  denseCP.H(rowOffset1+run1,run2) = tmp(run1,run2);
714  colOffset1 += getNX();
715  }
716  colOffset++;
717 
718  for( run3 = 0; run3 < N; run3++ ){
719  if( getNXA() != 0 ){
720  HDense.getSubBlock( rowOffset, colOffset, tmp, nn, getNXA() );
721  for( run1 = 0; run1 < nn; run1++ )
722  for( run2 = 0; run2 < getNXA(); run2++ )
723  denseCP.H(rowOffset1+run1,colOffset1+run2) = tmp(run1,run2);
724  colOffset1 += getNXA();
725  }
726  colOffset ++;
727  }
728 
729  if( getNP() != 0 ){
730  HDense.getSubBlock( rowOffset, colOffset, tmp, nn, getNP() );
731  for( run1 = 0; run1 < nn; run1++ )
732  for( run2 = 0; run2 < getNP(); run2++ )
733  denseCP.H(rowOffset1+run1,colOffset1+run2) = tmp(run1,run2);
734  colOffset1 += getNP();
735  }
736  colOffset ++;
737 
738  for( run3 = 0; run3 < N-1; run3++ ){
739  if( getNU() != 0 ){
740  HDense.getSubBlock( rowOffset, colOffset, tmp, nn, getNU() );
741  for( run1 = 0; run1 < nn; run1++ )
742  for( run2 = 0; run2 < getNU(); run2++ )
743  denseCP.H(rowOffset1+run1,colOffset1+run2) = tmp(run1,run2);
744  colOffset1 += getNU();
745  }
746  colOffset ++;
747  }
748 
749  for( run3 = 0; run3 < N-1; run3++ ){
750  if( getNW() != 0 ){
751  HDense.getSubBlock( rowOffset, colOffset, tmp, nn, getNW() );
752  for( run1 = 0; run1 < nn; run1++ )
753  for( run2 = 0; run2 < getNW(); run2++ )
754  denseCP.H(rowOffset1+run1,colOffset1+run2) = tmp(run1,run2);
755  colOffset1 += getNW();
756  }
757  colOffset ++;
758  }
759  rowOffset1 += nn;
760 
761  return SUCCESSFUL_RETURN;
762 }
763 
764 
766 
767 
768  uint run1, run2, run3;
769 
770  uint colOffset = 0;
771  uint colOffset1 = 0;
772 
773  uint N = getNumPoints();
774 
775  DMatrix tmp;
776 
777  if( getNX() != 0 ){
778  ADense.getSubBlock( rowOffset, 0, tmp, nn, getNX() );
779  for( run1 = 0; run1 < nn; run1++ )
780  for( run2 = 0; run2 < getNX(); run2++ )
781  denseCP.A(rowOffset1+run1,run2) = tmp(run1,run2);
782  colOffset1 += getNX();
783  }
784  colOffset++;
785 
786  for( run3 = 0; run3 < N; run3++ ){
787  if( getNXA() != 0 ){
788  ADense.getSubBlock( rowOffset, colOffset, tmp, nn, getNXA() );
789  for( run1 = 0; run1 < nn; run1++ )
790  for( run2 = 0; run2 < getNXA(); run2++ )
791  denseCP.A(rowOffset1+run1,colOffset1+run2) = tmp(run1,run2);
792  colOffset1 += getNXA();
793  }
794  colOffset ++;
795  }
796 
797  if( getNP() != 0 ){
798  ADense.getSubBlock( rowOffset, colOffset, tmp, nn, getNP() );
799  for( run1 = 0; run1 < nn; run1++ )
800  for( run2 = 0; run2 < getNP(); run2++ )
801  denseCP.A(rowOffset1+run1,colOffset1+run2) = tmp(run1,run2);
802  colOffset1 += getNP();
803  }
804  colOffset ++;
805 
806  for( run3 = 0; run3 < N-1; run3++ ){
807  if( getNU() != 0 ){
808  ADense.getSubBlock( rowOffset, colOffset, tmp, nn, getNU() );
809  for( run1 = 0; run1 < nn; run1++ )
810  for( run2 = 0; run2 < getNU(); run2++ )
811  denseCP.A(rowOffset1+run1,colOffset1+run2) = tmp(run1,run2);
812  colOffset1 += getNU();
813  }
814  colOffset ++;
815  }
816 
817  for( run3 = 0; run3 < N-1; run3++ ){
818  if( getNW() != 0 ){
819  ADense.getSubBlock( rowOffset, colOffset, tmp, nn, getNW() );
820  for( run1 = 0; run1 < nn; run1++ )
821  for( run2 = 0; run2 < getNW(); run2++ )
822  denseCP.A(rowOffset1+run1,colOffset1+run2) = tmp(run1,run2);
823  colOffset1 += getNW();
824  }
825  colOffset ++;
826  }
827  rowOffset1 += nn;
828 
829  return SUCCESSFUL_RETURN;
830 }
831 
832 
834 
835 
836  uint run1, run2, run3;
837 
838  uint colOffset = 0;
839  uint colOffset1 = 0;
840 
841  uint N = getNumPoints();
842 
843  DMatrix tmp;
844 
845  if( getNX() != 0 ){
846  T.getSubBlock( rowOffset, 0, tmp, nn, getNX() );
847  for( run1 = 0; run1 < nn; run1++ )
848  for( run2 = 0; run2 < getNX(); run2++ )
849  denseCP.A(rowOffset1+run1,run2) = tmp(run1,run2);
850  colOffset1 += getNX();
851  }
852  colOffset++;
853 
854  for( run3 = 0; run3 < N; run3++ ){
855  if( getNXA() != 0 ){
856  T.getSubBlock( rowOffset, colOffset, tmp, nn, getNXA() );
857  for( run1 = 0; run1 < nn; run1++ )
858  for( run2 = 0; run2 < getNXA(); run2++ )
859  denseCP.A(rowOffset1+run1,colOffset1+run2) = tmp(run1,run2);
860  colOffset1 += getNXA();
861  }
862  colOffset ++;
863  }
864 
865  if( getNP() != 0 ){
866  T.getSubBlock( rowOffset, colOffset, tmp, nn, getNP() );
867  for( run1 = 0; run1 < nn; run1++ )
868  for( run2 = 0; run2 < getNP(); run2++ )
869  denseCP.A(rowOffset1+run1,colOffset1+run2) = tmp(run1,run2);
870  colOffset1 += getNP();
871  }
872  colOffset ++;
873 
874  for( run3 = 0; run3 < N-1; run3++ ){
875  if( getNU() != 0 ){
876  T.getSubBlock( rowOffset, colOffset, tmp, nn, getNU() );
877  for( run1 = 0; run1 < nn; run1++ )
878  for( run2 = 0; run2 < getNU(); run2++ )
879  denseCP.A(rowOffset1+run1,colOffset1+run2) = tmp(run1,run2);
880  colOffset1 += getNU();
881  }
882  colOffset ++;
883  }
884 
885  for( run3 = 0; run3 < N-1; run3++ ){
886  if( getNW() != 0 ){
887  T.getSubBlock( rowOffset, colOffset, tmp, nn, getNW() );
888  for( run1 = 0; run1 < nn; run1++ )
889  for( run2 = 0; run2 < getNW(); run2++ )
890  denseCP.A(rowOffset1+run1,colOffset1+run2) = tmp(run1,run2);
891  colOffset1 += getNW();
892  }
893  colOffset ++;
894  }
895  rowOffset1 += nn;
896 
897  return SUCCESSFUL_RETURN;
898 }
899 
900 
902 
903  uint run1;
904 
905  DMatrix tmp;
906 
907  lbADense.getSubBlock( rowOffset, 0, tmp, nn, 1 );
908  for( run1 = 0; run1 < nn; run1++ )
909  denseCP.lbA(rowOffset1+run1) = tmp(run1,0);
910 
911  ubADense.getSubBlock( rowOffset, 0, tmp, nn, 1 );
912  for( run1 = 0; run1 < nn; run1++ )
913  denseCP.ubA(rowOffset1+run1) = tmp(run1,0);
914 
915  rowOffset1 += nn;
916 
917  return SUCCESSFUL_RETURN;
918 }
919 
920 
922 
923  uint run1;
924 
925  DMatrix tmp;
926 
927  lbDense.getSubBlock( rowOffset, 0, tmp, nn, 1 );
928  for( run1 = 0; run1 < nn; run1++ )
929  denseCP.lbA(rowOffset1+run1) = tmp(run1,0);
930 
931  ubDense.getSubBlock( rowOffset, 0, tmp, nn, 1 );
932  for( run1 = 0; run1 < nn; run1++ )
933  denseCP.ubA(rowOffset1+run1) = tmp(run1,0);
934 
935  rowOffset1 += nn;
936 
937  return SUCCESSFUL_RETURN;
938 }
939 
940 
942 
943  uint run1, run3;
944 
945  uint N = getNumPoints();
946 
947  uint rowOffset = N;
948  uint rowOffset1 = 0;
949 
950  DMatrix tmp;
951 
952  if( getNX() != 0 ){
953  lbDense.getSubBlock( 0, 0, tmp, getNX(), 1 );
954  for( run1 = 0; run1 < getNX(); run1++ )
955  denseCP.lb(rowOffset1+run1) = tmp(run1,0);
956  ubDense.getSubBlock( 0, 0, tmp, getNX(), 1 );
957  for( run1 = 0; run1 < getNX(); run1++ )
958  denseCP.ub(rowOffset1+run1) = tmp(run1,0);
959  rowOffset1 += getNX();
960  }
961 
962  for( run3 = 0; run3 < N; run3++ ){
963  if( getNXA() != 0 ){
964  lbDense.getSubBlock( rowOffset, 0, tmp, getNXA(), 1 );
965  for( run1 = 0; run1 < getNXA(); run1++ )
966  denseCP.lb(rowOffset1+run1) = tmp(run1,0);
967  ubDense.getSubBlock( rowOffset, 0, tmp, getNXA(), 1 );
968  for( run1 = 0; run1 < getNXA(); run1++ )
969  denseCP.ub(rowOffset1+run1) = tmp(run1,0);
970  rowOffset1 += getNXA();
971  }
972  rowOffset++;
973  }
974 
975  if( getNP() != 0 ){
976  lbDense.getSubBlock( rowOffset, 0, tmp, getNP(), 1 );
977  for( run1 = 0; run1 < getNP(); run1++ )
978  denseCP.lb(rowOffset1+run1) = tmp(run1,0);
979  ubDense.getSubBlock( rowOffset, 0, tmp, getNP(), 1 );
980  for( run1 = 0; run1 < getNP(); run1++ )
981  denseCP.ub(rowOffset1+run1) = tmp(run1,0);
982  rowOffset1 += getNP();
983  }
984  rowOffset++;
985 
986  for( run3 = 0; run3 < N-1; run3++ ){
987  if( getNU() != 0 ){
988  lbDense.getSubBlock( rowOffset, 0, tmp, getNU(), 1 );
989  for( run1 = 0; run1 < getNU(); run1++ )
990  denseCP.lb(rowOffset1+run1) = tmp(run1,0);
991  ubDense.getSubBlock( rowOffset, 0, tmp, getNU(), 1 );
992  for( run1 = 0; run1 < getNU(); run1++ )
993  denseCP.ub(rowOffset1+run1) = tmp(run1,0);
994  rowOffset1 += getNU();
995  }
996  rowOffset++;
997  }
998  rowOffset++;
999 
1000  for( run3 = 0; run3 < N-1; run3++ ){
1001  if( getNW() != 0 ){
1002  lbDense.getSubBlock( rowOffset, 0, tmp, getNW(), 1 );
1003  for( run1 = 0; run1 < getNW(); run1++ )
1004  denseCP.lb(rowOffset1+run1) = tmp(run1,0);
1005  ubDense.getSubBlock( rowOffset, 0, tmp, getNW(), 1 );
1006  for( run1 = 0; run1 < getNW(); run1++ )
1007  denseCP.ub(rowOffset1+run1) = tmp(run1,0);
1008  rowOffset1 += getNW();
1009  }
1010  rowOffset ++;
1011  }
1012 
1013  return SUCCESSFUL_RETURN;
1014 }
1015 
1016 
1018 
1019  uint run1, run3;
1020 
1021  uint rowOffset = 0;
1022  uint rowOffset1 = 0;
1023 
1024  uint N = getNumPoints();
1025 
1026  DMatrix tmp;
1027 
1028  if( getNX() != 0 ){
1029  gDense.getSubBlock( 0, 0, tmp, 1, getNX() );
1030  for( run1 = 0; run1 < getNX(); run1++ )
1031  denseCP.g(rowOffset1+run1) = tmp(0,run1);
1032  rowOffset1 += getNX();
1033  }
1034  rowOffset++;
1035 
1036  for( run3 = 0; run3 < N; run3++ ){
1037  if( getNXA() != 0 ){
1038  gDense.getSubBlock( 0, rowOffset, tmp, 1, getNXA() );
1039  for( run1 = 0; run1 < getNXA(); run1++ )
1040  denseCP.g(rowOffset1+run1) = tmp(0,run1);
1041  rowOffset1 += getNXA();
1042  }
1043  rowOffset++;
1044  }
1045 
1046  if( getNP() != 0 ){
1047  gDense.getSubBlock( 0, rowOffset, tmp, 1, getNP() );
1048  for( run1 = 0; run1 < getNP(); run1++ )
1049  denseCP.g(rowOffset1+run1) = tmp(0,run1);
1050  rowOffset1 += getNP();
1051  }
1052  rowOffset++;
1053 
1054  for( run3 = 0; run3 < N-1; run3++ ){
1055  if( getNU() != 0 ){
1056  gDense.getSubBlock( 0, rowOffset, tmp, 1, getNU() );
1057  for( run1 = 0; run1 < getNU(); run1++ )
1058  denseCP.g(rowOffset1+run1) = tmp(0,run1);
1059  rowOffset1 += getNU();
1060  }
1061  rowOffset++;
1062  }
1063 
1064  for( run3 = 0; run3 < N-1; run3++ ){
1065  if( getNW() != 0 ){
1066  gDense.getSubBlock( 0, rowOffset, tmp, 1, getNW() );
1067  for( run1 = 0; run1 < getNW(); run1++ )
1068  denseCP.g(rowOffset1+run1) = tmp(0,run1);
1069  rowOffset1 += getNW();
1070  }
1071  rowOffset ++;
1072  }
1073 
1074  return SUCCESSFUL_RETURN;
1075 }
1076 
1077 
1078 
1080  )
1081 {
1082  // do not expand if condensing is frozen
1083 // if ( condensingStatus == COS_FROZEN )
1084 // return SUCCESSFUL_RETURN;
1085 
1087  return ACADOERROR( RET_UNABLE_TO_EXPAND );
1088 
1089 
1090  uint run1,run2;
1091 
1092  int rowCount = 0;
1093  int rowCount1 = 0;
1094 
1095  uint nF = getNF();
1096  uint N = getNumPoints();
1097 
1098 
1099  DVector denseDualSolution( denseCP.getMergedDualSolution( ) );
1100 
1101  if( getNX() != 0 ){
1102 
1103  BlockMatrix primalDense;
1104  primalDense.init( 3*N, 1 );
1105 
1106  DMatrix tmp;
1107 
1108  rowCount = 0;
1109  rowCount1 = 0;
1110 
1111  tmp.init( getNX(), 1 );
1112  for( run1 = 0; run1 < getNX(); run1++ )
1113  tmp(run1,0) = (*denseCP.x)(run1);
1114  primalDense.setDense( 0, 0, tmp );
1115  rowCount1 += getNX();
1116  rowCount++;
1117 
1118  for( run2 = 0; run2 < N; run2++ ){
1119  if( getNXA() != 0 ){
1120  tmp.init( getNXA(), 1 );
1121  for( run1 = 0; run1 < getNXA(); run1++ )
1122  tmp(run1,0) = (*denseCP.x)(rowCount1+run1);
1123  primalDense.setDense( rowCount, 0, tmp );
1124  rowCount1 += getNXA();
1125  }
1126  rowCount++;
1127  }
1128 
1129  if( getNP() != 0 ){
1130  tmp.init( getNP(), 1 );
1131  for( run1 = 0; run1 < getNP(); run1++ )
1132  tmp(run1,0) = (*denseCP.x)(rowCount1+run1);
1133  primalDense.setDense( rowCount, 0, tmp );
1134  rowCount1 += getNP();
1135  }
1136  rowCount++;
1137 
1138  for( run2 = 0; run2 < N-1; run2++ ){
1139  if( getNU() != 0 ){
1140  tmp.init( getNU(), 1 );
1141  for( run1 = 0; run1 < getNU(); run1++ )
1142  tmp(run1,0) = (*denseCP.x)(rowCount1+run1);
1143  primalDense.setDense( rowCount, 0, tmp );
1144  rowCount1 += getNU();
1145  }
1146  rowCount++;
1147  }
1148 
1149  for( run2 = 0; run2 < N-1; run2++ ){
1150  if( getNW() != 0 ){
1151  tmp.init( getNW(), 1 );
1152  for( run1 = 0; run1 < getNW(); run1++ )
1153  tmp(run1,0) = (*denseCP.x)(rowCount1+run1);
1154  primalDense.setDense( rowCount, 0, tmp );
1155  rowCount1 += getNW();
1156  }
1157  rowCount++;
1158  }
1159 
1160  cp.deltaX = T*primalDense + d;
1161 
1162  BlockMatrix aux;
1163 
1164  aux = (cp.deltaX^cp.hessian) + cp.objectiveGradient;
1165 
1166  DVector aux2(N*getNX());
1167  aux2.setZero();
1168 
1169  int run = 0;
1170  uint run3, run4;
1171 
1172  for( run1 = 0; run1 < (uint)cp.constraintGradient.getNumRows(); run1++ ){
1173  for( run2 = 0; run2 < N; run2++ ){
1174  cp.constraintGradient.getSubBlock( run1, run2, tmp );
1175  for( run3 = 0; run3 < (uint)tmp.getNumRows(); run3++ ){
1176  for( run4 = 0; run4 < getNX(); run4++ ){
1177  aux2(run2*getNX()+run4) += denseDualSolution(nF+run+run3)*tmp(run3,run4);
1178  }
1179  }
1180  run += tmp.getNumRows();
1181  }
1182  }
1183 
1184  DVector aux3(N*getNX());
1185 
1186  for( run2 = 0; run2 < getNX(); run2++ )
1187  aux3(run2) = denseDualSolution(run2);
1188 
1189  run = nF;
1190  run += getNC();
1191 
1192  for( run1 = 0; run1 < N-1; run1++ )
1193  for( run2 = 0; run2 < getNX(); run2++ )
1194  aux3((run1+1)*getNX()+run2) = denseDualSolution(run+run1*getNX()+run2);
1195 
1196  DVector *aux4 = new DVector[N-1];
1197 
1198  // aux = x^T denseCP.H + denseCP.g (BlockMatrix)
1199  // aux2 = lambda^T denseCP.A (DVector )
1200  // aux3 = lambda_bound (DVector )
1201 
1202  for( run1 = 0; run1 < N-1; run1++ ){
1203  aux4[run1].init(getNX());
1204  aux.getSubBlock( 0, run1+1, tmp );
1205 
1206  for( run2 = 0; run2 < getNX(); run2++ ){
1207  aux4[run1](run2) = tmp(0,run2) - aux2( (run1+1)*getNX()+run2 ) - aux3( (run1+1)*getNX()+run2 );
1208  }
1209  }
1210 
1211  // aux4[...] = x^T denseCP.H + denseCP.g - lambda^T denseCP.A - lambda_bound
1212 
1213  DMatrix Gx;
1214  DVector *lambdaDyn;
1215  lambdaDyn = new DVector[N-1];
1216  lambdaDyn[N-2] = aux4[N-2];
1217 
1218  for( run1 = N-2; run1 >= 1; run1-- ){
1219  cp.dynGradient.getSubBlock( run1, 0, Gx );
1220  lambdaDyn[run1-1] = (Gx.transpose() * lambdaDyn[run1]) + aux4[run1-1];
1221  }
1222 
1223  cp.lambdaDynamic.init( N-1, 1 );
1224 
1225  for( run1 = 0; run1 < N-1; run1++ ){
1226  tmp.init( lambdaDyn[run1].getDim(), 1 );
1227  for( run2 = 0; run2 < lambdaDyn[run1].getDim(); run2++ )
1228  tmp(run2,0) = lambdaDyn[run1].operator()(run2);
1229  cp.lambdaDynamic.setDense( run1, 0, tmp );
1230  }
1231 
1232 
1233  int dynMode;
1234  get( DYNAMIC_SENSITIVITY, dynMode );
1235 
1236  if( dynMode == FORWARD_SENSITIVITY_LIFTED ){
1237  for( run1 = 0; run1 < N-1; run1++ ){
1238  tmp.init( lambdaDyn[run1].getDim(), 1 );
1239  tmp.setAll(1.0/((double) lambdaDyn[run1].getDim()));
1240  cp.lambdaDynamic.setDense( run1, 0, tmp );
1241  }
1242  }
1243 
1244  delete[] lambdaDyn;
1245  delete[] aux4;
1246  }
1247  else{
1248  DMatrix tmp ( getNP(),1 );
1249  cp.deltaX.init( 5, 1 );
1250 
1251  for( run1 = 0; run1 < getNP(); run1++ )
1252  tmp( run1, 0 ) = (*denseCP.x)(run1);
1253 
1254  cp.deltaX.setDense( 2, 0, tmp );
1255  }
1256 
1257 
1258  DMatrix tmp;
1259 
1261 
1262  int run = 0;
1263  for( run1 = 0; run1 < blockDims.getDim(); run1++ ){
1264  tmp.init( (int) blockDims(run1), 1 );
1265 
1266  for( run2 = 0; run2 < blockDims(run1); run2++ ){
1267  tmp(run2,0) = denseDualSolution(nF+run+run2);
1268  }
1269  cp.lambdaConstraint.setDense( run1, 0, tmp );
1270  run += (int) blockDims(run1);
1271  }
1272 
1273  cp.lambdaBound.init(4*N+1,1);
1274  rowCount = 0;
1275 
1276  if( getNX() != 0 ){
1277  tmp.init(getNX(),1);
1278  for( run1 = 0; run1 < getNX(); run1++ )
1279  tmp(run1,0) = denseDualSolution(run1);
1280  cp.lambdaBound.setDense( 0, 0, tmp );
1281  rowCount++;
1282  }
1283  rowCount1 = 0;
1284 
1285  for( run2 = 1; run2 < N; run2++ ){
1286  if( getNX() != 0 ){
1287  tmp.init( getNX(), 1 );
1288  for( run1 = 0; run1 < getNX(); run1++ )
1289  tmp(run1,0) = denseDualSolution(nF+getNC()+rowCount1+run1);
1290  cp.lambdaBound.setDense( rowCount, 0, tmp );
1291  rowCount1 += getNX();
1292  }
1293  rowCount++;
1294  }
1295  rowCount1 = getNX();
1296 
1297  for( run2 = 0; run2 < N; run2++ ){
1298  if( getNXA() != 0 ){
1299  tmp.init( getNXA(), 1 );
1300  for( run1 = 0; run1 < getNXA(); run1++ )
1301  tmp(run1,0) = denseDualSolution(rowCount1+run1);
1302  cp.lambdaBound.setDense( rowCount, 0, tmp );
1303  rowCount1 += getNXA();
1304  }
1305  rowCount++;
1306  }
1307 
1308  if( getNP() != 0 ){
1309  tmp.init( getNP(), 1 );
1310  for( run1 = 0; run1 < getNP(); run1++ )
1311  tmp(run1,0) = denseDualSolution(rowCount1+run1);
1312  cp.lambdaBound.setDense( rowCount, 0, tmp );
1313  rowCount1 += getNP();
1314  }
1315  rowCount++;
1316 
1317  for( run2 = 0; run2 < N-1; run2++ ){
1318  if( getNU() != 0 ){
1319  tmp.init( getNU(), 1 );
1320  for( run1 = 0; run1 < getNU(); run1++ )
1321  tmp(run1,0) = denseDualSolution(rowCount1+run1);
1322  cp.lambdaBound.setDense( rowCount, 0, tmp );
1323  rowCount1 += getNU();
1324  }
1325  rowCount++;
1326  }
1327  rowCount++;
1328 
1329  for( run2 = 0; run2 < N-1; run2++ ){
1330  if( getNW() != 0 ){
1331  tmp.init( getNW(), 1 );
1332  for( run1 = 0; run1 < getNW(); run1++ )
1333  tmp(run1,0) = denseDualSolution(rowCount1+run1);
1334  cp.lambdaBound.setDense( rowCount, 0, tmp );
1335  rowCount1 += getNW();
1336  }
1337  rowCount++;
1338  }
1339 
1340  if ( condensingStatus != COS_FROZEN )
1342 
1343  return SUCCESSFUL_RETURN;
1344 }
1345 
1346 
1347 
1349 {
1350  uint N = getNumPoints();
1351 
1352  T.init( 5*N, 3*N );
1353  d.init( 5*N, 1 );
1354 
1355  T.setZero();
1356  d.setZero();
1357 
1358  if( getNX() != 0 ) T.setIdentity( 0, 0, getNX() );
1359  for( uint run1 = 0; run1 < N; run1++ ){
1360  if( getNXA() != 0 ) T.setIdentity( N+run1, 1+run1, getNXA() );
1361  if( getNP() != 0 ) T.setIdentity( 2*N+run1, N+1 , getNP() );
1362  if( getNU() != 0 ){
1363  if( run1 != N-1 ) T.setIdentity( 3*N+run1, N+2+run1, getNU() );
1364  else T.setIdentity( 3*N+run1, N+1+run1, getNU() );
1365  }
1366  if( getNW() != 0 ){
1367  if( run1 != N-1 ) T.setIdentity( 4*N+run1, 2*N+1+run1, getNW() );
1368  else T.setIdentity( 4*N+run1, 2*N +run1, getNW() );
1369  }
1370  }
1371 
1373 
1374  return SUCCESSFUL_RETURN;
1375 }
1376 
1377 
1379  )
1380 {
1381  uint run1, run2;
1382  uint N = getNumPoints();
1383 
1384  DMatrix Gx;
1385  DMatrix G;
1386  DMatrix tmp;
1387 
1388  for( run1 = 0; run1 < N-1; run1++ )
1389  {
1390  // DIFFERENTIAL STATES:
1391  // --------------------
1392  cp.dynGradient.getSubBlock( run1, 0, Gx ); // Get the sensitivity G_x^i with respect to x
1393 
1394  if ( condensingStatus != COS_FROZEN )
1395  {
1396  T .getSubBlock( run1, 0, tmp ); // get the corresponding C_i .
1397 
1398  T.setDense( run1+1, 0, Gx*tmp ); // compute C_{i+1} := G_x^i * C_i
1399 
1400  // ALGEBRAIC STATES:
1401  // --------------------
1402 
1403  cp.dynGradient.getSubBlock( run1, 1, G );
1404 
1405  for( run2 = 0; run2 <= run1; run2++ ){
1406 
1407  T.getSubBlock( run1, run2+1, tmp );
1408 
1409  if( G.getDim() != 0 ){
1410 
1411  if( run1 == run2 ) T.setDense( run1+1, run2+1, G );
1412  else T.setDense( run1+1, run2+1, Gx*tmp );
1413  }
1414  }
1415 
1416  // PARAMETERS:
1417  // --------------------
1418 
1419  cp.dynGradient.getSubBlock( run1, 2 , G ); // Get the sensitivity G_p^i with respect to p
1420  T .getSubBlock( run1, N+1, tmp ); // get the corresponding D_p^i.
1421 
1422  if( tmp.getDim() != 0 ){
1423  if( G.getDim() != 0 )
1424  T.setDense( run1+1, N+1, Gx*tmp + G ); // compute D_p^{i+1} := G_x^i D_p^i + G_p^i
1425  }
1426  else{
1427  if( G.getDim() != 0 )
1428  T.setDense( run1+1, N+1, G ); // compute D_p^{i+1} := G_x^i D_p^i + G_p^i
1429  }
1430 
1431  // CONTROLS:
1432  // --------------------
1433 
1434  cp.dynGradient.getSubBlock( run1, 3, G );
1435  for( run2 = 0; run2 <= run1; run2++ ){
1436 
1437  T.getSubBlock( run1, run2+2+N, tmp );
1438 
1439  if( G.getDim() != 0 ){
1440 
1441  if( run1 == run2 ) T.setDense( run1+1, run2+2+N, G );
1442  else T.setDense( run1+1, run2+2+N, Gx*tmp );
1443  }
1444  }
1445 
1446  // DISTURBANCES:
1447  // --------------------
1448 
1449  cp.dynGradient.getSubBlock( run1, 4, G );
1450 
1451  for( run2 = 0; run2 <= run1; run2++ ){
1452 
1453  T.getSubBlock( run1, run2+1+2*N, tmp );
1454 
1455  if( G.getDim() != 0 ){
1456 
1457  if( run1 == run2 ) T.setDense( run1+1, run2+1+2*N, G );
1458  else T.setDense( run1+1, run2+1+2*N, Gx*tmp );
1459  }
1460  }
1461  }
1462 
1463  // RESIDUUM:
1464  // --------------------
1465 
1466  cp.dynResiduum.getSubBlock( run1, 0, G ); // Get the residuum b^i
1467  d .getSubBlock( run1, 0, tmp ); // get the corresponding d^i.
1468 
1469  if( tmp.getDim() != 0 ){
1470  if( G.getDim() != 0 )
1471  d.setDense( run1+1, 0, Gx*tmp + G ); // compute d^{i+1} := G_x^i d^i + b^i
1472  }
1473  else{
1474  if( G.getDim() != 0 )
1475  d.setDense( run1+1, 0, G ); // compute d^{i+1} := b^i
1476  }
1477  }
1478 
1479  return SUCCESSFUL_RETURN;
1480 }
1481 
1482 
1483 
1485  DenseCP& _denseCPrelaxed
1486  ) const
1487 {
1488  switch( infeasibleQPhandling )
1489  {
1490  case IQH_RELAX_L1:
1491  return setupRelaxedQPdataL2( _denseCPrelaxed );
1492 
1493  case IQH_RELAX_L2:
1494  return setupRelaxedQPdataL2( _denseCPrelaxed );
1495 
1496  default:
1498  }
1499 }
1500 
1501 
1503  ) const
1504 {
1506 
1507  uint nV = denseCP.getNV();
1508  uint nC = denseCP.getNC();
1509 
1510 
1511  /* allocate memory for relaxed QP data */
1512  _denseCPrelaxed.H.init( nV+nC,nV+nC );
1513  _denseCPrelaxed.A.init( nC,nV+nC );
1514  _denseCPrelaxed.g.init( nV+nC );
1515  _denseCPrelaxed.lb.init( nV+nC );
1516  _denseCPrelaxed.ub.init( nV+nC );
1517  _denseCPrelaxed.lbA = denseCP.lbA;
1518  _denseCPrelaxed.ubA = denseCP.ubA;
1519 
1520 // acadoPrintf( "frob: %e\n", H.getNorm( MN_FROBENIUS ) );
1521 // acadoPrintf( "trace: %e\n", H.getTrace() );
1522 // acadoPrintf( "min diag: %e\n", (H.getDiag()).getMin() );
1523 // acadoPrintf( "max diag: %e\n", (H.getDiag()).getMax() );
1524 // acadoPrintf( "max(abs(g)): %e\n", (g.getAbsolute()).getMax() );
1525 
1526  _denseCPrelaxed.H.setZero( );
1527  for( uint i=0; i<nV; ++i )
1528  for( uint j=0; j<nV; ++j )
1529  _denseCPrelaxed.H(i,j) = denseCP.H(i,j);
1530 
1531  _denseCPrelaxed.A.setZero( );
1532  for( uint i=0; i<nC; ++i )
1533  for( uint j=0; j<nV; ++j )
1534  _denseCPrelaxed.A(i,j) = denseCP.A(i,j);
1535  for( uint i=0; i<nC; ++i )
1536  _denseCPrelaxed.A(i,nV+i) = 1.0;
1537 
1538  _denseCPrelaxed.g.setZero( );
1539  for( uint i=0; i<nV; ++i )
1540  _denseCPrelaxed.g(i) = denseCP.g(i);
1541  double l1entry = 1.0e5 * ((denseCP.g).getAbsolute()).getMax( );
1542  for( uint i=nV; i<nV+nC; ++i )
1543  _denseCPrelaxed.g(i) = l1entry;
1544 
1545  _denseCPrelaxed.lb.setZero( );
1546  for( uint i=0; i<nV; ++i )
1547  _denseCPrelaxed.lb(i) = denseCP.lb(i);
1548  for( uint i=nV; i<nV+nC; ++i )
1549  _denseCPrelaxed.lb(i) = 0.0;
1550 
1551  _denseCPrelaxed.ub.setZero( );
1552  for( uint i=0; i<nV; ++i )
1553  _denseCPrelaxed.ub(i) = denseCP.ub(i);
1554  for( uint i=nV; i<nV+nC; ++i )
1555  _denseCPrelaxed.ub(i) = INFTY;
1556 
1557  return SUCCESSFUL_RETURN;
1558 }
1559 
1560 
1562  ) const
1563 {
1564  uint nV = denseCP.getNV();
1565  uint nC = denseCP.getNC();
1566 
1567 
1568  /* allocate memory for relaxed QP data */
1569  _denseCPrelaxed.H.init( nV+nC,nV+nC );
1570  _denseCPrelaxed.A.init( nC,nV+nC );
1571  _denseCPrelaxed.g.init( nV+nC );
1572  _denseCPrelaxed.lb.init( nV+nC );
1573  _denseCPrelaxed.ub.init( nV+nC );
1574  _denseCPrelaxed.lbA = denseCP.lbA;
1575  _denseCPrelaxed.ubA = denseCP.ubA;
1576 
1577 // acadoPrintf( "frob: %e\n", H.getNorm( MN_FROBENIUS ) );
1578 // acadoPrintf( "trace: %e\n", H.getTrace() );
1579 // acadoPrintf( "min diag: %e\n", (H.getDiag()).getMin() );
1580 // acadoPrintf( "max diag: %e\n", (H.getDiag()).getMax() );
1581 // acadoPrintf( "max(abs(g)): %e\n", (g.getAbsolute()).getMax() );
1582 
1583  _denseCPrelaxed.H.setZero( );
1584  for( uint i=0; i<nV; ++i )
1585  for( uint j=0; j<nV; ++j )
1586  _denseCPrelaxed.H(i,j) = denseCP.H(i,j);
1587  double diagonalEntry = 1.0e5 * acadoMax( (denseCP.H).getTrace( ),((denseCP.g).getAbsolute()).getMax( ) );
1588  for( uint i=nV; i<nV+nC; ++i )
1589  _denseCPrelaxed.H(i,i) = diagonalEntry;
1590 
1591  _denseCPrelaxed.A.setZero( );
1592  for( uint i=0; i<nC; ++i )
1593  for( uint j=0; j<nV; ++j )
1594  _denseCPrelaxed.A(i,j) = denseCP.A(i,j);
1595  for( uint i=0; i<nC; ++i )
1596  _denseCPrelaxed.A(i,nV+i) = 1.0;
1597 
1598  _denseCPrelaxed.g.setZero( );
1599  for( uint i=0; i<nV; ++i )
1600  _denseCPrelaxed.g(i) = denseCP.g(i);
1601 
1602  _denseCPrelaxed.lb.setZero( );
1603  for( uint i=0; i<nV; ++i )
1604  _denseCPrelaxed.lb(i) = denseCP.lb(i);
1605  for( uint i=nV; i<nV+nC; ++i )
1606  _denseCPrelaxed.lb(i) = -INFTY;
1607 
1608  _denseCPrelaxed.ub.setZero( );
1609  for( uint i=0; i<nV; ++i )
1610  _denseCPrelaxed.ub(i) = denseCP.ub(i);
1611  for( uint i=nV; i<nV+nC; ++i )
1612  _denseCPrelaxed.ub(i) = INFTY;
1613 
1614  return SUCCESSFUL_RETURN;
1615 }
1616 
1617 
1618 
1620  )
1621 {
1622  denseCP.init( getNF(),getNA() );
1623 
1624  denseCP.g.init( getNF() );
1625  denseCP.lb.init( getNF() );
1626  denseCP.ub.init( getNF() );
1627  denseCP.lbA.init( getNA() );
1628  denseCP.ubA.init( getNA() );
1629 
1630  cpSolver->init( &denseCP );
1631 
1632  // setup relaxed QP
1633 // uint nV_relaxed, nC_relaxed;
1634 // getRelaxedQPdimensions( infeasibleQPhandling,nV_relaxed,nC_relaxed );
1635 // cpSolverRelaxed->init( nV_relaxed,nC_relaxed );
1636 
1637  return SUCCESSFUL_RETURN;
1638 }
1639 
1640 
1641 
1643  InfeasibleQPhandling infeasibleQPhandling
1644  )
1645 {
1646 // denseCP.lb.print( "lb" );
1647 // // denseCP.ub.print( "ub" );
1648 // (denseCP.ub - denseCP.lb).print( "ub-lb" );
1649 
1650  if ( infeasibleQPhandling == IQH_UNDEFINED )
1651  return cpSolver->solve( &denseCP );
1652 
1653 
1654  /* relax QP... */
1655  DenseCP denseCPrelaxed;
1656 
1657  if ( setupRelaxedQPdata( infeasibleQPhandling,denseCPrelaxed ) != SUCCESSFUL_RETURN )
1659 
1660  /* ... and solve relaxed QP */
1661  if ( cpSolverRelaxed == 0 )
1663 
1664  if ( ( cpSolverRelaxed->getNumberOfVariables( ) != denseCPrelaxed.getNV() ) ||
1665  ( cpSolverRelaxed->getNumberOfConstraints( ) != denseCPrelaxed.getNC() ) )
1666  {
1667  cpSolverRelaxed->init( denseCPrelaxed.getNV(),denseCPrelaxed.getNC() );
1668  }
1669 
1670  returnValue returnvalue;
1671  returnvalue = cpSolverRelaxed->solve( &denseCPrelaxed );
1672 
1673  if ( ( returnvalue != SUCCESSFUL_RETURN ) && ( returnvalue != RET_QP_SOLUTION_REACHED_LIMIT ) )
1674  return returnvalue;
1675 
1676  const uint nV = denseCP.getNV();
1677  const uint nC = denseCP.getNC();
1678 
1679  DVector deltaDenseTmp;
1680  DVector lambdaDenseTmp;
1681 
1682  DVector deltaDense (nV);
1683  DVector lambdaDense(nV+nC);
1684 
1685  cpSolverRelaxed->getPrimalSolution( deltaDenseTmp );
1686  cpSolverRelaxed->getDualSolution ( lambdaDenseTmp );
1687 
1688  for( uint i=0; i < nV; ++i )
1689  deltaDense(i) = deltaDenseTmp(i);
1690 
1691 // double maxRelaxed = 0.0;
1692 //
1693 // for( uint i=nV; i < nV+nC; ++i )
1694 // if ( fabs(deltaDenseTmp(i)) > maxRelaxed )
1695 // maxRelaxed = fabs(deltaDenseTmp(i));
1696 // printf ("maxRelaxed: %e\n",maxRelaxed);
1697 
1698  for( uint i=0; i < nV; ++i )
1699  lambdaDense(i) = lambdaDenseTmp(i);
1700  for( uint i=0; i < nC; ++i )
1701  lambdaDense(nV+i) = lambdaDenseTmp(denseCPrelaxed.getNV()+i);
1702 
1703  // deltaDense.print( "xOpt" );
1704  denseCP.setQPsolution( deltaDense,lambdaDense );
1705 
1706  return returnvalue;
1707 }
1708 
1709 
1710 
1711 
1713 
1714 // end of file.
int getMax(int x, int y)
returnValue setLast(LogName _name, int lastValue, double time=-INFTY)
#define LOG(level)
Just define a handy macro for getting the logger.
DVector lb
Definition: dense_cp.hpp:384
#define N
HessianApproximationMode
T getMax() const
Definition: matrix.hpp:293
Data class for storing generic optimization variables.
Definition: ocp_iterate.hpp:57
Implements a very rudimentary block sparse matrix class.
returnValue projectHessian(DMatrix &H_, double dampingFactor)
void init(unsigned _nRows=0, unsigned _nCols=0)
Definition: matrix.hpp:135
virtual returnValue initializeCPsolver(InfeasibleQPhandling infeasibleQPhandling)
virtual returnValue finalizeSolve(BandedCP &cp)
virtual returnValue getVarianceCovariance(DMatrix &var)=0
int acadoMax(const int x, const int y)
Definition: acado_utils.cpp:64
Allows real time measurements based on the system&#39;s clock.
Definition: real_clock.hpp:53
BlockMatrix upperConstraintResiduum
Definition: banded_cp.hpp:114
const double INFTY
virtual returnValue getFirstControl(DVector &u0_) const
DVector ub
Definition: dense_cp.hpp:385
BooleanType areRealTimeParametersDefined() const
virtual returnValue reset()
Definition: clock.cpp:86
BlockMatrix deltaX
Definition: banded_cp.hpp:123
#define ACADOINFO(retval)
virtual returnValue unfreezeCondensing()
UserInteraction * userInteraction
BEGIN_NAMESPACE_ACADO const double EPS
returnValue setDense(uint rowIdx, uint colIdx, const DMatrix &value)
void init(unsigned _dim=0)
Definition: vector.hpp:155
BlockMatrix objectiveGradient
Definition: banded_cp.hpp:104
Computes eigenvalues and eigenvectors of selfadjoint matrices.
returnValue expand(BandedCP &cp)
Allows to pass back messages to the calling function.
iterative scaling algorithm to equilibrate rows and column norms in matrices
Definition: matrix.hpp:471
DVector g
Definition: dense_cp.hpp:382
BandedCPsolver & operator=(const BandedCPsolver &rhs)
returnValue generateHessianBlockLine(uint nn, uint rowOffset, uint &rowOffset1)
BEGIN_NAMESPACE_ACADO typedef unsigned int uint
Definition: acado_types.hpp:42
void setAll(const T &_value)
Definition: matrix.hpp:141
BlockMatrix lambdaConstraint
Definition: banded_cp.hpp:126
returnValue setQPsolution(const DVector &x_, const DVector &y_)
Definition: dense_cp.cpp:204
virtual returnValue setupRelaxedQPdata(InfeasibleQPhandling infeasibleQPhandling, DenseCP &_denseCPrelaxed) const
virtual returnValue freezeCondensing()
BooleanType isQP() const
Base class for algorithms solving banded conic programs arising in optimal control.
returnValue computeCondensingOperator(BandedCP &cp)
virtual returnValue solveQP(uint maxIter, InfeasibleQPhandling infeasibleQPhandling=IQH_UNDEFINED)
#define CLOSE_NAMESPACE_ACADO
returnValue generateConstraintBlockLine(uint nn, uint rowOffset, uint &rowOffset1)
virtual returnValue init(const DenseCP *cp_)=0
bool isEmpty() const
Definition: vector.hpp:176
BlockMatrix hessian
Definition: banded_cp.hpp:103
(not yet documented)
BlockMatrix lambdaDynamic
Definition: banded_cp.hpp:125
returnValue generateStateBoundVectors(uint nn, uint rowOffset, uint &rowOffset1)
virtual returnValue solve(DenseCP *cp_)=0
EIGEN_STRONG_INLINE Index rows() const
DMatrix H
Definition: dense_cp.hpp:381
virtual BandedCPsolver * clone() const
BlockMatrix lowerBoundResiduum
Definition: banded_cp.hpp:106
uint getNC() const
uint getNumPoints() const
DVector getMergedDualSolution() const
Definition: dense_cp.cpp:257
Returned value is a warning.
returnValue getSubBlock(uint rowIdx, uint colIdx, DMatrix &value) const
Data class for storing generic conic programs.
Definition: dense_cp.hpp:55
bool isSymmetric() const
Definition: matrix.cpp:139
CondensingBasedCPsolver & operator=(const CondensingBasedCPsolver &rhs)
virtual returnValue getVarianceCovariance(DMatrix &var)
printLevel
Definition: example1.py:55
unsigned getDim() const
Definition: matrix.hpp:181
T getConditionNumber() const
Definition: matrix.cpp:235
unsigned getDim() const
Definition: vector.hpp:172
returnValue init(uint _nRows, uint _nCols)
virtual returnValue getDualSolution(DVector &yOpt) const =0
static Logger & instance()
returnValueLevel getLogLevel()
const MatrixType & eigenvectors() const
Returns the eigenvectors of given matrix.
virtual returnValue prepareSolve(BandedCP &cp)
InfeasibleQPhandling
uint getNV() const
BlockMatrix lowerConstraintResiduum
Definition: banded_cp.hpp:113
BlockMatrix lambdaBound
Definition: banded_cp.hpp:124
BlockMatrix upperBoundResiduum
Definition: banded_cp.hpp:107
Derived & setZero(Index size)
virtual returnValue init(const DenseCP *cp)
DVector * x
Definition: dense_cp.hpp:398
virtual returnValue getParameters(DVector &p_) const
Encapsulates all user interaction for setting options, logging data and plotting results.
virtual returnValue solve(BandedCP &cp)
returnValue generateStateBoundBlockLine(uint nn, uint rowOffset, uint &rowOffset1)
void rhs(const real_t *x, real_t *f)
PrintLevel
virtual returnValue getPrimalSolution(DVector &xOpt) const =0
virtual returnValue start()
Definition: real_clock.cpp:78
virtual uint getNumberOfVariables() const =0
unsigned getNumRows() const
Definition: matrix.hpp:185
DVector ubA
Definition: dense_cp.hpp:389
#define ASSERT(x)
#define BT_TRUE
Definition: acado_types.hpp:47
T getMin() const
Definition: matrix.hpp:297
virtual returnValue setRealTimeParameters(const DVector &DeltaX, const DVector &DeltaP=emptyConstVector)
virtual DenseQPsolver * cloneDenseQPsolver() const =0
unsigned getNumCols() const
Definition: matrix.hpp:189
returnValue generateConstraintVectors(uint nn, uint rowOffset, uint &rowOffset1)
DMatrix A
Definition: dense_cp.hpp:387
BlockMatrix dynResiduum
Definition: banded_cp.hpp:110
virtual returnValue stop()
Definition: real_clock.cpp:107
virtual returnValue setupRelaxedQPdataL2(DenseCP &_denseCPrelaxed) const
BlockMatrix dynGradient
Definition: banded_cp.hpp:109
returnValue setZero(uint rowIdx, uint colIdx)
virtual returnValue solve(DenseCP *cp_)
virtual uint getNumberOfConstraints() const =0
#define ACADOWARNING(retval)
#define BEGIN_NAMESPACE_ACADO
DVector lbA
Definition: dense_cp.hpp:388
#define BT_FALSE
Definition: acado_types.hpp:49
returnValue condense(BandedCP &cp)
virtual returnValue solveCPsubproblem()
returnValue symmetrize()
Definition: matrix.cpp:156
returnValue init(uint nV_, uint nC_)
Definition: dense_cp.cpp:195
BlockMatrix constraintGradient
Definition: banded_cp.hpp:112
Solves banded conic programs arising in optimal control using condensing.
virtual returnValue init(const OCPiterate &iter_)
returnValue setIdentity(uint rowIdx, uint colIdx, uint dim)
virtual returnValue setupRelaxedQPdataL1(DenseCP &_denseCPrelaxed) const
virtual uint getNumberOfIterations() const =0
const RealVectorType & eigenvalues() const
Returns the eigenvalues of given matrix.
virtual DenseCPsolver * clone() const =0
#define ACADOERROR(retval)
virtual returnValue getTime(double &_elapsedTime)
Definition: clock.cpp:92
uint getNumRows() const
Data class for storing conic programs arising from optimal control.
Definition: banded_cp.hpp:56


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Mon Jun 10 2019 12:34:31