xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision 2b8d69ca7ea5fe9190df62c1dce3bbd66fce84dd)
1 
2 /*  --------------------------------------------------------------------
3 
4      This file implements a subclass of the SeqAIJ matrix class that uses
5      the SuperLU sparse solver.
6 */
7 
8 /*
9      Defines the data structure for the base matrix type (SeqAIJ)
10 */
11 #include <../src/mat/impls/aij/seq/aij.h>    /*I "petscmat.h" I*/
12 
13 /*
14      SuperLU include files
15 */
16 EXTERN_C_BEGIN
17 #if defined(PETSC_USE_COMPLEX)
18 #if defined(PETSC_USE_REAL_SINGLE)
19 #include <slu_cdefs.h>
20 #else
21 #include <slu_zdefs.h>
22 #endif
23 #else
24 #if defined(PETSC_USE_REAL_SINGLE)
25 #include <slu_sdefs.h>
26 #else
27 #include <slu_ddefs.h>
28 #endif
29 #endif
30 #include <slu_util.h>
31 EXTERN_C_END
32 
33 /*
34      This is the data we are "ADDING" to the SeqAIJ matrix type to get the SuperLU matrix type
35 */
36 typedef struct {
37   SuperMatrix       A,L,U,B,X;
38   superlu_options_t options;
39   PetscInt          *perm_c; /* column permutation vector */
40   PetscInt          *perm_r; /* row permutations from partial pivoting */
41   PetscInt          *etree;
42   PetscReal         *R, *C;
43   char              equed[1];
44   PetscInt          lwork;
45   void              *work;
46   PetscReal         rpg, rcond;
47   mem_usage_t       mem_usage;
48   MatStructure      flg;
49   SuperLUStat_t     stat;
50   Mat               A_dup;
51   PetscScalar       *rhs_dup;
52   GlobalLU_t        Glu;
53 
54   /* Flag to clean up (non-global) SuperLU objects during Destroy */
55   PetscBool CleanUpSuperLU;
56 } Mat_SuperLU;
57 
58 extern PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer);
59 extern PetscErrorCode MatLUFactorNumeric_SuperLU(Mat,Mat,const MatFactorInfo*);
60 extern PetscErrorCode MatDestroy_SuperLU(Mat);
61 extern PetscErrorCode MatView_SuperLU(Mat,PetscViewer);
62 extern PetscErrorCode MatAssemblyEnd_SuperLU(Mat,MatAssemblyType);
63 extern PetscErrorCode MatSolve_SuperLU(Mat,Vec,Vec);
64 extern PetscErrorCode MatMatSolve_SuperLU(Mat,Mat,Mat);
65 extern PetscErrorCode MatSolveTranspose_SuperLU(Mat,Vec,Vec);
66 extern PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,Mat,IS,IS,const MatFactorInfo*);
67 extern PetscErrorCode MatDuplicate_SuperLU(Mat, MatDuplicateOption, Mat*);
68 
69 /*
70     Utility function
71 */
72 #undef __FUNCT__
73 #define __FUNCT__ "MatFactorInfo_SuperLU"
74 PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
75 {
76   Mat_SuperLU       *lu= (Mat_SuperLU*)A->spptr;
77   PetscErrorCode    ierr;
78   superlu_options_t options;
79 
80   PetscFunctionBegin;
81   /* check if matrix is superlu_dist type */
82   if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0);
83 
84   options = lu->options;
85 
86   ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr);
87   ierr = PetscViewerASCIIPrintf(viewer,"  Equil: %s\n",(options.Equil != NO) ? "YES" : "NO");CHKERRQ(ierr);
88   ierr = PetscViewerASCIIPrintf(viewer,"  ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr);
89   ierr = PetscViewerASCIIPrintf(viewer,"  IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr);
90   ierr = PetscViewerASCIIPrintf(viewer,"  SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES" : "NO");CHKERRQ(ierr);
91   ierr = PetscViewerASCIIPrintf(viewer,"  DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr);
92   ierr = PetscViewerASCIIPrintf(viewer,"  PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES" : "NO");CHKERRQ(ierr);
93   ierr = PetscViewerASCIIPrintf(viewer,"  ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES" : "NO");CHKERRQ(ierr);
94   ierr = PetscViewerASCIIPrintf(viewer,"  RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr);
95   ierr = PetscViewerASCIIPrintf(viewer,"  ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES" : "NO");CHKERRQ(ierr);
96   ierr = PetscViewerASCIIPrintf(viewer,"  PrintStat: %s\n",(options.PrintStat != NO) ? "YES" : "NO");CHKERRQ(ierr);
97   ierr = PetscViewerASCIIPrintf(viewer,"  lwork: %D\n",lu->lwork);CHKERRQ(ierr);
98   if (A->factortype == MAT_FACTOR_ILU) {
99     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_DropTol: %g\n",options.ILU_DropTol);CHKERRQ(ierr);
100     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_FillTol: %g\n",options.ILU_FillTol);CHKERRQ(ierr);
101     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_FillFactor: %g\n",options.ILU_FillFactor);CHKERRQ(ierr);
102     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_DropRule: %D\n",options.ILU_DropRule);CHKERRQ(ierr);
103     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_Norm: %D\n",options.ILU_Norm);CHKERRQ(ierr);
104     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_MILU: %D\n",options.ILU_MILU);CHKERRQ(ierr);
105   }
106   PetscFunctionReturn(0);
107 }
108 
109 /*
110     These are the methods provided to REPLACE the corresponding methods of the
111    SeqAIJ matrix class. Other methods of SeqAIJ are not replaced
112 */
113 #undef __FUNCT__
114 #define __FUNCT__ "MatLUFactorNumeric_SuperLU"
115 PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info)
116 {
117   Mat_SuperLU    *lu = (Mat_SuperLU*)F->spptr;
118   Mat_SeqAIJ     *aa;
119   PetscErrorCode ierr;
120   PetscInt       sinfo;
121   PetscReal      ferr, berr;
122   NCformat       *Ustore;
123   SCformat       *Lstore;
124 
125   PetscFunctionBegin;
126   if (lu->flg == SAME_NONZERO_PATTERN) { /* successing numerical factorization */
127     lu->options.Fact = SamePattern;
128     /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
129     Destroy_SuperMatrix_Store(&lu->A);
130     if (lu->options.Equil) {
131       ierr = MatCopy_SeqAIJ(A,lu->A_dup,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
132     }
133     if (lu->lwork >= 0) {
134       PetscStackCall("SuperLU:Destroy_SuperNode_Matrix",Destroy_SuperNode_Matrix(&lu->L));
135       PetscStackCall("SuperLU:Destroy_CompCol_Matrix",Destroy_CompCol_Matrix(&lu->U));
136       lu->options.Fact = SamePattern;
137     }
138   }
139 
140   /* Create the SuperMatrix for lu->A=A^T:
141        Since SuperLU likes column-oriented matrices,we pass it the transpose,
142        and then solve A^T X = B in MatSolve(). */
143   if (lu->options.Equil) {
144     aa = (Mat_SeqAIJ*)(lu->A_dup)->data;
145   } else {
146     aa = (Mat_SeqAIJ*)(A)->data;
147   }
148 #if defined(PETSC_USE_COMPLEX)
149 #if defined(PETSC_USE_REAL_SINGLE)
150   PetscStackCall("SuperLU:cCreate_CompCol_Matrix",cCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(singlecomplex*)aa->a,aa->j,aa->i,SLU_NC,SLU_C,SLU_GE));
151 #else
152   PetscStackCall("SuperLU:zCreate_CompCol_Matrix",zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,SLU_NC,SLU_Z,SLU_GE));
153 #endif
154 #else
155 #if defined(PETSC_USE_REAL_SINGLE)
156   PetscStackCall("SuperLU:sCreate_CompCol_Matrix",sCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,SLU_NC,SLU_S,SLU_GE));
157 #else
158   PetscStackCall("SuperLU:dCreate_CompCol_Matrix",dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,SLU_NC,SLU_D,SLU_GE));
159 #endif
160 #endif
161 
162   /* Numerical factorization */
163   lu->B.ncol = 0;  /* Indicate not to solve the system */
164   if (F->factortype == MAT_FACTOR_LU) {
165 #if defined(PETSC_USE_COMPLEX)
166 #if defined(PETSC_USE_REAL_SINGLE)
167     PetscStackCall("SuperLU:cgssvx",cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
168                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
169                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
170 #else
171     PetscStackCall("SuperLU:zgssvx",zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
172                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
173                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
174 #endif
175 #else
176 #if defined(PETSC_USE_REAL_SINGLE)
177     PetscStackCall("SuperLU:sgssvx",sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
178                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
179                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
180 #else
181     PetscStackCall("SuperLU:dgssvx",dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
182                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
183                                      &lu->Glu,&lu->mem_usage, &lu->stat, &sinfo));
184 #endif
185 #endif
186   } else if (F->factortype == MAT_FACTOR_ILU) {
187     /* Compute the incomplete factorization, condition number and pivot growth */
188 #if defined(PETSC_USE_COMPLEX)
189 #if defined(PETSC_USE_REAL_SINGLE)
190     PetscStackCall("SuperLU:cgsisx",cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C,
191                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
192                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
193 #else
194     PetscStackCall("SuperLU:zgsisx",zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C,
195                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
196                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
197 #endif
198 #else
199 #if defined(PETSC_USE_REAL_SINGLE)
200     PetscStackCall("SuperLU:sgsisx",sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
201                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
202                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
203 #else
204     PetscStackCall("SuperLU:dgsisx",dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
205                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
206                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
207 #endif
208 #endif
209   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
210   if (!sinfo || sinfo == lu->A.ncol+1) {
211     if (lu->options.PivotGrowth) {
212       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. pivot growth = %e\n", lu->rpg);CHKERRQ(ierr);
213     }
214     if (lu->options.ConditionNumber) {
215       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. condition number = %e\n", lu->rcond);CHKERRQ(ierr);
216     }
217   } else if (sinfo > 0) {
218     if (A->erroriffailure) {
219       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo);
220     } else {
221       if (sinfo <= lu->A.ncol) {
222         if (lu->options.ILU_FillTol == 0.0) {
223           F->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
224         }
225         ierr = PetscInfo2(F,"Number of zero pivots %D, ILU_FillTol %g\n",sinfo,lu->options.ILU_FillTol);CHKERRQ(ierr);
226       } else if (sinfo == lu->A.ncol + 1) {
227         /*
228          U is nonsingular, but RCOND is less than machine
229  		      precision, meaning that the matrix is singular to
230  		      working precision. Nevertheless, the solution and
231  		      error bounds are computed because there are a number
232  		      of situations where the computed solution can be more
233  		      accurate than the value of RCOND would suggest.
234          */
235         ierr = PetscInfo1(F,"Matrix factor U is nonsingular, but is singular to working precision. The solution is computed. info %D",sinfo);CHKERRQ(ierr);
236       } else { /* sinfo > lu->A.ncol + 1 */
237         F->errortype = MAT_FACTOR_OUTMEMORY;
238         ierr = PetscInfo1(F,"Number of bytes allocated when memory allocation fails %D\n",sinfo);CHKERRQ(ierr);
239       }
240     }
241   } else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);
242 
243   if (lu->options.PrintStat) {
244     ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");CHKERRQ(ierr);
245     PetscStackCall("SuperLU:StatPrint",StatPrint(&lu->stat));
246     Lstore = (SCformat*) lu->L.Store;
247     Ustore = (NCformat*) lu->U.Store;
248     ierr   = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %D\n", Lstore->nnz);
249     ierr   = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %D\n", Ustore->nnz);
250     ierr   = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
251     ierr   = PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\n",
252                          lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6);
253   }
254 
255   lu->flg                = SAME_NONZERO_PATTERN;
256   F->ops->solve          = MatSolve_SuperLU;
257   F->ops->solvetranspose = MatSolveTranspose_SuperLU;
258   F->ops->matsolve       = NULL;
259   PetscFunctionReturn(0);
260 }
261 
262 #undef __FUNCT__
263 #define __FUNCT__ "MatGetDiagonal_SuperLU"
264 PetscErrorCode MatGetDiagonal_SuperLU(Mat A,Vec v)
265 {
266   PetscFunctionBegin;
267   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: SuperLU factor");
268   PetscFunctionReturn(0);
269 }
270 
271 #undef __FUNCT__
272 #define __FUNCT__ "MatDestroy_SuperLU"
273 PetscErrorCode MatDestroy_SuperLU(Mat A)
274 {
275   PetscErrorCode ierr;
276   Mat_SuperLU    *lu=(Mat_SuperLU*)A->spptr;
277 
278   PetscFunctionBegin;
279   if (lu && lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
280     PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->A));
281     PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->B));
282     PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->X));
283     PetscStackCall("SuperLU:StatFree",StatFree(&lu->stat));
284     if (lu->lwork >= 0) {
285       PetscStackCall("SuperLU:Destroy_SuperNode_Matrix",Destroy_SuperNode_Matrix(&lu->L));
286       PetscStackCall("SuperLU:Destroy_CompCol_Matrix",Destroy_CompCol_Matrix(&lu->U));
287     }
288   }
289   if (lu) {
290     ierr = PetscFree(lu->etree);CHKERRQ(ierr);
291     ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
292     ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
293     ierr = PetscFree(lu->R);CHKERRQ(ierr);
294     ierr = PetscFree(lu->C);CHKERRQ(ierr);
295     ierr = PetscFree(lu->rhs_dup);CHKERRQ(ierr);
296     ierr = MatDestroy(&lu->A_dup);CHKERRQ(ierr);
297   }
298   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
299 
300   /* clear composed functions */
301   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr);
302   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSuperluSetILUDropTol_C",NULL);CHKERRQ(ierr);
303 
304   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
305   PetscFunctionReturn(0);
306 }
307 
308 #undef __FUNCT__
309 #define __FUNCT__ "MatView_SuperLU"
310 PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
311 {
312   PetscErrorCode    ierr;
313   PetscBool         iascii;
314   PetscViewerFormat format;
315 
316   PetscFunctionBegin;
317   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
318   if (iascii) {
319     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
320     if (format == PETSC_VIEWER_ASCII_INFO) {
321       ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr);
322     }
323   }
324   PetscFunctionReturn(0);
325 }
326 
327 
328 #undef __FUNCT__
329 #define __FUNCT__ "MatSolve_SuperLU_Private"
330 PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x)
331 {
332   Mat_SuperLU       *lu = (Mat_SuperLU*)A->spptr;
333   const PetscScalar *barray;
334   PetscScalar       *xarray;
335   PetscErrorCode    ierr;
336   PetscInt          info,i,n;
337   PetscReal         ferr,berr;
338   static PetscBool  cite = PETSC_FALSE;
339 
340   PetscFunctionBegin;
341   if (lu->lwork == -1) PetscFunctionReturn(0);
342   ierr = PetscCitationsRegister("@article{superlu99,\n  author  = {James W. Demmel and Stanley C. Eisenstat and\n             John R. Gilbert and Xiaoye S. Li and Joseph W. H. Liu},\n  title = {A supernodal approach to sparse partial pivoting},\n  journal = {SIAM J. Matrix Analysis and Applications},\n  year = {1999},\n  volume  = {20},\n  number = {3},\n  pages = {720-755}\n}\n",&cite);CHKERRQ(ierr);
343 
344   ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr);
345   lu->B.ncol = 1;   /* Set the number of right-hand side */
346   if (lu->options.Equil && !lu->rhs_dup) {
347     /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */
348     ierr = PetscMalloc1(n,&lu->rhs_dup);CHKERRQ(ierr);
349   }
350   if (lu->options.Equil) {
351     /* Copy b into rsh_dup */
352     ierr   = VecGetArrayRead(b,&barray);CHKERRQ(ierr);
353     ierr   = PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));CHKERRQ(ierr);
354     ierr   = VecRestoreArrayRead(b,&barray);CHKERRQ(ierr);
355     barray = lu->rhs_dup;
356   } else {
357     ierr = VecGetArrayRead(b,&barray);CHKERRQ(ierr);
358   }
359   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
360 
361 #if defined(PETSC_USE_COMPLEX)
362 #if defined(PETSC_USE_REAL_SINGLE)
363   ((DNformat*)lu->B.Store)->nzval = (singlecomplex*)barray;
364   ((DNformat*)lu->X.Store)->nzval = (singlecomplex*)xarray;
365 #else
366   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
367   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
368 #endif
369 #else
370   ((DNformat*)lu->B.Store)->nzval = (void*)barray;
371   ((DNformat*)lu->X.Store)->nzval = xarray;
372 #endif
373 
374   lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
375   if (A->factortype == MAT_FACTOR_LU) {
376 #if defined(PETSC_USE_COMPLEX)
377 #if defined(PETSC_USE_REAL_SINGLE)
378     PetscStackCall("SuperLU:cgssvx",cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
379                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
380                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
381 #else
382     PetscStackCall("SuperLU:zgssvx",zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
383                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
384                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
385 #endif
386 #else
387 #if defined(PETSC_USE_REAL_SINGLE)
388     PetscStackCall("SuperLU:sgssvx",sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
389                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
390                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
391 #else
392     PetscStackCall("SuperLU:dgssvx",dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
393                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
394                                      &lu->Glu,&lu->mem_usage, &lu->stat, &info));
395 #endif
396 #endif
397   } else if (A->factortype == MAT_FACTOR_ILU) {
398 #if defined(PETSC_USE_COMPLEX)
399 #if defined(PETSC_USE_REAL_SINGLE)
400     PetscStackCall("SuperLU:cgsisx",cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
401                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
402                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
403 #else
404     PetscStackCall("SuperLU:zgsisx",zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
405                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
406                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
407 #endif
408 #else
409 #if defined(PETSC_USE_REAL_SINGLE)
410     PetscStackCall("SuperLU:sgsisx",sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
411                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
412                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
413 #else
414     PetscStackCall("SuperLU:dgsisx",dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
415                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
416                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
417 #endif
418 #endif
419   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
420   if (!lu->options.Equil) {
421     ierr = VecRestoreArrayRead(b,&barray);CHKERRQ(ierr);
422   }
423   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
424 
425   if (!info || info == lu->A.ncol+1) {
426     if (lu->options.IterRefine) {
427       ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
428       ierr = PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
429       for (i = 0; i < 1; ++i) {
430         ierr = PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr);
431       }
432     }
433   } else if (info > 0) {
434     if (lu->lwork == -1) {
435       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", info - lu->A.ncol);
436     } else {
437       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",info);
438     }
439   } else if (info < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info);
440 
441   if (lu->options.PrintStat) {
442     ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
443     PetscStackCall("SuperLU:StatPrint",StatPrint(&lu->stat));
444   }
445   PetscFunctionReturn(0);
446 }
447 
448 #undef __FUNCT__
449 #define __FUNCT__ "MatSolve_SuperLU"
450 PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
451 {
452   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
453   PetscErrorCode ierr;
454 
455   PetscFunctionBegin;
456   if (A->errortype) {
457     ierr = PetscInfo(A,"MatSolve is called with singular matrix factor, skip\n");CHKERRQ(ierr);
458     ierr = VecSetInf(x);CHKERRQ(ierr);
459     PetscFunctionReturn(0);
460   }
461 
462   lu->options.Trans = TRANS;
463   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
464   PetscFunctionReturn(0);
465 }
466 
467 #undef __FUNCT__
468 #define __FUNCT__ "MatSolveTranspose_SuperLU"
469 PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x)
470 {
471   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
472   PetscErrorCode ierr;
473 
474   PetscFunctionBegin;
475   if (A->errortype) {
476     ierr = PetscInfo(A,"MatSolve is called with singular matrix factor, skip\n");CHKERRQ(ierr);
477     ierr = VecSetInf(x);CHKERRQ(ierr);
478     PetscFunctionReturn(0);
479   }
480 
481   lu->options.Trans = NOTRANS;
482   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
483   PetscFunctionReturn(0);
484 }
485 
486 #undef __FUNCT__
487 #define __FUNCT__ "MatMatSolve_SuperLU"
488 PetscErrorCode MatMatSolve_SuperLU(Mat A,Mat B,Mat X)
489 {
490   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
491   PetscBool      flg;
492   PetscErrorCode ierr;
493 
494   PetscFunctionBegin;
495   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
496   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
497   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
498   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
499   lu->options.Trans = TRANS;
500   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_SuperLU() is not implemented yet");
501   PetscFunctionReturn(0);
502 }
503 
504 /*
505    Note the r permutation is ignored
506 */
507 #undef __FUNCT__
508 #define __FUNCT__ "MatLUFactorSymbolic_SuperLU"
509 PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
510 {
511   Mat_SuperLU *lu = (Mat_SuperLU*)(F->spptr);
512 
513   PetscFunctionBegin;
514   lu->flg                 = DIFFERENT_NONZERO_PATTERN;
515   lu->CleanUpSuperLU      = PETSC_TRUE;
516   F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
517   PetscFunctionReturn(0);
518 }
519 
520 #undef __FUNCT__
521 #define __FUNCT__ "MatSuperluSetILUDropTol_SuperLU"
522 static PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol)
523 {
524   Mat_SuperLU *lu= (Mat_SuperLU*)F->spptr;
525 
526   PetscFunctionBegin;
527   lu->options.ILU_DropTol = dtol;
528   PetscFunctionReturn(0);
529 }
530 
531 #undef __FUNCT__
532 #define __FUNCT__ "MatSuperluSetILUDropTol"
533 /*@
534   MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance
535    Logically Collective on Mat
536 
537    Input Parameters:
538 +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface
539 -  dtol - drop tolerance
540 
541   Options Database:
542 .   -mat_superlu_ilu_droptol <dtol>
543 
544    Level: beginner
545 
546    References:
547 .      SuperLU Users' Guide
548 
549 .seealso: MatGetFactor()
550 @*/
551 PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol)
552 {
553   PetscErrorCode ierr;
554 
555   PetscFunctionBegin;
556   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
557   PetscValidLogicalCollectiveReal(F,dtol,2);
558   ierr = PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol));CHKERRQ(ierr);
559   PetscFunctionReturn(0);
560 }
561 
562 #undef __FUNCT__
563 #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_superlu"
564 PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type)
565 {
566   PetscFunctionBegin;
567   *type = MATSOLVERSUPERLU;
568   PetscFunctionReturn(0);
569 }
570 
571 
572 /*MC
573   MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices
574   via the external package SuperLU.
575 
576   Use ./configure --download-superlu to have PETSc installed with SuperLU
577 
578   Use -pc_type lu -pc_factor_mat_solver_package superlu to us this direct solver
579 
580   Options Database Keys:
581 + -mat_superlu_equil <FALSE>            - Equil (None)
582 . -mat_superlu_colperm <COLAMD>         - (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD
583 . -mat_superlu_iterrefine <NOREFINE>    - (choose one of) NOREFINE SINGLE DOUBLE EXTRA
584 . -mat_superlu_symmetricmode: <FALSE>   - SymmetricMode (None)
585 . -mat_superlu_diagpivotthresh <1>      - DiagPivotThresh (None)
586 . -mat_superlu_pivotgrowth <FALSE>      - PivotGrowth (None)
587 . -mat_superlu_conditionnumber <FALSE>  - ConditionNumber (None)
588 . -mat_superlu_rowperm <NOROWPERM>      - (choose one of) NOROWPERM LargeDiag
589 . -mat_superlu_replacetinypivot <FALSE> - ReplaceTinyPivot (None)
590 . -mat_superlu_printstat <FALSE>        - PrintStat (None)
591 . -mat_superlu_lwork <0>                - size of work array in bytes used by factorization (None)
592 . -mat_superlu_ilu_droptol <0>          - ILU_DropTol (None)
593 . -mat_superlu_ilu_filltol <0>          - ILU_FillTol (None)
594 . -mat_superlu_ilu_fillfactor <0>       - ILU_FillFactor (None)
595 . -mat_superlu_ilu_droprull <0>         - ILU_DropRule (None)
596 . -mat_superlu_ilu_norm <0>             - ILU_Norm (None)
597 - -mat_superlu_ilu_milu <0>             - ILU_MILU (None)
598 
599    Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves
600 
601    Level: beginner
602 
603 .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage
604 M*/
605 
606 #undef __FUNCT__
607 #define __FUNCT__ "MatGetFactor_seqaij_superlu"
608 static PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F)
609 {
610   Mat            B;
611   Mat_SuperLU    *lu;
612   PetscErrorCode ierr;
613   PetscInt       indx,m=A->rmap->n,n=A->cmap->n;
614   PetscBool      flg,set;
615   PetscReal      real_input;
616   const char     *colperm[]   ={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
617   const char     *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
618   const char     *rowperm[]   ={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
619 
620   PetscFunctionBegin;
621   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
622   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
623   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
624   ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
625 
626   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
627     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_SuperLU;
628     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU;
629   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
630 
631   B->ops->destroy     = MatDestroy_SuperLU;
632   B->ops->view        = MatView_SuperLU;
633   B->ops->getdiagonal = MatGetDiagonal_SuperLU;
634   B->factortype       = ftype;
635   B->assembled        = PETSC_TRUE;           /* required by -ksp_view */
636   B->preallocated     = PETSC_TRUE;
637 
638   ierr = PetscNewLog(B,&lu);CHKERRQ(ierr);
639 
640   if (ftype == MAT_FACTOR_LU) {
641     set_default_options(&lu->options);
642     /* Comments from SuperLU_4.0/SRC/dgssvx.c:
643       "Whether or not the system will be equilibrated depends on the
644        scaling of the matrix A, but if equilibration is used, A is
645        overwritten by diag(R)*A*diag(C) and B by diag(R)*B
646        (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)."
647      We set 'options.Equil = NO' as default because additional space is needed for it.
648     */
649     lu->options.Equil = NO;
650   } else if (ftype == MAT_FACTOR_ILU) {
651     /* Set the default input options of ilu: */
652     PetscStackCall("SuperLU:ilu_set_default_options",ilu_set_default_options(&lu->options));
653   }
654   lu->options.PrintStat = NO;
655 
656   /* Initialize the statistics variables. */
657   PetscStackCall("SuperLU:StatInit",StatInit(&lu->stat));
658   lu->lwork = 0;   /* allocate space internally by system malloc */
659 
660   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr);
661   ierr = PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,NULL);CHKERRQ(ierr);
662   ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr);
663   if (flg) lu->options.ColPerm = (colperm_t)indx;
664   ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr);
665   if (flg) lu->options.IterRefine = (IterRefine_t)indx;
666   ierr = PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,&set);CHKERRQ(ierr);
667   if (set && flg) lu->options.SymmetricMode = YES;
668   ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&real_input,&flg);CHKERRQ(ierr);
669   if (flg) lu->options.DiagPivotThresh = (double) real_input;
670   ierr = PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,&set);CHKERRQ(ierr);
671   if (set && flg) lu->options.PivotGrowth = YES;
672   ierr = PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,&set);CHKERRQ(ierr);
673   if (set && flg) lu->options.ConditionNumber = YES;
674   ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[lu->options.RowPerm],&indx,&flg);CHKERRQ(ierr);
675   if (flg) lu->options.RowPerm = (rowperm_t)indx;
676   ierr = PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,&set);CHKERRQ(ierr);
677   if (set && flg) lu->options.ReplaceTinyPivot = YES;
678   ierr = PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,&set);CHKERRQ(ierr);
679   if (set && flg) lu->options.PrintStat = YES;
680   ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,NULL);CHKERRQ(ierr);
681   if (lu->lwork > 0) {
682     /* lwork is in bytes, hence PetscMalloc() is used here, not PetscMalloc1()*/
683     ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr);
684   } else if (lu->lwork != 0 && lu->lwork != -1) {
685     ierr      = PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
686     lu->lwork = 0;
687   }
688   /* ilu options */
689   ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&real_input,&flg);CHKERRQ(ierr);
690   if (flg) lu->options.ILU_DropTol = (double) real_input;
691   ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&real_input,&flg);CHKERRQ(ierr);
692   if (flg) lu->options.ILU_FillTol = (double) real_input;
693   ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&real_input,&flg);CHKERRQ(ierr);
694   if (flg) lu->options.ILU_FillFactor = (double) real_input;
695   ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,NULL);CHKERRQ(ierr);
696   ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr);
697   if (flg) lu->options.ILU_Norm = (norm_t)indx;
698   ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr);
699   if (flg) lu->options.ILU_MILU = (milu_t)indx;
700   PetscOptionsEnd();
701   if (lu->options.Equil == YES) {
702     /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */
703     ierr = MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);CHKERRQ(ierr);
704   }
705 
706   /* Allocate spaces (notice sizes are for the transpose) */
707   ierr = PetscMalloc1(m,&lu->etree);CHKERRQ(ierr);
708   ierr = PetscMalloc1(n,&lu->perm_r);CHKERRQ(ierr);
709   ierr = PetscMalloc1(m,&lu->perm_c);CHKERRQ(ierr);
710   ierr = PetscMalloc1(n,&lu->R);CHKERRQ(ierr);
711   ierr = PetscMalloc1(m,&lu->C);CHKERRQ(ierr);
712 
713   /* create rhs and solution x without allocate space for .Store */
714 #if defined(PETSC_USE_COMPLEX)
715 #if defined(PETSC_USE_REAL_SINGLE)
716   PetscStackCall("SuperLU:cCreate_Dense_Matrix(",cCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE));
717   PetscStackCall("SuperLU:cCreate_Dense_Matrix(",cCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE));
718 #else
719   PetscStackCall("SuperLU:zCreate_Dense_Matrix",zCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE));
720   PetscStackCall("SuperLU:zCreate_Dense_Matrix",zCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE));
721 #endif
722 #else
723 #if defined(PETSC_USE_REAL_SINGLE)
724   PetscStackCall("SuperLU:sCreate_Dense_Matrix",sCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE));
725   PetscStackCall("SuperLU:sCreate_Dense_Matrix",sCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE));
726 #else
727   PetscStackCall("SuperLU:dCreate_Dense_Matrix",dCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE));
728   PetscStackCall("SuperLU:dCreate_Dense_Matrix",dCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE));
729 #endif
730 #endif
731 
732   ierr     = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr);
733   ierr     = PetscObjectComposeFunction((PetscObject)B,"MatSuperluSetILUDropTol_C",MatSuperluSetILUDropTol_SuperLU);CHKERRQ(ierr);
734   B->spptr = lu;
735 
736   *F       = B;
737   PetscFunctionReturn(0);
738 }
739 
740 #undef __FUNCT__
741 #define __FUNCT__ "MatSolverPackageRegister_SuperLU"
742 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_SuperLU(void)
743 {
744   PetscErrorCode ierr;
745 
746   PetscFunctionBegin;
747   ierr = MatSolverPackageRegister(MATSOLVERSUPERLU,MATSEQAIJ,       MAT_FACTOR_LU,MatGetFactor_seqaij_superlu);CHKERRQ(ierr);
748   ierr = MatSolverPackageRegister(MATSOLVERSUPERLU,MATSEQAIJ,       MAT_FACTOR_ILU,MatGetFactor_seqaij_superlu);CHKERRQ(ierr);
749   PetscFunctionReturn(0);
750 }
751