MOOCHO (Single Doxygen Collection)
Version of the Day
|
Matrix subclass for banded symmetric positive definite matrices and their Cholesky factors. More...
#include <ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp>
Inherits MatrixSymWithOpFactorized.
Public Types | |
typedef Teuchos::RCP < MemMngPack::ReleaseResource > | release_resource_ptr_t |
Public Member Functions | |
MatrixSymPosDefBandedChol (size_type n=0, size_type kd=0, DMatrixSlice *MB=NULL, const release_resource_ptr_t &MB_release_resource_ptr=NULL, BLAS_Cpp::Uplo MB_uplo=BLAS_Cpp::lower, DMatrixSlice *UB=NULL, const release_resource_ptr_t &UB_release_resource_ptr=NULL, BLAS_Cpp::Uplo UB_uplo=BLAS_Cpp::lower, bool update_factor=false, value_type scale=1.0) | |
Construct and Initialize. More... | |
void | initialize (size_type n=0, size_type kd=0, DMatrixSlice *MB=NULL, const release_resource_ptr_t &MB_release_resource_ptr=NULL, BLAS_Cpp::Uplo MB_uplo=BLAS_Cpp::lower, DMatrixSlice *UB=NULL, const release_resource_ptr_t &UB_release_resource_ptr=NULL, BLAS_Cpp::Uplo UB_uplo=BLAS_Cpp::lower, bool update_factor=false, value_type scale=1.0) |
Initialize. More... | |
size_type | kd () const |
DMatrixSlice & | MB () |
Get view of MB. More... | |
const DMatrixSlice & | MB () const |
BLAS_Cpp::Uplo | MB_uplo () const |
DMatrixSlice & | UB () |
Get view of UB. More... | |
const DMatrixSlice & | UB () const |
BLAS_Cpp::Uplo | UB_uplo () const |
size_type | rows () const |
size_type | nz () const |
std::ostream & | output (std::ostream &out) const |
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 &vs_rhs3, 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 SpVectorSlice &sv_rhs3, value_type beta) const |
void | V_InvMtV (DVectorSlice *vs_lhs, BLAS_Cpp::Transp trans_rhs1, const DVectorSlice &vs_rhs2) const |
With throw exception if factorization is not allowed. More... | |
Private Member Functions | |
void | assert_initialized () const |
void | update_factorization () const |
Private Attributes | |
size_type | n_ |
size_type | kd_ |
DMatrixSlice | MB_ |
release_resource_ptr_t | MB_release_resource_ptr_ |
BLAS_Cpp::Uplo | MB_uplo_ |
DMatrixSlice | UB_ |
release_resource_ptr_t | UB_release_resource_ptr_ |
BLAS_Cpp::Uplo | UB_uplo_ |
bool | factor_updated_ |
value_type | scale_ |
Matrix subclass for banded symmetric positive definite matrices and their Cholesky factors.
This class is designed to support the LAPACK routines for banded symmetric positive definite matrices. The banded matrix and/or its cholesky factor are stored in simple flat rectangular matrices compatible with the LAPACK routines.
For example, for n = 8, kd = 3# the original matrix M# (if set) is stored in the following format MB#: {verbatim}
M MB
[ x x x x ] [ x x x x x x x x ] [ x x x x x ] lower triangle [ x x x x x x x o ] [ x x x x x x ] => [ x x x x x x o o ] [ x x x x x x x ] [ x x x x x o o o ] [ x x x x x x x ] [ x x x x x x ] [ o o o x x x x x ] [ x x x x x ] upper triangle [ o o x x x x x x ] [ x x x x ] => [ o x x x x x x x ] [ x x x x x x x x ] {verbatim} The Cholesky factor #U# is sorted in a similar format UB#. Technically, the matrix is M = scale * U'*U# so that M# may be negative definite as well.
Definition at line 78 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
typedef Teuchos::RCP< MemMngPack::ReleaseResource> ConstrainedOptPack::MatrixSymPosDefBandedChol::release_resource_ptr_t |
Definition at line 84 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
ConstrainedOptPack::MatrixSymPosDefBandedChol::MatrixSymPosDefBandedChol | ( | size_type | n = 0 , |
size_type | kd = 0 , |
||
DMatrixSlice * | MB = NULL , |
||
const release_resource_ptr_t & | MB_release_resource_ptr = NULL , |
||
BLAS_Cpp::Uplo | MB_uplo = BLAS_Cpp::lower , |
||
DMatrixSlice * | UB = NULL , |
||
const release_resource_ptr_t & | UB_release_resource_ptr = NULL , |
||
BLAS_Cpp::Uplo | UB_uplo = BLAS_Cpp::lower , |
||
bool | update_factor = false , |
||
value_type | scale = 1.0 |
||
) |
Construct and Initialize.
This constructor just calls #this->initialize(...)#.
void ConstrainedOptPack::MatrixSymPosDefBandedChol::initialize | ( | size_type | n = 0 , |
size_type | kd = 0 , |
||
DMatrixSlice * | MB = NULL , |
||
const release_resource_ptr_t & | MB_release_resource_ptr = NULL , |
||
BLAS_Cpp::Uplo | MB_uplo = BLAS_Cpp::lower , |
||
DMatrixSlice * | UB = NULL , |
||
const release_resource_ptr_t & | UB_release_resource_ptr = NULL , |
||
BLAS_Cpp::Uplo | UB_uplo = BLAS_Cpp::lower , |
||
bool | update_factor = false , |
||
value_type | scale = 1.0 |
||
) |
Initialize.
If called with all of the default arguments then #this# will become uninitialized.
ToDo: Finish pre and post conditions!
n | [in] Determines the size of the banded matrix (n x n). If n == 0# then all of the following arguments should be left at their defaults and #this# will become uninitialized. |
kd | [in] Determines the band width of the matrix as defined by xPBTRF(...). |
MB | [in/state] If MB != NULL# then this matrix (size (kd+1) x n) is used to store the original banded matrix MB# in the format of xPBTRF(...). This matrix must be initialized on input. |
MB_release_resource_ptr | [in] Only significant if MB != NULL#. Points to a resource to be released when MB# is no longer needed. |
MB_uplo | [in] Determines if MB# is stores the upper or lower triangular elements. |
UB | [in/state] If UB != NULL# then this matrix (size (kd+1) x n) is used to store the Cholesky factor of the banded matrix in the format of xPBTRF(...). This matrix may or may not be initialized on input. If #update_factor == false# this this matrix must already be initialized. If #update_factor == true# then this matrix will be computed. 2 * If UB == NULL# then storage for the Cholesky factor will be computed on the fly and will be factored. |
UB_release_resource_ptr | [in] Only significant if UB != NULL#. Points to a resource to be released when UB# is no longer needed. |
UB_uplo | [in] Determines if UB# is stores the upper or lower triangular elements. |
update_factor | [in] If true then the factor will be computed within this function call. |
scale | [in] Only significant if MB != NULL# or UB != NULL# (see intro). |
|
inline |
Definition at line 232 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
|
inline |
Get view of MB.
Definition at line 238 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
|
inline |
Definition at line 244 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
|
inline |
Definition at line 250 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
|
inline |
Get view of UB.
Definition at line 256 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
|
inline |
Definition at line 262 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
|
inline |
Definition at line 268 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
size_type ConstrainedOptPack::MatrixSymPosDefBandedChol::rows | ( | ) | const |
size_type ConstrainedOptPack::MatrixSymPosDefBandedChol::nz | ( | ) | const |
std::ostream& ConstrainedOptPack::MatrixSymPosDefBandedChol::output | ( | std::ostream & | out | ) | const |
void ConstrainedOptPack::MatrixSymPosDefBandedChol::Vp_StMtV | ( | DVectorSlice * | vs_lhs, |
value_type | alpha, | ||
BLAS_Cpp::Transp | trans_rhs1, | ||
const DVectorSlice & | vs_rhs2, | ||
value_type | beta | ||
) | const |
void ConstrainedOptPack::MatrixSymPosDefBandedChol::Vp_StMtV | ( | DVectorSlice * | vs_lhs, |
value_type | alpha, | ||
BLAS_Cpp::Transp | trans_rhs1, | ||
const SpVectorSlice & | sv_rhs2, | ||
value_type | beta | ||
) | const |
void ConstrainedOptPack::MatrixSymPosDefBandedChol::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 & | vs_rhs3, | ||
value_type | beta | ||
) | const |
void ConstrainedOptPack::MatrixSymPosDefBandedChol::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 SpVectorSlice & | sv_rhs3, | ||
value_type | beta | ||
) | const |
void ConstrainedOptPack::MatrixSymPosDefBandedChol::V_InvMtV | ( | DVectorSlice * | vs_lhs, |
BLAS_Cpp::Transp | trans_rhs1, | ||
const DVectorSlice & | vs_rhs2 | ||
) | const |
With throw exception if factorization is not allowed.
|
private |
|
private |
|
private |
Definition at line 209 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
|
private |
Definition at line 210 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
|
private |
Definition at line 211 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
|
private |
Definition at line 212 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
|
private |
Definition at line 213 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
|
mutableprivate |
Definition at line 214 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
|
mutableprivate |
Definition at line 215 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
|
mutableprivate |
Definition at line 216 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
|
mutableprivate |
Definition at line 217 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.
|
private |
Definition at line 218 of file ConstrainedOptPack_MatrixSymPosDefBandedChol.hpp.