ConstrainedOptPack: C++ Tools for Constrained (and Unconstrained) Optimization
Version of the Day
|
Matrix subclass for variable reduction orthogonal matrix R = Gc(:,con_decomp)'*Y. More...
#include <ConstrainedOptPack_MatrixDecompRangeOrthog.hpp>
Constructors/initializers | |
MatrixDecompRangeOrthog () | |
Constructs uninitialized. More... | |
MatrixDecompRangeOrthog (const C_ptr_t &C_ptr, const D_ptr_t &D_ptr, const S_ptr_t &S_ptr) | |
Calls this->initialize() . More... | |
void | initialize (const C_ptr_t &C_ptr, const D_ptr_t &D_ptr, const S_ptr_t &S_ptr) |
Initialize the matrix object. More... | |
void | set_uninitialized () |
Make uninitialized. More... | |
Access | |
const C_ptr_t & | C_ptr () const |
const D_ptr_t & | D_ptr () const |
const S_ptr_t & | S_ptr () const |
Overridden from MatrixOp | |
size_type | rows () const |
size_type | cols () const |
const VectorSpace & | space_cols () const |
const VectorSpace & | space_rows () const |
std::ostream & | output (std::ostream &out) const |
void | Vp_StMtV (VectorMutable *v_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2, value_type beta) const |
Overridden from MatrixOpNonsing | |
void | V_InvMtV (VectorMutable *v_lhs, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2) const |
Matrix subclass for variable reduction orthogonal matrix R = Gc(:,con_decomp)'*Y.
This matrix class is used to represent the matrix:
R = C*(I + D*D') inv(R) = inv(I + D*D') * inv(C) = (I - D * inv(I + D'*D) * D') * inv(C) \______/ S
Above, the expresion for inv(R)
is derived using the Sherman-Morrison-Woodbury formula. The nonsingular matrix S = I + D'*D
is setup by the client, along with the basis matrix C
and the direct sensitivity matrix D
.
Definition at line 67 of file ConstrainedOptPack_MatrixDecompRangeOrthog.hpp.
ConstrainedOptPack::MatrixDecompRangeOrthog::MatrixDecompRangeOrthog | ( | ) |
Constructs uninitialized.
Postconditions:
this->set_uninitialized()
. Definition at line 56 of file ConstrainedOptPack_MatrixDecompRangeOrthog.cpp.
ConstrainedOptPack::MatrixDecompRangeOrthog::MatrixDecompRangeOrthog | ( | const C_ptr_t & | C_ptr, |
const D_ptr_t & | D_ptr, | ||
const S_ptr_t & | S_ptr | ||
) |
Calls this->initialize()
.
Definition at line 59 of file ConstrainedOptPack_MatrixDecompRangeOrthog.cpp.
void ConstrainedOptPack::MatrixDecompRangeOrthog::initialize | ( | const C_ptr_t & | C_ptr, |
const D_ptr_t & | D_ptr, | ||
const S_ptr_t & | S_ptr | ||
) |
Initialize the matrix object.
C_ptr | [in] |
D_ptr | [in] |
S_ptr | [in] |
Preconditions:
C_ptr.get() != NULL
(throw std::invalid_argument
) D_ptr.get() != NULL
(throw std::invalid_argument
) S_ptr.get() != NULL
(throw std::invalid_argument
) C_ptr->space_rows().is_compatible(D_ptr->space_cols()) == true
(throw VectorSpace::IncompatibleVectorSpaces
) S_ptr->space_cols().is_compatible(D_ptr->space_rows()) == true
(throw VectorSpace::IncompatibleVectorSpaces
) Postconditions:
this->C_ptr().get() == C_ptr.get()
this->D_ptr().get() == D_ptr.get()
this->S_ptr().get() == S_ptr.get()
this->space_cols().is_compatible(C_ptr->space_cols())
this->space_rows().is_compatible(C_ptr->space_rows())
Definition at line 68 of file ConstrainedOptPack_MatrixDecompRangeOrthog.cpp.
void ConstrainedOptPack::MatrixDecompRangeOrthog::set_uninitialized | ( | ) |
|
inline |
Definition at line 214 of file ConstrainedOptPack_MatrixDecompRangeOrthog.hpp.
|
inline |
Definition at line 221 of file ConstrainedOptPack_MatrixDecompRangeOrthog.hpp.
|
inline |
Definition at line 228 of file ConstrainedOptPack_MatrixDecompRangeOrthog.hpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixBase.
Definition at line 109 of file ConstrainedOptPack_MatrixDecompRangeOrthog.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixBase.
Definition at line 114 of file ConstrainedOptPack_MatrixDecompRangeOrthog.cpp.
|
virtual |
Implements AbstractLinAlgPack::MatrixBase.
Definition at line 119 of file ConstrainedOptPack_MatrixDecompRangeOrthog.cpp.
|
virtual |
Implements AbstractLinAlgPack::MatrixBase.
Definition at line 124 of file ConstrainedOptPack_MatrixDecompRangeOrthog.cpp.
|
virtual |
Reimplemented from AbstractLinAlgPack::MatrixOp.
Definition at line 129 of file ConstrainedOptPack_MatrixDecompRangeOrthog.cpp.
void ConstrainedOptPack::MatrixDecompRangeOrthog::Vp_StMtV | ( | VectorMutable * | v_lhs, |
value_type | alpha, | ||
BLAS_Cpp::Transp | trans_rhs1, | ||
const Vector & | v_rhs2, | ||
value_type | beta | ||
) | const |
Definition at line 143 of file ConstrainedOptPack_MatrixDecompRangeOrthog.cpp.
void ConstrainedOptPack::MatrixDecompRangeOrthog::V_InvMtV | ( | VectorMutable * | v_lhs, |
BLAS_Cpp::Transp | trans_rhs1, | ||
const Vector & | v_rhs2 | ||
) | const |
Definition at line 205 of file ConstrainedOptPack_MatrixDecompRangeOrthog.cpp.