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