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