MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Types | Public Member Functions | Private Types | Private Member Functions | Private Attributes | List of all members
ConstrainedOptPack::QPSchur::ActiveSet Class Reference

Represents and manages the active set for the QPSchur algorithm. More...

#include <ConstrainedOptPack_QPSchur.hpp>

Public Types

typedef QPSchurPack::QP QP
 
typedef MatrixSymAddDelUpdateable MSADU
 

Public Member Functions

 STANDARD_COMPOSITION_MEMBERS (MatrixSymAddDelUpdateableWithOpNonsingular, schur_comp)
 «std comp» members for schur complement matrix S_hat. More...
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (MSADU::PivotTolerances, pivot_tols)
 Set the tolerances to use when updating the schur complement. More...
 
 ActiveSet (const schur_comp_ptr_t &schur_comp, MSADU::PivotTolerances pivot_tols=MSADU::PivotTolerances(1e-6, 1e-8, 1e-8))
 

Private Types

typedef std::vector< int > ij_map_t
 
typedef std::map< int, size_type > s_map_t
 
typedef std::vector< EBoundsbnds_t
 
typedef std::vector< int > l_fxfx_t
 
typedef std::vector< size_type > P_row_t
 
typedef std::vector< size_type > P_col_t
 

Private Member Functions

void assert_initialized () const
 
void assert_s (size_type s) const
 
void reinitialize_matrices (bool test)
 
bool remove_augmented_element (size_type sd, bool force_refactorization, MatrixSymAddDelUpdateable::EEigenValType eigen_val_drop, std::ostream *out, EOutputLevel output_level, bool allow_any_cond)
 
 ActiveSet ()
 

Private Attributes

bool initialized_
 
bool test_
 
QPqp_
 
const QP::x_init_tx_init_
 
size_type n_
 
size_type n_R_
 
size_type m_
 
size_type m_breve_
 
size_type q_plus_hat_
 
size_type q_F_hat_
 
size_type q_C_hat_
 
ij_map_t ij_map_
 
DVector constr_norm_
 
bnds_t bnds_
 
l_fxfx_t l_fxfx_
 
U_hat_t U_hat_
 
GenPermMatrixSlice P_XF_hat_
 
P_row_t P_XF_hat_row_
 
P_row_t P_XF_hat_col_
 
GenPermMatrixSlice P_FC_hat_
 
P_row_t P_FC_hat_row_
 
P_row_t P_FC_hat_col_
 
GenPermMatrixSlice P_plus_hat_
 
P_row_t P_plus_hat_row_
 
P_row_t P_plus_hat_col_
 
GenPermMatrixSlice Q_XD_hat_
 
P_row_t Q_XD_hat_row_
 
P_row_t Q_XD_hat_col_
 
DVector d_hat_
 
DVector z_hat_
 
DVector p_z_hat_
 
DVector mu_D_hat_
 
DVector p_mu_D_hat_
 

Update the active set.

void initialize (QP &qp, size_type num_act_change, const int ij_act_change[], const EBounds bnds[], bool test, bool salvage_init_schur_comp, std::ostream *out, EOutputLevel output_level)
 Initialize with an additional active set. More...
 
void refactorize_schur_comp ()
 Reinitialize the schur complement factorization for the current active set. More...
 
bool add_constraint (size_type ja, EBounds bnd_ja, bool update_steps, std::ostream *out, EOutputLevel output_level, bool force_refactorization=true, bool allow_any_cond=false)
 Add a constraint to the active set then refactorize the schur complemnt (if forced). More...
 
bool drop_constraint (int jd, std::ostream *out, EOutputLevel output_level, bool force_refactorization=true, bool allow_any_cond=false)
 Drop a constraint from the active set then refactorize the schur complement (if forced). More...
 
bool drop_add_constraints (int jd, size_type ja, EBounds bnd_ja, bool update_steps, std::ostream *out, EOutputLevel output_level)
 Drop a constraint from, then add a constraint to the active set and refactorize the schur complement. More...
 

access the QP

QPqp ()
 
const QPqp () const
 

Access the active sets quantities.

size_type q_hat () const
 Return the total size of the schur complement. More...
 
size_type q_plus_hat () const
 Return the number of constraints from A_bar added to the active set. More...
 
size_type q_F_hat () const
 Return the number of variables that where initially fixed but are currently free or fixed to another bound. More...
 
size_type q_C_hat () const
 Return the number of variables that where initially fixed but are currently fixed to another bound. More...
 
size_type q_D_hat () const
 Return the number of variables that where initially fixed and are still currently fixed to their intial bounds. More...
 
int ij_map (size_type s) const
 Returns -i for row & column of S_bar for an initially fixed variable left out of Ko that became free and returns j for the constraint a(j)'*x that was added to the active set. More...
 
size_type s_map (int ij) const
 Map from a constraint or initially fixed variable to a row and column in the schur complement S_bar. More...
 
value_type constr_norm (size_type s) const
 Returns ||a(j)||2 where j = ij_map(s). More...
 
EBounds bnd (size_type s) const
 Return which bound is active for the active constraint. More...
 
size_type l_fxfx (size_type k) const
 Returns the indice of x_X(l) of the initially fixed variables that are still fixed at their original bounds. More...
 
const U_hat_tU_hat () const
 
const MatrixSymOpNonsing & S_hat () const
 
const GenPermMatrixSlice & P_XF_hat () const
 
const GenPermMatrixSlice & P_FC_hat () const
 
const GenPermMatrixSlice & P_plus_hat () const
 
const GenPermMatrixSlice & Q_XD_hat () const
 
const DVectorSlice d_hat () const
 
DVectorSlice z_hat ()
 
const DVectorSlice z_hat () const
 
DVectorSlice p_z_hat ()
 
const DVectorSlice p_z_hat () const
 
DVectorSlice mu_D_hat ()
 
const DVectorSlice mu_D_hat () const
 
DVectorSlice p_mu_D_hat ()
 
const DVectorSlice p_mu_D_hat () const
 
bool is_init_fixed (size_type j) const
 Determine if a constriant was an initially fixed variable. More...
 
bool all_dof_used_up () const
 Returns true if all the degrees of freedom of the QP are used up. More...
 

Detailed Description

Represents and manages the active set for the QPSchur algorithm.

This is a concrete type that encapsulates the maintaince of the active set and abstracts quantities associated with it.

At each iteration the dual active-set QP algorithm must solve the system:

 [ Ko       U_hat ]   [    v  ]   [  fo   ]
 [ U_hat'   V_hat ] * [ z_hat ] = [ d_hat ]

Above, U_hat contains the updates to the KKT system for adding constraints to the active set and freeing variables that where initially fixed and therefore left out of Ko.

This object maintains references to objects that represent the current augmented KKT system:

MatrixOp                      : U_hat        <: R^((n_R+m) x q_hat)
MatrixSymOp                   : V_hat        <: R^(q_hat x q_hat)
MatrixSymOpNonsing        : S_hat        <: R^(q_hat x q_hat)
GenPermMatrixSlice                : P_XF_hat     <: R^(n x q_hat)             (q_F_hat nonzeros)
GenPermMatrixSlice                : P_FC_hat     <: R^(q_hat x q_hat)         (q_C_hat nonzeros)
GenPermMatrixSlice                : P_plus_hat   <: R^((n+m_breve) x q_hat)   (q_plus_hat nonzeros)
GenPermMatrixSlice                : Q_XD_hat     <: R^(n x q_D_hat)           (q_D_hat nonzeros)
DVector                            : d_hat        <: R^(q_hat)
DVector                            : z_hat        <: R^(q_hat)

Definition at line 727 of file ConstrainedOptPack_QPSchur.hpp.

Member Typedef Documentation

Definition at line 734 of file ConstrainedOptPack_QPSchur.hpp.

typedef MatrixSymAddDelUpdateable ConstrainedOptPack::QPSchur::ActiveSet::MSADU

Definition at line 736 of file ConstrainedOptPack_QPSchur.hpp.

typedef std::vector<int> ConstrainedOptPack::QPSchur::ActiveSet::ij_map_t
private

Definition at line 965 of file ConstrainedOptPack_QPSchur.hpp.

typedef std::map<int,size_type> ConstrainedOptPack::QPSchur::ActiveSet::s_map_t
private

Definition at line 967 of file ConstrainedOptPack_QPSchur.hpp.

Definition at line 969 of file ConstrainedOptPack_QPSchur.hpp.

typedef std::vector<int> ConstrainedOptPack::QPSchur::ActiveSet::l_fxfx_t
private

Definition at line 971 of file ConstrainedOptPack_QPSchur.hpp.

typedef std::vector<size_type> ConstrainedOptPack::QPSchur::ActiveSet::P_row_t
private

Definition at line 973 of file ConstrainedOptPack_QPSchur.hpp.

typedef std::vector<size_type> ConstrainedOptPack::QPSchur::ActiveSet::P_col_t
private

Definition at line 975 of file ConstrainedOptPack_QPSchur.hpp.

Constructor & Destructor Documentation

ConstrainedOptPack::QPSchur::ActiveSet::ActiveSet ( const schur_comp_ptr_t &  schur_comp,
MSADU::PivotTolerances  pivot_tols = MSADU::PivotTolerances( 1e-6,1e-8,1e-8 ) 
)

Definition at line 1021 of file ConstrainedOptPack_QPSchur.cpp.

ConstrainedOptPack::QPSchur::ActiveSet::ActiveSet ( )
private

Member Function Documentation

ConstrainedOptPack::QPSchur::ActiveSet::STANDARD_COMPOSITION_MEMBERS ( MatrixSymAddDelUpdateableWithOpNonsingular  ,
schur_comp   
)

«std comp» members for schur complement matrix S_hat.

Warning: Resetting schur_comp will cause a reinitialization to an empty active set.

ConstrainedOptPack::QPSchur::ActiveSet::STANDARD_MEMBER_COMPOSITION_MEMBERS ( MSADU::PivotTolerances  ,
pivot_tols   
)

Set the tolerances to use when updating the schur complement.

void ConstrainedOptPack::QPSchur::ActiveSet::initialize ( QP qp,
size_type  num_act_change,
const int  ij_act_change[],
const EBounds  bnds[],
bool  test,
bool  salvage_init_schur_comp,
std::ostream *  out,
EOutputLevel  output_level 
)

Initialize with an additional active set.

If the initial schur complement is not full rank then an LDConstraintException exception will be thrown. The active set will contain all of the constraints it can such that the schur complement is nonsingular.

Definition at line 1039 of file ConstrainedOptPack_QPSchur.cpp.

void ConstrainedOptPack::QPSchur::ActiveSet::refactorize_schur_comp ( )

Reinitialize the schur complement factorization for the current active set.

ToDo: Finish documentation

Definition at line 1471 of file ConstrainedOptPack_QPSchur.cpp.

bool ConstrainedOptPack::QPSchur::ActiveSet::add_constraint ( size_type  ja,
EBounds  bnd_ja,
bool  update_steps,
std::ostream *  out,
EOutputLevel  output_level,
bool  force_refactorization = true,
bool  allow_any_cond = false 
)

Add a constraint to the active set then refactorize the schur complemnt (if forced).

ToDo: Finish documentation

If the new KKT system is singular then the exeption MatrixSymAddDelUpdateable::SingularUpdateException will be thrown but the old KKT system will be kept intact.

If the reduced Hessian for the new KKT system does not have the correct inertia then the exception MatrixSymAddDelUpdateable::WrongInertiaUpdateException will be thrown but the old KKT system will be kept intact.

Returns
Returns true if any output was sent to *out.

Definition at line 1477 of file ConstrainedOptPack_QPSchur.cpp.

bool ConstrainedOptPack::QPSchur::ActiveSet::drop_constraint ( int  jd,
std::ostream *  out,
EOutputLevel  output_level,
bool  force_refactorization = true,
bool  allow_any_cond = false 
)

Drop a constraint from the active set then refactorize the schur complement (if forced).

ToDo: Finish documentation

Returns true if any output was sent to *out.

Definition at line 1716 of file ConstrainedOptPack_QPSchur.cpp.

bool ConstrainedOptPack::QPSchur::ActiveSet::drop_add_constraints ( int  jd,
size_type  ja,
EBounds  bnd_ja,
bool  update_steps,
std::ostream *  out,
EOutputLevel  output_level 
)

Drop a constraint from, then add a constraint to the active set and refactorize the schur complement.

ToDo: Finish documentation

Returns true if any output was sent to *out.

Definition at line 1940 of file ConstrainedOptPack_QPSchur.cpp.

QPSchur::ActiveSet::QP & ConstrainedOptPack::QPSchur::ActiveSet::qp ( )

Definition at line 1954 of file ConstrainedOptPack_QPSchur.cpp.

const QPSchur::ActiveSet::QP & ConstrainedOptPack::QPSchur::ActiveSet::qp ( ) const

Definition at line 1961 of file ConstrainedOptPack_QPSchur.cpp.

size_type ConstrainedOptPack::QPSchur::ActiveSet::q_hat ( ) const

Return the total size of the schur complement.

q_hat = q_plus_hat + q_F_hat + q_C_hat.

Definition at line 1967 of file ConstrainedOptPack_QPSchur.cpp.

size_type ConstrainedOptPack::QPSchur::ActiveSet::q_plus_hat ( ) const

Return the number of constraints from A_bar added to the active set.

Definition at line 1973 of file ConstrainedOptPack_QPSchur.cpp.

size_type ConstrainedOptPack::QPSchur::ActiveSet::q_F_hat ( ) const

Return the number of variables that where initially fixed but are currently free or fixed to another bound.

Definition at line 1979 of file ConstrainedOptPack_QPSchur.cpp.

size_type ConstrainedOptPack::QPSchur::ActiveSet::q_C_hat ( ) const

Return the number of variables that where initially fixed but are currently fixed to another bound.

Definition at line 1985 of file ConstrainedOptPack_QPSchur.cpp.

size_type ConstrainedOptPack::QPSchur::ActiveSet::q_D_hat ( ) const

Return the number of variables that where initially fixed and are still currently fixed to their intial bounds.

Definition at line 1991 of file ConstrainedOptPack_QPSchur.cpp.

int ConstrainedOptPack::QPSchur::ActiveSet::ij_map ( size_type  s) const

Returns -i for row & column of S_bar for an initially fixed variable left out of Ko that became free and returns j for the constraint a(j)'*x that was added to the active set.

1 <= s <= q_hat

Definition at line 1997 of file ConstrainedOptPack_QPSchur.cpp.

size_type ConstrainedOptPack::QPSchur::ActiveSet::s_map ( int  ij) const

Map from a constraint or initially fixed variable to a row and column in the schur complement S_bar.

To determine if an initially fixed variable x(i) is now free call s_map(-i). If s_map(-i) returns zero then x(i) is still fixed. Otherwise s_map(-i) returns the row and column in S_bar for this change in the active set.

To determine if a constraint a(j)'*x is part of the active set call s_map(j). If s_map(j) returns zero then a(j)'*x is not part of the active set. Otherwise s_map(j) returns the row and column in S_bar for this change in the active set.

Definition at line 2003 of file ConstrainedOptPack_QPSchur.cpp.

value_type ConstrainedOptPack::QPSchur::ActiveSet::constr_norm ( size_type  s) const

Returns ||a(j)||2 where j = ij_map(s).

If ij_map(s) < 0, the this function returns zero.

1 <= s <= q_hat

Definition at line 2012 of file ConstrainedOptPack_QPSchur.cpp.

EBounds ConstrainedOptPack::QPSchur::ActiveSet::bnd ( size_type  s) const

Return which bound is active for the active constraint.

Definition at line 2018 of file ConstrainedOptPack_QPSchur.cpp.

size_type ConstrainedOptPack::QPSchur::ActiveSet::l_fxfx ( size_type  k) const

Returns the indice of x_X(l) of the initially fixed variables that are still fixed at their original bounds.

i <= k <= q_D_hat

Definition at line 2024 of file ConstrainedOptPack_QPSchur.cpp.

const QPSchur::U_hat_t & ConstrainedOptPack::QPSchur::ActiveSet::U_hat ( ) const

Definition at line 2030 of file ConstrainedOptPack_QPSchur.cpp.

const MatrixSymOpNonsing & ConstrainedOptPack::QPSchur::ActiveSet::S_hat ( ) const

Definition at line 2036 of file ConstrainedOptPack_QPSchur.cpp.

const GenPermMatrixSlice & ConstrainedOptPack::QPSchur::ActiveSet::P_XF_hat ( ) const

Definition at line 2042 of file ConstrainedOptPack_QPSchur.cpp.

const GenPermMatrixSlice & ConstrainedOptPack::QPSchur::ActiveSet::P_FC_hat ( ) const

Definition at line 2048 of file ConstrainedOptPack_QPSchur.cpp.

const GenPermMatrixSlice & ConstrainedOptPack::QPSchur::ActiveSet::P_plus_hat ( ) const

Definition at line 2054 of file ConstrainedOptPack_QPSchur.cpp.

const GenPermMatrixSlice & ConstrainedOptPack::QPSchur::ActiveSet::Q_XD_hat ( ) const

Definition at line 2060 of file ConstrainedOptPack_QPSchur.cpp.

const DVectorSlice ConstrainedOptPack::QPSchur::ActiveSet::d_hat ( ) const

Definition at line 2066 of file ConstrainedOptPack_QPSchur.cpp.

DVectorSlice ConstrainedOptPack::QPSchur::ActiveSet::z_hat ( )

Definition at line 2072 of file ConstrainedOptPack_QPSchur.cpp.

const DVectorSlice ConstrainedOptPack::QPSchur::ActiveSet::z_hat ( ) const

Definition at line 2078 of file ConstrainedOptPack_QPSchur.cpp.

DVectorSlice ConstrainedOptPack::QPSchur::ActiveSet::p_z_hat ( )

Definition at line 2084 of file ConstrainedOptPack_QPSchur.cpp.

const DVectorSlice ConstrainedOptPack::QPSchur::ActiveSet::p_z_hat ( ) const

Definition at line 2090 of file ConstrainedOptPack_QPSchur.cpp.

DVectorSlice ConstrainedOptPack::QPSchur::ActiveSet::mu_D_hat ( )

Definition at line 2096 of file ConstrainedOptPack_QPSchur.cpp.

const DVectorSlice ConstrainedOptPack::QPSchur::ActiveSet::mu_D_hat ( ) const

Definition at line 2102 of file ConstrainedOptPack_QPSchur.cpp.

DVectorSlice ConstrainedOptPack::QPSchur::ActiveSet::p_mu_D_hat ( )

Definition at line 2108 of file ConstrainedOptPack_QPSchur.cpp.

const DVectorSlice ConstrainedOptPack::QPSchur::ActiveSet::p_mu_D_hat ( ) const

Definition at line 2114 of file ConstrainedOptPack_QPSchur.cpp.

bool ConstrainedOptPack::QPSchur::ActiveSet::is_init_fixed ( size_type  j) const

Determine if a constriant was an initially fixed variable.

This function will return true if:

j <= n && x_init(j) != FREE

This is just a function of convienience

Definition at line 2120 of file ConstrainedOptPack_QPSchur.cpp.

bool ConstrainedOptPack::QPSchur::ActiveSet::all_dof_used_up ( ) const

Returns true if all the degrees of freedom of the QP are used up.

Definition at line 2126 of file ConstrainedOptPack_QPSchur.cpp.

void ConstrainedOptPack::QPSchur::ActiveSet::assert_initialized ( ) const
private

Definition at line 2133 of file ConstrainedOptPack_QPSchur.cpp.

void ConstrainedOptPack::QPSchur::ActiveSet::assert_s ( size_type  s) const
private

Definition at line 2141 of file ConstrainedOptPack_QPSchur.cpp.

void ConstrainedOptPack::QPSchur::ActiveSet::reinitialize_matrices ( bool  test)
private

Definition at line 2146 of file ConstrainedOptPack_QPSchur.cpp.

bool ConstrainedOptPack::QPSchur::ActiveSet::remove_augmented_element ( size_type  sd,
bool  force_refactorization,
MatrixSymAddDelUpdateable::EEigenValType  eigen_val_drop,
std::ostream *  out,
EOutputLevel  output_level,
bool  allow_any_cond 
)
private

Definition at line 2182 of file ConstrainedOptPack_QPSchur.cpp.

Member Data Documentation

bool ConstrainedOptPack::QPSchur::ActiveSet::initialized_
private

Definition at line 980 of file ConstrainedOptPack_QPSchur.hpp.

bool ConstrainedOptPack::QPSchur::ActiveSet::test_
private

Definition at line 981 of file ConstrainedOptPack_QPSchur.hpp.

QP* ConstrainedOptPack::QPSchur::ActiveSet::qp_
private

Definition at line 982 of file ConstrainedOptPack_QPSchur.hpp.

const QP::x_init_t* ConstrainedOptPack::QPSchur::ActiveSet::x_init_
private

Definition at line 983 of file ConstrainedOptPack_QPSchur.hpp.

size_type ConstrainedOptPack::QPSchur::ActiveSet::n_
private

Definition at line 984 of file ConstrainedOptPack_QPSchur.hpp.

size_type ConstrainedOptPack::QPSchur::ActiveSet::n_R_
private

Definition at line 985 of file ConstrainedOptPack_QPSchur.hpp.

size_type ConstrainedOptPack::QPSchur::ActiveSet::m_
private

Definition at line 986 of file ConstrainedOptPack_QPSchur.hpp.

size_type ConstrainedOptPack::QPSchur::ActiveSet::m_breve_
private

Definition at line 987 of file ConstrainedOptPack_QPSchur.hpp.

size_type ConstrainedOptPack::QPSchur::ActiveSet::q_plus_hat_
private

Definition at line 988 of file ConstrainedOptPack_QPSchur.hpp.

size_type ConstrainedOptPack::QPSchur::ActiveSet::q_F_hat_
private

Definition at line 989 of file ConstrainedOptPack_QPSchur.hpp.

size_type ConstrainedOptPack::QPSchur::ActiveSet::q_C_hat_
private

Definition at line 990 of file ConstrainedOptPack_QPSchur.hpp.

ij_map_t ConstrainedOptPack::QPSchur::ActiveSet::ij_map_
private

Definition at line 991 of file ConstrainedOptPack_QPSchur.hpp.

DVector ConstrainedOptPack::QPSchur::ActiveSet::constr_norm_
private

Definition at line 993 of file ConstrainedOptPack_QPSchur.hpp.

bnds_t ConstrainedOptPack::QPSchur::ActiveSet::bnds_
private

Definition at line 994 of file ConstrainedOptPack_QPSchur.hpp.

l_fxfx_t ConstrainedOptPack::QPSchur::ActiveSet::l_fxfx_
private

Definition at line 995 of file ConstrainedOptPack_QPSchur.hpp.

U_hat_t ConstrainedOptPack::QPSchur::ActiveSet::U_hat_
private

Definition at line 996 of file ConstrainedOptPack_QPSchur.hpp.

GenPermMatrixSlice ConstrainedOptPack::QPSchur::ActiveSet::P_XF_hat_
private

Definition at line 1004 of file ConstrainedOptPack_QPSchur.hpp.

P_row_t ConstrainedOptPack::QPSchur::ActiveSet::P_XF_hat_row_
private

Definition at line 1005 of file ConstrainedOptPack_QPSchur.hpp.

P_row_t ConstrainedOptPack::QPSchur::ActiveSet::P_XF_hat_col_
private

Definition at line 1006 of file ConstrainedOptPack_QPSchur.hpp.

GenPermMatrixSlice ConstrainedOptPack::QPSchur::ActiveSet::P_FC_hat_
private

Definition at line 1015 of file ConstrainedOptPack_QPSchur.hpp.

P_row_t ConstrainedOptPack::QPSchur::ActiveSet::P_FC_hat_row_
private

Definition at line 1016 of file ConstrainedOptPack_QPSchur.hpp.

P_row_t ConstrainedOptPack::QPSchur::ActiveSet::P_FC_hat_col_
private

Definition at line 1017 of file ConstrainedOptPack_QPSchur.hpp.

GenPermMatrixSlice ConstrainedOptPack::QPSchur::ActiveSet::P_plus_hat_
private

Definition at line 1025 of file ConstrainedOptPack_QPSchur.hpp.

P_row_t ConstrainedOptPack::QPSchur::ActiveSet::P_plus_hat_row_
private

Definition at line 1026 of file ConstrainedOptPack_QPSchur.hpp.

P_row_t ConstrainedOptPack::QPSchur::ActiveSet::P_plus_hat_col_
private

Definition at line 1027 of file ConstrainedOptPack_QPSchur.hpp.

GenPermMatrixSlice ConstrainedOptPack::QPSchur::ActiveSet::Q_XD_hat_
private

Definition at line 1033 of file ConstrainedOptPack_QPSchur.hpp.

P_row_t ConstrainedOptPack::QPSchur::ActiveSet::Q_XD_hat_row_
private

Definition at line 1034 of file ConstrainedOptPack_QPSchur.hpp.

P_row_t ConstrainedOptPack::QPSchur::ActiveSet::Q_XD_hat_col_
private

Definition at line 1035 of file ConstrainedOptPack_QPSchur.hpp.

DVector ConstrainedOptPack::QPSchur::ActiveSet::d_hat_
private

Definition at line 1037 of file ConstrainedOptPack_QPSchur.hpp.

DVector ConstrainedOptPack::QPSchur::ActiveSet::z_hat_
private

Definition at line 1038 of file ConstrainedOptPack_QPSchur.hpp.

DVector ConstrainedOptPack::QPSchur::ActiveSet::p_z_hat_
private

Definition at line 1039 of file ConstrainedOptPack_QPSchur.hpp.

DVector ConstrainedOptPack::QPSchur::ActiveSet::mu_D_hat_
private

Definition at line 1040 of file ConstrainedOptPack_QPSchur.hpp.

DVector ConstrainedOptPack::QPSchur::ActiveSet::p_mu_D_hat_
private

Definition at line 1041 of file ConstrainedOptPack_QPSchur.hpp.


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