PaStiX Handbook  6.3.0
pastix_zcores.h
Go to the documentation of this file.
1 /**
2  * @file pastix_zcores.h
3  *
4  * PaStiX kernel header.
5  *
6  * @copyright 2011-2023 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
7  * Univ. Bordeaux. All rights reserved.
8  *
9  * @version 6.3.0
10  * @author Mathieu Faverge
11  * @author Pierre Ramet
12  * @author Xavier Lacoste
13  * @author Esragul Korkmaz
14  * @author Gregoire Pichon
15  * @author Tony Delarue
16  * @author Alycia Lisito
17  * @date 2023-03-30
18  * @generated from /builds/solverstack/pastix/kernels/pastix_zcores.h, normal z -> z, Mon Aug 28 13:40:23 2023
19  *
20  */
21 #ifndef _pastix_zcores_h_
22 #define _pastix_zcores_h_
23 
24 #ifndef DOXYGEN_SHOULD_SKIP_THIS
25 #define pastix_cblk_lock( cblk_ ) pastix_atomic_lock( &((cblk_)->lock) )
26 #define pastix_cblk_unlock( cblk_ ) pastix_atomic_unlock( &((cblk_)->lock) )
27 #endif /* DOXYGEN_SHOULD_SKIP_THIS */
28 
29 /**
30  * @addtogroup kernel_blas_lapack
31  * @{
32  * This module contains all the BLAS and LAPACK-like kernels that are working
33  * on lapack layout matrices.
34  *
35  * @name PastixComplex64 BLAS kernels
36  * @{
37  */
38 void core_zplrnt( int m,
39  int n,
40  pastix_complex64_t *A,
41  int lda,
42  int gM,
43  int m0,
44  int n0,
45  unsigned long long int seed );
46 void core_zgetmo( int m,
47  int n,
48  const pastix_complex64_t *A,
49  int lda,
50  pastix_complex64_t *B,
51  int ldb );
52 int core_zgeadd( pastix_trans_t trans,
53  pastix_int_t M,
54  pastix_int_t N,
55  pastix_complex64_t alpha,
56  const pastix_complex64_t *A,
57  pastix_int_t LDA,
58  pastix_complex64_t beta,
59  pastix_complex64_t *B,
60  pastix_int_t LDB );
61 int core_zgemdm( pastix_trans_t transA,
62  pastix_trans_t transB,
63  int M,
64  int N,
65  int K,
66  pastix_complex64_t alpha,
67  const pastix_complex64_t *A,
68  int LDA,
69  const pastix_complex64_t *B,
70  int LDB,
71  pastix_complex64_t beta,
72  pastix_complex64_t *C,
73  int LDC,
74  const pastix_complex64_t *D,
75  int incD,
76  pastix_complex64_t *WORK,
77  int LWORK );
78 int core_zpqrcp( double tol,
79  pastix_int_t maxrank,
80  int full_update,
81  pastix_int_t nb,
82  pastix_int_t m,
83  pastix_int_t n,
84  pastix_complex64_t *A,
85  pastix_int_t lda,
86  pastix_int_t *jpvt,
87  pastix_complex64_t *tau,
88  pastix_complex64_t *work,
89  pastix_int_t lwork,
90  double *rwork );
91 int core_zrqrcp( double tol,
92  pastix_int_t maxrank,
93  int refine,
94  pastix_int_t nb,
95  pastix_int_t m,
96  pastix_int_t n,
97  pastix_complex64_t *A,
98  pastix_int_t lda,
99  pastix_int_t *jpvt,
100  pastix_complex64_t *tau,
101  pastix_complex64_t *work,
102  pastix_int_t lwork,
103  double *rwork );
104 int core_zrqrrt( double tol,
105  pastix_int_t maxrank,
106  pastix_int_t nb,
107  pastix_int_t m,
108  pastix_int_t n,
109  pastix_complex64_t *A,
110  pastix_int_t lda,
111  pastix_complex64_t *tau,
112  pastix_complex64_t *B,
113  pastix_int_t ldb,
114  pastix_complex64_t *tau_b,
115  pastix_complex64_t *work,
116  pastix_int_t lwork,
117  double normA );
118 int core_ztqrcp( double tol,
119  pastix_int_t maxrank,
120  int unused,
121  pastix_int_t nb,
122  pastix_int_t m,
123  pastix_int_t n,
124  pastix_complex64_t *A,
125  pastix_int_t lda,
126  pastix_int_t *jpvt,
127  pastix_complex64_t *tau,
128  pastix_complex64_t *work,
129  pastix_int_t lwork,
130  double *rwork );
131 int core_ztradd( pastix_uplo_t uplo,
132  pastix_trans_t trans,
133  pastix_int_t M,
134  pastix_int_t N,
135  pastix_complex64_t alpha,
136  const pastix_complex64_t *A,
137  pastix_int_t LDA,
138  pastix_complex64_t beta,
139  pastix_complex64_t *B,
140  pastix_int_t LDB);
141 int core_zscalo( pastix_trans_t trans,
142  pastix_int_t M,
143  pastix_int_t N,
144  const pastix_complex64_t *A,
145  pastix_int_t lda,
146  const pastix_complex64_t *D,
147  pastix_int_t ldd,
148  pastix_complex64_t *B,
149  pastix_int_t ldb );
150 
151 /**
152  * @}
153  * @name PastixComplex64 Othogonalization kernels for low-rank updates
154  * @{
155  */
156 pastix_fixdbl_t core_zlrorthu_fullqr( pastix_int_t M,
157  pastix_int_t N,
158  pastix_int_t rank,
159  pastix_complex64_t *U,
160  pastix_int_t ldu,
161  pastix_complex64_t *V,
162  pastix_int_t ldv );
163 pastix_fixdbl_t core_zlrorthu_partialqr( pastix_int_t M,
164  pastix_int_t N,
165  pastix_int_t r1,
166  pastix_int_t *r2ptr,
167  pastix_int_t offx,
168  pastix_int_t offy,
169  pastix_complex64_t *U,
170  pastix_int_t ldu,
171  pastix_complex64_t *V,
172  pastix_int_t ldv );
173 pastix_fixdbl_t core_zlrorthu_cgs( pastix_int_t M1,
174  pastix_int_t N1,
175  pastix_int_t M2,
176  pastix_int_t N2,
177  pastix_int_t r1,
178  pastix_int_t *r2ptr,
179  pastix_int_t offx,
180  pastix_int_t offy,
181  pastix_complex64_t *U,
182  pastix_int_t ldu,
183  pastix_complex64_t *V,
184  pastix_int_t ldv );
185 
186 /**
187  * @}
188  * @name PastixComplex64 LAPACK kernels
189  * @{
190  */
191 void core_zpotrfsp( pastix_int_t n,
192  pastix_complex64_t *A,
193  pastix_int_t lda,
194  pastix_int_t *nbpivots,
195  double criterion );
196 void core_zpxtrfsp( pastix_int_t n,
197  pastix_complex64_t *A,
198  pastix_int_t lda,
199  pastix_int_t *nbpivots,
200  double criterion );
201 void core_zgetrfsp( pastix_int_t n,
202  pastix_complex64_t *A,
203  pastix_int_t lda,
204  pastix_int_t *nbpivots,
205  double criterion );
206 #if defined(PRECISION_z) || defined(PRECISION_c)
207 void core_zhetrfsp( pastix_int_t n,
208  pastix_complex64_t *A,
209  pastix_int_t lda,
210  pastix_int_t *nbpivots,
211  double criterion );
212 #endif
213 void core_zsytrfsp( pastix_int_t n,
214  pastix_complex64_t *A,
215  pastix_int_t lda,
216  pastix_int_t *nbpivots,
217  double criterion );
218 
219 /**
220  * @}
221  * @}
222  *
223  * @addtogroup kernel_fact
224  * @{
225  * This module contains all the kernel working at the solver matrix structure
226  * level for the numerical factorization step.
227  *
228  * @name PastixComplex64 cblk-BLAS CPU kernels
229  * @{
230  */
231 
232 int cpucblk_zgeaddsp1d( const SolverCblk *cblk1,
233  SolverCblk *cblk2,
234  const pastix_complex64_t *L1,
235  pastix_complex64_t *L2,
236  const pastix_complex64_t *U1,
237  pastix_complex64_t *U2 );
238 
239 pastix_fixdbl_t cpucblk_zgemmsp( pastix_coefside_t sideA,
240  pastix_trans_t trans,
241  const SolverCblk *cblk,
242  const SolverBlok *blok,
243  SolverCblk *fcblk,
244  const void *A,
245  const void *B,
246  void *C,
247  pastix_complex64_t *work,
248  pastix_int_t lwork,
249  const pastix_lr_t *lowrank );
250 void cpucblk_ztrsmsp( pastix_side_t side,
251  pastix_uplo_t uplo,
252  pastix_trans_t trans,
253  pastix_diag_t diag,
254  const SolverCblk *cblk,
255  const void *A,
256  void *C,
257  const pastix_lr_t *lowrank );
258 void cpucblk_zscalo ( pastix_trans_t trans,
259  SolverCblk *cblk,
260  void *dataL,
261  void *dataLD );
262 
263 pastix_fixdbl_t cpublok_zgemmsp( pastix_trans_t trans,
264  const SolverCblk *cblk,
265  SolverCblk *fcblk,
266  pastix_int_t blok_mk,
267  pastix_int_t blok_nk,
268  pastix_int_t blok_mn,
269  const void *A,
270  const void *B,
271  void *C,
272  const pastix_lr_t *lowrank );
273 pastix_fixdbl_t cpublok_ztrsmsp( pastix_side_t side,
274  pastix_uplo_t uplo,
275  pastix_trans_t trans,
276  pastix_diag_t diag,
277  const SolverCblk *cblk,
278  pastix_int_t blok_m,
279  const void *A,
280  void *C,
281  const pastix_lr_t *lowrank );
282 void cpublok_zscalo ( pastix_trans_t trans,
283  SolverCblk *cblk,
284  pastix_int_t blok_m,
285  const void *A,
286  const void *dataD,
287  void *dataB );
288 
289 /**
290  * @}
291  * @name PastixComplex64 cblk LU kernels
292  * @{
293  */
294 int cpucblk_zgetrfsp1d_getrf( SolverMatrix *solvmtx,
295  SolverCblk *cblk,
296  void *L,
297  void *U );
298 int cpucblk_zgetrfsp1d_panel( SolverMatrix *solvmtx,
299  SolverCblk *cblk,
300  void *L,
301  void *U );
302 int cpucblk_zgetrfsp1d ( SolverMatrix *solvmtx,
303  SolverCblk *cblk,
304  pastix_complex64_t *work,
305  pastix_int_t lwork );
306 
307 /**
308  * @}
309  * @name PastixComplex64 cblk Cholesky kernels
310  * @{
311  */
312 int cpucblk_zpotrfsp1d_potrf( SolverMatrix *solvmtx,
313  SolverCblk *cblk,
314  void *dataL );
315 int cpucblk_zpotrfsp1d_panel( SolverMatrix *solvmtx,
316  SolverCblk *cblk,
317  void *L );
318 int cpucblk_zpotrfsp1d( SolverMatrix *solvmtx,
319  SolverCblk *cblk,
320  pastix_complex64_t *work,
321  pastix_int_t lwork );
322 
323 /**
324  * @}
325  */
326 
327 #if defined(PRECISION_z) || defined(PRECISION_c)
328  /**
329  * @name PastixComplex64 cblk LDL^h kernels
330  * @{
331  */
332 int cpucblk_zhetrfsp1d_hetrf( SolverMatrix *solvmtx,
333  SolverCblk *cblk,
334  void *L );
335 int cpucblk_zhetrfsp1d_panel( SolverMatrix *solvmtx,
336  SolverCblk *cblk,
337  void *L,
338  void *DLh );
339 int cpucblk_zhetrfsp1d( SolverMatrix *solvmtx,
340  SolverCblk *cblk,
341  pastix_complex64_t *work1,
342  pastix_complex64_t *work2,
343  pastix_int_t lwork );
344 
345 /**
346  * @}
347  * @name PastixComplex64 cblk LL^t kernels
348  * @{
349  */
350 int cpucblk_zpxtrfsp1d_pxtrf( SolverMatrix *solvmtx,
351  SolverCblk *cblk,
352  void *dataL );
353 int cpucblk_zpxtrfsp1d_panel( SolverMatrix *solvmtx,
354  SolverCblk *cblk,
355  void *dataL );
356 int cpucblk_zpxtrfsp1d( SolverMatrix *solvmtx,
357  SolverCblk *cblk,
358  pastix_complex64_t *work,
359  pastix_int_t lwork );
360 
361 /**
362  * @}
363  */
364 #endif
365 
366  /**
367  * @name PastixComplex64 cblk LDL^t kernels
368  * @{
369  */
370 int cpucblk_zsytrfsp1d_sytrf( SolverMatrix *solvmtx,
371  SolverCblk *cblk,
372  void *dataL );
373 int cpucblk_zsytrfsp1d_panel( SolverMatrix *solvmtx,
374  SolverCblk *cblk,
375  void *L,
376  void *DLt );
377 int cpucblk_zsytrfsp1d( SolverMatrix *solvmtx,
378  SolverCblk *cblk,
379  pastix_complex64_t *Dlt,
380  pastix_complex64_t *work,
381  pastix_int_t lwork );
382 
383 /**
384  * @}
385  * @name PastixComplex64 initialization and additionnal routines
386  * @{
387  */
388 void cpucblk_zalloc_lrws( const SolverCblk *cblk,
389  pastix_lrblock_t *lrblok,
390  pastix_complex64_t *ws );
392  SolverCblk *cblk,
393  int rkmax );
395  SolverCblk *cblk );
397  SolverCblk *cblk );
399  SolverCblk *cblk );
401  const SolverMatrix *solvmtx,
402  const pastix_bcsc_t *bcsc,
403  pastix_int_t itercblk );
405  const SolverMatrix *solvmtx,
406  const pastix_bcsc_t *bcsc,
407  pastix_int_t itercblk,
408  const char *directory );
409 void cpucblk_zgetschur( const SolverCblk *cblk,
410  int upper_part,
411  pastix_complex64_t *S,
412  pastix_int_t lds );
414  const SolverCblk *cblk,
415  FILE *stream );
417  const SolverCblk *cblkA,
418  SolverCblk *cblkB );
420  double alpha,
421  const SolverCblk *cblkA,
422  SolverCblk *cblkB,
423  const pastix_lr_t *lowrank );
424 
425 /**
426  * @}
427  * @name PastixComplex64 MPI routines
428  * @{
429  */
430 int cpucblk_zincoming_deps( int mt_flag,
431  pastix_coefside_t side,
432  SolverMatrix *solvmtx,
433  SolverCblk *cblk );
435  SolverMatrix *solvmtx,
436  const SolverCblk *cblk,
437  SolverCblk *fcbk );
439  pastix_int_t sched,
440  SolverMatrix *solvmtx );
441 void cpucblk_zupdate_reqtab( SolverMatrix *solvmtx );
442 #if defined( PASTIX_WITH_MPI )
443 void cpucblk_zmpi_progress( pastix_coefside_t side,
444  SolverMatrix *solvmtx,
445  int threadid );
446 void cpucblk_zisend_rhs_bwd( SolverMatrix *solvmtx,
447  pastix_rhs_t rhsb,
448  SolverCblk *cblk );
449 #endif
450 void cpucblk_zmpi_rhs_fwd_progress( const args_solve_t *enums,
451  SolverMatrix *solvmtx,
452  pastix_rhs_t rhsb,
453  int threadid );
455  SolverMatrix *solvmtx,
456  pastix_rhs_t rhsb,
457  const SolverCblk *cblk,
458  SolverCblk *fcbk );
459 int cpucblk_zincoming_rhs_fwd_deps( int rank,
460  const args_solve_t *enums,
461  SolverMatrix *solvmtx,
462  SolverCblk *cblk,
463  pastix_rhs_t rhsb );
465  pastix_int_t sched,
466  SolverMatrix *solvmtx,
467  pastix_rhs_t rhsb );
468 
469 void cpucblk_zmpi_rhs_bwd_progress( const args_solve_t *enums,
470  SolverMatrix *solvmtx,
471  pastix_rhs_t rhsb,
472  int threadid );
474  SolverMatrix *solvmtx,
475  pastix_rhs_t rhsb,
476  const SolverCblk *cblk,
477  SolverCblk *fcbk );
478 int cpucblk_zincoming_rhs_bwd_deps( int rank,
479  const args_solve_t *enums,
480  SolverMatrix *solvmtx,
481  SolverCblk *cblk,
482  pastix_rhs_t rhsb );
484  pastix_int_t sched,
485  SolverMatrix *solvmtx,
486  pastix_rhs_t rhsb );
487 void cpucblk_zsend_rhs_forward( const SolverMatrix *solvmtx,
488  SolverCblk *cblk,
489  pastix_rhs_t b );
490 void cpucblk_zrecv_rhs_forward( const SolverMatrix *solvmtx,
491  SolverCblk *cblk,
492  pastix_complex64_t *work,
493  pastix_rhs_t b );
494 void cpucblk_zsend_rhs_backward( const SolverMatrix *solvmtx,
495  SolverCblk *cblk,
496  pastix_rhs_t b );
497 void cpucblk_zrecv_rhs_backward( const SolverMatrix *solvmtx,
498  SolverCblk *cblk,
499  pastix_rhs_t b );
500 
501 /**
502  * @}
503  * @name PastixComplex64 compression/uncompression routines
504  * @{
505  */
506 pastix_fixdbl_t cpublok_zcompress( const pastix_lr_t *lowrank,
507  pastix_int_t M,
508  pastix_int_t N,
509  pastix_lrblock_t *blok );
510 pastix_int_t cpucblk_zcompress( const SolverMatrix *solvmtx,
511  pastix_coefside_t side,
512  int max_ilulvl,
513  SolverCblk *cblk );
515  SolverCblk *cblk );
517  const SolverMatrix *solvmtx,
518  SolverCblk *cblk,
519  pastix_int_t *orig,
520  pastix_int_t *gain );
521 
522 /**
523  * @}
524  * @}
525  *
526  * @addtogroup kernel_solve
527  * @{
528  * This module contains all the kernel working on the solver matrix structure
529  * for the solve step.
530  *
531  */
532 
534  pastix_uplo_t uplo,
535  pastix_trans_t trans,
536  pastix_diag_t diag,
537  const SolverCblk *cblk,
538  int nrhs,
539  const void *dataA,
540  pastix_complex64_t *b,
541  int ldb );
543  pastix_trans_t trans,
544  pastix_int_t nrhs,
545  const SolverCblk *cblk,
546  const SolverBlok *blok,
547  SolverCblk *fcbk,
548  const void *dataA,
549  const pastix_complex64_t *B,
550  pastix_int_t ldb,
551  pastix_complex64_t *C,
552  pastix_int_t ldc );
553 
554 void solve_cblk_ztrsmsp_forward( const args_solve_t *enums,
555  SolverMatrix *datacode,
556  const SolverCblk *cblk,
557  pastix_rhs_t b );
558 void solve_cblk_ztrsmsp_backward( const args_solve_t *enums,
559  SolverMatrix *datacode,
560  SolverCblk *cblk,
561  pastix_rhs_t b );
562 
563 void solve_cblk_zdiag( const SolverCblk *cblk,
564  int nrhs,
565  pastix_complex64_t *b,
566  int ldb,
567  pastix_complex64_t *work );
568 /**
569  * @}
570  *
571  * @addtogroup kernel_fact_null
572  * @{
573  * This module contains the three terms update functions for the LDL^t and
574  * LDL^h factorizations.
575  *
576  */
577 #if defined(PRECISION_z) || defined(PRECISION_c)
578 void core_zhetrfsp1d_gemm( const SolverCblk *cblk,
579  const SolverBlok *blok,
580  SolverCblk *fcblk,
581  const pastix_complex64_t *L,
582  pastix_complex64_t *C,
583  pastix_complex64_t *work );
584 #endif
585 void core_zsytrfsp1d_gemm( const SolverCblk *cblk,
586  const SolverBlok *blok,
587  SolverCblk *fcblk,
588  const pastix_complex64_t *L,
589  pastix_complex64_t *C,
590  pastix_complex64_t *work );
591 
592 /**
593  * @}
594  */
595 
596 #endif /* _pastix_zcores_h_ */
int cpucblk_zhetrfsp1d(SolverMatrix *solvmtx, SolverCblk *cblk, pastix_complex64_t *DLh, pastix_complex64_t *work, pastix_int_t lwork)
Perform the LDL^h factorization of a given panel and apply all its updates.
int cpucblk_zhetrfsp1d_hetrf(SolverMatrix *solvmtx, SolverCblk *cblk, void *dataL)
Computes the LDL^h factorization of the diagonal block in a panel.
int cpucblk_zhetrfsp1d_panel(SolverMatrix *solvmtx, SolverCblk *cblk, void *L, void *DLh)
Compute the LDL^h factorization of one panel.
void core_zhetrfsp(pastix_int_t n, pastix_complex64_t *A, pastix_int_t lda, pastix_int_t *nbpivots, double criterion)
Compute the block static pivoting factorization of the hermitian matrix n-by-n A such that A = L * D ...
void core_zhetrfsp1d_gemm(const SolverCblk *cblk, const SolverBlok *blok, SolverCblk *fcblk, const pastix_complex64_t *L, pastix_complex64_t *C, pastix_complex64_t *work)
int cpucblk_zpxtrfsp1d_panel(SolverMatrix *solvmtx, SolverCblk *cblk, void *L)
Compute the LL^t factorization of one panel.
int cpucblk_zpxtrfsp1d(SolverMatrix *solvmtx, SolverCblk *cblk, pastix_complex64_t *work, pastix_int_t lwork)
Perform the LL^t factorization of a given panel and apply all its updates.
int cpucblk_zpxtrfsp1d_pxtrf(SolverMatrix *solvmtx, SolverCblk *cblk, void *dataL)
Compute the LL^t factorization of the diagonal block in a panel.
int core_zrqrcp(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)
Compute a randomized QR factorization.
Definition: core_zrqrcp.c:114
void core_zgetrfsp(pastix_int_t n, pastix_complex64_t *A, pastix_int_t lda, pastix_int_t *nbpivots, double criterion)
Compute the block static pivoting LU factorization of the matrix m-by-n A = L * U.
pastix_fixdbl_t core_zlrorthu_fullqr(pastix_int_t M, pastix_int_t N, pastix_int_t rank, pastix_complex64_t *U, pastix_int_t ldu, pastix_complex64_t *V, pastix_int_t ldv)
Try to orthognalize the u part of the low-rank form, and update the v part accordingly using full QR.
Definition: core_zlrothu.c:83
void core_zgetmo(int m, int n, const pastix_complex64_t *A, int lda, pastix_complex64_t *B, int ldb)
Transposes a m-by-n matrix out of place using an extra workspace of size m-by-n.
Definition: core_zgetmo.c:49
int core_ztqrcp(double tol, pastix_int_t maxrank, int unused, 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)
Compute a randomized QR factorization with truncated updates.
Definition: core_ztqrcp.c:99
int core_zgeadd(pastix_trans_t trans, pastix_int_t M, pastix_int_t N, pastix_complex64_t alpha, const pastix_complex64_t *A, pastix_int_t LDA, pastix_complex64_t beta, pastix_complex64_t *B, pastix_int_t LDB)
Add two matrices together.
Definition: core_zgeadd.c:78
void core_zsytrfsp(pastix_int_t n, pastix_complex64_t *A, pastix_int_t lda, pastix_int_t *nbpivots, double criterion)
Compute the block static pivoting factorization of the symmetric matrix n-by-n A such that A = L * D ...
void core_zplrnt(int m, int n, pastix_complex64_t *A, int lda, int gM, int m0, int n0, unsigned long long int seed)
Generate a random tile.
Definition: core_zplrnt.c:91
pastix_fixdbl_t core_zlrorthu_cgs(pastix_int_t M1, pastix_int_t N1, pastix_int_t M2, pastix_int_t N2, pastix_int_t r1, pastix_int_t *r2ptr, pastix_int_t offx, pastix_int_t offy, pastix_complex64_t *U, pastix_int_t ldu, pastix_complex64_t *V, pastix_int_t ldv)
Try to orthognalize the U part of the low-rank form, and update the V part accordingly using CGS.
Definition: core_zlrothu.c:430
int core_zpqrcp(double tol, pastix_int_t maxrank, int full_update, 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)
Compute a rank-reavealing QR factorization.
Definition: core_zpqrcp.c:105
int core_zrqrrt(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)
Compute a randomized QR factorization with rotation technique.
Definition: core_zrqrrt.c:126
int core_zscalo(pastix_trans_t trans, pastix_int_t M, pastix_int_t N, const pastix_complex64_t *A, pastix_int_t lda, const pastix_complex64_t *D, pastix_int_t ldd, pastix_complex64_t *B, pastix_int_t ldb)
Scale a matrix by a diagonal out of place.
Definition: core_zscalo.c:73
int core_ztradd(pastix_uplo_t uplo, pastix_trans_t trans, pastix_int_t M, pastix_int_t N, pastix_complex64_t alpha, const pastix_complex64_t *A, pastix_int_t LDA, pastix_complex64_t beta, pastix_complex64_t *B, pastix_int_t LDB)
Add two triangular matrices together as in PBLAS pztradd.
Definition: core_ztradd.c:79
int core_zgemdm(pastix_trans_t transA, pastix_trans_t transB, int M, int N, int K, pastix_complex64_t alpha, const pastix_complex64_t *A, int LDA, const pastix_complex64_t *B, int LDB, pastix_complex64_t beta, pastix_complex64_t *C, int LDC, const pastix_complex64_t *D, int incD, pastix_complex64_t *WORK, int LWORK)
Perform one of the following matrix-matrix operations.
Definition: core_zgemdm.c:139
pastix_fixdbl_t core_zlrorthu_partialqr(pastix_int_t M, pastix_int_t N, pastix_int_t r1, pastix_int_t *r2ptr, pastix_int_t offx, pastix_int_t offy, pastix_complex64_t *U, pastix_int_t ldu, pastix_complex64_t *V, pastix_int_t ldv)
Try to orthognalize the U part of the low-rank form, and update the V part accordingly using partial ...
Definition: core_zlrothu.c:197
void core_zpotrfsp(pastix_int_t n, pastix_complex64_t *A, pastix_int_t lda, pastix_int_t *nbpivots, double criterion)
Compute the block static pivoting Cholesky factorization of the matrix n-by-n A = L * L^t .
void core_zpxtrfsp(pastix_int_t n, pastix_complex64_t *A, pastix_int_t lda, pastix_int_t *nbpivots, double criterion)
Compute the block static pivoting LL^t factorization of the matrix n-by-n A = L * L^t .
void core_zsytrfsp1d_gemm(const SolverCblk *cblk, const SolverBlok *blok, SolverCblk *fcblk, const pastix_complex64_t *L, pastix_complex64_t *C, pastix_complex64_t *work)
void cpucblk_zrelease_rhs_bwd_deps(const args_solve_t *enums, SolverMatrix *solvmtx, pastix_rhs_t rhsb, const SolverCblk *cblk, SolverCblk *fcbk)
Release the dependencies of the given cblk after an update.
void cpucblk_ztrsmsp(pastix_side_t side, pastix_uplo_t uplo, pastix_trans_t trans, pastix_diag_t diag, const SolverCblk *cblk, const void *A, void *C, const pastix_lr_t *lowrank)
Compute the updates associated to a column of off-diagonal blocks.
Definition: core_ztrsmsp.c:356
void cpucblk_zrequest_rhs_bwd_cleanup(const args_solve_t *enums, pastix_int_t sched, SolverMatrix *solvmtx, pastix_rhs_t rhsb)
Waitall routine for current cblk request.
void cpucblk_zsend_rhs_backward(const SolverMatrix *solvmtx, SolverCblk *cblk, pastix_rhs_t b)
Send the rhs associated to a cblk->lcolidx to the remote node.
void cpucblk_zinit(pastix_coefside_t side, const SolverMatrix *solvmtx, const pastix_bcsc_t *bcsc, pastix_int_t itercblk, const char *directory)
Fully initialize a single cblk.
void cpucblk_zuncompress(pastix_coefside_t side, SolverCblk *cblk)
Uncompress a single column block from low-rank format to full-rank format.
int cpucblk_zincoming_rhs_fwd_deps(int rank, const args_solve_t *enums, SolverMatrix *solvmtx, SolverCblk *cblk, pastix_rhs_t rhsb)
Wait for incoming dependencies, and return when cblk->ctrbcnt has reached 0.
int cpucblk_zpotrfsp1d_potrf(SolverMatrix *solvmtx, SolverCblk *cblk, void *dataL)
Compute the Cholesky factorization of the diagonal block in a panel.
int cpucblk_zgetrfsp1d_getrf(SolverMatrix *solvmtx, SolverCblk *cblk, void *L, void *U)
Compute the LU factorization of the diagonal block in a panel.
int cpucblk_zsytrfsp1d_panel(SolverMatrix *solvmtx, SolverCblk *cblk, void *L, void *DLt)
Compute the LDL^t factorization of one panel.
void cpucblk_zalloc(pastix_coefside_t side, SolverCblk *cblk)
Allocate the cblk structure to store the coefficient.
void cpucblk_zgetschur(const SolverCblk *cblk, int upper_part, pastix_complex64_t *S, pastix_int_t lds)
Extract a cblk panel of the Schur complement to a dense lapack form.
int cpucblk_zdiff(pastix_coefside_t side, const SolverCblk *cblkA, SolverCblk *cblkB)
Compare two column blocks in full-rank format.
Definition: cpucblk_zdiff.c:63
void cpucblk_zalloc_fr(pastix_coefside_t side, SolverCblk *cblk)
Allocate the cblk structure to store the coefficient.
void cpucblk_zadd(pastix_coefside_t side, double alpha, const SolverCblk *cblkA, SolverCblk *cblkB, const pastix_lr_t *lowrank)
Add two column bloks in full rank format.
Definition: cpucblk_zadd.c:379
void cpucblk_zalloc_lr(pastix_coefside_t side, SolverCblk *cblk, int rkmax)
Allocate the cblk structure to store the coefficient.
pastix_fixdbl_t cpucblk_zgemmsp(pastix_coefside_t sideA, pastix_trans_t trans, const SolverCblk *cblk, const SolverBlok *blok, SolverCblk *fcblk, const void *A, const void *B, void *C, pastix_complex64_t *work, pastix_int_t lwork, const pastix_lr_t *lowrank)
Compute the updates associated to one off-diagonal block.
int cpucblk_zincoming_deps(int mt_flag, pastix_coefside_t side, SolverMatrix *solvmtx, SolverCblk *cblk)
Wait for incoming dependencies, and return when cblk->ctrbcnt has reached 0.
void cpucblk_zdump(pastix_coefside_t side, const SolverCblk *cblk, FILE *stream)
Dump a single column block into a FILE in a human readale format.
Definition: coeftab_zinit.c:54
void cpucblk_zrequest_rhs_fwd_cleanup(const args_solve_t *enums, pastix_int_t sched, SolverMatrix *solvmtx, pastix_rhs_t rhsb)
Waitall routine for current cblk request.
int cpucblk_zgeaddsp1d(const SolverCblk *cblk1, SolverCblk *cblk2, const pastix_complex64_t *L1, pastix_complex64_t *L2, const pastix_complex64_t *U1, pastix_complex64_t *U2)
Add two column blocks together.
Definition: core_zgeadd.c:245
int cpucblk_zgetrfsp1d(SolverMatrix *solvmtx, SolverCblk *cblk, pastix_complex64_t *work, pastix_int_t lwork)
Perform the LU factorization of a given panel and apply all its updates.
void cpucblk_zsend_rhs_forward(const SolverMatrix *solvmtx, SolverCblk *cblk, pastix_rhs_t b)
Send the rhs associated to a cblk->lcolidx to the remote node.
pastix_fixdbl_t cpublok_zcompress(const pastix_lr_t *lowrank, pastix_int_t M, pastix_int_t N, pastix_lrblock_t *blok)
Compress a single block from full-rank to low-rank format.
void cpucblk_zrelease_deps(pastix_coefside_t side, SolverMatrix *solvmtx, const SolverCblk *cblk, SolverCblk *fcbk)
Release the dependencies of the given cblk after an update.
int cpucblk_zincoming_rhs_bwd_deps(int rank, const args_solve_t *enums, SolverMatrix *solvmtx, SolverCblk *cblk, pastix_rhs_t rhsb)
Wait for incoming dependencies, and return when cblk->ctrbcnt has reached 0.
int cpucblk_zpotrfsp1d(SolverMatrix *solvmtx, SolverCblk *cblk, pastix_complex64_t *work, pastix_int_t lwork)
Perform the Cholesky factorization of a given panel and apply all its updates.
int cpucblk_zpotrfsp1d_panel(SolverMatrix *solvmtx, SolverCblk *cblk, void *L)
Compute the Cholesky factorization of one panel.
void cpucblk_zfree(pastix_coefside_t side, SolverCblk *cblk)
Free the cblk structure that store the coefficient.
pastix_fixdbl_t cpublok_zgemmsp(pastix_trans_t trans, const SolverCblk *cblk, SolverCblk *fcblk, pastix_int_t blok_mk, pastix_int_t blok_nk, pastix_int_t blok_mn, const void *A, const void *B, void *C, const pastix_lr_t *lowrank)
Compute the CPU gemm associated to a couple of off-diagonal blocks.
pastix_int_t cpucblk_zcompress(const SolverMatrix *solvmtx, pastix_coefside_t side, int max_ilulvl, SolverCblk *cblk)
Compress a single column block from full-rank to low-rank format.
void cpucblk_zrecv_rhs_backward(const SolverMatrix *solvmtx, SolverCblk *cblk, pastix_rhs_t b)
Receive the rhs associated to a cblk->lcolidx to the remote node.
pastix_fixdbl_t cpublok_ztrsmsp(pastix_side_t side, pastix_uplo_t uplo, pastix_trans_t trans, pastix_diag_t diag, const SolverCblk *cblk, pastix_int_t blok_m, const void *A, void *C, const pastix_lr_t *lowrank)
Compute the updates associated to one off-diagonal block.
Definition: core_ztrsmsp.c:690
void cpucblk_zmemory(pastix_coefside_t side, const SolverMatrix *solvmtx, SolverCblk *cblk, pastix_int_t *orig, pastix_int_t *gain)
Return the memory gain of the low-rank form over the full-rank form for a single column-block.
void cpublok_zscalo(pastix_trans_t trans, SolverCblk *cblk, pastix_int_t blok_m, const void *A, const void *dataD, void *dataB)
Copy the lower terms of the block with scaling for the two-terms algorithm.
Definition: core_zscalo.c:315
int cpucblk_zgetrfsp1d_panel(SolverMatrix *solvmtx, SolverCblk *cblk, void *L, void *U)
Compute the LU factorization of one panel.
void cpucblk_zrequest_cleanup(pastix_coefside_t side, pastix_int_t sched, SolverMatrix *solvmtx)
Waitall routine for current cblk request.
void cpucblk_zrelease_rhs_fwd_deps(const args_solve_t *enums, SolverMatrix *solvmtx, pastix_rhs_t rhsb, const SolverCblk *cblk, SolverCblk *fcbk)
Release the dependencies of the given cblk after an update.
int cpucblk_zsytrfsp1d_sytrf(SolverMatrix *solvmtx, SolverCblk *cblk, void *dataL)
Computes the LDL^t factorization of the diagonal block in a panel.
int cpucblk_zsytrfsp1d(SolverMatrix *solvmtx, SolverCblk *cblk, pastix_complex64_t *Dlt, pastix_complex64_t *work, pastix_int_t lwork)
Perform the LDL^t factorization of a given panel and apply all its updates.
void cpucblk_zalloc_lrws(const SolverCblk *cblk, pastix_lrblock_t *lrblok, pastix_complex64_t *ws)
Initialize lrblock structure from a workspace from all blocks of the cblk associated.
Definition: cpucblk_zinit.c:96
void cpucblk_zrecv_rhs_forward(const SolverMatrix *solvmtx, SolverCblk *cblk, pastix_complex64_t *work, pastix_rhs_t b)
Receive the rhs associated to a cblk->lcolidx to the remote node.
void cpucblk_zscalo(pastix_trans_t trans, SolverCblk *cblk, void *dataL, void *dataLD)
Copy the L term with scaling for the two-terms algorithm.
Definition: core_zscalo.c:170
void cpucblk_zfillin(pastix_coefside_t side, const SolverMatrix *solvmtx, const pastix_bcsc_t *bcsc, pastix_int_t itercblk)
Initialize the coeftab structure from the internal bcsc.
Structure to define the type of function to use for the low-rank kernels and their parameters.
The block low-rank structure to hold a matrix in low-rank form.
void solve_blok_zgemm(pastix_side_t side, pastix_trans_t trans, pastix_int_t nrhs, const SolverCblk *cblk, const SolverBlok *blok, SolverCblk *fcbk, const void *dataA, const pastix_complex64_t *B, pastix_int_t ldb, pastix_complex64_t *C, pastix_int_t ldc)
Apply a solve gemm update related to a single block of the matrix A.
void solve_cblk_zdiag(const SolverCblk *cblk, int nrhs, pastix_complex64_t *b, int ldb, pastix_complex64_t *work)
Apply the diagonal solve related to one cblk to all the right hand side.
void solve_cblk_ztrsmsp_forward(const args_solve_t *enums, SolverMatrix *datacode, const SolverCblk *cblk, pastix_rhs_t b)
Apply a forward solve related to one cblk to all the right hand side.
void solve_blok_ztrsm(pastix_side_t side, pastix_uplo_t uplo, pastix_trans_t trans, pastix_diag_t diag, const SolverCblk *cblk, int nrhs, const void *dataA, pastix_complex64_t *b, int ldb)
Apply a solve trsm update related to a diagonal block of the matrix A.
Definition: solve_ztrsmsp.c:75
void solve_cblk_ztrsmsp_backward(const args_solve_t *enums, SolverMatrix *datacode, SolverCblk *cblk, pastix_rhs_t b)
Apply a backward solve related to one cblk to all the right hand side.
enum pastix_diag_e pastix_diag_t
Diagonal.
enum pastix_uplo_e pastix_uplo_t
Upper/Lower part.
enum pastix_side_e pastix_side_t
Side of the operation.
enum pastix_trans_e pastix_trans_t
Transpostion.
enum pastix_coefside_e pastix_coefside_t
Data blocks used in the kernel.
Arguments for the solve.
Definition: solver.h:85
Solver block structure.
Definition: solver.h:137
Solver column block structure.
Definition: solver.h:156