xref: /petsc/src/mat/impls/aij/mpi/superlu_dist/superlu_dist.c (revision f1869c40995ebfee8175f23a5f765bd6679fe3a8)
1 
2 /*
3         Provides an interface to the SuperLU_DIST sparse solver
4 */
5 
6 #include <../src/mat/impls/aij/seq/aij.h>
7 #include <../src/mat/impls/aij/mpi/mpiaij.h>
8 #if defined(PETSC_HAVE_STDLIB_H) /* This is to get around weird problem with SuperLU on cray */
9 #include <stdlib.h>
10 #endif
11 
12 #if defined(PETSC_USE_64BIT_INDICES)
13 /* ugly SuperLU_Dist variable telling it to use long long int */
14 #define _LONGINT
15 #endif
16 
17 EXTERN_C_BEGIN
18 #if defined(PETSC_USE_COMPLEX)
19 #include <superlu_zdefs.h>
20 #else
21 #include <superlu_ddefs.h>
22 #endif
23 EXTERN_C_END
24 
25 /*
26     GLOBAL - The sparse matrix and right hand side are all stored initially on process 0. Should be called centralized
27     DISTRIBUTED - The sparse matrix and right hand size are initially stored across the entire MPI communicator.
28 */
29 typedef enum {GLOBAL,DISTRIBUTED} SuperLU_MatInputMode;
30 const char *SuperLU_MatInputModes[] = {"GLOBAL","DISTRIBUTED","SuperLU_MatInputMode","PETSC_",0};
31 
32 typedef struct {
33   int_t                  nprow,npcol,*row,*col;
34   gridinfo_t             grid;
35   superlu_dist_options_t options;
36   SuperMatrix            A_sup;
37   ScalePermstruct_t      ScalePermstruct;
38   LUstruct_t             LUstruct;
39   int                    StatPrint;
40   SuperLU_MatInputMode   MatInputMode;
41   SOLVEstruct_t          SOLVEstruct;
42   fact_t                 FactPattern;
43   MPI_Comm               comm_superlu;
44 #if defined(PETSC_USE_COMPLEX)
45   doublecomplex          *val;
46 #else
47   double                 *val;
48 #endif
49   PetscBool              matsolve_iscalled,matmatsolve_iscalled;
50   PetscBool              CleanUpSuperLU_Dist;  /* Flag to clean up (non-global) SuperLU objects during Destroy */
51 } Mat_SuperLU_DIST;
52 
53 
54 #undef __FUNCT__
55 #define __FUNCT__ "MatSuperluDistGetDiagU_SuperLU_DIST"
56 PetscErrorCode MatSuperluDistGetDiagU_SuperLU_DIST(Mat F,PetscScalar *diagU)
57 {
58   Mat_SuperLU_DIST  *lu= (Mat_SuperLU_DIST*)F->data;
59 
60   PetscFunctionBegin;
61 #if defined(PETSC_USE_COMPLEX)
62   PetscStackCall("SuperLU_DIST:pzGetDiagU",pzGetDiagU(F->rmap->N,&lu->LUstruct,&lu->grid,(doublecomplex*)diagU));
63 #else
64   PetscStackCall("SuperLU_DIST:pdGetDiagU",pdGetDiagU(F->rmap->N,&lu->LUstruct,&lu->grid,diagU));
65 #endif
66   PetscFunctionReturn(0);
67 }
68 
69 #undef __FUNCT__
70 #define __FUNCT__ "MatSuperluDistGetDiagU"
71 PetscErrorCode MatSuperluDistGetDiagU(Mat F,PetscScalar *diagU)
72 {
73   PetscErrorCode    ierr;
74 
75   PetscFunctionBegin;
76   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
77   ierr = PetscTryMethod(F,"MatSuperluDistGetDiagU_C",(Mat,PetscScalar*),(F,diagU));CHKERRQ(ierr);
78   PetscFunctionReturn(0);
79 }
80 
81 #undef __FUNCT__
82 #define __FUNCT__ "MatDestroy_SuperLU_DIST"
83 static PetscErrorCode MatDestroy_SuperLU_DIST(Mat A)
84 {
85   PetscErrorCode   ierr;
86   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data;
87 
88   PetscFunctionBegin;
89   if (lu->CleanUpSuperLU_Dist) {
90     /* Deallocate SuperLU_DIST storage */
91     if (lu->MatInputMode == GLOBAL) {
92       PetscStackCall("SuperLU_DIST:Destroy_CompCol_Matrix_dist",Destroy_CompCol_Matrix_dist(&lu->A_sup));
93     } else {
94       PetscStackCall("SuperLU_DIST:Destroy_CompRowLoc_Matrix_dist",Destroy_CompRowLoc_Matrix_dist(&lu->A_sup));
95       if (lu->options.SolveInitialized) {
96 #if defined(PETSC_USE_COMPLEX)
97         PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
98 #else
99         PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
100 #endif
101       }
102     }
103     PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(A->cmap->N, &lu->grid, &lu->LUstruct));
104     PetscStackCall("SuperLU_DIST:ScalePermstructFree",ScalePermstructFree(&lu->ScalePermstruct));
105     PetscStackCall("SuperLU_DIST:LUstructFree",LUstructFree(&lu->LUstruct));
106 
107     /* Release the SuperLU_DIST process grid. */
108     PetscStackCall("SuperLU_DIST:superlu_gridexit",superlu_gridexit(&lu->grid));
109     ierr = MPI_Comm_free(&(lu->comm_superlu));CHKERRQ(ierr);
110   }
111   ierr = PetscFree(A->data);CHKERRQ(ierr);
112   /* clear composed functions */
113   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr);
114   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSuperluDistGetDiagU_C",NULL);CHKERRQ(ierr);
115 
116   PetscFunctionReturn(0);
117 }
118 
119 #undef __FUNCT__
120 #define __FUNCT__ "MatSolve_SuperLU_DIST"
121 static PetscErrorCode MatSolve_SuperLU_DIST(Mat A,Vec b_mpi,Vec x)
122 {
123   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data;
124   PetscErrorCode   ierr;
125   PetscMPIInt      size;
126   PetscInt         m=A->rmap->n,M=A->rmap->N,N=A->cmap->N;
127   SuperLUStat_t    stat;
128   double           berr[1];
129   PetscScalar      *bptr;
130   PetscInt         nrhs=1;
131   Vec              x_seq;
132   IS               iden;
133   VecScatter       scat;
134   int              info; /* SuperLU_Dist info code is ALWAYS an int, even with long long indices */
135   static PetscBool cite = PETSC_FALSE;
136 
137   PetscFunctionBegin;
138   ierr = PetscCitationsRegister("@article{lidemmel03,\n  author = {Xiaoye S. Li and James W. Demmel},\n  title = {{SuperLU_DIST}: A Scalable Distributed-Memory Sparse Direct\n           Solver for Unsymmetric Linear Systems},\n  journal = {ACM Trans. Mathematical Software},\n  volume = {29},\n  number = {2},\n  pages = {110-140},\n  year = 2003\n}\n",&cite);CHKERRQ(ierr);
139   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
140   if (size > 1 && lu->MatInputMode == GLOBAL) {
141     /* global mat input, convert b to x_seq */
142     ierr = VecCreateSeq(PETSC_COMM_SELF,N,&x_seq);CHKERRQ(ierr);
143     ierr = ISCreateStride(PETSC_COMM_SELF,N,0,1,&iden);CHKERRQ(ierr);
144     ierr = VecScatterCreate(b_mpi,iden,x_seq,iden,&scat);CHKERRQ(ierr);
145     ierr = ISDestroy(&iden);CHKERRQ(ierr);
146 
147     ierr = VecScatterBegin(scat,b_mpi,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
148     ierr = VecScatterEnd(scat,b_mpi,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
149     ierr = VecGetArray(x_seq,&bptr);CHKERRQ(ierr);
150   } else { /* size==1 || distributed mat input */
151     if (lu->options.SolveInitialized && !lu->matsolve_iscalled) {
152       /* see comments in MatMatSolve() */
153 #if defined(PETSC_USE_COMPLEX)
154       PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
155 #else
156       PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
157 #endif
158       lu->options.SolveInitialized = NO;
159     }
160     ierr = VecCopy(b_mpi,x);CHKERRQ(ierr);
161     ierr = VecGetArray(x,&bptr);CHKERRQ(ierr);
162   }
163 
164   if (lu->options.Fact != FACTORED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"SuperLU_DIST options.Fact mush equal FACTORED");
165 
166   PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat));        /* Initialize the statistics variables. */
167   if (lu->MatInputMode == GLOBAL) {
168 #if defined(PETSC_USE_COMPLEX)
169     PetscStackCall("SuperLU_DIST:pzgssvx_ABglobal",pzgssvx_ABglobal(&lu->options,&lu->A_sup,&lu->ScalePermstruct,(doublecomplex*)bptr,M,nrhs,&lu->grid,&lu->LUstruct,berr,&stat,&info));
170 #else
171     PetscStackCall("SuperLU_DIST:pdgssvx_ABglobal",pdgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct,bptr,M,nrhs,&lu->grid,&lu->LUstruct,berr,&stat,&info));
172 #endif
173   } else { /* distributed mat input */
174 #if defined(PETSC_USE_COMPLEX)
175     PetscStackCall("SuperLU_DIST:pzgssvx",pzgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,(doublecomplex*)bptr,m,nrhs,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
176 #else
177     PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,bptr,m,nrhs,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
178 #endif
179   }
180   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"pdgssvx fails, info: %d\n",info);
181 
182   if (lu->options.PrintStat) PStatPrint(&lu->options, &stat, &lu->grid);      /* Print the statistics. */
183   PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat));
184 
185   if (size > 1 && lu->MatInputMode == GLOBAL) {
186     /* convert seq x to mpi x */
187     ierr = VecRestoreArray(x_seq,&bptr);CHKERRQ(ierr);
188     ierr = VecScatterBegin(scat,x_seq,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
189     ierr = VecScatterEnd(scat,x_seq,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
190     ierr = VecScatterDestroy(&scat);CHKERRQ(ierr);
191     ierr = VecDestroy(&x_seq);CHKERRQ(ierr);
192   } else {
193     ierr = VecRestoreArray(x,&bptr);CHKERRQ(ierr);
194 
195     lu->matsolve_iscalled    = PETSC_TRUE;
196     lu->matmatsolve_iscalled = PETSC_FALSE;
197   }
198   PetscFunctionReturn(0);
199 }
200 
201 #undef __FUNCT__
202 #define __FUNCT__ "MatMatSolve_SuperLU_DIST"
203 static PetscErrorCode MatMatSolve_SuperLU_DIST(Mat A,Mat B_mpi,Mat X)
204 {
205   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data;
206   PetscErrorCode   ierr;
207   PetscMPIInt      size;
208   PetscInt         M=A->rmap->N,m=A->rmap->n,nrhs;
209   SuperLUStat_t    stat;
210   double           berr[1];
211   PetscScalar      *bptr;
212   int              info; /* SuperLU_Dist info code is ALWAYS an int, even with long long indices */
213   PetscBool        flg;
214 
215   PetscFunctionBegin;
216   if (lu->options.Fact != FACTORED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"SuperLU_DIST options.Fact mush equal FACTORED");
217   ierr = PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
218   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
219   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
220   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
221 
222   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
223   if (size > 1 && lu->MatInputMode == GLOBAL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatInputMode=GLOBAL for nproc %d>1 is not supported",size);
224   /* size==1 or distributed mat input */
225   if (lu->options.SolveInitialized && !lu->matmatsolve_iscalled) {
226     /* communication pattern of SOLVEstruct is unlikely created for matmatsolve,
227        thus destroy it and create a new SOLVEstruct.
228        Otherwise it may result in memory corruption or incorrect solution
229        See src/mat/examples/tests/ex125.c */
230 #if defined(PETSC_USE_COMPLEX)
231     PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
232 #else
233     PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
234 #endif
235     lu->options.SolveInitialized = NO;
236   }
237   ierr = MatCopy(B_mpi,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
238 
239   ierr = MatGetSize(B_mpi,NULL,&nrhs);CHKERRQ(ierr);
240 
241   PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat));        /* Initialize the statistics variables. */
242   ierr = MatDenseGetArray(X,&bptr);CHKERRQ(ierr);
243   if (lu->MatInputMode == GLOBAL) { /* size == 1 */
244 #if defined(PETSC_USE_COMPLEX)
245     PetscStackCall("SuperLU_DIST:pzgssvx_ABglobal",pzgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct,(doublecomplex*)bptr, M, nrhs,&lu->grid, &lu->LUstruct, berr, &stat, &info));
246 #else
247     PetscStackCall("SuperLU_DIST:pdgssvx_ABglobal",pdgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct,bptr, M, nrhs, &lu->grid, &lu->LUstruct, berr, &stat, &info));
248 #endif
249   } else { /* distributed mat input */
250 #if defined(PETSC_USE_COMPLEX)
251     PetscStackCall("SuperLU_DIST:pzgssvx",pzgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,(doublecomplex*)bptr,m,nrhs,&lu->grid, &lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
252 #else
253     PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,bptr,m,nrhs,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
254 #endif
255   }
256   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"pdgssvx fails, info: %d\n",info);
257   ierr = MatDenseRestoreArray(X,&bptr);CHKERRQ(ierr);
258 
259   if (lu->options.PrintStat) PStatPrint(&lu->options, &stat, &lu->grid); /* Print the statistics. */
260   PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat));
261   lu->matsolve_iscalled    = PETSC_FALSE;
262   lu->matmatsolve_iscalled = PETSC_TRUE;
263   PetscFunctionReturn(0);
264 }
265 
266 
267 #undef __FUNCT__
268 #define __FUNCT__ "MatLUFactorNumeric_SuperLU_DIST"
269 static PetscErrorCode MatLUFactorNumeric_SuperLU_DIST(Mat F,Mat A,const MatFactorInfo *info)
270 {
271   Mat              *tseq,A_seq = NULL;
272   Mat_SeqAIJ       *aa,*bb;
273   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)(F)->data;
274   PetscErrorCode   ierr;
275   PetscInt         M=A->rmap->N,N=A->cmap->N,i,*ai,*aj,*bi,*bj,nz,rstart,*garray,
276                    m=A->rmap->n, colA_start,j,jcol,jB,countA,countB,*bjj,*ajj;
277   int              sinfo;   /* SuperLU_Dist info flag is always an int even with long long indices */
278   PetscMPIInt      size;
279   SuperLUStat_t    stat;
280   double           *berr=0;
281   IS               isrow;
282 #if defined(PETSC_USE_COMPLEX)
283   doublecomplex    *av, *bv;
284 #else
285   double           *av, *bv;
286 #endif
287 
288   PetscFunctionBegin;
289   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
290 
291   if (lu->MatInputMode == GLOBAL) { /* global mat input */
292     if (size > 1) { /* convert mpi A to seq mat A */
293       ierr = ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow);CHKERRQ(ierr);
294       ierr = MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&tseq);CHKERRQ(ierr);
295       ierr = ISDestroy(&isrow);CHKERRQ(ierr);
296 
297       A_seq = *tseq;
298       ierr  = PetscFree(tseq);CHKERRQ(ierr);
299       aa    = (Mat_SeqAIJ*)A_seq->data;
300     } else {
301       PetscBool flg;
302       ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
303       if (flg) {
304         Mat_MPIAIJ *At = (Mat_MPIAIJ*)A->data;
305         A = At->A;
306       }
307       aa =  (Mat_SeqAIJ*)A->data;
308     }
309 
310     /* Convert Petsc NR matrix to SuperLU_DIST NC.
311        Note: memories of lu->val, col and row are allocated by CompRow_to_CompCol_dist()! */
312     if (lu->options.Fact != DOFACT) {/* successive numeric factorization, sparsity pattern is reused. */
313       PetscStackCall("SuperLU_DIST:Destroy_CompCol_Matrix_dist",Destroy_CompCol_Matrix_dist(&lu->A_sup));
314       if (lu->FactPattern == SamePattern_SameRowPerm) {
315         lu->options.Fact = SamePattern_SameRowPerm; /* matrix has similar numerical values */
316       } else { /* lu->FactPattern == SamePattern */
317         PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(N, &lu->grid, &lu->LUstruct));
318         lu->options.Fact = SamePattern;
319       }
320     }
321 #if defined(PETSC_USE_COMPLEX)
322     PetscStackCall("SuperLU_DIST:zCompRow_to_CompCol_dist",zCompRow_to_CompCol_dist(M,N,aa->nz,(doublecomplex*)aa->a,(int_t*)aa->j,(int_t*)aa->i,&lu->val,&lu->col, &lu->row));
323 #else
324     PetscStackCall("SuperLU_DIST:dCompRow_to_CompCol_dist",dCompRow_to_CompCol_dist(M,N,aa->nz,aa->a,(int_t*)aa->j,(int_t*)aa->i,&lu->val, &lu->col, &lu->row));
325 #endif
326 
327     /* Create compressed column matrix A_sup. */
328 #if defined(PETSC_USE_COMPLEX)
329     PetscStackCall("SuperLU_DIST:zCreate_CompCol_Matrix_dist",zCreate_CompCol_Matrix_dist(&lu->A_sup, M, N, aa->nz, lu->val, lu->col, lu->row, SLU_NC, SLU_Z, SLU_GE));
330 #else
331     PetscStackCall("SuperLU_DIST:dCreate_CompCol_Matrix_dist",dCreate_CompCol_Matrix_dist(&lu->A_sup, M, N, aa->nz, lu->val, lu->col, lu->row, SLU_NC, SLU_D, SLU_GE));
332 #endif
333   } else { /* distributed mat input */
334     Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data;
335     aa=(Mat_SeqAIJ*)(mat->A)->data;
336     bb=(Mat_SeqAIJ*)(mat->B)->data;
337     ai=aa->i; aj=aa->j;
338     bi=bb->i; bj=bb->j;
339 #if defined(PETSC_USE_COMPLEX)
340     av=(doublecomplex*)aa->a;
341     bv=(doublecomplex*)bb->a;
342 #else
343     av=aa->a;
344     bv=bb->a;
345 #endif
346     rstart = A->rmap->rstart;
347     nz     = aa->nz + bb->nz;
348     garray = mat->garray;
349 
350     if (lu->options.Fact == DOFACT) { /* first numeric factorization */
351 #if defined(PETSC_USE_COMPLEX)
352       PetscStackCall("SuperLU_DIST:zallocateA_dist",zallocateA_dist(m, nz, &lu->val, &lu->col, &lu->row));
353 #else
354       PetscStackCall("SuperLU_DIST:dallocateA_dist",dallocateA_dist(m, nz, &lu->val, &lu->col, &lu->row));
355 #endif
356     } else { /* successive numeric factorization, sparsity pattern and perm_c are reused. */
357       if (lu->FactPattern == SamePattern_SameRowPerm) {
358         lu->options.Fact = SamePattern_SameRowPerm; /* matrix has similar numerical values */
359       } else if (lu->FactPattern == SamePattern) {
360         PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(N, &lu->grid, &lu->LUstruct)); /* Deallocate L and U matrices. */
361         lu->options.Fact = SamePattern;
362       } else if (lu->FactPattern == DOFACT) {
363         SUPERLU_FREE( lu->A_sup.Store );
364         PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(N, &lu->grid, &lu->LUstruct));
365         lu->options.Fact = DOFACT;
366       } else {
367         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"options.Fact must be one of SamePattern SamePattern_SameRowPerm");
368       }
369     }
370     nz = 0;
371     for (i=0; i<m; i++) {
372       lu->row[i] = nz;
373       countA     = ai[i+1] - ai[i];
374       countB     = bi[i+1] - bi[i];
375       if (aj) {
376         ajj = aj + ai[i]; /* ptr to the beginning of this row */
377       } else {
378         ajj = NULL;
379       }
380       bjj = bj + bi[i];
381 
382       /* B part, smaller col index */
383       if (aj) {
384         colA_start = rstart + ajj[0]; /* the smallest global col index of A */
385       } else { /* superlu_dist does not require matrix has diagonal entries, thus aj=NULL would work */
386         colA_start = rstart;
387       }
388       jB         = 0;
389       for (j=0; j<countB; j++) {
390         jcol = garray[bjj[j]];
391         if (jcol > colA_start) {
392           jB = j;
393           break;
394         }
395         lu->col[nz]   = jcol;
396         lu->val[nz++] = *bv++;
397         if (j==countB-1) jB = countB;
398       }
399 
400       /* A part */
401       for (j=0; j<countA; j++) {
402         lu->col[nz]   = rstart + ajj[j];
403         lu->val[nz++] = *av++;
404       }
405 
406       /* B part, larger col index */
407       for (j=jB; j<countB; j++) {
408         lu->col[nz]   = garray[bjj[j]];
409         lu->val[nz++] = *bv++;
410       }
411     }
412     lu->row[m] = nz;
413 
414     if (lu->options.Fact == DOFACT) {
415 #if defined(PETSC_USE_COMPLEX)
416       PetscStackCall("SuperLU_DIST:zCreate_CompRowLoc_Matrix_dist",zCreate_CompRowLoc_Matrix_dist(&lu->A_sup, M, N, nz, m, rstart,lu->val, lu->col, lu->row, SLU_NR_loc, SLU_Z, SLU_GE));
417 #else
418       PetscStackCall("SuperLU_DIST:dCreate_CompRowLoc_Matrix_dist",dCreate_CompRowLoc_Matrix_dist(&lu->A_sup, M, N, nz, m, rstart,lu->val, lu->col, lu->row, SLU_NR_loc, SLU_D, SLU_GE));
419 #endif
420     }
421   }
422 
423   /* Factor the matrix. */
424   PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat));   /* Initialize the statistics variables. */
425   if (lu->MatInputMode == GLOBAL) { /* global mat input */
426 #if defined(PETSC_USE_COMPLEX)
427     PetscStackCall("SuperLU_DIST:pzgssvx_ABglobal",pzgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, M, 0,&lu->grid, &lu->LUstruct, berr, &stat, &sinfo));
428 #else
429     PetscStackCall("SuperLU_DIST:pdgssvx_ABglobal",pdgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, M, 0,&lu->grid, &lu->LUstruct, berr, &stat, &sinfo));
430 #endif
431   } else { /* distributed mat input */
432 #if defined(PETSC_USE_COMPLEX)
433     PetscStackCall("SuperLU_DIST:pzgssvx",pzgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, m, 0, &lu->grid,&lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &sinfo));
434 #else
435     PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, m, 0, &lu->grid,&lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &sinfo));
436 #endif
437   }
438 
439   if (sinfo > 0) {
440     if (A->erroriffailure) {
441       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo);
442     } else {
443       if (sinfo <= lu->A_sup.ncol) {
444         F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
445         ierr = PetscInfo1(F,"U(i,i) is exactly zero, i= %D\n",sinfo);CHKERRQ(ierr);
446       } else if (sinfo > lu->A_sup.ncol) {
447         /*
448          number of bytes allocated when memory allocation
449          failure occurred, plus A->ncol.
450          */
451         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
452         ierr = PetscInfo1(F,"Number of bytes allocated when memory allocation fails %D\n",sinfo);CHKERRQ(ierr);
453       }
454     }
455   } else if (sinfo < 0) {
456     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, argument in p*gssvx() had an illegal value", sinfo);
457   }
458 
459   if (lu->MatInputMode == GLOBAL && size > 1) {
460     ierr = MatDestroy(&A_seq);CHKERRQ(ierr);
461   }
462 
463   if (lu->options.PrintStat) {
464     PStatPrint(&lu->options, &stat, &lu->grid);  /* Print the statistics. */
465   }
466   PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat));
467   (F)->assembled    = PETSC_TRUE;
468   (F)->preallocated = PETSC_TRUE;
469   lu->options.Fact  = FACTORED; /* The factored form of A is supplied. Local option used by this func. only */
470   PetscFunctionReturn(0);
471 }
472 
473 /* Note the Petsc r and c permutations are ignored */
474 #undef __FUNCT__
475 #define __FUNCT__ "MatLUFactorSymbolic_SuperLU_DIST"
476 static PetscErrorCode MatLUFactorSymbolic_SuperLU_DIST(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
477 {
478   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)F->data;
479   PetscInt         M   = A->rmap->N,N=A->cmap->N;
480 
481   PetscFunctionBegin;
482   /* Initialize the SuperLU process grid. */
483   PetscStackCall("SuperLU_DIST:superlu_gridinit",superlu_gridinit(lu->comm_superlu, lu->nprow, lu->npcol, &lu->grid));
484 
485   /* Initialize ScalePermstruct and LUstruct. */
486   PetscStackCall("SuperLU_DIST:ScalePermstructInit",ScalePermstructInit(M, N, &lu->ScalePermstruct));
487   PetscStackCall("SuperLU_DIST:LUstructInit",LUstructInit(N, &lu->LUstruct));
488   F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU_DIST;
489   F->ops->solve           = MatSolve_SuperLU_DIST;
490   F->ops->matsolve        = MatMatSolve_SuperLU_DIST;
491   lu->CleanUpSuperLU_Dist = PETSC_TRUE;
492   PetscFunctionReturn(0);
493 }
494 
495 #undef __FUNCT__
496 #define __FUNCT__ "MatFactorGetSolverPackage_aij_superlu_dist"
497 static PetscErrorCode MatFactorGetSolverPackage_aij_superlu_dist(Mat A,const MatSolverPackage *type)
498 {
499   PetscFunctionBegin;
500   *type = MATSOLVERSUPERLU_DIST;
501   PetscFunctionReturn(0);
502 }
503 
504 #undef __FUNCT__
505 #define __FUNCT__ "MatFactorInfo_SuperLU_DIST"
506 static PetscErrorCode MatFactorInfo_SuperLU_DIST(Mat A,PetscViewer viewer)
507 {
508   Mat_SuperLU_DIST       *lu=(Mat_SuperLU_DIST*)A->data;
509   superlu_dist_options_t options;
510   PetscErrorCode         ierr;
511 
512   PetscFunctionBegin;
513   /* check if matrix is superlu_dist type */
514   if (A->ops->solve != MatSolve_SuperLU_DIST) PetscFunctionReturn(0);
515 
516   options = lu->options;
517   ierr    = PetscViewerASCIIPrintf(viewer,"SuperLU_DIST run parameters:\n");CHKERRQ(ierr);
518   ierr    = PetscViewerASCIIPrintf(viewer,"  Process grid nprow %D x npcol %D \n",lu->nprow,lu->npcol);CHKERRQ(ierr);
519   ierr    = PetscViewerASCIIPrintf(viewer,"  Equilibrate matrix %s \n",PetscBools[options.Equil != NO]);CHKERRQ(ierr);
520   ierr    = PetscViewerASCIIPrintf(viewer,"  Matrix input mode %d \n",lu->MatInputMode);CHKERRQ(ierr);
521   ierr    = PetscViewerASCIIPrintf(viewer,"  Replace tiny pivots %s \n",PetscBools[options.ReplaceTinyPivot != NO]);CHKERRQ(ierr);
522   ierr    = PetscViewerASCIIPrintf(viewer,"  Use iterative refinement %s \n",PetscBools[options.IterRefine == SLU_DOUBLE]);CHKERRQ(ierr);
523   ierr    = PetscViewerASCIIPrintf(viewer,"  Processors in row %d col partition %d \n",lu->nprow,lu->npcol);CHKERRQ(ierr);
524   ierr    = PetscViewerASCIIPrintf(viewer,"  Row permutation %s \n",(options.RowPerm == NOROWPERM) ? "NATURAL" : "LargeDiag");CHKERRQ(ierr);
525 
526   switch (options.ColPerm) {
527   case NATURAL:
528     ierr = PetscViewerASCIIPrintf(viewer,"  Column permutation NATURAL\n");CHKERRQ(ierr);
529     break;
530   case MMD_AT_PLUS_A:
531     ierr = PetscViewerASCIIPrintf(viewer,"  Column permutation MMD_AT_PLUS_A\n");CHKERRQ(ierr);
532     break;
533   case MMD_ATA:
534     ierr = PetscViewerASCIIPrintf(viewer,"  Column permutation MMD_ATA\n");CHKERRQ(ierr);
535     break;
536   case METIS_AT_PLUS_A:
537     ierr = PetscViewerASCIIPrintf(viewer,"  Column permutation METIS_AT_PLUS_A\n");CHKERRQ(ierr);
538     break;
539   case PARMETIS:
540     ierr = PetscViewerASCIIPrintf(viewer,"  Column permutation PARMETIS\n");CHKERRQ(ierr);
541     break;
542   default:
543     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown column permutation");
544   }
545 
546   ierr = PetscViewerASCIIPrintf(viewer,"  Parallel symbolic factorization %s \n",PetscBools[options.ParSymbFact != NO]);CHKERRQ(ierr);
547 
548   if (lu->FactPattern == SamePattern) {
549     ierr = PetscViewerASCIIPrintf(viewer,"  Repeated factorization SamePattern\n");CHKERRQ(ierr);
550   } else {
551     ierr = PetscViewerASCIIPrintf(viewer,"  Repeated factorization SamePattern_SameRowPerm\n");CHKERRQ(ierr);
552   }
553   PetscFunctionReturn(0);
554 }
555 
556 #undef __FUNCT__
557 #define __FUNCT__ "MatView_SuperLU_DIST"
558 static PetscErrorCode MatView_SuperLU_DIST(Mat A,PetscViewer viewer)
559 {
560   PetscErrorCode    ierr;
561   PetscBool         iascii;
562   PetscViewerFormat format;
563 
564   PetscFunctionBegin;
565   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
566   if (iascii) {
567     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
568     if (format == PETSC_VIEWER_ASCII_INFO) {
569       ierr = MatFactorInfo_SuperLU_DIST(A,viewer);CHKERRQ(ierr);
570     }
571   }
572   PetscFunctionReturn(0);
573 }
574 
575 #undef __FUNCT__
576 #define __FUNCT__ "MatGetFactor_aij_superlu_dist"
577 static PetscErrorCode MatGetFactor_aij_superlu_dist(Mat A,MatFactorType ftype,Mat *F)
578 {
579   Mat                    B;
580   Mat_SuperLU_DIST       *lu;
581   PetscErrorCode         ierr;
582   PetscInt               M=A->rmap->N,N=A->cmap->N,indx;
583   PetscMPIInt            size;
584   superlu_dist_options_t options;
585   PetscBool              flg;
586   const char             *colperm[]     = {"NATURAL","MMD_AT_PLUS_A","MMD_ATA","METIS_AT_PLUS_A","PARMETIS"};
587   const char             *rowperm[]     = {"LargeDiag","NATURAL"};
588   const char             *factPattern[] = {"SamePattern","SamePattern_SameRowPerm","DOFACT"};
589   PetscBool              set;
590 
591   PetscFunctionBegin;
592   /* Create the factorization matrix */
593   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
594   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);CHKERRQ(ierr);
595   ierr = PetscStrallocpy("superlu_dist",&((PetscObject)B)->type_name);CHKERRQ(ierr);
596   ierr = MatSetUp(B);CHKERRQ(ierr);
597   B->ops->getinfo          = MatGetInfo_External;
598   B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU_DIST;
599   B->ops->view             = MatView_SuperLU_DIST;
600   B->ops->destroy          = MatDestroy_SuperLU_DIST;
601 
602   B->factortype = MAT_FACTOR_LU;
603 
604   /* set solvertype */
605   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
606   ierr = PetscStrallocpy(MATSOLVERSUPERLU_DIST,&B->solvertype);CHKERRQ(ierr);
607 
608   ierr    = PetscNewLog(B,&lu);CHKERRQ(ierr);
609   B->data = lu;
610 
611   /* Set the default input options:
612      options.Fact              = DOFACT;
613      options.Equil             = YES;
614      options.ParSymbFact       = NO;
615      options.ColPerm           = METIS_AT_PLUS_A;
616      options.RowPerm           = LargeDiag;
617      options.ReplaceTinyPivot  = YES;
618      options.IterRefine        = DOUBLE;
619      options.Trans             = NOTRANS;
620      options.SolveInitialized  = NO; -hold the communication pattern used MatSolve() and MatMatSolve()
621      options.RefineInitialized = NO;
622      options.PrintStat         = YES;
623   */
624   set_default_options_dist(&options);
625 
626   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(lu->comm_superlu));CHKERRQ(ierr);
627   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
628   /* Default num of process columns and rows */
629   lu->npcol = (int_t) (0.5 + PetscSqrtReal((PetscReal)size));
630   if (!lu->npcol) lu->npcol = 1;
631   while (lu->npcol > 0) {
632     lu->nprow = (int_t) (size/lu->npcol);
633     if (size == lu->nprow * lu->npcol) break;
634     lu->npcol--;
635   }
636 
637   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"SuperLU_Dist Options","Mat");CHKERRQ(ierr);
638   ierr = PetscOptionsInt("-mat_superlu_dist_r","Number rows in processor partition","None",lu->nprow,(PetscInt*)&lu->nprow,NULL);CHKERRQ(ierr);
639   ierr = PetscOptionsInt("-mat_superlu_dist_c","Number columns in processor partition","None",lu->npcol,(PetscInt*)&lu->npcol,NULL);CHKERRQ(ierr);
640   if (size != lu->nprow * lu->npcol) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number of processes %d must equal to nprow %d * npcol %d",size,lu->nprow,lu->npcol);
641 
642   lu->MatInputMode = DISTRIBUTED;
643 
644   ierr = PetscOptionsEnum("-mat_superlu_dist_matinput","Matrix input mode (global or distributed)","None",SuperLU_MatInputModes,(PetscEnum)lu->MatInputMode,(PetscEnum*)&lu->MatInputMode,NULL);CHKERRQ(ierr);
645   if (lu->MatInputMode == DISTRIBUTED && size == 1) lu->MatInputMode = GLOBAL;
646 
647   ierr = PetscOptionsBool("-mat_superlu_dist_equil","Equilibrate matrix","None",options.Equil ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
648   if (set && !flg) options.Equil = NO;
649 
650   ierr = PetscOptionsEList("-mat_superlu_dist_rowperm","Row permutation","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr);
651   if (flg) {
652     switch (indx) {
653     case 0:
654       options.RowPerm = LargeDiag;
655       break;
656     case 1:
657       options.RowPerm = NOROWPERM;
658       break;
659     }
660   }
661 
662   ierr = PetscOptionsEList("-mat_superlu_dist_colperm","Column permutation","None",colperm,5,colperm[3],&indx,&flg);CHKERRQ(ierr);
663   if (flg) {
664     switch (indx) {
665     case 0:
666       options.ColPerm = NATURAL;
667       break;
668     case 1:
669       options.ColPerm = MMD_AT_PLUS_A;
670       break;
671     case 2:
672       options.ColPerm = MMD_ATA;
673       break;
674     case 3:
675       options.ColPerm = METIS_AT_PLUS_A;
676       break;
677     case 4:
678       options.ColPerm = PARMETIS;   /* only works for np>1 */
679       break;
680     default:
681       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown column permutation");
682     }
683   }
684 
685   options.ReplaceTinyPivot = NO;
686   ierr = PetscOptionsBool("-mat_superlu_dist_replacetinypivot","Replace tiny pivots","None",options.ReplaceTinyPivot ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
687   if (set && flg) options.ReplaceTinyPivot = YES;
688 
689   options.ParSymbFact = NO;
690   ierr = PetscOptionsBool("-mat_superlu_dist_parsymbfact","Parallel symbolic factorization","None",PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
691   if (set && flg && size>1) {
692     if (lu->MatInputMode == GLOBAL) {
693 #if defined(PETSC_USE_INFO)
694       ierr = PetscInfo(A,"Warning: '-mat_superlu_dist_parsymbfact' is ignored because MatInputMode=GLOBAL\n");CHKERRQ(ierr);
695 #endif
696     } else {
697 #if defined(PETSC_HAVE_PARMETIS)
698       options.ParSymbFact = YES;
699       options.ColPerm     = PARMETIS;   /* in v2.2, PARMETIS is forced for ParSymbFact regardless of user ordering setting */
700 #else
701       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"parsymbfact needs PARMETIS");
702 #endif
703     }
704   }
705 
706   lu->FactPattern = SamePattern_SameRowPerm;
707   ierr = PetscOptionsEList("-mat_superlu_dist_fact","Sparsity pattern for repeated matrix factorization","None",factPattern,3,factPattern[0],&indx,&flg);CHKERRQ(ierr);
708   if (flg) {
709     switch (indx) {
710     case 0:
711       lu->FactPattern = SamePattern;
712       break;
713     case 1:
714       lu->FactPattern = SamePattern_SameRowPerm;
715       break;
716     case 2:
717       lu->FactPattern = DOFACT;
718       break;
719     }
720   }
721 
722   options.IterRefine = NOREFINE;
723   ierr               = PetscOptionsBool("-mat_superlu_dist_iterrefine","Use iterative refinement","None",options.IterRefine == NOREFINE ? PETSC_FALSE : PETSC_TRUE ,&flg,&set);CHKERRQ(ierr);
724   if (set) {
725     if (flg) options.IterRefine = SLU_DOUBLE;
726     else options.IterRefine = NOREFINE;
727   }
728 
729   if (PetscLogPrintInfo) options.PrintStat = YES;
730   else options.PrintStat = NO;
731   ierr = PetscOptionsBool("-mat_superlu_dist_statprint","Print factorization information","None",(PetscBool)options.PrintStat,(PetscBool*)&options.PrintStat,NULL);CHKERRQ(ierr);
732   ierr = PetscOptionsEnd();CHKERRQ(ierr);
733 
734   lu->options              = options;
735   lu->options.Fact         = DOFACT;
736   lu->matsolve_iscalled    = PETSC_FALSE;
737   lu->matmatsolve_iscalled = PETSC_FALSE;
738 
739   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_aij_superlu_dist);CHKERRQ(ierr);
740   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSuperluDistGetDiagU_C",MatSuperluDistGetDiagU_SuperLU_DIST);CHKERRQ(ierr);
741 
742   *F = B;
743   PetscFunctionReturn(0);
744 }
745 
746 #undef __FUNCT__
747 #define __FUNCT__ "MatSolverPackageRegister_SuperLU_DIST"
748 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_SuperLU_DIST(void)
749 {
750   PetscErrorCode ierr;
751   PetscFunctionBegin;
752   ierr = MatSolverPackageRegister(MATSOLVERSUPERLU_DIST,MATMPIAIJ,  MAT_FACTOR_LU,MatGetFactor_aij_superlu_dist);CHKERRQ(ierr);
753   ierr = MatSolverPackageRegister(MATSOLVERSUPERLU_DIST,MATSEQAIJ,  MAT_FACTOR_LU,MatGetFactor_aij_superlu_dist);CHKERRQ(ierr);
754   PetscFunctionReturn(0);
755 }
756 
757 /*MC
758   MATSOLVERSUPERLU_DIST - Parallel direct solver package for LU factorization
759 
760   Use ./configure --download-superlu_dist --download-parmetis --download-metis --download-ptscotch  to have PETSc installed with SuperLU_DIST
761 
762   Use -pc_type lu -pc_factor_mat_solver_package superlu_dist to us this direct solver
763 
764    Works with AIJ matrices
765 
766   Options Database Keys:
767 + -mat_superlu_dist_r <n> - number of rows in processor partition
768 . -mat_superlu_dist_c <n> - number of columns in processor partition
769 . -mat_superlu_dist_matinput <0,1> - matrix input mode; 0=global, 1=distributed
770 . -mat_superlu_dist_equil - equilibrate the matrix
771 . -mat_superlu_dist_rowperm <LargeDiag,NATURAL> - row permutation
772 . -mat_superlu_dist_colperm <MMD_AT_PLUS_A,MMD_ATA,NATURAL> - column permutation
773 . -mat_superlu_dist_replacetinypivot - replace tiny pivots
774 . -mat_superlu_dist_fact <SamePattern> - (choose one of) SamePattern SamePattern_SameRowPerm
775 . -mat_superlu_dist_iterrefine - use iterative refinement
776 - -mat_superlu_dist_statprint - print factorization information
777 
778    Level: beginner
779 
780 .seealso: PCLU
781 
782 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
783 
784 M*/
785 
786