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