MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Typedefs | Functions
DenseLinAlgLAPack Namespace Reference

Classes

class  FactorizationException
 Exception for factorization error. More...
 

Typedefs

typedef DenseLinAlgPack::value_type value_type
 
typedef DenseLinAlgPack::index_type index_type
 
typedef DenseLinAlgPack::size_type size_type
 

Functions

void potrf (DMatrixSliceTriEle *A)
 Calls xPOTRF to compute the cholesky factorization of a symmetric positive definte matrix. More...
 
void geqrf (DMatrixSlice *A, DVectorSlice *tau, DVectorSlice *work)
 Calls xGEQRF to compute the QR factorization of a matrix A. More...
 
void ormrq (BLAS_Cpp::Side side, BLAS_Cpp::Transp trans, const DMatrixSlice &A, const DVectorSlice &tau, DMatrixSlice *C, DVectorSlice *work)
 Calls xORMRQ to compute a matrix matrix op(C)=op(Q)*op(C) where Q is stored in A and tau computed from geqrf(..). More...
 
void sytrf (DMatrixSliceTriEle *A, FortranTypes::f_int ipiv[], DVectorSlice *work)
 Calls xSYTRF to compute the P*A*P' = L'*D*L factorization of a symmetric possibly indefinite matrix. More...
 
void sytrs (const DMatrixSliceTriEle &A, FortranTypes::f_int ipiv[], DMatrixSlice *B, DVectorSlice *work)
 Calls xSYTRS(...) to compute the solution of the factorized system A * X = B where A was factorized by xSYTRF(...). More...
 
void getrf (DMatrixSlice *A, FortranTypes::f_int ipiv[], FortranTypes::f_int *rank)
 Calls xGETRF to compute the P'*A = L*U factorization of a general rectuangular matrix. More...
 
void getrs (const DMatrixSlice &LU, const FortranTypes::f_int ipiv[], BLAS_Cpp::Transp transp, DMatrixSlice *B)
 Calls xGETRS to solve linear systems with the factors of P'*A = L*U generated by xGETRF. More...
 

Typedef Documentation

Definition at line 16 of file DenseLinAlgLAPack_Types.hpp.

Definition at line 17 of file DenseLinAlgLAPack_Types.hpp.

Definition at line 18 of file DenseLinAlgLAPack_Types.hpp.

Function Documentation

void DenseLinAlgLAPack::potrf ( DMatrixSliceTriEle *  A)

Calls xPOTRF to compute the cholesky factorization of a symmetric positive definte matrix.

On input A contains the upper or lower trianagular elements of a symmetric positive definite matrix. On output it contians the cholesky factor U (A = U'*U) if #A.uplo() == upper# and L (A=L*L') if #A.uplo() == lower# of the matrix.

If the matrix is not positive definite a #FactorizationException# will be thrown with an error message.

Definition at line 53 of file DenseLinAlgLAPack.cpp.

void DenseLinAlgLAPack::geqrf ( DMatrixSlice *  A,
DVectorSlice *  tau,
DVectorSlice *  work 
)

Calls xGEQRF to compute the QR factorization of a matrix A.

The factors A = Q*R are stored in the storage space for A and in an extra vector tau of length k = min(A.rows(),A.cols()). The matrix R is stored in the upper diagonal of A on exit and Q is stored as a set of elementary reflectors in the elements of A below the diagonal and in the extra vector tau.

See LAPACK documentation of xGEQRF for more detail of these parameters.

If there is a problem with one of the arguments a std::invalid_argument exception will be thrown.

Parameters
A[in/out] The mxn matrix to be factorized on input. On output contains R and part of Q.
tau[out] DVector of length min(m,n) for the other part of Q.
work[work] Workspace vector of length > N.

Definition at line 73 of file DenseLinAlgLAPack.cpp.

void DenseLinAlgLAPack::ormrq ( BLAS_Cpp::Side  side,
BLAS_Cpp::Transp  trans,
const DMatrixSlice &  A,
const DVectorSlice &  tau,
DMatrixSlice *  C,
DVectorSlice *  work 
)

Calls xORMRQ to compute a matrix matrix op(C)=op(Q)*op(C) where Q is stored in A and tau computed from geqrf(..).

See the LAPACK documentation for xQRMRQ to see a description of the parameters.

Parameters
side[in] Determines if op(Q) is on the right (C=op(Q)*C) or on the left (C=C*op(Q)).
trans[in] Determines if op(Q) = Q or Q'.
A[in] Contains part of Q stored in the lower diagonal elements as set by geqrf(...).
tau[in] Contains the other part of Q as set by geqrf(...).
C[in/out] On input contains C and on output contains op(Q)*C (left) or C*op(Q) (right).
work[work] Workspace vector of length > N if side==left or > M if side==right.

Definition at line 91 of file DenseLinAlgLAPack.cpp.

void DenseLinAlgLAPack::sytrf ( DMatrixSliceTriEle *  A,
FortranTypes::f_int  ipiv[],
DVectorSlice *  work 
)

Calls xSYTRF to compute the P*A*P' = L'*D*L factorization of a symmetric possibly indefinite matrix.

See LAPACK documentation for xSYTRF so see a description of the function.

This factorization can be used to solve linear systems by calling sytrs(...). This function will throw #FactorizationException# if the matrix is singular.

Parameters
A[in/out] On input contains the elements of the symmetric matrix to be factorized. On output contains the factorization.
ipiv[out] Array of length >= A.rows(). On output contains permutation matrix P.
work[work] Workspace used to compute factoization. The length must be greater than 0 but optimal is n * NB. NB is the optimal block size (64?), n = A.rows().

Definition at line 111 of file DenseLinAlgLAPack.cpp.

void DenseLinAlgLAPack::sytrs ( const DMatrixSliceTriEle &  A,
FortranTypes::f_int  ipiv[],
DMatrixSlice *  B,
DVectorSlice *  work 
)

Calls xSYTRS(...) to compute the solution of the factorized system A * X = B where A was factorized by xSYTRF(...).

See LAPACK documentation for xSYTRS so see a description of the function.

This factorization can be used to solve linear systems by calling sytrs(...).

Parameters
A[in] Factorization of A computed by sytrf(...).
ipiv[out] Array of length >= A.rows() computed by sytrf(...).
B[in/out] On input contains the rhs matrix B. On output contiains the solution matrix X.
work[work] Workspace used to compute factoization. The length must be greater than 0 but optimal is n * NB. NB is the optimal block size (64?), n = A.rows().

Definition at line 131 of file DenseLinAlgLAPack.cpp.

void DenseLinAlgLAPack::getrf ( DMatrixSlice *  A,
FortranTypes::f_int  ipiv[],
FortranTypes::f_int *  rank 
)

Calls xGETRF to compute the P'*A = L*U factorization of a general rectuangular matrix.

See LAPACK documentation for xGETRF so see a description of the function.

This factorization can be used to solve linear systems by calling getrs(...). This function will throw #FactorizationException# if the matrix is singular.

Parameters
A[in/out] On input contains the elements of the general matrix to be factorized. On output contains the factorization.
ipiv[out] Array of length >= A.rows(). On output contains permutation matrix P.
rank[out] Returns the rank of the matrix.

Definition at line 151 of file DenseLinAlgLAPack.cpp.

void DenseLinAlgLAPack::getrs ( const DMatrixSlice &  LU,
const FortranTypes::f_int  ipiv[],
BLAS_Cpp::Transp  transp,
DMatrixSlice *  B 
)

Calls xGETRS to solve linear systems with the factors of P'*A = L*U generated by xGETRF.

See LAPACK documentation for xGETRS so see a description of the function.

This factorization can be used to solve linear systems by calling getrs(...). This function will throw #FactorizationException# if the matrix is singular.

Parameters
LU[in] On input contains the elements of the LU factors.
ipiv[out] Array of length >= LU.rows() that contains permutation matrix P.
tranp[in] Determines of (P*L*U) * X = Y (no_trans) is solved or (P*L*U)' * X = Y (trans) is solved.
B[in/out] On input, contains the rhs matrix Y, on output, contains the solution matix X.

Definition at line 169 of file DenseLinAlgLAPack.cpp.