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