xref: /petsc/src/mat/impls/aij/mpi/superlu_dist/superlu_dist.c (revision 201f5f94dcee918e9f0ce123f34f0acb7c97f68b)
1 /*$Id: superlu_DIST.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/
2 /*
3         Provides an interface to the SuperLU_DIST_2.0 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 arround weird problem with SuperLU on cray */
9 #include "stdlib.h"
10 #endif
11 
12 /*MC
13   SuperLU - an external package providing director solvers for MPIAIJ
14   matrices.
15 
16   Notes: PETSc does not natively support PCLU for MPIAIJ
17   matrices. However, if SuperLU_dist is installed (see the manual for
18   instructions on how to declare the existence of external packages),
19   it can be specified as the solver to be used by a call to
20   MatUseSuperLU_DIST_MPIAIJ(). This is only supported for
21   double precision real.
22 
23   Options Database Keys:
24 + -mat_aij_superlu_dist_r <n> : number of rows in processor partition
25 . -mat_aij_superlu_dist_c <n> : number of columns in processor partition
26 . -mat_aij_superlu_dist_matinput 0|1 : matrix input mode; 0=global, 1=distributed
27 . -mat_aij_superlu_dist_equil :, equilibrate the matrix
28 . -mat_aij_superlu_dist_rowperm LargeDiag|NATURAL : row permutation
29 . -mat_aij_superlu_dist_colperm MMD_AT_PLUS_A|MMD_ATA|COLAMD|NATURAL : column permutation
30 . -mat_aij_superlu_dist_replacetinypivot : replace tiny pivots
31 . -mat_aij_superlu_dist_iterrefine : use iterative refinement
32 - -mat_aij_superlu_dist_statprint : print factorization information
33 
34 .seealso: MatUseSuperLU_DIST_MPIAIJ(), PCLU
35 M*/
36 
37 EXTERN_C_BEGIN
38 #if defined(PETSC_USE_COMPLEX)
39 #include "superlu_zdefs.h"
40 #else
41 #include "superlu_ddefs.h"
42 #endif
43 EXTERN_C_END
44 
45 typedef enum { GLOBAL,DISTRIBUTED
46 } SuperLU_MatInputMode;
47 
48 typedef struct {
49   int_t                   nprow,npcol,*row,*col;
50   gridinfo_t              grid;
51   superlu_options_t       options;
52   SuperMatrix             A_sup;
53   ScalePermstruct_t       ScalePermstruct;
54   LUstruct_t              LUstruct;
55   int                     StatPrint;
56   int                     MatInputMode;
57   SOLVEstruct_t           SOLVEstruct;
58   MatStructure            flg;
59   MPI_Comm                comm_superlu;
60 #if defined(PETSC_USE_COMPLEX)
61   doublecomplex           *val;
62 #else
63   double                  *val;
64 #endif
65 
66   /* A few function pointers for inheritance */
67   int (*MatView)(Mat,PetscViewer);
68   int (*MatAssemblyEnd)(Mat,MatAssemblyType);
69   int (*MatDestroy)(Mat);
70 
71   /* Flag to clean up (non-global) SuperLU objects during Destroy */
72   PetscTruth CleanUpSuperLUDist;
73 } Mat_MPIAIJ_SuperLU_DIST;
74 
75 #undef __FUNCT__
76 #define __FUNCT__ "MatDestroy_MPIAIJ_SuperLU_DIST"
77 int MatDestroy_MPIAIJ_SuperLU_DIST(Mat A)
78 {
79   Mat_MPIAIJ              *a  = (Mat_MPIAIJ*)A->data;
80   Mat_MPIAIJ_SuperLU_DIST *lu = (Mat_MPIAIJ_SuperLU_DIST*)A->spptr;
81   int                     ierr, size=a->size,(*destroy)(Mat);
82 
83   PetscFunctionBegin;
84   if (lu->CleanUpSuperLUDist) {
85     /* Deallocate SuperLU_DIST storage */
86     if (lu->MatInputMode == GLOBAL) {
87       Destroy_CompCol_Matrix_dist(&lu->A_sup);
88     } else {
89       Destroy_CompRowLoc_Matrix_dist(&lu->A_sup);
90       if ( lu->options.SolveInitialized ) {
91 #if defined(PETSC_USE_COMPLEX)
92         zSolveFinalize(&lu->options, &lu->SOLVEstruct);
93 #else
94         dSolveFinalize(&lu->options, &lu->SOLVEstruct);
95 #endif
96       }
97     }
98     Destroy_LU(A->N, &lu->grid, &lu->LUstruct);
99     ScalePermstructFree(&lu->ScalePermstruct);
100     LUstructFree(&lu->LUstruct);
101 
102     /* Release the SuperLU_DIST process grid. */
103     superlu_gridexit(&lu->grid);
104 
105     ierr = MPI_Comm_free(&(lu->comm_superlu));CHKERRQ(ierr);
106   }
107   destroy = lu->MatDestroy;
108   ierr = PetscFree(lu);CHKERRQ(ierr);
109   ierr = (*destroy)(A);CHKERRQ(ierr);
110 
111   PetscFunctionReturn(0);
112 }
113 
114 extern int MatMPIAIJFactorInfo_SuperLu(Mat A,PetscViewer viewer);
115 
116 #undef __FUNCT__
117 #define __FUNCT__ "MatView_MPIAIJ_Spooles_DIST"
118 int MatView_MPIAIJ_SuperLU_DIST(Mat A,PetscViewer viewer)
119 {
120   int                     ierr;
121   PetscTruth              isascii;
122   PetscViewerFormat       format;
123   Mat_MPIAIJ_SuperLU_DIST *lu=(Mat_MPIAIJ_SuperLU_DIST*)(A->spptr);
124 
125   PetscFunctionBegin;
126   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
127 
128   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
129   if (isascii) {
130     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
131     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
132       ierr = MatMPIAIJFactorInfo_SuperLu(A,viewer);CHKERRQ(ierr);
133     }
134   }
135   PetscFunctionReturn(0);
136 }
137 
138 #undef __FUNCT__
139 #define __FUNCT__ "MatAssemblyEnd_MPIAIJ_SuperLU_DIST"
140 int MatAssemblyEnd_MPIAIJ_SuperLU_DIST(Mat A,MatAssemblyType mode) {
141   int                     ierr;
142   Mat_MPIAIJ_SuperLU_DIST *lu=(Mat_MPIAIJ_SuperLU_DIST*)(A->spptr);
143 
144   PetscFunctionBegin;
145   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
146   ierr = MatUseSuperLU_DIST_MPIAIJ(A);CHKERRQ(ierr);
147   PetscFunctionReturn(0);
148 }
149 
150 #undef __FUNCT__
151 #define __FUNCT__ "MatSolve_MPIAIJ_SuperLU_DIST"
152 int MatSolve_MPIAIJ_SuperLU_DIST(Mat A,Vec b_mpi,Vec x)
153 {
154   Mat_MPIAIJ              *aa = (Mat_MPIAIJ*)A->data;
155   Mat_MPIAIJ_SuperLU_DIST *lu = (Mat_MPIAIJ_SuperLU_DIST*)A->spptr;
156   int                     ierr, size=aa->size;
157   int                     m=A->M, N=A->N;
158   SuperLUStat_t           stat;
159   double                  berr[1];
160   PetscScalar             *bptr;
161   int                     info, nrhs=1;
162   Vec                     x_seq;
163   IS                      iden;
164   VecScatter              scat;
165   PetscLogDouble          time0,time,time_min,time_max;
166 
167   PetscFunctionBegin;
168   if (size > 1) {
169     if (lu->MatInputMode == GLOBAL) { /* global mat input, convert b to x_seq */
170       ierr = VecCreateSeq(PETSC_COMM_SELF,N,&x_seq);CHKERRQ(ierr);
171       ierr = ISCreateStride(PETSC_COMM_SELF,N,0,1,&iden);CHKERRQ(ierr);
172       ierr = VecScatterCreate(b_mpi,iden,x_seq,iden,&scat);CHKERRQ(ierr);
173       ierr = ISDestroy(iden);CHKERRQ(ierr);
174 
175       ierr = VecScatterBegin(b_mpi,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
176       ierr = VecScatterEnd(b_mpi,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
177       ierr = VecGetArray(x_seq,&bptr);CHKERRQ(ierr);
178     } else { /* distributed mat input */
179       ierr = VecCopy(b_mpi,x);CHKERRQ(ierr);
180       ierr = VecGetArray(x,&bptr);CHKERRQ(ierr);
181     }
182   } else { /* size == 1 */
183     ierr = VecCopy(b_mpi,x);CHKERRQ(ierr);
184     ierr = VecGetArray(x,&bptr);CHKERRQ(ierr);
185   }
186 
187   lu->options.Fact = FACTORED; /* The factored form of A is supplied. Local option used by this func. only.*/
188 
189   PStatInit(&stat);        /* Initialize the statistics variables. */
190   if (lu->StatPrint) {
191     ierr = MPI_Barrier(A->comm);CHKERRQ(ierr); /* to be removed */
192     ierr = PetscGetTime(&time0);CHKERRQ(ierr);  /* to be removed */
193   }
194   if (lu->MatInputMode == GLOBAL) {
195 #if defined(PETSC_USE_COMPLEX)
196     pzgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct,(doublecomplex*)bptr, m, nrhs,
197                    &lu->grid, &lu->LUstruct, berr, &stat, &info);
198 #else
199     pdgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct,bptr, m, nrhs,
200                    &lu->grid, &lu->LUstruct, berr, &stat, &info);
201 #endif
202   } else { /* distributed mat input */
203 #if defined(PETSC_USE_COMPLEX)
204     pzgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, (doublecomplex*)bptr, A->M, nrhs, &lu->grid,
205 	    &lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &info);
206     if (info) SETERRQ1(1,"pzgssvx fails, info: %d\n",info);
207 #else
208     pdgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, bptr, A->M, nrhs, &lu->grid,
209 	    &lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &info);
210     if (info) SETERRQ1(1,"pdgssvx fails, info: %d\n",info);
211 #endif
212   }
213   if (lu->StatPrint) {
214     ierr = PetscGetTime(&time);CHKERRQ(ierr);  /* to be removed */
215      PStatPrint(&lu->options, &stat, &lu->grid);     /* Print the statistics. */
216   }
217   PStatFree(&stat);
218 
219   if (size > 1) {
220     if (lu->MatInputMode == GLOBAL){ /* convert seq x to mpi x */
221       ierr = VecRestoreArray(x_seq,&bptr);CHKERRQ(ierr);
222       ierr = VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
223       ierr = VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
224       ierr = VecScatterDestroy(scat);CHKERRQ(ierr);
225       ierr = VecDestroy(x_seq);CHKERRQ(ierr);
226     } else {
227       ierr = VecRestoreArray(x,&bptr);CHKERRQ(ierr);
228     }
229   } else {
230     ierr = VecRestoreArray(x,&bptr);CHKERRQ(ierr);
231   }
232   if (lu->StatPrint) {
233     time0 = time - time0;
234     ierr = MPI_Reduce(&time0,&time_max,1,MPI_DOUBLE,MPI_MAX,0,A->comm);CHKERRQ(ierr);
235     ierr = MPI_Reduce(&time0,&time_min,1,MPI_DOUBLE,MPI_MIN,0,A->comm);CHKERRQ(ierr);
236     ierr = MPI_Reduce(&time0,&time,1,MPI_DOUBLE,MPI_SUM,0,A->comm);CHKERRQ(ierr);
237     time = time/size; /* average time */
238     ierr = PetscPrintf(A->comm, "  Time for superlu_dist solve (max/min/avg): %g / %g / %g\n\n",time_max,time_min,time);CHKERRQ(ierr);
239   }
240   PetscFunctionReturn(0);
241 }
242 
243 #undef __FUNCT__
244 #define __FUNCT__ "MatLUFactorNumeric_MPIAIJ_SuperLU_DIST"
245 int MatLUFactorNumeric_MPIAIJ_SuperLU_DIST(Mat A,Mat *F)
246 {
247   Mat_MPIAIJ              *fac = (Mat_MPIAIJ*)(*F)->data,*mat;
248   Mat                     *tseq,A_seq = PETSC_NULL;
249   Mat_SeqAIJ              *aa,*bb;
250   Mat_MPIAIJ_SuperLU_DIST *lu = (Mat_MPIAIJ_SuperLU_DIST*)(*F)->spptr;
251   int                     M=A->M,N=A->N,info,ierr,size=fac->size,i,*ai,*aj,*bi,*bj,nz,rstart,*garray,
252                           m=A->m, irow,colA_start,j,jcol,jB,countA,countB,*bjj,*ajj;
253   SuperLUStat_t           stat;
254   double                  *berr=0;
255   IS                      isrow;
256   PetscLogDouble          time0[2],time[2],time_min[2],time_max[2];
257 #if defined(PETSC_USE_COMPLEX)
258   doublecomplex           *av, *bv;
259 #else
260   double                  *av, *bv;
261 #endif
262 
263   PetscFunctionBegin;
264   if (lu->StatPrint) {
265     ierr = MPI_Barrier(A->comm);CHKERRQ(ierr);
266     ierr = PetscGetTime(&time0[0]);CHKERRQ(ierr);
267   }
268 
269   if (lu->MatInputMode == GLOBAL) { /* global mat input */
270     if (size > 1) { /* convert mpi A to seq mat A */
271       ierr = ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow); CHKERRQ(ierr);
272       ierr = MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&tseq); CHKERRQ(ierr);
273       ierr = ISDestroy(isrow);CHKERRQ(ierr);
274 
275       A_seq = *tseq;
276       ierr = PetscFree(tseq);CHKERRQ(ierr);
277       aa =  (Mat_SeqAIJ*)A_seq->data;
278     } else {
279       aa =  (Mat_SeqAIJ*)A->data;
280     }
281 
282     /* Allocate storage, then convert Petsc NR matrix to SuperLU_DIST NC */
283     if (lu->flg == DIFFERENT_NONZERO_PATTERN) {/* first numeric factorization */
284 #if defined(PETSC_USE_COMPLEX)
285       zallocateA_dist(N, aa->nz, &lu->val, &lu->col, &lu->row);
286 #else
287       dallocateA_dist(N, aa->nz, &lu->val, &lu->col, &lu->row);
288 #endif
289     } else { /* successive numeric factorization, sparsity pattern is reused. */
290       Destroy_CompCol_Matrix_dist(&lu->A_sup);
291       Destroy_LU(N, &lu->grid, &lu->LUstruct);
292       lu->options.Fact = SamePattern;
293     }
294 #if defined(PETSC_USE_COMPLEX)
295     zCompRow_to_CompCol_dist(M,N,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,&lu->val,&lu->col, &lu->row);
296 #else
297     dCompRow_to_CompCol_dist(M,N,aa->nz,aa->a,aa->j,aa->i,&lu->val, &lu->col, &lu->row);
298 #endif
299 
300     /* Create compressed column matrix A_sup. */
301 #if defined(PETSC_USE_COMPLEX)
302     zCreate_CompCol_Matrix_dist(&lu->A_sup, M, N, aa->nz, lu->val, lu->col, lu->row, SLU_NC, SLU_Z, SLU_GE);
303 #else
304     dCreate_CompCol_Matrix_dist(&lu->A_sup, M, N, aa->nz, lu->val, lu->col, lu->row, SLU_NC, SLU_D, SLU_GE);
305 #endif
306   } else { /* distributed mat input */
307     mat =  (Mat_MPIAIJ*)A->data;
308     aa=(Mat_SeqAIJ*)(mat->A)->data;
309     bb=(Mat_SeqAIJ*)(mat->B)->data;
310     ai=aa->i; aj=aa->j;
311     bi=bb->i; bj=bb->j;
312 #if defined(PETSC_USE_COMPLEX)
313     av=(doublecomplex*)aa->a;
314     bv=(doublecomplex*)bb->a;
315 #else
316     av=aa->a;
317     bv=bb->a;
318 #endif
319     rstart = mat->rstart;
320     nz     = aa->nz + bb->nz;
321     garray = mat->garray;
322     rstart = mat->rstart;
323 
324     if (lu->flg == DIFFERENT_NONZERO_PATTERN) {/* first numeric factorization */
325 #if defined(PETSC_USE_COMPLEX)
326       zallocateA_dist(m, nz, &lu->val, &lu->col, &lu->row);
327 #else
328       dallocateA_dist(m, nz, &lu->val, &lu->col, &lu->row);
329 #endif
330     } else { /* successive numeric factorization, sparsity pattern and perm_c are reused. */
331       /* Destroy_CompRowLoc_Matrix_dist(&lu->A_sup);  */ /* crash! */
332       Destroy_LU(N, &lu->grid, &lu->LUstruct);
333       lu->options.Fact = SamePattern;
334     }
335     nz = 0; jB = 0; irow = mat->rstart;
336     for ( i=0; i<m; i++ ) {
337       lu->row[i] = nz;
338       countA = ai[i+1] - ai[i];
339       countB = bi[i+1] - bi[i];
340       ajj = aj + ai[i];  /* ptr to the beginning of this row */
341       bjj = bj + bi[i];
342 
343       /* B part, smaller col index */
344       colA_start = mat->rstart + ajj[0]; /* the smallest global col index of A */
345       for (j=0; j<countB; j++){
346         jcol = garray[bjj[j]];
347         if (jcol > colA_start) {
348           jB = j;
349           break;
350         }
351         lu->col[nz] = jcol;
352         lu->val[nz++] = *bv++;
353         if (j==countB-1) jB = countB;
354       }
355 
356       /* A part */
357       for (j=0; j<countA; j++){
358         lu->col[nz] = mat->rstart + ajj[j];
359         lu->val[nz++] = *av++;
360       }
361 
362       /* B part, larger col index */
363       for (j=jB; j<countB; j++){
364         lu->col[nz] = garray[bjj[j]];
365         lu->val[nz++] = *bv++;
366       }
367     }
368     lu->row[m] = nz;
369 #if defined(PETSC_USE_COMPLEX)
370     zCreate_CompRowLoc_Matrix_dist(&lu->A_sup, M, N, nz, m, rstart,
371 				   lu->val, lu->col, lu->row, SLU_NR_loc, SLU_Z, SLU_GE);
372 #else
373     dCreate_CompRowLoc_Matrix_dist(&lu->A_sup, M, N, nz, m, rstart,
374 				   lu->val, lu->col, lu->row, SLU_NR_loc, SLU_D, SLU_GE);
375 #endif
376   }
377   if (lu->StatPrint) {
378     ierr = PetscGetTime(&time[0]);CHKERRQ(ierr);
379     time0[0] = time[0] - time0[0];
380   }
381 
382   /* Factor the matrix. */
383   PStatInit(&stat);   /* Initialize the statistics variables. */
384 
385   if (lu->StatPrint) {
386     ierr = MPI_Barrier(A->comm);CHKERRQ(ierr);
387     ierr = PetscGetTime(&time0[1]);CHKERRQ(ierr);
388   }
389 
390   if (lu->MatInputMode == GLOBAL) { /* global mat input */
391 #if defined(PETSC_USE_COMPLEX)
392     pzgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, M, 0,
393                    &lu->grid, &lu->LUstruct, berr, &stat, &info);
394 #else
395     pdgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, M, 0,
396                    &lu->grid, &lu->LUstruct, berr, &stat, &info);
397 #endif
398   } else { /* distributed mat input */
399 #if defined(PETSC_USE_COMPLEX)
400     pzgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, M, 0, &lu->grid,
401 	    &lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &info);
402     if (info) SETERRQ1(1,"pzgssvx fails, info: %d\n",info);
403 #else
404     pdgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, M, 0, &lu->grid,
405 	    &lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &info);
406     if (info) SETERRQ1(1,"pdgssvx fails, info: %d\n",info);
407 #endif
408   }
409   if (lu->StatPrint) {
410     ierr = PetscGetTime(&time[1]);CHKERRQ(ierr);  /* to be removed */
411     time0[1] = time[1] - time0[1];
412     if (lu->StatPrint) PStatPrint(&lu->options, &stat, &lu->grid);  /* Print the statistics. */
413   }
414   PStatFree(&stat);
415 
416   if (lu->MatInputMode == GLOBAL && size > 1){
417     ierr = MatDestroy(A_seq);CHKERRQ(ierr);
418   }
419 
420   if (lu->StatPrint) {
421     ierr = MPI_Reduce(time0,time_max,2,MPI_DOUBLE,MPI_MAX,0,A->comm);
422     ierr = MPI_Reduce(time0,time_min,2,MPI_DOUBLE,MPI_MIN,0,A->comm);
423     ierr = MPI_Reduce(time0,time,2,MPI_DOUBLE,MPI_SUM,0,A->comm);
424     for (i=0; i<2; i++) time[i] = time[i]/size; /* average time */
425     ierr = PetscPrintf(A->comm, "  Time for mat conversion (max/min/avg):    %g / %g / %g\n",time_max[0],time_min[0],time[0]);
426     ierr = PetscPrintf(A->comm, "  Time for superlu_dist fact (max/min/avg): %g / %g / %g\n\n",time_max[1],time_min[1],time[1]);
427   }
428   (*F)->assembled = PETSC_TRUE;
429   lu->flg         = SAME_NONZERO_PATTERN;
430   PetscFunctionReturn(0);
431 }
432 
433 /* Note the Petsc r and c permutations are ignored */
434 #undef __FUNCT__
435 #define __FUNCT__ "MatLUFactorSymbolic_MPIAIJ_SuperLU_DIST"
436 int MatLUFactorSymbolic_MPIAIJ_SuperLU_DIST(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
437 {
438   Mat                     B;
439   Mat_MPIAIJ_SuperLU_DIST *lu;
440   int                     ierr,M=A->M,N=A->N,size;
441   superlu_options_t       options;
442   char                    buff[32];
443   PetscTruth              flg;
444   char                    *ptype[] = {"MMD_AT_PLUS_A","NATURAL","MMD_ATA","COLAMD"};
445   char                    *prtype[] = {"LargeDiag","NATURAL"};
446   PetscFunctionBegin;
447 
448   /* Create the factorization matrix */
449   ierr = MatCreate(A->comm,A->m,A->n,M,N,&B);CHKERRQ(ierr);
450   ierr = MatSetType(B,MATSUPERLUDIST);CHKERRQ(ierr);
451   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);
452   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
453 
454   B->ops->solve            = MatSolve_MPIAIJ_SuperLU_DIST;
455   B->factor                = FACTOR_LU;
456 
457   lu = (Mat_MPIAIJ_SuperLU_DIST*)(B->spptr);
458 
459   /* Set the input options */
460   set_default_options(&options);
461   lu->MatInputMode = GLOBAL;
462   ierr = MPI_Comm_dup(A->comm,&(lu->comm_superlu));CHKERRQ(ierr);
463 
464   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
465   lu->nprow = size/2;               /* Default process rows.      */
466   if (lu->nprow == 0) lu->nprow = 1;
467   lu->npcol = size/lu->nprow;           /* Default process columns.   */
468 
469   ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU_Dist Options","Mat");CHKERRQ(ierr);
470 
471     ierr = PetscOptionsInt("-mat_aij_superlu_dist_r","Number rows in processor partition","None",lu->nprow,&lu->nprow,PETSC_NULL);CHKERRQ(ierr);
472     ierr = PetscOptionsInt("-mat_aij_superlu_dist_c","Number columns in processor partition","None",lu->npcol,&lu->npcol,PETSC_NULL);CHKERRQ(ierr);
473     if (size != lu->nprow * lu->npcol) SETERRQ(1,"Number of processes should be equal to nprow*npcol");
474 
475     ierr = PetscOptionsInt("-mat_aij_superlu_dist_matinput","Matrix input mode (0: GLOBAL; 1: DISTRIBUTED)","None",lu->MatInputMode,&lu->MatInputMode,PETSC_NULL);CHKERRQ(ierr);
476     if(lu->MatInputMode == DISTRIBUTED && size == 1) lu->MatInputMode = GLOBAL;
477 
478     ierr = PetscOptionsLogical("-mat_aij_superlu_dist_equil","Equilibrate matrix","None",PETSC_TRUE,&flg,0);CHKERRQ(ierr);
479     if (!flg) {
480       options.Equil = NO;
481     }
482 
483     ierr = PetscOptionsEList("-mat_aij_superlu_dist_rowperm","Row permutation","None",prtype,2,prtype[0],buff,32,&flg);CHKERRQ(ierr);
484     while (flg) {
485       ierr = PetscStrcmp(buff,"LargeDiag",&flg);CHKERRQ(ierr);
486       if (flg) {
487         options.RowPerm = LargeDiag;
488         break;
489       }
490       ierr = PetscStrcmp(buff,"NATURAL",&flg);CHKERRQ(ierr);
491       if (flg) {
492         options.RowPerm = NOROWPERM;
493         break;
494       }
495       SETERRQ1(1,"Unknown row permutation %s",buff);
496     }
497 
498     ierr = PetscOptionsEList("-mat_aij_superlu_dist_colperm","Column permutation","None",ptype,4,ptype[0],buff,32,&flg);CHKERRQ(ierr);
499     while (flg) {
500       ierr = PetscStrcmp(buff,"MMD_AT_PLUS_A",&flg);CHKERRQ(ierr);
501       if (flg) {
502         options.ColPerm = MMD_AT_PLUS_A;
503         break;
504       }
505       ierr = PetscStrcmp(buff,"NATURAL",&flg);CHKERRQ(ierr);
506       if (flg) {
507         options.ColPerm = NATURAL;
508         break;
509       }
510       ierr = PetscStrcmp(buff,"MMD_ATA",&flg);CHKERRQ(ierr);
511       if (flg) {
512         options.ColPerm = MMD_ATA;
513         break;
514       }
515       ierr = PetscStrcmp(buff,"COLAMD",&flg);CHKERRQ(ierr);
516       if (flg) {
517         options.ColPerm = COLAMD;
518         break;
519       }
520       SETERRQ1(1,"Unknown column permutation %s",buff);
521     }
522 
523     ierr = PetscOptionsLogical("-mat_aij_superlu_dist_replacetinypivot","Replace tiny pivots","None",PETSC_TRUE,&flg,0);CHKERRQ(ierr);
524     if (!flg) {
525       options.ReplaceTinyPivot = NO;
526     }
527 
528     options.IterRefine = NOREFINE;
529     ierr = PetscOptionsLogical("-mat_aij_superlu_dist_iterrefine","Use iterative refinement","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
530     if (flg) {
531       options.IterRefine = DOUBLE;
532     }
533 
534     if (PetscLogPrintInfo) {
535       lu->StatPrint = (int)PETSC_TRUE;
536     } else {
537       lu->StatPrint = (int)PETSC_FALSE;
538     }
539     ierr = PetscOptionsLogical("-mat_aij_superlu_dist_statprint","Print factorization information","None",
540                               (PetscTruth)lu->StatPrint,(PetscTruth*)&lu->StatPrint,0);CHKERRQ(ierr);
541   PetscOptionsEnd();
542 
543   /* Initialize the SuperLU process grid. */
544   superlu_gridinit(lu->comm_superlu, lu->nprow, lu->npcol, &lu->grid);
545 
546   /* Initialize ScalePermstruct and LUstruct. */
547   ScalePermstructInit(M, N, &lu->ScalePermstruct);
548   LUstructInit(M, N, &lu->LUstruct);
549 
550   lu->options            = options;
551   lu->flg                = DIFFERENT_NONZERO_PATTERN;
552   lu->CleanUpSuperLUDist = PETSC_TRUE;
553   *F = B;
554   PetscFunctionReturn(0);
555 }
556 
557 #undef __FUNCT__
558 #define __FUNCT__ "MatUseSuperLU_DIST_MPIAIJ"
559 /*@
560   MatUseSuperLU_DIST_MPIAIJ -- Let the distributed SuperLU package
561   be the solver used by PCLU
562 
563   Input Parameter:
564 . A - the matrix
565 
566   Level: intermediate
567 
568 .seealso: SuperLU, PCSetType(), PCLU
569 @*/
570 int MatUseSuperLU_DIST_MPIAIJ(Mat A)
571 {
572   PetscFunctionBegin;
573   A->ops->lufactorsymbolic = MatLUFactorSymbolic_MPIAIJ_SuperLU_DIST;
574   A->ops->lufactornumeric  = MatLUFactorNumeric_MPIAIJ_SuperLU_DIST;
575   PetscFunctionReturn(0);
576 }
577 
578 #undef __FUNCT__
579 #define __FUNCT__ "MatMPIAIJFactorInfo_SuperLu"
580 int MatMPIAIJFactorInfo_SuperLu(Mat A,PetscViewer viewer)
581 {
582   Mat_MPIAIJ_SuperLU_DIST *lu= (Mat_MPIAIJ_SuperLU_DIST*)A->spptr;
583   superlu_options_t       options;
584   int                     ierr;
585   char                    *colperm;
586 
587   PetscFunctionBegin;
588   /* check if matrix is superlu_dist type */
589   if (A->ops->solve != MatSolve_MPIAIJ_SuperLU_DIST) PetscFunctionReturn(0);
590 
591   options = lu->options;
592   ierr = PetscViewerASCIIPrintf(viewer,"SuperLU_DIST run parameters:\n");CHKERRQ(ierr);
593   ierr = PetscViewerASCIIPrintf(viewer,"  Equilibrate matrix %s \n",(options.Equil != NO) ? "true": "false");CHKERRQ(ierr);
594   ierr = PetscViewerASCIIPrintf(viewer,"  Replace tiny pivots %s \n",(options.ReplaceTinyPivot != NO) ? "true": "false");CHKERRQ(ierr);
595   ierr = PetscViewerASCIIPrintf(viewer,"  Use iterative refinement %s \n",(options.IterRefine == DOUBLE) ? "true": "false");CHKERRQ(ierr);
596   ierr = PetscViewerASCIIPrintf(viewer,"  Processors in row %d col partition %d \n",lu->nprow,lu->npcol);CHKERRQ(ierr);
597   ierr = PetscViewerASCIIPrintf(viewer,"  Row permutation %s \n",(options.RowPerm == NOROWPERM) ? "NATURAL": "LargeDiag");CHKERRQ(ierr);
598   if (options.ColPerm == NATURAL) {
599     colperm = "NATURAL";
600   } else if (options.ColPerm == MMD_AT_PLUS_A) {
601     colperm = "MMD_AT_PLUS_A";
602   } else if (options.ColPerm == MMD_ATA) {
603     colperm = "MMD_ATA";
604   } else if (options.ColPerm == COLAMD) {
605     colperm = "COLAMD";
606   } else {
607     SETERRQ(1,"Unknown column permutation");
608   }
609   ierr = PetscViewerASCIIPrintf(viewer,"  Column permutation %s \n",colperm);CHKERRQ(ierr);
610   PetscFunctionReturn(0);
611 }
612 
613 EXTERN_C_BEGIN
614 #undef __FUNCT__
615 #define __FUNCT__ "MatCreate_MPIAIJ_SuperLU_DIST"
616 int MatCreate_MPIAIJ_SuperLU_DIST(Mat A) {
617   int                     ierr,size;
618   MPI_Comm                comm;
619   Mat_MPIAIJ_SuperLU_DIST *lu;
620 
621   PetscFunctionBegin;
622   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
623   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
624   if (size == 1) {
625     ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
626   } else {
627     ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
628   }
629   ierr = MatUseSuperLU_DIST_MPIAIJ(A);CHKERRQ(ierr);
630 
631   ierr                   = PetscNew(Mat_MPIAIJ_SuperLU_DIST,&lu);CHKERRQ(ierr);
632   lu->MatView            = A->ops->view;
633   lu->MatAssemblyEnd     = A->ops->assemblyend;
634   lu->MatDestroy         = A->ops->destroy;
635   lu->CleanUpSuperLUDist = PETSC_FALSE;
636   A->spptr               = (void*)lu;
637   A->ops->view           = MatView_MPIAIJ_SuperLU_DIST;
638   A->ops->assemblyend    = MatAssemblyEnd_MPIAIJ_SuperLU_DIST;
639   A->ops->destroy        = MatDestroy_MPIAIJ_SuperLU_DIST;
640   PetscFunctionReturn(0);
641 }
642 EXTERN_C_END
643 
644 EXTERN_C_BEGIN
645 #undef __FUNCT__
646 #define __FUNCT__ "MatLoad_MPIAIJ_SuperLU_DIST"
647 int MatLoad_MPIAIJ_SuperLU_DIST(PetscViewer viewer,MatType type,Mat *A) {
648   int      ierr,size,(*r)(PetscViewer,MatType,Mat*);
649   MPI_Comm comm;
650 
651   PetscFunctionBegin;
652   ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
653   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
654   if (size == 1) {
655     ierr = PetscFListFind(comm,MatLoadList,MATSEQAIJ,(void(**)(void))&r);CHKERRQ(ierr);
656   } else {
657     ierr = PetscFListFind(comm,MatLoadList,MATMPIAIJ,(void(**)(void))&r);CHKERRQ(ierr);
658   }
659   ierr = (*r)(viewer,type,A);CHKERRQ(ierr);
660   PetscFunctionReturn(0);
661 }
662 EXTERN_C_END
663