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