Go to the documentation of this file.
19 #ifndef _pastix_slrcores_h_
20 #define _pastix_slrcores_h_
36 int core_slrsze (
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 );
43 pastix_int_t offx, pastix_int_t offy );
57 pastix_int_t M, pastix_int_t N,
99 pastix_atomic_lock_t *
lock;
105 #define PASTE_CORE_SLRMM_PARAMS(_a_) \
106 const pastix_lr_t *lowrank = (_a_)->lowrank; \
107 pastix_trans_t transA = (_a_)->transA; \
108 pastix_trans_t transB = (_a_)->transB; \
109 pastix_int_t M = (_a_)->M; \
110 pastix_int_t N = (_a_)->N; \
111 pastix_int_t K = (_a_)->K; \
112 pastix_int_t Cm = (_a_)->Cm; \
113 pastix_int_t Cn = (_a_)->Cn; \
114 pastix_int_t offx = (_a_)->offx; \
115 pastix_int_t offy = (_a_)->offy; \
116 float alpha = (_a_)->alpha; \
117 const pastix_lrblock_t *A = (_a_)->A; \
118 const pastix_lrblock_t *B = (_a_)->B; \
119 float beta = (_a_)->beta; \
120 pastix_lrblock_t *C = (_a_)->C; \
121 float *work = (_a_)->work; \
122 pastix_int_t lwork = (_a_)->lwork; \
123 pastix_atomic_lock_t *lock = (_a_)->lock;
128 #define PASTE_CORE_SLRMM_VOID \
154 static inline float *
162 params->
lwused += newsize;
198 pastix_int_t Brkmin );
202 pastix_int_t Arkmin );
236 pastix_fixdbl_t
core_sge2lr_svd(
int use_reltol, pastix_fixdbl_t tol, pastix_int_t rklimit,
237 pastix_int_t m, pastix_int_t n,
242 pastix_int_t offx, pastix_int_t offy);
257 typedef int (*core_srrqr_cp_t)(
float tol, pastix_int_t maxrank,
int refine, pastix_int_t nb,
258 pastix_int_t m, pastix_int_t n,
259 float *A, pastix_int_t lda,
260 pastix_int_t *jpvt,
float *tau,
261 float *work, pastix_int_t lwork,
float *rwork );
263 typedef int (*core_srrqr_rt_t)(
float tol, pastix_int_t maxrank, pastix_int_t nb,
264 pastix_int_t m, pastix_int_t n,
265 float *A, pastix_int_t lda,
float *tau,
266 float *B, pastix_int_t ldb,
float *tau_b,
267 float *work, pastix_int_t lwork,
float normA );
275 pastix_fixdbl_t
core_sge2lr_pqrcp(
int use_reltol, pastix_fixdbl_t tol, pastix_int_t rklimit,
276 pastix_int_t m, pastix_int_t n,
281 pastix_int_t offx, pastix_int_t offy );
283 pastix_fixdbl_t
core_sge2lr_rqrcp(
int use_reltol, pastix_fixdbl_t tol, pastix_int_t rklimit,
284 pastix_int_t m, pastix_int_t n,
289 pastix_int_t offx, pastix_int_t offy );
291 pastix_fixdbl_t
core_sge2lr_tqrcp(
int use_reltol, pastix_fixdbl_t tol, pastix_int_t rklimit,
292 pastix_int_t m, pastix_int_t n,
297 pastix_int_t offx, pastix_int_t offy );
299 pastix_fixdbl_t
core_sge2lr_rqrrt(
int use_reltol, pastix_fixdbl_t tol, pastix_int_t rklimit,
300 pastix_int_t m, pastix_int_t n,
305 int use_reltol, pastix_fixdbl_t tol, pastix_int_t rklimit,
306 pastix_int_t m, pastix_int_t n,
307 const void *Avoid, pastix_int_t lda,
310 int use_reltol, pastix_fixdbl_t tol, pastix_int_t rklimit,
311 pastix_int_t m, pastix_int_t n,
312 const void *Avoid, pastix_int_t lda,
319 pastix_int_t offx, pastix_int_t offy );
344 const float *A, pastix_int_t lda,
345 const float *B, pastix_int_t ldb );
Structure to define the type of function to use for the low-rank kernels and their parameters.
const pastix_lrblock_t * B
pastix_fixdbl_t core_sge2lr_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.
void core_slrconcatenate_u(float 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, float *u1u2)
Concatenate left parts of two low-rank matrices.
void core_slralloc(pastix_int_t M, pastix_int_t N, pastix_int_t rkmax, pastix_lrblock_t *A)
Allocate a low-rank matrix.
enum pastix_trans_e pastix_trans_t
Transpostion.
pastix_fixdbl_t core_srradd_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.
int core_slrdbg_check_orthogonality(pastix_int_t M, pastix_int_t N, const float *A, pastix_int_t lda)
Check the orthogonality of the matrix A.
pastix_fixdbl_t core_sge2lr_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_sfrfr2lr(core_slrmm_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.
The block low-rank structure to hold a matrix in low-rank form.
pastix_fixdbl_t core_sfrfr2fr(core_slrmm_t *params)
Perform the full-rank operation C = alpha * op(A) * op(B) + beta C.
pastix_fixdbl_t core_sge2lr_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.
void core_slrfree(pastix_lrblock_t *A)
Free a low-rank matrix.
void core_slrcpy(const pastix_lr_t *lowrank, pastix_trans_t transA, float 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.
Structure to store all the parameters of the core_slrmm family functions.
pastix_fixdbl_t core_sge2lr_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_sge2lr_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_slradd(core_slrmm_t *params, const pastix_lrblock_t *AB, pastix_trans_t transV, int infomask)
Perform the addition of two low-rank matrices.
enum pastix_normtype_e pastix_normtype_t
Norms.
const pastix_lrblock_t * A
struct core_slrmm_s core_slrmm_t
Structure to store all the parameters of the core_slrmm family functions.
pastix_fixdbl_t core_slrmm(core_slrmm_t *params)
Compute the matrix matrix product when involved matrices are stored in a low-rank structure.
char * core_slrpack(pastix_int_t M, pastix_int_t N, const pastix_lrblock_t *A, char *buffer)
Pack low-rank data by side.
int core_slrsze(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.
void core_slrdbg_printsvd(pastix_int_t M, pastix_int_t N, const float *A, pastix_int_t lda)
Print the svd values of the given matrix.
pastix_fixdbl_t core_slrlr2fr(core_slrmm_t *params)
Perform the operation C = alpha * op(A) * op(B) + beta C, with A and B low-rank and C full-rank.
pastix_fixdbl_t core_srradd_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.
pastix_fixdbl_t core_slrfr2lr(core_slrmm_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.
int core_slrdbg_check_orthogonality_AB(pastix_int_t M, pastix_int_t NA, pastix_int_t NB, const float *A, pastix_int_t lda, const float *B, pastix_int_t ldb)
Check the orthogonality of the matrix A relatively to the matrix B.
pastix_fixdbl_t core_srradd_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.
const char * core_slrunpack2(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.
void core_slrconcatenate_v(pastix_trans_t transA1, float 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, float *v1v2)
Concatenate right parts of two low-rank matrices.
pastix_fixdbl_t core_srradd_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_sge2lr_qrcp(core_srrqr_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_sfrlr2lr(core_slrmm_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.
size_t core_slrgetsize(pastix_int_t M, pastix_int_t N, pastix_lrblock_t *A)
Compute the size of a block to send in LR.
const pastix_lr_t * lowrank
char * core_slrunpack(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.
static float * core_slrmm_getws(core_slrmm_t *params, ssize_t newsize)
Function to get a workspace pointer if space is available in the one provided.
pastix_fixdbl_t core_srradd_qr(core_srrqr_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 decom...
pastix_fixdbl_t core_slrfr2fr(core_slrmm_t *params)
Perform the operation C = alpha * op(A) * op(B) + beta C, with B and C full-rank and A low-rank.
int core_slr2ge(pastix_trans_t trans, pastix_int_t M, pastix_int_t N, const pastix_lrblock_t *Alr, float *A, pastix_int_t lda)
Convert a low rank matrix into a dense matrix.
pastix_fixdbl_t core_sfrlr2fr(core_slrmm_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_sge2lr_qrrt(core_srrqr_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_slrlr2lr(core_slrmm_t *params, pastix_lrblock_t *AB, int *infomask)
Perform the operation AB = op(A) * op(B), with A, B, and AB low-rank.
float core_slrnrm(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.
pastix_atomic_lock_t * lock