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