ConstrainedOptPack: C++ Tools for Constrained (and Unconstrained) Optimization
Version of the Day
|
Matrix class for non-singular Hessian matrix augmented with a terms for "Big M" relaxation variables. More...
#include <ConstrainedOptPack_MatrixSymHessianRelaxNonSing.hpp>
Public Types | |
typedef Teuchos::RCP< const MatrixSymOpNonsing > | G_ptr_t |
typedef Teuchos::RCP < VectorMutable > | vec_mut_ptr_t |
typedef Teuchos::RCP< const VectorSpace > | space_ptr_t |
Constructors/initializers | |
MatrixSymHessianRelaxNonSing () | |
Construct uninitialized. More... | |
MatrixSymHessianRelaxNonSing (const G_ptr_t &G_ptr, const vec_mut_ptr_t &M_diag_ptr, const space_ptr_t &space=Teuchos::null) | |
Constructor (calls initialize() ). More... | |
void | initialize (const G_ptr_t &G_ptr, const vec_mut_ptr_t &M_diag_ptr, const space_ptr_t &space=Teuchos::null) |
Initialize the Hessian and the relaxation terms. More... | |
void | set_uninitialized () |
Set uninitalized. More... | |
const G_ptr_t & | G_ptr () const |
const vec_mut_ptr_t & | M_diag_ptr () const |
const MatrixSymOpNonsing & | G () const |
const AbstractLinAlgPack::MatrixSymDiagStd & | M () const |
Overridden from MatrixOp | |
const VectorSpace & | space_cols () const |
bool | Mp_StM (MatrixOp *mwo_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs) const |
void | Vp_StMtV (VectorMutable *v_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2, value_type beta) const |
void | Vp_StMtV (VectorMutable *v_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice &sv_rhs2, value_type beta) const |
void | Vp_StPtMtV (VectorMutable *v_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, BLAS_Cpp::Transp M_rhs2_trans, const Vector &v_rhs3, value_type beta) const |
void | Vp_StPtMtV (VectorMutable *v_lhs, value_type alpha, const GenPermMatrixSlice &P_rhs1, BLAS_Cpp::Transp P_rhs1_trans, BLAS_Cpp::Transp M_rhs2_trans, const SpVectorSlice &sv_rhs3, value_type beta) const |
Overridden form MatrixSymOp | |
void | Mp_StPtMtP (MatrixSymOp *sym_lhs, value_type alpha, EMatRhsPlaceHolder dummy_place_holder, const GenPermMatrixSlice &gpms_rhs, BLAS_Cpp::Transp gpms_rhs_trans, value_type beta) const |
Overridden from MatrixNonsing | |
void | V_InvMtV (VectorMutable *v_lhs, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2) const |
void | V_InvMtV (VectorMutable *v_lhs, BLAS_Cpp::Transp trans_rhs1, const SpVectorSlice &sv_rhs2) const |
Matrix class for non-singular Hessian matrix augmented with a terms for "Big M" relaxation variables.
The matrix that is formed is: {verbatim}
H = [ G ] [ M ] {verbatim} where M
is a diagonal matrix made up of entries M_diag
.
Definition at line 62 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.hpp.
typedef Teuchos::RCP<const MatrixSymOpNonsing> ConstrainedOptPack::MatrixSymHessianRelaxNonSing::G_ptr_t |
Definition at line 68 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.hpp.
typedef Teuchos::RCP<VectorMutable> ConstrainedOptPack::MatrixSymHessianRelaxNonSing::vec_mut_ptr_t |
Definition at line 70 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.hpp.
typedef Teuchos::RCP<const VectorSpace> ConstrainedOptPack::MatrixSymHessianRelaxNonSing::space_ptr_t |
Definition at line 72 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.hpp.
ConstrainedOptPack::MatrixSymHessianRelaxNonSing::MatrixSymHessianRelaxNonSing | ( | ) |
Construct uninitialized.
Postconditions:
this->rows() == 0
&this->G_ptr().get() == NULL
&this->M_diag_ptr().get() == NULL
Definition at line 149 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.cpp.
ConstrainedOptPack::MatrixSymHessianRelaxNonSing::MatrixSymHessianRelaxNonSing | ( | const G_ptr_t & | G_ptr, |
const vec_mut_ptr_t & | M_diag_ptr, | ||
const space_ptr_t & | space = Teuchos::null |
||
) |
Constructor (calls initialize()
).
Definition at line 153 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.cpp.
void ConstrainedOptPack::MatrixSymHessianRelaxNonSing::initialize | ( | const G_ptr_t & | G_ptr, |
const vec_mut_ptr_t & | M_diag_ptr, | ||
const space_ptr_t & | space = Teuchos::null |
||
) |
Initialize the Hessian and the relaxation terms.
Preconditions:
G_ptr.get() != NULL
(throw std::invalid_argument
) G_ptr->rows() > 0
(throw std::invalid_argument
) M_diag_ptr.get() != NULL
(throw std::invalid_argument
) M_diag_ptr->dim() > 0
(throw std::invalid_argument
) Postconditions:
this->rows() == G_ptr->rows() + M_diag->dim()
&this->G() == G_ptr.get()
&this->M().diag() == M_diag_ptr.get()
G_ptr | [in] Smart pointer to matrix that this object will represent. The underlying matrix object *G_ptr.get() should not be modified without calling this->initialize() again. |
M_diag_ptr | [in] Smart pointer to the diagonal for M . All of the elements in this vector must be nonzero! This vector is used to initalize the diagonal matrix M and this vector should not be modified externally without calling this->initialize() again. |
space | [in] Smart pointer to the space used for this->space_cols() and this->space_rows() . If space.get() == NULL then VectorSpaceCompoisteStd will be used which is constructed from the spaces G_ptr->space_cols() and M_diag_ptr->space() . |
Definition at line 163 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.cpp.
void ConstrainedOptPack::MatrixSymHessianRelaxNonSing::set_uninitialized | ( | ) |
Set uninitalized.
|
inline |
Definition at line 229 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.hpp.
|
inline |
Definition at line 236 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.hpp.
|
inline |
Definition at line 243 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.hpp.
|
inline |
Definition at line 251 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.hpp.
|
virtual |
Implements AbstractLinAlgPack::MatrixBase.
Definition at line 196 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.cpp.
bool ConstrainedOptPack::MatrixSymHessianRelaxNonSing::Mp_StM | ( | MatrixOp * | mwo_lhs, |
value_type | alpha, | ||
BLAS_Cpp::Transp | trans_rhs | ||
) | const |
Definition at line 202 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.cpp.
void ConstrainedOptPack::MatrixSymHessianRelaxNonSing::Vp_StMtV | ( | VectorMutable * | v_lhs, |
value_type | alpha, | ||
BLAS_Cpp::Transp | trans_rhs1, | ||
const Vector & | v_rhs2, | ||
value_type | beta | ||
) | const |
Definition at line 220 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.cpp.
void ConstrainedOptPack::MatrixSymHessianRelaxNonSing::Vp_StMtV | ( | VectorMutable * | v_lhs, |
value_type | alpha, | ||
BLAS_Cpp::Transp | trans_rhs1, | ||
const SpVectorSlice & | sv_rhs2, | ||
value_type | beta | ||
) | const |
Definition at line 237 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.cpp.
void ConstrainedOptPack::MatrixSymHessianRelaxNonSing::Vp_StPtMtV | ( | VectorMutable * | v_lhs, |
value_type | alpha, | ||
const GenPermMatrixSlice & | P_rhs1, | ||
BLAS_Cpp::Transp | P_rhs1_trans, | ||
BLAS_Cpp::Transp | M_rhs2_trans, | ||
const Vector & | v_rhs3, | ||
value_type | beta | ||
) | const |
Definition at line 254 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.cpp.
void ConstrainedOptPack::MatrixSymHessianRelaxNonSing::Vp_StPtMtV | ( | VectorMutable * | v_lhs, |
value_type | alpha, | ||
const GenPermMatrixSlice & | P_rhs1, | ||
BLAS_Cpp::Transp | P_rhs1_trans, | ||
BLAS_Cpp::Transp | M_rhs2_trans, | ||
const SpVectorSlice & | sv_rhs3, | ||
value_type | beta | ||
) | const |
Definition at line 269 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.cpp.
void ConstrainedOptPack::MatrixSymHessianRelaxNonSing::Mp_StPtMtP | ( | MatrixSymOp * | sym_lhs, |
value_type | alpha, | ||
EMatRhsPlaceHolder | dummy_place_holder, | ||
const GenPermMatrixSlice & | gpms_rhs, | ||
BLAS_Cpp::Transp | gpms_rhs_trans, | ||
value_type | beta | ||
) | const |
Definition at line 285 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.cpp.
void ConstrainedOptPack::MatrixSymHessianRelaxNonSing::V_InvMtV | ( | VectorMutable * | v_lhs, |
BLAS_Cpp::Transp | trans_rhs1, | ||
const Vector & | v_rhs2 | ||
) | const |
Definition at line 372 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.cpp.
void ConstrainedOptPack::MatrixSymHessianRelaxNonSing::V_InvMtV | ( | VectorMutable * | v_lhs, |
BLAS_Cpp::Transp | trans_rhs1, | ||
const SpVectorSlice & | sv_rhs2 | ||
) | const |
Definition at line 387 of file ConstrainedOptPack_MatrixSymHessianRelaxNonSing.cpp.