SuperLU Distributed 8.2.1
Distributed memory sparse direct solver
superlu_sdefs.h
Go to the documentation of this file.
1
28#ifndef __SUPERLU_SDEFS /* allow multiple inclusions */
29#define __SUPERLU_SDEFS
30
31/*
32 * File name: superlu_sdefs.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 float*, 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 float*, 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 float *R;
79 float *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 float **Lnzval_bc_ptr; /* size ceil(NSUPERS/Pc) */
104 float *Lnzval_bc_dat; /* size sum of sizes of Lnzval_bc_ptr[lk]) */
105 long int *Lnzval_bc_offset; /* size ceil(NSUPERS/Pc) */
107
108 float **Linv_bc_ptr; /* size ceil(NSUPERS/Pc) */
109 float *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 float **Uinv_bc_ptr; /* size ceil(NSUPERS/Pc) */
120 float *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 float **Unzval_br_ptr; /* size ceil(NSUPERS/Pr) */
130 float *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 float *Lval_buf; /* Buffer for the remote nonzeros of L */
142 int_t *Usub_buf; /* Buffer for the remote subscripts of U */
143 float *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 float *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 float *Uval_buf_2[MAX_LOOKAHEADS]; /* Buffer for the remote nonzeros of U */
149 float *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 float *d_Lnzval_bc_dat;
224 long int *d_Lnzval_bc_offset;
225 float *d_Linv_bc_dat ;
226 float *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 float *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} sLocalLU_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 float *val_tosend; /* X values to be sent to other processes */
276 float *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 psgsmv_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 psgsmv_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 psgsmv() */
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{
309 int_t nub;
310 int_t klst;
311 int_t ldu;
312 int_t* usub;
313 float* uval;
315
316typedef struct
317{
318 int_t *lsub;
319 float *lusup;
320 int_t luptr0;
321 int_t nlb; //number of l blocks
322 int_t nsupr;
324
325
326
327/* HyP_t is the data structure to assist HALO offload of Schur-complement. */
328typedef struct
329{
330 Remain_info_t *lookAhead_info, *Remain_info;
331 Ublock_info_t *Ublock_info, *Ublock_info_Phi;
332
333 int_t first_l_block_acc , first_u_block_acc;
334 int_t last_offload ;
335 int_t *Lblock_dirty_bit, * Ublock_dirty_bit;
336 float *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
343 int_t jj_cpu;
344 float *bigU_Phi;
345 float *bigU_host;
346 int_t Lnbrow;
347 int_t Rnbrow;
348
349 int_t buffer_size;
350 int_t bigu_size;
351 int_t offloadCondition;
352 int_t superlu_acc_offload;
353 int_t nGPUStreams;
354} HyP_t;
355
356typedef struct
357{
359 float * Lval_buf ;
361 float * Uval_buf ;
363
365 int_t knsupc,
366 HyP_t* HyP,
367 SCT_t* SCT,
368 SuperLUStat_t *stat
369 );
370
371typedef struct
372{
373 gEtreeInfo_t gEtreeInfo;
374 int_t* iperm_c_supno;
375 int_t* myNodeCount;
376 int_t* myTreeIdxs;
377 int_t* myZeroTrIdxs;
378 int_t** treePerm;
379 sForest_t** sForests;
380 int_t* supernode2treeMap;
383
384typedef struct
385{
386 float *bigU;
387 float *bigV;
388} sscuBufs_t;
389
390typedef struct
391{
395
396typedef struct
397{
398 Ublock_info_t* Ublock_info;
399 Remain_info_t* Remain_info;
400 uPanelInfo_t* uPanelInfo;
401 lPanelInfo_t* lPanelInfo;
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, float *, int_t *, int_t *,
424extern void
426 float **, 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 float *, int_t);
442
443extern void sallocateA_dist (int_t, int_t, float **, int_t **, int_t **);
444extern void sGenXtrue_dist (int_t, int_t, float *, int_t);
445extern void sFillRHS_dist (char *, int_t, float *, int_t,
446 SuperMatrix *, float *, int_t);
447extern int screate_matrix(SuperMatrix *, int, float **, int *,
448 float **, int *, FILE *, gridinfo_t *);
449extern int screate_matrix_rb(SuperMatrix *, int, float **, int *,
450 float **, int *, FILE *, gridinfo_t *);
451extern int screate_matrix_dat(SuperMatrix *, int, float **, int *,
452 float **, int *, FILE *, gridinfo_t *);
453extern int screate_matrix_postfix(SuperMatrix *, int, float **, int *,
454 float **, int *, FILE *, char *, gridinfo_t *);
455
456extern void sScalePermstructInit(const int_t, const int_t,
459
460/* Driver related */
461extern void sgsequ_dist (SuperMatrix *, float *, float *, float *,
462 float *, float *, int_t *);
463extern float slangs_dist (char *, SuperMatrix *);
464extern void slaqgs_dist (SuperMatrix *, float *, float *, float,
465 float, float, char *);
466extern void psgsequ (SuperMatrix *, float *, float *, float *,
467 float *, float *, int_t *, gridinfo_t *);
468extern float pslangs (char *, SuperMatrix *, gridinfo_t *);
469extern void pslaqgs (SuperMatrix *, float *, float *, float,
470 float, float, char *);
471extern int psPermute_Dense_Matrix(int_t, int_t, int_t [], int_t[],
472 float [], int, float [], int, int,
473 gridinfo_t *);
474
475extern int sp_strsv_dist (char *, char *, char *, SuperMatrix *,
476 SuperMatrix *, float *, int *);
477extern int sp_sgemv_dist (char *, float, SuperMatrix *, float *,
478 int, float, float *, int);
479extern int sp_sgemm_dist (char *, int, float, SuperMatrix *,
480 float *, int, float, float *, int);
481
485 sScalePermstruct_t *, float *,
486 int, int, gridinfo_t *, sLUstruct_t *, float *,
487 SuperLUStat_t *, int *);
488extern float psdistribute(fact_t, int_t, SuperMatrix *,
492 sScalePermstruct_t *, float *,
493 int, int, gridinfo_t *, sLUstruct_t *,
494 sSOLVEstruct_t *, float *, SuperLUStat_t *, int *);
501 int_t [], int_t [], gridinfo_t *grid,
503extern void pxgstrs_finalize(pxgstrs_comm_t *);
504extern int sldperm_dist(int, int, int_t, int_t [], int_t [],
505 float [], int_t *, float [], float []);
506extern int sstatic_schedule(superlu_dist_options_t *, int, int,
508 int_t *, int_t *, int *);
509extern void sLUstructInit(const int_t, sLUstruct_t *);
510extern void sLUstructFree(sLUstruct_t *);
511extern void sDestroy_LU(int_t, gridinfo_t *, sLUstruct_t *);
512extern void sDestroy_Tree(int_t, gridinfo_t *, sLUstruct_t *);
513extern void sscatter_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, float *tempv,
516 int* indirect_thread, int* indirect2,
517 int_t ** Lrowind_bc_ptr, float **Lnzval_bc_ptr,
518 gridinfo_t * grid);
519extern void sscatter_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, float* tempv,
522 int_t ** Ufstnz_br_ptr, float **Unzval_br_ptr,
523 gridinfo_t * grid);
524extern int_t psgstrf(superlu_dist_options_t *, int, int, float anorm,
526
527/* #define GPU_PROF
528#define IPM_PROF */
529
530/* Solve related */
532 float *, int_t, int, SuperLUStat_t *, int *);
534 float *, int_t, int_t, int_t, int, sSOLVEstruct_t *,
535 SuperLUStat_t *, int *);
536extern void psgstrf2_trsm(superlu_dist_options_t * options, int_t k0, int_t k,
537 double thresh, Glu_persist_t *, gridinfo_t *,
538 sLocalLU_t *, MPI_Request *, int tag_ub,
539 SuperLUStat_t *, int *info);
540extern void psgstrs2_omp(int_t k0, int_t k, Glu_persist_t *, gridinfo_t *,
542extern int_t psReDistribute_B_to_X(float *B, int_t m_loc, int nrhs, int_t ldb,
543 int_t fst_row, int_t *ilsum, float *x,
546extern void slsum_fmod(float *, float *, float *, float *,
547 int, int, int_t , int *fmod, int_t, int_t, int_t,
549 MPI_Request [], SuperLUStat_t *);
550extern void slsum_bmod(float *, float *, float *,
551 int, int_t, int *bmod, int_t *, Ucb_indptr_t **,
552 int_t **, int_t *, gridinfo_t *, sLocalLU_t *,
553 MPI_Request [], SuperLUStat_t *);
554
555extern void slsum_fmod_inv(float *, float *, float *, float *,
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 slsum_fmod_inv_master(float *, float *, float *, float *,
560 int, int, int_t , int *fmod, int_t,
562 SuperLUStat_t **, int_t, int_t, int_t, int_t, int, int);
563extern void slsum_bmod_inv(float *, float *, float *, float *,
564 int, int_t, int *bmod, int_t *, Ucb_indptr_t **,
565 int_t **, int_t *, gridinfo_t *, sLocalLU_t *,
566 SuperLUStat_t **, int_t *, int_t *, int_t, int_t, int, int);
567extern void slsum_bmod_inv_master(float *, float *, float *, float *,
568 int, int_t, int *bmod, int_t *, Ucb_indptr_t **,
569 int_t **, int_t *, gridinfo_t *, sLocalLU_t *,
570 SuperLUStat_t **, int_t, int_t, int, int);
571
572extern void sComputeLevelsets(int , int_t , gridinfo_t *,
574
575#ifdef GPU_ACC
576extern void slsum_fmod_inv_gpu_wrap(int_t, int_t, int_t, int_t, float *, float *, int, int, int_t , int *fmod, C_Tree *, C_Tree *, int_t *, int_t *, int64_t *, float *, int64_t *, float *, int64_t *, int_t *, int64_t *, int_t *, gridinfo_t *, float * , float * , int_t );
577extern void dlsum_bmod_inv_gpu_wrap(int_t, int_t, int_t, int_t, float *, float *,int,int, int_t , int *bmod, C_Tree *, C_Tree *, int_t *, int_t *,int_t *, int64_t *, float *, int64_t *, int_t *, int64_t *, Ucb_indptr_t *, int64_t *, float *, int64_t *,int_t *,gridinfo_t *);
578#endif
579
580extern void psgsrfs(int_t, SuperMatrix *, float, sLUstruct_t *,
582 float [], int_t, float [], int_t, int,
583 sSOLVEstruct_t *, float *, SuperLUStat_t *, int *);
584extern void psgsrfs_ABXglobal(int_t, SuperMatrix *, float, sLUstruct_t *,
585 gridinfo_t *, float *, int_t, float *, int_t,
586 int, float *, SuperLUStat_t *, int *);
588 gridinfo_t *, int_t *, int_t *[],
589 float *[], int_t *[], int_t []);
590extern int psgsmv_AXglobal(int_t, int_t [], float [], int_t [],
591 float [], float []);
592extern int psgsmv_AXglobal_abs(int_t, int_t [], float [], int_t [],
593 float [], float []);
594extern void psgsmv_init(SuperMatrix *, int_t *, gridinfo_t *,
595 psgsmv_comm_t *);
596extern void psgsmv(int_t, SuperMatrix *, gridinfo_t *, psgsmv_comm_t *,
597 float x[], float ax[]);
598extern void psgsmv_finalize(psgsmv_comm_t *);
599
600/* Memory-related */
601extern float *floatMalloc_dist(int_t);
602extern float *floatCalloc_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 sZeroLblocks(int, int, gridinfo_t *, sLUstruct_t *);
616extern void sZeroUblocks(int iam, int n, gridinfo_t *, sLUstruct_t *);
617extern void sfill_dist (float *, int_t, float);
618extern void sinf_norm_error_dist (int_t, int_t, float*, int_t,
619 float*, int_t, gridinfo_t*);
620extern void psinf_norm_error(int, int_t, int_t, float [], int_t,
621 float [], int_t , MPI_Comm);
622extern void sreadhb_dist (int, FILE *, int_t *, int_t *, int_t *,
623 float **, int_t **, int_t **);
624extern void sreadtriple_dist(FILE *, int_t *, int_t *, int_t *,
625 float **, int_t **, int_t **);
626extern void sreadtriple_noheader(FILE *, int_t *, int_t *, int_t *,
627 float **, int_t **, int_t **);
628extern void sreadrb_dist(int, FILE *, int_t *, int_t *, int_t *,
629 float **, int_t **, int_t **);
630extern void sreadMM_dist(FILE *, int_t *, int_t *, int_t *,
631 float **, int_t **, int_t **);
632extern int sread_binary(FILE *, int_t *, int_t *, int_t *,
633 float **, int_t **, int_t **);
634
635/* Distribute the data for numerical factorization */
639extern void psGetDiagU(int_t, sLUstruct_t *, gridinfo_t *, float *);
640
642
643/* Routines for debugging */
644extern void sPrintLblocks(int, int_t, gridinfo_t *, Glu_persist_t *,
645 sLocalLU_t *);
646extern void sPrintUblocks(int, int_t, gridinfo_t *, Glu_persist_t *,
647 sLocalLU_t *);
652extern void Printfloat5(char *, int_t, float *);
653extern int file_Printfloat5(FILE *, char *, int_t, float *);
654
655extern void sGenCOOLblocks(int, int_t, gridinfo_t*,
656 Glu_persist_t*, sLocalLU_t *, int_t** , int_t** , float ** , int_t* , int_t* );
657extern void sGenCSCLblocks(int, int_t, gridinfo_t*,
658 Glu_persist_t*, sLocalLU_t *, float **, int_t **, int_t **, int_t*, int_t*);
659extern void sGenCSRLblocks(int, int_t, gridinfo_t*,
660 Glu_persist_t*, sLocalLU_t *, float **, int_t **, int_t **, int_t*, int_t*);
661
662
663/* BLAS */
664
665#ifdef USE_VENDOR_BLAS
666extern void sgemm_(const char*, const char*, const int*, const int*, const int*,
667 const float*, const float*, const int*, const float*,
668 const int*, const float*, float*, const int*, int, int);
669extern void strsv_(char*, char*, char*, int*, float*, int*,
670 float*, int*, int, int, int);
671extern void strsm_(const char*, const char*, const char*, const char*,
672 const int*, const int*, const float*, const float*, const int*,
673 float*, const int*, int, int, int, int);
674extern void sgemv_(const char *, const int *, const int *, const float *,
675 const float *a, const int *, const float *, const int *,
676 const float *, float *, const int *, int);
677
678#else
679extern int sgemm_(const char*, const char*, const int*, const int*, const int*,
680 const float*, const float*, const int*, const float*,
681 const int*, const float*, float*, const int*);
682extern int strsv_(char*, char*, char*, int*, float*, int*,
683 float*, int*);
684extern int strsm_(const char*, const char*, const char*, const char*,
685 const int*, const int*, const float*, const float*, const int*,
686 float*, const int*);
687extern void sgemv_(const char *, const int *, const int *, const float *,
688 const float *a, const int *, const float *, const int *,
689 const float *, float *, const int *);
690#endif
691
692extern void sger_(const int*, const int*, const float*,
693 const float*, const int*, const float*, const int*,
694 float*, const int*);
695
696extern int sscal_(const int *n, const float *alpha, float *dx, const int *incx);
697extern int saxpy_(const int *n, const float *alpha, const float *x,
698 const int *incx, float *y, const int *incy);
699
700/* SuperLU BLAS interface: ssuperlu_blas.c */
701extern int superlu_sgemm(const char *transa, const char *transb,
702 int m, int n, int k, float alpha, float *a,
703 int lda, float *b, int ldb, float beta, float *c, int ldc);
704extern int superlu_strsm(const char *sideRL, const char *uplo,
705 const char *transa, const char *diag, const int m, const int n,
706 const float alpha, const float *a,
707 const int lda, float *b, const int ldb);
708extern int superlu_sger(const int m, const int n, const float alpha,
709 const float *x, const int incx, const float *y,
710 const int incy, float *a, const int lda);
711extern int superlu_sscal(const int n, const float alpha, float *x, const int incx);
712extern int superlu_saxpy(const int n, const float alpha,
713 const float *x, const int incx, float *y, const int incy);
714extern int superlu_sgemv(const char *trans, const int m,
715 const int n, const float alpha, const float *a,
716 const int lda, const float *x, const int incx,
717 const float beta, float *y, const int incy);
718extern int superlu_strsv(char *uplo, char *trans, char *diag,
719 int n, float *a, int lda, float *x, int incx);
720
721#ifdef SLU_HAVE_LAPACK
722extern void strtri_(char*, char*, int*, float*, int*, int*);
723#endif
724
725/*==== For 3D code ====*/
726extern int screate_matrix3d(SuperMatrix *A, int nrhs, float **rhs,
727 int *ldb, float **x, int *ldx,
728 FILE *fp, gridinfo3d_t *grid3d);
729extern int screate_matrix_postfix3d(SuperMatrix *A, int nrhs, float **rhs,
730 int *ldb, float **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 sGatherNRformat_loc3d(fact_t Fact, NRformat_loc *A, float *B,
736 int ldb, int nrhs, gridinfo3d_t *grid3d,
737 NRformat_loc3d **);
738extern int sScatter_B3d(NRformat_loc3d *A3d, gridinfo3d_t *grid3d);
739
741 sScalePermstruct_t *, float B[], int ldb, int nrhs,
743 float *berr, SuperLUStat_t *, int *info);
744extern int_t psgstrf3d(superlu_dist_options_t *, int m, int n, float anorm,
746 gridinfo3d_t *, SuperLUStat_t *, int *);
747extern void sInit_HyP(HyP_t* HyP, sLocalLU_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, float *L_mat, int ldl,
755 float *U_mat, int ldu, float *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, float **Lnzval_bc_ptr,
762 int_t **Ufstnz_br_ptr, float **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
772sblock_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 float *L_mat, int_t ldl, float *U_mat, int_t ldu,
775 float *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, float **Lnzval_bc_ptr,
782 int_t **Ufstnz_br_ptr, float **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, sLUstruct_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, sLUstruct_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, sLUstruct_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, sLUstruct_t *, gridinfo_t*,
818 SCT_t*SCT, SuperLUStat_t * );
819
820 /* from gather.h */
821extern void sgather_u(int_t num_u_blks,
822 Ublock_info_t *Ublock_info, int_t * usub,
823 float *uval, float *bigU, int_t ldu,
824 int_t *xsup, int_t klst /* for SuperSize */
825 );
826
827extern void sgather_l( int_t num_LBlk, int_t knsupc,
828 Remain_info_t *L_info,
829 float * lval, int_t LD_lval,
830 float * L_buff );
831
832extern void sRgather_L(int_t k, int_t *lsub, float *lusup, gEtreeInfo_t*,
834 int_t *myIperm, int_t *iperm_c_supno );
835extern void sRgather_U(int_t k, int_t jj0, int_t *usub, float *uval,
836 float *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 sLUstruct_t *LUstruct, gridinfo3d_t * grid3d);
844extern void sDestroy_trf3Dpartition(trf3Dpartition_t *trf3Dpartition, gridinfo3d_t *grid3d);
845
846extern void s3D_printMemUse(trf3Dpartition_t* trf3Dpartition,
847 sLUstruct_t *LUstruct, gridinfo3d_t * grid3d);
848
849//extern int* getLastDep(gridinfo_t *grid, SuperLUStat_t *stat,
850// superlu_dist_options_t *options, sLocalLU_t *Llu,
851// int_t* xsup, int_t num_look_aheads, int_t nsupers,
852// int_t * iperm_c_supno);
853
854extern void sinit3DLUstructForest( int_t* myTreeIdxs, int_t* myZeroTrIdxs,
855 sForest_t** sForests, sLUstruct_t* LUstruct,
856 gridinfo3d_t* grid3d);
857
858extern int_t sgatherAllFactoredLUFr(int_t* myZeroTrIdxs, sForest_t* sForests,
859 sLUstruct_t* LUstruct, gridinfo3d_t* grid3d,
860 SCT_t* SCT );
861
862 /* The following are from pdgstrf2.h */
863extern int_t sLpanelUpdate(int_t off0, int_t nsupc, float* ublk_ptr,
864 int_t ld_ujrow, float* lusup, int_t nsupr, SCT_t*);
866 double thresh, float *BlockUFactor, Glu_persist_t *,
868 SuperLUStat_t *, int *info, SCT_t*);
869extern int_t sTrs2_GatherU(int_t iukp, int_t rukp, int_t klst,
870 int_t nsupc, int_t ldu, int_t *usub,
871 float* uval, float *tempv);
872extern int_t sTrs2_ScatterU(int_t iukp, int_t rukp, int_t klst,
873 int_t nsupc, int_t ldu, int_t *usub,
874 float* uval, float *tempv);
875extern int_t sTrs2_GatherTrsmScatter(int_t klst, int_t iukp, int_t rukp,
876 int_t *usub, float* uval, float *tempv,
877 int_t knsupc, int nsupr, float* lusup,
878 Glu_persist_t *Glu_persist) ;
879extern void psgstrs2
880#ifdef _CRAY
881(
882 int_t m, int_t k0, int_t k, Glu_persist_t *Glu_persist, gridinfo_t *grid,
883 sLocalLU_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 sLocalLU_t *Llu, SuperLUStat_t *stat
889);
890#endif
891
892extern void psgstrf2(superlu_dist_options_t *, int_t nsupers, int_t k0,
893 int_t k, double thresh, Glu_persist_t *, gridinfo_t *,
894 sLocalLU_t *, MPI_Request *, int, SuperLUStat_t *, int *);
895
896 /* from p3dcomm.h */
897extern int_t sAllocLlu_3d(int_t nsupers, sLUstruct_t * LUstruct, gridinfo3d_t* grid3d);
898extern int_t sp3dScatter(int_t n, sLUstruct_t * LUstruct, gridinfo3d_t* grid3d);
900 sLUstruct_t * LUstruct, gridinfo3d_t* grid3d);
902 sLUstruct_t * LUstruct, gridinfo3d_t* grid3d);
903extern int_t scollect3dLpanels(int_t layer, int_t nsupers, sLUstruct_t * LUstruct, gridinfo3d_t* grid3d);
904extern int_t scollect3dUpanels(int_t layer, int_t nsupers, sLUstruct_t * LUstruct, gridinfo3d_t* grid3d);
905extern int_t sp3dCollect(int_t layer, int_t n, sLUstruct_t * LUstruct, gridinfo3d_t* grid3d);
906/*zero out LU non zero entries*/
907extern int_t szeroSetLU(int_t nnodes, int_t* nodeList , sLUstruct_t *, gridinfo3d_t*);
908extern int sAllocGlu_3d(int_t n, int_t nsupers, sLUstruct_t *);
910extern int sDeAllocGlu_3d(sLUstruct_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 float* Lval_buf, float* Uval_buf,
919 sLUstruct_t* LUstruct, gridinfo3d_t* grid3d, SCT_t* SCT);
920/*reduces all nodelists required in a level*/
921extern int sreduceAllAncestors3d(int_t ilvl, int_t* myNodeCount,
922 int_t** treePerm,
923 sLUValSubBuf_t* LUvsb,
924 sLUstruct_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, sLUValSubBuf_t* LUvsb,
934 sLUstruct_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 sinit3DLUstruct( int_t* myTreeIdxs, int_t* myZeroTrIdxs,
944 int_t* nodeCount, int_t** nodeList,
945 sLUstruct_t* LUstruct, gridinfo3d_t* grid3d);
946
948 sLUstruct_t* LUstruct, gridinfo3d_t* grid3d, SCT_t* SCT);
949int_t szRecvLPanel(int_t k, int_t sender, float alpha,
950 float beta, float* Lval_buf,
951 sLUstruct_t* LUstruct, gridinfo3d_t* grid3d, SCT_t* SCT);
953 sLUstruct_t* LUstruct, gridinfo3d_t* grid3d, SCT_t* SCT);
954int_t szRecvUPanel(int_t k, int_t sender, float alpha,
955 float beta, float* Uval_buf,
956 sLUstruct_t* LUstruct, gridinfo3d_t* grid3d, SCT_t* SCT);
957
958 /* from communication_aux.h */
960 gridinfo_t *, int* msgcnt, MPI_Request *,
961 int **ToSendR, int_t *xsup, int );
963 gridinfo_t *, int* msgcnt, int **ToSendR,
964 int_t *xsup , SCT_t*, int);
965extern int_t sIBcast_UPanel(int_t k, int_t k0, int_t* usub, float* uval,
966 gridinfo_t *, int* msgcnt, MPI_Request *,
967 int *ToSendD, int );
968extern int_t sBcast_UPanel(int_t k, int_t k0, int_t* usub, float* uval,
969 gridinfo_t *, int* msgcnt, int *ToSendD, SCT_t*, int);
970extern int_t sIrecv_LPanel (int_t k, int_t k0, int_t* Lsub_buf,
971 float* Lval_buf, gridinfo_t *,
972 MPI_Request *, sLocalLU_t *, int);
973extern int_t sIrecv_UPanel(int_t k, int_t k0, int_t* Usub_buf, float*,
974 sLocalLU_t *, gridinfo_t*, MPI_Request *, int);
975extern int_t sWait_URecv(MPI_Request *, int* msgcnt, SCT_t *);
976extern int_t sWait_LRecv(MPI_Request*, int* msgcnt, int* msgcntsU,
977 gridinfo_t *, SCT_t*);
978extern int_t sISend_UDiagBlock(int_t k0, float *ublk_ptr, int_t size,
979 MPI_Request *, gridinfo_t *, int);
980extern int_t sRecv_UDiagBlock(int_t k0, float *ublk_ptr, int_t size,
981 int_t src, gridinfo_t *, SCT_t*, int);
982extern int_t sPackLBlock(int_t k, float* Dest, Glu_persist_t *,
984extern int_t sISend_LDiagBlock(int_t k0, float *lblk_ptr, int_t size,
985 MPI_Request *, gridinfo_t *, int);
986extern int_t sIRecv_UDiagBlock(int_t k0, float *ublk_ptr, int_t size,
987 int_t src, MPI_Request *, gridinfo_t *,
988 SCT_t*, int);
989extern int_t sIRecv_LDiagBlock(int_t k0, float *L_blk_ptr, int_t size,
990 int_t src, MPI_Request *, gridinfo_t*, SCT_t*, int);
991extern int_t sUDiagBlockRecvWait( int_t k, int_t* IrecvPlcd_D, int_t* factored_L,
992 MPI_Request *, gridinfo_t *, sLUstruct_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 sIBcast_UDiagBlock(int_t k, float *ublk_ptr, int_t size,
997 MPI_Request *, gridinfo_t *);
998extern int_t sIBcast_LDiagBlock(int_t k, float *lblk_ptr, int_t size,
999 MPI_Request *, gridinfo_t *);
1000#endif
1001
1002 /* from trfCommWrapper.h */
1004 float *BlockUFactor, float *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 sLUstruct_t *LUstruct, SuperLUStat_t *, int *info,
1009 SCT_t *, int tag_ub);
1010extern int_t sUPanelTrSolve( int_t k, float* BlockLFactor, float* bigV,
1013extern int_t sLPanelUpdate(int_t k, int_t* IrecvPlcd_D, int_t* factored_L,
1014 MPI_Request *, float* BlockUFactor, gridinfo_t *,
1015 sLUstruct_t *, SCT_t *);
1016extern int_t sUPanelUpdate(int_t k, int_t* factored_U, MPI_Request *,
1017 float* BlockLFactor, float* bigV,
1020extern int_t sIBcastRecvLPanel(int_t k, int_t k0, int* msgcnt,
1021 MPI_Request *, MPI_Request *,
1022 int_t* Lsub_buf, float* Lval_buf,
1024 SCT_t *, int tag_ub);
1025extern int_t sIBcastRecvUPanel(int_t k, int_t k0, int* msgcnt, MPI_Request *,
1026 MPI_Request *, int_t* Usub_buf, float* Uval_buf,
1027 gridinfo_t *, sLUstruct_t *, SCT_t *, int tag_ub);
1028extern int_t sWaitL(int_t k, int* msgcnt, int* msgcntU, MPI_Request *,
1029 MPI_Request *, gridinfo_t *, sLUstruct_t *, SCT_t *);
1030extern int_t sWaitU(int_t k, int* msgcnt, MPI_Request *, MPI_Request *,
1031 gridinfo_t *, sLUstruct_t *, SCT_t *);
1032extern int_t sLPanelTrSolve(int_t k, int_t* factored_L, float* BlockUFactor,
1033 gridinfo_t *, sLUstruct_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 float *bigU, int_t* Lsub_buf,
1043 float* Lval_buf, int_t* Usub_buf,
1044 float* Uval_buf, gridinfo_t *, sLUstruct_t *);
1048 sLUValSubBuf_t* LUvsb, gridinfo_t *,
1049 sLUstruct_t *, HyP_t*);
1050extern float* sgetBigV(int_t, int_t);
1052// permutation from superLU default
1053
1054 /* from treeFactorization.h */
1056extern int_t sinitScuBufs(int_t ldt, int_t num_threads, int_t nsupers,
1058extern int sfreeScuBufs(sscuBufs_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 sscuBufs_t *scuBufs, // contains buffers for schur complement update
1066 packLUInfo_t*packLUInfo,
1067 msgs_t*msgs,
1068 sLUValSubBuf_t* LUvsb,
1069 sdiagFactBufs_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 sLUstruct_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 sscuBufs_t *scuBufs, // contains buffers for schur complement update
1086 packLUInfo_t*packLUInfo,
1087 msgs_t*msgs,
1088 sLUValSubBuf_t* LUvsb,
1089 sdiagFactBufs_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 sLUstruct_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 sscuBufs_t *scuBufs, // contains buffers for schur complement update
1105 packLUInfo_t*packLUInfo,
1106 msgs_t*msgs,
1107 sLUValSubBuf_t* LUvsb,
1108 sdiagFactBufs_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 sLUstruct_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 sscuBufs_t *scuBufs, // contains buffers for schur complement update
1123 packLUInfo_t*packLUInfo,
1124 msgs_t**msgss, // size=num Look ahead
1125 sLUValSubBuf_t** LUvsbs, // size=num Look ahead
1126 sdiagFactBufs_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 sLUstruct_t *LUstruct, gridinfo3d_t * grid3d, SuperLUStat_t *stat,
1135 double thresh, SCT_t *SCT, int tag_ub,
1136 int *info
1137);
1139extern int sLluBufFreeArr(int_t numLA, sLUValSubBuf_t **LUvsbs);
1141extern int sfreeDiagFactBufsArr(int_t mxLeafNode, sdiagFactBufs_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 sscuBufs_t *scuBufs, // contains buffers for schur complement update
1153 packLUInfo_t*packLUInfo,
1154 msgs_t**msgss, // size=num Look ahead
1155 sLUValSubBuf_t** LUvsbs, // size=num Look ahead
1156 sdiagFactBufs_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 sLUstruct_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 strtri_
Definition: superlu_FCnames.h:160
#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
int screate_matrix_postfix(SuperMatrix *, int, float **, int *, float **, int *, FILE *, char *, gridinfo_t *)
Definition: screate_matrix.c:76
int_t szRecvLPanel(int_t k, int_t sender, float alpha, float beta, float *Lval_buf, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
void Printfloat5(char *, int_t, float *)
Definition: sutil_dist.c:543
void psgssvx(superlu_dist_options_t *, SuperMatrix *, sScalePermstruct_t *, float *, int, int, gridinfo_t *, sLUstruct_t *, sSOLVEstruct_t *, float *, SuperLUStat_t *, int *)
int_t initPackLUInfo(int_t nsupers, packLUInfo_t *packLUInfo)
Definition: treeFactorization.c:367
int screate_matrix(SuperMatrix *, int, float **, int *, float **, int *, FILE *, gridinfo_t *)
Definition: screate_matrix.c:348
int screate_matrix_postfix3d(SuperMatrix *A, int nrhs, float **rhs, int *ldb, float **x, int *ldx, FILE *fp, char *postfix, gridinfo3d_t *grid3d)
Definition: screate_matrix3d.c:72
int_t sTrs2_GatherTrsmScatter(int_t klst, int_t iukp, int_t rukp, int_t *usub, float *uval, float *tempv, int_t knsupc, int nsupr, float *lusup, Glu_persist_t *Glu_persist)
Definition: psgstrf2.c:716
void sDestroy_A3d_gathered_on_2d(sSOLVEstruct_t *, gridinfo3d_t *)
Definition: psutil.c:798
int updateDirtyBit(int_t k0, HyP_t *HyP, gridinfo_t *grid)
Definition: sec_structs.c:618
int superlu_sgemv(const char *trans, const int m, const int n, const float alpha, const float *a, const int lda, const float *x, const int incx, const float beta, float *y, const int incy)
void sCreate_CompCol_Matrix_dist(SuperMatrix *, int_t, int_t, int_t, float *, int_t *, int_t *, Stype_t, Dtype_t, Mtype_t)
void sDestroy_trf3Dpartition(trf3Dpartition_t *trf3Dpartition, gridinfo3d_t *grid3d)
void psgstrs_Bglobal(int_t, sLUstruct_t *, gridinfo_t *, float *, int_t, int, SuperLUStat_t *, int *)
Definition: psgstrs_Bglobal.c:104
int saxpy_(const int *n, const float *alpha, const float *x, const int *incx, float *y, const int *incy)
int sPrint_CompRowLoc_Matrix_dist(SuperMatrix *)
int_t sWait_LRecv(MPI_Request *, int *msgcnt, int *msgcntsU, gridinfo_t *, SCT_t *)
int superlu_sscal(const int n, const float alpha, float *x, const int incx)
void sGatherNRformat_loc3d(fact_t Fact, NRformat_loc *A, float *B, int ldb, int nrhs, gridinfo3d_t *grid3d, NRformat_loc3d **)
int psCompRow_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: psutil.c:29
void pslaqgs(SuperMatrix *, float *, float *, float, float, float, char *)
Definition: pslaqgs.c:85
int psgsmv_AXglobal_abs(int_t, int_t[], float[], int_t[], float[], float[])
Definition: psgsmv_AXglobal.c:285
void sCreate_SuperNode_Matrix_dist(SuperMatrix *, int_t, int_t, int_t, float *, int_t *, int_t *, int_t *, int_t *, int_t *, Stype_t, Dtype_t, Mtype_t)
int_t szRecvUPanel(int_t k, int_t sender, float alpha, float beta, float *Uval_buf, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
void psgssvx3d(superlu_dist_options_t *, SuperMatrix *, sScalePermstruct_t *, float B[], int ldb, int nrhs, gridinfo3d_t *, sLUstruct_t *, sSOLVEstruct_t *, float *berr, SuperLUStat_t *, int *info)
Definition: psgssvx3d.c:501
void psgsrfs_ABXglobal(int_t, SuperMatrix *, float, sLUstruct_t *, gridinfo_t *, float *, int_t, float *, int_t, int, float *, SuperLUStat_t *, int *)
Definition: psgsrfs_ABXglobal.c:125
int_t sAllocLlu_3d(int_t nsupers, sLUstruct_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
void sCompRow_to_CompCol_dist(int_t, int_t, int_t, float *, int_t *, int_t *, float **, int_t **, int_t **)
void psgstrs2(int_t m, int_t k0, int_t k, Glu_persist_t *Glu_persist, gridinfo_t *grid, sLocalLU_t *Llu, SuperLUStat_t *stat)
int_t sscatter3dLPanels(int_t nsupers, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
void psinf_norm_error(int, int_t, int_t, float[], int_t, float[], int_t, MPI_Comm)
Check the inf-norm of the error vector.
Definition: psutil.c:825
int_t sscatter3dUPanels(int_t nsupers, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
int_t sp3dCollect(int_t layer, int_t n, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
void sreadrb_dist(int, FILE *, int_t *, int_t *, int_t *, float **, int_t **, int_t **)
Definition: sreadrb.c:275
void sPrint_Dense_Matrix_dist(SuperMatrix *)
int_t sDiagFactIBCast(int_t k, int_t k0, float *BlockUFactor, float *BlockLFactor, int_t *IrecvPlcd_D, MPI_Request *, MPI_Request *, MPI_Request *, MPI_Request *, gridinfo_t *, superlu_dist_options_t *, double thresh, sLUstruct_t *LUstruct, SuperLUStat_t *, int *info, SCT_t *, int tag_ub)
int sAllocGlu_3d(int_t n, int_t nsupers, sLUstruct_t *)
Definition: sutil_dist.c:432
void sGenXtrue_dist(int_t, int_t, float *, int_t)
Definition: sutil_dist.c:486
void sClone_CompRowLoc_Matrix_dist(SuperMatrix *, SuperMatrix *)
int getNsupers(int, Glu_persist_t *)
Definition: trfAux.c:42
void slsum_fmod_inv(float *, float *, float *, float *, int, int_t, int *fmod, int_t *, gridinfo_t *, sLocalLU_t *, SuperLUStat_t **, int_t *, int_t *, int_t, int_t, int_t, int_t, int, int)
Definition: psgstrs_lsum.c:416
float psdistribute(fact_t, int_t, SuperMatrix *, sScalePermstruct_t *, Glu_freeable_t *, sLUstruct_t *, gridinfo_t *)
Definition: psdistribute.c:326
int sDeAllocLlu_3d(int_t n, sLUstruct_t *, gridinfo3d_t *)
Definition: sutil_dist.c:450
void psgstrs(int_t, sLUstruct_t *, sScalePermstruct_t *, gridinfo_t *, float *, int_t, int_t, int_t, int, sSOLVEstruct_t *, SuperLUStat_t *, int *)
Definition: psgstrs.c:840
int sread_binary(FILE *, int_t *, int_t *, int_t *, float **, int_t **, int_t **)
Definition: sbinary_io.c:4
void sComputeLevelsets(int, int_t, gridinfo_t *, Glu_persist_t *, sLocalLU_t *, int_t *)
void sInit_HyP(HyP_t *HyP, sLocalLU_t *Llu, int_t mcb, int_t mrb)
void sPrint_CompCol_Matrix_dist(SuperMatrix *)
int_t szeroSetLU(int_t nnodes, int_t *nodeList, sLUstruct_t *, gridinfo3d_t *)
void sPrintLblocks(int, int_t, gridinfo_t *, Glu_persist_t *, sLocalLU_t *)
Print the blocks in the factored matrix L.
Definition: sutil_dist.c:570
int_t sblock_gemm_scatterTopLeft(int_t lb, int_t j, float *bigV, int_t knsupc, int_t klst, int_t *lsub, int_t *usub, int_t ldt, int *indirect, int *indirect2, HyP_t *HyP, sLUstruct_t *, gridinfo_t *, SCT_t *SCT, SuperLUStat_t *)
int_t sIBcastRecvUPanel(int_t k, int_t k0, int *msgcnt, MPI_Request *, MPI_Request *, int_t *Usub_buf, float *Uval_buf, gridinfo_t *, sLUstruct_t *, SCT_t *, int tag_ub)
int file_sPrint_CompRowLoc_Matrix_dist(FILE *fp, SuperMatrix *A)
void Free_HyP(HyP_t *HyP)
Definition: sec_structs.c:594
void psGetDiagU(int_t, sLUstruct_t *, gridinfo_t *, float *)
Definition: psGetDiagU.c:66
int_t ancestorFactor(int_t ilvl, sForest_t *sforest, commRequests_t **comReqss, sscuBufs_t *scuBufs, packLUInfo_t *packLUInfo, msgs_t **msgss, sLUValSubBuf_t **LUvsbs, sdiagFactBufs_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, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SuperLUStat_t *stat, double thresh, SCT_t *SCT, int tag_ub, int *info)
int_t sISend_LDiagBlock(int_t k0, float *lblk_ptr, int_t size, MPI_Request *, gridinfo_t *, int)
void sger_(const int *, const int *, const float *, const float *, const int *, const float *, const int *, float *, const int *)
int_t sreduceAncestors3d(int_t sender, int_t receiver, int_t nnodes, int_t *nodeList, float *Lval_buf, float *Uval_buf, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
int_t scollect3dUpanels(int_t layer, int_t nsupers, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
void sgemv_(const char *, const int *, const int *, const float *, const float *a, const int *, const float *, const int *, const float *, float *, const int *)
trf3Dpartition_t * sinitTrf3Dpartition(int_t nsupers, superlu_dist_options_t *options, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
void sinf_norm_error_dist(int_t, int_t, float *, int_t, float *, int_t, gridinfo_t *)
Check the inf-norm of the error vector.
Definition: sutil_dist.c:522
void sScalePermstructFree(sScalePermstruct_t *)
Deallocate ScalePermstruct.
Definition: sutil_dist.c:409
void sCreate_CompRowLoc_Matrix_dist(SuperMatrix *, int_t, int_t, int_t, int_t, int_t, float *, int_t *, int_t *, Stype_t, Dtype_t, Mtype_t)
int sstatic_schedule(superlu_dist_options_t *, int, int, sLUstruct_t *, gridinfo_t *, SuperLUStat_t *, int_t *, int_t *, int *)
Definition: sstatic_schedule.c:46
int_t checkRecvUDiag(int_t k, commRequests_t *comReqs, gridinfo_t *grid, SCT_t *SCT)
Definition: treeFactorization.c:401
float * floatCalloc_dist(int_t)
Definition: smemory_dist.c:162
int sSolveInit(superlu_dist_options_t *, SuperMatrix *, int_t[], int_t[], int_t, sLUstruct_t *, gridinfo_t *, sSOLVEstruct_t *)
Initialize the data structure for the solution phase.
Definition: psutil.c:690
int superlu_saxpy(const int n, const float alpha, const float *x, const int incx, float *y, const int incy)
int_t sblock_gemm_scatterBottomRight(int_t lb, int_t j, float *bigV, int_t knsupc, int_t klst, int_t *lsub, int_t *usub, int_t ldt, int *indirect, int *indirect2, HyP_t *HyP, sLUstruct_t *, gridinfo_t *, SCT_t *SCT, SuperLUStat_t *)
int superlu_sgemm(const char *transa, const char *transb, int m, int n, int k, float alpha, float *a, int lda, float *b, int ldb, float beta, float *c, int ldc)
void sCopy_CompRowLoc_Matrix_dist(SuperMatrix *, SuperMatrix *)
Definition: sutil_dist.c:324
void sgsequ_dist(SuperMatrix *, float *, float *, float *, float *, float *, int_t *)
Definition: sgsequ_dist.c:94
void sGenCSRLblocks(int, int_t, gridinfo_t *, Glu_persist_t *, sLocalLU_t *, float **, int_t **, int_t **, int_t *, int_t *)
int_t sblock_gemm_scatterTopRight(int_t lb, int_t j, float *bigV, int_t knsupc, int_t klst, int_t *lsub, int_t *usub, int_t ldt, int *indirect, int *indirect2, HyP_t *HyP, sLUstruct_t *, gridinfo_t *, SCT_t *SCT, SuperLUStat_t *)
int sDeAllocGlu_3d(sLUstruct_t *)
Definition: sutil_dist.c:442
void psgssvx_ABglobal(superlu_dist_options_t *, SuperMatrix *, sScalePermstruct_t *, float *, int, int, gridinfo_t *, sLUstruct_t *, float *, SuperLUStat_t *, int *)
void psgstrf2_trsm(superlu_dist_options_t *options, int_t k0, int_t k, double thresh, Glu_persist_t *, gridinfo_t *, sLocalLU_t *, MPI_Request *, int tag_ub, SuperLUStat_t *, int *info)
Definition: psgstrf2.c:143
int_t sIBcast_UPanel(int_t k, int_t k0, int_t *usub, float *uval, gridinfo_t *, int *msgcnt, MPI_Request *, int *ToSendD, int)
void sinit3DLUstructForest(int_t *myTreeIdxs, int_t *myZeroTrIdxs, sForest_t **sForests, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
int_t ssparseTreeFactor(int_t nnodes, int_t *perm_c_supno, treeTopoInfo_t *treeTopoInfo, commRequests_t *comReqs, sscuBufs_t *scuBufs, packLUInfo_t *packLUInfo, msgs_t *msgs, sLUValSubBuf_t *LUvsb, sdiagFactBufs_t *dFBuf, factStat_t *factStat, factNodelists_t *fNlists, superlu_dist_options_t *options, int_t *gIperm_c_supno, int_t ldt, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SuperLUStat_t *stat, double thresh, SCT_t *SCT, int *info)
int strsm_(const char *, const char *, const char *, const char *, const int *, const int *, const float *, const float *, const int *, float *, const int *)
int sreduceAllAncestors3d(int_t ilvl, int_t *myNodeCount, int_t **treePerm, sLUValSubBuf_t *LUvsb, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
void sreadtriple_dist(FILE *, int_t *, int_t *, int_t *, float **, int_t **, int_t **)
Definition: sreadtriple.c:35
int s_c2cpp_GetHWPM(SuperMatrix *, gridinfo_t *, sScalePermstruct_t *)
Definition: s_c2cpp_GetHWPM.cpp:54
void sScalePermstructInit(const int_t, const int_t, sScalePermstruct_t *)
Allocate storage in ScalePermstruct.
Definition: sutil_dist.c:398
void sscatter_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, float *tempv, int *indirect_thread, int *indirect2, int_t **Lrowind_bc_ptr, float **Lnzval_bc_ptr, gridinfo_t *grid)
int screate_matrix_rb(SuperMatrix *, int, float **, int *, float **, int *, FILE *, gridinfo_t *)
void psgstrf2(superlu_dist_options_t *, int_t nsupers, int_t k0, int_t k, double thresh, Glu_persist_t *, gridinfo_t *, sLocalLU_t *, MPI_Request *, int, SuperLUStat_t *, int *)
int sfreeScuBufs(sscuBufs_t *scuBufs)
int_t sdenseTreeFactor(int_t nnnodes, int_t *perm_c_supno, commRequests_t *comReqs, sscuBufs_t *scuBufs, packLUInfo_t *packLUInfo, msgs_t *msgs, sLUValSubBuf_t *LUvsb, sdiagFactBufs_t *dFBuf, factStat_t *factStat, factNodelists_t *fNlists, superlu_dist_options_t *options, int_t *gIperm_c_supno, int_t ldt, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SuperLUStat_t *stat, double thresh, SCT_t *SCT, int tag_ub, int *info)
int_t sIBcast_LPanel(int_t k, int_t k0, int_t *lsub, float *lusup, gridinfo_t *, int *msgcnt, MPI_Request *, int **ToSendR, int_t *xsup, int)
void sreadhb_dist(int, FILE *, int_t *, int_t *, int_t *, float **, int_t **, int_t **)
Definition: sreadhb.c:107
int freePackLUInfo(packLUInfo_t *packLUInfo)
Definition: treeFactorization.c:376
int_t psgstrs_init(int_t, int_t, int_t, int_t, int_t[], int_t[], gridinfo_t *grid, Glu_persist_t *, sSOLVEstruct_t *)
Definition: psutil.c:561
int_t sWait_URecv(MPI_Request *, int *msgcnt, SCT_t *)
int_t sIRecv_UDiagBlock(int_t k0, float *ublk_ptr, int_t size, int_t src, MPI_Request *, gridinfo_t *, SCT_t *, int)
int_t sinitScuBufs(int_t ldt, int_t num_threads, int_t nsupers, sscuBufs_t *, sLUstruct_t *, gridinfo_t *)
int_t sIrecv_LPanel(int_t k, int_t k0, int_t *Lsub_buf, float *Lval_buf, gridinfo_t *, MPI_Request *, sLocalLU_t *, int)
void sGenCSCLblocks(int, int_t, gridinfo_t *, Glu_persist_t *, sLocalLU_t *, float **, int_t **, int_t **, int_t *, int_t *)
int superlu_strsv(char *uplo, char *trans, char *diag, int n, float *a, int lda, float *x, int incx)
void sZeroLblocks(int, int, gridinfo_t *, sLUstruct_t *)
Sets all entries of matrix L to zero.
Definition: sutil_dist.c:619
float sdistribute(fact_t, int_t, SuperMatrix *, Glu_freeable_t *, sLUstruct_t *, gridinfo_t *)
Definition: sdistribute.c:63
int_t psgstrf(superlu_dist_options_t *, int, int, float anorm, sLUstruct_t *, gridinfo_t *, SuperLUStat_t *, int *)
Definition: psgstrf.c:242
void sGenCOOLblocks(int, int_t, gridinfo_t *, Glu_persist_t *, sLocalLU_t *, int_t **, int_t **, float **, int_t *, int_t *)
int_t sPackLBlock(int_t k, float *Dest, Glu_persist_t *, gridinfo_t *, sLocalLU_t *)
int_t sLPanelUpdate(int_t k, int_t *IrecvPlcd_D, int_t *factored_L, MPI_Request *, float *BlockUFactor, gridinfo_t *, sLUstruct_t *, SCT_t *)
int_t sUPanelUpdate(int_t k, int_t *factored_U, MPI_Request *, float *BlockLFactor, float *bigV, int_t ldt, Ublock_info_t *, gridinfo_t *, sLUstruct_t *, SuperLUStat_t *, SCT_t *)
void Local_Sgstrf2(superlu_dist_options_t *options, int_t k, double thresh, float *BlockUFactor, Glu_persist_t *, gridinfo_t *, sLocalLU_t *, SuperLUStat_t *, int *info, SCT_t *)
int_t treeFactor(int_t nnnodes, int_t *perm_c_supno, commRequests_t *comReqs, sscuBufs_t *scuBufs, packLUInfo_t *packLUInfo, msgs_t *msgs, sLUValSubBuf_t *LUvsb, sdiagFactBufs_t *dFBuf, factStat_t *factStat, factNodelists_t *fNlists, superlu_dist_options_t *options, int_t *gIperm_c_supno, int_t ldt, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SuperLUStat_t *stat, double thresh, SCT_t *SCT, int *info)
#define MAX_LOOKAHEADS
Definition: superlu_sdefs.h:96
int_t szSendLPanel(int_t k, int_t receiver, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
void * duser_malloc_dist(int_t, int_t)
Definition: dmemory_dist.c:30
int_t LDiagBlockRecvWait(int_t k, int_t *factored_U, MPI_Request *, gridinfo_t *)
Definition: communication_aux.c:218
int_t sIRecv_LDiagBlock(int_t k0, float *L_blk_ptr, int_t size, int_t src, MPI_Request *, gridinfo_t *, SCT_t *, int)
int superlu_strsm(const char *sideRL, const char *uplo, const char *transa, const char *diag, const int m, const int n, const float alpha, const float *a, const int lda, float *b, const int ldb)
void sDestroy_LU(int_t, gridinfo_t *, sLUstruct_t *)
Destroy distributed L & U matrices.
Definition: psutil.c:435
int psgsmv_AXglobal(int_t, int_t[], float[], int_t[], float[], float[])
Definition: psgsmv_AXglobal.c:259
int_t sUPanelTrSolve(int_t k, float *BlockLFactor, float *bigV, int_t ldt, Ublock_info_t *, gridinfo_t *, sLUstruct_t *, SuperLUStat_t *, SCT_t *)
int psgsmv_AXglobal_setup(SuperMatrix *, Glu_persist_t *, gridinfo_t *, int_t *, int_t *[], float *[], int_t *[], int_t[])
void slsum_fmod_inv_master(float *, float *, float *, float *, int, int, int_t, int *fmod, int_t, int_t *, gridinfo_t *, sLocalLU_t *, SuperLUStat_t **, int_t, int_t, int_t, int_t, int, int)
Definition: psgstrs_lsum.c:956
void sZero_CompRowLoc_Matrix_dist(SuperMatrix *)
Sets all entries of a matrix to zero, A_{i,j}=0, for i,j=1,..,n.
Definition: sutil_dist.c:339
void slsum_fmod(float *, float *, float *, float *, int, int, int_t, int *fmod, int_t, int_t, int_t, int_t *, gridinfo_t *, sLocalLU_t *, MPI_Request[], SuperLUStat_t *)
Definition: psgstrs_lsum.c:62
void sallocateA_dist(int_t, int_t, float **, int_t **, int_t **)
Definition: smemory_dist.c:147
int_t sIrecv_UPanel(int_t k, int_t k0, int_t *Usub_buf, float *, sLocalLU_t *, gridinfo_t *, MPI_Request *, int)
void psgsequ(SuperMatrix *, float *, float *, float *, float *, float *, int_t *, gridinfo_t *)
Definition: psgsequ.c:86
void sgather_u(int_t num_u_blks, Ublock_info_t *Ublock_info, int_t *usub, float *uval, float *bigU, int_t ldu, int_t *xsup, int_t klst)
void sLUstructInit(const int_t, sLUstruct_t *)
Allocate storage in LUstruct.
Definition: psutil.c:402
int_t sgatherFactoredLU(int_t sender, int_t receiver, int_t nnodes, int_t *nodeList, sLUValSubBuf_t *LUvsb, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
void sSolveFinalize(superlu_dist_options_t *, sSOLVEstruct_t *)
Release the resources used for the solution phase.
Definition: psutil.c:778
int sscal_(const int *n, const float *alpha, float *dx, const int *incx)
int_t sSchurComplementSetupGPU(int_t k, msgs_t *msgs, packLUInfo_t *, int_t *, int_t *, int_t *, gEtreeInfo_t *, factNodelists_t *, sscuBufs_t *, sLUValSubBuf_t *LUvsb, gridinfo_t *, sLUstruct_t *, HyP_t *)
void psCompute_Diag_Inv(int_t, sLUstruct_t *, gridinfo_t *, SuperLUStat_t *, int *)
Definition: psgstrs.c:649
void sreadtriple_noheader(FILE *, int_t *, int_t *, int_t *, float **, int_t **, int_t **)
Definition: sreadtriple_noheader.c:35
int_t sTrs2_GatherU(int_t iukp, int_t rukp, int_t klst, int_t nsupc, int_t ldu, int_t *usub, float *uval, float *tempv)
Definition: psgstrf2.c:669
void slaqgs_dist(SuperMatrix *, float *, float *, float, float, float, char *)
Definition: slaqgs_dist.c:93
int_t sblock_gemm_scatterBottomLeft(int_t lb, int_t j, float *bigV, int_t knsupc, int_t klst, int_t *lsub, int_t *usub, int_t ldt, int *indirect, int *indirect2, HyP_t *HyP, sLUstruct_t *, gridinfo_t *, SCT_t *SCT, SuperLUStat_t *)
int file_Printfloat5(FILE *, char *, int_t, float *)
Definition: sutil_dist.c:555
void sCreate_Dense_Matrix_dist(SuperMatrix *, int_t, int_t, float *, int_t, Stype_t, Dtype_t, Mtype_t)
void sblock_gemm_scatter(int_t lb, int_t j, Ublock_info_t *Ublock_info, Remain_info_t *Remain_info, float *L_mat, int ldl, float *U_mat, int ldu, float *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, float **Lnzval_bc_ptr, int_t **Ufstnz_br_ptr, float **Unzval_br_ptr, int_t *xsup, gridinfo_t *, SuperLUStat_t *)
int_t sUDiagBlockRecvWait(int_t k, int_t *IrecvPlcd_D, int_t *factored_L, MPI_Request *, gridinfo_t *, sLUstruct_t *, SCT_t *)
int strsv_(char *, char *, char *, int *, float *, int *, float *, int *)
void pxgstrs_finalize(pxgstrs_comm_t *)
Definition: util.c:266
float sdist_psymbtonum(fact_t, int_t, SuperMatrix *, sScalePermstruct_t *, Pslu_freeable_t *, sLUstruct_t *, gridinfo_t *)
Definition: pssymbfact_distdata.c:1188
int_t sLpanelUpdate(int_t off0, int_t nsupc, float *ublk_ptr, int_t ld_ujrow, float *lusup, int_t nsupr, SCT_t *)
int sldperm_dist(int, int, int_t, int_t[], int_t[], float[], int_t *, float[], float[])
Definition: sldperm_dist.c:96
int_t sinitDiagFactBufs(int_t ldt, sdiagFactBufs_t *dFBuf)
int sLluBufFreeArr(int_t numLA, sLUValSubBuf_t **LUvsbs)
sLUValSubBuf_t ** sLluBufInitArr(int_t numLA, sLUstruct_t *LUstruct)
float * sgetBigV(int_t, int_t)
int_t psReDistribute_B_to_X(float *B, int_t m_loc, int nrhs, int_t ldb, int_t fst_row, int_t *ilsum, float *x, sScalePermstruct_t *, Glu_persist_t *, gridinfo_t *, sSOLVEstruct_t *)
Definition: psgstrs.c:155
void sRgather_U(int_t k, int_t jj0, int_t *usub, float *uval, float *bigU, gEtreeInfo_t *, Glu_persist_t *, gridinfo_t *, HyP_t *, int_t *myIperm, int_t *iperm_c_supno, int_t *perm_u)
int_t psgstrf3d(superlu_dist_options_t *, int m, int n, float anorm, trf3Dpartition_t *, SCT_t *, sLUstruct_t *, gridinfo3d_t *, SuperLUStat_t *, int *)
Definition: psgstrf3d.c:121
void psgsmv_finalize(psgsmv_comm_t *)
Definition: psgsmv.c:371
int sp_strsv_dist(char *, char *, char *, SuperMatrix *, SuperMatrix *, float *, int *)
Definition: ssp_blas2_dist.c:95
void sLUstructFree(sLUstruct_t *)
Deallocate LUstruct.
Definition: psutil.c:416
int sgemm_(const char *, const char *, const int *, const int *, const int *, const float *, const float *, const int *, const float *, const int *, const float *, float *, const int *)
int_t sTrs2_ScatterU(int_t iukp, int_t rukp, int_t klst, int_t nsupc, int_t ldu, int_t *usub, float *uval, float *tempv)
Definition: psgstrf2.c:694
void sCopy_CompCol_Matrix_dist(SuperMatrix *, SuperMatrix *)
int sp_sgemv_dist(char *, float, SuperMatrix *, float *, int, float, float *, int)
SpGEMV.
Definition: ssp_blas2_dist.c:391
void slsum_bmod_inv(float *, float *, float *, float *, int, int_t, int *bmod, int_t *, Ucb_indptr_t **, int_t **, int_t *, gridinfo_t *, sLocalLU_t *, SuperLUStat_t **, int_t *, int_t *, int_t, int_t, int, int)
Definition: psgstrs_lsum.c:1355
float pslangs(char *, SuperMatrix *, gridinfo_t *)
Definition: pslangs.c:65
void sfill_dist(float *, int_t, float)
Fills a float precision array with a given value.
Definition: sutil_dist.c:512
void slsum_bmod(float *, float *, float *, int, int_t, int *bmod, int_t *, Ucb_indptr_t **, int_t **, int_t *, gridinfo_t *, sLocalLU_t *, MPI_Request[], SuperLUStat_t *)
Definition: psgstrs_lsum.c:246
int_t sRecv_UDiagBlock(int_t k0, float *ublk_ptr, int_t size, int_t src, gridinfo_t *, SCT_t *, int)
int_t ssparseTreeFactor_ASYNC(sForest_t *sforest, commRequests_t **comReqss, sscuBufs_t *scuBufs, packLUInfo_t *packLUInfo, msgs_t **msgss, sLUValSubBuf_t **LUvsbs, sdiagFactBufs_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, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SuperLUStat_t *stat, double thresh, SCT_t *SCT, int tag_ub, int *info)
int_t sSchurComplementSetup(int_t k, int *msgcnt, Ublock_info_t *, Remain_info_t *, uPanelInfo_t *, lPanelInfo_t *, int_t *, int_t *, int_t *, float *bigU, int_t *Lsub_buf, float *Lval_buf, int_t *Usub_buf, float *Uval_buf, gridinfo_t *, sLUstruct_t *)
void sRgather_L(int_t k, int_t *lsub, float *lusup, gEtreeInfo_t *, Glu_persist_t *, gridinfo_t *, HyP_t *, int_t *myIperm, int_t *iperm_c_supno)
int psPermute_Dense_Matrix(int_t, int_t, int_t[], int_t[], float[], int, float[], 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: psutil.c:291
void sscatter_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, float *tempv, int_t **Ufstnz_br_ptr, float **Unzval_br_ptr, gridinfo_t *grid)
void sScaleAdd_CompRowLoc_Matrix_dist(SuperMatrix *, SuperMatrix *, float)
Scale and add: adds a scalar multiple of one matrix to another. A_{i,j} = c * A_{i,...
Definition: sutil_dist.c:381
void sDestroy_Tree(int_t, gridinfo_t *, sLUstruct_t *)
Destroy broadcast and reduction trees used in triangular solve.
Definition: psutil.c:854
int_t sLluBufInit(sLUValSubBuf_t *, sLUstruct_t *)
void slsum_bmod_inv_master(float *, float *, float *, float *, int, int_t, int *bmod, int_t *, Ucb_indptr_t **, int_t **, int_t *, gridinfo_t *, sLocalLU_t *, SuperLUStat_t **, int_t, int_t, int, int)
Definition: psgstrs_lsum.c:1819
int_t sISend_UDiagBlock(int_t k0, float *ublk_ptr, int_t size, MPI_Request *, gridinfo_t *, int)
void duser_free_dist(int_t, int_t)
Definition: dmemory_dist.c:49
int_t sBcast_UPanel(int_t k, int_t k0, int_t *usub, float *uval, gridinfo_t *, int *msgcnt, int *ToSendD, SCT_t *, int)
int_t szSendUPanel(int_t k, int_t receiver, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
int_t checkRecvLDiag(int_t k, commRequests_t *comReqs, gridinfo_t *, SCT_t *)
Definition: treeFactorization.c:422
void sScaleAddId_CompRowLoc_Matrix_dist(SuperMatrix *, float)
Scale and add I: scales a matrix and adds an identity. A_{i,j} = c * A_{i,j} + \delta_{i,...
Definition: sutil_dist.c:356
int_t sLPanelTrSolve(int_t k, int_t *factored_L, float *BlockUFactor, gridinfo_t *, sLUstruct_t *)
int_t sWaitL(int_t k, int *msgcnt, int *msgcntU, MPI_Request *, MPI_Request *, gridinfo_t *, sLUstruct_t *, SCT_t *)
float slangs_dist(char *, SuperMatrix *)
Definition: slangs_dist.c:72
void sCopy_Dense_Matrix_dist(int_t, int_t, float *, int_t, float *, int_t)
int_t sWaitU(int_t k, int *msgcnt, MPI_Request *, MPI_Request *, gridinfo_t *, sLUstruct_t *, SCT_t *)
int_t sQuerySpace_dist(int_t, sLUstruct_t *, gridinfo_t *, SuperLUStat_t *, superlu_dist_mem_usage_t *)
Definition: smemory_dist.c:73
int_t sinit3DLUstruct(int_t *myTreeIdxs, int_t *myZeroTrIdxs, int_t *nodeCount, int_t **nodeList, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
int sp_sgemm_dist(char *, int, float, SuperMatrix *, float *, int, float, float *, int)
Definition: ssp_blas3_dist.c:126
int screate_matrix3d(SuperMatrix *A, int nrhs, float **rhs, int *ldb, float **x, int *ldx, FILE *fp, gridinfo3d_t *grid3d)
int_t sgatherAllFactoredLU(trf3Dpartition_t *trf3Dpartition, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
int_t scollect3dLpanels(int_t layer, int_t nsupers, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
sdiagFactBufs_t ** sinitDiagFactBufsArr(int_t mxLeafNode, int_t ldt, gridinfo_t *grid)
void s3D_printMemUse(trf3Dpartition_t *trf3Dpartition, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
Definition: smemory_dist.c:243
float * sgetBigU(int_t, gridinfo_t *, sLUstruct_t *)
float * floatMalloc_dist(int_t)
Definition: smemory_dist.c:155
void sgather_l(int_t num_LBlk, int_t knsupc, Remain_info_t *L_info, float *lval, int_t LD_lval, float *L_buff)
void sPrintUblocks(int, int_t, gridinfo_t *, Glu_persist_t *, sLocalLU_t *)
Print the blocks in the factored matrix U.
Definition: sutil_dist.c:750
void sFillRHS_dist(char *, int_t, float *, int_t, SuperMatrix *, float *, int_t)
Let rhs[i] = sum of i-th row of A, so the solution vector is all 1's.
Definition: sutil_dist.c:499
void sZeroUblocks(int iam, int n, gridinfo_t *, sLUstruct_t *)
Sets all entries of matrix U to zero.
Definition: sutil_dist.c:791
void psgstrs2_omp(int_t k0, int_t k, Glu_persist_t *, gridinfo_t *, sLocalLU_t *, Ublock_info_t *, SuperLUStat_t *)
Definition: psgstrf2.c:762
int sfreeDiagFactBufsArr(int_t mxLeafNode, sdiagFactBufs_t **dFBufs)
int_t sBcast_LPanel(int_t k, int_t k0, int_t *lsub, float *lusup, gridinfo_t *, int *msgcnt, int **ToSendR, int_t *xsup, SCT_t *, int)
int_t sp3dScatter(int_t n, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d)
int_t sIBcastRecvLPanel(int_t k, int_t k0, int *msgcnt, MPI_Request *, MPI_Request *, int_t *Lsub_buf, float *Lval_buf, int_t *factored, gridinfo_t *, sLUstruct_t *, SCT_t *, int tag_ub)
void psgsmv_init(SuperMatrix *, int_t *, gridinfo_t *, psgsmv_comm_t *)
Definition: psgsmv.c:27
int sScatter_B3d(NRformat_loc3d *A3d, gridinfo3d_t *grid3d)
void psgsmv(int_t, SuperMatrix *, gridinfo_t *, psgsmv_comm_t *, float x[], float ax[])
Definition: psgsmv.c:235
int_t sgatherAllFactoredLUFr(int_t *myZeroTrIdxs, sForest_t *sForests, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d, SCT_t *SCT)
int screate_matrix_dat(SuperMatrix *, int, float **, int *, float **, int *, FILE *, gridinfo_t *)
void psgsrfs(int_t, SuperMatrix *, float, sLUstruct_t *, sScalePermstruct_t *, gridinfo_t *, float[], int_t, float[], int_t, int, sSOLVEstruct_t *, float *, SuperLUStat_t *, int *)
int superlu_sger(const int m, const int n, const float alpha, const float *x, const int incx, const float *y, const int incy, float *a, const int lda)
void sreadMM_dist(FILE *, int_t *, int_t *, int_t *, float **, int_t **, int_t **)
Definition: sreadMM.c:38
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
float * lookAhead_L_buff
Definition: superlu_sdefs.h:336
float * bigU_host
Definition: superlu_sdefs.h:345
float * bigU_Phi
Definition: superlu_sdefs.h:344
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_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
float * lusup
Definition: superlu_sdefs.h:319
Definition: superlu_defs.h:947
Definition: superlu_ddefs.h:397
Definition: superlu_sdefs.h:263
int_t * ind_torecv
Definition: superlu_sdefs.h:266
float * val_tosend
Definition: superlu_sdefs.h:275
int_t * ptr_ind_torecv
Definition: superlu_sdefs.h:269
int * RecvCounts
Definition: superlu_sdefs.h:273
int_t TotalValSend
Definition: superlu_sdefs.h:279
int_t TotalIndSend
Definition: superlu_sdefs.h:277
int_t * ind_tosend
Definition: superlu_sdefs.h:265
int_t * ptr_ind_tosend
Definition: superlu_sdefs.h:267
int_t * extern_start
Definition: superlu_sdefs.h:264
float * val_torecv
Definition: superlu_sdefs.h:276
int * SendCounts
Definition: superlu_sdefs.h:271
Definition: superlu_defs.h:551
Definition: superlu_defs.h:901
Definition: superlu_sdefs.h:357
float * Uval_buf
Definition: superlu_sdefs.h:361
float * Lval_buf
Definition: superlu_sdefs.h:359
int_t * Usub_buf
Definition: superlu_sdefs.h:360
int_t * Lsub_buf
Definition: superlu_sdefs.h:358
Definition: superlu_sdefs.h:254
Glu_persist_t * Glu_persist
Definition: superlu_sdefs.h:256
int_t * etree
Definition: superlu_sdefs.h:255
sLocalLU_t * Llu
Definition: superlu_sdefs.h:257
char dt
Definition: superlu_sdefs.h:258
Definition: superlu_sdefs.h:97
float ** Lnzval_bc_ptr
Definition: superlu_sdefs.h:103
long int * Ufstnz_br_offset
Definition: superlu_sdefs.h:126
int_t ** Ucb_valptr
Definition: superlu_sdefs.h:208
int ** ToSendR
Definition: superlu_sdefs.h:161
int ** fsendx_plist
Definition: superlu_sdefs.h:166
int_t ** ut_sendx_plist
Definition: superlu_sdefs.h:195
int_t SolveMsgVol
Definition: superlu_sdefs.h:182
int nfrecvx
Definition: superlu_sdefs.h:168
int_t SolveMsgSent
Definition: superlu_sdefs.h:181
int_t L_SOLVE
Definition: superlu_sdefs.h:190
long int * Linv_bc_offset
Definition: superlu_sdefs.h:110
float ** Unzval_br_ptr
Definition: superlu_sdefs.h:129
int_t * Ufstnz_br_dat
Definition: superlu_sdefs.h:125
int_t ldalsum
Definition: superlu_sdefs.h:180
long int Unzval_br_cnt
Definition: superlu_sdefs.h:132
int nfsendx
Definition: superlu_sdefs.h:169
int_t * ilsum
Definition: superlu_sdefs.h:178
int_t ** Lrowind_bc_ptr
Definition: superlu_sdefs.h:98
long int * Lrowind_bc_offset
Definition: superlu_sdefs.h:100
float * Unzval_br_dat
Definition: superlu_sdefs.h:130
long int Ucb_valcnt
Definition: superlu_sdefs.h:211
long int Uinv_bc_cnt
Definition: superlu_sdefs.h:122
int_t ** Ufstnz_br_ptr
Definition: superlu_sdefs.h:124
int * brecv
Definition: superlu_sdefs.h:172
long int * Lindval_loc_bc_offset
Definition: superlu_sdefs.h:115
int_t inv
Definition: superlu_sdefs.h:217
C_Tree * LBtree_ptr
Definition: superlu_sdefs.h:135
int_t n_utsendx
Definition: superlu_sdefs.h:197
C_Tree * LRtree_ptr
Definition: superlu_sdefs.h:136
int_t * ut_ilsum
Definition: superlu_sdefs.h:193
Ucb_indptr_t * Ucb_inddat
Definition: superlu_sdefs.h:204
int_t n
Definition: superlu_sdefs.h:214
int_t * Lrowind_bc_dat
Definition: superlu_sdefs.h:99
int nbsendx
Definition: superlu_sdefs.h:174
long int * Uinv_bc_offset
Definition: superlu_sdefs.h:121
int * ToRecv
Definition: superlu_sdefs.h:159
int * frecv
Definition: superlu_sdefs.h:167
Ucb_indptr_t ** Ucb_indptr
Definition: superlu_sdefs.h:203
int_t ut_ldalsum
Definition: superlu_sdefs.h:192
int_t n_utrecvx
Definition: superlu_sdefs.h:198
int_t nroot
Definition: superlu_sdefs.h:200
int_t nleaf
Definition: superlu_sdefs.h:215
C_Tree * URtree_ptr
Definition: superlu_sdefs.h:138
int * bmod
Definition: superlu_sdefs.h:170
int_t n_utrecvmod
Definition: superlu_sdefs.h:199
int_t * Ucb_valdat
Definition: superlu_sdefs.h:209
int_t * utmod
Definition: superlu_sdefs.h:194
int * ToSendD
Definition: superlu_sdefs.h:160
int_t * Lindval_loc_bc_dat
Definition: superlu_sdefs.h:114
long int Lnzval_bc_cnt
Definition: superlu_sdefs.h:106
int_t ** Lrowind_bc_2_lsum
Definition: superlu_sdefs.h:118
int * mod_bit
Definition: superlu_sdefs.h:175
int_t * Urbs
Definition: superlu_sdefs.h:202
int_t FRECV
Definition: superlu_sdefs.h:191
long int * Ucb_indoffset
Definition: superlu_sdefs.h:205
int_t UT_SOLVE
Definition: superlu_sdefs.h:189
long int Ucb_indcnt
Definition: superlu_sdefs.h:206
float * ujrow
Definition: superlu_sdefs.h:149
int * fmod
Definition: superlu_sdefs.h:165
long int Ufstnz_br_cnt
Definition: superlu_sdefs.h:127
int nbrecvx
Definition: superlu_sdefs.h:173
int_t * ut_modbit
Definition: superlu_sdefs.h:201
float ** Uinv_bc_ptr
Definition: superlu_sdefs.h:119
int_t nfrecvmod
Definition: superlu_sdefs.h:216
float * Lnzval_bc_dat
Definition: superlu_sdefs.h:104
int ** bsendx_plist
Definition: superlu_sdefs.h:171
C_Tree * UBtree_ptr
Definition: superlu_sdefs.h:137
long int Linv_bc_cnt
Definition: superlu_sdefs.h:111
int_t ** Lindval_loc_bc_ptr
Definition: superlu_sdefs.h:113
long int * Lnzval_bc_offset
Definition: superlu_sdefs.h:105
float ** Linv_bc_ptr
Definition: superlu_sdefs.h:108
long int Lrowind_bc_cnt
Definition: superlu_sdefs.h:101
long int * Ucb_valoffset
Definition: superlu_sdefs.h:210
long int Lindval_loc_bc_cnt
Definition: superlu_sdefs.h:116
float * Linv_bc_dat
Definition: superlu_sdefs.h:109
int_t * utrecv
Definition: superlu_sdefs.h:196
float * Uinv_bc_dat
Definition: superlu_sdefs.h:120
long int * Unzval_br_offset
Definition: superlu_sdefs.h:131
int_t * Unnz
Definition: superlu_sdefs.h:117
Definition: superlu_sdefs.h:284
int_t * xrow_to_proc
Definition: superlu_sdefs.h:295
int_t * row_to_proc
Definition: superlu_sdefs.h:285
psgsmv_comm_t * gsmv_comm
Definition: superlu_sdefs.h:288
int_t * diag_len
Definition: superlu_sdefs.h:287
NRformat_loc3d * A3d
Definition: superlu_sdefs.h:296
int_t * inv_perm_c
Definition: superlu_sdefs.h:286
pxgstrs_comm_t * gstrs_comm
Definition: superlu_sdefs.h:290
int_t * A_colind_gsmv
Definition: superlu_sdefs.h:291
Definition: superlu_sdefs.h:76
int_t * perm_r
Definition: superlu_sdefs.h:80
float * R
Definition: superlu_sdefs.h:78
DiagScale_t DiagScale
Definition: superlu_sdefs.h:77
int_t * perm_c
Definition: superlu_sdefs.h:81
float * C
Definition: superlu_sdefs.h:79
Definition: superlu_sdefs.h:391
float * BlockUFactor
Definition: superlu_sdefs.h:393
float * BlockLFactor
Definition: superlu_sdefs.h:392
Definition: superlu_sdefs.h:385
float * bigV
Definition: superlu_sdefs.h:387
float * bigU
Definition: superlu_sdefs.h:386
Definition: superlu_defs.h:744
Definition: superlu_defs.h:712
Definition: superlu_defs.h:882
Definition: superlu_ddefs.h:372
sLUValSubBuf_t * LUvsb
Definition: superlu_sdefs.h:381
Definition: superlu_ddefs.h:308
float * uval
Definition: superlu_sdefs.h:313
Definitions which are precision-neutral.