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