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 | List of all members
ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd Class Reference

General (and flexible) implementation class for a QPSchur QP problem. More...

#include <ConstrainedOptPack_QPInitFixedFreeStd.hpp>

Inheritance diagram for ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd:
Inheritance graph
[legend]

Public Member Functions

 QPInitFixedFreeStd ()
 Construct uninitialized. More...
 
void initialize (const DVectorSlice &g, const MatrixSymOp &G, const MatrixOp *A, size_type n_R, const size_type i_x_free[], const size_type i_x_fixed[], const EBounds bnd_fixed[], const DVectorSlice &b_X, const MatrixSymOpNonsing &Ko, const DVectorSlice &fo, Constraints *constraints, std::ostream *out=NULL, bool test_setup=false, value_type warning_tol=1e-10, value_type error_tol=1e-5, bool print_all_warnings=false)
 Initialize. More...
 
- Public Member Functions inherited from ConstrainedOptPack::QPSchurPack::QP
virtual ~QP ()
 
virtual void dump_qp (std::ostream &out)
 Dump the definition of the QP to a stream. More...
 

Overridden from QP

size_type n () const
 
size_type m () const
 
const DVectorSlice g () const
 
const MatrixSymOp & G () const
 
const MatrixOp & A () const
 
size_type n_R () const
 
const x_init_tx_init () const
 
const l_x_X_map_tl_x_X_map () const
 
const i_x_X_map_ti_x_X_map () const
 
const DVectorSlice b_X () const
 
const GenPermMatrixSlice & Q_R () const
 
const GenPermMatrixSlice & Q_X () const
 
const MatrixSymOpNonsing & Ko () const
 
const DVectorSlice fo () const
 
Constraintsconstraints ()
 
const Constraintsconstraints () const
 

Additional Inherited Members

- Public Types inherited from ConstrainedOptPack::QPSchurPack::QP
typedef
vector_one_based_checked
< EBounds > 
x_init_t
 
typedef
vector_one_based_checked
< size_type
l_x_X_map_t
 
typedef
vector_one_based_checked
< size_type
i_x_X_map_t
 
typedef QPSchurPack::Constraints Constraints
 

Detailed Description

General (and flexible) implementation class for a QPSchur QP problem.

The basic idea of this class is to just build the QP from its various components in a way that is easy and flexible for the client. The class will also do consistency testing if asked to.

Definition at line 58 of file ConstrainedOptPack_QPInitFixedFreeStd.hpp.

Constructor & Destructor Documentation

ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::QPInitFixedFreeStd ( )

Construct uninitialized.

Definition at line 52 of file ConstrainedOptPack_QPInitFixedFreeStd.cpp.

Member Function Documentation

void ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::initialize ( const DVectorSlice &  g,
const MatrixSymOp &  G,
const MatrixOp *  A,
size_type  n_R,
const size_type  i_x_free[],
const size_type  i_x_fixed[],
const EBounds  bnd_fixed[],
const DVectorSlice &  b_X,
const MatrixSymOpNonsing &  Ko,
const DVectorSlice &  fo,
Constraints constraints,
std::ostream *  out = NULL,
bool  test_setup = false,
value_type  warning_tol = 1e-10,
value_type  error_tol = 1e-5,
bool  print_all_warnings = false 
)

Initialize.

The pointers and references to the objects pointed to by the arguments to this function must not be modified by the caller. Copies of these objects are not made internally so these objects must remain valid while this object is in use.

If the sizes of the arguments do not match up or some consistency test fails then exceptions may be thrown with (hopefully) helpful messages.

Parameters
g[in] vector (size n): objective gradient
G[in] matrix (size n x n): objective Hessian
A[in] matrix (size n x m): full rank equality constraints in Ko. If A==NULL then there are no equality constraints in Ko and m will be zero.
n_R[in] number of initially free variables
i_x_free[in] array (size n_R): i_x_free[l-1], l = 1...n_R defines the matrix Q_R as:
Q_R(:,l) = e(i_x_free[l-1]), l = 1...n_R
The ordering of these indices is significant. It is allowed for i_x_free == NULL in which case it will be considered to be identity.
i_x_fixed[in] array (size n_X = n - n_R): i_x_fixed[l-1], l = 1...n_X defines the matrix Q_X as:
Q_X(:,l) = e(i_x_fixed[l-1]), l = 1...n_X
The ordering of these indices is significant.
bnd_fixed[in] array (size n_X = n - n_R): bnd_fixed[l-1], l = 1...n_X defines the initial active set as:
{verbatim} / LOWER : b_X(l) = xL(i_x_fixed[l-1]) bnd_fixed[l-1] = | UPPER : b_X(l) = xU(i_x_fixed[l-1]) \ EQUALITY : b_X(l) = xL(i) = xU(i) (i = i_x_fixed[l-1]) {verbatim}
b_X[in] vector (size n_X = n - n_R): Initial varaible bounds (see bnd_fixed)
Ko[in] matrix (size (n_R+m) x (n_R+m)): Initial KKT matrix
fo[in] vector (size n_R + m): Initial KKT system rhs vector
constraints[in] Constraints object for the extra constraints cL_bar <= A_bar'*x <= cU_bar
out[out] If out!=NULL, then any warning or error messages will be printed here.
test_setup[in] If set to true, then consistency checks will be made on all the input arguments. The cost of the tests will not be too excessive in runtime or storge costs and do not completly validate everything
waring_tol[in] Warning tolerance for tests.
error_tol[in] Error tolerance for tests. If the relative error of any test exceeds this limit, then an error message will be printed to out (if out!=NULL) and then a runtime exception will be thrown.
print_all_warnings[in] If set to true, then any relative errors for tests that are above warning_tol will be printed to out (if out!= NULL) (O(n) output). Otherwise, if false, then only the number of violations and the maximum violation will be printed (O(1) output).

Definition at line 62 of file ConstrainedOptPack_QPInitFixedFreeStd.cpp.

size_type ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::n ( ) const
virtual
size_type ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::m ( ) const
virtual
const DVectorSlice ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::g ( ) const
virtual
const MatrixSymOp & ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::G ( ) const
virtual
const MatrixOp & ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::A ( ) const
virtual
size_type ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::n_R ( ) const
virtual
const QP::x_init_t & ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::x_init ( ) const
virtual
const QP::l_x_X_map_t & ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::l_x_X_map ( ) const
virtual
const QP::i_x_X_map_t & ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::i_x_X_map ( ) const
virtual
const DVectorSlice ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::b_X ( ) const
virtual
const GenPermMatrixSlice & ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::Q_R ( ) const
virtual
const GenPermMatrixSlice & ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::Q_X ( ) const
virtual
const MatrixSymOpNonsing & ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::Ko ( ) const
virtual
const DVectorSlice ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::fo ( ) const
virtual
Constraints & ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::constraints ( )
virtual
const Constraints & ConstrainedOptPack::QPSchurPack::QPInitFixedFreeStd::constraints ( ) const
virtual

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