ConstrainedOptPack: C++ Tools for Constrained (and Unconstrained) Optimization  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Classes | Public Member Functions | Static Public Member Functions | Protected Types | Protected Member Functions | List of all members
ConstrainedOptPack::QPSchur Class Reference

Solves a Quadratic Program with a dual QP method using a schur complement factorization. More...

#include <ConstrainedOptPack_QPSchur.hpp>

Classes

class  ActiveSet
 Represents and manages the active set for the QPSchur algorithm. More...
 
class  DualInfeasibleException
 Thrown if during the course of the primal-dual iteration a non-dual feasible point if found. More...
 
class  InconsistantConstraintsException
 Thrown if constraints are inconsistant (no feasible region) More...
 
class  NumericalInstabilityException
 Thrown if there is some numerical instability. More...
 
class  TestFailed
 Thrown if a test failed. More...
 
class  U_hat_t
 Represents the matrix U_hat. More...
 

Public Member Functions

const ActiveSetact_set () const
 Return a reference to the active set object. More...
 

Static Public Member Functions

static void dump_act_set_quantities (const ActiveSet &act_set, std::ostream &out, bool print_S_hat=true)
 Dump all the active set quantities for debugging. More...
 

Protected Types

enum  EPDSteps
 
enum  EIterRefineReturn
 

Protected Member Functions

virtual ESolveReturn qp_algo (EPDSteps first_step, std::ostream *out, EOutputLevel output_level, ERunTests test_what, const DVectorSlice &vo, ActiveSet *act_set, DVectorSlice *v, DVectorSlice *x, size_type *iter, size_type *num_adds, size_type *num_drops, size_type *iter_refine_num_resid, size_type *iter_refine_num_solves, StopWatchPack::stopwatch *timer)
 Run the algorithm from a dual feasible point. More...
 
virtual void set_x (const ActiveSet &act_set, const DVectorSlice &v, DVectorSlice *x)
 Set the values in x for all the variables. More...
 
virtual void set_multipliers (const ActiveSet &act_set, const DVectorSlice &v, SpVector *mu, DVectorSlice *lambda, SpVector *lambda_breve)
 Map from the active set to the sparse multipliers for the inequality constraints. More...
 
bool timeout_return (StopWatchPack::stopwatch *timer, std::ostream *out, EOutputLevel output_level) const
 Determine if time has run out and if we should return. More...
 
EIterRefineReturn iter_refine (const ActiveSet &act_set, std::ostream *out, EOutputLevel output_level, const value_type ao, const DVectorSlice *bo, const value_type aa, const DVectorSlice *ba, DVectorSlice *v, DVectorSlice *z, size_type *iter_refine_num_resid, size_type *iter_refine_num_solves)
 Perform iterative refinement on the augmented KKT system for the current active set. More...
 

Public Types

enum  ERunTests
 Enumeration for if to run internal tests or not. More...
 
enum  ESolveReturn
 solve_qp return values More...
 
enum  EOutputLevel
 Output level. More...
 
typedef QPSchurPack::QP QP
 
typedef MatrixSymAddDelUpdateable MSADU
 
static value_type DEGENERATE_MULT = std::numeric_limits<value_type>::min()
 Value for near degenerate lagrange multipliers. More...
 

Public Member functions

 STANDARD_COMPOSITION_MEMBERS (MatrixSymAddDelUpdateableWithOpNonsingular, schur_comp)
 Schur complement matrix object S_hat. More...
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (size_type, max_iter)
 Set the maximum number of primal-dual QP iterations to take. More...
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, max_real_runtime)
 Set the maximum wall clock runtime (in minutes). More...
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, feas_tol)
 Set the feasibility tolerance for the constriants. More...
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, loose_feas_tol)
 Set a looser feasibility tolerance ( > feas_tol ) More...
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, dual_infeas_tol)
 Set the tolerence where a scaled Langrange multiplier is considered degenerate. More...
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, huge_primal_step)
 Set the tolerence for the size of the step in the primal space that is considered to be a near infinite step. This is used to determine if the KKT system is near singular. More...
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, huge_dual_step)
 Set the tolerence for the size of the step in the dual space that is considered to be a near infinite step. This is used to determine if the constriants are infeasible. More...
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, warning_tol)
 <<std member="" comp>="">> members for the warning tolerance for tests. More...
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, error_tol)
 <<std member="" comp>="">> members for the error tolerance for tests. More...
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (size_type, iter_refine_min_iter)
 Set the minimum number of refinement iterations to perform when using iterative refinement. More...
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (size_type, iter_refine_max_iter)
 Set the maximum number of refinement iterations to perform when using iterative refinement. More...
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, iter_refine_opt_tol)
 Set the maxinum scaled tolerance the residual of the optimality conditions must be before terminating iterative refinement. More...
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, iter_refine_feas_tol)
 Set the maxinum scaled tolerance the residual of the feasibility conditions must be before terminating iterative refinement. More...
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (bool, iter_refine_at_solution)
 Set whether iterative refinement is automatically used once the solution is found. More...
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (bool, salvage_init_schur_comp)
 Set whether a singular initial schur complement will attempted to be salvaged by adding as many nonsingular rows/cols as possible. More...
 
void pivot_tols (MSADU::PivotTolerances pivot_tols)
 Set the tolerances to use when updating the schur complement. More...
 
MSADU::PivotTolerances pivot_tols () const
 
virtual ~QPSchur ()
 
 QPSchur (const schur_comp_ptr_t &schur_comp=Teuchos::null, size_type max_iter=100, value_type max_real_runtime=1e+20, value_type feas_tol=1e-8, value_type loose_feas_tol=1e-6, value_type dual_infeas_tol=1e-12, value_type huge_primal_step=1e+20, value_type huge_dual_step=1e+20, value_type warning_tol=1e-10, value_type error_tol=1e-5, size_type iter_refine_min_iter=1, size_type iter_refine_max_iter=3, value_type iter_refine_opt_tol=1e-12, value_type iter_refine_feas_tol=1e-12, bool iter_refine_at_solution=true, bool salvage_init_schur_comp=true, MSADU::PivotTolerances pivot_tols=MSADU::PivotTolerances(1e-8, 1e-11, 1e-11))
 
virtual ESolveReturn solve_qp (QP &qp, size_type num_act_change, const int ij_act_change[], const EBounds bnds[], std::ostream *out, EOutputLevel output_level, ERunTests test_what, DVectorSlice *x, SpVector *mu, DVectorSlice *lambda, SpVector *lambda_breve, size_type *iter, size_type *num_adds, size_type *num_drops)
 Solve a QP. More...
 

Detailed Description

Solves a Quadratic Program with a dual QP method using a schur complement factorization.

See the paper "QPSchur: A Primal-Dual Active-Set Quadratic Programming Algorithm Using a Schur Complement Factorization Method" for a description of what this class does.

Definition at line 395 of file ConstrainedOptPack_QPSchur.hpp.

Member Typedef Documentation

Definition at line 402 of file ConstrainedOptPack_QPSchur.hpp.

typedef MatrixSymAddDelUpdateable ConstrainedOptPack::QPSchur::MSADU

Definition at line 404 of file ConstrainedOptPack_QPSchur.hpp.

Member Enumeration Documentation

Enumeration for if to run internal tests or not.

Definition at line 419 of file ConstrainedOptPack_QPSchur.hpp.

solve_qp return values

Definition at line 421 of file ConstrainedOptPack_QPSchur.hpp.

Output level.

Definition at line 433 of file ConstrainedOptPack_QPSchur.hpp.

Definition at line 1083 of file ConstrainedOptPack_QPSchur.hpp.

Definition at line 1118 of file ConstrainedOptPack_QPSchur.hpp.

Constructor & Destructor Documentation

virtual ConstrainedOptPack::QPSchur::~QPSchur ( )
inlinevirtual

Definition at line 530 of file ConstrainedOptPack_QPSchur.hpp.

ConstrainedOptPack::QPSchur::QPSchur ( const schur_comp_ptr_t &  schur_comp = Teuchos::null,
size_type  max_iter = 100,
value_type  max_real_runtime = 1e+20,
value_type  feas_tol = 1e-8,
value_type  loose_feas_tol = 1e-6,
value_type  dual_infeas_tol = 1e-12,
value_type  huge_primal_step = 1e+20,
value_type  huge_dual_step = 1e+20,
value_type  warning_tol = 1e-10,
value_type  error_tol = 1e-5,
size_type  iter_refine_min_iter = 1,
size_type  iter_refine_max_iter = 3,
value_type  iter_refine_opt_tol = 1e-12,
value_type  iter_refine_feas_tol = 1e-12,
bool  iter_refine_at_solution = true,
bool  salvage_init_schur_comp = true,
MSADU::PivotTolerances  pivot_tols = MSADU::PivotTolerances( 1e-8,1e-11,1e-11 ) 
)

Definition at line 2243 of file ConstrainedOptPack_QPSchur.cpp.

Member Function Documentation

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

Schur complement matrix object S_hat.

ConstrainedOptPack::QPSchur::STANDARD_MEMBER_COMPOSITION_MEMBERS ( size_type  ,
max_iter   
)

Set the maximum number of primal-dual QP iterations to take.

ConstrainedOptPack::QPSchur::STANDARD_MEMBER_COMPOSITION_MEMBERS ( value_type  ,
max_real_runtime   
)

Set the maximum wall clock runtime (in minutes).

ConstrainedOptPack::QPSchur::STANDARD_MEMBER_COMPOSITION_MEMBERS ( value_type  ,
feas_tol   
)

Set the feasibility tolerance for the constriants.

ConstrainedOptPack::QPSchur::STANDARD_MEMBER_COMPOSITION_MEMBERS ( value_type  ,
loose_feas_tol   
)

Set a looser feasibility tolerance ( > feas_tol )

ConstrainedOptPack::QPSchur::STANDARD_MEMBER_COMPOSITION_MEMBERS ( value_type  ,
dual_infeas_tol   
)

Set the tolerence where a scaled Langrange multiplier is considered degenerate.

ConstrainedOptPack::QPSchur::STANDARD_MEMBER_COMPOSITION_MEMBERS ( value_type  ,
huge_primal_step   
)

Set the tolerence for the size of the step in the primal space that is considered to be a near infinite step. This is used to determine if the KKT system is near singular.

ConstrainedOptPack::QPSchur::STANDARD_MEMBER_COMPOSITION_MEMBERS ( value_type  ,
huge_dual_step   
)

Set the tolerence for the size of the step in the dual space that is considered to be a near infinite step. This is used to determine if the constriants are infeasible.

ConstrainedOptPack::QPSchur::STANDARD_MEMBER_COMPOSITION_MEMBERS ( value_type  ,
warning_tol   
)

<<std member="" comp>="">> members for the warning tolerance for tests.

ConstrainedOptPack::QPSchur::STANDARD_MEMBER_COMPOSITION_MEMBERS ( value_type  ,
error_tol   
)

<<std member="" comp>="">> members for the error tolerance for tests.

ConstrainedOptPack::QPSchur::STANDARD_MEMBER_COMPOSITION_MEMBERS ( size_type  ,
iter_refine_min_iter   
)

Set the minimum number of refinement iterations to perform when using iterative refinement.

ConstrainedOptPack::QPSchur::STANDARD_MEMBER_COMPOSITION_MEMBERS ( size_type  ,
iter_refine_max_iter   
)

Set the maximum number of refinement iterations to perform when using iterative refinement.

ConstrainedOptPack::QPSchur::STANDARD_MEMBER_COMPOSITION_MEMBERS ( value_type  ,
iter_refine_opt_tol   
)

Set the maxinum scaled tolerance the residual of the optimality conditions must be before terminating iterative refinement.

ConstrainedOptPack::QPSchur::STANDARD_MEMBER_COMPOSITION_MEMBERS ( value_type  ,
iter_refine_feas_tol   
)

Set the maxinum scaled tolerance the residual of the feasibility conditions must be before terminating iterative refinement.

ConstrainedOptPack::QPSchur::STANDARD_MEMBER_COMPOSITION_MEMBERS ( bool  ,
iter_refine_at_solution   
)

Set whether iterative refinement is automatically used once the solution is found.

ConstrainedOptPack::QPSchur::STANDARD_MEMBER_COMPOSITION_MEMBERS ( bool  ,
salvage_init_schur_comp   
)

Set whether a singular initial schur complement will attempted to be salvaged by adding as many nonsingular rows/cols as possible.

void ConstrainedOptPack::QPSchur::pivot_tols ( MSADU::PivotTolerances  pivot_tols)

Set the tolerances to use when updating the schur complement.

Definition at line 2233 of file ConstrainedOptPack_QPSchur.cpp.

QPSchur::MSADU::PivotTolerances ConstrainedOptPack::QPSchur::pivot_tols ( ) const

Definition at line 2238 of file ConstrainedOptPack_QPSchur.cpp.

QPSchur::ESolveReturn ConstrainedOptPack::QPSchur::solve_qp ( QP qp,
size_type  num_act_change,
const int  ij_act_change[],
const EBounds  bnds[],
std::ostream *  out,
EOutputLevel  output_level,
ERunTests  test_what,
DVectorSlice *  x,
SpVector *  mu,
DVectorSlice *  lambda,
SpVector *  lambda_breve,
size_type iter,
size_type num_adds,
size_type num_drops 
)
virtual

Solve a QP.

If the initial schur complement turns out to have the wrong inertia then the QP is nonconvex, and the exception WrongInteriaUpdateExecption will be thrown. Otherwise, unless some other strange exception is thrown, this function will return normally (see return).

Parameters
qp[in] The abstraction for the QP being solved
num_act_change[in] The number of changes to the active set before the primal-dual QP algorithm starts.
ij_act_change[in] Array (size num_act_change): specifying how to initialize the active set. If i = -ij_act_change(s) > 0 then the initially fixed variable x(i) is to be freed. If j = ij_act_change(s) > 0 then the constraint a_bar(j)'*x is to be added to the active set to the bound bnd(s). The order of these changes can significantly effect the performance of the algorithm if these changes are not part of the optimal active set. Put changes that you are sure about earlier in the list and those that you are not a sure about later.
bnd[in] Array (size num_act_change): bnd(s) gives which bound to make active. If ij_act_change(s) < 0 then this is ignored.
out[out] output stream. Iteration information is printed according to output_level. If output_level == NO_OUTPUT then out may be NULL. If out==NULL, then output_level is forced to NO_OUTPUT
output_level[in] Specifies the level of output (see EOutputLevel).
test_what[in] Determines if internal validation tests are performed. The optimality conditions for the QP are not checked internally, since this is something that client can (and should) do independently. RUN_TESTS : As many validation/consistency tests are performed internally as possible. If a test fails then a TestFailed execption will be thrown. NO_TEST : No tests are performed internally. This is to allow the fastest possible execution.
x[out] vector (size qp.n()): Solution or current iteration value
mu[out] sparse vector (size qp.n()): Optimal lagrange multipliers for bound constraints. On output mu->is_sorted() == true.
lambda[out] vector (size q.m()): Optimal lagrange multipliers for equality constraints.
lambda_breve[out] sparse vector (size qp.constraints().m_breve()) for the active constraints in A_breve. On output lambda_breve->is_sorted() == true.
iter[out] The number of warm start drops and primal-dual iterations.
num_adds[out] The number of updates to the active set where a constraint was added. These do not include initially fixed variables.
num_drops[out] The number of updates to the active set where a constraint was dropped. These include constraints dropped during a warm start as well as during the primal-dual QP iterations.
Returns
Returns the type of point that x, mu , lambda and lambda_breve represents.

Definition at line 2281 of file ConstrainedOptPack_QPSchur.cpp.

const QPSchur::ActiveSet & ConstrainedOptPack::QPSchur::act_set ( ) const

Return a reference to the active set object.

Definition at line 2842 of file ConstrainedOptPack_QPSchur.cpp.

void ConstrainedOptPack::QPSchur::dump_act_set_quantities ( const ActiveSet act_set,
std::ostream &  out,
bool  print_S_hat = true 
)
static

Dump all the active set quantities for debugging.

Definition at line 4919 of file ConstrainedOptPack_QPSchur.cpp.

QPSchur::ESolveReturn ConstrainedOptPack::QPSchur::qp_algo ( EPDSteps  first_step,
std::ostream *  out,
EOutputLevel  output_level,
ERunTests  test_what,
const DVectorSlice &  vo,
ActiveSet act_set,
DVectorSlice *  v,
DVectorSlice *  x,
size_type iter,
size_type num_adds,
size_type num_drops,
size_type iter_refine_num_resid,
size_type iter_refine_num_solves,
StopWatchPack::stopwatch timer 
)
protectedvirtual

Run the algorithm from a dual feasible point.

By default, the algorithm should start with first_step = PICK_VIOLATED_CONSTRAINT if we are starting with a dual feasible point.

Definition at line 2849 of file ConstrainedOptPack_QPSchur.cpp.

void ConstrainedOptPack::QPSchur::set_x ( const ActiveSet act_set,
const DVectorSlice &  v,
DVectorSlice *  x 
)
protectedvirtual

Set the values in x for all the variables.

Definition at line 4515 of file ConstrainedOptPack_QPSchur.cpp.

void ConstrainedOptPack::QPSchur::set_multipliers ( const ActiveSet act_set,
const DVectorSlice &  v,
SpVector *  mu,
DVectorSlice *  lambda,
SpVector *  lambda_breve 
)
protectedvirtual

Map from the active set to the sparse multipliers for the inequality constraints.

Definition at line 4529 of file ConstrainedOptPack_QPSchur.cpp.

bool ConstrainedOptPack::QPSchur::timeout_return ( StopWatchPack::stopwatch timer,
std::ostream *  out,
EOutputLevel  output_level 
) const
protected

Determine if time has run out and if we should return.

Definition at line 4590 of file ConstrainedOptPack_QPSchur.cpp.

QPSchur::EIterRefineReturn ConstrainedOptPack::QPSchur::iter_refine ( const ActiveSet act_set,
std::ostream *  out,
EOutputLevel  output_level,
const value_type  ao,
const DVectorSlice *  bo,
const value_type  aa,
const DVectorSlice *  ba,
DVectorSlice *  v,
DVectorSlice *  z,
size_type iter_refine_num_resid,
size_type iter_refine_num_solves 
)
protected

Perform iterative refinement on the augmented KKT system for the current active set.

[   Ko     U_hat ] [ v ] + [ ao * bo ]
[ U_hat'   V_hat ] [ z ]   [ aa * ba ]

Returns true if iterative refinement satisfied the convergence criteria.

Definition at line 4605 of file ConstrainedOptPack_QPSchur.cpp.

Member Data Documentation

value_type ConstrainedOptPack::QPSchur::DEGENERATE_MULT = std::numeric_limits<value_type>::min()
static

Value for near degenerate lagrange multipliers.

Definition at line 442 of file ConstrainedOptPack_QPSchur.hpp.


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