42 #ifndef MATRIX_SYM_POS_DEF_BANDED_CHOL_H
43 #define MATRIX_SYM_POS_DEF_BANDED_CHOL_H
45 #include "ConstrainedOptPack_Types.hpp"
46 #include "AbstractLinAlgPack/src/MatrixSymWithOpFactorized.hpp"
47 #include "DenseLinAlgPack_DMatrixClass.hpp"
48 #include "Miref_count_ptr.h"
49 #include "MiReleaseResource.h"
51 namespace ConstrainedOptPack {
96 ,DMatrixSlice *
MB = NULL
99 ,DMatrixSlice *
UB = NULL
102 ,
bool update_factor =
false
103 ,value_type scale = 1.0
144 ,DMatrixSlice *
MB = NULL
147 ,DMatrixSlice *
UB = NULL
150 ,
bool update_factor =
false
151 ,value_type scale = 1.0
160 const DMatrixSlice&
MB()
const;
167 const DMatrixSlice&
UB()
const;
179 std::ostream&
output(std::ostream& out)
const;
182 ,
const DVectorSlice& vs_rhs2, value_type beta)
const;
185 ,
const SpVectorSlice& sv_rhs2, value_type beta)
const;
187 void Vp_StPtMtV(DVectorSlice* vs_lhs, value_type alpha
190 ,
const DVectorSlice& vs_rhs3, value_type beta)
const;
192 void Vp_StPtMtV(DVectorSlice* vs_lhs, value_type alpha
195 ,
const SpVectorSlice& sv_rhs3, value_type beta)
const;
202 ,
const DVectorSlice& vs_rhs2)
const;
214 mutable DMatrixSlice UB_;
217 mutable bool factor_updated_;
223 void assert_initialized()
const;
224 void update_factorization()
const;
275 #endif // MATRIX_SYM_POS_DEF_BANDED_CHOL_H
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
DMatrixSlice & UB()
Get view of UB.
BLAS_Cpp::Uplo MB_uplo() 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
DMatrixSlice & MB()
Get view of MB.
void V_InvMtV(DVectorSlice *vs_lhs, BLAS_Cpp::Transp trans_rhs1, const DVectorSlice &vs_rhs2) const
With throw exception if factorization is not allowed.
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.
BLAS_Cpp::Uplo UB_uplo() const
Matrix subclass for banded symmetric positive definite matrices and their Cholesky factors...
Teuchos::RCP< MemMngPack::ReleaseResource > release_resource_ptr_t
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.