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