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