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