qpOASES  3.2.1
An Implementation of the Online Active Set Strategy
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Friends
QProblemB Class Reference

Implements the online active set strategy for box-constrained QPs. More...

#include <QProblemB.hpp>

Inheritance diagram for QProblemB:
QProblem SQProblem SQProblemSchur

List of all members.

Public Member Functions

 QProblemB ()
 QProblemB (int_t _nV, HessianType _hessianType=HST_UNKNOWN, BooleanType allocDenseMats=BT_TRUE)
 QProblemB (const QProblemB &rhs)
virtual ~QProblemB ()
virtual QProblemBoperator= (const QProblemB &rhs)
virtual returnValue reset ()
returnValue init (SymmetricMatrix *_H, const real_t *const _g, const real_t *const _lb, const real_t *const _ub, int_t &nWSR, real_t *const cputime=0, const real_t *const xOpt=0, const real_t *const yOpt=0, const Bounds *const guessedBounds=0, const real_t *const _R=0)
returnValue init (const real_t *const _H, const real_t *const _g, const real_t *const _lb, const real_t *const _ub, int_t &nWSR, real_t *const cputime=0, const real_t *const xOpt=0, const real_t *const yOpt=0, const Bounds *const guessedBounds=0, const real_t *const _R=0)
returnValue init (const char *const H_file, const char *const g_file, const char *const lb_file, const char *const ub_file, int_t &nWSR, real_t *const cputime=0, const real_t *const xOpt=0, const real_t *const yOpt=0, const Bounds *const guessedBounds=0, const char *const R_file=0)
returnValue hotstart (const real_t *const g_new, const real_t *const lb_new, const real_t *const ub_new, int_t &nWSR, real_t *const cputime=0, const Bounds *const guessedBounds=0)
returnValue hotstart (const char *const g_file, const char *const lb_file, const char *const ub_file, int_t &nWSR, real_t *const cputime=0, const Bounds *const guessedBounds=0)
virtual returnValue getWorkingSet (real_t *workingSet)
virtual returnValue getWorkingSetBounds (real_t *workingSetB)
virtual returnValue getWorkingSetConstraints (real_t *workingSetC)
returnValue getBounds (Bounds &_bounds) const
int_t getNV () const
int_t getNFR () const
int_t getNFX () const
int_t getNFV () const
virtual int_t getNZ () const
real_t getObjVal () const
real_t getObjVal (const real_t *const _x) const
returnValue getPrimalSolution (real_t *const xOpt) const
virtual returnValue getDualSolution (real_t *const yOpt) const
QProblemStatus getStatus () const
BooleanType isInitialised () const
BooleanType isSolved () const
BooleanType isInfeasible () const
BooleanType isUnbounded () const
HessianType getHessianType () const
returnValue setHessianType (HessianType _hessianType)
BooleanType usingRegularisation () const
Options getOptions () const
returnValue setOptions (const Options &_options)
PrintLevel getPrintLevel () const
returnValue setPrintLevel (PrintLevel _printlevel)
uint_t getCount () const
returnValue resetCounter ()
virtual returnValue printProperties ()
returnValue printOptions () const

Protected Member Functions

returnValue clear ()
returnValue copy (const QProblemB &rhs)
returnValue determineHessianType ()
virtual returnValue setupSubjectToType ()
virtual returnValue setupSubjectToType (const real_t *const lb_new, const real_t *const ub_new)
virtual returnValue computeCholesky ()
virtual returnValue setupInitialCholesky ()
returnValue obtainAuxiliaryWorkingSet (const real_t *const xOpt, const real_t *const yOpt, const Bounds *const guessedBounds, Bounds *auxiliaryBounds) const
returnValue areBoundsConsistent (const real_t *const lb, const real_t *const ub) const
virtual returnValue backsolveR (const real_t *const b, BooleanType transposed, real_t *const a) const
virtual returnValue backsolveR (const real_t *const b, BooleanType transposed, BooleanType removingBound, real_t *const a) const
returnValue determineDataShift (const real_t *const g_new, const real_t *const lb_new, const real_t *const ub_new, real_t *const delta_g, real_t *const delta_lb, real_t *const delta_ub, BooleanType &Delta_bB_isZero)
returnValue setupQPdata (SymmetricMatrix *_H, const real_t *const _g, const real_t *const _lb, const real_t *const _ub)
returnValue setupQPdata (const real_t *const _H, const real_t *const _g, const real_t *const _lb, const real_t *const _ub)
returnValue setupQPdataFromFile (const char *const H_file, const char *const g_file, const char *const lb_file, const char *const ub_file)
returnValue loadQPvectorsFromFile (const char *const g_file, const char *const lb_file, const char *const ub_file, real_t *const g_new, real_t *const lb_new, real_t *const ub_new) const
returnValue setInfeasibilityFlag (returnValue returnvalue, BooleanType doThrowError=BT_FALSE)
BooleanType isCPUtimeLimitExceeded (const real_t *const cputime, real_t starttime, int_t nWSR) const
returnValue regulariseHessian ()
returnValue setH (SymmetricMatrix *H_new)
returnValue setH (const real_t *const H_new)
returnValue setG (const real_t *const g_new)
returnValue setLB (const real_t *const lb_new)
returnValue setLB (int_t number, real_t value)
returnValue setUB (const real_t *const ub_new)
returnValue setUB (int_t number, real_t value)
void computeGivens (real_t xold, real_t yold, real_t &xnew, real_t &ynew, real_t &c, real_t &s) const
void applyGivens (real_t c, real_t s, real_t nu, real_t xold, real_t yold, real_t &xnew, real_t &ynew) const
real_t getRelativeHomotopyLength (const real_t *const g_new, const real_t *const lb_new, const real_t *const ub_new)
virtual returnValue performRamping ()
returnValue updateFarBounds (real_t curFarBound, int_t nRamp, const real_t *const lb_new, real_t *const lb_new_far, const real_t *const ub_new, real_t *const ub_new_far) const
returnValue performRatioTest (int_t nIdx, const int_t *const idxList, const SubjectTo *const subjectTo, const real_t *const num, const real_t *const den, real_t epsNum, real_t epsDen, real_t &t, int_t &BC_idx) const
BooleanType isBlocking (real_t num, real_t den, real_t epsNum, real_t epsDen, real_t &t) const
SymSparseMatcreateDiagSparseMat (int_t n, real_t diagVal=1.0)
virtual returnValue setupAuxiliaryQP (const Bounds *const guessedBounds)

Protected Attributes

BooleanType freeHessian
SymmetricMatrixH
real_tg
real_tlb
real_tub
Bounds bounds
real_tR
BooleanType haveCholesky
real_tx
real_ty
real_t tau
QProblemStatus status
BooleanType infeasible
BooleanType unbounded
HessianType hessianType
real_t regVal
uint_t count
real_tdelta_xFR_TMP
real_t ramp0
real_t ramp1
int_t rampOffset
Options options
Flipper flipper
TabularOutput tabularOutput

Private Member Functions

returnValue solveInitialQP (const real_t *const xOpt, const real_t *const yOpt, const Bounds *const guessedBounds, const real_t *const _R, int_t &nWSR, real_t *const cputime)
returnValue solveQP (const real_t *const g_new, const real_t *const lb_new, const real_t *const ub_new, int_t &nWSR, real_t *const cputime, int_t nWSRperformed=0, BooleanType isFirstCall=BT_TRUE)
returnValue solveRegularisedQP (const real_t *const g_new, const real_t *const lb_new, const real_t *const ub_new, int_t &nWSR, real_t *const cputime, int_t nWSRperformed=0, BooleanType isFirstCall=BT_TRUE)
returnValue setupAuxiliaryWorkingSet (const Bounds *const auxiliaryBounds, BooleanType setupAfresh)
returnValue setupAuxiliaryQPsolution (const real_t *const xOpt, const real_t *const yOpt)
returnValue setupAuxiliaryQPgradient ()
returnValue setupAuxiliaryQPbounds (BooleanType useRelaxation)
returnValue determineStepDirection (const real_t *const delta_g, const real_t *const delta_lb, const real_t *const delta_ub, BooleanType Delta_bB_isZero, real_t *const delta_xFX, real_t *const delta_xFR, real_t *const delta_yFX)
returnValue performStep (const real_t *const delta_g, const real_t *const delta_lb, const real_t *const delta_ub, const real_t *const delta_xFX, const real_t *const delta_xFR, const real_t *const delta_yFX, int_t &BC_idx, SubjectToStatus &BC_status)
returnValue changeActiveSet (int_t BC_idx, SubjectToStatus BC_status)
virtual returnValue performDriftCorrection ()
BooleanType shallRefactorise (const Bounds *const guessedBounds) const
returnValue addBound (int_t number, SubjectToStatus B_status, BooleanType updateCholesky)
returnValue removeBound (int_t number, BooleanType updateCholesky)
returnValue printIteration (int_t iter, int_t BC_idx, SubjectToStatus BC_status, real_t homotopyLength, BooleanType isFirstCall=BT_TRUE)

Friends

class SolutionAnalysis

Detailed Description

Class for setting up and solving quadratic programs with bounds (= box constraints) only. The main feature is the possibily to use the newly developed online active set strategy for parametric quadratic programming.

Author:
Hans Joachim Ferreau, Andreas Potschka, Christian Kirches
Version:
3.2
Date:
2007-2017

Constructor & Destructor Documentation

QProblemB::QProblemB ( int_t  _nV,
HessianType  _hessianType = HST_UNKNOWN,
BooleanType  allocDenseMats = BT_TRUE 
)

Constructor which takes the QP dimension and Hessian type information. If the Hessian is the zero (i.e. HST_ZERO) or the identity matrix (i.e. HST_IDENTITY), respectively, no memory is allocated for it and a NULL pointer can be passed for it to the init() functions.

Parameters:
_nVNumber of variables.
_hessianTypeType of Hessian matrix.
allocDenseMatsEnable allocation of dense matrices.

References bounds, BT_FALSE, BT_TRUE, count, delta_xFR_TMP, Options::finalRamping, flipper, freeHessian, g, getGlobalMessageHandler(), H, haveCholesky, hessianType, infeasible, Bounds::init(), Flipper::init(), Options::initialRamping, lb, options, PL_NONE, printCopyrightNotice(), Options::printLevel, QPS_NOTINITIALISED, R, ramp0, ramp1, rampOffset, real_t, regVal, MessageHandling::reset(), RET_INVALID_ARGUMENTS, setPrintLevel(), status, tau, THROWERROR, ub, unbounded, x, and y.

QProblemB::QProblemB ( const QProblemB rhs)

Copy constructor (deep copy).

Parameters:
rhsRhs object.

References BT_FALSE, copy(), freeHessian, and H.

QProblemB::~QProblemB ( ) [virtual]

Member Function Documentation

returnValue QProblemB::addBound ( int_t  number,
SubjectToStatus  B_status,
BooleanType  updateCholesky 
) [private]

Adds a bound to active set (specialised version for the case where no constraints exist).

Returns:
SUCCESSFUL_RETURN
RET_ADDBOUND_FAILED
Parameters:
numberNumber of bound to be added to active set.
B_statusStatus of new active bound.
updateCholeskyFlag indicating if Cholesky decomposition shall be updated.

References applyGivens(), bounds, BT_TRUE, computeGivens(), Bounds::getFree(), Indexlist::getIndex(), getNFR(), getNV(), getStatus(), hessianType, HST_IDENTITY, HST_ZERO, TabularOutput::idxAddB, Bounds::moveFreeToFixed(), QPS_AUXILIARYQPSOLVED, QPS_HOMOTOPYQPSOLVED, QPS_NOTINITIALISED, QPS_PREPARINGAUXILIARYQP, QPS_SOLVED, real_t, RET_ADDBOUND_FAILED, RET_UNKNOWN_BUG, RR, SUCCESSFUL_RETURN, tabularOutput, and THROWERROR.

Referenced by changeActiveSet(), and setupAuxiliaryWorkingSet().

void QProblemB::applyGivens ( real_t  c,
real_t  s,
real_t  nu,
real_t  xold,
real_t  yold,
real_t xnew,
real_t ynew 
) const [inline, protected]

Applies Givens matrix determined by c and s (cf. computeGivens).

Returns:
SUCCESSFUL_RETURN
Parameters:
cCosine entry of Givens matrix.
sSine entry of Givens matrix.
nuFurther factor: s/(1+c).
xoldMatrix entry to be transformed corresponding to the normalised entry of the original matrix.
yoldMatrix entry to be transformed corresponding to the annihilated entry of the original matrix.
xnewOutput: Transformed matrix entry corresponding to the normalised entry of the original matrix.
ynewOutput: Transformed matrix entry corresponding to the annihilated entry of the original matrix.

Referenced by QProblem::addBound(), addBound(), QProblem::addConstraint(), SQProblemSchur::calcDetSchur(), QProblem::removeBound(), QProblem::removeConstraint(), and SQProblemSchur::updateSchurQR().

returnValue QProblemB::areBoundsConsistent ( const real_t *const  lb,
const real_t *const  ub 
) const [protected]

Decides if lower bounds are smaller than upper bounds

Returns:
SUCCESSFUL_RETURN
RET_QP_INFEASIBLE
Parameters:
lbVector of lower bounds
ubVector of upper bounds

References EPS, getNV(), RET_QP_INFEASIBLE, and SUCCESSFUL_RETURN.

Referenced by QProblem::areBoundsConsistent(), and hotstart().

returnValue QProblemB::backsolveR ( const real_t *const  b,
BooleanType  transposed,
real_t *const  a 
) const [protected, virtual]

Solves the system Ra = b or R^Ta = b where R is an upper triangular matrix.

Returns:
SUCCESSFUL_RETURN
RET_DIV_BY_ZERO
Parameters:
bRight hand side vector.
transposedIndicates if the transposed system shall be solved.
aOutput: Solution vector

Reimplemented in SQProblemSchur.

References BT_FALSE.

Referenced by QProblem::determineStepDirection(), determineStepDirection(), QProblem::removeBound(), removeBound(), and QProblem::removeConstraint().

returnValue QProblemB::backsolveR ( const real_t *const  b,
BooleanType  transposed,
BooleanType  removingBound,
real_t *const  a 
) const [protected, virtual]

Solves the system Ra = b or R^Ta = b where R is an upper triangular matrix.
Special variant for the case that this function is called from within "removeBound()".

Returns:
SUCCESSFUL_RETURN
RET_DIV_BY_ZERO
Parameters:
bRight hand side vector.
transposedIndicates if the transposed system shall be solved.
removingBoundIndicates if function is called from "removeBound()".
aOutput: Solution vector

Reimplemented in SQProblemSchur.

References BT_FALSE, BT_TRUE, getAbs(), getNV(), getNZ(), real_t, RET_DIV_BY_ZERO, RR, SUCCESSFUL_RETURN, THROWERROR, and ZERO.

returnValue QProblemB::changeActiveSet ( int_t  BC_idx,
SubjectToStatus  BC_status 
) [private]

Updates active set.

Returns:
SUCCESSFUL_RETURN
RET_REMOVE_FROM_ACTIVESET_FAILED
RET_ADD_TO_ACTIVESET_FAILED
Parameters:
BC_idxIndex of blocking constraint.
BC_statusStatus of blocking constraint.

References __FILE__, __FUNC__, __LINE__, addBound(), BT_TRUE, getGlobalMessageHandler(), MAX_STRING_LENGTH, removeBound(), RET_ADD_TO_ACTIVESET, RET_ADD_TO_ACTIVESET_FAILED, RET_OPTIMAL_SOLUTION_FOUND, RET_REMOVE_FROM_ACTIVESET, RET_REMOVE_FROM_ACTIVESET_FAILED, ST_INACTIVE, ST_LOWER, ST_UNDEFINED, SUCCESSFUL_RETURN, THROWERROR, MessageHandling::throwInfo(), VS_VISIBLE, and y.

Referenced by solveQP().

returnValue QProblemB::clear ( ) [protected]

Frees all allocated memory.

Returns:
SUCCESSFUL_RETURN

Reimplemented in QProblem, and SQProblemSchur.

References BT_TRUE, delta_xFR_TMP, freeHessian, g, H, lb, R, SUCCESSFUL_RETURN, ub, x, and y.

Referenced by operator=(), and ~QProblemB().

returnValue QProblemB::computeCholesky ( ) [protected, virtual]

Computes the Cholesky decomposition of the (simply projected) Hessian (i.e. R^T*R = Z^T*H*Z). It only works in the case where Z is a simple projection matrix! Note: If Hessian turns out not to be positive definite, the Hessian type is set to HST_SEMIDEF accordingly.

Returns:
SUCCESSFUL_RETURN
RET_HESSIAN_NOT_SPD
RET_INDEXLIST_CORRUPTED

References bounds, BT_TRUE, Options::epsRegularisation, getAbs(), Matrix::getCol(), Bounds::getFree(), getMin(), getNFR(), Indexlist::getNumberArray(), getNV(), getSqrt(), H, hessianType, HST_IDENTITY, HST_SEMIDEF, HST_ZERO, options, POTRF, R, regVal, RET_CHOLESKY_OF_ZERO_HESSIAN, RET_HESSIAN_NOT_SPD, RR, SUCCESSFUL_RETURN, THROWERROR, and usingRegularisation().

Referenced by QProblem::computeProjectedCholesky(), setupAuxiliaryQP(), setupInitialCholesky(), and solveQP().

void QProblemB::computeGivens ( real_t  xold,
real_t  yold,
real_t xnew,
real_t ynew,
real_t c,
real_t s 
) const [inline, protected]

Computes parameters for the Givens matrix G for which [x,y]*G = [z,0]

Returns:
SUCCESSFUL_RETURN
Parameters:
xoldMatrix entry to be normalised.
yoldMatrix entry to be annihilated.
xnewOutput: Normalised matrix entry.
ynewOutput: Annihilated matrix entry.
cOutput: Cosine entry of Givens matrix.
sOutput: Sine entry of Givens matrix.

References BT_TRUE, getAbs(), getSqrt(), isZero(), and real_t.

Referenced by QProblem::addBound(), addBound(), QProblem::addConstraint(), SQProblemSchur::calcDetSchur(), QProblem::removeBound(), QProblem::removeConstraint(), and SQProblemSchur::updateSchurQR().

returnValue QProblemB::copy ( const QProblemB rhs) [protected]

Copies all members from given rhs object.

Returns:
SUCCESSFUL_RETURN
Parameters:
rhsRhs object.

References bounds, BT_TRUE, count, delta_xFR_TMP, SymmetricMatrix::duplicateSym(), flipper, freeHessian, g, getNV(), H, haveCholesky, hessianType, infeasible, lb, options, Options::printLevel, R, ramp0, ramp1, rampOffset, real_t, regVal, setG(), setLB(), setPrintLevel(), setUB(), status, SUCCESSFUL_RETURN, tau, ub, unbounded, x, and y.

Referenced by operator=(), and QProblemB().

SymSparseMat * QProblemB::createDiagSparseMat ( int_t  n,
real_t  diagVal = 1.0 
) [protected]

Creates a sparse diagonal (square-)matrix which is a given multiple of the identity matrix.

Returns:
Diagonal matrix
Parameters:
nRow/column dimension of matrix to be created.
diagValValue of all diagonal entries.

References SparseMatrix::createDiagInfo(), Matrix::doFreeMemory(), and real_t.

Referenced by QProblem::computeProjectedCholesky().

returnValue QProblemB::determineDataShift ( const real_t *const  g_new,
const real_t *const  lb_new,
const real_t *const  ub_new,
real_t *const  delta_g,
real_t *const  delta_lb,
real_t *const  delta_ub,
BooleanType Delta_bB_isZero 
) [protected]

Determines step direction of the shift of the QP data.

Returns:
SUCCESSFUL_RETURN
Parameters:
g_newNew gradient vector.
lb_newNew lower bounds.
ub_newNew upper bounds.
delta_gOutput: Step direction of gradient vector.
delta_lbOutput: Step direction of lower bounds.
delta_ubOutput: Step direction of upper bounds.
Delta_bB_isZeroOutput: Indicates if active bounds are to be shifted.

References bounds, BT_FALSE, BT_TRUE, EPS, g, getAbs(), Bounds::getFixed(), getNFX(), Indexlist::getNumberArray(), getNV(), INFTY, lb, SUCCESSFUL_RETURN, and ub.

Referenced by solveQP().

returnValue QProblemB::determineStepDirection ( const real_t *const  delta_g,
const real_t *const  delta_lb,
const real_t *const  delta_ub,
BooleanType  Delta_bB_isZero,
real_t *const  delta_xFX,
real_t *const  delta_xFR,
real_t *const  delta_yFX 
) [private]

Determines step direction of the homotopy path.

Returns:
SUCCESSFUL_RETURN
RET_STEPDIRECTION_FAILED_CHOLESKY
Parameters:
delta_gStep direction of gradient vector.
delta_lbStep direction of lower bounds.
delta_ubStep direction of upper bounds.
Delta_bB_isZeroIndicates if active bounds are to be shifted.
delta_xFXOutput: Primal homotopy step direction of fixed variables.
delta_xFROutput: Primal homotopy step direction of free variables.
delta_yFXOutput: Dual homotopy step direction of fixed variables' multiplier.

References backsolveR(), bounds, BT_FALSE, BT_TRUE, delta_xFR_TMP, Options::epsIterRef, getAbs(), Bounds::getFixed(), Bounds::getFree(), getNFR(), getNFX(), Indexlist::getNumberArray(), SubjectTo::getStatus(), H, hessianType, HST_IDENTITY, HST_ZERO, Options::numRefinementSteps, options, real_t, regVal, RET_STEPDIRECTION_FAILED_CHOLESKY, ST_LOWER, SUCCESSFUL_RETURN, THROWERROR, Matrix::times(), and usingRegularisation().

Referenced by solveQP().

Returns current bounds object of the QP (deep copy).

Returns:
SUCCESSFUL_RETURN
RET_QPOBJECT_NOT_SETUP
Parameters:
_boundsOutput: Bounds object.

References bounds, getNV(), RET_QPOBJECT_NOT_SETUP, SUCCESSFUL_RETURN, and THROWERROR.

Referenced by SolutionAnalysis::checkCurvatureOnStronglyActiveConstraints().

uint_t QProblemB::getCount ( ) const [inline]

Returns the current number of QP problems solved.

Returns:
Number of QP problems solved.

References count.

returnValue QProblemB::getDualSolution ( real_t *const  yOpt) const [virtual]

Returns the dual solution vector.

Returns:
SUCCESSFUL_RETURN
RET_QP_NOT_SOLVED
Parameters:
yOptOutput: Dual solution vector (if QP has been solved).

Reimplemented in QProblem.

References getNV(), getStatus(), QPS_AUXILIARYQPSOLVED, QPS_HOMOTOPYQPSOLVED, QPS_SOLVED, RET_QP_NOT_SOLVED, SUCCESSFUL_RETURN, and y.

Referenced by solveOqpBenchmark().

Returns Hessian type flag (type is not determined due to this call!).

Returns:
Hessian type.

References hessianType.

Referenced by SolutionAnalysis::getKktViolation().

int_t QProblemB::getNFR ( ) const [inline]
int_t QProblemB::getNFV ( ) const [inline]

Returns the number of implicitly fixed variables.

Returns:
Number of implicitly fixed variables.

References bounds, and Bounds::getNFV().

Referenced by setupInitialCholesky(), and QProblem::setupInitialCholesky().

int_t QProblemB::getNFX ( ) const [inline]
int_t QProblemB::getNV ( ) const [inline]

Returns the number of variables.

Returns:
Number of variables.

References bounds, and Bounds::getNV().

Referenced by QProblem::addBound(), addBound(), QProblem::addBound_checkLI(), SQProblemSchur::addBound_checkLISchur(), SQProblemSchur::addBound_ensureLI(), QProblem::addBound_ensureLI(), QProblem::addConstraint(), QProblem::addConstraint_checkLI(), SQProblemSchur::addConstraint_checkLISchur(), SQProblemSchur::addConstraint_ensureLI(), QProblem::addConstraint_ensureLI(), areBoundsConsistent(), backsolveR(), QProblem::changeActiveSet(), computeCholesky(), QProblem::computeProjectedCholesky(), QProblem::copy(), copy(), determineDataShift(), determineHessianType(), QProblem::determineStepDirection(), QProblem::ensureNonzeroCurvature(), getBounds(), QProblem::getConstraints(), getDualSolution(), QProblem::getDualSolution(), QProblem::getFreeVariablesFlags(), SolutionAnalysis::getKktViolation(), getObjVal(), getPrimalSolution(), getRelativeHomotopyLength(), SolutionAnalysis::getVarianceCovariance(), QProblem::getWorkingSet(), getWorkingSetBounds(), SQProblem::hotstart(), hotstart(), QProblem::hotstart(), init(), QProblem::init(), loadQPvectorsFromFile(), obtainAuxiliaryWorkingSet(), QProblem::obtainAuxiliaryWorkingSet(), QProblem::performDriftCorrection(), performDriftCorrection(), performRamping(), QProblem::performRamping(), QProblem::performStep(), performStep(), QProblem::printIteration(), printIteration(), QProblem::printProperties(), printProperties(), regulariseHessian(), QProblem::removeBound(), removeBound(), QProblem::removeConstraint(), QProblem::reset(), reset(), QProblem::setA(), setG(), setH(), setLB(), QProblem::setLBA(), setUB(), QProblem::setUBA(), SQProblemSchur::setupAuxiliaryQP(), QProblem::setupAuxiliaryQP(), setupAuxiliaryQP(), QProblem::setupAuxiliaryQPbounds(), setupAuxiliaryQPbounds(), QProblem::setupAuxiliaryQPgradient(), setupAuxiliaryQPgradient(), QProblem::setupAuxiliaryQPsolution(), setupAuxiliaryQPsolution(), SQProblemSchur::setupAuxiliaryWorkingSet(), QProblem::setupAuxiliaryWorkingSet(), setupAuxiliaryWorkingSet(), setupInitialCholesky(), QProblem::setupInitialCholesky(), SQProblem::setupNewAuxiliaryQP(), QProblem::setupQPdata(), setupQPdataFromFile(), QProblem::setupQPdataFromFile(), setupSubjectToType(), QProblem::setupTQfactorisation(), QProblem::shallRefactorise(), shallRefactorise(), QProblem::solveCurrentEQP(), QProblem::solveInitialQP(), solveInitialQP(), QProblem::solveQP(), solveQP(), QProblem::solveRegularisedQP(), solveRegularisedQP(), QProblem::updateActivitiesForHotstart(), updateFarBounds(), QProblem::updateFarBounds(), QProblem::writeQpDataIntoMatFile(), and QProblem::writeQpWorkspaceIntoMatFile().

int_t QProblemB::getNZ ( ) const [virtual]

Returns the dimension of null space.

Returns:
Dimension of null space.

Reimplemented in QProblem.

References getNFR().

Referenced by backsolveR().

Returns the optimal objective function value.

Returns:
finite value: Optimal objective function value (QP was solved)
+infinity: QP was not yet solved

References getStatus(), INFTY, QPS_AUXILIARYQPSOLVED, QPS_HOMOTOPYQPSOLVED, QPS_SOLVED, real_t, and x.

Referenced by main().

real_t QProblemB::getObjVal ( const real_t *const  _x) const

Returns the objective function value at an arbitrary point x.

Returns:
Objective function value at point x
Parameters:
_xPoint at which the objective function shall be evaluated.

References BT_TRUE, g, getNV(), H, hessianType, HST_IDENTITY, HST_ZERO, real_t, regVal, Matrix::times(), and usingRegularisation().

Options QProblemB::getOptions ( ) const [inline]

Returns current options struct.

Returns:
Current options struct.

References options.

Returns the primal solution vector.

Returns:
SUCCESSFUL_RETURN
RET_QP_NOT_SOLVED
Parameters:
xOptOutput: Primal solution vector (if QP has been solved).

References getNV(), getStatus(), QPS_AUXILIARYQPSOLVED, QPS_HOMOTOPYQPSOLVED, QPS_SOLVED, RET_QP_NOT_SOLVED, SUCCESSFUL_RETURN, and x.

Referenced by main(), and solveOqpBenchmark().

Returns the print level.

Returns:
Print level.

References options, and Options::printLevel.

real_t QProblemB::getRelativeHomotopyLength ( const real_t *const  g_new,
const real_t *const  lb_new,
const real_t *const  ub_new 
) [protected]

Compute relative length of homotopy in data space for termination criterion.

Returns:
Relative length in data space.
Parameters:
g_newFinal gradient.
lb_newFinal lower variable bounds.
ub_newFinal upper variable bounds.

References g, getAbs(), getNV(), lb, real_t, and ub.

Referenced by solveQP().

returnValue QProblemB::getWorkingSet ( real_t workingSet) [virtual]

Writes a vector with the state of the working set

Returns:
SUCCESSFUL_RETURN
RET_INVALID_ARGUMENTS
Parameters:
workingSetOutput: array containing state of the working set.

Reimplemented in QProblem.

References getWorkingSetBounds().

returnValue QProblemB::getWorkingSetBounds ( real_t workingSetB) [virtual]

Writes a vector with the state of the working set of bounds

Returns:
SUCCESSFUL_RETURN
RET_INVALID_ARGUMENTS
Parameters:
workingSetBOutput: array containing state of the working set of bounds.

Reimplemented in QProblem.

References bounds, getNV(), SubjectTo::getStatus(), RET_INVALID_ARGUMENTS, ST_LOWER, ST_UPPER, SUCCESSFUL_RETURN, and THROWERROR.

Referenced by SolutionAnalysis::getKktViolation(), and getWorkingSet().

Writes a vector with the state of the working set of constraints

Returns:
SUCCESSFUL_RETURN
RET_INVALID_ARGUMENTS
Parameters:
workingSetCOutput: array containing state of the working set of constraints.

Reimplemented in QProblem.

References RET_INVALID_ARGUMENTS, SUCCESSFUL_RETURN, and THROWERROR.

returnValue QProblemB::hotstart ( const real_t *const  g_new,
const real_t *const  lb_new,
const real_t *const  ub_new,
int_t nWSR,
real_t *const  cputime = 0,
const Bounds *const  guessedBounds = 0 
)

Solves an initialised QP sequence using the online active set strategy. By default, QP solution is started from previous solution. If a guess for the working set is provided, an initialised homotopy is performed.

Note: This function internally calls solveQP/solveRegularisedQP for solving an initialised QP!

Returns:
SUCCESSFUL_RETURN
RET_MAX_NWSR_REACHED
RET_HOTSTART_FAILED_AS_QP_NOT_INITIALISED
RET_HOTSTART_FAILED
RET_SHIFT_DETERMINATION_FAILED
RET_STEPDIRECTION_DETERMINATION_FAILED
RET_STEPLENGTH_DETERMINATION_FAILED
RET_HOMOTOPY_STEP_FAILED
RET_HOTSTART_STOPPED_INFEASIBILITY
RET_HOTSTART_STOPPED_UNBOUNDEDNESS
RET_SETUP_AUXILIARYQP_FAILED
Parameters:
g_newGradient of neighbouring QP to be solved.
lb_newLower bounds of neighbouring QP to be solved.
If no lower bounds exist, a NULL pointer can be passed.
ub_newUpper bounds of neighbouring QP to be solved.
If no upper bounds exist, a NULL pointer can be passed.
nWSRInput: Maximum number of working set recalculations;
Output: Number of performed working set recalculations.
cputimeInput: Maximum CPU time allowed for QP solution.
Output: CPU time spent for QP solution (or to perform nWSR iterations).
guessedBoundsOptimal working set of bounds for solution (xOpt,yOpt).
(If a null pointer is passed, the previous working set is kept!)

References areBoundsConsistent(), Options::boundTolerance, BT_FALSE, BT_TRUE, count, Options::enableFarBounds, getAbs(), getCPUtime(), getNV(), Options::growFarBounds, haveCholesky, infeasible, INFTY, Options::initialFarBounds, options, QPS_AUXILIARYQPSOLVED, QPS_HOMOTOPYQPSOLVED, QPS_SOLVED, rampOffset, real_t, RET_HOTSTART_STOPPED_INFEASIBILITY, RET_HOTSTART_STOPPED_UNBOUNDEDNESS, RET_QPOBJECT_NOT_SETUP, RET_SETUP_AUXILIARYQP_FAILED, setInfeasibilityFlag(), setupAuxiliaryQP(), setupInitialCholesky(), solveRegularisedQP(), status, SUCCESSFUL_RETURN, THROWERROR, unbounded, updateFarBounds(), and x.

Referenced by hotstart(), main(), solveInitialQP(), and solveOqpBenchmark().

returnValue QProblemB::hotstart ( const char *const  g_file,
const char *const  lb_file,
const char *const  ub_file,
int_t nWSR,
real_t *const  cputime = 0,
const Bounds *const  guessedBounds = 0 
)

Solves an initialised QP sequence using the online active set strategy, where QP data is read from files. By default, QP solution is started from previous solution. If a guess for the working set is provided, an initialised homotopy is performed.

Note: This function internally calls solveQP/solveRegularisedQP for solving an initialised QP!

Returns:
SUCCESSFUL_RETURN
RET_MAX_NWSR_REACHED
RET_HOTSTART_FAILED_AS_QP_NOT_INITIALISED
RET_HOTSTART_FAILED
RET_SHIFT_DETERMINATION_FAILED
RET_STEPDIRECTION_DETERMINATION_FAILED
RET_STEPLENGTH_DETERMINATION_FAILED
RET_HOMOTOPY_STEP_FAILED
RET_HOTSTART_STOPPED_INFEASIBILITY
RET_HOTSTART_STOPPED_UNBOUNDEDNESS
RET_UNABLE_TO_READ_FILE
RET_SETUP_AUXILIARYQP_FAILED
RET_INVALID_ARGUMENTS
Parameters:
g_fileName of file where gradient, of neighbouring QP to be solved, is stored.
lb_fileName of file where lower bounds, of neighbouring QP to be solved, is stored.
If no lower bounds exist, a NULL pointer can be passed.
ub_fileName of file where upper bounds, of neighbouring QP to be solved, is stored.
If no upper bounds exist, a NULL pointer can be passed.
nWSRInput: Maximum number of working set recalculations;
Output: Number of performed working set recalculations.
cputimeInput: Maximum CPU time allowed for QP solution.
Output: CPU time spent for QP solution (or to perform nWSR iterations).
guessedBoundsOptimal working set of bounds for solution (xOpt,yOpt).
(If a null pointer is passed, the previous working set is kept!)

References getNV(), hotstart(), loadQPvectorsFromFile(), real_t, RET_INVALID_ARGUMENTS, RET_QPOBJECT_NOT_SETUP, RET_UNABLE_TO_READ_FILE, SUCCESSFUL_RETURN, and THROWERROR.

returnValue QProblemB::init ( SymmetricMatrix _H,
const real_t *const  _g,
const real_t *const  _lb,
const real_t *const  _ub,
int_t nWSR,
real_t *const  cputime = 0,
const real_t *const  xOpt = 0,
const real_t *const  yOpt = 0,
const Bounds *const  guessedBounds = 0,
const real_t *const  _R = 0 
)

Initialises a simply bounded QP problem with given QP data and tries to solve it using at most nWSR iterations. Depending on the parameter constellation it:
1. 0, 0, 0 : starts with xOpt = 0, yOpt = 0 and gB empty (or all implicit equality bounds),
2. xOpt, 0, 0 : starts with xOpt, yOpt = 0 and obtain gB by "clipping",
3. 0, yOpt, 0 : starts with xOpt = 0, yOpt and obtain gB from yOpt != 0,
4. 0, 0, gB: starts with xOpt = 0, yOpt = 0 and gB,
5. xOpt, yOpt, 0 : starts with xOpt, yOpt and obtain gB from yOpt != 0,
6. xOpt, 0, gB: starts with xOpt, yOpt = 0 and gB,
7. xOpt, yOpt, gB: starts with xOpt, yOpt and gB (assume them to be consistent!)

Note: This function internally calls solveInitialQP for initialisation!

Returns:
SUCCESSFUL_RETURN
RET_INIT_FAILED
RET_INIT_FAILED_CHOLESKY
RET_INIT_FAILED_HOTSTART
RET_INIT_FAILED_INFEASIBILITY
RET_INIT_FAILED_UNBOUNDEDNESS
RET_MAX_NWSR_REACHED
RET_INVALID_ARGUMENTS
Parameters:
_HHessian matrix (a shallow copy is made).
_gGradient vector.
_lbLower bounds (on variables).
If no lower bounds exist, a NULL pointer can be passed.
_ubUpper bounds (on variables).
If no upper bounds exist, a NULL pointer can be passed.
nWSRInput: Maximum number of working set recalculations when using initial homotopy.
Output: Number of performed working set recalculations.
cputimeInput: Maximum CPU time allowed for QP initialisation.
Output: CPU time spent for QP initialisation (if pointer passed).
xOptOptimal primal solution vector. A NULL pointer can be passed.
(If a null pointer is passed, the old primal solution is kept!)
yOptOptimal dual solution vector. A NULL pointer can be passed.
(If a null pointer is passed, the old dual solution is kept!)
guessedBoundsOptimal working set of bounds for solution (xOpt,yOpt).
(If a null pointer is passed, all bounds are assumed inactive!)
_RPre-computed (upper triangular) Cholesky factor of Hessian matrix. The Cholesky factor must be stored in a real_t array of size nV*nV in row-major format. Note: Only used if xOpt/yOpt and gB are NULL!
(If a null pointer is passed, Cholesky decomposition is computed internally!)

References BT_TRUE, getNV(), SubjectTo::getStatus(), isInitialised(), reset(), RET_INVALID_ARGUMENTS, RET_NO_CHOLESKY_WITH_INITIAL_GUESS, RET_QP_ALREADY_INITIALISED, RET_QPOBJECT_NOT_SETUP, setupQPdata(), solveInitialQP(), ST_UNDEFINED, SUCCESSFUL_RETURN, THROWERROR, and THROWWARNING.

Referenced by main(), and solveOqpBenchmark().

returnValue QProblemB::init ( const real_t *const  _H,
const real_t *const  _g,
const real_t *const  _lb,
const real_t *const  _ub,
int_t nWSR,
real_t *const  cputime = 0,
const real_t *const  xOpt = 0,
const real_t *const  yOpt = 0,
const Bounds *const  guessedBounds = 0,
const real_t *const  _R = 0 
)

Initialises a simply bounded QP problem with given QP data and tries to solve it using at most nWSR iterations. Depending on the parameter constellation it:
1. 0, 0, 0 : starts with xOpt = 0, yOpt = 0 and gB empty (or all implicit equality bounds),
2. xOpt, 0, 0 : starts with xOpt, yOpt = 0 and obtain gB by "clipping",
3. 0, yOpt, 0 : starts with xOpt = 0, yOpt and obtain gB from yOpt != 0,
4. 0, 0, gB: starts with xOpt = 0, yOpt = 0 and gB,
5. xOpt, yOpt, 0 : starts with xOpt, yOpt and obtain gB from yOpt != 0,
6. xOpt, 0, gB: starts with xOpt, yOpt = 0 and gB,
7. xOpt, yOpt, gB: starts with xOpt, yOpt and gB (assume them to be consistent!)

Note: This function internally calls solveInitialQP for initialisation!

Returns:
SUCCESSFUL_RETURN
RET_INIT_FAILED
RET_INIT_FAILED_CHOLESKY
RET_INIT_FAILED_HOTSTART
RET_INIT_FAILED_INFEASIBILITY
RET_INIT_FAILED_UNBOUNDEDNESS
RET_MAX_NWSR_REACHED
RET_INVALID_ARGUMENTS
Parameters:
_HHessian matrix (a shallow copy is made).
If Hessian matrix is trivial, a NULL pointer can be passed.
_gGradient vector.
_lbLower bounds (on variables).
If no lower bounds exist, a NULL pointer can be passed.
_ubUpper bounds (on variables).
If no upper bounds exist, a NULL pointer can be passed.
nWSRInput: Maximum number of working set recalculations when using initial homotopy.
Output: Number of performed working set recalculations.
cputimeInput: Maximum CPU time allowed for QP initialisation.
Output: CPU time spent for QP initialisation (if pointer passed).
xOptOptimal primal solution vector. A NULL pointer can be passed.
(If a null pointer is passed, the old primal solution is kept!)
yOptOptimal dual solution vector. A NULL pointer can be passed.
(If a null pointer is passed, the old dual solution is kept!)
guessedBoundsOptimal working set of bounds for solution (xOpt,yOpt).
(If a null pointer is passed, all bounds are assumed inactive!)
_RPre-computed (upper triangular) Cholesky factor of Hessian matrix. The Cholesky factor must be stored in a real_t array of size nV*nV in row-major format. Note: Only used if xOpt/yOpt and gB are NULL!
(If a null pointer is passed, Cholesky decomposition is computed internally!)

References BT_TRUE, getNV(), SubjectTo::getStatus(), isInitialised(), reset(), RET_INVALID_ARGUMENTS, RET_NO_CHOLESKY_WITH_INITIAL_GUESS, RET_QP_ALREADY_INITIALISED, RET_QPOBJECT_NOT_SETUP, setupQPdata(), solveInitialQP(), ST_UNDEFINED, SUCCESSFUL_RETURN, THROWERROR, and THROWWARNING.

returnValue QProblemB::init ( const char *const  H_file,
const char *const  g_file,
const char *const  lb_file,
const char *const  ub_file,
int_t nWSR,
real_t *const  cputime = 0,
const real_t *const  xOpt = 0,
const real_t *const  yOpt = 0,
const Bounds *const  guessedBounds = 0,
const char *const  R_file = 0 
)

Initialises a simply bounded QP problem with given QP data to be read from files and solves it using at most nWSR iterations. Depending on the parameter constellation it:
1. 0, 0, 0 : starts with xOpt = 0, yOpt = 0 and gB empty (or all implicit equality bounds),
2. xOpt, 0, 0 : starts with xOpt, yOpt = 0 and obtain gB by "clipping",
3. 0, yOpt, 0 : starts with xOpt = 0, yOpt and obtain gB from yOpt != 0,
4. 0, 0, gB: starts with xOpt = 0, yOpt = 0 and gB,
5. xOpt, yOpt, 0 : starts with xOpt, yOpt and obtain gB from yOpt != 0,
6. xOpt, 0, gB: starts with xOpt, yOpt = 0 and gB,
7. xOpt, yOpt, gB: starts with xOpt, yOpt and gB (assume them to be consistent!)

Note: This function internally calls solveInitialQP for initialisation!

Returns:
SUCCESSFUL_RETURN
RET_INIT_FAILED
RET_INIT_FAILED_CHOLESKY
RET_INIT_FAILED_HOTSTART
RET_INIT_FAILED_INFEASIBILITY
RET_INIT_FAILED_UNBOUNDEDNESS
RET_MAX_NWSR_REACHED
RET_UNABLE_TO_READ_FILE
Parameters:
H_fileName of file where Hessian matrix is stored.
If Hessian matrix is trivial, a NULL pointer can be passed.
g_fileName of file where gradient vector is stored.
lb_fileName of file where lower bound vector.
If no lower bounds exist, a NULL pointer can be passed.
ub_fileName of file where upper bound vector.
If no upper bounds exist, a NULL pointer can be passed.
nWSRInput: Maximum number of working set recalculations when using initial homotopy.
Output: Number of performed working set recalculations.
cputimeInput: Maximum CPU time allowed for QP initialisation.
Output: CPU time spent for QP initialisation (if pointer passed).
xOptOptimal primal solution vector. A NULL pointer can be passed.
(If a null pointer is passed, the old primal solution is kept!)
yOptOptimal dual solution vector. A NULL pointer can be passed.
(If a null pointer is passed, the old dual solution is kept!)
guessedBoundsOptimal working set of bounds for solution (xOpt,yOpt).
(If a null pointer is passed, all bounds are assumed inactive!)
R_fileName of the file where a pre-computed (upper triangular) Cholesky factor of the Hessian matrix is stored.
(If a null pointer is passed, Cholesky decomposition is computed internally!)

References BT_TRUE, getNV(), SubjectTo::getStatus(), isInitialised(), R, readFromFile(), reset(), RET_INVALID_ARGUMENTS, RET_NO_CHOLESKY_WITH_INITIAL_GUESS, RET_QP_ALREADY_INITIALISED, RET_QPOBJECT_NOT_SETUP, RET_UNABLE_TO_READ_FILE, setupQPdataFromFile(), solveInitialQP(), ST_UNDEFINED, SUCCESSFUL_RETURN, THROWERROR, and THROWWARNING.

BooleanType QProblemB::isBlocking ( real_t  num,
real_t  den,
real_t  epsNum,
real_t  epsDen,
real_t t 
) const [inline, protected]

Checks whether given ratio is blocking, i.e. limits the maximum step length along the homotopy path to a value lower than given one.

Returns:
SUCCESSFUL_RETURN
Parameters:
numNumerator for performing the ratio test.
denDenominator for performing the ratio test.
epsNumNumerator tolerance.
epsDenDenominator tolerance.
tInput: Current maximum step length along the homotopy path, Output: Updated maximum possible step length along the homotopy path.

References BT_FALSE, and BT_TRUE.

Referenced by performRatioTest().

BooleanType QProblemB::isCPUtimeLimitExceeded ( const real_t *const  cputime,
real_t  starttime,
int_t  nWSR 
) const [protected]

Determines if next QP iteration can be performed within given CPU time limit.

Returns:
BT_TRUE: CPU time limit is exceeded, stop QP solution.
BT_FALSE: Sufficient CPU time for next QP iteration.
Parameters:
cputimeMaximum CPU time allowed for QP solution.
starttimeStart time of current QP solution.
nWSRNumber of working set recalculations performed so far.

References BT_FALSE, BT_TRUE, getCPUtime(), and real_t.

Referenced by QProblem::solveQP(), and solveQP().

Returns if the QP is infeasible.

Returns:
BT_TRUE: QP infeasible
BT_FALSE: QP feasible (or not known to be infeasible!)

References infeasible.

Referenced by QProblem::solveInitialQP(), solveInitialQP(), and QProblem::solveQP().

Returns if the QProblem object is initialised.

Returns:
BT_TRUE: QProblemB initialised
BT_FALSE: QProblemB not initialised

References BT_FALSE, BT_TRUE, QPS_NOTINITIALISED, and status.

Referenced by init(), and QProblem::init().

BooleanType QProblemB::isSolved ( ) const [inline]

Returns if the QP has been solved.

Returns:
BT_TRUE: QProblemB solved
BT_FALSE: QProblemB not solved

References BT_FALSE, BT_TRUE, QPS_SOLVED, and status.

BooleanType QProblemB::isUnbounded ( ) const [inline]

Returns if the QP is unbounded.

Returns:
BT_TRUE: QP unbounded
BT_FALSE: QP unbounded (or not known to be unbounded!)

References unbounded.

Referenced by QProblem::solveInitialQP(), and solveInitialQP().

returnValue QProblemB::loadQPvectorsFromFile ( const char *const  g_file,
const char *const  lb_file,
const char *const  ub_file,
real_t *const  g_new,
real_t *const  lb_new,
real_t *const  ub_new 
) const [protected]

Loads new QP vectors from files (internal members are not affected!).

Returns:
SUCCESSFUL_RETURN
RET_UNABLE_TO_OPEN_FILE
RET_UNABLE_TO_READ_FILE
RET_INVALID_ARGUMENTS
Parameters:
g_fileName of file where gradient, of neighbouring QP to be solved, is stored.
lb_fileName of file where lower bounds, of neighbouring QP to be solved, is stored.
If no lower bounds exist, a NULL pointer can be passed.
ub_fileName of file where upper bounds, of neighbouring QP to be solved, is stored.
If no upper bounds exist, a NULL pointer can be passed.
g_newOutput: Gradient of neighbouring QP to be solved.
lb_newOutput: Lower bounds of neighbouring QP to be solved
ub_newOutput: Upper bounds of neighbouring QP to be solved

References getNV(), readFromFile(), RET_INVALID_ARGUMENTS, SUCCESSFUL_RETURN, and THROWERROR.

Referenced by hotstart().

returnValue QProblemB::obtainAuxiliaryWorkingSet ( const real_t *const  xOpt,
const real_t *const  yOpt,
const Bounds *const  guessedBounds,
Bounds auxiliaryBounds 
) const [protected]

Obtains the desired working set for the auxiliary initial QP in accordance with the user specifications

Returns:
SUCCESSFUL_RETURN
RET_OBTAINING_WORKINGSET_FAILED
RET_INVALID_ARGUMENTS
Parameters:
xOptOptimal primal solution vector. If a NULL pointer is passed, all entries are assumed to be zero.
yOptOptimal dual solution vector. If a NULL pointer is passed, all entries are assumed to be zero.
guessedBoundsGuessed working set for solution (xOpt,yOpt).
auxiliaryBoundsInput: Allocated bound object.
Output: Working set for auxiliary QP.

References bounds, Options::boundTolerance, EPS, getNV(), SubjectTo::getStatus(), SubjectTo::getType(), Options::initialStatusBounds, lb, options, RET_INVALID_ARGUMENTS, RET_OBTAINING_WORKINGSET_FAILED, Bounds::setupBound(), ST_EQUALITY, ST_INACTIVE, ST_LOWER, ST_UNBOUNDED, ST_UPPER, SUCCESSFUL_RETURN, THROWERROR, and ub.

Referenced by QProblem::obtainAuxiliaryWorkingSet(), and solveInitialQP().

QProblemB & QProblemB::operator= ( const QProblemB rhs) [virtual]

Assignment operator (deep copy).

Parameters:
rhsRhs object.

References clear(), and copy().

Drift correction at end of each active set iteration

Returns:
SUCCESSFUL_RETURN

Reimplemented in QProblem.

References bounds, getMax(), getMin(), getNV(), SubjectTo::getStatus(), SubjectTo::getType(), lb, setupAuxiliaryQPgradient(), ST_BOUNDED, ST_DISABLED, ST_EQUALITY, ST_INACTIVE, ST_INFEASIBLE_LOWER, ST_INFEASIBLE_UPPER, ST_LOWER, ST_UNBOUNDED, ST_UNDEFINED, ST_UNKNOWN, ST_UPPER, ub, x, and y.

Referenced by solveQP().

returnValue QProblemB::performRamping ( ) [protected, virtual]

Ramping Strategy to avoid ties. Modifies homotopy start without changing current active set.

Returns:
SUCCESSFUL_RETURN

Reimplemented in QProblem.

References bounds, getNV(), SubjectTo::getStatus(), SubjectTo::getType(), lb, ramp0, ramp1, rampOffset, real_t, setupAuxiliaryQPgradient(), ST_DISABLED, ST_EQUALITY, ST_INACTIVE, ST_LOWER, ST_UNBOUNDED, ST_UPPER, SUCCESSFUL_RETURN, ub, x, and y.

Referenced by solveQP().

returnValue QProblemB::performRatioTest ( int_t  nIdx,
const int_t *const  idxList,
const SubjectTo *const  subjectTo,
const real_t *const  num,
const real_t *const  den,
real_t  epsNum,
real_t  epsDen,
real_t t,
int_t BC_idx 
) const [protected]

Performs robustified ratio test yield the maximum possible step length along the homotopy path.

Returns:
SUCCESSFUL_RETURN
Parameters:
nIdxNumber of ratios to be checked.
idxListArray containing the indices of all ratios to be checked.
subjectToBound/Constraint object corresponding to ratios to be checked.
numArray containing all numerators for performing the ratio test.
denArray containing all denominators for performing the ratio test.
epsNumNumerator tolerance.
epsDenDenominator tolerance.
tOutput: Maximum possible step length along the homotopy path.
BC_idxOutput: Index of blocking constraint.

References BT_TRUE, SubjectTo::getStatus(), SubjectTo::getType(), isBlocking(), ST_EQUALITY, ST_INACTIVE, ST_LOWER, ST_UPPER, and SUCCESSFUL_RETURN.

Referenced by SQProblemSchur::addBound_ensureLI(), QProblem::addBound_ensureLI(), SQProblemSchur::addConstraint_ensureLI(), QProblem::addConstraint_ensureLI(), QProblem::performStep(), and performStep().

returnValue QProblemB::performStep ( const real_t *const  delta_g,
const real_t *const  delta_lb,
const real_t *const  delta_ub,
const real_t *const  delta_xFX,
const real_t *const  delta_xFR,
const real_t *const  delta_yFX,
int_t BC_idx,
SubjectToStatus BC_status 
) [private]

Determines the maximum possible step length along the homotopy path and performs this step (without changing working set).

Returns:
SUCCESSFUL_RETURN
RET_QP_INFEASIBLE
Parameters:
delta_gStep direction of gradient.
delta_lbStep direction of lower bounds.
delta_ubStep direction of upper bounds.
delta_xFXPrimal homotopy step direction of fixed variables.
delta_xFRPrimal homotopy step direction of free variables.
delta_yFXDual homotopy step direction of fixed variables' multiplier.
BC_idxOutput: Index of blocking constraint.
BC_statusOutput: Status of blocking constraint.

References __FILE__, __FUNC__, __LINE__, bounds, BT_FALSE, Options::epsDen, Options::epsNum, g, Bounds::getFixed(), Bounds::getFree(), getGlobalMessageHandler(), getMax(), getNFR(), getNFX(), Indexlist::getNumberArray(), getNV(), SubjectTo::hasNoLower(), SubjectTo::hasNoUpper(), lb, MAX_STRING_LENGTH, options, performRatioTest(), real_t, RET_STEPSIZE, RET_STEPSIZE_NONPOSITIVE, ST_INACTIVE, ST_LOWER, ST_UNDEFINED, ST_UPPER, SUCCESSFUL_RETURN, tau, MessageHandling::throwInfo(), MessageHandling::throwWarning(), ub, VS_VISIBLE, x, y, and ZERO.

Referenced by solveQP().

returnValue QProblemB::printIteration ( int_t  iter,
int_t  BC_idx,
SubjectToStatus  BC_status,
real_t  homotopyLength,
BooleanType  isFirstCall = BT_TRUE 
) [private]

Prints concise information on the current iteration.

Returns:
SUCCESSFUL_RETURN
Parameters:
iterNumber of current iteration.
BC_idxIndex of blocking bound.
BC_statusStatus of blocking bound.
homotopyLengthCurrent homotopy distance.
isFirstCallIndicating whether this is the first call for current QP.

References BT_TRUE, count, EPS, TabularOutput::excAddB, TabularOutput::excRemB, g, getAbs(), getNFX(), getNV(), H, hessianType, HST_ZERO, TabularOutput::idxAddB, TabularOutput::idxRemB, lb, MAX_STRING_LENGTH, myPrintf(), options, PL_DEBUG_ITER, PL_MEDIUM, PL_TABULAR, Options::printLevel, real_t, RET_INVALID_ARGUMENTS, ST_INACTIVE, ST_UNDEFINED, SUCCESSFUL_RETURN, tabularOutput, tau, THROWERROR, Matrix::times(), ub, x, and y.

Referenced by solveQP().

Prints a list of all options and their current values.

Returns:
SUCCESSFUL_RETURN

References options, and Options::print().

Referenced by main().

returnValue QProblemB::removeBound ( int_t  number,
BooleanType  updateCholesky 
) [private]

Resets QP problem counter (to zero).

Returns:
SUCCESSFUL_RETURN.

References count, and SUCCESSFUL_RETURN.

returnValue QProblemB::setG ( const real_t *const  g_new) [inline, protected]

Changes gradient vector of the QP.

Returns:
SUCCESSFUL_RETURN
RET_INVALID_ARGUMENTS
Parameters:
g_newNew gradient vector (with correct dimension!).

References g, getNV(), real_t, RET_INVALID_ARGUMENTS, RET_QPOBJECT_NOT_SETUP, SUCCESSFUL_RETURN, and THROWERROR.

Referenced by copy(), and setupQPdata().

returnValue QProblemB::setH ( SymmetricMatrix H_new) [inline, protected]

Sets Hessian matrix of the QP.

Returns:
SUCCESSFUL_RETURN
Parameters:
H_newNew Hessian matrix (a shallow copy is made).

References BT_FALSE, BT_TRUE, freeHessian, H, and SUCCESSFUL_RETURN.

Referenced by SQProblemSchur::setupAuxiliaryQP(), SQProblem::setupNewAuxiliaryQP(), setupQPdata(), and setupQPdataFromFile().

returnValue QProblemB::setH ( const real_t *const  H_new) [inline, protected]

Sets dense Hessian matrix of the QP. If a null pointer is passed and a) hessianType is HST_IDENTITY, nothing is done, b) hessianType is not HST_IDENTITY, Hessian matrix is set to zero.

Returns:
SUCCESSFUL_RETURN
Parameters:
H_newNew dense Hessian matrix (with correct dimension!), a shallow copy is made.

References BT_FALSE, BT_TRUE, freeHessian, getNV(), H, hessianType, HST_IDENTITY, HST_ZERO, real_t, and SUCCESSFUL_RETURN.

Changes the print level.

Returns:
SUCCESSFUL_RETURN
Parameters:
_hessianTypeNew Hessian type.

References hessianType, and SUCCESSFUL_RETURN.

returnValue QProblemB::setInfeasibilityFlag ( returnValue  returnvalue,
BooleanType  doThrowError = BT_FALSE 
) [protected]

Sets internal infeasibility flag and throws given error in case the far bound strategy is not enabled (as QP might actually not be infeasible in this case).

Returns:
RET_HOTSTART_STOPPED_INFEASIBILITY
RET_ENSURELI_FAILED_CYCLING
RET_ENSURELI_FAILED_NOINDEX
Parameters:
returnvalueReturnvalue to be tunneled.
doThrowErrorFlag forcing to throw an error.

References BT_FALSE, BT_TRUE, Options::enableFarBounds, infeasible, options, and THROWERROR.

Referenced by SQProblemSchur::addBound_ensureLI(), QProblem::addBound_ensureLI(), SQProblemSchur::addConstraint_ensureLI(), QProblem::addConstraint_ensureLI(), hotstart(), QProblem::hotstart(), QProblem::solveQP(), and solveQP().

returnValue QProblemB::setLB ( const real_t *const  lb_new) [inline, protected]

Changes lower bound vector of the QP.

Returns:
SUCCESSFUL_RETURN
RET_QPOBJECT_NOT_SETUP
Parameters:
lb_newNew lower bound vector (with correct dimension!).

References getNV(), INFTY, lb, real_t, RET_QPOBJECT_NOT_SETUP, SUCCESSFUL_RETURN, and THROWERROR.

Referenced by copy(), and setupQPdata().

returnValue QProblemB::setLB ( int_t  number,
real_t  value 
) [inline, protected]

Changes single entry of lower bound vector of the QP.

Returns:
SUCCESSFUL_RETURN
RET_QPOBJECT_NOT_SETUP
RET_INDEX_OUT_OF_BOUNDS
Parameters:
numberNumber of entry to be changed.
valueNew value for entry of lower bound vector.

References getNV(), lb, RET_INDEX_OUT_OF_BOUNDS, RET_QPOBJECT_NOT_SETUP, SUCCESSFUL_RETURN, and THROWERROR.

returnValue QProblemB::setOptions ( const Options _options) [inline]

Overrides current options with given ones.

Returns:
SUCCESSFUL_RETURN
Parameters:
_optionsNew options.

References Options::ensureConsistency(), options, Options::printLevel, setPrintLevel(), and SUCCESSFUL_RETURN.

Referenced by main(), and solveOqpBenchmark().

returnValue QProblemB::setUB ( const real_t *const  ub_new) [inline, protected]

Changes upper bound vector of the QP.

Returns:
SUCCESSFUL_RETURN
RET_QPOBJECT_NOT_SETUP
Parameters:
ub_newNew upper bound vector (with correct dimension!).

References getNV(), INFTY, real_t, RET_QPOBJECT_NOT_SETUP, SUCCESSFUL_RETURN, THROWERROR, and ub.

Referenced by copy(), and setupQPdata().

returnValue QProblemB::setUB ( int_t  number,
real_t  value 
) [inline, protected]

Changes single entry of upper bound vector of the QP.

Returns:
SUCCESSFUL_RETURN
RET_QPOBJECT_NOT_SETUP
RET_INDEX_OUT_OF_BOUNDS
Parameters:
numberNumber of entry to be changed.
valueNew value for entry of upper bound vector.

References getNV(), RET_INDEX_OUT_OF_BOUNDS, RET_QPOBJECT_NOT_SETUP, SUCCESSFUL_RETURN, THROWERROR, and ub.

returnValue QProblemB::setupAuxiliaryQP ( const Bounds *const  guessedBounds) [protected, virtual]

Updates QP vectors, working sets and internal data structures in order to start from an optimal solution corresponding to initial guesses of the working set for bounds

Returns:
SUCCESSFUL_RETURN
RET_SETUP_AUXILIARYQP_FAILED
Parameters:
guessedBoundsInitial guess for working set of bounds.

References bounds, BT_FALSE, BT_TRUE, computeCholesky(), getNV(), SubjectTo::getStatus(), Bounds::init(), QPS_PREPARINGAUXILIARYQP, RET_SETUP_AUXILIARYQP_FAILED, Bounds::setupAllFree(), setupAuxiliaryQPbounds(), setupAuxiliaryQPgradient(), setupAuxiliaryWorkingSet(), setupSubjectToType(), shallRefactorise(), ST_INACTIVE, status, SUCCESSFUL_RETURN, THROWERROR, and y.

Referenced by hotstart().

Sets up bounds of the auxiliary initial QP for given optimal primal/dual solution and given initial working set (assumes that members X, Y and BOUNDS have already been initialised!).

Returns:
SUCCESSFUL_RETURN
RET_UNKNOWN_BUG
Parameters:
useRelaxationFlag indicating if inactive bounds shall be relaxed.

References Options::boundRelaxation, bounds, BT_TRUE, getNV(), SubjectTo::getStatus(), SubjectTo::getType(), lb, options, RET_UNKNOWN_BUG, ST_EQUALITY, ST_INACTIVE, ST_INFEASIBLE_LOWER, ST_INFEASIBLE_UPPER, ST_LOWER, ST_UPPER, SUCCESSFUL_RETURN, THROWERROR, ub, and x.

Referenced by setupAuxiliaryQP(), and solveInitialQP().

Sets up gradient of the auxiliary initial QP for given optimal primal/dual solution and given initial working set (assumes that members X, Y and BOUNDS have already been initialised!).

Returns:
SUCCESSFUL_RETURN

Reimplemented in QProblem.

References BT_FALSE, g, getNV(), H, hessianType, HST_IDENTITY, HST_ZERO, regVal, SUCCESSFUL_RETURN, Matrix::times(), usingRegularisation(), x, and y.

Referenced by performDriftCorrection(), performRamping(), setupAuxiliaryQP(), and solveInitialQP().

returnValue QProblemB::setupAuxiliaryQPsolution ( const real_t *const  xOpt,
const real_t *const  yOpt 
) [private]

Sets up the optimal primal/dual solution of the auxiliary initial QP.

Returns:
SUCCESSFUL_RETURN
Parameters:
xOptOptimal primal solution vector. If a NULL pointer is passed, all entries are set to zero.
yOptOptimal dual solution vector. If a NULL pointer is passed, all entries are set to zero.

Reimplemented in QProblem.

References getNV(), SUCCESSFUL_RETURN, x, and y.

Referenced by solveInitialQP().

returnValue QProblemB::setupAuxiliaryWorkingSet ( const Bounds *const  auxiliaryBounds,
BooleanType  setupAfresh 
) [private]

Sets up bound data structure according to auxiliaryBounds. (If the working set shall be setup afresh, make sure that bounds data structure has been resetted!)

Returns:
SUCCESSFUL_RETURN
RET_SETUP_WORKINGSET_FAILED
RET_INVALID_ARGUMENTS
RET_UNKNOWN_BUG
Parameters:
auxiliaryBoundsWorking set for auxiliary QP.
setupAfreshFlag indicating if given working set shall be setup afresh or by updating the current one.

References addBound(), bounds, BT_FALSE, BT_TRUE, getNV(), SubjectTo::getStatus(), removeBound(), RET_INVALID_ARGUMENTS, RET_SETUP_WORKINGSET_FAILED, RET_UNKNOWN_BUG, ST_INACTIVE, ST_LOWER, ST_UNDEFINED, ST_UPPER, SUCCESSFUL_RETURN, and THROWERROR.

Referenced by setupAuxiliaryQP(), and solveInitialQP().

returnValue QProblemB::setupInitialCholesky ( ) [protected, virtual]

Computes initial Cholesky decomposition of the (simply projected) Hessian making use of the function computeCholesky().

Returns:
SUCCESSFUL_RETURN
RET_HESSIAN_NOT_SPD
RET_INDEXLIST_CORRUPTED

Reimplemented in QProblem.

References BT_TRUE, computeCholesky(), Options::enableRegularisation, getNFR(), getNFV(), getNV(), haveCholesky, options, regulariseHessian(), RET_HESSIAN_NOT_SPD, RET_INIT_FAILED_CHOLESKY, RET_INIT_FAILED_REGULARISATION, and SUCCESSFUL_RETURN.

Referenced by hotstart().

returnValue QProblemB::setupQPdata ( SymmetricMatrix _H,
const real_t *const  _g,
const real_t *const  _lb,
const real_t *const  _ub 
) [protected]

Sets up internal QP data.

Returns:
SUCCESSFUL_RETURN
RET_INVALID_ARGUMENTS
Parameters:
_HHessian matrix.
_gGradient vector.
_lbLower bounds (on variables).
If no lower bounds exist, a NULL pointer can be passed.
_ubUpper bounds (on variables).
If no upper bounds exist, a NULL pointer can be passed.

References RET_INVALID_ARGUMENTS, setG(), setH(), setLB(), setUB(), SUCCESSFUL_RETURN, and THROWERROR.

Referenced by init(), and QProblem::setupQPdata().

returnValue QProblemB::setupQPdata ( const real_t *const  _H,
const real_t *const  _g,
const real_t *const  _lb,
const real_t *const  _ub 
) [protected]

Sets up internal QP data. If the current Hessian is trivial (i.e. HST_ZERO or HST_IDENTITY) but a non-trivial one is given, memory for Hessian is allocated and it is set to the given one.

Returns:
SUCCESSFUL_RETURN
RET_INVALID_ARGUMENTS
RET_NO_HESSIAN_SPECIFIED
Parameters:
_HHessian matrix.
If Hessian matrix is trivial,a NULL pointer can be passed.
_gGradient vector.
_lbLower bounds (on variables).
If no lower bounds exist, a NULL pointer can be passed.
_ubUpper bounds (on variables).
If no upper bounds exist, a NULL pointer can be passed.

References RET_INVALID_ARGUMENTS, setG(), setH(), setLB(), setUB(), SUCCESSFUL_RETURN, and THROWERROR.

returnValue QProblemB::setupQPdataFromFile ( const char *const  H_file,
const char *const  g_file,
const char *const  lb_file,
const char *const  ub_file 
) [protected]

Sets up internal QP data by loading it from files. If the current Hessian is trivial (i.e. HST_ZERO or HST_IDENTITY) but a non-trivial one is given, memory for Hessian is allocated and it is set to the given one.

Returns:
SUCCESSFUL_RETURN
RET_UNABLE_TO_OPEN_FILE
RET_UNABLE_TO_READ_FILE
RET_INVALID_ARGUMENTS
RET_NO_HESSIAN_SPECIFIED
Parameters:
H_fileName of file where Hessian matrix, of neighbouring QP to be solved, is stored.
If Hessian matrix is trivial,a NULL pointer can be passed.
g_fileName of file where gradient, of neighbouring QP to be solved, is stored.
lb_fileName of file where lower bounds, of neighbouring QP to be solved, is stored.
If no lower bounds exist, a NULL pointer can be passed.
ub_fileName of file where upper bounds, of neighbouring QP to be solved, is stored.
If no upper bounds exist, a NULL pointer can be passed.

References Matrix::doFreeMemory(), g, getNV(), H, INFTY, lb, readFromFile(), real_t, RET_INVALID_ARGUMENTS, setH(), SUCCESSFUL_RETURN, THROWERROR, and ub.

Referenced by init().

returnValue QProblemB::setupSubjectToType ( ) [protected, virtual]

Determines type of existing constraints and bounds (i.e. implicitly fixed, unbounded etc.).

Returns:
SUCCESSFUL_RETURN
RET_SETUPSUBJECTTOTYPE_FAILED

Reimplemented in QProblem.

References lb, and ub.

Referenced by setupAuxiliaryQP(), QProblem::setupSubjectToType(), solveInitialQP(), solveQP(), and QProblem::updateActivitiesForHotstart().

returnValue QProblemB::setupSubjectToType ( const real_t *const  lb_new,
const real_t *const  ub_new 
) [protected, virtual]

Determines type of new constraints and bounds (i.e. implicitly fixed, unbounded etc.).

Returns:
SUCCESSFUL_RETURN
RET_SETUPSUBJECTTOTYPE_FAILED
Parameters:
lb_newNew lower bounds.
ub_newNew upper bounds.

References bounds, Options::boundTolerance, BT_FALSE, BT_TRUE, Options::enableEqualities, Options::enableFarBounds, getNV(), INFTY, lb, options, SubjectTo::setNoLower(), SubjectTo::setNoUpper(), SubjectTo::setType(), ST_BOUNDED, ST_EQUALITY, ST_UNBOUNDED, SUCCESSFUL_RETURN, and ub.

BooleanType QProblemB::shallRefactorise ( const Bounds *const  guessedBounds) const [private]

Determines if it is more efficient to refactorise the matrices when hotstarting or not (i.e. better to update the existing factorisations).

Returns:
BT_TRUE iff matrices shall be refactorised afresh
Parameters:
guessedBoundsGuessed new working set.

References bounds, BT_FALSE, BT_TRUE, Bounds::getNFX(), getNV(), SubjectTo::getStatus(), hessianType, HST_INDEF, and HST_SEMIDEF.

Referenced by setupAuxiliaryQP().

returnValue QProblemB::solveInitialQP ( const real_t *const  xOpt,
const real_t *const  yOpt,
const Bounds *const  guessedBounds,
const real_t *const  _R,
int_t nWSR,
real_t *const  cputime 
) [private]

Solves a QProblemB whose QP data is assumed to be stored in the member variables. A guess for its primal/dual optimal solution vectors and the corresponding optimal working set can be provided. Note: This function is internally called by all init functions!

Returns:
SUCCESSFUL_RETURN
RET_INIT_FAILED
RET_INIT_FAILED_CHOLESKY
RET_INIT_FAILED_HOTSTART
RET_INIT_FAILED_INFEASIBILITY
RET_INIT_FAILED_UNBOUNDEDNESS
RET_MAX_NWSR_REACHED
Parameters:
xOptOptimal primal solution vector.
yOptOptimal dual solution vector.
guessedBoundsOptimal working set of bounds for solution (xOpt,yOpt).
_RPre-computed (upper triangular) Cholesky factor of Hessian matrix.
nWSRInput: Maximum number of working set recalculations;
Output: Number of performed working set recalculations.
cputimeInput: Maximum CPU time allowed for QP solution.
Output: CPU time spent for QP solution (or to perform nWSR iterations).

References bounds, BT_FALSE, BT_TRUE, determineHessianType(), g, getCPUtime(), getNV(), haveCholesky, hessianType, hotstart(), HST_SEMIDEF, HST_ZERO, Options::initialStatusBounds, isInfeasible(), isUnbounded(), lb, obtainAuxiliaryWorkingSet(), options, QPS_AUXILIARYQPSOLVED, QPS_NOTINITIALISED, QPS_PREPARINGAUXILIARYQP, R, real_t, regulariseHessian(), RET_INIT_FAILED, RET_INIT_FAILED_HOTSTART, RET_INIT_FAILED_INFEASIBILITY, RET_INIT_FAILED_REGULARISATION, RET_INIT_FAILED_UNBOUNDEDNESS, RET_INIT_SUCCESSFUL, RET_MAX_NWSR_REACHED, RET_NO_CHOLESKY_WITH_INITIAL_GUESS, RR, Bounds::setupAllFree(), setupAuxiliaryQPbounds(), setupAuxiliaryQPgradient(), setupAuxiliaryQPsolution(), setupAuxiliaryWorkingSet(), setupSubjectToType(), ST_INACTIVE, status, SUCCESSFUL_RETURN, THROWERROR, THROWINFO, THROWWARNING, and ub.

Referenced by init().

returnValue QProblemB::solveQP ( const real_t *const  g_new,
const real_t *const  lb_new,
const real_t *const  ub_new,
int_t nWSR,
real_t *const  cputime,
int_t  nWSRperformed = 0,
BooleanType  isFirstCall = BT_TRUE 
) [private]

Solves an initialised QProblemB using online active set strategy. Note: This function is internally called by all hotstart functions!

Returns:
SUCCESSFUL_RETURN
RET_MAX_NWSR_REACHED
RET_HOTSTART_FAILED_AS_QP_NOT_INITIALISED
RET_HOTSTART_FAILED
RET_SHIFT_DETERMINATION_FAILED
RET_STEPDIRECTION_DETERMINATION_FAILED
RET_STEPLENGTH_DETERMINATION_FAILED
RET_HOMOTOPY_STEP_FAILED
RET_HOTSTART_STOPPED_INFEASIBILITY
RET_HOTSTART_STOPPED_UNBOUNDEDNESS
Parameters:
g_newGradient of neighbouring QP to be solved.
lb_newLower bounds of neighbouring QP to be solved.
If no lower bounds exist, a NULL pointer can be passed.
ub_newUpper bounds of neighbouring QP to be solved.
If no upper bounds exist, a NULL pointer can be passed.
nWSRInput: Maximum number of working set recalculations;
Output: Number of performed working set recalculations.
cputimeInput: Maximum CPU time allowed for QP solution.
Output: CPU time spent for QP solution (or to perform nWSR iterations).
nWSRperformedNumber of working set recalculations already performed to solve this QP within previous solveQP() calls. This number is always zero, except for successive calls from solveRegularisedQP() or when using the far bound strategy.
isFirstCallIndicating whether this is the first call for current QP.

References __FILE__, __FUNC__, __LINE__, BT_FALSE, BT_TRUE, changeActiveSet(), computeCholesky(), determineDataShift(), determineStepDirection(), Options::enableCholeskyRefactorisation, Options::enableDriftCorrection, Options::enableRamping, EPS, TabularOutput::excAddB, TabularOutput::excAddC, TabularOutput::excRemB, TabularOutput::excRemC, getCPUtime(), getGlobalMessageHandler(), getNV(), getRelativeHomotopyLength(), getStatus(), TabularOutput::idxAddB, TabularOutput::idxAddC, TabularOutput::idxRemB, TabularOutput::idxRemC, infeasible, isCPUtimeLimitExceeded(), MAX_STRING_LENGTH, options, performDriftCorrection(), performRamping(), performStep(), PL_HIGH, printIteration(), Options::printLevel, QPS_HOMOTOPYQPSOLVED, QPS_NOTINITIALISED, QPS_PERFORMINGHOMOTOPY, QPS_PREPARINGAUXILIARYQP, QPS_SOLVED, real_t, RET_HOMOTOPY_STEP_FAILED, RET_HOTSTART_FAILED, RET_HOTSTART_FAILED_AS_QP_NOT_INITIALISED, RET_HOTSTART_STOPPED_INFEASIBILITY, RET_HOTSTART_STOPPED_UNBOUNDEDNESS, RET_ITERATION_STARTED, RET_MAX_NWSR_REACHED, RET_OPTIMAL_SOLUTION_FOUND, RET_PRINT_ITERATION_FAILED, RET_SHIFT_DETERMINATION_FAILED, RET_STEPDIRECTION_DETERMINATION_FAILED, RET_STEPLENGTH_DETERMINATION_FAILED, setInfeasibilityFlag(), setupSubjectToType(), status, SUCCESSFUL_RETURN, tabularOutput, tau, Options::terminationTolerance, THROWERROR, MessageHandling::throwInfo(), THROWINFO, MessageHandling::throwWarning(), unbounded, and VS_VISIBLE.

Referenced by solveRegularisedQP().

returnValue QProblemB::solveRegularisedQP ( const real_t *const  g_new,
const real_t *const  lb_new,
const real_t *const  ub_new,
int_t nWSR,
real_t *const  cputime,
int_t  nWSRperformed = 0,
BooleanType  isFirstCall = BT_TRUE 
) [private]

Solves an initialised QProblemB using online active set strategy. Note: This function is internally called by all hotstart functions!

Returns:
SUCCESSFUL_RETURN
RET_MAX_NWSR_REACHED
RET_HOTSTART_FAILED_AS_QP_NOT_INITIALISED
RET_HOTSTART_FAILED
RET_SHIFT_DETERMINATION_FAILED
RET_STEPDIRECTION_DETERMINATION_FAILED
RET_STEPLENGTH_DETERMINATION_FAILED
RET_HOMOTOPY_STEP_FAILED
RET_HOTSTART_STOPPED_INFEASIBILITY
RET_HOTSTART_STOPPED_UNBOUNDEDNESS
Parameters:
g_newGradient of neighbouring QP to be solved.
lb_newLower bounds of neighbouring QP to be solved.
If no lower bounds exist, a NULL pointer can be passed.
ub_newUpper bounds of neighbouring QP to be solved.
If no upper bounds exist, a NULL pointer can be passed.
nWSRInput: Maximum number of working set recalculations;
Output: Number of performed working set recalculations.
cputimeInput: Maximum CPU time allowed for QP solution.
Output: CPU time spent for QP solution (or to perform nWSR iterations).
nWSRperformedNumber of working set recalculations already performed to solve this QP within previous solveRegularisedQP() calls. This number is always zero, except for successive calls when using the far bound strategy.
isFirstCallIndicating whether this is the first call for current QP.

References BT_FALSE, g, getNV(), Options::numRegularisationSteps, options, real_t, regVal, RET_FEWER_REGSTEPS_NWSR, RET_MAX_NWSR_REACHED, RET_NO_REGSTEP_NWSR, solveQP(), SUCCESSFUL_RETURN, THROWWARNING, usingRegularisation(), and x.

Referenced by hotstart().

returnValue QProblemB::updateFarBounds ( real_t  curFarBound,
int_t  nRamp,
const real_t *const  lb_new,
real_t *const  lb_new_far,
const real_t *const  ub_new,
real_t *const  ub_new_far 
) const [protected]

...

Parameters:
curFarBound...
nRamp...
lb_new...
lb_new_far...
ub_new...
ub_new_far...

References BT_TRUE, Options::enableRamping, getMax(), getMin(), getNV(), options, ramp0, ramp1, rampOffset, real_t, and SUCCESSFUL_RETURN.

Referenced by hotstart().


Member Data Documentation

Data structure for problem's bounds.

Referenced by SQProblemSchur::addBound(), QProblem::addBound(), addBound(), SQProblemSchur::addBound_checkLISchur(), SQProblemSchur::addBound_ensureLI(), QProblem::addBound_ensureLI(), QProblem::addConstraint(), QProblem::addConstraint_checkLI(), SQProblemSchur::addConstraint_checkLISchur(), SQProblemSchur::addConstraint_ensureLI(), QProblem::addConstraint_ensureLI(), SolutionAnalysis::checkCurvatureOnStronglyActiveConstraints(), computeCholesky(), QProblem::computeProjectedCholesky(), copy(), SQProblemSchur::correctInertia(), determineDataShift(), QProblem::determineDataShift(), QProblem::determineStepDirection(), determineStepDirection(), SQProblemSchur::determineStepDirection2(), QProblem::dropInfeasibles(), QProblem::ensureNonzeroCurvature(), getBounds(), QProblem::getFreeVariablesFlags(), getNFR(), getNFV(), getNFX(), getNV(), SolutionAnalysis::getVarianceCovariance(), getWorkingSetBounds(), QProblem::hotstart(), obtainAuxiliaryWorkingSet(), QProblem::performDriftCorrection(), performDriftCorrection(), performRamping(), QProblem::performRamping(), QProblem::performStep(), performStep(), QProblem::printProperties(), printProperties(), QProblemB(), SQProblemSchur::removeBound(), QProblem::removeBound(), removeBound(), QProblem::removeConstraint(), SQProblemSchur::repairSingularWorkingSet(), reset(), SQProblemSchur::resetSchurComplement(), SQProblemSchur::setupAuxiliaryQP(), QProblem::setupAuxiliaryQP(), setupAuxiliaryQP(), QProblem::setupAuxiliaryQPbounds(), setupAuxiliaryQPbounds(), SQProblemSchur::setupAuxiliaryWorkingSet(), QProblem::setupAuxiliaryWorkingSet(), setupAuxiliaryWorkingSet(), SQProblem::setupNewAuxiliaryQP(), setupSubjectToType(), QProblem::setupTQfactorisation(), QProblem::shallRefactorise(), shallRefactorise(), QProblem::solveCurrentEQP(), QProblem::solveInitialQP(), solveInitialQP(), QProblem::updateActivitiesForHotstart(), and QProblem::writeQpWorkspaceIntoMatFile().

uint_t QProblemB::count [protected]

Counts the number of hotstart function calls.

Referenced by copy(), getCount(), hotstart(), QProblem::hotstart(), QProblem::printIteration(), printIteration(), QProblemB(), and resetCounter().

Temporary for determineStepDirection

Referenced by clear(), copy(), QProblem::determineStepDirection(), determineStepDirection(), and QProblemB().

Struct for making a temporary copy of the matrix factorisations.

Referenced by copy(), QProblem::QProblem(), QProblemB(), QProblem::removeBound(), removeBound(), QProblem::removeConstraint(), QProblem::reset(), and reset().

Flag indicating whether the Hessian matrix needs to be de-allocated.

Referenced by clear(), copy(), QProblemB(), setH(), and SQProblem::setupNewAuxiliaryQP().

real_t* QProblemB::g [protected]

Flag indicating whether Cholesky decomposition has already been setup.

Referenced by copy(), hotstart(), QProblem::hotstart(), QProblemB(), reset(), setupInitialCholesky(), QProblem::setupInitialCholesky(), QProblem::solveInitialQP(), and solveInitialQP().

real_t* QProblemB::lb [protected]

Struct containing all user-defined options for solving QPs.

Referenced by SQProblemSchur::addBound(), QProblem::addBound_checkLI(), SQProblemSchur::addBound_checkLISchur(), SQProblemSchur::addBound_ensureLI(), QProblem::addBound_ensureLI(), SQProblemSchur::addConstraint(), QProblem::addConstraint_checkLI(), SQProblemSchur::addConstraint_checkLISchur(), SQProblemSchur::addConstraint_ensureLI(), QProblem::addConstraint_ensureLI(), SQProblemSchur::addToSchurComplement(), QProblem::changeActiveSet(), computeCholesky(), QProblem::computeProjectedCholesky(), copy(), SQProblemSchur::correctInertia(), SQProblemSchur::deleteFromSchurComplement(), determineHessianType(), QProblem::determineStepDirection(), determineStepDirection(), SQProblemSchur::determineStepDirection2(), QProblem::dropInfeasibles(), QProblem::ensureNonzeroCurvature(), getOptions(), getPrintLevel(), hotstart(), QProblem::hotstart(), obtainAuxiliaryWorkingSet(), QProblem::obtainAuxiliaryWorkingSet(), QProblem::performStep(), performStep(), QProblem::printIteration(), printIteration(), printOptions(), QProblem::printProperties(), printProperties(), QProblemB(), regulariseHessian(), SQProblemSchur::removeBound(), QProblem::removeBound(), removeBound(), SQProblemSchur::removeConstraint(), QProblem::removeConstraint(), SQProblemSchur::repairSingularWorkingSet(), reset(), SQProblemSchur::resetSchurComplement(), setInfeasibilityFlag(), setOptions(), setPrintLevel(), SQProblemSchur::setupAuxiliaryQP(), QProblem::setupAuxiliaryQPbounds(), setupAuxiliaryQPbounds(), SQProblemSchur::setupAuxiliaryWorkingSet(), QProblem::setupAuxiliaryWorkingSet(), setupInitialCholesky(), QProblem::setupInitialCholesky(), SQProblem::setupNewAuxiliaryQP(), setupSubjectToType(), QProblem::setupSubjectToType(), QProblem::solveInitialQP(), solveInitialQP(), QProblem::solveQP(), solveQP(), QProblem::solveRegularisedQP(), solveRegularisedQP(), SQProblemSchur::undoDeleteFromSchurComplement(), QProblem::updateActivitiesForHotstart(), updateFarBounds(), QProblem::updateFarBounds(), and SQProblemSchur::updateSchurQR().

real_t* QProblemB::R [protected]
real_t QProblemB::ramp0 [protected]
real_t QProblemB::ramp1 [protected]
real_t QProblemB::tau [protected]
real_t* QProblemB::ub [protected]
real_t* QProblemB::x [protected]
real_t* QProblemB::y [protected]

The documentation for this class was generated from the following files: