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 Types | Public Member Functions | Protected Member Functions | List of all members
ConstrainedOptPack::MatrixHessianSuperBasic Class Reference

Matrix class that represents a hessian matrix where only the super submatrix for the super basic variables need be nonsingular. More...

#include <ConstrainedOptPack_MatrixHessianSuperBasic.hpp>

Inheritance diagram for ConstrainedOptPack::MatrixHessianSuperBasic:
Inheritance graph
[legend]

Public Types

typedef Teuchos::RCP< const
MatrixSymWithOpFactorized > 
B_RR_ptr_t
 
typedef Teuchos::RCP< const
MatrixOp > 
B_RX_ptr_t
 
typedef Teuchos::RCP< const
MatrixSymOp > 
B_XX_ptr_t
 
typedef std::vector< EBounds > bnd_fixed_t
 

Public Member Functions

 MatrixHessianSuperBasic ()
 Constructs to uninitialized. More...
 
virtual void initialize (size_type n, size_type n_R, const size_type i_x_free[], const size_type i_x_fixed[], const EBounds bnd_fixed[], const B_RR_ptr_t &B_RR_ptr, const B_RX_ptr_t &B_RX_ptr, BLAS_Cpp::Transp B_RX_trans, const B_XX_ptr_t &B_XX_ptr)
 Initialize the matrix. More...
 

Protected Member Functions

void assert_initialized () const
 

Provide access to constituent members

const GenPermMatrixSlice & Q_R () const
 
const GenPermMatrixSlice & Q_X () const
 
const bnd_fixed_tbnd_fixed () const
 
const B_RR_ptr_tB_RR_ptr () const
 
const B_RX_ptr_tB_RX_ptr () const
 
BLAS_Cpp::Transp B_RX_trans () const
 
const B_XX_ptr_tB_XX_ptr () const
 

Overridden from Matrix

size_type rows () const
 

Overridden from MatrixOp

void Vp_StMtV (DVectorSlice *vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const DVectorSlice &vs_rhs2, value_type beta) const
 
void Vp_StMtV (DVectorSlice *vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice &sv_rhs2, value_type beta) const
 
void Vp_StPtMtV (DVectorSlice *vs_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, BLAS_Cpp::Transp M_rhs2_trans, const DVectorSlice &sv_rhs3, value_type beta) const
 
value_type transVtMtV (const SpVectorSlice &sv_rhs1, BLAS_Cpp::Transp trans_rhs2, const SpVectorSlice &sv_rhs3) const
 

Detailed Description

Matrix class that represents a hessian matrix where only the super submatrix for the super basic variables need be nonsingular.

Given a Hessian matrix #B# and a partitioning of the variables #Q = [ Q_R Q_X ]# into free (superbasic) Q_R'*x# and fixed (nonbasic) Q_X'*x# variables, this class represents the following matrix: {verbatim} [n,n] == size(B) [n,n] == size(Q), Q is an othogonal permutation matrix (i.e. Q*Q' = Q'*Q = I) [n,n_R] == size(Q_R) [n,n_X] == size(Q_X)

B = Q*Q'*B*Q*Q' = [ Q_R Q_X ] * [ Q_R'*B*Q_R Q_R'*B*Q_X ] * [ Q_R' ] [ Q_X'*B*Q_R Q_X'*B*Q_X ] [ Q_X' ]

= [ Q_R Q_X ] * [ B_RR op(B_RX) ] * [ Q_R' ] [ op(B_RX') B_XX ] [ Q_X' ]

= Q_R*B_RR*Q_R' + Q_R*op(B_RX)*Q_X' + Q_X*op(B_RX')*Q_R + Q_X*B_XX*Q_X' {verbatim} Above, we allow the prepresentation of #op(B_RX) = Q_R'*B*Q_X# to be transposed to allow for more flexibility. Since #B_RR# and #B_XX# are symmetric, we do not need to worry about transpose or not. For this class the matrix #B_RR# is required to be symmetric and nonsingular (#MatrixSymWithOpFactorized# interface), but not necessarily positive definite. This is the condition necessary for the Hessian when projected into the active constraints at the solution for an NLP. The other matrices hold no other special properties other than #B_XX# being symmetric of course.

The default compiler generated constructors are allowed and initialize the matrix to uninitialized by default.

Definition at line 87 of file ConstrainedOptPack_MatrixHessianSuperBasic.hpp.

Member Typedef Documentation

Constructor & Destructor Documentation

ConstrainedOptPack::MatrixHessianSuperBasic::MatrixHessianSuperBasic ( )

Constructs to uninitialized.

Member Function Documentation

virtual void ConstrainedOptPack::MatrixHessianSuperBasic::initialize ( size_type  n,
size_type  n_R,
const size_type  i_x_free[],
const size_type  i_x_fixed[],
const EBounds  bnd_fixed[],
const B_RR_ptr_t B_RR_ptr,
const B_RX_ptr_t B_RX_ptr,
BLAS_Cpp::Transp  B_RX_trans,
const B_XX_ptr_t B_XX_ptr 
)
virtual

Initialize the matrix.

Preconditions:{itemize} [#i_x_free != NULL#] #0 < i_x_free[l-1] <= n, l = 1...n_R# (throw ???) [#i_x_fixed != NULL#]#0 < i_x_fixed[l-1] <= n, l = 1...n-n_R# (throw ???) [#i_x_free != NULL && i_x_fixed != NULL#] #i_x_free[l-1] != i_x_fixed[p-1], l = 1...n_R, p = 1...n-n_X# (throw ???) [#n_R > 0#] #B_RR_ptr.get() != NULL && B_RR_ptr->rows() == n_R && B_RR_ptr->cols() == n_R# (throw #std::invalid_argument#) [#n_R == 0#] #B_RR_ptr.get() == NULL# (throw #std::invalid_argument#) [#n_R < n && B_RX_ptr.get() != NULL && B_RX_trans == no_trans#] B_RX_ptr->rows() == n_R && B_RX_ptr->cols() == n-n_R# (throw #std::invalid_argument#) [#n_R < n && B_RX_ptr.get() != NULL && B_RX_trans == trans#] B_RX_ptr->rows() == n-n_R && B_RX_ptr->cols() == n_R# (throw #std::invalid_argument#) [#n_R == n#] #B_RX_ptr.get() == NULL# (throw ##std::invalid_argument#) [#n_R < n#] #B_XX_ptr.get() != NULL && B_XX_ptr->rows() == n-n_R && B_XX_ptr->cols() == n-n_R# (throw #std::invalid_argument#) [#n_R == n#] #B_XX_ptr.get() == NULL# (throw ##std::invalid_argument#) {itemize}

Parameters
n[in] number of variables (used for consistency checking)
n_R[in] number of initially free variables (used for consistency checking)
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. If #n == n_R# then #i_x_free == NULL# is allowed in which case it is assumed to be identity. If #n_R == 0# then of course #i_x_free == NULL# is allowed.
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. If #n_R==0# then #i_x_fixed == NULL# is allowed in which case it is assumed to be identity. If #n_R == n# then of course #i_x_fixed == NULL# is allowed.
bnd_fixed[in] array (size #n_X#): bnd_fixed[l-1], l = 1...n_X# defines what bound the variables in #i_x_fixed[l-1], l = 1...n_X# are fixed at: LOWER#, UPPER# or #EQUALITY#. If #n_R == n# then of course bnd_fixed == NULL# is allowed.
B_RR_ptr[in] Smart pointer to matrix #B_RR# (size #n_R x n_R#) for the free (super basic) variables. #B_RR_ptr.get() != NULL# must be true if #n_R > # or an exception will be thrown. if #n_R == 0# then #B_RR_ptr.get() == NULL# may be true.
B_RX_ptr[in] Smart pointer to matrix #B_RX# (size #n_R x n_X# if B_RX_trans==no_trans# or #n_X x n_R# if B_RX_trans==trans#) for the cross terms of free (super basic) and fixed (nonbasic) variables. It is allowed for #B_RX_ptr.get() == NULL#.
B_RX_trans[in] Determines if op(B_RX) = B_RX (no_trans#) or op(B_RX) = B_RX' (trans#). Ignored if #n_R == n#.
B_XX_ptr[in] Smart pointer to matrix B_XX (size #n_X x n_X#) for the fixed (nonbasic) variables. #B_XX_ptr.get() != NULL# must be true if #n_R < n# or an exception will be thrown. If #n_R == n# then #B_XX_ptr.get() == NULL# may be true.

Reimplemented in ConstrainedOptPack::MatrixHessianSuperBasicInitDiagonal.

const GenPermMatrixSlice & ConstrainedOptPack::MatrixHessianSuperBasic::Q_R ( ) const
inline
const GenPermMatrixSlice & ConstrainedOptPack::MatrixHessianSuperBasic::Q_X ( ) const
inline
const MatrixHessianSuperBasic::bnd_fixed_t & ConstrainedOptPack::MatrixHessianSuperBasic::bnd_fixed ( ) const
inline
const MatrixHessianSuperBasic::B_RR_ptr_t & ConstrainedOptPack::MatrixHessianSuperBasic::B_RR_ptr ( ) const
inline
const MatrixHessianSuperBasic::B_RX_ptr_t & ConstrainedOptPack::MatrixHessianSuperBasic::B_RX_ptr ( ) const
inline
BLAS_Cpp::Transp ConstrainedOptPack::MatrixHessianSuperBasic::B_RX_trans ( ) const
inline
const MatrixHessianSuperBasic::B_XX_ptr_t & ConstrainedOptPack::MatrixHessianSuperBasic::B_XX_ptr ( ) const
inline
void ConstrainedOptPack::MatrixHessianSuperBasic::Vp_StMtV ( DVectorSlice *  vs_lhs,
value_type  alpha,
BLAS_Cpp::Transp  trans_rhs1,
const DVectorSlice &  vs_rhs2,
value_type  beta 
) const

void ConstrainedOptPack::MatrixHessianSuperBasic::Vp_StMtV ( DVectorSlice *  vs_lhs,
value_type  alpha,
BLAS_Cpp::Transp  trans_rhs1,
const SpVectorSlice &  sv_rhs2,
value_type  beta 
) const

void ConstrainedOptPack::MatrixHessianSuperBasic::Vp_StPtMtV ( DVectorSlice *  vs_lhs,
value_type  alpha,
const GenPermMatrixSlice &  P_rhs1,
BLAS_Cpp::Transp  P_rhs1_trans,
BLAS_Cpp::Transp  M_rhs2_trans,
const DVectorSlice &  sv_rhs3,
value_type  beta 
) const

value_type ConstrainedOptPack::MatrixHessianSuperBasic::transVtMtV ( const SpVectorSlice &  sv_rhs1,
BLAS_Cpp::Transp  trans_rhs2,
const SpVectorSlice &  sv_rhs3 
) const

void ConstrainedOptPack::MatrixHessianSuperBasic::assert_initialized ( ) const
protected


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