qpOASES  3.2.2
An Implementation of the Online Active Set Strategy
SQProblemSchur.hpp
Go to the documentation of this file.
1 /*
2  * This file is part of qpOASES.
3  *
4  * qpOASES -- An Implementation of the Online Active Set Strategy.
5  * Copyright (C) 2007-2014 by Hans Joachim Ferreau, Andreas Potschka,
6  * Christian Kirches et al. All rights reserved.
7  *
8  * qpOASES is free software; you can redistribute it and/or
9  * modify it under the terms of the GNU Lesser General Public
10  * License as published by the Free Software Foundation; either
11  * version 2.1 of the License, or (at your option) any later version.
12  *
13  * qpOASES is distributed in the hope that it will be useful,
14  * but WITHOUT ANY WARRANTY; without even the implied warranty of
15  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
16  * See the GNU Lesser General Public License for more details.
17  *
18  * You should have received a copy of the GNU Lesser General Public
19  * License along with qpOASES; if not, write to the Free Software
20  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
21  *
22  */
23 
24 
38 #ifndef QPOASES_SQPROBLEMSCHUR_HPP
39 #define QPOASES_SQPROBLEMSCHUR_HPP
40 
41 
42 #include <qpOASES/SQProblem.hpp>
43 #include <qpOASES/SparseSolver.hpp>
44 
45 
47 
48 
60 class SQProblemSchur : public SQProblem
61 {
62  /* allow SolutionAnalysis class to access private members */
63  friend class SolutionAnalysis;
64 
65  /*
66  * PUBLIC MEMBER FUNCTIONS
67  */
68  public:
70  SQProblemSchur( );
71 
77  SQProblemSchur( int_t _nV,
78  int_t _nC,
79  HessianType _hessianType = HST_UNKNOWN,
80  int_t maxSchurUpdates = 75
81  );
82 
84  SQProblemSchur( const SQProblemSchur& rhs
85  );
86 
88  virtual ~SQProblemSchur( );
89 
91  virtual SQProblemSchur& operator=( const SQProblemSchur& rhs
92  );
93 
97  virtual returnValue reset( );
98 
101  returnValue resetSchurComplement( BooleanType allowInertiaCorrection );
102 
104  inline int_t getNumFactorizations( ) const;
105 
106  /*
107  * PROTECTED MEMBER FUNCTIONS
108  */
109  protected:
112  returnValue clear( );
113 
116  returnValue copy( const SQProblemSchur& rhs
117  );
118 
122 
127 
131 
143  Matrix *A_new,
145  const real_t *lb_new,
146  const real_t *ub_new,
147  const real_t *lbA_new,
148  const real_t *ubA_new
149  );
150 
157  virtual returnValue setupAuxiliaryWorkingSet( const Bounds* const auxiliaryBounds,
158  const Constraints* const auxiliaryConstraints,
159  BooleanType setupAfresh
160  );
161 
162 
168  virtual returnValue addConstraint( int_t number,
169  SubjectToStatus C_status,
170  BooleanType updateCholesky,
171  BooleanType ensureLI = BT_TRUE
172  );
173 
179  virtual returnValue addConstraint_checkLI( int_t number
180  );
181 
190  virtual returnValue addConstraint_ensureLI( int_t number,
191  SubjectToStatus C_status
192  );
193 
199  virtual returnValue addBound( int_t number,
200  SubjectToStatus B_status,
201  BooleanType updateCholesky,
202  BooleanType ensureLI = BT_TRUE
203  );
204 
209  virtual returnValue addBound_checkLI( int_t number
210  );
211 
220  virtual returnValue addBound_ensureLI( int_t number,
221  SubjectToStatus B_status
222  );
223 
229  virtual returnValue removeConstraint( int_t number,
230  BooleanType updateCholesky,
231  BooleanType allowFlipping = BT_FALSE,
232  BooleanType ensureNZC = BT_FALSE
233  );
234 
240  virtual returnValue removeBound( int_t number,
241  BooleanType updateCholesky,
242  BooleanType allowFlipping = BT_FALSE,
243  BooleanType ensureNZC = BT_FALSE
244  );
245 
248  virtual returnValue backsolveT( const real_t* const b,
249  BooleanType transposed,
250  real_t* const a
251  ) const;
252 
255  virtual returnValue backsolveR( const real_t* const b,
256  BooleanType transposed,
257  real_t* const a
258  ) const;
259 
262  virtual returnValue backsolveR( const real_t* const b,
263  BooleanType transposed,
264  BooleanType removingBound,
265  real_t* const a
266  ) const;
267 
268 
273  virtual returnValue determineStepDirection( const real_t* const delta_g,
274  const real_t* const delta_lbA,
275  const real_t* const delta_ubA,
276  const real_t* const delta_lb,
277  const real_t* const delta_ub,
278  BooleanType Delta_bC_isZero,
279  BooleanType Delta_bB_isZero,
280  real_t* const delta_xFX,
281  real_t* const delta_xFR,
282  real_t* const delta_yAC,
283  real_t* const delta_yFX
284  );
285 
286  virtual returnValue determineStepDirection2( const real_t* const delta_g,
287  const real_t* const delta_lbA,
288  const real_t* const delta_ubA,
289  const real_t* const delta_lb,
290  const real_t* const delta_ub,
291  BooleanType Delta_bC_isZero,
292  BooleanType Delta_bB_isZero,
293  real_t* const delta_xFX,
294  real_t* const delta_xFR,
295  real_t* const delta_yAC,
296  real_t* const delta_yFX
297  );
298 
299  /*
300  * PRIVATE MEMBER FUNCTION
301  */
302  private:
309  real_t* const xiC,
310  real_t* const xiX
311  );
312 
319  real_t* const xiC,
320  real_t* const xiX
321  );
322 
325  returnValue computeMTimes( real_t alpha, const real_t* const x, real_t beta, real_t* const y );
326 
329  returnValue computeMTransTimes( real_t alpha, const real_t* const x, real_t beta, real_t* const y );
330 
332  returnValue addToSchurComplement( int_t number, SchurUpdateType update, int_t numNonzerosM, const sparse_int_t* M_pos, const real_t* const M_vals, int_t numNonzerosN, const sparse_int_t* Npos, const real_t* const Nvals, real_t N_diag );
333 
336 
339 
341  real_t calcDetSchur( int_t idxDel );
342 
344  returnValue updateSchurQR( int_t idxDel );
345 
347  returnValue backsolveSchurQR( int_t dimS, const real_t* const rhs, int_t dimRhs, real_t* const sol );
348 
351 
354 
355  returnValue stepCalcRhs( int_t nFR, int_t nFX, int_t nAC, int_t* FR_idx, int_t* FX_idx, int_t* AC_idx, real_t& rhs_max, const real_t* const delta_g,
356  const real_t* const delta_lbA, const real_t* const delta_ubA,
357  const real_t* const delta_lb, const real_t* const delta_ub,
358  BooleanType Delta_bC_isZero, BooleanType Delta_bB_isZero,
359  real_t* const delta_xFX, real_t* const delta_xFR,
360  real_t* const delta_yAC, real_t* const delta_yFX
361  );
362 
363  returnValue stepCalcReorder(int_t nFR, int_t nAC, int_t* FR_idx, int_t* AC_idx, int_t nFRStart, int_t nACStart,
364  int_t* FR_idxStart, int_t* AC_idxStart, int_t* FR_iSort, int_t* FR_iSortStart,
365  int_t* AC_iSort, int_t* AC_iSortStart, real_t* rhs
366  );
367 
368  returnValue stepCalcBacksolveSchur( int_t nFR, int_t nFX, int_t nAC, int_t* FR_idx, int_t* FX_idx, int_t* AC_idx,
369  int_t dim, real_t* rhs, real_t* sol
370  );
371 
372  returnValue stepCalcReorder2( int_t nFR, int_t nAC, int_t* FR_idx, int_t* AC_idx, int_t nFRStart, int_t nACStart,
373  int_t* FR_idxStart, int_t* AC_idxStart, int_t* FR_iSort, int_t* FR_iSortStart,
374  int_t* AC_iSort, int_t* AC_iSortStart, real_t* sol, real_t* const delta_xFR, real_t* const delta_yAC
375  );
376 
377  returnValue stepCalcResid( int_t nFR, int_t nFX, int_t nAC, int_t* FR_idx, int_t* FX_idx, int_t* AC_idx,
378  BooleanType Delta_bC_isZero, real_t* const delta_xFX, real_t* const delta_xFR,
379  real_t* const delta_yAC, const real_t* const delta_g,
380  const real_t* const delta_lbA, const real_t* const delta_ubA, real_t& rnrm
381  );
382 
383  returnValue stepCalcDeltayFx( int_t nFR, int_t nFX, int_t nAC, int_t* FX_idx, const real_t* const delta_g,
384  real_t* const delta_xFX, real_t* const delta_xFR, real_t* const delta_yAC, real_t* const delta_yFX
385  );
386 
387  /*
388  * PROTECTED MEMBER VARIABLES
389  */
390  protected:
413 };
414 
415 
417 
419 
420 #endif /* QPOASES_QPROBLEMSCHUR_HPP */
421 
422 
423 /*
424  * end of file
425  */
returnValue
Defines all symbols for global return values.
Definition: MessageHandling.hpp:65
BooleanType
Definition: Types.hpp:204
@ BT_TRUE
Definition: Types.hpp:206
@ BT_FALSE
Definition: Types.hpp:205
HessianType
Definition: Types.hpp:249
@ HST_UNKNOWN
Definition: Types.hpp:256
SubjectToStatus
Definition: Types.hpp:273
int_t sparse_int_t
Definition: Types.hpp:199
SchurUpdateType
Definition: Types.hpp:284
int int_t
Definition: Types.hpp:180
BEGIN_NAMESPACE_QPOASES typedef double real_t
Definition: Types.hpp:171
#define END_NAMESPACE_QPOASES
Definition: Types.hpp:110
#define BEGIN_NAMESPACE_QPOASES
Definition: Types.hpp:107
Manages working sets of bounds (i.e. box constraints).
Definition: Bounds.hpp:57
Manages working sets of constraints.
Definition: Constraints.hpp:57
Stores and manages index lists.
Definition: Indexlist.hpp:56
Abstract base class for interfacing tailored matrix-vector operations.
Definition: Matrices.hpp:60
real_t * y
Definition: QProblemB.hpp:984
real_t * x
Definition: QProblemB.hpp:983
virtual returnValue setupAuxiliaryQP(const Bounds *const guessedBounds)
Definition: QProblemB.cpp:2980
Implements the online active set strategy for QPs with varying, sparse matrices.
Definition: SQProblemSchur.hpp:61
returnValue deleteFromSchurComplement(int_t idx, BooleanType allowUndo=BT_FALSE)
Definition: SQProblemSchur.cpp:3255
virtual returnValue addConstraint_checkLI(int_t number)
Definition: SQProblemSchur.cpp:841
real_t calcDetSchur(int_t idxDel)
Definition: SQProblemSchur.cpp:2198
virtual returnValue removeBound(int_t number, BooleanType updateCholesky, BooleanType allowFlipping=BT_FALSE, BooleanType ensureNZC=BT_FALSE)
Definition: SQProblemSchur.cpp:1805
int_t * schurUpdateIndex
Definition: SQProblemSchur.hpp:403
returnValue correctInertia()
Definition: SQProblemSchur.cpp:3439
virtual SQProblemSchur & operator=(const SQProblemSchur &rhs)
Definition: SQProblemSchur.cpp:183
returnValue clear()
Definition: SQProblemSchur.cpp:218
real_t * M_vals
Definition: SQProblemSchur.hpp:407
virtual returnValue addBound_ensureLI(int_t number, SubjectToStatus B_status)
Definition: SQProblemSchur.cpp:1364
returnValue copy(const SQProblemSchur &rhs)
Definition: SQProblemSchur.cpp:242
returnValue addConstraint_checkLISchur(int_t number, real_t *const xiC, real_t *const xiX)
Definition: SQProblemSchur.cpp:862
virtual returnValue setupTQfactorisation()
Definition: SQProblemSchur.cpp:646
SparseSolver * sparseSolver
Definition: SQProblemSchur.hpp:391
int_t numFactorizations
Definition: SQProblemSchur.hpp:401
real_t rcondS
Definition: SQProblemSchur.hpp:400
virtual returnValue backsolveR(const real_t *const b, BooleanType transposed, real_t *const a) const
Definition: SQProblemSchur.cpp:2180
returnValue addToSchurComplement(int_t number, SchurUpdateType update, int_t numNonzerosM, const sparse_int_t *M_pos, const real_t *const M_vals, int_t numNonzerosN, const sparse_int_t *Npos, const real_t *const Nvals, real_t N_diag)
Definition: SQProblemSchur.cpp:3173
returnValue backsolveSchurQR(int_t dimS, const real_t *const rhs, int_t dimRhs, real_t *const sol)
Definition: SQProblemSchur.cpp:2430
sparse_int_t * M_ir
Definition: SQProblemSchur.hpp:408
virtual returnValue backsolveT(const real_t *const b, BooleanType transposed, real_t *const a) const
Definition: SQProblemSchur.cpp:2171
sparse_int_t * M_jc
Definition: SQProblemSchur.hpp:409
int_t getNumFactorizations() const
Definition: SQProblemSchur.ipp:47
int_t nS
Definition: SQProblemSchur.hpp:394
returnValue resetSchurComplement(BooleanType allowInertiaCorrection)
Definition: SQProblemSchur.cpp:2995
Indexlist boundsFreeStart
Definition: SQProblemSchur.hpp:411
virtual returnValue addBound(int_t number, SubjectToStatus B_status, BooleanType updateCholesky, BooleanType ensureLI=BT_TRUE)
Definition: SQProblemSchur.cpp:1122
real_t * Q_
Definition: SQProblemSchur.hpp:397
real_t * S
Definition: SQProblemSchur.hpp:393
virtual returnValue computeInitialCholesky()
Definition: SQProblemSchur.cpp:637
virtual returnValue addConstraint(int_t number, SubjectToStatus C_status, BooleanType updateCholesky, BooleanType ensureLI=BT_TRUE)
Definition: SQProblemSchur.cpp:655
virtual ~SQProblemSchur()
Definition: SQProblemSchur.cpp:172
real_t detS
Definition: SQProblemSchur.hpp:399
virtual returnValue computeProjectedCholesky()
Definition: SQProblemSchur.cpp:628
SchurUpdateType * schurUpdate
Definition: SQProblemSchur.hpp:404
returnValue addBound_checkLISchur(int_t number, real_t *const xiC, real_t *const xiX)
Definition: SQProblemSchur.cpp:1277
returnValue undoDeleteFromSchurComplement(int_t idx)
Definition: SQProblemSchur.cpp:3354
virtual returnValue determineStepDirection2(const real_t *const delta_g, const real_t *const delta_lbA, const real_t *const delta_ubA, const real_t *const delta_lb, const real_t *const delta_ub, BooleanType Delta_bC_isZero, BooleanType Delta_bB_isZero, real_t *const delta_xFX, real_t *const delta_xFR, real_t *const delta_yAC, real_t *const delta_yFX)
Definition: SQProblemSchur.cpp:2852
returnValue computeMTransTimes(real_t alpha, const real_t *const x, real_t beta, real_t *const y)
Definition: SQProblemSchur.cpp:3144
virtual returnValue removeConstraint(int_t number, BooleanType updateCholesky, BooleanType allowFlipping=BT_FALSE, BooleanType ensureNZC=BT_FALSE)
Definition: SQProblemSchur.cpp:1530
virtual returnValue determineStepDirection(const real_t *const delta_g, const real_t *const delta_lbA, const real_t *const delta_ubA, const real_t *const delta_lb, const real_t *const delta_ub, BooleanType Delta_bC_isZero, BooleanType Delta_bB_isZero, real_t *const delta_xFX, real_t *const delta_xFR, real_t *const delta_yAC, real_t *const delta_yFX)
Definition: SQProblemSchur.cpp:2821
returnValue computeMTimes(real_t alpha, const real_t *const x, real_t beta, real_t *const y)
Definition: SQProblemSchur.cpp:3127
virtual returnValue setupAuxiliaryWorkingSet(const Bounds *const auxiliaryBounds, const Constraints *const auxiliaryConstraints, BooleanType setupAfresh)
Definition: SQProblemSchur.cpp:519
int_t nSmax
Definition: SQProblemSchur.hpp:395
returnValue repairSingularWorkingSet()
Definition: SQProblemSchur.cpp:3536
virtual returnValue addConstraint_ensureLI(int_t number, SubjectToStatus C_status)
Definition: SQProblemSchur.cpp:954
returnValue updateSchurQR(int_t idxDel)
Definition: SQProblemSchur.cpp:2291
real_t * R_
Definition: SQProblemSchur.hpp:398
Indexlist constraintsActiveStart
Definition: SQProblemSchur.hpp:412
virtual returnValue addBound_checkLI(int_t number)
Definition: SQProblemSchur.cpp:1257
SQProblemSchur()
Definition: SQProblemSchur.cpp:69
virtual returnValue reset()
Definition: SQProblemSchur.cpp:198
int_t M_physicallength
Definition: SQProblemSchur.hpp:406
Implements the online active set strategy for QPs with varying matrices.
Definition: SQProblem.hpp:60
Provides additional tools for analysing QP solutions.
Definition: SolutionAnalysis.hpp:58
Base class for linear solvers that are used in a Schur-complement implementation in qpOASES.
Definition: SparseSolver.hpp:55
Abstract base class for interfacing matrix-vector operations tailored to symmetric matrices.
Definition: Matrices.hpp:293