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