SuperLU Distributed 8.2.1
Distributed memory sparse direct solver
superlu_ddefs.h
Go to the documentation of this file.
1
28#ifndef __SUPERLU_DDEFS /* allow multiple inclusions */
29#define __SUPERLU_DDEFS
30
31/*
32 * File name: superlu_ddefs.h
33 * Purpose: Distributed SuperLU data types and function prototypes
34 * History:
35 */
36
37#include "superlu_defs.h"
38
39/*
40 *-- The structure used to store matrix A of the linear system and
41 * several vectors describing the transformations done to matrix A.
42 *
43 * A (SuperMatrix*)
44 * Matrix A in A*X=B, of dimension (A->nrow, A->ncol).
45 * The number of linear equations is A->nrow. The type of A can be:
46 * Stype = SLU_NC; Dtype = SLU_D; Mtype = SLU_GE.
47 *
48 * DiagScale (DiagScale_t)
49 * Specifies the form of equilibration that was done.
50 * = NOEQUIL: No equilibration.
51 * = ROW: Row equilibration, i.e., A was premultiplied by diag(R).
52 * = COL: Column equilibration, i.e., A was postmultiplied by diag(C).
53 * = BOTH: Both row and column equilibration, i.e., A was replaced
54 * by diag(R)*A*diag(C).
55 *
56 * R double*, dimension (A->nrow)
57 * The row scale factors for A.
58 * If DiagScale = ROW or BOTH, A is multiplied on the left by diag(R).
59 * If DiagScale = NOEQUIL or COL, R is not defined.
60 *
61 * C double*, dimension (A->ncol)
62 * The column scale factors for A.
63 * If DiagScale = COL or BOTH, A is multiplied on the right by diag(C).
64 * If DiagScale = NOEQUIL or ROW, C is not defined.
65 *
66 * perm_r (int*) dimension (A->nrow)
67 * Row permutation vector which defines the permutation matrix Pr,
68 * perm_r[i] = j means row i of A is in position j in Pr*A.
69 *
70 * perm_c (int*) dimension (A->ncol)
71 * Column permutation vector, which defines the
72 * permutation matrix Pc; perm_c[i] = j means column i of A is
73 * in position j in A*Pc.
74 *
75 */
76typedef struct {
78 double *R;
79 double *C;
83
84#if 0 // Sherry: move to superlu_defs.h
85/*-- Auxiliary data type used in PxGSTRS/PxGSTRS1. */
86typedef struct {
87 int_t lbnum; /* Row block number (local). */
88 int_t indpos; /* Starting position in Uindex[]. */
90#endif
91
92/*
93 * On each processor, the blocks in L are stored in compressed block
94 * column format, the blocks in U are stored in compressed block row format.
95 */
96#define MAX_LOOKAHEADS 50
97typedef struct {
98 int_t **Lrowind_bc_ptr; /* size ceil(NSUPERS/Pc) */
99 int_t *Lrowind_bc_dat; /* size sum of sizes of Lrowind_bc_ptr[lk]) */
100 long int *Lrowind_bc_offset; /* size ceil(NSUPERS/Pc) */
102
103 double **Lnzval_bc_ptr; /* size ceil(NSUPERS/Pc) */
104 double *Lnzval_bc_dat; /* size sum of sizes of Lnzval_bc_ptr[lk]) */
105 long int *Lnzval_bc_offset; /* size ceil(NSUPERS/Pc) */
107
108 double **Linv_bc_ptr; /* size ceil(NSUPERS/Pc) */
109 double *Linv_bc_dat; /* size sum of sizes of Linv_bc_ptr[lk]) */
110 long int *Linv_bc_offset; /* size ceil(NSUPERS/Pc) */
111 long int Linv_bc_cnt;
112
113 int_t **Lindval_loc_bc_ptr; /* size ceil(NSUPERS/Pc) pointers to locations in Lrowind_bc_ptr and Lnzval_bc_ptr */
114 int_t *Lindval_loc_bc_dat; /* size sum of sizes of Lindval_loc_bc_ptr[lk]) */
115 long int *Lindval_loc_bc_offset; /* size ceil(NSUPERS/Pc) */
117 int_t *Unnz; /* number of nonzeros per block column in U*/
118 int_t **Lrowind_bc_2_lsum; /* size ceil(NSUPERS/Pc) map indices of Lrowind_bc_ptr to indices of lsum */
119 double **Uinv_bc_ptr; /* size ceil(NSUPERS/Pc) */
120 double *Uinv_bc_dat; /* size sum of sizes of Linv_bc_ptr[lk]) */
121 long int *Uinv_bc_offset; /* size ceil(NSUPERS/Pc) */
122 long int Uinv_bc_cnt;
123
124 int_t **Ufstnz_br_ptr; /* size ceil(NSUPERS/Pr) */
125 int_t *Ufstnz_br_dat; /* size sum of sizes of Ufstnz_br_ptr[lk]) */
126 long int *Ufstnz_br_offset; /* size ceil(NSUPERS/Pr) */
128
129 double **Unzval_br_ptr; /* size ceil(NSUPERS/Pr) */
130 double *Unzval_br_dat; /* size sum of sizes of Unzval_br_ptr[lk]) */
131 long int *Unzval_br_offset; /* size ceil(NSUPERS/Pr) */
133
134 /*-- Data structures used for broadcast and reduction trees. --*/
135 C_Tree *LBtree_ptr; /* size ceil(NSUPERS/Pc) */
136 C_Tree *LRtree_ptr; /* size ceil(NSUPERS/Pr) */
137 C_Tree *UBtree_ptr; /* size ceil(NSUPERS/Pc) */
138 C_Tree *URtree_ptr; /* size ceil(NSUPERS/Pr) */
139#if 0
140 int_t *Lsub_buf; /* Buffer for the remote subscripts of L */
141 double *Lval_buf; /* Buffer for the remote nonzeros of L */
142 int_t *Usub_buf; /* Buffer for the remote subscripts of U */
143 double *Uval_buf; /* Buffer for the remote nonzeros of U */
144#endif
145 int_t *Lsub_buf_2[MAX_LOOKAHEADS]; /* Buffers for the remote subscripts of L*/
146 double *Lval_buf_2[MAX_LOOKAHEADS]; /* Buffers for the remote nonzeros of L */
147 int_t *Usub_buf_2[MAX_LOOKAHEADS]; /* Buffer for the remote subscripts of U */
148 double *Uval_buf_2[MAX_LOOKAHEADS]; /* Buffer for the remote nonzeros of U */
149 double *ujrow; /* used in panel factorization. */
150 int_t bufmax[NBUFFERS]; /* Maximum buffer size across all MPI ranks:
151 * 0 : maximum size of Lsub_buf[]
152 * 1 : maximum size of Lval_buf[]
153 * 2 : maximum size of Usub_buf[]
154 * 3 : maximum size of Uval_buf[]
155 * 4 : maximum size of tempv[LDA]
156 */
157
158 /*-- Record communication schedule for factorization. --*/
159 int *ToRecv; /* Recv from no one (0), left (1), and up (2).*/
160 int *ToSendD; /* Whether need to send down block row. */
161 int **ToSendR; /* List of processes to send right block col. */
162
163 /*-- Record communication schedule for forward/back solves. --*/
164 /* 1/15/22 Sherry: changed int_t to int type */
165 int *fmod; /* Modification count for L-solve */
166 int **fsendx_plist; /* Column process list to send down Xk */
167 int *frecv; /* Modifications to be recv'd in proc row */
168 int nfrecvx; /* Number of Xk I will receive in L-solve */
169 int nfsendx; /* Number of Xk I will send in L-solve */
170 int *bmod; /* Modification count for U-solve */
171 int **bsendx_plist; /* Column process list to send down Xk */
172 int *brecv; /* Modifications to be recv'd in proc row */
173 int nbrecvx; /* Number of Xk I will receive in U-solve */
174 int nbsendx; /* Number of Xk I will send in U-solve */
175 int *mod_bit; /* Flag contribution from each row blocks */
176
177 /*-- Auxiliary arrays used for forward/back solves. --*/
178 int_t *ilsum; /* Starting position of each supernode in lsum
179 (local) */
180 int_t ldalsum; /* LDA of lsum (local) */
181 int_t SolveMsgSent; /* Number of actual messages sent in LU-solve */
182 int_t SolveMsgVol; /* Volume of messages sent in the solve phase */
183
184
185 /*********************/
186 /* The following variables are used in the hybrid solver */
187
188 /*-- Counts to be used in U^{-T} triangular solve. -- */
192 int_t ut_ldalsum; /* LDA of lsum (local) */
193 int_t *ut_ilsum; /* ilsum in column-wise */
194 int_t *utmod; /* Modification count for Ut-solve. */
195 int_t **ut_sendx_plist; /* Row process list to send down Xk */
196 int_t *utrecv; /* Modifications to be recev'd in proc column. */
197 int_t n_utsendx; /* Number of Xk I will receive */
198 int_t n_utrecvx; /* Number of Xk I will send */
203 Ucb_indptr_t **Ucb_indptr;/* Vertical linked list pointing to Uindex[] */
205 long int *Ucb_indoffset;
206 long int Ucb_indcnt;
207
208 int_t **Ucb_valptr; /* Vertical linked list pointing to Unzval[] */
210 long int *Ucb_valoffset;
211 long int Ucb_valcnt;
212
213 /* some additional counters for L solve */
217 int_t inv; /* whether the diagonal block is inverted*/
218
219 /* The following variables are used in GPU trisolve*/
220#ifdef GPU_ACC
221 int_t *d_Lrowind_bc_dat;
222 long int *d_Lrowind_bc_offset;
223 double *d_Lnzval_bc_dat;
224 long int *d_Lnzval_bc_offset;
225 double *d_Linv_bc_dat ;
226 double *d_Uinv_bc_dat ;
227 long int *d_Linv_bc_offset ;
228 long int *d_Uinv_bc_offset ;
229 int_t *d_Lindval_loc_bc_dat ;
230 long int *d_Lindval_loc_bc_offset ;
231
232 int_t *d_Urbs;
233 int_t *d_Ufstnz_br_dat;
234 long int *d_Ufstnz_br_offset;
235 double *d_Unzval_br_dat;
236 long int *d_Unzval_br_offset;
237
238 int_t *d_Ucb_valdat;
239 long int *d_Ucb_valoffset;
240 Ucb_indptr_t *d_Ucb_inddat;
241 long int *d_Ucb_indoffset;
242
243 int_t *d_ilsum ;
244 int_t *d_xsup ;
245 C_Tree *d_LBtree_ptr ;
246 C_Tree *d_LRtree_ptr ;
247 C_Tree *d_UBtree_ptr ;
248 C_Tree *d_URtree_ptr ;
249#endif
250
251} dLocalLU_t;
252
253
254typedef struct {
258 char dt;
260
261
262/*-- Data structure for communication during matrix-vector multiplication. */
263typedef struct {
265 int_t *ind_tosend; /* X indeices to be sent to other processes */
266 int_t *ind_torecv; /* X indeices to be received from other processes */
267 int_t *ptr_ind_tosend;/* Printers to ind_tosend[] (Size procs)
268 (also point to val_torecv) */
269 int_t *ptr_ind_torecv;/* Printers to ind_torecv[] (Size procs)
270 (also point to val_tosend) */
271 int *SendCounts; /* Numbers of X indices to be sent
272 (also numbers of X values to be received) */
273 int *RecvCounts; /* Numbers of X indices to be received
274 (also numbers of X values to be sent) */
275 double *val_tosend; /* X values to be sent to other processes */
276 double *val_torecv; /* X values to be received from other processes */
277 int_t TotalIndSend; /* Total number of indices to be sent
278 (also total number of values to be received) */
279 int_t TotalValSend; /* Total number of values to be sent.
280 (also total number of indices to be received) */
282
283/*-- Data structure holding the information for the solution phase --*/
284typedef struct {
287 int_t num_diag_procs, *diag_procs, *diag_len;
288 pdgsmv_comm_t *gsmv_comm; /* communication metadata for SpMV,
289 required by IterRefine. */
290 pxgstrs_comm_t *gstrs_comm; /* communication metadata for SpTRSV. */
291 int_t *A_colind_gsmv; /* After pdgsmv_init(), the global column
292 indices of A are translated into the relative
293 positions in the gathered x-vector.
294 This is re-used in repeated calls to pdgsmv() */
295 int_t *xrow_to_proc; /* used by PDSLin */
296 NRformat_loc3d* A3d; /* Point to 3D {A, B} gathered on 2D layer 0.
297 This needs to be peresistent between
298 3D factorization and solve. */
300
301
302
303/*==== For 3D code ====*/
304
305// new structures for pdgstrf_4_8
306
307typedef struct
308{
313 double* uval;
315
316typedef struct
317{
319 double *lusup;
321 int_t nlb; //number of l blocks
324
325
326
327/* HyP_t is the data structure to assist HALO offload of Schur-complement. */
328typedef struct
329{
331 Ublock_info_t *Ublock_info, *Ublock_info_Phi;
332
333 int_t first_l_block_acc , first_u_block_acc;
335 int_t *Lblock_dirty_bit, * Ublock_dirty_bit;
336 double *lookAhead_L_buff, *Remain_L_buff;
337 int_t lookAheadBlk; /* number of blocks in look-ahead window */
338 int_t RemainBlk ; /* number of blocks outside look-ahead window */
339 int_t num_look_aheads, nsupers;
340 int_t ldu, ldu_Phi;
341 int_t num_u_blks, num_u_blks_Phi;
342
344 double *bigU_Phi;
345 double *bigU_host;
348
354} HyP_t;
355
356typedef struct
357{
359 double * Lval_buf ;
361 double * Uval_buf ;
363
365 int_t knsupc,
366 HyP_t* HyP,
367 SCT_t* SCT,
368 SuperLUStat_t *stat
369 );
370
371typedef struct
372{
383
384typedef struct
385{
386 double *bigU;
387 double *bigV;
388} dscuBufs_t;
389
390typedef struct
391{
395
396typedef struct
397{
403
404//#endif
405/*=====================*/
406
407/***********************************************************************
408 * Function prototypes
409 ***********************************************************************/
410
411#ifdef __cplusplus
412extern "C" {
413#endif
414
415
416/* Supernodal LU factor related */
417extern void
420extern void
422 int_t, double *, int_t *, int_t *,
424extern void
426 double **, int_t **, int_t **);
427extern int
429 SuperMatrix *);
430extern void
432extern void
435extern void
437 int_t *, int_t *, int_t *, int_t *, int_t *,
439extern void
441 double *, int_t);
442
443extern void dallocateA_dist (int_t, int_t, double **, int_t **, int_t **);
444extern void dGenXtrue_dist (int_t, int_t, double *, int_t);
445extern void dFillRHS_dist (char *, int_t, double *, int_t,
446 SuperMatrix *, double *, int_t);
447extern int dcreate_matrix(SuperMatrix *, int, double **, int *,
448 double **, int *, FILE *, gridinfo_t *);
449extern int dcreate_matrix_rb(SuperMatrix *, int, double **, int *,
450 double **, int *, FILE *, gridinfo_t *);
451extern int dcreate_matrix_dat(SuperMatrix *, int, double **, int *,
452 double **, int *, FILE *, gridinfo_t *);
453extern int dcreate_matrix_postfix(SuperMatrix *, int, double **, int *,
454 double **, int *, FILE *, char *, gridinfo_t *);
455
456extern void dScalePermstructInit(const int_t, const int_t,
459
460/* Driver related */
461extern void dgsequ_dist (SuperMatrix *, double *, double *, double *,
462 double *, double *, int_t *);
463extern double dlangs_dist (char *, SuperMatrix *);
464extern void dlaqgs_dist (SuperMatrix *, double *, double *, double,
465 double, double, char *);
466extern void pdgsequ (SuperMatrix *, double *, double *, double *,
467 double *, double *, int_t *, gridinfo_t *);
468extern double pdlangs (char *, SuperMatrix *, gridinfo_t *);
469extern void pdlaqgs (SuperMatrix *, double *, double *, double,
470 double, double, char *);
471extern int pdPermute_Dense_Matrix(int_t, int_t, int_t [], int_t[],
472 double [], int, double [], int, int,
473 gridinfo_t *);
474
475extern int sp_dtrsv_dist (char *, char *, char *, SuperMatrix *,
476 SuperMatrix *, double *, int *);
477extern int sp_dgemv_dist (char *, double, SuperMatrix *, double *,
478 int, double, double *, int);
479extern int sp_dgemm_dist (char *, int, double, SuperMatrix *,
480 double *, int, double, double *, int);
481
485 dScalePermstruct_t *, double *,
486 int, int, gridinfo_t *, dLUstruct_t *, double *,
487 SuperLUStat_t *, int *);
488extern float pddistribute(fact_t, int_t, SuperMatrix *,
492 dScalePermstruct_t *, double *,
493 int, int, gridinfo_t *, dLUstruct_t *,
494 dSOLVEstruct_t *, double *, SuperLUStat_t *, int *);
501 int_t [], int_t [], gridinfo_t *grid,
503extern void pxgstrs_finalize(pxgstrs_comm_t *);
504extern int dldperm_dist(int, int, int_t, int_t [], int_t [],
505 double [], int_t *, double [], double []);
506extern int dstatic_schedule(superlu_dist_options_t *, int, int,
508 int_t *, int_t *, int *);
509extern void dLUstructInit(const int_t, dLUstruct_t *);
510extern void dLUstructFree(dLUstruct_t *);
511extern void dDestroy_LU(int_t, gridinfo_t *, dLUstruct_t *);
512extern void dDestroy_Tree(int_t, gridinfo_t *, dLUstruct_t *);
513extern void dscatter_l (int ib, int ljb, int nsupc, int_t iukp, int_t* xsup,
514 int klst, int nbrow, int_t lptr, int temp_nbrow,
515 int_t* usub, int_t* lsub, double *tempv,
516 int* indirect_thread, int* indirect2,
517 int_t ** Lrowind_bc_ptr, double **Lnzval_bc_ptr,
518 gridinfo_t * grid);
519extern void dscatter_u (int ib, int jb, int nsupc, int_t iukp, int_t * xsup,
520 int klst, int nbrow, int_t lptr, int temp_nbrow,
521 int_t* lsub, int_t* usub, double* tempv,
522 int_t ** Ufstnz_br_ptr, double **Unzval_br_ptr,
523 gridinfo_t * grid);
524extern int_t pdgstrf(superlu_dist_options_t *, int, int, double anorm,
526
527/* #define GPU_PROF
528#define IPM_PROF */
529
530/* Solve related */
532 double *, int_t, int, SuperLUStat_t *, int *);
534 double *, int_t, int_t, int_t, int, dSOLVEstruct_t *,
535 SuperLUStat_t *, int *);
536extern void pdgstrf2_trsm(superlu_dist_options_t * options, int_t k0, int_t k,
537 double thresh, Glu_persist_t *, gridinfo_t *,
538 dLocalLU_t *, MPI_Request *, int tag_ub,
539 SuperLUStat_t *, int *info);
540extern void pdgstrs2_omp(int_t k0, int_t k, Glu_persist_t *, gridinfo_t *,
542extern int_t pdReDistribute_B_to_X(double *B, int_t m_loc, int nrhs, int_t ldb,
543 int_t fst_row, int_t *ilsum, double *x,
546extern void dlsum_fmod(double *, double *, double *, double *,
547 int, int, int_t , int *fmod, int_t, int_t, int_t,
549 MPI_Request [], SuperLUStat_t *);
550extern void dlsum_bmod(double *, double *, double *,
551 int, int_t, int *bmod, int_t *, Ucb_indptr_t **,
552 int_t **, int_t *, gridinfo_t *, dLocalLU_t *,
553 MPI_Request [], SuperLUStat_t *);
554
555extern void dlsum_fmod_inv(double *, double *, double *, double *,
556 int, int_t , int *fmod,
558 SuperLUStat_t **, int_t *, int_t *, int_t, int_t, int_t, int_t, int, int);
559extern void dlsum_fmod_inv_master(double *, double *, double *, double *,
560 int, int, int_t , int *fmod, int_t,
562 SuperLUStat_t **, int_t, int_t, int_t, int_t, int, int);
563extern void dlsum_bmod_inv(double *, double *, double *, double *,
564 int, int_t, int *bmod, int_t *, Ucb_indptr_t **,
565 int_t **, int_t *, gridinfo_t *, dLocalLU_t *,
566 SuperLUStat_t **, int_t *, int_t *, int_t, int_t, int, int);
567extern void dlsum_bmod_inv_master(double *, double *, double *, double *,
568 int, int_t, int *bmod, int_t *, Ucb_indptr_t **,
569 int_t **, int_t *, gridinfo_t *, dLocalLU_t *,
570 SuperLUStat_t **, int_t, int_t, int, int);
571
572extern void dComputeLevelsets(int , int_t , gridinfo_t *,
574
575#ifdef GPU_ACC
576extern void dlsum_fmod_inv_gpu_wrap(int_t, int_t, int_t, int_t, double *, double *, int, int, int_t , int *fmod, C_Tree *, C_Tree *, int_t *, int_t *, int64_t *, double *, int64_t *, double *, int64_t *, int_t *, int64_t *, int_t *, gridinfo_t *, double * , double * , int_t );
577extern void dlsum_bmod_inv_gpu_wrap(int_t, int_t, int_t, int_t, double *, double *,int,int, int_t , int *bmod, C_Tree *, C_Tree *, int_t *, int_t *,int_t *, int64_t *, double *, int64_t *, int_t *, int64_t *, Ucb_indptr_t *, int64_t *, double *, int64_t *,int_t *,gridinfo_t *);
578#endif
579
580extern void pdgsrfs(int_t, SuperMatrix *, double, dLUstruct_t *,
582 double [], int_t, double [], int_t, int,
583 dSOLVEstruct_t *, double *, SuperLUStat_t *, int *);
584extern void pdgsrfs_ABXglobal(int_t, SuperMatrix *, double, dLUstruct_t *,
585 gridinfo_t *, double *, int_t, double *, int_t,
586 int, double *, SuperLUStat_t *, int *);
588 gridinfo_t *, int_t *, int_t *[],
589 double *[], int_t *[], int_t []);
590extern int pdgsmv_AXglobal(int_t, int_t [], double [], int_t [],
591 double [], double []);
592extern int pdgsmv_AXglobal_abs(int_t, int_t [], double [], int_t [],
593 double [], double []);
594extern void pdgsmv_init(SuperMatrix *, int_t *, gridinfo_t *,
595 pdgsmv_comm_t *);
596extern void pdgsmv(int_t, SuperMatrix *, gridinfo_t *, pdgsmv_comm_t *,
597 double x[], double ax[]);
598extern void pdgsmv_finalize(pdgsmv_comm_t *);
599
600/* Memory-related */
601extern double *doubleMalloc_dist(int_t);
602extern double *doubleCalloc_dist(int_t);
603extern void *duser_malloc_dist (int_t, int_t);
604extern void duser_free_dist (int_t, int_t);
607
608/* Auxiliary routines */
609
615extern void dZeroLblocks(int, int, gridinfo_t *, dLUstruct_t *);
616extern void dZeroUblocks(int iam, int n, gridinfo_t *, dLUstruct_t *);
617extern void dfill_dist (double *, int_t, double);
618extern void dinf_norm_error_dist (int_t, int_t, double*, int_t,
619 double*, int_t, gridinfo_t*);
620extern void pdinf_norm_error(int, int_t, int_t, double [], int_t,
621 double [], int_t , MPI_Comm);
622extern void dreadhb_dist (int, FILE *, int_t *, int_t *, int_t *,
623 double **, int_t **, int_t **);
624extern void dreadtriple_dist(FILE *, int_t *, int_t *, int_t *,
625 double **, int_t **, int_t **);
626extern void dreadtriple_noheader(FILE *, int_t *, int_t *, int_t *,
627 double **, int_t **, int_t **);
628extern void dreadrb_dist(int, FILE *, int_t *, int_t *, int_t *,
629 double **, int_t **, int_t **);
630extern void dreadMM_dist(FILE *, int_t *, int_t *, int_t *,
631 double **, int_t **, int_t **);
632extern int dread_binary(FILE *, int_t *, int_t *, int_t *,
633 double **, int_t **, int_t **);
634
635/* Distribute the data for numerical factorization */
639extern void pdGetDiagU(int_t, dLUstruct_t *, gridinfo_t *, double *);
640
642
643/* Routines for debugging */
644extern void dPrintLblocks(int, int_t, gridinfo_t *, Glu_persist_t *,
645 dLocalLU_t *);
646extern void dPrintUblocks(int, int_t, gridinfo_t *, Glu_persist_t *,
647 dLocalLU_t *);
652extern void Printdouble5(char *, int_t, double *);
653extern int file_Printdouble5(FILE *, char *, int_t, double *);
654
655extern void dGenCOOLblocks(int, int_t, gridinfo_t*,
656 Glu_persist_t*, dLocalLU_t *, int_t** , int_t** , double ** , int_t* , int_t* );
657extern void dGenCSCLblocks(int, int_t, gridinfo_t*,
658 Glu_persist_t*, dLocalLU_t *, double **, int_t **, int_t **, int_t*, int_t*);
659extern void dGenCSRLblocks(int, int_t, gridinfo_t*,
660 Glu_persist_t*, dLocalLU_t *, double **, int_t **, int_t **, int_t*, int_t*);
661
662
663/* BLAS */
664
665#ifdef USE_VENDOR_BLAS
666extern void dgemm_(const char*, const char*, const int*, const int*, const int*,
667 const double*, const double*, const int*, const double*,
668 const int*, const double*, double*, const int*, int, int);
669extern void dtrsv_(char*, char*, char*, int*, double*, int*,
670 double*, int*, int, int, int);
671extern void dtrsm_(const char*, const char*, const char*, const char*,
672 const int*, const int*, const double*, const double*, const int*,
673 double*, const int*, int, int, int, int);
674extern void dgemv_(const char *, const int *, const int *, const double *,
675 const double *a, const int *, const double *, const int *,
676 const double *, double *, const int *, int);
677
678#else
679extern int dgemm_(const char*, const char*, const int*, const int*, const int*,
680 const double*, const double*, const int*, const double*,
681 const int*, const double*, double*, const int*);
682extern int dtrsv_(char*, char*, char*, int*, double*, int*,
683 double*, int*);
684extern int dtrsm_(const char*, const char*, const char*, const char*,
685 const int*, const int*, const double*, const double*, const int*,
686 double*, const int*);
687extern void dgemv_(const char *, const int *, const int *, const double *,
688 const double *a, const int *, const double *, const int *,
689 const double *, double *, const int *);
690#endif
691
692extern void dger_(const int*, const int*, const double*,
693 const double*, const int*, const double*, const int*,
694 double*, const int*);
695
696extern int dscal_(const int *n, const double *alpha, double *dx, const int *incx);
697extern int daxpy_(const int *n, const double *alpha, const double *x,
698 const int *incx, double *y, const int *incy);
699
700/* SuperLU BLAS interface: dsuperlu_blas.c */
701extern int superlu_dgemm(const char *transa, const char *transb,
702 int m, int n, int k, double alpha, double *a,
703 int lda, double *b, int ldb, double beta, double *c, int ldc);
704extern int superlu_dtrsm(const char *sideRL, const char *uplo,
705 const char *transa, const char *diag, const int m, const int n,
706 const double alpha, const double *a,
707 const int lda, double *b, const int ldb);
708extern int superlu_dger(const int m, const int n, const double alpha,
709 const double *x, const int incx, const double *y,
710 const int incy, double *a, const int lda);
711extern int superlu_dscal(const int n, const double alpha, double *x, const int incx);
712extern int superlu_daxpy(const int n, const double alpha,
713 const double *x, const int incx, double *y, const int incy);
714extern int superlu_dgemv(const char *trans, const int m,
715 const int n, const double alpha, const double *a,
716 const int lda, const double *x, const int incx,
717 const double beta, double *y, const int incy);
718extern int superlu_dtrsv(char *uplo, char *trans, char *diag,
719 int n, double *a, int lda, double *x, int incx);
720
721#ifdef SLU_HAVE_LAPACK
722extern void dtrtri_(char*, char*, int*, double*, int*, int*);
723#endif
724
725/*==== For 3D code ====*/
726extern int dcreate_matrix3d(SuperMatrix *A, int nrhs, double **rhs,
727 int *ldb, double **x, int *ldx,
728 FILE *fp, gridinfo3d_t *grid3d);
729extern int dcreate_matrix_postfix3d(SuperMatrix *A, int nrhs, double **rhs,
730 int *ldb, double **x, int *ldx,
731 FILE *fp, char * postfix, gridinfo3d_t *grid3d);
732
733/* Matrix distributed in NRformat_loc in 3D process grid. It converts
734 it to a NRformat_loc distributed in 2D grid in grid-0 */
735extern void dGatherNRformat_loc3d(fact_t Fact, NRformat_loc *A, double *B,
736 int ldb, int nrhs, gridinfo3d_t *grid3d,
737 NRformat_loc3d **);
738extern int dScatter_B3d(NRformat_loc3d *A3d, gridinfo3d_t *grid3d);
739
741 dScalePermstruct_t *, double B[], int ldb, int nrhs,
743 double *berr, SuperLUStat_t *, int *info);
744extern int_t pdgstrf3d(superlu_dist_options_t *, int m, int n, double anorm,
746 gridinfo3d_t *, SuperLUStat_t *, int *);
747extern void dInit_HyP(HyP_t* HyP, dLocalLU_t *Llu, int_t mcb, int_t mrb );
748extern void Free_HyP(HyP_t* HyP);
749extern int updateDirtyBit(int_t k0, HyP_t* HyP, gridinfo_t* grid);
750
751 /* from scatter.h */
752extern void
754 Remain_info_t *Remain_info, double *L_mat, int ldl,
755 double *U_mat, int ldu, double *bigV,
756 // int_t jj0,
757 int_t knsupc, int_t klst,
758 int_t *lsub, int_t *usub, int_t ldt,
759 int_t thread_id,
760 int *indirect, int *indirect2,
761 int_t **Lrowind_bc_ptr, double **Lnzval_bc_ptr,
762 int_t **Ufstnz_br_ptr, double **Unzval_br_ptr,
763 int_t *xsup, gridinfo_t *, SuperLUStat_t *
764#ifdef SCATTER_PROFILE
765 , double *Host_TheadScatterMOP, double *Host_TheadScatterTimer
766#endif
767 );
768
769#ifdef _OPENMP
770/*this version uses a lock to prevent multiple thread updating the same block*/
771extern void
772dblock_gemm_scatter_lock( int_t lb, int_t j, omp_lock_t* lock,
773 Ublock_info_t *Ublock_info, Remain_info_t *Remain_info,
774 double *L_mat, int_t ldl, double *U_mat, int_t ldu,
775 double *bigV,
776 // int_t jj0,
777 int_t knsupc, int_t klst,
778 int_t *lsub, int_t *usub, int_t ldt,
779 int_t thread_id,
780 int *indirect, int *indirect2,
781 int_t **Lrowind_bc_ptr, double **Lnzval_bc_ptr,
782 int_t **Ufstnz_br_ptr, double **Unzval_br_ptr,
783 int_t *xsup, gridinfo_t *
784#ifdef SCATTER_PROFILE
785 , double *Host_TheadScatterMOP, double *Host_TheadScatterTimer
786#endif
787 );
788#endif
789
790extern int_t
792 int_t knsupc, int_t klst, int_t* lsub,
793 int_t * usub, int_t ldt,
794 int* indirect, int* indirect2,
795 HyP_t* HyP, dLUstruct_t *, gridinfo_t*,
796 SCT_t*SCT, SuperLUStat_t *
797 );
798extern int_t
800 int_t knsupc, int_t klst, int_t* lsub,
801 int_t * usub, int_t ldt,
802 int* indirect, int* indirect2,
803 HyP_t* HyP, dLUstruct_t *, gridinfo_t*,
804 SCT_t*SCT, SuperLUStat_t * );
805extern int_t
807 int_t knsupc, int_t klst, int_t* lsub,
808 int_t * usub, int_t ldt,
809 int* indirect, int* indirect2,
810 HyP_t* HyP, dLUstruct_t *, gridinfo_t*,
811 SCT_t*SCT, SuperLUStat_t * );
812extern int_t
814 int_t knsupc, int_t klst, int_t* lsub,
815 int_t * usub, int_t ldt,
816 int* indirect, int* indirect2,
817 HyP_t* HyP, dLUstruct_t *, gridinfo_t*,
818 SCT_t*SCT, SuperLUStat_t * );
819
820 /* from gather.h */
821extern void dgather_u(int_t num_u_blks,
822 Ublock_info_t *Ublock_info, int_t * usub,
823 double *uval, double *bigU, int_t ldu,
824 int_t *xsup, int_t klst /* for SuperSize */
825 );
826
827extern void dgather_l( int_t num_LBlk, int_t knsupc,
828 Remain_info_t *L_info,
829 double * lval, int_t LD_lval,
830 double * L_buff );
831
832extern void dRgather_L(int_t k, int_t *lsub, double *lusup, gEtreeInfo_t*,
834 int_t *myIperm, int_t *iperm_c_supno );
835extern void dRgather_U(int_t k, int_t jj0, int_t *usub, double *uval,
836 double *bigU, gEtreeInfo_t*, Glu_persist_t *,
837 gridinfo_t *, HyP_t *, int_t *myIperm,
838 int_t *iperm_c_supno, int_t *perm_u);
839
840 /* from xtrf3Dpartition.h */
842 superlu_dist_options_t *options,
843 dLUstruct_t *LUstruct, gridinfo3d_t * grid3d);
844extern void dDestroy_trf3Dpartition(trf3Dpartition_t *trf3Dpartition, gridinfo3d_t *grid3d);
845
846extern void d3D_printMemUse(trf3Dpartition_t* trf3Dpartition,
847 dLUstruct_t *LUstruct, gridinfo3d_t * grid3d);
848
849//extern int* getLastDep(gridinfo_t *grid, SuperLUStat_t *stat,
850// superlu_dist_options_t *options, dLocalLU_t *Llu,
851// int_t* xsup, int_t num_look_aheads, int_t nsupers,
852// int_t * iperm_c_supno);
853
854extern void dinit3DLUstructForest( int_t* myTreeIdxs, int_t* myZeroTrIdxs,
855 sForest_t** sForests, dLUstruct_t* LUstruct,
856 gridinfo3d_t* grid3d);
857
858extern int_t dgatherAllFactoredLUFr(int_t* myZeroTrIdxs, sForest_t* sForests,
859 dLUstruct_t* LUstruct, gridinfo3d_t* grid3d,
860 SCT_t* SCT );
861
862 /* The following are from pdgstrf2.h */
863extern int_t dLpanelUpdate(int_t off0, int_t nsupc, double* ublk_ptr,
864 int_t ld_ujrow, double* lusup, int_t nsupr, SCT_t*);
865extern void Local_Dgstrf2(superlu_dist_options_t *options, int_t k,
866 double thresh, double *BlockUFactor, Glu_persist_t *,
868 SuperLUStat_t *, int *info, SCT_t*);
869extern int_t dTrs2_GatherU(int_t iukp, int_t rukp, int_t klst,
870 int_t nsupc, int_t ldu, int_t *usub,
871 double* uval, double *tempv);
872extern int_t dTrs2_ScatterU(int_t iukp, int_t rukp, int_t klst,
873 int_t nsupc, int_t ldu, int_t *usub,
874 double* uval, double *tempv);
875extern int_t dTrs2_GatherTrsmScatter(int_t klst, int_t iukp, int_t rukp,
876 int_t *usub, double* uval, double *tempv,
877 int_t knsupc, int nsupr, double* lusup,
878 Glu_persist_t *Glu_persist) ;
879extern void pdgstrs2
880#ifdef _CRAY
881(
882 int_t m, int_t k0, int_t k, Glu_persist_t *Glu_persist, gridinfo_t *grid,
883 dLocalLU_t *Llu, SuperLUStat_t *stat, _fcd ftcs1, _fcd ftcs2, _fcd ftcs3
884);
885#else
886(
887 int_t m, int_t k0, int_t k, Glu_persist_t *Glu_persist, gridinfo_t *grid,
888 dLocalLU_t *Llu, SuperLUStat_t *stat
889);
890#endif
891
892extern void pdgstrf2(superlu_dist_options_t *, int_t nsupers, int_t k0,
893 int_t k, double thresh, Glu_persist_t *, gridinfo_t *,
894 dLocalLU_t *, MPI_Request *, int, SuperLUStat_t *, int *);
895
896 /* from p3dcomm.h */
897extern int_t dAllocLlu_3d(int_t nsupers, dLUstruct_t * LUstruct, gridinfo3d_t* grid3d);
898extern int_t dp3dScatter(int_t n, dLUstruct_t * LUstruct, gridinfo3d_t* grid3d);
900 dLUstruct_t * LUstruct, gridinfo3d_t* grid3d);
902 dLUstruct_t * LUstruct, gridinfo3d_t* grid3d);
903extern int_t dcollect3dLpanels(int_t layer, int_t nsupers, dLUstruct_t * LUstruct, gridinfo3d_t* grid3d);
904extern int_t dcollect3dUpanels(int_t layer, int_t nsupers, dLUstruct_t * LUstruct, gridinfo3d_t* grid3d);
905extern int_t dp3dCollect(int_t layer, int_t n, dLUstruct_t * LUstruct, gridinfo3d_t* grid3d);
906/*zero out LU non zero entries*/
907extern int_t dzeroSetLU(int_t nnodes, int_t* nodeList , dLUstruct_t *, gridinfo3d_t*);
908extern int dAllocGlu_3d(int_t n, int_t nsupers, dLUstruct_t *);
910extern int dDeAllocGlu_3d(dLUstruct_t *);
911
912/* Reduces L and U panels of nodes in the List nodeList (size=nnnodes)
913receiver[L(nodelist)] =sender[L(nodelist)] +receiver[L(nodelist)]
914receiver[U(nodelist)] =sender[U(nodelist)] +receiver[U(nodelist)]
915*/
917 int_t nnodes, int_t* nodeList,
918 double* Lval_buf, double* Uval_buf,
919 dLUstruct_t* LUstruct, gridinfo3d_t* grid3d, SCT_t* SCT);
920/*reduces all nodelists required in a level*/
921extern int dreduceAllAncestors3d(int_t ilvl, int_t* myNodeCount,
922 int_t** treePerm,
923 dLUValSubBuf_t* LUvsb,
924 dLUstruct_t* LUstruct,
925 gridinfo3d_t* grid3d,
926 SCT_t* SCT );
927/*
928 Copies factored L and U panels from sender grid to receiver grid
929 receiver[L(nodelist)] <-- sender[L(nodelist)];
930 receiver[U(nodelist)] <-- sender[U(nodelist)];
931*/
933 int_t nnodes, int_t *nodeList, dLUValSubBuf_t* LUvsb,
934 dLUstruct_t* LUstruct, gridinfo3d_t* grid3d,SCT_t* SCT );
935
936/*Gathers all the L and U factors to grid 0 for solve stage
937 By repeatidly calling above function*/
939 gridinfo3d_t* grid3d, SCT_t* SCT );
940
941/*Distributes data in each layer and initilizes ancestors
942 as zero in required nodes*/
943int_t dinit3DLUstruct( int_t* myTreeIdxs, int_t* myZeroTrIdxs,
944 int_t* nodeCount, int_t** nodeList,
945 dLUstruct_t* LUstruct, gridinfo3d_t* grid3d);
946
948 dLUstruct_t* LUstruct, gridinfo3d_t* grid3d, SCT_t* SCT);
949int_t dzRecvLPanel(int_t k, int_t sender, double alpha,
950 double beta, double* Lval_buf,
951 dLUstruct_t* LUstruct, gridinfo3d_t* grid3d, SCT_t* SCT);
953 dLUstruct_t* LUstruct, gridinfo3d_t* grid3d, SCT_t* SCT);
954int_t dzRecvUPanel(int_t k, int_t sender, double alpha,
955 double beta, double* Uval_buf,
956 dLUstruct_t* LUstruct, gridinfo3d_t* grid3d, SCT_t* SCT);
957
958 /* from communication_aux.h */
959extern int_t dIBcast_LPanel (int_t k, int_t k0, int_t* lsub, double* lusup,
960 gridinfo_t *, int* msgcnt, MPI_Request *,
961 int **ToSendR, int_t *xsup, int );
962extern int_t dBcast_LPanel(int_t k, int_t k0, int_t* lsub, double* lusup,
963 gridinfo_t *, int* msgcnt, int **ToSendR,
964 int_t *xsup , SCT_t*, int);
965extern int_t dIBcast_UPanel(int_t k, int_t k0, int_t* usub, double* uval,
966 gridinfo_t *, int* msgcnt, MPI_Request *,
967 int *ToSendD, int );
968extern int_t dBcast_UPanel(int_t k, int_t k0, int_t* usub, double* uval,
969 gridinfo_t *, int* msgcnt, int *ToSendD, SCT_t*, int);
970extern int_t dIrecv_LPanel (int_t k, int_t k0, int_t* Lsub_buf,
971 double* Lval_buf, gridinfo_t *,
972 MPI_Request *, dLocalLU_t *, int);
973extern int_t dIrecv_UPanel(int_t k, int_t k0, int_t* Usub_buf, double*,
974 dLocalLU_t *, gridinfo_t*, MPI_Request *, int);
975extern int_t dWait_URecv(MPI_Request *, int* msgcnt, SCT_t *);
976extern int_t dWait_LRecv(MPI_Request*, int* msgcnt, int* msgcntsU,
977 gridinfo_t *, SCT_t*);
978extern int_t dISend_UDiagBlock(int_t k0, double *ublk_ptr, int_t size,
979 MPI_Request *, gridinfo_t *, int);
980extern int_t dRecv_UDiagBlock(int_t k0, double *ublk_ptr, int_t size,
981 int_t src, gridinfo_t *, SCT_t*, int);
982extern int_t dPackLBlock(int_t k, double* Dest, Glu_persist_t *,
984extern int_t dISend_LDiagBlock(int_t k0, double *lblk_ptr, int_t size,
985 MPI_Request *, gridinfo_t *, int);
986extern int_t dIRecv_UDiagBlock(int_t k0, double *ublk_ptr, int_t size,
987 int_t src, MPI_Request *, gridinfo_t *,
988 SCT_t*, int);
989extern int_t dIRecv_LDiagBlock(int_t k0, double *L_blk_ptr, int_t size,
990 int_t src, MPI_Request *, gridinfo_t*, SCT_t*, int);
991extern int_t dUDiagBlockRecvWait( int_t k, int_t* IrecvPlcd_D, int_t* factored_L,
992 MPI_Request *, gridinfo_t *, dLUstruct_t *, SCT_t *);
993extern int_t LDiagBlockRecvWait( int_t k, int_t* factored_U, MPI_Request *, gridinfo_t *);
994
995#if (MPI_VERSION>2)
996extern int_t dIBcast_UDiagBlock(int_t k, double *ublk_ptr, int_t size,
997 MPI_Request *, gridinfo_t *);
998extern int_t dIBcast_LDiagBlock(int_t k, double *lblk_ptr, int_t size,
999 MPI_Request *, gridinfo_t *);
1000#endif
1001
1002 /* from trfCommWrapper.h */
1004 double *BlockUFactor, double *BlockLFactor,
1005 int_t* IrecvPlcd_D, MPI_Request *, MPI_Request *,
1006 MPI_Request *, MPI_Request *, gridinfo_t *,
1007 superlu_dist_options_t *, double thresh,
1008 dLUstruct_t *LUstruct, SuperLUStat_t *, int *info,
1009 SCT_t *, int tag_ub);
1010extern int_t dUPanelTrSolve( int_t k, double* BlockLFactor, double* bigV,
1013extern int_t dLPanelUpdate(int_t k, int_t* IrecvPlcd_D, int_t* factored_L,
1014 MPI_Request *, double* BlockUFactor, gridinfo_t *,
1015 dLUstruct_t *, SCT_t *);
1016extern int_t dUPanelUpdate(int_t k, int_t* factored_U, MPI_Request *,
1017 double* BlockLFactor, double* bigV,
1020extern int_t dIBcastRecvLPanel(int_t k, int_t k0, int* msgcnt,
1021 MPI_Request *, MPI_Request *,
1022 int_t* Lsub_buf, double* Lval_buf,
1024 SCT_t *, int tag_ub);
1025extern int_t dIBcastRecvUPanel(int_t k, int_t k0, int* msgcnt, MPI_Request *,
1026 MPI_Request *, int_t* Usub_buf, double* Uval_buf,
1027 gridinfo_t *, dLUstruct_t *, SCT_t *, int tag_ub);
1028extern int_t dWaitL(int_t k, int* msgcnt, int* msgcntU, MPI_Request *,
1029 MPI_Request *, gridinfo_t *, dLUstruct_t *, SCT_t *);
1030extern int_t dWaitU(int_t k, int* msgcnt, MPI_Request *, MPI_Request *,
1031 gridinfo_t *, dLUstruct_t *, SCT_t *);
1032extern int_t dLPanelTrSolve(int_t k, int_t* factored_L, double* BlockUFactor,
1033 gridinfo_t *, dLUstruct_t *);
1034
1035 /* from trfAux.h */
1036extern int getNsupers(int, Glu_persist_t *);
1037extern int_t initPackLUInfo(int_t nsupers, packLUInfo_t* packLUInfo);
1038extern int freePackLUInfo(packLUInfo_t* packLUInfo);
1041 lPanelInfo_t *, int_t*, int_t *, int_t *,
1042 double *bigU, int_t* Lsub_buf,
1043 double* Lval_buf, int_t* Usub_buf,
1044 double* Uval_buf, gridinfo_t *, dLUstruct_t *);
1048 dLUValSubBuf_t* LUvsb, gridinfo_t *,
1049 dLUstruct_t *, HyP_t*);
1050extern double* dgetBigV(int_t, int_t);
1051extern double* dgetBigU(int_t, gridinfo_t *, dLUstruct_t *);
1052// permutation from superLU default
1053
1054 /* from treeFactorization.h */
1056extern int_t dinitScuBufs(int_t ldt, int_t num_threads, int_t nsupers,
1058extern int dfreeScuBufs(dscuBufs_t* scuBufs);
1059
1060// the generic tree factoring code
1062 int_t nnnodes, // number of nodes in the tree
1063 int_t *perm_c_supno, // list of nodes in the order of factorization
1064 commRequests_t *comReqs, // lists of communication requests
1065 dscuBufs_t *scuBufs, // contains buffers for schur complement update
1066 packLUInfo_t*packLUInfo,
1067 msgs_t*msgs,
1068 dLUValSubBuf_t* LUvsb,
1069 ddiagFactBufs_t *dFBuf,
1070 factStat_t *factStat,
1071 factNodelists_t *fNlists,
1072 superlu_dist_options_t *options,
1073 int_t * gIperm_c_supno,
1074 int_t ldt,
1075 dLUstruct_t *LUstruct, gridinfo3d_t * grid3d, SuperLUStat_t *stat,
1076 double thresh, SCT_t *SCT,
1077 int *info
1078);
1079
1081 int_t nnodes, // number of nodes in the tree
1082 int_t *perm_c_supno, // list of nodes in the order of factorization
1083 treeTopoInfo_t* treeTopoInfo,
1084 commRequests_t *comReqs, // lists of communication requests
1085 dscuBufs_t *scuBufs, // contains buffers for schur complement update
1086 packLUInfo_t*packLUInfo,
1087 msgs_t*msgs,
1088 dLUValSubBuf_t* LUvsb,
1089 ddiagFactBufs_t *dFBuf,
1090 factStat_t *factStat,
1091 factNodelists_t *fNlists,
1092 superlu_dist_options_t *options,
1093 int_t * gIperm_c_supno,
1094 int_t ldt,
1095 dLUstruct_t *LUstruct, gridinfo3d_t * grid3d, SuperLUStat_t *stat,
1096 double thresh, SCT_t *SCT,
1097 int *info
1098);
1099
1101 int_t nnnodes, // number of nodes in the tree
1102 int_t *perm_c_supno, // list of nodes in the order of factorization
1103 commRequests_t *comReqs, // lists of communication requests
1104 dscuBufs_t *scuBufs, // contains buffers for schur complement update
1105 packLUInfo_t*packLUInfo,
1106 msgs_t*msgs,
1107 dLUValSubBuf_t* LUvsb,
1108 ddiagFactBufs_t *dFBuf,
1109 factStat_t *factStat,
1110 factNodelists_t *fNlists,
1111 superlu_dist_options_t *options,
1112 int_t * gIperm_c_supno,
1113 int_t ldt,
1114 dLUstruct_t *LUstruct, gridinfo3d_t * grid3d, SuperLUStat_t *stat,
1115 double thresh, SCT_t *SCT, int tag_ub,
1116 int *info
1117);
1118
1120 sForest_t* sforest,
1121 commRequests_t **comReqss, // lists of communication requests // size maxEtree level
1122 dscuBufs_t *scuBufs, // contains buffers for schur complement update
1123 packLUInfo_t*packLUInfo,
1124 msgs_t**msgss, // size=num Look ahead
1125 dLUValSubBuf_t** LUvsbs, // size=num Look ahead
1126 ddiagFactBufs_t **dFBufs, // size maxEtree level
1127 factStat_t *factStat,
1128 factNodelists_t *fNlists,
1129 gEtreeInfo_t* gEtreeInfo, // global etree info
1130 superlu_dist_options_t *options,
1131 int_t * gIperm_c_supno,
1132 int_t ldt,
1133 HyP_t* HyP,
1134 dLUstruct_t *LUstruct, gridinfo3d_t * grid3d, SuperLUStat_t *stat,
1135 double thresh, SCT_t *SCT, int tag_ub,
1136 int *info
1137);
1139extern int dLluBufFreeArr(int_t numLA, dLUValSubBuf_t **LUvsbs);
1141extern int dfreeDiagFactBufsArr(int_t mxLeafNode, ddiagFactBufs_t** dFBufs);
1143extern int_t checkRecvUDiag(int_t k, commRequests_t *comReqs,
1144 gridinfo_t *grid, SCT_t *SCT);
1145extern int_t checkRecvLDiag(int_t k, commRequests_t *comReqs, gridinfo_t *, SCT_t *);
1146
1147 /* from ancFactorization.h (not called) */
1149 int_t ilvl, // level of factorization
1150 sForest_t* sforest,
1151 commRequests_t **comReqss, // lists of communication requests // size maxEtree level
1152 dscuBufs_t *scuBufs, // contains buffers for schur complement update
1153 packLUInfo_t*packLUInfo,
1154 msgs_t**msgss, // size=num Look ahead
1155 dLUValSubBuf_t** LUvsbs, // size=num Look ahead
1156 ddiagFactBufs_t **dFBufs, // size maxEtree level
1157 factStat_t *factStat,
1158 factNodelists_t *fNlists,
1159 gEtreeInfo_t* gEtreeInfo, // global etree info
1160 superlu_dist_options_t *options,
1161 int_t * gIperm_c_supno,
1162 int_t ldt,
1163 HyP_t* HyP,
1164 dLUstruct_t *LUstruct, gridinfo3d_t * grid3d, SuperLUStat_t *stat,
1165 double thresh, SCT_t *SCT, int tag_ub, int *info
1166);
1167
1168/*== end 3D prototypes ===================*/
1169
1170
1171#ifdef __cplusplus
1172 }
1173#endif
1174
1175#endif /* __SUPERLU_dDEFS */
1176
int j
Definition: dutil_dist.c:248
#define dtrtri_
Definition: superlu_FCnames.h:161
int superlu_dger(const int m, const int n, const double alpha, const double *x, const int incx, const double *y, const int incy, double *a, const int lda)
int_t dLluBufInit(dLUValSubBuf_t *, dLUstruct_t *)
int_t dscatter3dLPanels(int_t nsupers, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
int dstatic_schedule(superlu_dist_options_t *, int, int, dLUstruct_t *, gridinfo_t *, SuperLUStat_t *, int_t *, int_t *, int *)
Definition: dstatic_schedule.c:46
void dSolveFinalize(superlu_dist_options_t *, dSOLVEstruct_t *)
Release the resources used for the solution phase.
Definition: pdutil.c:874
void dScaleAddId_CompRowLoc_Matrix_dist(SuperMatrix *, double)
Scale and add I: scales a matrix and adds an identity. A_{i,j} = c * A_{i,j} + \delta_{i,...
Definition: dutil_dist.c:356
int_t initPackLUInfo(int_t nsupers, packLUInfo_t *packLUInfo)
Definition: treeFactorization.c:367
void dScaleAdd_CompRowLoc_Matrix_dist(SuperMatrix *, SuperMatrix *, double)
Scale and add: adds a scalar multiple of one matrix to another. A_{i,j} = c * A_{i,...
Definition: dutil_dist.c:381
void dscatter_l(int ib, int ljb, int nsupc, int_t iukp, int_t *xsup, int klst, int nbrow, int_t lptr, int temp_nbrow, int_t *usub, int_t *lsub, double *tempv, int *indirect_thread, int *indirect2, int_t **Lrowind_bc_ptr, double **Lnzval_bc_ptr, gridinfo_t *grid)
int dDeAllocLlu_3d(int_t n, dLUstruct_t *, gridinfo3d_t *)
Definition: dutil_dist.c:450
int_t pdReDistribute_B_to_X(double *B, int_t m_loc, int nrhs, int_t ldb, int_t fst_row, int_t *ilsum, double *x, dScalePermstruct_t *, Glu_persist_t *, gridinfo_t *, dSOLVEstruct_t *)
Definition: pdgstrs.c:347
int_t dPackLBlock(int_t k, double *Dest, Glu_persist_t *, gridinfo_t *, dLocalLU_t *)
int updateDirtyBit(int_t k0, HyP_t *HyP, gridinfo_t *grid)
Definition: sec_structs.c:618
void pdGetDiagU(int_t, dLUstruct_t *, gridinfo_t *, double *)
Definition: pdGetDiagU.c:66
int superlu_daxpy(const int n, const double alpha, const double *x, const int incx, double *y, const int incy)
int dcreate_matrix(SuperMatrix *, int, double **, int *, double **, int *, FILE *, gridinfo_t *)
Definition: dcreate_matrix.c:349
void pdgsrfs(int_t, SuperMatrix *, double, dLUstruct_t *, dScalePermstruct_t *, gridinfo_t *, double[], int_t, double[], int_t, int, dSOLVEstruct_t *, double *, SuperLUStat_t *, int *)
int_t dQuerySpace_dist(int_t, dLUstruct_t *, gridinfo_t *, SuperLUStat_t *, superlu_dist_mem_usage_t *)
Definition: dmemory_dist.c:73
int dcreate_matrix3d(SuperMatrix *A, int nrhs, double **rhs, int *ldb, double **x, int *ldx, FILE *fp, gridinfo3d_t *grid3d)
void dGenCSCLblocks(int, int_t, gridinfo_t *, Glu_persist_t *, dLocalLU_t *, double **, int_t **, int_t **, int_t *, int_t *)
Definition: dutil_dist.c:901
void dZeroUblocks(int iam, int n, gridinfo_t *, dLUstruct_t *)
Sets all entries of matrix U to zero.
Definition: dutil_dist.c:1249
void pdgsmv_init(SuperMatrix *, int_t *, gridinfo_t *, pdgsmv_comm_t *)
Definition: pdgsmv.c:27
int_t dISend_UDiagBlock(int_t k0, double *ublk_ptr, int_t size, MPI_Request *, gridinfo_t *, int)
int dfreeScuBufs(dscuBufs_t *scuBufs)
int_t dUDiagBlockRecvWait(int_t k, int_t *IrecvPlcd_D, int_t *factored_L, MPI_Request *, gridinfo_t *, dLUstruct_t *, SCT_t *)
void pdgsmv(int_t, SuperMatrix *, gridinfo_t *, pdgsmv_comm_t *, double x[], double ax[])
Definition: pdgsmv.c:235
int_t dreduceAncestors3d(int_t sender, int_t receiver, int_t nnodes, int_t *nodeList, double *Lval_buf, double *Uval_buf, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
void dGenCSRLblocks(int, int_t, gridinfo_t *, Glu_persist_t *, dLocalLU_t *, double **, int_t **, int_t **, int_t *, int_t *)
Definition: dutil_dist.c:1054
int_t dTrs2_GatherTrsmScatter(int_t klst, int_t iukp, int_t rukp, int_t *usub, double *uval, double *tempv, int_t knsupc, int nsupr, double *lusup, Glu_persist_t *Glu_persist)
Definition: pdgstrf2.c:716
int_t dcollect3dUpanels(int_t layer, int_t nsupers, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
int_t dIRecv_UDiagBlock(int_t k0, double *ublk_ptr, int_t size, int_t src, MPI_Request *, gridinfo_t *, SCT_t *, int)
int_t dinitScuBufs(int_t ldt, int_t num_threads, int_t nsupers, dscuBufs_t *, dLUstruct_t *, gridinfo_t *)
int_t ddenseTreeFactor(int_t nnnodes, int_t *perm_c_supno, commRequests_t *comReqs, dscuBufs_t *scuBufs, packLUInfo_t *packLUInfo, msgs_t *msgs, dLUValSubBuf_t *LUvsb, ddiagFactBufs_t *dFBuf, factStat_t *factStat, factNodelists_t *fNlists, superlu_dist_options_t *options, int_t *gIperm_c_supno, int_t ldt, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SuperLUStat_t *stat, double thresh, SCT_t *SCT, int tag_ub, int *info)
void dinit3DLUstructForest(int_t *myTreeIdxs, int_t *myZeroTrIdxs, sForest_t **sForests, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
int_t scuStatUpdate(int_t knsupc, HyP_t *HyP, SCT_t *SCT, SuperLUStat_t *stat)
Definition: sec_structs.c:635
int pdgsmv_AXglobal_abs(int_t, int_t[], double[], int_t[], double[], double[])
Definition: pdgsmv_AXglobal.c:285
int_t dgatherAllFactoredLU(trf3Dpartition_t *trf3Dpartition, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
void dreadrb_dist(int, FILE *, int_t *, int_t *, int_t *, double **, int_t **, int_t **)
Definition: dreadrb.c:275
void pdCompute_Diag_Inv(int_t, dLUstruct_t *, gridinfo_t *, SuperLUStat_t *, int *)
Definition: pdgstrs.c:839
void pdgsequ(SuperMatrix *, double *, double *, double *, double *, double *, int_t *, gridinfo_t *)
Definition: pdgsequ.c:86
int dscal_(const int *n, const double *alpha, double *dx, const int *incx)
int_t dWaitL(int_t k, int *msgcnt, int *msgcntU, MPI_Request *, MPI_Request *, gridinfo_t *, dLUstruct_t *, SCT_t *)
int dtrsm_(const char *, const char *, const char *, const char *, const int *, const int *, const double *, const double *, const int *, double *, const int *)
void dGenXtrue_dist(int_t, int_t, double *, int_t)
Definition: dutil_dist.c:486
void d3D_printMemUse(trf3Dpartition_t *trf3Dpartition, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
Definition: dmemory_dist.c:243
void dgather_l(int_t num_LBlk, int_t knsupc, Remain_info_t *L_info, double *lval, int_t LD_lval, double *L_buff)
int sp_dgemv_dist(char *, double, SuperMatrix *, double *, int, double, double *, int)
SpGEMV.
Definition: dsp_blas2_dist.c:391
int_t dzRecvLPanel(int_t k, int_t sender, double alpha, double beta, double *Lval_buf, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
int_t dIBcast_UPanel(int_t k, int_t k0, int_t *usub, double *uval, gridinfo_t *, int *msgcnt, MPI_Request *, int *ToSendD, int)
int getNsupers(int, Glu_persist_t *)
Definition: trfAux.c:42
int dfreeDiagFactBufsArr(int_t mxLeafNode, ddiagFactBufs_t **dFBufs)
double pdlangs(char *, SuperMatrix *, gridinfo_t *)
Definition: pdlangs.c:65
int daxpy_(const int *n, const double *alpha, const double *x, const int *incx, double *y, const int *incy)
trf3Dpartition_t * dinitTrf3Dpartition(int_t nsupers, superlu_dist_options_t *options, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
void dGenCOOLblocks(int, int_t, gridinfo_t *, Glu_persist_t *, dLocalLU_t *, int_t **, int_t **, double **, int_t *, int_t *)
Definition: dutil_dist.c:794
int_t dSchurComplementSetup(int_t k, int *msgcnt, Ublock_info_t *, Remain_info_t *, uPanelInfo_t *, lPanelInfo_t *, int_t *, int_t *, int_t *, double *bigU, int_t *Lsub_buf, double *Lval_buf, int_t *Usub_buf, double *Uval_buf, gridinfo_t *, dLUstruct_t *)
void dDestroy_A3d_gathered_on_2d(dSOLVEstruct_t *, gridinfo3d_t *)
Definition: pdutil.c:894
int_t dzSendLPanel(int_t k, int_t receiver, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
int_t dp3dScatter(int_t n, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
void dCreate_CompRowLoc_Matrix_dist(SuperMatrix *, int_t, int_t, int_t, int_t, int_t, double *, int_t *, int_t *, Stype_t, Dtype_t, Mtype_t)
void Free_HyP(HyP_t *HyP)
Definition: sec_structs.c:594
void dlsum_bmod_inv(double *, double *, double *, double *, int, int_t, int *bmod, int_t *, Ucb_indptr_t **, int_t **, int_t *, gridinfo_t *, dLocalLU_t *, SuperLUStat_t **, int_t *, int_t *, int_t, int_t, int, int)
Definition: pdgstrs_lsum.c:1355
int_t dblock_gemm_scatterTopLeft(int_t lb, int_t j, double *bigV, int_t knsupc, int_t klst, int_t *lsub, int_t *usub, int_t ldt, int *indirect, int *indirect2, HyP_t *HyP, dLUstruct_t *, gridinfo_t *, SCT_t *SCT, SuperLUStat_t *)
void dlsum_fmod(double *, double *, double *, double *, int, int, int_t, int *fmod, int_t, int_t, int_t, int_t *, gridinfo_t *, dLocalLU_t *, MPI_Request[], SuperLUStat_t *)
Definition: pdgstrs_lsum.c:62
int_t dIRecv_LDiagBlock(int_t k0, double *L_blk_ptr, int_t size, int_t src, MPI_Request *, gridinfo_t *, SCT_t *, int)
void dScalePermstructFree(dScalePermstruct_t *)
Deallocate ScalePermstruct.
Definition: dutil_dist.c:409
int_t dzSendUPanel(int_t k, int_t receiver, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
int pdPermute_Dense_Matrix(int_t, int_t, int_t[], int_t[], double[], int, double[], int, int, gridinfo_t *)
Permute the distributed dense matrix: B <= perm(X). perm[i] = j means the i-th row of X is in the j-t...
Definition: pdutil.c:294
void dRgather_U(int_t k, int_t jj0, int_t *usub, double *uval, double *bigU, gEtreeInfo_t *, Glu_persist_t *, gridinfo_t *, HyP_t *, int_t *myIperm, int_t *iperm_c_supno, int_t *perm_u)
void dlsum_bmod_inv_master(double *, double *, double *, double *, int, int_t, int *bmod, int_t *, Ucb_indptr_t **, int_t **, int_t *, gridinfo_t *, dLocalLU_t *, SuperLUStat_t **, int_t, int_t, int, int)
Definition: pdgstrs_lsum.c:1819
int_t checkRecvUDiag(int_t k, commRequests_t *comReqs, gridinfo_t *grid, SCT_t *SCT)
Definition: treeFactorization.c:401
int pdgsmv_AXglobal(int_t, int_t[], double[], int_t[], double[], double[])
Definition: pdgsmv_AXglobal.c:259
int_t dAllocLlu_3d(int_t nsupers, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
void pdgstrs2(int_t m, int_t k0, int_t k, Glu_persist_t *Glu_persist, gridinfo_t *grid, dLocalLU_t *Llu, SuperLUStat_t *stat)
void dgemv_(const char *, const int *, const int *, const double *, const double *a, const int *, const double *, const int *, const double *, double *, const int *)
ddiagFactBufs_t ** dinitDiagFactBufsArr(int_t mxLeafNode, int_t ldt, gridinfo_t *grid)
int_t dIBcastRecvLPanel(int_t k, int_t k0, int *msgcnt, MPI_Request *, MPI_Request *, int_t *Lsub_buf, double *Lval_buf, int_t *factored, gridinfo_t *, dLUstruct_t *, SCT_t *, int tag_ub)
void dinf_norm_error_dist(int_t, int_t, double *, int_t, double *, int_t, gridinfo_t *)
Check the inf-norm of the error vector.
Definition: dutil_dist.c:522
float ddistribute(fact_t, int_t, SuperMatrix *, Glu_freeable_t *, dLUstruct_t *, gridinfo_t *)
Definition: ddistribute.c:65
int pdgsmv_AXglobal_setup(SuperMatrix *, Glu_persist_t *, gridinfo_t *, int_t *, int_t *[], double *[], int_t *[], int_t[])
void pdgsrfs_ABXglobal(int_t, SuperMatrix *, double, dLUstruct_t *, gridinfo_t *, double *, int_t, double *, int_t, int, double *, SuperLUStat_t *, int *)
Definition: pdgsrfs_ABXglobal.c:125
void dGatherNRformat_loc3d(fact_t Fact, NRformat_loc *A, double *B, int ldb, int nrhs, gridinfo3d_t *grid3d, NRformat_loc3d **)
int_t dcollect3dLpanels(int_t layer, int_t nsupers, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
int dcreate_matrix_postfix3d(SuperMatrix *A, int nrhs, double **rhs, int *ldb, double **x, int *ldx, FILE *fp, char *postfix, gridinfo3d_t *grid3d)
Definition: dcreate_matrix3d.c:72
void dreadMM_dist(FILE *, int_t *, int_t *, int_t *, double **, int_t **, int_t **)
Definition: dreadMM.c:38
void dInit_HyP(HyP_t *HyP, dLocalLU_t *Llu, int_t mcb, int_t mrb)
void dlsum_fmod_inv_master(double *, double *, double *, double *, int, int, int_t, int *fmod, int_t, int_t *, gridinfo_t *, dLocalLU_t *, SuperLUStat_t **, int_t, int_t, int_t, int_t, int, int)
Definition: pdgstrs_lsum.c:956
int freePackLUInfo(packLUInfo_t *packLUInfo)
Definition: treeFactorization.c:376
void dDestroy_Tree(int_t, gridinfo_t *, dLUstruct_t *)
Destroy broadcast and reduction trees used in triangular solve.
Definition: pdutil.c:437
void dblock_gemm_scatter(int_t lb, int_t j, Ublock_info_t *Ublock_info, Remain_info_t *Remain_info, double *L_mat, int ldl, double *U_mat, int ldu, double *bigV, int_t knsupc, int_t klst, int_t *lsub, int_t *usub, int_t ldt, int_t thread_id, int *indirect, int *indirect2, int_t **Lrowind_bc_ptr, double **Lnzval_bc_ptr, int_t **Ufstnz_br_ptr, double **Unzval_br_ptr, int_t *xsup, gridinfo_t *, SuperLUStat_t *)
int_t dscatter3dUPanels(int_t nsupers, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
void dDestroy_LU(int_t, gridinfo_t *, dLUstruct_t *)
Destroy distributed L & U matrices.
Definition: pdutil.c:485
void pdgssvx3d(superlu_dist_options_t *, SuperMatrix *, dScalePermstruct_t *, double B[], int ldb, int nrhs, gridinfo3d_t *, dLUstruct_t *, dSOLVEstruct_t *, double *berr, SuperLUStat_t *, int *info)
Definition: pdgssvx3d.c:501
int_t dWait_LRecv(MPI_Request *, int *msgcnt, int *msgcntsU, gridinfo_t *, SCT_t *)
void dCopy_CompRowLoc_Matrix_dist(SuperMatrix *, SuperMatrix *)
Definition: dutil_dist.c:324
void dCreate_CompCol_Matrix_dist(SuperMatrix *, int_t, int_t, int_t, double *, int_t *, int_t *, Stype_t, Dtype_t, Mtype_t)
int_t dgatherAllFactoredLUFr(int_t *myZeroTrIdxs, sForest_t *sForests, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
int_t dTrs2_ScatterU(int_t iukp, int_t rukp, int_t klst, int_t nsupc, int_t ldu, int_t *usub, double *uval, double *tempv)
Definition: pdgstrf2.c:694
int_t dDiagFactIBCast(int_t k, int_t k0, double *BlockUFactor, double *BlockLFactor, int_t *IrecvPlcd_D, MPI_Request *, MPI_Request *, MPI_Request *, MPI_Request *, gridinfo_t *, superlu_dist_options_t *, double thresh, dLUstruct_t *LUstruct, SuperLUStat_t *, int *info, SCT_t *, int tag_ub)
int_t dsparseTreeFactor(int_t nnodes, int_t *perm_c_supno, treeTopoInfo_t *treeTopoInfo, commRequests_t *comReqs, dscuBufs_t *scuBufs, packLUInfo_t *packLUInfo, msgs_t *msgs, dLUValSubBuf_t *LUvsb, ddiagFactBufs_t *dFBuf, factStat_t *factStat, factNodelists_t *fNlists, superlu_dist_options_t *options, int_t *gIperm_c_supno, int_t ldt, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SuperLUStat_t *stat, double thresh, SCT_t *SCT, int *info)
void dPrint_Dense_Matrix_dist(SuperMatrix *)
int_t dblock_gemm_scatterTopRight(int_t lb, int_t j, double *bigV, int_t knsupc, int_t klst, int_t *lsub, int_t *usub, int_t ldt, int *indirect, int *indirect2, HyP_t *HyP, dLUstruct_t *, gridinfo_t *, SCT_t *SCT, SuperLUStat_t *)
void dRgather_L(int_t k, int_t *lsub, double *lusup, gEtreeInfo_t *, Glu_persist_t *, gridinfo_t *, HyP_t *, int_t *myIperm, int_t *iperm_c_supno)
int_t dblock_gemm_scatterBottomLeft(int_t lb, int_t j, double *bigV, int_t knsupc, int_t klst, int_t *lsub, int_t *usub, int_t ldt, int *indirect, int *indirect2, HyP_t *HyP, dLUstruct_t *, gridinfo_t *, SCT_t *SCT, SuperLUStat_t *)
void pdgstrs(int_t, dLUstruct_t *, dScalePermstruct_t *, gridinfo_t *, double *, int_t, int_t, int_t, int, dSOLVEstruct_t *, SuperLUStat_t *, int *)
Definition: pdgstrs.c:1030
int dScatter_B3d(NRformat_loc3d *A3d, gridinfo3d_t *grid3d)
void dgather_u(int_t num_u_blks, Ublock_info_t *Ublock_info, int_t *usub, double *uval, double *bigU, int_t ldu, int_t *xsup, int_t klst)
int file_dPrint_CompRowLoc_Matrix_dist(FILE *fp, SuperMatrix *A)
int_t dIBcastRecvUPanel(int_t k, int_t k0, int *msgcnt, MPI_Request *, MPI_Request *, int_t *Usub_buf, double *Uval_buf, gridinfo_t *, dLUstruct_t *, SCT_t *, int tag_ub)
void dZeroLblocks(int, int, gridinfo_t *, dLUstruct_t *)
Sets all entries of matrix L to zero.
Definition: dutil_dist.c:619
int sp_dtrsv_dist(char *, char *, char *, SuperMatrix *, SuperMatrix *, double *, int *)
Definition: dsp_blas2_dist.c:95
#define MAX_LOOKAHEADS
Definition: superlu_ddefs.h:96
int dldperm_dist(int, int, int_t, int_t[], int_t[], double[], int_t *, double[], double[])
Definition: dldperm_dist.c:96
int_t dLpanelUpdate(int_t off0, int_t nsupc, double *ublk_ptr, int_t ld_ujrow, double *lusup, int_t nsupr, SCT_t *)
void * duser_malloc_dist(int_t, int_t)
Definition: dmemory_dist.c:30
void Local_Dgstrf2(superlu_dist_options_t *options, int_t k, double thresh, double *BlockUFactor, Glu_persist_t *, gridinfo_t *, dLocalLU_t *, SuperLUStat_t *, int *info, SCT_t *)
Definition: pdgstrf2.c:402
int_t LDiagBlockRecvWait(int_t k, int_t *factored_U, MPI_Request *, gridinfo_t *)
Definition: communication_aux.c:218
void pdgstrs_Bglobal(int_t, dLUstruct_t *, gridinfo_t *, double *, int_t, int, SuperLUStat_t *, int *)
Definition: pdgstrs_Bglobal.c:104
int_t dWait_URecv(MPI_Request *, int *msgcnt, SCT_t *)
void dLUstructFree(dLUstruct_t *)
Deallocate LUstruct.
Definition: pdutil.c:419
int dPrint_CompRowLoc_Matrix_dist(SuperMatrix *)
void dCreate_SuperNode_Matrix_dist(SuperMatrix *, int_t, int_t, int_t, double *, int_t *, int_t *, int_t *, int_t *, int_t *, Stype_t, Dtype_t, Mtype_t)
void pdgssvx_ABglobal(superlu_dist_options_t *, SuperMatrix *, dScalePermstruct_t *, double *, int, int, gridinfo_t *, dLUstruct_t *, double *, SuperLUStat_t *, int *)
void dPrint_CompCol_Matrix_dist(SuperMatrix *)
int dtrsv_(char *, char *, char *, int *, double *, int *, double *, int *)
int_t dIrecv_LPanel(int_t k, int_t k0, int_t *Lsub_buf, double *Lval_buf, gridinfo_t *, MPI_Request *, dLocalLU_t *, int)
void pdgstrf2_trsm(superlu_dist_options_t *options, int_t k0, int_t k, double thresh, Glu_persist_t *, gridinfo_t *, dLocalLU_t *, MPI_Request *, int tag_ub, SuperLUStat_t *, int *info)
Definition: pdgstrf2.c:143
int_t pdgstrs_init(int_t, int_t, int_t, int_t, int_t[], int_t[], gridinfo_t *grid, Glu_persist_t *, dSOLVEstruct_t *)
Definition: pdutil.c:657
int_t treeFactor(int_t nnnodes, int_t *perm_c_supno, commRequests_t *comReqs, dscuBufs_t *scuBufs, packLUInfo_t *packLUInfo, msgs_t *msgs, dLUValSubBuf_t *LUvsb, ddiagFactBufs_t *dFBuf, factStat_t *factStat, factNodelists_t *fNlists, superlu_dist_options_t *options, int_t *gIperm_c_supno, int_t ldt, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SuperLUStat_t *stat, double thresh, SCT_t *SCT, int *info)
int_t dISend_LDiagBlock(int_t k0, double *lblk_ptr, int_t size, MPI_Request *, gridinfo_t *, int)
int_t dgatherFactoredLU(int_t sender, int_t receiver, int_t nnodes, int_t *nodeList, dLUValSubBuf_t *LUvsb, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
int_t dinitDiagFactBufs(int_t ldt, ddiagFactBufs_t *dFBuf)
int_t pdgstrf(superlu_dist_options_t *, int, int, double anorm, dLUstruct_t *, gridinfo_t *, SuperLUStat_t *, int *)
Definition: pdgstrf.c:242
double * dgetBigU(int_t, gridinfo_t *, dLUstruct_t *)
int_t dSchurComplementSetupGPU(int_t k, msgs_t *msgs, packLUInfo_t *, int_t *, int_t *, int_t *, gEtreeInfo_t *, factNodelists_t *, dscuBufs_t *, dLUValSubBuf_t *LUvsb, gridinfo_t *, dLUstruct_t *, HyP_t *)
int dAllocGlu_3d(int_t n, int_t nsupers, dLUstruct_t *)
Definition: dutil_dist.c:432
int_t dBcast_UPanel(int_t k, int_t k0, int_t *usub, double *uval, gridinfo_t *, int *msgcnt, int *ToSendD, SCT_t *, int)
int file_Printdouble5(FILE *, char *, int_t, double *)
Definition: dutil_dist.c:555
void dreadhb_dist(int, FILE *, int_t *, int_t *, int_t *, double **, int_t **, int_t **)
Definition: dreadhb.c:107
int_t dTrs2_GatherU(int_t iukp, int_t rukp, int_t klst, int_t nsupc, int_t ldu, int_t *usub, double *uval, double *tempv)
Definition: pdgstrf2.c:669
void pdgstrs2_omp(int_t k0, int_t k, Glu_persist_t *, gridinfo_t *, dLocalLU_t *, Ublock_info_t *, SuperLUStat_t *)
Definition: pdgstrf2.c:762
int superlu_dgemv(const char *trans, const int m, const int n, const double alpha, const double *a, const int lda, const double *x, const int incx, const double beta, double *y, const int incy)
void dCompRow_to_CompCol_dist(int_t, int_t, int_t, double *, int_t *, int_t *, double **, int_t **, int_t **)
void dscatter_u(int ib, int jb, int nsupc, int_t iukp, int_t *xsup, int klst, int nbrow, int_t lptr, int temp_nbrow, int_t *lsub, int_t *usub, double *tempv, int_t **Ufstnz_br_ptr, double **Unzval_br_ptr, gridinfo_t *grid)
void dgsequ_dist(SuperMatrix *, double *, double *, double *, double *, double *, int_t *)
Definition: dgsequ_dist.c:84
int_t dRecv_UDiagBlock(int_t k0, double *ublk_ptr, int_t size, int_t src, gridinfo_t *, SCT_t *, int)
void dLUstructInit(const int_t, dLUstruct_t *)
Allocate storage in LUstruct.
Definition: pdutil.c:405
void dCreate_Dense_Matrix_dist(SuperMatrix *, int_t, int_t, double *, int_t, Stype_t, Dtype_t, Mtype_t)
void dallocateA_dist(int_t, int_t, double **, int_t **, int_t **)
Definition: dmemory_dist.c:147
void dlsum_fmod_inv(double *, double *, double *, double *, int, int_t, int *fmod, int_t *, gridinfo_t *, dLocalLU_t *, SuperLUStat_t **, int_t *, int_t *, int_t, int_t, int_t, int_t, int, int)
Definition: pdgstrs_lsum.c:416
int dcreate_matrix_rb(SuperMatrix *, int, double **, int *, double **, int *, FILE *, gridinfo_t *)
float pddistribute(fact_t, int_t, SuperMatrix *, dScalePermstruct_t *, Glu_freeable_t *, dLUstruct_t *, gridinfo_t *)
Definition: pddistribute.c:328
int superlu_dtrsm(const char *sideRL, const char *uplo, const char *transa, const char *diag, const int m, const int n, const double alpha, const double *a, const int lda, double *b, const int ldb)
void pxgstrs_finalize(pxgstrs_comm_t *)
Definition: util.c:266
void dlsum_bmod(double *, double *, double *, int, int_t, int *bmod, int_t *, Ucb_indptr_t **, int_t **, int_t *, gridinfo_t *, dLocalLU_t *, MPI_Request[], SuperLUStat_t *)
Definition: pdgstrs_lsum.c:246
void pdinf_norm_error(int, int_t, int_t, double[], int_t, double[], int_t, MPI_Comm)
Check the inf-norm of the error vector.
Definition: pdutil.c:921
int dread_binary(FILE *, int_t *, int_t *, int_t *, double **, int_t **, int_t **)
Definition: dbinary_io.c:4
double * doubleMalloc_dist(int_t)
Definition: dmemory_dist.c:155
int pdCompRow_loc_to_CompCol_global(int_t, SuperMatrix *, gridinfo_t *, SuperMatrix *)
Gather A from the distributed compressed row format to global A in compressed column format.
Definition: pdutil.c:32
int dSolveInit(superlu_dist_options_t *, SuperMatrix *, int_t[], int_t[], int_t, dLUstruct_t *, gridinfo_t *, dSOLVEstruct_t *)
Initialize the data structure for the solution phase.
Definition: pdutil.c:786
int_t dIBcast_LPanel(int_t k, int_t k0, int_t *lsub, double *lusup, gridinfo_t *, int *msgcnt, MPI_Request *, int **ToSendR, int_t *xsup, int)
void pdgstrf2(superlu_dist_options_t *, int_t nsupers, int_t k0, int_t k, double thresh, Glu_persist_t *, gridinfo_t *, dLocalLU_t *, MPI_Request *, int, SuperLUStat_t *, int *)
void dPrintUblocks(int, int_t, gridinfo_t *, Glu_persist_t *, dLocalLU_t *)
Print the blocks in the factored matrix U.
Definition: dutil_dist.c:1208
int dreduceAllAncestors3d(int_t ilvl, int_t *myNodeCount, int_t **treePerm, dLUValSubBuf_t *LUvsb, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
int_t dzeroSetLU(int_t nnodes, int_t *nodeList, dLUstruct_t *, gridinfo3d_t *)
int_t dWaitU(int_t k, int *msgcnt, MPI_Request *, MPI_Request *, gridinfo_t *, dLUstruct_t *, SCT_t *)
int dcreate_matrix_postfix(SuperMatrix *, int, double **, int *, double **, int *, FILE *, char *, gridinfo_t *)
Definition: dcreate_matrix.c:76
int_t dUPanelTrSolve(int_t k, double *BlockLFactor, double *bigV, int_t ldt, Ublock_info_t *, gridinfo_t *, dLUstruct_t *, SuperLUStat_t *, SCT_t *)
int superlu_dgemm(const char *transa, const char *transb, int m, int n, int k, double alpha, double *a, int lda, double *b, int ldb, double beta, double *c, int ldc)
int superlu_dscal(const int n, const double alpha, double *x, const int incx)
int dDeAllocGlu_3d(dLUstruct_t *)
Definition: dutil_dist.c:442
void dClone_CompRowLoc_Matrix_dist(SuperMatrix *, SuperMatrix *)
double dlangs_dist(char *, SuperMatrix *)
Definition: dlangs_dist.c:62
int dLluBufFreeArr(int_t numLA, dLUValSubBuf_t **LUvsbs)
int_t ancestorFactor(int_t ilvl, sForest_t *sforest, commRequests_t **comReqss, dscuBufs_t *scuBufs, packLUInfo_t *packLUInfo, msgs_t **msgss, dLUValSubBuf_t **LUvsbs, ddiagFactBufs_t **dFBufs, factStat_t *factStat, factNodelists_t *fNlists, gEtreeInfo_t *gEtreeInfo, superlu_dist_options_t *options, int_t *gIperm_c_supno, int_t ldt, HyP_t *HyP, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SuperLUStat_t *stat, double thresh, SCT_t *SCT, int tag_ub, int *info)
void pdgsmv_finalize(pdgsmv_comm_t *)
Definition: pdgsmv.c:371
void dScalePermstructInit(const int_t, const int_t, dScalePermstruct_t *)
Allocate storage in ScalePermstruct.
Definition: dutil_dist.c:398
void dlaqgs_dist(SuperMatrix *, double *, double *, double, double, double, char *)
Definition: dlaqgs_dist.c:81
float ddist_psymbtonum(fact_t, int_t, SuperMatrix *, dScalePermstruct_t *, Pslu_freeable_t *, dLUstruct_t *, gridinfo_t *)
Definition: pdsymbfact_distdata.c:1188
int d_c2cpp_GetHWPM(SuperMatrix *, gridinfo_t *, dScalePermstruct_t *)
Definition: d_c2cpp_GetHWPM.cpp:53
int dgemm_(const char *, const char *, const int *, const int *, const int *, const double *, const double *, const int *, const double *, const int *, const double *, double *, const int *)
void dger_(const int *, const int *, const double *, const double *, const int *, const double *, const int *, double *, const int *)
int superlu_dtrsv(char *uplo, char *trans, char *diag, int n, double *a, int lda, double *x, int incx)
void dCopy_Dense_Matrix_dist(int_t, int_t, double *, int_t, double *, int_t)
int dcreate_matrix_dat(SuperMatrix *, int, double **, int *, double **, int *, FILE *, gridinfo_t *)
int_t dBcast_LPanel(int_t k, int_t k0, int_t *lsub, double *lusup, gridinfo_t *, int *msgcnt, int **ToSendR, int_t *xsup, SCT_t *, int)
int_t dLPanelTrSolve(int_t k, int_t *factored_L, double *BlockUFactor, gridinfo_t *, dLUstruct_t *)
void dreadtriple_dist(FILE *, int_t *, int_t *, int_t *, double **, int_t **, int_t **)
Definition: dreadtriple.c:35
int_t pdgstrf3d(superlu_dist_options_t *, int m, int n, double anorm, trf3Dpartition_t *, SCT_t *, dLUstruct_t *, gridinfo3d_t *, SuperLUStat_t *, int *)
Definition: pdgstrf3d.c:121
void duser_free_dist(int_t, int_t)
Definition: dmemory_dist.c:49
void dZero_CompRowLoc_Matrix_dist(SuperMatrix *)
Sets all entries of a matrix to zero, A_{i,j}=0, for i,j=1,..,n.
Definition: dutil_dist.c:339
void dDestroy_trf3Dpartition(trf3Dpartition_t *trf3Dpartition, gridinfo3d_t *grid3d)
int_t dzRecvUPanel(int_t k, int_t sender, double alpha, double beta, double *Uval_buf, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
int_t dinit3DLUstruct(int_t *myTreeIdxs, int_t *myZeroTrIdxs, int_t *nodeCount, int_t **nodeList, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
int_t checkRecvLDiag(int_t k, commRequests_t *comReqs, gridinfo_t *, SCT_t *)
Definition: treeFactorization.c:422
double * dgetBigV(int_t, int_t)
void dFillRHS_dist(char *, int_t, double *, int_t, SuperMatrix *, double *, int_t)
Let rhs[i] = sum of i-th row of A, so the solution vector is all 1's.
Definition: dutil_dist.c:499
void dreadtriple_noheader(FILE *, int_t *, int_t *, int_t *, double **, int_t **, int_t **)
Definition: dreadtriple_noheader.c:35
void dComputeLevelsets(int, int_t, gridinfo_t *, Glu_persist_t *, dLocalLU_t *, int_t *)
Definition: dutil_dist.c:752
int_t dsparseTreeFactor_ASYNC(sForest_t *sforest, commRequests_t **comReqss, dscuBufs_t *scuBufs, packLUInfo_t *packLUInfo, msgs_t **msgss, dLUValSubBuf_t **LUvsbs, ddiagFactBufs_t **dFBufs, factStat_t *factStat, factNodelists_t *fNlists, gEtreeInfo_t *gEtreeInfo, superlu_dist_options_t *options, int_t *gIperm_c_supno, int_t ldt, HyP_t *HyP, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SuperLUStat_t *stat, double thresh, SCT_t *SCT, int tag_ub, int *info)
double * doubleCalloc_dist(int_t)
Definition: dmemory_dist.c:162
int_t dUPanelUpdate(int_t k, int_t *factored_U, MPI_Request *, double *BlockLFactor, double *bigV, int_t ldt, Ublock_info_t *, gridinfo_t *, dLUstruct_t *, SuperLUStat_t *, SCT_t *)
void Printdouble5(char *, int_t, double *)
Definition: dutil_dist.c:543
int_t dp3dCollect(int_t layer, int_t n, dLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
void dPrintLblocks(int, int_t, gridinfo_t *, Glu_persist_t *, dLocalLU_t *)
Print the blocks in the factored matrix L.
Definition: dutil_dist.c:570
dLUValSubBuf_t ** dLluBufInitArr(int_t numLA, dLUstruct_t *LUstruct)
int sp_dgemm_dist(char *, int, double, SuperMatrix *, double *, int, double, double *, int)
Definition: dsp_blas3_dist.c:124
void pdlaqgs(SuperMatrix *, double *, double *, double, double, double, char *)
Definition: pdlaqgs.c:85
int_t dblock_gemm_scatterBottomRight(int_t lb, int_t j, double *bigV, int_t knsupc, int_t klst, int_t *lsub, int_t *usub, int_t ldt, int *indirect, int *indirect2, HyP_t *HyP, dLUstruct_t *, gridinfo_t *, SCT_t *SCT, SuperLUStat_t *)
int_t dIrecv_UPanel(int_t k, int_t k0, int_t *Usub_buf, double *, dLocalLU_t *, gridinfo_t *, MPI_Request *, int)
int_t dLPanelUpdate(int_t k, int_t *IrecvPlcd_D, int_t *factored_L, MPI_Request *, double *BlockUFactor, gridinfo_t *, dLUstruct_t *, SCT_t *)
void dCopy_CompCol_Matrix_dist(SuperMatrix *, SuperMatrix *)
void dfill_dist(double *, int_t, double)
Fills a double precision array with a given value.
Definition: dutil_dist.c:512
void pdgssvx(superlu_dist_options_t *, SuperMatrix *, dScalePermstruct_t *, double *, int, int, gridinfo_t *, dLUstruct_t *, dSOLVEstruct_t *, double *, SuperLUStat_t *, int *)
#define NBUFFERS
Definition: superlu_defs.h:194
int int_t
Definition: superlu_defs.h:114
DiagScale_t
Definition: superlu_enum_consts.h:35
fact_t
Definition: superlu_enum_consts.h:30
Mtype_t
Definition: supermatrix.h:42
Dtype_t
Definition: supermatrix.h:35
Stype_t
Definition: supermatrix.h:22
integer, parameter, public lsub
Definition: superlupara.f90:35
integer, parameter, public trans
Definition: superlupara.f90:35
integer, parameter, public factored
Definition: superlupara.f90:35
integer, parameter, public lusup
Definition: superlupara.f90:35
integer, parameter, public usub
Definition: superlupara.f90:35
Definition: superlu_defs.h:1142
Definition: superlu_defs.h:490
Definition: superlu_defs.h:435
Definition: superlu_ddefs.h:329
int_t * Lblock_dirty_bit
Definition: superlu_ddefs.h:335
int_t nGPUStreams
Definition: superlu_ddefs.h:353
double * lookAhead_L_buff
Definition: superlu_ddefs.h:336
int_t buffer_size
Definition: superlu_ddefs.h:349
int_t superlu_acc_offload
Definition: superlu_ddefs.h:352
double * bigU_host
Definition: superlu_ddefs.h:345
int_t bigu_size
Definition: superlu_ddefs.h:350
int_t Rnbrow
Definition: superlu_ddefs.h:347
Ublock_info_t * Ublock_info
Definition: superlu_ddefs.h:331
int_t RemainBlk
Definition: superlu_ddefs.h:338
int_t first_l_block_acc
Definition: superlu_ddefs.h:333
double * bigU_Phi
Definition: superlu_ddefs.h:344
int_t num_u_blks
Definition: superlu_ddefs.h:341
int_t offloadCondition
Definition: superlu_ddefs.h:351
int_t last_offload
Definition: superlu_ddefs.h:334
int_t Lnbrow
Definition: superlu_ddefs.h:346
int_t ldu
Definition: superlu_ddefs.h:340
Remain_info_t * lookAhead_info
Definition: superlu_ddefs.h:330
int_t nsupers
Definition: superlu_ddefs.h:339
int_t lookAheadBlk
Definition: superlu_ddefs.h:337
int_t jj_cpu
Definition: superlu_ddefs.h:343
Definition: supermatrix.h:194
Definition: supermatrix.h:176
Definition: psymbfact.h:57
Definition: superlu_defs.h:770
Definition: util_dist.h:172
Definition: util_dist.h:95
Definition: supermatrix.h:54
Definition: superlu_defs.h:760
Definition: superlu_defs.h:752
Definition: superlu_defs.h:924
Definition: superlu_ddefs.h:357
int_t * Usub_buf
Definition: superlu_ddefs.h:360
int_t * Lsub_buf
Definition: superlu_ddefs.h:358
double * Lval_buf
Definition: superlu_ddefs.h:359
double * Uval_buf
Definition: superlu_ddefs.h:361
Definition: superlu_ddefs.h:254
Glu_persist_t * Glu_persist
Definition: superlu_ddefs.h:256
int_t * etree
Definition: superlu_ddefs.h:255
char dt
Definition: superlu_ddefs.h:258
dLocalLU_t * Llu
Definition: superlu_ddefs.h:257
Definition: superlu_ddefs.h:97
int * ToRecv
Definition: superlu_ddefs.h:159
int * frecv
Definition: superlu_ddefs.h:167
int ** ToSendR
Definition: superlu_ddefs.h:161
Ucb_indptr_t ** Ucb_indptr
Definition: superlu_ddefs.h:203
long int * Lnzval_bc_offset
Definition: superlu_ddefs.h:105
long int * Lrowind_bc_offset
Definition: superlu_ddefs.h:100
int_t nroot
Definition: superlu_ddefs.h:200
long int Ucb_indcnt
Definition: superlu_ddefs.h:206
long int Unzval_br_cnt
Definition: superlu_ddefs.h:132
int_t * Ucb_valdat
Definition: superlu_ddefs.h:209
int_t * Lindval_loc_bc_dat
Definition: superlu_ddefs.h:114
double * Linv_bc_dat
Definition: superlu_ddefs.h:109
long int * Lindval_loc_bc_offset
Definition: superlu_ddefs.h:115
long int Uinv_bc_cnt
Definition: superlu_ddefs.h:122
int_t ldalsum
Definition: superlu_ddefs.h:180
long int Lindval_loc_bc_cnt
Definition: superlu_ddefs.h:116
int_t ** Lrowind_bc_ptr
Definition: superlu_ddefs.h:98
int_t FRECV
Definition: superlu_ddefs.h:191
long int Ucb_valcnt
Definition: superlu_ddefs.h:211
long int * Uinv_bc_offset
Definition: superlu_ddefs.h:121
int_t n
Definition: superlu_ddefs.h:214
double * Uinv_bc_dat
Definition: superlu_ddefs.h:120
int nfrecvx
Definition: superlu_ddefs.h:168
long int Linv_bc_cnt
Definition: superlu_ddefs.h:111
long int Lnzval_bc_cnt
Definition: superlu_ddefs.h:106
int_t inv
Definition: superlu_ddefs.h:217
int nbrecvx
Definition: superlu_ddefs.h:173
long int * Unzval_br_offset
Definition: superlu_ddefs.h:131
int_t nfrecvmod
Definition: superlu_ddefs.h:216
int * bmod
Definition: superlu_ddefs.h:170
int ** fsendx_plist
Definition: superlu_ddefs.h:166
C_Tree * URtree_ptr
Definition: superlu_ddefs.h:138
int_t * ilsum
Definition: superlu_ddefs.h:178
int nfsendx
Definition: superlu_ddefs.h:169
int_t * ut_modbit
Definition: superlu_ddefs.h:201
int * ToSendD
Definition: superlu_ddefs.h:160
double * ujrow
Definition: superlu_ddefs.h:149
int_t SolveMsgVol
Definition: superlu_ddefs.h:182
int_t * Lrowind_bc_dat
Definition: superlu_ddefs.h:99
C_Tree * LRtree_ptr
Definition: superlu_ddefs.h:136
int_t n_utrecvmod
Definition: superlu_ddefs.h:199
int_t * Unnz
Definition: superlu_ddefs.h:117
C_Tree * LBtree_ptr
Definition: superlu_ddefs.h:135
long int Lrowind_bc_cnt
Definition: superlu_ddefs.h:101
int_t nleaf
Definition: superlu_ddefs.h:215
long int Ufstnz_br_cnt
Definition: superlu_ddefs.h:127
int nbsendx
Definition: superlu_ddefs.h:174
double ** Uinv_bc_ptr
Definition: superlu_ddefs.h:119
int_t ** Ucb_valptr
Definition: superlu_ddefs.h:208
int_t ut_ldalsum
Definition: superlu_ddefs.h:192
double ** Unzval_br_ptr
Definition: superlu_ddefs.h:129
int_t UT_SOLVE
Definition: superlu_ddefs.h:189
int_t ** Ufstnz_br_ptr
Definition: superlu_ddefs.h:124
int_t ** ut_sendx_plist
Definition: superlu_ddefs.h:195
int_t * utrecv
Definition: superlu_ddefs.h:196
int * fmod
Definition: superlu_ddefs.h:165
int_t * utmod
Definition: superlu_ddefs.h:194
int_t * Urbs
Definition: superlu_ddefs.h:202
int_t ** Lrowind_bc_2_lsum
Definition: superlu_ddefs.h:118
int_t * ut_ilsum
Definition: superlu_ddefs.h:193
double * Lnzval_bc_dat
Definition: superlu_ddefs.h:104
long int * Ucb_indoffset
Definition: superlu_ddefs.h:205
long int * Linv_bc_offset
Definition: superlu_ddefs.h:110
int_t n_utrecvx
Definition: superlu_ddefs.h:198
double * Unzval_br_dat
Definition: superlu_ddefs.h:130
double ** Lnzval_bc_ptr
Definition: superlu_ddefs.h:103
C_Tree * UBtree_ptr
Definition: superlu_ddefs.h:137
double ** Linv_bc_ptr
Definition: superlu_ddefs.h:108
int * mod_bit
Definition: superlu_ddefs.h:175
int ** bsendx_plist
Definition: superlu_ddefs.h:171
int_t * Ufstnz_br_dat
Definition: superlu_ddefs.h:125
int_t SolveMsgSent
Definition: superlu_ddefs.h:181
int_t n_utsendx
Definition: superlu_ddefs.h:197
long int * Ufstnz_br_offset
Definition: superlu_ddefs.h:126
long int * Ucb_valoffset
Definition: superlu_ddefs.h:210
int * brecv
Definition: superlu_ddefs.h:172
int_t ** Lindval_loc_bc_ptr
Definition: superlu_ddefs.h:113
Ucb_indptr_t * Ucb_inddat
Definition: superlu_ddefs.h:204
int_t L_SOLVE
Definition: superlu_ddefs.h:190
Definition: superlu_ddefs.h:284
int_t * xrow_to_proc
Definition: superlu_ddefs.h:295
pdgsmv_comm_t * gsmv_comm
Definition: superlu_ddefs.h:288
pxgstrs_comm_t * gstrs_comm
Definition: superlu_ddefs.h:290
int_t * A_colind_gsmv
Definition: superlu_ddefs.h:291
NRformat_loc3d * A3d
Definition: superlu_ddefs.h:296
int_t * row_to_proc
Definition: superlu_ddefs.h:285
int_t * inv_perm_c
Definition: superlu_ddefs.h:286
int_t * diag_len
Definition: superlu_ddefs.h:287
Definition: superlu_ddefs.h:76
int_t * perm_r
Definition: superlu_ddefs.h:80
DiagScale_t DiagScale
Definition: superlu_ddefs.h:77
double * R
Definition: superlu_ddefs.h:78
int_t * perm_c
Definition: superlu_ddefs.h:81
double * C
Definition: superlu_ddefs.h:79
Definition: superlu_ddefs.h:391
double * BlockLFactor
Definition: superlu_ddefs.h:392
double * BlockUFactor
Definition: superlu_ddefs.h:393
Definition: superlu_ddefs.h:385
double * bigU
Definition: superlu_ddefs.h:386
double * bigV
Definition: superlu_ddefs.h:387
Definition: superlu_defs.h:937
Definition: superlu_defs.h:839
Definition: superlu_defs.h:890
Definition: superlu_defs.h:398
Definition: superlu_defs.h:388
Definition: superlu_ddefs.h:317
int_t * lsub
Definition: superlu_ddefs.h:318
int_t nsupr
Definition: superlu_ddefs.h:322
int_t luptr0
Definition: superlu_ddefs.h:320
int_t nlb
Definition: superlu_ddefs.h:321
double * lusup
Definition: superlu_ddefs.h:319
Definition: superlu_defs.h:947
Definition: superlu_ddefs.h:397
uPanelInfo_t * uPanelInfo
Definition: superlu_ddefs.h:400
Ublock_info_t * Ublock_info
Definition: superlu_ddefs.h:398
Remain_info_t * Remain_info
Definition: superlu_ddefs.h:399
lPanelInfo_t * lPanelInfo
Definition: superlu_ddefs.h:401
Definition: superlu_ddefs.h:263
int_t TotalIndSend
Definition: superlu_ddefs.h:277
int_t * extern_start
Definition: superlu_ddefs.h:264
double * val_torecv
Definition: superlu_ddefs.h:276
int_t * ind_torecv
Definition: superlu_ddefs.h:266
double * val_tosend
Definition: superlu_ddefs.h:275
int_t * ptr_ind_tosend
Definition: superlu_ddefs.h:267
int * SendCounts
Definition: superlu_ddefs.h:271
int_t * ptr_ind_torecv
Definition: superlu_ddefs.h:269
int_t TotalValSend
Definition: superlu_ddefs.h:279
int_t * ind_tosend
Definition: superlu_ddefs.h:265
int * RecvCounts
Definition: superlu_ddefs.h:273
Definition: superlu_defs.h:551
Definition: superlu_defs.h:901
Definition: superlu_defs.h:744
Definition: superlu_defs.h:712
Definition: superlu_defs.h:882
Definition: superlu_ddefs.h:372
int_t ** treePerm
Definition: superlu_ddefs.h:378
gEtreeInfo_t gEtreeInfo
Definition: superlu_ddefs.h:373
int_t * myZeroTrIdxs
Definition: superlu_ddefs.h:377
dLUValSubBuf_t * LUvsb
Definition: superlu_ddefs.h:381
sForest_t ** sForests
Definition: superlu_ddefs.h:379
int_t * iperm_c_supno
Definition: superlu_ddefs.h:374
int_t * myNodeCount
Definition: superlu_ddefs.h:375
int_t * supernode2treeMap
Definition: superlu_ddefs.h:380
int_t * myTreeIdxs
Definition: superlu_ddefs.h:376
Definition: superlu_ddefs.h:308
int_t * usub
Definition: superlu_ddefs.h:312
int_t nub
Definition: superlu_ddefs.h:309
double * uval
Definition: superlu_ddefs.h:313
int_t klst
Definition: superlu_ddefs.h:310
int_t ldu
Definition: superlu_ddefs.h:311
Definitions which are precision-neutral.