PaStiX Handbook 6.4.0
|
#include "pastix_lowrank.h"
Go to the source code of this file.
Data Structures | |
struct | core_zlrmm_s |
Structure to store all the parameters of the core_zlrmm family functions. More... | |
Typedefs | |
PastixComplex64 main template to convert a full rank matrix to low-rank | |
typedef int(* | core_zrrqr_cp_t) (double tol, pastix_int_t maxrank, int refine, pastix_int_t nb, pastix_int_t m, pastix_int_t n, pastix_complex64_t *A, pastix_int_t lda, pastix_int_t *jpvt, pastix_complex64_t *tau, pastix_complex64_t *work, pastix_int_t lwork, double *rwork) |
TODO. | |
typedef int(* | core_zrrqr_rt_t) (double tol, pastix_int_t maxrank, pastix_int_t nb, pastix_int_t m, pastix_int_t n, pastix_complex64_t *A, pastix_int_t lda, pastix_complex64_t *tau, pastix_complex64_t *B, pastix_int_t ldb, pastix_complex64_t *tau_b, pastix_complex64_t *work, pastix_int_t lwork, double normA) |
TODO. | |
Functions | |
pastix_fixdbl_t | core_zlrmm (core_zlrmm_t *params) |
Compute the matrix matrix product when involved matrices are stored in a low-rank structure. | |
PastixComplex64 low-rank kernels | |
void | core_zlralloc (pastix_int_t M, pastix_int_t N, pastix_int_t rkmax, pastix_lrblock_t *A) |
Allocate a low-rank matrix. | |
void | core_zlrfree (pastix_lrblock_t *A) |
Free a low-rank matrix. | |
int | core_zlrsze (int copy, pastix_int_t M, pastix_int_t N, pastix_lrblock_t *A, pastix_int_t newrk, pastix_int_t newrkmax, pastix_int_t rklimit) |
Resize a low-rank matrix. | |
int | core_zlr2ge (pastix_trans_t trans, pastix_int_t M, pastix_int_t N, const pastix_lrblock_t *Alr, pastix_complex64_t *A, pastix_int_t lda) |
Convert a low rank matrix into a dense matrix. | |
void | core_zlrcpy (const pastix_lr_t *lowrank, pastix_trans_t transA, pastix_complex64_t alpha, pastix_int_t M1, pastix_int_t N1, const pastix_lrblock_t *A, pastix_int_t M2, pastix_int_t N2, pastix_lrblock_t *B, pastix_int_t offx, pastix_int_t offy) |
Copy a small low-rank structure into a large one. | |
void | core_zlrconcatenate_u (pastix_complex64_t alpha, pastix_int_t M1, pastix_int_t N1, const pastix_lrblock_t *A, pastix_int_t M2, pastix_lrblock_t *B, pastix_int_t offx, pastix_complex64_t *u1u2) |
Concatenate left parts of two low-rank matrices. | |
void | core_zlrconcatenate_v (pastix_trans_t transA1, pastix_complex64_t alpha, pastix_int_t M1, pastix_int_t N1, const pastix_lrblock_t *A, pastix_int_t N2, pastix_lrblock_t *B, pastix_int_t offy, pastix_complex64_t *v1v2) |
Concatenate right parts of two low-rank matrices. | |
double | core_zlrnrm (pastix_normtype_t ntype, int transV, pastix_int_t M, pastix_int_t N, const pastix_lrblock_t *A) |
Compute the norm of a low-rank matrix. | |
size_t | core_zlrgetsize (pastix_int_t M, pastix_int_t N, pastix_lrblock_t *A) |
Compute the size of a block to send in LR. | |
char * | core_zlrpack (pastix_int_t M, pastix_int_t N, const pastix_lrblock_t *A, char *buffer) |
Pack low-rank data by side. | |
char * | core_zlrunpack (pastix_int_t M, pastix_int_t N, pastix_lrblock_t *A, char *buffer) |
Unpack low rank data and fill the cblk concerned by the computation. | |
const char * | core_zlrunpack2 (pastix_int_t M, pastix_int_t N, pastix_lrblock_t *A, const char *input, char **outptr) |
Unpack low rank data and fill the cblk concerned by the computation. | |
update_fr Functions to perform the update on a full-rank matrix | |
pastix_fixdbl_t | core_zfrfr2fr (core_zlrmm_t *params) |
Perform the full-rank operation C = alpha * op(A) * op(B) + beta C. | |
pastix_fixdbl_t | core_zfrlr2fr (core_zlrmm_t *params) |
Perform the operation C = alpha * op(A) * op(B) + beta C, with A and C full-rank and B low-rank. | |
pastix_fixdbl_t | core_zlrfr2fr (core_zlrmm_t *params) |
Perform the operation C = alpha * op(A) * op(B) + beta C, with B and C full-rank and A low-rank. | |
pastix_fixdbl_t | core_zlrlr2fr (core_zlrmm_t *params) |
Perform the operation C = alpha * op(A) * op(B) + beta C, with A and B low-rank and C full-rank. | |
update_lr Functions to prepare the AB product for an update on a low-rank matrix | |
pastix_fixdbl_t | core_zfrfr2lr (core_zlrmm_t *params, pastix_lrblock_t *AB, int *infomask, pastix_int_t Kmax) |
Perform the operation AB = op(A) * op(B), with A and B full-rank and AB low-rank. | |
pastix_fixdbl_t | core_zfrlr2lr (core_zlrmm_t *params, pastix_lrblock_t *AB, int *infomask, pastix_int_t Brkmin) |
Perform the operation AB = op(A) * op(B), with A full-rank and B and AB low-rank. | |
pastix_fixdbl_t | core_zlrfr2lr (core_zlrmm_t *params, pastix_lrblock_t *AB, int *infomask, pastix_int_t Arkmin) |
Perform the operation AB = op(A) * op(B), with B full-rank and A and AB low-rank. | |
pastix_fixdbl_t | core_zlrlr2lr (core_zlrmm_t *params, pastix_lrblock_t *AB, int *infomask) |
Perform the operation AB = op(A) * op(B), with A, B, and AB low-rank. | |
add_lr Functions to add the AB contribution in a low-rank format to any C matrix | |
pastix_fixdbl_t | core_zlradd (core_zlrmm_t *params, const pastix_lrblock_t *AB, pastix_trans_t transV, int infomask) |
Perform the addition of two low-rank matrices. | |
PastixComplex64 SVD low-rank kernels | |
pastix_fixdbl_t | core_zge2lr_svd (int use_reltol, pastix_fixdbl_t tol, pastix_int_t rklimit, pastix_int_t m, pastix_int_t n, const void *Avoid, pastix_int_t lda, pastix_lrblock_t *Alr) |
Convert a full rank matrix in a low rank matrix, using SVD. | |
pastix_fixdbl_t | core_zrradd_svd (const pastix_lr_t *lowrank, pastix_trans_t transA1, const void *alphaptr, pastix_int_t M1, pastix_int_t N1, const pastix_lrblock_t *A, pastix_int_t M2, pastix_int_t N2, pastix_lrblock_t *B, pastix_int_t offx, pastix_int_t offy) |
Add two LR structures A=(-u1) v1^T and B=u2 v2^T into u2 v2^T. | |
PastixComplex64 Rank Revealing QR kernels for low-rank | |
pastix_fixdbl_t | core_zge2lr_pqrcp (int use_reltol, pastix_fixdbl_t tol, pastix_int_t rklimit, pastix_int_t m, pastix_int_t n, const void *Avoid, pastix_int_t lda, pastix_lrblock_t *Alr) |
Convert a full rank matrix in a low rank matrix, using PQRCP. | |
pastix_fixdbl_t | core_zrradd_pqrcp (const pastix_lr_t *lowrank, pastix_trans_t transA1, const void *alphaptr, pastix_int_t M1, pastix_int_t N1, const pastix_lrblock_t *A, pastix_int_t M2, pastix_int_t N2, pastix_lrblock_t *B, pastix_int_t offx, pastix_int_t offy) |
Add two LR structures A=(-u1) v1^T and B=u2 v2^T into u2 v2^T. | |
pastix_fixdbl_t | core_zge2lr_rqrcp (int use_reltol, pastix_fixdbl_t tol, pastix_int_t rklimit, pastix_int_t m, pastix_int_t n, const void *Avoid, pastix_int_t lda, pastix_lrblock_t *Alr) |
Convert a full rank matrix in a low rank matrix, using RQRCP. | |
pastix_fixdbl_t | core_zrradd_rqrcp (const pastix_lr_t *lowrank, pastix_trans_t transA1, const void *alphaptr, pastix_int_t M1, pastix_int_t N1, const pastix_lrblock_t *A, pastix_int_t M2, pastix_int_t N2, pastix_lrblock_t *B, pastix_int_t offx, pastix_int_t offy) |
Add two LR structures A=(-u1) v1^T and B=u2 v2^T into u2 v2^T. | |
pastix_fixdbl_t | core_zge2lr_tqrcp (int use_reltol, pastix_fixdbl_t tol, pastix_int_t rklimit, pastix_int_t m, pastix_int_t n, const void *Avoid, pastix_int_t lda, pastix_lrblock_t *Alr) |
Convert a full rank matrix in a low rank matrix, using TQRCP. | |
pastix_fixdbl_t | core_zrradd_tqrcp (const pastix_lr_t *lowrank, pastix_trans_t transA1, const void *alphaptr, pastix_int_t M1, pastix_int_t N1, const pastix_lrblock_t *A, pastix_int_t M2, pastix_int_t N2, pastix_lrblock_t *B, pastix_int_t offx, pastix_int_t offy) |
Add two LR structures A=(-u1) v1^T and B=u2 v2^T into u2 v2^T. | |
pastix_fixdbl_t | core_zge2lr_rqrrt (int use_reltol, pastix_fixdbl_t tol, pastix_int_t rklimit, pastix_int_t m, pastix_int_t n, const void *Avoid, pastix_int_t lda, pastix_lrblock_t *Alr) |
Convert a full rank matrix in a low rank matrix, using RQRRT. | |
pastix_fixdbl_t | core_zge2lr_qrcp (core_zrrqr_cp_t rrqrfct, int use_reltol, pastix_fixdbl_t tol, pastix_int_t rklimit, pastix_int_t m, pastix_int_t n, const void *Avoid, pastix_int_t lda, pastix_lrblock_t *Alr) |
Template to convert a full rank matrix into a low rank matrix through QR decompositions. | |
pastix_fixdbl_t | core_zge2lr_qrrt (core_zrrqr_rt_t rrqrfct, int use_reltol, pastix_fixdbl_t tol, pastix_int_t rklimit, pastix_int_t m, pastix_int_t n, const void *Avoid, pastix_int_t lda, pastix_lrblock_t *Alr) |
Template to convert a full rank matrix into a low rank matrix through QR decompositions. | |
pastix_fixdbl_t | core_zrradd_qr (core_zrrqr_cp_t rrqrfct, const pastix_lr_t *lowrank, pastix_trans_t transA1, const void *alphaptr, pastix_int_t M1, pastix_int_t N1, const pastix_lrblock_t *A, pastix_int_t M2, pastix_int_t N2, pastix_lrblock_t *B, pastix_int_t offx, pastix_int_t offy) |
Template to perform the addition of two low-rank structures with compression kernel based on QR decomposition. | |
PastixComplex64 low-rank debug functions | |
void | core_zlrdbg_printsvd (pastix_int_t M, pastix_int_t N, const pastix_complex64_t *A, pastix_int_t lda) |
Print the svd values of the given matrix. | |
int | core_zlrdbg_check_orthogonality (pastix_int_t M, pastix_int_t N, const pastix_complex64_t *A, pastix_int_t lda) |
Check the orthogonality of the matrix A. | |
int | core_zlrdbg_check_orthogonality_AB (pastix_int_t M, pastix_int_t NA, pastix_int_t NB, const pastix_complex64_t *A, pastix_int_t lda, const pastix_complex64_t *B, pastix_int_t ldb) |
Check the orthogonality of the matrix A relatively to the matrix B. | |
PastixComplex64 LRMM low-rank kernels | |
#define | PASTE_CORE_ZLRMM_PARAMS(_a_) |
Initialize all the parameters of the core_zlrmm family functions to ease the access. | |
#define | PASTE_CORE_ZLRMM_VOID |
Void all the parameters of the core_zlrmm family functions to silent warnings. | |
typedef struct core_zlrmm_s | core_zlrmm_t |
Structure to store all the parameters of the core_zlrmm family functions. | |
static pastix_complex64_t * | core_zlrmm_getws (core_zlrmm_t *params, ssize_t newsize) |
Function to get a workspace pointer if space is available in the one provided. | |
PaStiX kernel header.
Definition in file pastix_zlrcores.h.