xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
114b396bbSKris Buschelman 
22ef1f0ffSBarry Smith /*
35a46d3faSBarry Smith      This file implements a subclass of the SeqAIJ matrix class that uses
4c2b89b5dSBarry Smith      the SuperLU sparse solver.
514b396bbSKris Buschelman */
614b396bbSKris Buschelman 
75a46d3faSBarry Smith /*
85a46d3faSBarry Smith      Defines the data structure for the base matrix type (SeqAIJ)
95a46d3faSBarry Smith */
10c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
1114b396bbSKris Buschelman 
125a46d3faSBarry Smith /*
135a46d3faSBarry Smith      SuperLU include files
145a46d3faSBarry Smith */
1514b396bbSKris Buschelman EXTERN_C_BEGIN
1614b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
173cb270beSHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
183cb270beSHong Zhang     #include <slu_cdefs.h>
193cb270beSHong Zhang   #else
20c6db04a5SJed Brown     #include <slu_zdefs.h>
213cb270beSHong Zhang   #endif
223cb270beSHong Zhang #else
233cb270beSHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
243cb270beSHong Zhang     #include <slu_sdefs.h>
2514b396bbSKris Buschelman   #else
26c6db04a5SJed Brown     #include <slu_ddefs.h>
2714b396bbSKris Buschelman   #endif
283cb270beSHong Zhang #endif
29c6db04a5SJed Brown #include <slu_util.h>
3014b396bbSKris Buschelman EXTERN_C_END
3114b396bbSKris Buschelman 
325a46d3faSBarry Smith /*
33245fece9SBarry Smith      This is the data that defines the SuperLU factored matrix type
345a46d3faSBarry Smith */
3514b396bbSKris Buschelman typedef struct {
3675af56d4SHong Zhang   SuperMatrix       A, L, U, B, X;
3775af56d4SHong Zhang   superlu_options_t options;
38da7a1d00SHong Zhang   PetscInt         *perm_c; /* column permutation vector */
39da7a1d00SHong Zhang   PetscInt         *perm_r; /* row permutations from partial pivoting */
40da7a1d00SHong Zhang   PetscInt         *etree;
41da7a1d00SHong Zhang   PetscReal        *R, *C;
4275af56d4SHong Zhang   char              equed[1];
43da7a1d00SHong Zhang   PetscInt          lwork;
4475af56d4SHong Zhang   void             *work;
45da7a1d00SHong Zhang   PetscReal         rpg, rcond;
4675af56d4SHong Zhang   mem_usage_t       mem_usage;
4775af56d4SHong Zhang   MatStructure      flg;
485d8b2efaSHong Zhang   SuperLUStat_t     stat;
49cae5a23dSHong Zhang   Mat               A_dup;
50cae5a23dSHong Zhang   PetscScalar      *rhs_dup;
51c147c726SHong Zhang   GlobalLU_t        Glu;
522366e350SStefano Zampini   PetscBool         needconversion;
5314b396bbSKris Buschelman 
5414b396bbSKris Buschelman   /* Flag to clean up (non-global) SuperLU objects during Destroy */
55ace3abfcSBarry Smith   PetscBool CleanUpSuperLU;
56f0c56d0fSKris Buschelman } Mat_SuperLU;
5714b396bbSKris Buschelman 
585a46d3faSBarry Smith /*
595a46d3faSBarry Smith     Utility function
605a46d3faSBarry Smith */
61d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Info_SuperLU(Mat A, PetscViewer viewer)
62d71ae5a4SJacob Faibussowitsch {
63245fece9SBarry Smith   Mat_SuperLU      *lu = (Mat_SuperLU *)A->data;
645a46d3faSBarry Smith   superlu_options_t options;
655a46d3faSBarry Smith 
665a46d3faSBarry Smith   PetscFunctionBegin;
675a46d3faSBarry Smith   options = lu->options;
682205254eSKarl Rupp 
699566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "SuperLU run parameters:\n"));
709566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Equil: %s\n", (options.Equil != NO) ? "YES" : "NO"));
719566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  ColPerm: %" PetscInt_FMT "\n", options.ColPerm));
729566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  IterRefine: %" PetscInt_FMT "\n", options.IterRefine));
739566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  SymmetricMode: %s\n", (options.SymmetricMode != NO) ? "YES" : "NO"));
749566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  DiagPivotThresh: %g\n", options.DiagPivotThresh));
759566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  PivotGrowth: %s\n", (options.PivotGrowth != NO) ? "YES" : "NO"));
769566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  ConditionNumber: %s\n", (options.ConditionNumber != NO) ? "YES" : "NO"));
779566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  RowPerm: %" PetscInt_FMT "\n", options.RowPerm));
789566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  ReplaceTinyPivot: %s\n", (options.ReplaceTinyPivot != NO) ? "YES" : "NO"));
799566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  PrintStat: %s\n", (options.PrintStat != NO) ? "YES" : "NO"));
809566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  lwork: %" PetscInt_FMT "\n", lu->lwork));
81d5f3da31SBarry Smith   if (A->factortype == MAT_FACTOR_ILU) {
829566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  ILU_DropTol: %g\n", options.ILU_DropTol));
839566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  ILU_FillTol: %g\n", options.ILU_FillTol));
849566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  ILU_FillFactor: %g\n", options.ILU_FillFactor));
859566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  ILU_DropRule: %" PetscInt_FMT "\n", options.ILU_DropRule));
869566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  ILU_Norm: %" PetscInt_FMT "\n", options.ILU_Norm));
879566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  ILU_MILU: %" PetscInt_FMT "\n", options.ILU_MILU));
88cffbb591SHong Zhang   }
893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
905a46d3faSBarry Smith }
915a46d3faSBarry Smith 
92d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SuperLU_Private(Mat A, Vec b, Vec x)
93d71ae5a4SJacob Faibussowitsch {
94245fece9SBarry Smith   Mat_SuperLU       *lu = (Mat_SuperLU *)A->data;
95245fece9SBarry Smith   const PetscScalar *barray;
96245fece9SBarry Smith   PetscScalar       *xarray;
97245fece9SBarry Smith   PetscInt           info, i, n;
98245fece9SBarry Smith   PetscReal          ferr, berr;
99245fece9SBarry Smith   static PetscBool   cite = PETSC_FALSE;
100245fece9SBarry Smith 
101245fece9SBarry Smith   PetscFunctionBegin;
1023ba16761SJacob Faibussowitsch   if (lu->lwork == -1) PetscFunctionReturn(PETSC_SUCCESS);
1039371c9d4SSatish Balay   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 "
1049371c9d4SSatish Balay                                    "pivoting},\n  journal = {SIAM J. Matrix Analysis and Applications},\n  year = {1999},\n  volume  = {20},\n  number = {3},\n  pages = {720-755}\n}\n",
1059371c9d4SSatish Balay                                    &cite));
106245fece9SBarry Smith 
1079566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x, &n));
108245fece9SBarry Smith   lu->B.ncol = 1; /* Set the number of right-hand side */
109245fece9SBarry Smith   if (lu->options.Equil && !lu->rhs_dup) {
110245fece9SBarry Smith     /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */
1119566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &lu->rhs_dup));
112245fece9SBarry Smith   }
113245fece9SBarry Smith   if (lu->options.Equil) {
114245fece9SBarry Smith     /* Copy b into rsh_dup */
1159566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(b, &barray));
1169566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(lu->rhs_dup, barray, n));
1179566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(b, &barray));
118245fece9SBarry Smith     barray = lu->rhs_dup;
119245fece9SBarry Smith   } else {
1209566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(b, &barray));
121245fece9SBarry Smith   }
1229566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &xarray));
123245fece9SBarry Smith 
124245fece9SBarry Smith #if defined(PETSC_USE_COMPLEX)
125245fece9SBarry Smith   #if defined(PETSC_USE_REAL_SINGLE)
126245fece9SBarry Smith   ((DNformat *)lu->B.Store)->nzval = (singlecomplex *)barray;
127245fece9SBarry Smith   ((DNformat *)lu->X.Store)->nzval = (singlecomplex *)xarray;
128245fece9SBarry Smith   #else
129245fece9SBarry Smith   ((DNformat *)lu->B.Store)->nzval = (doublecomplex *)barray;
130245fece9SBarry Smith   ((DNformat *)lu->X.Store)->nzval = (doublecomplex *)xarray;
131245fece9SBarry Smith   #endif
132245fece9SBarry Smith #else
133245fece9SBarry Smith   ((DNformat *)lu->B.Store)->nzval = (void *)barray;
134245fece9SBarry Smith   ((DNformat *)lu->X.Store)->nzval = xarray;
135245fece9SBarry Smith #endif
136245fece9SBarry Smith 
137245fece9SBarry Smith   lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
138245fece9SBarry Smith   if (A->factortype == MAT_FACTOR_LU) {
139245fece9SBarry Smith #if defined(PETSC_USE_COMPLEX)
140245fece9SBarry Smith   #if defined(PETSC_USE_REAL_SINGLE)
1419371c9d4SSatish Balay     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));
142245fece9SBarry Smith   #else
1439371c9d4SSatish Balay     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));
144245fece9SBarry Smith   #endif
145245fece9SBarry Smith #else
146245fece9SBarry Smith   #if defined(PETSC_USE_REAL_SINGLE)
1479371c9d4SSatish Balay     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));
148245fece9SBarry Smith   #else
1499371c9d4SSatish Balay     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));
150245fece9SBarry Smith   #endif
151245fece9SBarry Smith #endif
152245fece9SBarry Smith   } else if (A->factortype == MAT_FACTOR_ILU) {
153245fece9SBarry Smith #if defined(PETSC_USE_COMPLEX)
154245fece9SBarry Smith   #if defined(PETSC_USE_REAL_SINGLE)
1559371c9d4SSatish Balay     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));
156245fece9SBarry Smith   #else
1579371c9d4SSatish Balay     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));
158245fece9SBarry Smith   #endif
159245fece9SBarry Smith #else
160245fece9SBarry Smith   #if defined(PETSC_USE_REAL_SINGLE)
1619371c9d4SSatish Balay     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));
162245fece9SBarry Smith   #else
1639371c9d4SSatish Balay     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));
164245fece9SBarry Smith   #endif
165245fece9SBarry Smith #endif
166245fece9SBarry Smith   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported");
16748a46eb9SPierre Jolivet   if (!lu->options.Equil) PetscCall(VecRestoreArrayRead(b, &barray));
1689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &xarray));
169245fece9SBarry Smith 
170245fece9SBarry Smith   if (!info || info == lu->A.ncol + 1) {
171245fece9SBarry Smith     if (lu->options.IterRefine) {
1729566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Iterative Refinement:\n"));
1739566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"));
1749566063dSJacob Faibussowitsch       for (i = 0; i < 1; ++i) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  %8d%8d%16e%16e\n", i + 1, lu->stat.RefineSteps, ferr, berr));
175245fece9SBarry Smith     }
176245fece9SBarry Smith   } else if (info > 0) {
177245fece9SBarry Smith     if (lu->lwork == -1) {
1789566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "  ** Estimated memory: %" PetscInt_FMT " bytes\n", info - lu->A.ncol));
179245fece9SBarry Smith     } else {
1809566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Warning: gssvx() returns info %" PetscInt_FMT "\n", info));
181245fece9SBarry Smith     }
1825f80ce2aSJacob Faibussowitsch   } 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);
183245fece9SBarry Smith 
184245fece9SBarry Smith   if (lu->options.PrintStat) {
1859566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "MatSolve__SuperLU():\n"));
186e77caa6dSBarry Smith     PetscStackCallExternalVoid("SuperLU:StatPrint", StatPrint(&lu->stat));
187245fece9SBarry Smith   }
1883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
189245fece9SBarry Smith }
190245fece9SBarry Smith 
191d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SuperLU(Mat A, Vec b, Vec x)
192d71ae5a4SJacob Faibussowitsch {
193245fece9SBarry Smith   Mat_SuperLU *lu = (Mat_SuperLU *)A->data;
19478bc9606SHong Zhang   trans_t      oldOption;
195245fece9SBarry Smith 
196245fece9SBarry Smith   PetscFunctionBegin;
197603e8f96SBarry Smith   if (A->factorerrortype) {
1989566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "MatSolve is called with singular matrix factor, skip\n"));
1999566063dSJacob Faibussowitsch     PetscCall(VecSetInf(x));
2003ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
201245fece9SBarry Smith   }
202245fece9SBarry Smith 
20378bc9606SHong Zhang   oldOption         = lu->options.Trans;
204245fece9SBarry Smith   lu->options.Trans = TRANS;
2059566063dSJacob Faibussowitsch   PetscCall(MatSolve_SuperLU_Private(A, b, x));
20678bc9606SHong Zhang   lu->options.Trans = oldOption;
2073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
208245fece9SBarry Smith }
209245fece9SBarry Smith 
210d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SuperLU(Mat A, Vec b, Vec x)
211d71ae5a4SJacob Faibussowitsch {
212245fece9SBarry Smith   Mat_SuperLU *lu = (Mat_SuperLU *)A->data;
21378bc9606SHong Zhang   trans_t      oldOption;
214245fece9SBarry Smith 
215245fece9SBarry Smith   PetscFunctionBegin;
216603e8f96SBarry Smith   if (A->factorerrortype) {
2179566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "MatSolve is called with singular matrix factor, skip\n"));
2189566063dSJacob Faibussowitsch     PetscCall(VecSetInf(x));
2193ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
220245fece9SBarry Smith   }
221245fece9SBarry Smith 
22278bc9606SHong Zhang   oldOption         = lu->options.Trans;
223245fece9SBarry Smith   lu->options.Trans = NOTRANS;
2249566063dSJacob Faibussowitsch   PetscCall(MatSolve_SuperLU_Private(A, b, x));
22578bc9606SHong Zhang   lu->options.Trans = oldOption;
2263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
227245fece9SBarry Smith }
228245fece9SBarry Smith 
229d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F, Mat A, const MatFactorInfo *info)
230d71ae5a4SJacob Faibussowitsch {
231245fece9SBarry Smith   Mat_SuperLU *lu = (Mat_SuperLU *)F->data;
232cae5a23dSHong Zhang   Mat_SeqAIJ  *aa;
2335a46d3faSBarry Smith   PetscInt     sinfo;
2345a46d3faSBarry Smith   PetscReal    ferr, berr;
2355a46d3faSBarry Smith   NCformat    *Ustore;
2365a46d3faSBarry Smith   SCformat    *Lstore;
2375a46d3faSBarry Smith 
2385a46d3faSBarry Smith   PetscFunctionBegin;
239da81f932SPierre Jolivet   if (lu->flg == SAME_NONZERO_PATTERN) { /* successive numerical factorization */
2405a46d3faSBarry Smith     lu->options.Fact = SamePattern;
2415a46d3faSBarry Smith     /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
2425a46d3faSBarry Smith     Destroy_SuperMatrix_Store(&lu->A);
2431baa6e33SBarry Smith     if (lu->A_dup) PetscCall(MatCopy_SeqAIJ(A, lu->A_dup, SAME_NONZERO_PATTERN));
2442366e350SStefano Zampini 
2455a46d3faSBarry Smith     if (lu->lwork >= 0) {
246e77caa6dSBarry Smith       PetscStackCallExternalVoid("SuperLU:Destroy_SuperNode_Matrix", Destroy_SuperNode_Matrix(&lu->L));
247e77caa6dSBarry Smith       PetscStackCallExternalVoid("SuperLU:Destroy_CompCol_Matrix", Destroy_CompCol_Matrix(&lu->U));
2485a46d3faSBarry Smith       lu->options.Fact = SamePattern;
2495a46d3faSBarry Smith     }
2505a46d3faSBarry Smith   }
2515a46d3faSBarry Smith 
2525a46d3faSBarry Smith   /* Create the SuperMatrix for lu->A=A^T:
2535a46d3faSBarry Smith        Since SuperLU likes column-oriented matrices,we pass it the transpose,
2545a46d3faSBarry Smith        and then solve A^T X = B in MatSolve(). */
2552366e350SStefano Zampini   if (lu->A_dup) {
256cae5a23dSHong Zhang     aa = (Mat_SeqAIJ *)(lu->A_dup)->data;
257cae5a23dSHong Zhang   } else {
258cae5a23dSHong Zhang     aa = (Mat_SeqAIJ *)(A)->data;
259cae5a23dSHong Zhang   }
2605a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX)
2613cb270beSHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
262e77caa6dSBarry Smith   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));
2633cb270beSHong Zhang   #else
264e77caa6dSBarry Smith   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));
2653cb270beSHong Zhang   #endif
2663cb270beSHong Zhang #else
2673cb270beSHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
268e77caa6dSBarry Smith   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));
2695a46d3faSBarry Smith   #else
270e77caa6dSBarry Smith   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));
2715a46d3faSBarry Smith   #endif
2723cb270beSHong Zhang #endif
2735a46d3faSBarry Smith 
2745a46d3faSBarry Smith   /* Numerical factorization */
2755a46d3faSBarry Smith   lu->B.ncol = 0; /* Indicate not to solve the system */
276d5f3da31SBarry Smith   if (F->factortype == MAT_FACTOR_LU) {
2775a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX)
2783cb270beSHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
2799371c9d4SSatish Balay     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));
2803cb270beSHong Zhang   #else
2819371c9d4SSatish Balay     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));
2823cb270beSHong Zhang   #endif
2833cb270beSHong Zhang #else
2843cb270beSHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
2859371c9d4SSatish Balay     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));
2865a46d3faSBarry Smith   #else
2879371c9d4SSatish Balay     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));
2885a46d3faSBarry Smith   #endif
2893cb270beSHong Zhang #endif
290d5f3da31SBarry Smith   } else if (F->factortype == MAT_FACTOR_ILU) {
291cffbb591SHong Zhang     /* Compute the incomplete factorization, condition number and pivot growth */
292cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX)
2933cb270beSHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
2949371c9d4SSatish Balay     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));
2953cb270beSHong Zhang   #else
2969371c9d4SSatish Balay     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));
2973cb270beSHong Zhang   #endif
2983cb270beSHong Zhang #else
2993cb270beSHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
3009371c9d4SSatish Balay     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));
301cffbb591SHong Zhang   #else
3029371c9d4SSatish Balay     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));
303cffbb591SHong Zhang   #endif
3043cb270beSHong Zhang #endif
305f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported");
3065a46d3faSBarry Smith   if (!sinfo || sinfo == lu->A.ncol + 1) {
30748a46eb9SPierre Jolivet     if (lu->options.PivotGrowth) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Recip. pivot growth = %e\n", lu->rpg));
30848a46eb9SPierre Jolivet     if (lu->options.ConditionNumber) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Recip. condition number = %e\n", lu->rcond));
3095a46d3faSBarry Smith   } else if (sinfo > 0) {
310675d1226SHong Zhang     if (A->erroriffailure) {
31198921bdaSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot in row %" PetscInt_FMT, sinfo);
312675d1226SHong Zhang     } else {
313675d1226SHong Zhang       if (sinfo <= lu->A.ncol) {
314ad540459SPierre Jolivet         if (lu->options.ILU_FillTol == 0.0) F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3159566063dSJacob Faibussowitsch         PetscCall(PetscInfo(F, "Number of zero pivots %" PetscInt_FMT ", ILU_FillTol %g\n", sinfo, lu->options.ILU_FillTol));
316675d1226SHong Zhang       } else if (sinfo == lu->A.ncol + 1) {
317675d1226SHong Zhang         /*
318675d1226SHong Zhang          U is nonsingular, but RCOND is less than machine
319675d1226SHong Zhang                       precision, meaning that the matrix is singular to
320675d1226SHong Zhang                       working precision. Nevertheless, the solution and
321675d1226SHong Zhang                       error bounds are computed because there are a number
322675d1226SHong Zhang                       of situations where the computed solution can be more
323675d1226SHong Zhang                       accurate than the value of RCOND would suggest.
324675d1226SHong Zhang          */
3259566063dSJacob Faibussowitsch         PetscCall(PetscInfo(F, "Matrix factor U is nonsingular, but is singular to working precision. The solution is computed. info %" PetscInt_FMT, sinfo));
326675d1226SHong Zhang       } else { /* sinfo > lu->A.ncol + 1 */
327603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
3289566063dSJacob Faibussowitsch         PetscCall(PetscInfo(F, "Number of bytes allocated when memory allocation fails %" PetscInt_FMT "\n", sinfo));
329675d1226SHong Zhang       }
330675d1226SHong Zhang     }
33198921bdaSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "info = %" PetscInt_FMT ", the %" PetscInt_FMT "-th argument in gssvx() had an illegal value", sinfo, -sinfo);
3325a46d3faSBarry Smith 
3335a46d3faSBarry Smith   if (lu->options.PrintStat) {
3349566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "MatLUFactorNumeric_SuperLU():\n"));
335e77caa6dSBarry Smith     PetscStackCallExternalVoid("SuperLU:StatPrint", StatPrint(&lu->stat));
3365a46d3faSBarry Smith     Lstore = (SCformat *)lu->L.Store;
3375a46d3faSBarry Smith     Ustore = (NCformat *)lu->U.Store;
3389566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "  No of nonzeros in factor L = %" PetscInt_FMT "\n", Lstore->nnz));
3399566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "  No of nonzeros in factor U = %" PetscInt_FMT "\n", Ustore->nnz));
3409566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "  No of nonzeros in L+U = %" PetscInt_FMT "\n", Lstore->nnz + Ustore->nnz - lu->A.ncol));
3419566063dSJacob Faibussowitsch     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));
3425a46d3faSBarry Smith   }
3435a46d3faSBarry Smith 
3445a46d3faSBarry Smith   lu->flg                = SAME_NONZERO_PATTERN;
3451d5ca7f3SHong Zhang   F->ops->solve          = MatSolve_SuperLU;
3461d5ca7f3SHong Zhang   F->ops->solvetranspose = MatSolveTranspose_SuperLU;
3471aef8b4cSStefano Zampini   F->ops->matsolve       = NULL;
3483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3495a46d3faSBarry Smith }
3505a46d3faSBarry Smith 
351d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_SuperLU(Mat A)
352d71ae5a4SJacob Faibussowitsch {
353245fece9SBarry Smith   Mat_SuperLU *lu = (Mat_SuperLU *)A->data;
35414b396bbSKris Buschelman 
35514b396bbSKris Buschelman   PetscFunctionBegin;
356245fece9SBarry Smith   if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
357e77caa6dSBarry Smith     PetscStackCallExternalVoid("SuperLU:Destroy_SuperMatrix_Store", Destroy_SuperMatrix_Store(&lu->A));
358e77caa6dSBarry Smith     PetscStackCallExternalVoid("SuperLU:Destroy_SuperMatrix_Store", Destroy_SuperMatrix_Store(&lu->B));
359e77caa6dSBarry Smith     PetscStackCallExternalVoid("SuperLU:Destroy_SuperMatrix_Store", Destroy_SuperMatrix_Store(&lu->X));
360e77caa6dSBarry Smith     PetscStackCallExternalVoid("SuperLU:StatFree", StatFree(&lu->stat));
3610e742b69SHong Zhang     if (lu->lwork >= 0) {
362e77caa6dSBarry Smith       PetscStackCallExternalVoid("SuperLU:Destroy_SuperNode_Matrix", Destroy_SuperNode_Matrix(&lu->L));
363e77caa6dSBarry Smith       PetscStackCallExternalVoid("SuperLU:Destroy_CompCol_Matrix", Destroy_CompCol_Matrix(&lu->U));
3640e742b69SHong Zhang     }
3650e742b69SHong Zhang   }
3669566063dSJacob Faibussowitsch   PetscCall(PetscFree(lu->etree));
3679566063dSJacob Faibussowitsch   PetscCall(PetscFree(lu->perm_r));
3689566063dSJacob Faibussowitsch   PetscCall(PetscFree(lu->perm_c));
3699566063dSJacob Faibussowitsch   PetscCall(PetscFree(lu->R));
3709566063dSJacob Faibussowitsch   PetscCall(PetscFree(lu->C));
3719566063dSJacob Faibussowitsch   PetscCall(PetscFree(lu->rhs_dup));
3729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lu->A_dup));
3739566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
3740e742b69SHong Zhang 
375d954db57SHong Zhang   /* clear composed functions */
3769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL));
3779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSuperluSetILUDropTol_C", NULL));
3783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
37914b396bbSKris Buschelman }
38014b396bbSKris Buschelman 
381d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_SuperLU(Mat A, PetscViewer viewer)
382d71ae5a4SJacob Faibussowitsch {
383ace3abfcSBarry Smith   PetscBool         iascii;
38414b396bbSKris Buschelman   PetscViewerFormat format;
38514b396bbSKris Buschelman 
38614b396bbSKris Buschelman   PetscFunctionBegin;
3879566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
38832077d6dSBarry Smith   if (iascii) {
3899566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
39048a46eb9SPierre Jolivet     if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_SuperLU(A, viewer));
39114b396bbSKris Buschelman   }
3923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39314b396bbSKris Buschelman }
39414b396bbSKris Buschelman 
395d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatSolve_SuperLU(Mat A, Mat B, Mat X)
396d71ae5a4SJacob Faibussowitsch {
397245fece9SBarry Smith   Mat_SuperLU *lu = (Mat_SuperLU *)A->data;
398cd723cd1SBarry Smith   PetscBool    flg;
399e0b74bf9SHong Zhang 
400e0b74bf9SHong Zhang   PetscFunctionBegin;
4019566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATSEQDENSE, MATMPIDENSE, NULL));
4025f80ce2aSJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix B must be MATDENSE matrix");
4039566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &flg, MATSEQDENSE, MATMPIDENSE, NULL));
4045f80ce2aSJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix X must be MATDENSE matrix");
4052205254eSKarl Rupp   lu->options.Trans = TRANS;
406e0b74bf9SHong Zhang   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatMatSolve_SuperLU() is not implemented yet");
4073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
408e0b74bf9SHong Zhang }
409e0b74bf9SHong Zhang 
410d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info)
411d71ae5a4SJacob Faibussowitsch {
412245fece9SBarry Smith   Mat_SuperLU *lu = (Mat_SuperLU *)(F->data);
41326cc229bSBarry Smith   PetscInt     indx;
41426cc229bSBarry Smith   PetscBool    flg, set;
41526cc229bSBarry Smith   PetscReal    real_input;
41626cc229bSBarry Smith   const char  *colperm[]    = {"NATURAL", "MMD_ATA", "MMD_AT_PLUS_A", "COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
41726cc229bSBarry Smith   const char  *iterrefine[] = {"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
41826cc229bSBarry Smith   const char  *rowperm[]    = {"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
419b24902e0SBarry Smith 
420b24902e0SBarry Smith   PetscFunctionBegin;
42126cc229bSBarry Smith   /* Set options to F */
42226cc229bSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "SuperLU Options", "Mat");
42326cc229bSBarry Smith   PetscCall(PetscOptionsBool("-mat_superlu_equil", "Equil", "None", (PetscBool)lu->options.Equil, (PetscBool *)&lu->options.Equil, NULL));
42426cc229bSBarry Smith   PetscCall(PetscOptionsEList("-mat_superlu_colperm", "ColPerm", "None", colperm, 4, colperm[3], &indx, &flg));
42526cc229bSBarry Smith   if (flg) lu->options.ColPerm = (colperm_t)indx;
42626cc229bSBarry Smith   PetscCall(PetscOptionsEList("-mat_superlu_iterrefine", "IterRefine", "None", iterrefine, 4, iterrefine[0], &indx, &flg));
42726cc229bSBarry Smith   if (flg) lu->options.IterRefine = (IterRefine_t)indx;
42826cc229bSBarry Smith   PetscCall(PetscOptionsBool("-mat_superlu_symmetricmode", "SymmetricMode", "None", (PetscBool)lu->options.SymmetricMode, &flg, &set));
42926cc229bSBarry Smith   if (set && flg) lu->options.SymmetricMode = YES;
43026cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_superlu_diagpivotthresh", "DiagPivotThresh", "None", lu->options.DiagPivotThresh, &real_input, &flg));
43126cc229bSBarry Smith   if (flg) lu->options.DiagPivotThresh = (double)real_input;
43226cc229bSBarry Smith   PetscCall(PetscOptionsBool("-mat_superlu_pivotgrowth", "PivotGrowth", "None", (PetscBool)lu->options.PivotGrowth, &flg, &set));
43326cc229bSBarry Smith   if (set && flg) lu->options.PivotGrowth = YES;
43426cc229bSBarry Smith   PetscCall(PetscOptionsBool("-mat_superlu_conditionnumber", "ConditionNumber", "None", (PetscBool)lu->options.ConditionNumber, &flg, &set));
43526cc229bSBarry Smith   if (set && flg) lu->options.ConditionNumber = YES;
43626cc229bSBarry Smith   PetscCall(PetscOptionsEList("-mat_superlu_rowperm", "rowperm", "None", rowperm, 2, rowperm[lu->options.RowPerm], &indx, &flg));
43726cc229bSBarry Smith   if (flg) lu->options.RowPerm = (rowperm_t)indx;
43826cc229bSBarry Smith   PetscCall(PetscOptionsBool("-mat_superlu_replacetinypivot", "ReplaceTinyPivot", "None", (PetscBool)lu->options.ReplaceTinyPivot, &flg, &set));
43926cc229bSBarry Smith   if (set && flg) lu->options.ReplaceTinyPivot = YES;
44026cc229bSBarry Smith   PetscCall(PetscOptionsBool("-mat_superlu_printstat", "PrintStat", "None", (PetscBool)lu->options.PrintStat, &flg, &set));
44126cc229bSBarry Smith   if (set && flg) lu->options.PrintStat = YES;
44226cc229bSBarry Smith   PetscCall(PetscOptionsInt("-mat_superlu_lwork", "size of work array in bytes used by factorization", "None", lu->lwork, &lu->lwork, NULL));
44326cc229bSBarry Smith   if (lu->lwork > 0) {
44426cc229bSBarry Smith     /* lwork is in bytes, hence PetscMalloc() is used here, not PetscMalloc1()*/
44526cc229bSBarry Smith     PetscCall(PetscMalloc(lu->lwork, &lu->work));
44626cc229bSBarry Smith   } else if (lu->lwork != 0 && lu->lwork != -1) {
44726cc229bSBarry Smith     PetscCall(PetscPrintf(PETSC_COMM_SELF, "   Warning: lwork %" PetscInt_FMT " is not supported by SUPERLU. The default lwork=0 is used.\n", lu->lwork));
44826cc229bSBarry Smith     lu->lwork = 0;
44926cc229bSBarry Smith   }
45026cc229bSBarry Smith   /* ilu options */
45126cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_superlu_ilu_droptol", "ILU_DropTol", "None", lu->options.ILU_DropTol, &real_input, &flg));
45226cc229bSBarry Smith   if (flg) lu->options.ILU_DropTol = (double)real_input;
45326cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_superlu_ilu_filltol", "ILU_FillTol", "None", lu->options.ILU_FillTol, &real_input, &flg));
45426cc229bSBarry Smith   if (flg) lu->options.ILU_FillTol = (double)real_input;
45526cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_superlu_ilu_fillfactor", "ILU_FillFactor", "None", lu->options.ILU_FillFactor, &real_input, &flg));
45626cc229bSBarry Smith   if (flg) lu->options.ILU_FillFactor = (double)real_input;
45726cc229bSBarry Smith   PetscCall(PetscOptionsInt("-mat_superlu_ilu_droprull", "ILU_DropRule", "None", lu->options.ILU_DropRule, &lu->options.ILU_DropRule, NULL));
45826cc229bSBarry Smith   PetscCall(PetscOptionsInt("-mat_superlu_ilu_norm", "ILU_Norm", "None", lu->options.ILU_Norm, &indx, &flg));
45926cc229bSBarry Smith   if (flg) lu->options.ILU_Norm = (norm_t)indx;
46026cc229bSBarry Smith   PetscCall(PetscOptionsInt("-mat_superlu_ilu_milu", "ILU_MILU", "None", lu->options.ILU_MILU, &indx, &flg));
46126cc229bSBarry Smith   if (flg) lu->options.ILU_MILU = (milu_t)indx;
46226cc229bSBarry Smith   PetscOptionsEnd();
46326cc229bSBarry Smith 
464b24902e0SBarry Smith   lu->flg                 = DIFFERENT_NONZERO_PATTERN;
465b24902e0SBarry Smith   lu->CleanUpSuperLU      = PETSC_TRUE;
4661d5ca7f3SHong Zhang   F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
4672366e350SStefano Zampini 
4682366e350SStefano Zampini   /* if we are here, the nonzero pattern has changed unless the user explicitly called MatLUFactorSymbolic */
4699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lu->A_dup));
47048a46eb9SPierre Jolivet   if (lu->needconversion) PetscCall(MatConvert(A, MATSEQAIJ, MAT_INITIAL_MATRIX, &lu->A_dup));
4712366e350SStefano Zampini   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 */
4729566063dSJacob Faibussowitsch     PetscCall(MatDuplicate_SeqAIJ(A, MAT_COPY_VALUES, &lu->A_dup));
4732366e350SStefano Zampini   }
4743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
475b24902e0SBarry Smith }
476b24902e0SBarry Smith 
477d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F, PetscReal dtol)
478d71ae5a4SJacob Faibussowitsch {
479245fece9SBarry Smith   Mat_SuperLU *lu = (Mat_SuperLU *)F->data;
4805ccb76cbSHong Zhang 
4815ccb76cbSHong Zhang   PetscFunctionBegin;
4825ccb76cbSHong Zhang   lu->options.ILU_DropTol = dtol;
4833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4845ccb76cbSHong Zhang }
4855ccb76cbSHong Zhang 
4865ccb76cbSHong Zhang /*@
4875ccb76cbSHong Zhang   MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance
488147403d9SBarry Smith 
489c3339decSBarry Smith    Logically Collective
4905ccb76cbSHong Zhang 
4915ccb76cbSHong Zhang    Input Parameters:
4922ef1f0ffSBarry Smith +  F - the factored matrix obtained by calling `MatGetFactor()`
4935ccb76cbSHong Zhang -  dtol - drop tolerance
4945ccb76cbSHong Zhang 
4953c7db156SBarry Smith   Options Database Key:
496147403d9SBarry Smith .   -mat_superlu_ilu_droptol <dtol> - the drop tolerance
4975ccb76cbSHong Zhang 
4985ccb76cbSHong Zhang    Level: beginner
4995ccb76cbSHong Zhang 
50096a0c994SBarry Smith    References:
501606c0280SSatish Balay .  * - SuperLU Users' Guide
5025ccb76cbSHong Zhang 
503*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MATSOLVERSUPERLU`
5045ccb76cbSHong Zhang @*/
505d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSuperluSetILUDropTol(Mat F, PetscReal dtol)
506d71ae5a4SJacob Faibussowitsch {
5075ccb76cbSHong Zhang   PetscFunctionBegin;
5085ccb76cbSHong Zhang   PetscValidHeaderSpecific(F, MAT_CLASSID, 1);
50969fbec6eSBarry Smith   PetscValidLogicalCollectiveReal(F, dtol, 2);
510cac4c232SBarry Smith   PetscTryMethod(F, "MatSuperluSetILUDropTol_C", (Mat, PetscReal), (F, dtol));
5113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5125ccb76cbSHong Zhang }
5135ccb76cbSHong Zhang 
514d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFactorGetSolverType_seqaij_superlu(Mat A, MatSolverType *type)
515d71ae5a4SJacob Faibussowitsch {
51635bd34faSBarry Smith   PetscFunctionBegin;
5172692d6eeSBarry Smith   *type = MATSOLVERSUPERLU;
5183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
51935bd34faSBarry Smith }
52035bd34faSBarry Smith 
521b24902e0SBarry Smith /*MC
522ba992d64SSatish Balay   MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices
523b24902e0SBarry Smith   via the external package SuperLU.
524b24902e0SBarry Smith 
5252ef1f0ffSBarry Smith   Use `./configure --download-superlu` to have PETSc installed with SuperLU
526b24902e0SBarry Smith 
5272ef1f0ffSBarry Smith   Use `-pc_type lu` `-pc_factor_mat_solver_type superlu` to use this direct solver
528c2b89b5dSBarry Smith 
529b24902e0SBarry Smith   Options Database Keys:
530e08999f5SMatthew G Knepley + -mat_superlu_equil <FALSE>            - Equil (None)
5312ef1f0ffSBarry Smith . -mat_superlu_colperm <COLAMD>         - (choose one of) `NATURAL`, `MMD_ATA MMD_AT_PLUS_A`, `COLAMD`
5322ef1f0ffSBarry Smith . -mat_superlu_iterrefine <NOREFINE>    - (choose one of) `NOREFINE`, `SINGLE`, `DOUBLE`, `EXTRA`
533e08999f5SMatthew G Knepley . -mat_superlu_symmetricmode: <FALSE>   - SymmetricMode (None)
534e08999f5SMatthew G Knepley . -mat_superlu_diagpivotthresh <1>      - DiagPivotThresh (None)
535e08999f5SMatthew G Knepley . -mat_superlu_pivotgrowth <FALSE>      - PivotGrowth (None)
536e08999f5SMatthew G Knepley . -mat_superlu_conditionnumber <FALSE>  - ConditionNumber (None)
5372ef1f0ffSBarry Smith . -mat_superlu_rowperm <NOROWPERM>      - (choose one of) `NOROWPERM`, `LargeDiag`
538e08999f5SMatthew G Knepley . -mat_superlu_replacetinypivot <FALSE> - ReplaceTinyPivot (None)
539e08999f5SMatthew G Knepley . -mat_superlu_printstat <FALSE>        - PrintStat (None)
540e08999f5SMatthew G Knepley . -mat_superlu_lwork <0>                - size of work array in bytes used by factorization (None)
541e08999f5SMatthew G Knepley . -mat_superlu_ilu_droptol <0>          - ILU_DropTol (None)
542e08999f5SMatthew G Knepley . -mat_superlu_ilu_filltol <0>          - ILU_FillTol (None)
543e08999f5SMatthew G Knepley . -mat_superlu_ilu_fillfactor <0>       - ILU_FillFactor (None)
544e08999f5SMatthew G Knepley . -mat_superlu_ilu_droprull <0>         - ILU_DropRule (None)
545e08999f5SMatthew G Knepley . -mat_superlu_ilu_norm <0>             - ILU_Norm (None)
546e08999f5SMatthew G Knepley - -mat_superlu_ilu_milu <0>             - ILU_MILU (None)
547b24902e0SBarry Smith 
5482ef1f0ffSBarry Smith    Level: beginner
5492ef1f0ffSBarry Smith 
55095452b02SPatrick Sanan    Notes:
55111a5261eSBarry Smith     Do not confuse this with `MATSOLVERSUPERLU_DIST` which is for parallel sparse solves
5525c9eb25fSBarry Smith 
55311a5261eSBarry Smith     Cannot use ordering provided by PETSc, provides its own.
5542c7c0729SBarry Smith 
555*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCLU`, `PCILU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`, `PCFactorSetMatSolverType()`, `MatSolverType`
556b24902e0SBarry Smith M*/
557b24902e0SBarry Smith 
558d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_seqaij_superlu(Mat A, MatFactorType ftype, Mat *F)
559d71ae5a4SJacob Faibussowitsch {
56014b396bbSKris Buschelman   Mat          B;
561f0c56d0fSKris Buschelman   Mat_SuperLU *lu;
56226cc229bSBarry Smith   PetscInt     m = A->rmap->n, n = A->cmap->n;
56314b396bbSKris Buschelman 
56414b396bbSKris Buschelman   PetscFunctionBegin;
5659566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
5669566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, PETSC_DETERMINE, PETSC_DETERMINE));
5679566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("superlu", &((PetscObject)B)->type_name));
5689566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
56966e17bc3SBarry Smith   B->trivialsymbolic = PETSC_TRUE;
570cffbb591SHong Zhang   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
571b24902e0SBarry Smith     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_SuperLU;
572cffbb591SHong Zhang     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU;
573b3a44c85SBarry Smith   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported");
574cffbb591SHong Zhang 
5759566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
5769566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERSUPERLU, &B->solvertype));
57700c67f3bSHong Zhang 
578b9c12af5SBarry Smith   B->ops->getinfo = MatGetInfo_External;
579b24902e0SBarry Smith   B->ops->destroy = MatDestroy_SuperLU;
5803519fcdcSHong Zhang   B->ops->view    = MatView_SuperLU;
581d5f3da31SBarry Smith   B->factortype   = ftype;
58294b7f48cSBarry Smith   B->assembled    = PETSC_TRUE; /* required by -ksp_view */
5835c9eb25fSBarry Smith   B->preallocated = PETSC_TRUE;
58414b396bbSKris Buschelman 
5854dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&lu));
586cae5a23dSHong Zhang 
587cffbb591SHong Zhang   if (ftype == MAT_FACTOR_LU) {
5889ce50919SHong Zhang     set_default_options(&lu->options);
5893d6581fbSHong Zhang     /* Comments from SuperLU_4.0/SRC/dgssvx.c:
5903d6581fbSHong Zhang       "Whether or not the system will be equilibrated depends on the
5913d6581fbSHong Zhang        scaling of the matrix A, but if equilibration is used, A is
5923d6581fbSHong Zhang        overwritten by diag(R)*A*diag(C) and B by diag(R)*B
5933d6581fbSHong Zhang        (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)."
5943d6581fbSHong Zhang      We set 'options.Equil = NO' as default because additional space is needed for it.
5953d6581fbSHong Zhang     */
5963d6581fbSHong Zhang     lu->options.Equil = NO;
597cffbb591SHong Zhang   } else if (ftype == MAT_FACTOR_ILU) {
5980c74a584SJed Brown     /* Set the default input options of ilu: */
599e77caa6dSBarry Smith     PetscStackCallExternalVoid("SuperLU:ilu_set_default_options", ilu_set_default_options(&lu->options));
600cffbb591SHong Zhang   }
6019ce50919SHong Zhang   lu->options.PrintStat = NO;
6021a2e2f44SHong Zhang 
6035d8b2efaSHong Zhang   /* Initialize the statistics variables. */
604e77caa6dSBarry Smith   PetscStackCallExternalVoid("SuperLU:StatInit", StatInit(&lu->stat));
6058914a3f7SHong Zhang   lu->lwork = 0; /* allocate space internally by system malloc */
6069ce50919SHong Zhang 
6075d8b2efaSHong Zhang   /* Allocate spaces (notice sizes are for the transpose) */
6089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &lu->etree));
6099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &lu->perm_r));
6109566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &lu->perm_c));
6119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &lu->R));
6129566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &lu->C));
6135d8b2efaSHong Zhang 
6145d8b2efaSHong Zhang   /* create rhs and solution x without allocate space for .Store */
6155d8b2efaSHong Zhang #if defined(PETSC_USE_COMPLEX)
6163cb270beSHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
617e77caa6dSBarry Smith   PetscStackCallExternalVoid("SuperLU:cCreate_Dense_Matrix(", cCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE));
618e77caa6dSBarry Smith   PetscStackCallExternalVoid("SuperLU:cCreate_Dense_Matrix(", cCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE));
6193cb270beSHong Zhang   #else
620e77caa6dSBarry Smith   PetscStackCallExternalVoid("SuperLU:zCreate_Dense_Matrix", zCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE));
621e77caa6dSBarry Smith   PetscStackCallExternalVoid("SuperLU:zCreate_Dense_Matrix", zCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE));
6223cb270beSHong Zhang   #endif
6233cb270beSHong Zhang #else
6243cb270beSHong Zhang   #if defined(PETSC_USE_REAL_SINGLE)
625e77caa6dSBarry Smith   PetscStackCallExternalVoid("SuperLU:sCreate_Dense_Matrix", sCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE));
626e77caa6dSBarry Smith   PetscStackCallExternalVoid("SuperLU:sCreate_Dense_Matrix", sCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE));
6275d8b2efaSHong Zhang   #else
628e77caa6dSBarry Smith   PetscStackCallExternalVoid("SuperLU:dCreate_Dense_Matrix", dCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE));
629e77caa6dSBarry Smith   PetscStackCallExternalVoid("SuperLU:dCreate_Dense_Matrix", dCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE));
6305d8b2efaSHong Zhang   #endif
6313cb270beSHong Zhang #endif
6325d8b2efaSHong Zhang 
6339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_superlu));
6349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSuperluSetILUDropTol_C", MatSuperluSetILUDropTol_SuperLU));
635245fece9SBarry Smith   B->data = lu;
63620be8e61SHong Zhang 
637899d7b4fSKris Buschelman   *F = B;
6383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63914b396bbSKris Buschelman }
640d954db57SHong Zhang 
641d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_seqsell_superlu(Mat A, MatFactorType ftype, Mat *F)
642d71ae5a4SJacob Faibussowitsch {
6432366e350SStefano Zampini   Mat_SuperLU *lu;
6442366e350SStefano Zampini 
6452366e350SStefano Zampini   PetscFunctionBegin;
6469566063dSJacob Faibussowitsch   PetscCall(MatGetFactor_seqaij_superlu(A, ftype, F));
6472366e350SStefano Zampini   lu                 = (Mat_SuperLU *)((*F)->data);
6482366e350SStefano Zampini   lu->needconversion = PETSC_TRUE;
6493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6502366e350SStefano Zampini }
6512366e350SStefano Zampini 
652d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_SuperLU(void)
653d71ae5a4SJacob Faibussowitsch {
65442c9c57cSBarry Smith   PetscFunctionBegin;
6559566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSUPERLU, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_seqaij_superlu));
6569566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSUPERLU, MATSEQAIJ, MAT_FACTOR_ILU, MatGetFactor_seqaij_superlu));
6579566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSUPERLU, MATSEQSELL, MAT_FACTOR_LU, MatGetFactor_seqsell_superlu));
6589566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSUPERLU, MATSEQSELL, MAT_FACTOR_ILU, MatGetFactor_seqsell_superlu));
6593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
66042c9c57cSBarry Smith }
661