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

Tests the optimality conditions of the output from a QPSolverRelaxed object. More...

#include <ConstrainedOptPack_QPSolverRelaxedTester.hpp>

Public Member Functions

 STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, opt_warning_tol)
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, opt_error_tol)
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, feas_warning_tol)
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, feas_error_tol)
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, comp_warning_tol)
 
 STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, comp_error_tol)
 
 QPSolverRelaxedTester (value_type opt_warning_tol=1e-10, value_type opt_error_tol=1e-5, value_type feas_warning_tol=1e-10, value_type feas_error_tol=1e-5, value_type comp_warning_tol=1e-10, value_type comp_error_tol=1e-5)
 
virtual ~QPSolverRelaxedTester ()
 
virtual bool check_optimality_conditions (QPSolverStats::ESolutionType solution_type, const value_type infinite_bound, std::ostream *out, bool print_all_warnings, bool print_vectors, const Vector &g, const MatrixSymOp &G, value_type etaL, const Vector &dL, const Vector &dU, const MatrixOp &E, BLAS_Cpp::Transp trans_E, const Vector &b, const Vector &eL, const Vector &eU, const MatrixOp &F, BLAS_Cpp::Transp trans_F, const Vector &f, const value_type *obj_d, const value_type *eta, const Vector *d, const Vector *nu, const Vector *mu, const Vector *Ed, const Vector *lambda, const Vector *Fd)
 Check the optimality conditions for the solved (or partially solved) QP. More...
 
virtual bool check_optimality_conditions (QPSolverStats::ESolutionType solution_type, const value_type infinite_bound, std::ostream *out, bool print_all_warnings, bool print_vectors, const Vector &g, const MatrixSymOp &G, value_type etaL, const Vector &dL, const Vector &dU, const MatrixOp &E, BLAS_Cpp::Transp trans_E, const Vector &b, const Vector &eL, const Vector &eU, const value_type *obj_d, const value_type *eta, const Vector *d, const Vector *nu, const Vector *mu, const Vector *Ed)
 Check the optimality conditions without general equality constrants. More...
 
virtual bool check_optimality_conditions (QPSolverStats::ESolutionType solution_type, const value_type infinite_bound, std::ostream *out, bool print_all_warnings, bool print_vectors, const Vector &g, const MatrixSymOp &G, value_type etaL, const Vector &dL, const Vector &dU, const MatrixOp &F, BLAS_Cpp::Transp trans_F, const Vector &f, const value_type *obj_d, const value_type *eta, const Vector *d, const Vector *nu, const Vector *lambda, const Vector *Fd)
 Check the optimality conditions general inequality constrants. More...
 
virtual bool check_optimality_conditions (QPSolverStats::ESolutionType solution_type, const value_type infinite_bound, std::ostream *out, bool print_all_warnings, bool print_vectors, const Vector &g, const MatrixSymOp &G, const Vector &dL, const Vector &dU, const value_type *obj_d, const Vector *d, const Vector *nu)
 Check the optimality conditions without general equality or inequality constrants (no relaxation needed). More...
 
virtual bool check_optimality_conditions (QPSolverStats::ESolutionType solution_type, const value_type infinite_bound, std::ostream *out, bool print_all_warnings, bool print_vectors, const Vector &g, const MatrixSymOp &G, value_type etaL, const Vector *dL, const Vector *dU, const MatrixOp *E, BLAS_Cpp::Transp trans_E, const Vector *b, const Vector *eL, const Vector *eU, const MatrixOp *F, BLAS_Cpp::Transp trans_F, const Vector *f, const value_type *obj_d, const value_type *eta, const Vector *d, const Vector *nu, const Vector *mu, const Vector *Ed, const Vector *lambda, const Vector *Fd)
 This is a more flexible function where the client can set different constraints to be included. More...
 

Protected Member Functions

virtual bool imp_check_optimality_conditions (QPSolverStats::ESolutionType solution_type, const value_type infinite_bound, std::ostream *out, bool print_all_warnings, bool print_vectors, const Vector &g, const MatrixSymOp &G, value_type etaL, const Vector *dL, const Vector *dU, const MatrixOp *E, BLAS_Cpp::Transp trans_E, const Vector *b, const Vector *eL, const Vector *eU, const MatrixOp *F, BLAS_Cpp::Transp trans_F, const Vector *f, const value_type *obj_d, const value_type *eta, const Vector *d, const Vector *nu, const Vector *mu, const Vector *Ed, const Vector *lambda, const Vector *Fd)
 Subclasses are to override this to implement the testing code. More...
 

Detailed Description

Tests the optimality conditions of the output from a QPSolverRelaxed object.

For the given QP and its solution (if solved) this class tests the optimality conditions.

The optimality conditions checked are:

 Linear dependence of gradients:
   
 (2)  d(L)/d(d) = g + G*d - nuL + nuU + op(E)'*(- muL + muU) + op(F)'*lambda
                = g + G*d + nu + op(E)'*mu + op(F)'*lambda = 0
   
      where: nu = nuU - nuL, mu = muU - muL

 Feasibility:
   
 (4.1)  etaL <=  eta
 (4.2)  dL   <=  d                       <= dU
 (4.3)  eL   <=  op(E)*d - b*eta         <= eU
 (4.4)  op(F)*d + (1 - eta) * f  = 0

 Complementarity:

 (5.1)  nu(i) * (dL - d)(i), if nu(i) <= 0, i = 1...n
 (5.2)  nu(i) * (d - dU)(i), if nu(i) >= 0, i = 1...n
 (5.3)  mu(j) * (eL - op(E)*d + b*eta)(j), if mu(j) <= 0, j = 1...m_in
 (5.4)  mu(j) * (op(E)*d - b*eta - eU)(j), if mu(j) >= 0, j = 1...m_in

The realtive error of each of these conditions is checked. Specifically, here is how the errors are computed which are compared to the error and warning tolerances:

   opt_err = || g + G*d + nu + op(E)'*mu + op(F)'*lambda ||inf / (1 + opt_scale)
                 
   feas_err = ||max(op(A)*x-b,0)||inf / ( 1 + ||op(A)*x||inf )

   comp_err(i) = gamma(i) * (op(A)*x - b)(i) / ( 1 + opt_scale + ||op(A).row(i)'*x||inf )
       ,for gamma(i) != 0

    where:
     op(A)*x <= b
     opt_scale = ||g||inf + ||G*d||inf + ||nu||inf + ||op(E)'*mu||inf + ||op(F)'*lambda||inf
                   + |(eta - etaL) * (b'*mu + f'*lambda)|

Above, op(A)*x <= b can represent any of the constraints in (4.1)-(4.4).

Any elements of opt_err(i) >= opt_warning_tol will result in an error message printed to *out. Any elements of opt_err(i) >= opt_error_tol will cause the checks to stop and false to be returned from the function check_optimiality_conditions().

Any elements of feas_err(i) >= feas_warning_tol will result in an error message printed to *out. Any elements of feas_err(i) >= feas_error_tol will cause the checks to stop and false to be returned from the function check_optimiality_conditions().

Any elements of comp_err(i) >= comp_warning_tol will result in an error message printed to *out. Any elements of comp_err(i) >= comp_error_tol will cause the checks to stop and false to be returned from the function check_optimiality_conditions().

The goal of these tests is to first and foremost to catch gross programming errors. These tests can also be used to help flag and catch illconditioning in the problem or instability in the QP solver. The importance of such tests can not be overstated. The scalings above are done to try to adjust for the scaling of the problem. Note that we are accounting for very big numbers but not for very small numbers very well and therefore tests may be conserative in some cases. At the very least we account for loss of precision due to catastrophic cancelation that occurs when subtracting large numbers and expecting to get zero. The purpose of including the term |b'*mu + f'*lambda| is to account for the situation where the relaxation is needed and kappa != 0 and therefore, form the condition d(M)/d(eta) - kappa - b'*mu - f'*lambda = 0 (with d(M)/d(eta) very large), the multipliers mu and lambda will be very large and will contribute to much roundoff errors.

As shown above, the complementarity conditions (5.1)-(5.4) are specifically checked. These should be satisfied for any solution type other than a SUBOPTIMAL_POINT return value from QPSolverRelaxed::solve_qp(). Checking the complementarity conditions for an active-set QP solver is just checking that the active constraints are satisfied. Checking them for an iterior-point solver is critical to ensure that the system was solved to satisfactory tolerance. By scaling the active-constraint violation by the Langrange multiplier we emphasis the feasibility of those constraints that have the greatest impact on the objective function. In other words, all things being equal, we are more concerned with a tight feasibility tolerance for constraints with larger lagrange multipliers than for those with smaller multipliers. The complementarity error is also scaled by the inverse of the sums of the optimality scaling opt_scale and the size of the constraint residual. By scaling by the max term opt_scale in the linear dependence of gradients we are trying to adjust the effect of the lagrange multiplier. Therefore, if the gradient of the objective g+G*d is large then opt_scale will account for this.

Definition at line 146 of file ConstrainedOptPack_QPSolverRelaxedTester.hpp.

Constructor & Destructor Documentation

ConstrainedOptPack::QPSolverRelaxedTester::QPSolverRelaxedTester ( value_type  opt_warning_tol = 1e-10,
value_type  opt_error_tol = 1e-5,
value_type  feas_warning_tol = 1e-10,
value_type  feas_error_tol = 1e-5,
value_type  comp_warning_tol = 1e-10,
value_type  comp_error_tol = 1e-5 
)
virtual ConstrainedOptPack::QPSolverRelaxedTester::~QPSolverRelaxedTester ( )
inlinevirtual

Member Function Documentation

ConstrainedOptPack::QPSolverRelaxedTester::STANDARD_MEMBER_COMPOSITION_MEMBERS ( value_type  ,
opt_warning_tol   
)

ConstrainedOptPack::QPSolverRelaxedTester::STANDARD_MEMBER_COMPOSITION_MEMBERS ( value_type  ,
opt_error_tol   
)

ConstrainedOptPack::QPSolverRelaxedTester::STANDARD_MEMBER_COMPOSITION_MEMBERS ( value_type  ,
feas_warning_tol   
)

ConstrainedOptPack::QPSolverRelaxedTester::STANDARD_MEMBER_COMPOSITION_MEMBERS ( value_type  ,
feas_error_tol   
)

ConstrainedOptPack::QPSolverRelaxedTester::STANDARD_MEMBER_COMPOSITION_MEMBERS ( value_type  ,
comp_warning_tol   
)

ConstrainedOptPack::QPSolverRelaxedTester::STANDARD_MEMBER_COMPOSITION_MEMBERS ( value_type  ,
comp_error_tol   
)

bool ConstrainedOptPack::QPSolverRelaxedTester::check_optimality_conditions ( QPSolverStats::ESolutionType  solution_type,
const value_type  infinite_bound,
std::ostream *  out,
bool  print_all_warnings,
bool  print_vectors,
const Vector &  g,
const MatrixSymOp &  G,
value_type  etaL,
const Vector &  dL,
const Vector &  dU,
const MatrixOp &  E,
BLAS_Cpp::Transp  trans_E,
const Vector &  b,
const Vector &  eL,
const Vector &  eU,
const MatrixOp &  F,
BLAS_Cpp::Transp  trans_F,
const Vector &  f,
const value_type *  obj_d,
const value_type *  eta,
const Vector *  d,
const Vector *  nu,
const Vector *  mu,
const Vector *  Ed,
const Vector *  lambda,
const Vector *  Fd 
)
virtual

Check the optimality conditions for the solved (or partially solved) QP.

The default implementation calls the function check_optimality_conditions() which accepts various sets of constraints.

Parameters
solution_type[in] Value returned from QPSolverRelaxed::solve_qp(). Even though all of the optimality conditions are checked the optimality conditions that are actually enforced is determined by this argument. OPTIMAL_SOLUTION : All of the optimality conditions are enforced. PRIMAL_FEASIBLE_POINT : Only the optimality conditions (4.1)-(4.4) are enforced. DUAL_FEASIBLE_POINT: Only the optimality condtions (2) and (6.1)-(6.4) are enforced. SUBOPTIMAL_POINT : None of the optimality conditions are enforced.
out[out] If !=NULL, the output is sent to this stream.
print_all_warnings[in] If true, then all errors greater than warning_tol will be printed.
g[in] Input to QPSolverRelaxed::solve_qp().
G[in] Input to QPSolverRelaxed::solve_qp().
etaL[in] Input to QPSolverRelaxed::solve_qp().
dL[in] Input to QPSolverRelaxed::solve_qp().
dU[in] Input to QPSolverRelaxed::solve_qp().
E[in] Input to QPSolverRelaxed::solve_qp().
trans_E[in] Input to QPSolverRelaxed::solve_qp().
b[in] Input to QPSolverRelaxed::solve_qp().
eL[in] Input to QPSolverRelaxed::solve_qp().
eU[in] Input to QPSolverRelaxed::solve_qp().
F[in] Input to QPSolverRelaxed::solve_qp().
trans_F[in] Input to QPSolverRelaxed::solve_qp().
f[in] Input to QPSolverRelaxed::solve_qp().
obj_d[in] Output from QPSolverRelaxed::solve_qp().
eta[in] Output from QPSolverRelaxed::solve_qp().
d[in] Output from QPSolverRelaxed::solve_qp().
nu[in] Output from QPSolverRelaxed::solve_qp().
mu[in] Output from QPSolverRelaxed::solve_qp().
Ed[in] Output from QPSolverRelaxed::solve_qp().
Fd[in] Output from QPSolverRelaxed::solve_qp().
Returns
true if all of the errors are greater than the error tolerances , otherwise it returns false

Definition at line 164 of file ConstrainedOptPack_QPSolverRelaxedTester.cpp.

bool ConstrainedOptPack::QPSolverRelaxedTester::check_optimality_conditions ( QPSolverStats::ESolutionType  solution_type,
const value_type  infinite_bound,
std::ostream *  out,
bool  print_all_warnings,
bool  print_vectors,
const Vector &  g,
const MatrixSymOp &  G,
value_type  etaL,
const Vector &  dL,
const Vector &  dU,
const MatrixOp &  E,
BLAS_Cpp::Transp  trans_E,
const Vector &  b,
const Vector &  eL,
const Vector &  eU,
const value_type *  obj_d,
const value_type *  eta,
const Vector *  d,
const Vector *  nu,
const Vector *  mu,
const Vector *  Ed 
)
virtual

Check the optimality conditions without general equality constrants.

Definition at line 187 of file ConstrainedOptPack_QPSolverRelaxedTester.cpp.

bool ConstrainedOptPack::QPSolverRelaxedTester::check_optimality_conditions ( QPSolverStats::ESolutionType  solution_type,
const value_type  infinite_bound,
std::ostream *  out,
bool  print_all_warnings,
bool  print_vectors,
const Vector &  g,
const MatrixSymOp &  G,
value_type  etaL,
const Vector &  dL,
const Vector &  dU,
const MatrixOp &  F,
BLAS_Cpp::Transp  trans_F,
const Vector &  f,
const value_type *  obj_d,
const value_type *  eta,
const Vector *  d,
const Vector *  nu,
const Vector *  lambda,
const Vector *  Fd 
)
virtual

Check the optimality conditions general inequality constrants.

Definition at line 208 of file ConstrainedOptPack_QPSolverRelaxedTester.cpp.

bool ConstrainedOptPack::QPSolverRelaxedTester::check_optimality_conditions ( QPSolverStats::ESolutionType  solution_type,
const value_type  infinite_bound,
std::ostream *  out,
bool  print_all_warnings,
bool  print_vectors,
const Vector &  g,
const MatrixSymOp &  G,
const Vector &  dL,
const Vector &  dU,
const value_type *  obj_d,
const Vector *  d,
const Vector *  nu 
)
virtual

Check the optimality conditions without general equality or inequality constrants (no relaxation needed).

Definition at line 228 of file ConstrainedOptPack_QPSolverRelaxedTester.cpp.

bool ConstrainedOptPack::QPSolverRelaxedTester::check_optimality_conditions ( QPSolverStats::ESolutionType  solution_type,
const value_type  infinite_bound,
std::ostream *  out,
bool  print_all_warnings,
bool  print_vectors,
const Vector &  g,
const MatrixSymOp &  G,
value_type  etaL,
const Vector *  dL,
const Vector *  dU,
const MatrixOp *  E,
BLAS_Cpp::Transp  trans_E,
const Vector *  b,
const Vector *  eL,
const Vector *  eU,
const MatrixOp *  F,
BLAS_Cpp::Transp  trans_F,
const Vector *  f,
const value_type *  obj_d,
const value_type *  eta,
const Vector *  d,
const Vector *  nu,
const Vector *  mu,
const Vector *  Ed,
const Vector *  lambda,
const Vector *  Fd 
)
virtual

This is a more flexible function where the client can set different constraints to be included.

Definition at line 245 of file ConstrainedOptPack_QPSolverRelaxedTester.cpp.

bool ConstrainedOptPack::QPSolverRelaxedTester::imp_check_optimality_conditions ( QPSolverStats::ESolutionType  solution_type,
const value_type  infinite_bound,
std::ostream *  out,
bool  print_all_warnings,
bool  print_vectors,
const Vector &  g,
const MatrixSymOp &  G,
value_type  etaL,
const Vector *  dL,
const Vector *  dU,
const MatrixOp *  E,
BLAS_Cpp::Transp  trans_E,
const Vector *  b,
const Vector *  eL,
const Vector *  eU,
const MatrixOp *  F,
BLAS_Cpp::Transp  trans_F,
const Vector *  f,
const value_type *  obj_d,
const value_type *  eta,
const Vector *  d,
const Vector *  nu,
const Vector *  mu,
const Vector *  Ed,
const Vector *  lambda,
const Vector *  Fd 
)
protectedvirtual

Subclasses are to override this to implement the testing code.

There is a default implementation that is very general and should be considered good enough for most applications.

Definition at line 276 of file ConstrainedOptPack_QPSolverRelaxedTester.cpp.


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