#include <ConstrainedOptPack_MatrixHessianSuperBasic.hpp>
Inheritance diagram for ConstrainedOptPack::MatrixHessianSuperBasic:
Provide access to constituent members | |
| const GenPermMatrixSlice & | Q_R () const |
| | |
| const GenPermMatrixSlice & | Q_X () const |
| | |
| const bnd_fixed_t & | bnd_fixed () const |
| | |
| const B_RR_ptr_t & | B_RR_ptr () const |
| | |
| const B_RX_ptr_t & | B_RX_ptr () const |
| | |
| BLAS_Cpp::Transp | B_RX_trans () const |
| | |
| const B_XX_ptr_t & | B_XX_ptr () 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 |
| | |
Public Types | |
| typedef Teuchos::RefCountPtr< const MatrixSymWithOpFactorized > | B_RR_ptr_t |
| | |
| typedef Teuchos::RefCountPtr< const MatrixOp > | B_RX_ptr_t |
| | |
| typedef Teuchos::RefCountPtr< const MatrixSymOp > | B_XX_ptr_t |
| | |
| typedef std::vector< EBounds > | bnd_fixed_t |
| | |
Public Member Functions | |
| MatrixHessianSuperBasic () | |
| Constructs to uninitialized. | |
| 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. | |
Protected Member Functions | |
| void | assert_initialized () const |
| | |
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 74 of file ConstrainedOptPack_MatrixHessianSuperBasic.hpp.
|
|
Definition at line 81 of file ConstrainedOptPack_MatrixHessianSuperBasic.hpp. |
|
|
Definition at line 84 of file ConstrainedOptPack_MatrixHessianSuperBasic.hpp. |
|
|
Definition at line 87 of file ConstrainedOptPack_MatrixHessianSuperBasic.hpp. |
|
|
Definition at line 90 of file ConstrainedOptPack_MatrixHessianSuperBasic.hpp. |
|
|
Constructs to uninitialized.
Definition at line 45 of file ConstrainedOptPack_MatrixHessianSuperBasic.cpp. |
|
||||||||||||||||||||||||||||||||||||||||
|
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}
Reimplemented in ConstrainedOptPack::MatrixHessianSuperBasicInitDiagonal. Definition at line 49 of file ConstrainedOptPack_MatrixHessianSuperBasic.cpp. |
|
|
Definition at line 255 of file ConstrainedOptPack_MatrixHessianSuperBasic.hpp. |
|
|
Definition at line 262 of file ConstrainedOptPack_MatrixHessianSuperBasic.hpp. |
|
|
Definition at line 270 of file ConstrainedOptPack_MatrixHessianSuperBasic.hpp. |
|
|
Definition at line 277 of file ConstrainedOptPack_MatrixHessianSuperBasic.hpp. |
|
|
Definition at line 285 of file ConstrainedOptPack_MatrixHessianSuperBasic.hpp. |
|
|
Definition at line 291 of file ConstrainedOptPack_MatrixHessianSuperBasic.hpp. |
|
|
Definition at line 299 of file ConstrainedOptPack_MatrixHessianSuperBasic.hpp. |
|
||||||||||||||||||||||||
|
Definition at line 156 of file ConstrainedOptPack_MatrixHessianSuperBasic.cpp. |
|
||||||||||||||||||||||||
|
Definition at line 228 of file ConstrainedOptPack_MatrixHessianSuperBasic.cpp. |
|
||||||||||||||||||||||||||||||||
|
Definition at line 300 of file ConstrainedOptPack_MatrixHessianSuperBasic.cpp. |
|
||||||||||||||||
|
Definition at line 384 of file ConstrainedOptPack_MatrixHessianSuperBasic.cpp. |
|
|
Definition at line 474 of file ConstrainedOptPack_MatrixHessianSuperBasic.cpp. |
1.3.9.1