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