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