xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 1d27aa22b2f6148b2c4e3f06a75e0638d6493e09)
11677d0b8SKris Buschelman /*
2d515b9b4SShri Abhyankar    Provides an interface to the UMFPACK sparse solver available through SuiteSparse version 4.2.1
32d4e2982SHong Zhang 
49e475b0dSSatish Balay    When build with PETSC_USE_64BIT_INDICES this will use Suitesparse_long as the
58592901bSBarry Smith    integer type in UMFPACK, otherwise it will use int. This means
68592901bSBarry Smith    all integers in this file as simply declared as PetscInt. Also it means
77de69702SBarry Smith    that one cannot use 64BIT_INDICES on 32-bit pointer systems [as Suitesparse_long is 32-bit only]
82d4e2982SHong Zhang 
91677d0b8SKris Buschelman */
10c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
111677d0b8SKris Buschelman 
128592901bSBarry Smith #if defined(PETSC_USE_64BIT_INDICES)
138592901bSBarry Smith   #if defined(PETSC_USE_COMPLEX)
148592901bSBarry Smith     #define umfpack_UMF_free_symbolic umfpack_zl_free_symbolic
158592901bSBarry Smith     #define umfpack_UMF_free_numeric  umfpack_zl_free_numeric
16d755ee67SBarry Smith     /* the type casts are needed because PetscInt is long long while SuiteSparse_long is long and compilers warn even when they are identical */
17d755ee67SBarry Smith     #define umfpack_UMF_wsolve(a, b, c, d, e, f, g, h, i, j, k, l, m, n) umfpack_zl_wsolve(a, (SuiteSparse_long *)b, (SuiteSparse_long *)c, d, e, f, g, h, i, (SuiteSparse_long *)j, k, l, (SuiteSparse_long *)m, n)
18d755ee67SBarry Smith     #define umfpack_UMF_numeric(a, b, c, d, e, f, g, h)                  umfpack_zl_numeric((SuiteSparse_long *)a, (SuiteSparse_long *)b, c, d, e, f, g, h)
198592901bSBarry Smith     #define umfpack_UMF_report_numeric                                   umfpack_zl_report_numeric
208592901bSBarry Smith     #define umfpack_UMF_report_control                                   umfpack_zl_report_control
218592901bSBarry Smith     #define umfpack_UMF_report_status                                    umfpack_zl_report_status
228592901bSBarry Smith     #define umfpack_UMF_report_info                                      umfpack_zl_report_info
238592901bSBarry Smith     #define umfpack_UMF_report_symbolic                                  umfpack_zl_report_symbolic
24d755ee67SBarry Smith     #define umfpack_UMF_qsymbolic(a, b, c, d, e, f, g, h, i, j)          umfpack_zl_qsymbolic(a, b, (SuiteSparse_long *)c, (SuiteSparse_long *)d, e, f, (SuiteSparse_long *)g, h, i, j)
25d755ee67SBarry Smith     #define umfpack_UMF_symbolic(a, b, c, d, e, f, g, h, i)              umfpack_zl_symbolic(a, b, (SuiteSparse_long *)c, (SuiteSparse_long *)d, e, f, g, h, i)
268592901bSBarry Smith     #define umfpack_UMF_defaults                                         umfpack_zl_defaults
278592901bSBarry Smith 
288592901bSBarry Smith   #else
298592901bSBarry Smith     #define umfpack_UMF_free_symbolic                           umfpack_dl_free_symbolic
308592901bSBarry Smith     #define umfpack_UMF_free_numeric                            umfpack_dl_free_numeric
31d755ee67SBarry Smith     #define umfpack_UMF_wsolve(a, b, c, d, e, f, g, h, i, j, k) umfpack_dl_wsolve(a, (SuiteSparse_long *)b, (SuiteSparse_long *)c, d, e, f, g, h, i, (SuiteSparse_long *)j, k)
32d755ee67SBarry Smith     #define umfpack_UMF_numeric(a, b, c, d, e, f, g)            umfpack_dl_numeric((SuiteSparse_long *)a, (SuiteSparse_long *)b, c, d, e, f, g)
338592901bSBarry Smith     #define umfpack_UMF_report_numeric                          umfpack_dl_report_numeric
348592901bSBarry Smith     #define umfpack_UMF_report_control                          umfpack_dl_report_control
358592901bSBarry Smith     #define umfpack_UMF_report_status                           umfpack_dl_report_status
368592901bSBarry Smith     #define umfpack_UMF_report_info                             umfpack_dl_report_info
378592901bSBarry Smith     #define umfpack_UMF_report_symbolic                         umfpack_dl_report_symbolic
38d755ee67SBarry Smith     #define umfpack_UMF_qsymbolic(a, b, c, d, e, f, g, h, i)    umfpack_dl_qsymbolic(a, b, (SuiteSparse_long *)c, (SuiteSparse_long *)d, e, (SuiteSparse_long *)f, g, h, i)
39d755ee67SBarry Smith     #define umfpack_UMF_symbolic(a, b, c, d, e, f, g, h)        umfpack_dl_symbolic(a, b, (SuiteSparse_long *)c, (SuiteSparse_long *)d, e, f, g, h)
408592901bSBarry Smith     #define umfpack_UMF_defaults                                umfpack_dl_defaults
418592901bSBarry Smith   #endif
428592901bSBarry Smith 
438592901bSBarry Smith #else
448592901bSBarry Smith   #if defined(PETSC_USE_COMPLEX)
458592901bSBarry Smith     #define umfpack_UMF_free_symbolic   umfpack_zi_free_symbolic
468592901bSBarry Smith     #define umfpack_UMF_free_numeric    umfpack_zi_free_numeric
478592901bSBarry Smith     #define umfpack_UMF_wsolve          umfpack_zi_wsolve
488592901bSBarry Smith     #define umfpack_UMF_numeric         umfpack_zi_numeric
498592901bSBarry Smith     #define umfpack_UMF_report_numeric  umfpack_zi_report_numeric
508592901bSBarry Smith     #define umfpack_UMF_report_control  umfpack_zi_report_control
518592901bSBarry Smith     #define umfpack_UMF_report_status   umfpack_zi_report_status
528592901bSBarry Smith     #define umfpack_UMF_report_info     umfpack_zi_report_info
538592901bSBarry Smith     #define umfpack_UMF_report_symbolic umfpack_zi_report_symbolic
548592901bSBarry Smith     #define umfpack_UMF_qsymbolic       umfpack_zi_qsymbolic
558592901bSBarry Smith     #define umfpack_UMF_symbolic        umfpack_zi_symbolic
568592901bSBarry Smith     #define umfpack_UMF_defaults        umfpack_zi_defaults
578592901bSBarry Smith 
588592901bSBarry Smith   #else
598592901bSBarry Smith     #define umfpack_UMF_free_symbolic   umfpack_di_free_symbolic
608592901bSBarry Smith     #define umfpack_UMF_free_numeric    umfpack_di_free_numeric
618592901bSBarry Smith     #define umfpack_UMF_wsolve          umfpack_di_wsolve
628592901bSBarry Smith     #define umfpack_UMF_numeric         umfpack_di_numeric
638592901bSBarry Smith     #define umfpack_UMF_report_numeric  umfpack_di_report_numeric
648592901bSBarry Smith     #define umfpack_UMF_report_control  umfpack_di_report_control
658592901bSBarry Smith     #define umfpack_UMF_report_status   umfpack_di_report_status
668592901bSBarry Smith     #define umfpack_UMF_report_info     umfpack_di_report_info
678592901bSBarry Smith     #define umfpack_UMF_report_symbolic umfpack_di_report_symbolic
688592901bSBarry Smith     #define umfpack_UMF_qsymbolic       umfpack_di_qsymbolic
698592901bSBarry Smith     #define umfpack_UMF_symbolic        umfpack_di_symbolic
708592901bSBarry Smith     #define umfpack_UMF_defaults        umfpack_di_defaults
718592901bSBarry Smith   #endif
728592901bSBarry Smith #endif
738592901bSBarry Smith 
741677d0b8SKris Buschelman EXTERN_C_BEGIN
75c6db04a5SJed Brown #include <umfpack.h>
761677d0b8SKris Buschelman EXTERN_C_END
771677d0b8SKris Buschelman 
78fcfd50ebSBarry Smith static const char *const UmfpackOrderingTypes[] = {"CHOLMOD", "AMD", "GIVEN", "METIS", "BEST", "NONE", "USER", "UmfpackOrderingTypes", "UMFPACK_ORDERING_", 0};
79fcfd50ebSBarry Smith 
801677d0b8SKris Buschelman typedef struct {
811677d0b8SKris Buschelman   void        *Symbolic, *Numeric;
821677d0b8SKris Buschelman   double       Info[UMFPACK_INFO], Control[UMFPACK_CONTROL], *W;
83ae26a541SJed Brown   PetscInt    *Wi, *perm_c;
84ae26a541SJed Brown   Mat          A; /* Matrix used for factorization */
851677d0b8SKris Buschelman   MatStructure flg;
861677d0b8SKris Buschelman 
871677d0b8SKris Buschelman   /* Flag to clean up UMFPACK objects during Destroy */
88ace3abfcSBarry Smith   PetscBool CleanUpUMFPACK;
89f0c56d0fSKris Buschelman } Mat_UMFPACK;
90f0c56d0fSKris Buschelman 
91d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_UMFPACK(Mat A)
92d71ae5a4SJacob Faibussowitsch {
93b9c12af5SBarry Smith   Mat_UMFPACK *lu = (Mat_UMFPACK *)A->data;
941677d0b8SKris Buschelman 
951677d0b8SKris Buschelman   PetscFunctionBegin;
96b9c12af5SBarry Smith   if (lu->CleanUpUMFPACK) {
978592901bSBarry Smith     umfpack_UMF_free_symbolic(&lu->Symbolic);
988592901bSBarry Smith     umfpack_UMF_free_numeric(&lu->Numeric);
999566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->Wi));
1009566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->W));
1019566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->perm_c));
1021677d0b8SKris Buschelman   }
1039566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lu->A));
1042e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL));
1059566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
1063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1071677d0b8SKris Buschelman }
1081677d0b8SKris Buschelman 
109d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_UMFPACK_Private(Mat A, Vec b, Vec x, int uflag)
110d71ae5a4SJacob Faibussowitsch {
111b9c12af5SBarry Smith   Mat_UMFPACK       *lu = (Mat_UMFPACK *)A->data;
112ae26a541SJed Brown   Mat_SeqAIJ        *a  = (Mat_SeqAIJ *)lu->A->data;
113d9ca1df4SBarry Smith   PetscScalar       *av = a->a, *xa;
114d9ca1df4SBarry Smith   const PetscScalar *ba;
115ae26a541SJed Brown   PetscInt          *ai = a->i, *aj = a->j, status;
116fcf6e9abSJed Brown   static PetscBool   cite = PETSC_FALSE;
1171677d0b8SKris Buschelman 
1181677d0b8SKris Buschelman   PetscFunctionBegin;
1193ba16761SJacob Faibussowitsch   if (!A->rmap->n) PetscFunctionReturn(PETSC_SUCCESS);
1209371c9d4SSatish Balay   PetscCall(PetscCitationsRegister("@article{davis2004algorithm,\n  title={Algorithm 832: {UMFPACK} V4.3---An Unsymmetric-Pattern Multifrontal Method},\n  author={Davis, Timothy A},\n  journal={ACM Transactions on Mathematical Software (TOMS)},\n  "
1219371c9d4SSatish Balay                                    "volume={30},\n  number={2},\n  pages={196--199},\n  year={2004},\n  publisher={ACM}\n}\n",
1229371c9d4SSatish Balay                                    &cite));
1232d4e2982SHong Zhang   /* solve Ax = b by umfpack_*_wsolve */
1242d4e2982SHong Zhang 
125ae26a541SJed Brown   if (!lu->Wi) { /* first time, allocate working space for wsolve */
1269566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(A->rmap->n, &lu->Wi));
1279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5 * A->rmap->n, &lu->W));
128ae26a541SJed Brown   }
129ae26a541SJed Brown 
1309566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(b, &ba));
1319566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &xa));
13279c34000SHong Zhang #if defined(PETSC_USE_COMPLEX)
133b0043f70SBarry Smith   status = umfpack_UMF_wsolve(uflag, ai, aj, (PetscReal *)av, NULL, (PetscReal *)xa, NULL, (PetscReal *)ba, NULL, lu->Numeric, lu->Control, lu->Info, lu->Wi, lu->W);
13479c34000SHong Zhang #else
1354b019bd2SJed Brown   status = umfpack_UMF_wsolve(uflag, ai, aj, av, xa, ba, lu->Numeric, lu->Control, lu->Info, lu->Wi, lu->W);
13679c34000SHong Zhang #endif
1378592901bSBarry Smith   umfpack_UMF_report_info(lu->Control, lu->Info);
1381677d0b8SKris Buschelman   if (status < 0) {
1398592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
140e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "umfpack_UMF_wsolve failed");
1411677d0b8SKris Buschelman   }
1421677d0b8SKris Buschelman 
1439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(b, &ba));
1449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &xa));
1453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1464b019bd2SJed Brown }
1472a325a84SHong Zhang 
148d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_UMFPACK(Mat A, Vec b, Vec x)
149d71ae5a4SJacob Faibussowitsch {
1504b019bd2SJed Brown   PetscFunctionBegin;
1514b019bd2SJed Brown   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
1529566063dSJacob Faibussowitsch   PetscCall(MatSolve_UMFPACK_Private(A, b, x, UMFPACK_Aat));
1533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1544b019bd2SJed Brown }
1554b019bd2SJed Brown 
156d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolveTranspose_UMFPACK(Mat A, Vec b, Vec x)
157d71ae5a4SJacob Faibussowitsch {
1584b019bd2SJed Brown   PetscFunctionBegin;
1594b019bd2SJed Brown   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
1609566063dSJacob Faibussowitsch   PetscCall(MatSolve_UMFPACK_Private(A, b, x, UMFPACK_A));
1613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1621677d0b8SKris Buschelman }
1631677d0b8SKris Buschelman 
164d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat F, Mat A, const MatFactorInfo *info)
165d71ae5a4SJacob Faibussowitsch {
166b9c12af5SBarry Smith   Mat_UMFPACK *lu = (Mat_UMFPACK *)(F)->data;
167ae26a541SJed Brown   Mat_SeqAIJ  *a  = (Mat_SeqAIJ *)A->data;
168ae26a541SJed Brown   PetscInt    *ai = a->i, *aj = a->j, status;
169ae26a541SJed Brown   PetscScalar *av = a->a;
1701677d0b8SKris Buschelman 
1711677d0b8SKris Buschelman   PetscFunctionBegin;
1723ba16761SJacob Faibussowitsch   if (!A->rmap->n) PetscFunctionReturn(PETSC_SUCCESS);
1731677d0b8SKris Buschelman   /* numeric factorization of A' */
1742d4e2982SHong Zhang 
175ad540459SPierre Jolivet   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) umfpack_UMF_free_numeric(&lu->Numeric);
1762d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1778592901bSBarry Smith   status = umfpack_UMF_numeric(ai, aj, (double *)av, NULL, lu->Symbolic, &lu->Numeric, lu->Control, lu->Info);
1782d4e2982SHong Zhang #else
1798592901bSBarry Smith   status = umfpack_UMF_numeric(ai, aj, av, lu->Symbolic, &lu->Numeric, lu->Control, lu->Info);
1808592901bSBarry Smith #endif
1819f42a82aSMatthew Knepley   if (status < 0) {
1828592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
183e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "umfpack_UMF_numeric failed");
1849f42a82aSMatthew Knepley   }
1852d4e2982SHong Zhang   /* report numeric factorization of A' when Control[PRL] > 3 */
1868592901bSBarry Smith   (void)umfpack_UMF_report_numeric(lu->Numeric, lu->Control);
1871677d0b8SKris Buschelman 
1889566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)A));
1899566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lu->A));
1902205254eSKarl Rupp 
191ae26a541SJed Brown   lu->A                  = A;
1922a325a84SHong Zhang   lu->flg                = SAME_NONZERO_PATTERN;
1932a325a84SHong Zhang   lu->CleanUpUMFPACK     = PETSC_TRUE;
1944b019bd2SJed Brown   F->ops->solve          = MatSolve_UMFPACK;
1954b019bd2SJed Brown   F->ops->solvetranspose = MatSolveTranspose_UMFPACK;
1963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1971677d0b8SKris Buschelman }
1981677d0b8SKris Buschelman 
199d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info)
200d71ae5a4SJacob Faibussowitsch {
201ae26a541SJed Brown   Mat_SeqAIJ  *a  = (Mat_SeqAIJ *)A->data;
202b9c12af5SBarry Smith   Mat_UMFPACK *lu = (Mat_UMFPACK *)(F->data);
20326cc229bSBarry Smith   PetscInt     i, *ai = a->i, *aj = a->j, m = A->rmap->n, n = A->cmap->n, status, idx;
204989213a4SSatish Balay #if !defined(PETSC_USE_COMPLEX)
205ae26a541SJed Brown   PetscScalar *av = a->a;
206989213a4SSatish Balay #endif
2075d0c19d7SBarry Smith   const PetscInt *ra;
20826cc229bSBarry Smith   const char     *strategy[] = {"AUTO", "UNSYMMETRIC", "SYMMETRIC"};
20926cc229bSBarry Smith   const char     *scale[]    = {"NONE", "SUM", "MAX"};
21026cc229bSBarry Smith   PetscBool       flg;
2111677d0b8SKris Buschelman 
2121677d0b8SKris Buschelman   PetscFunctionBegin;
213a643c07aSBarry Smith   (F)->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
2143ba16761SJacob Faibussowitsch   if (!n) PetscFunctionReturn(PETSC_SUCCESS);
21526cc229bSBarry Smith 
21626cc229bSBarry Smith   /* Set options to F */
21726cc229bSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "UMFPACK Options", "Mat");
21826cc229bSBarry Smith   /* Control parameters used by reporting routiones */
21926cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_prl", "Control[UMFPACK_PRL]", "None", lu->Control[UMFPACK_PRL], &lu->Control[UMFPACK_PRL], NULL));
22026cc229bSBarry Smith 
22126cc229bSBarry Smith   /* Control parameters for symbolic factorization */
22226cc229bSBarry Smith   PetscCall(PetscOptionsEList("-mat_umfpack_strategy", "ordering and pivoting strategy", "None", strategy, 3, strategy[0], &idx, &flg));
22326cc229bSBarry Smith   if (flg) {
22426cc229bSBarry Smith     switch (idx) {
225d71ae5a4SJacob Faibussowitsch     case 0:
226d71ae5a4SJacob Faibussowitsch       lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO;
227d71ae5a4SJacob Faibussowitsch       break;
228d71ae5a4SJacob Faibussowitsch     case 1:
229d71ae5a4SJacob Faibussowitsch       lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC;
230d71ae5a4SJacob Faibussowitsch       break;
231d71ae5a4SJacob Faibussowitsch     case 2:
232d71ae5a4SJacob Faibussowitsch       lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC;
233d71ae5a4SJacob Faibussowitsch       break;
23426cc229bSBarry Smith     }
23526cc229bSBarry Smith   }
23626cc229bSBarry Smith   PetscCall(PetscOptionsEList("-mat_umfpack_ordering", "Internal ordering method", "None", UmfpackOrderingTypes, PETSC_STATIC_ARRAY_LENGTH(UmfpackOrderingTypes), UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]], &idx, &flg));
23726cc229bSBarry Smith   if (flg) lu->Control[UMFPACK_ORDERING] = (int)idx;
23826cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_dense_col", "Control[UMFPACK_DENSE_COL]", "None", lu->Control[UMFPACK_DENSE_COL], &lu->Control[UMFPACK_DENSE_COL], NULL));
23926cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_dense_row", "Control[UMFPACK_DENSE_ROW]", "None", lu->Control[UMFPACK_DENSE_ROW], &lu->Control[UMFPACK_DENSE_ROW], NULL));
24026cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_amd_dense", "Control[UMFPACK_AMD_DENSE]", "None", lu->Control[UMFPACK_AMD_DENSE], &lu->Control[UMFPACK_AMD_DENSE], NULL));
24126cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_block_size", "Control[UMFPACK_BLOCK_SIZE]", "None", lu->Control[UMFPACK_BLOCK_SIZE], &lu->Control[UMFPACK_BLOCK_SIZE], NULL));
24226cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_fixq", "Control[UMFPACK_FIXQ]", "None", lu->Control[UMFPACK_FIXQ], &lu->Control[UMFPACK_FIXQ], NULL));
24326cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_aggressive", "Control[UMFPACK_AGGRESSIVE]", "None", lu->Control[UMFPACK_AGGRESSIVE], &lu->Control[UMFPACK_AGGRESSIVE], NULL));
24426cc229bSBarry Smith 
24526cc229bSBarry Smith   /* Control parameters used by numeric factorization */
24626cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_pivot_tolerance", "Control[UMFPACK_PIVOT_TOLERANCE]", "None", lu->Control[UMFPACK_PIVOT_TOLERANCE], &lu->Control[UMFPACK_PIVOT_TOLERANCE], NULL));
24726cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_sym_pivot_tolerance", "Control[UMFPACK_SYM_PIVOT_TOLERANCE]", "None", lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE], &lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE], NULL));
24826cc229bSBarry Smith   PetscCall(PetscOptionsEList("-mat_umfpack_scale", "Control[UMFPACK_SCALE]", "None", scale, 3, scale[0], &idx, &flg));
24926cc229bSBarry Smith   if (flg) {
25026cc229bSBarry Smith     switch (idx) {
251d71ae5a4SJacob Faibussowitsch     case 0:
252d71ae5a4SJacob Faibussowitsch       lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE;
253d71ae5a4SJacob Faibussowitsch       break;
254d71ae5a4SJacob Faibussowitsch     case 1:
255d71ae5a4SJacob Faibussowitsch       lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM;
256d71ae5a4SJacob Faibussowitsch       break;
257d71ae5a4SJacob Faibussowitsch     case 2:
258d71ae5a4SJacob Faibussowitsch       lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX;
259d71ae5a4SJacob Faibussowitsch       break;
26026cc229bSBarry Smith     }
26126cc229bSBarry Smith   }
26226cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_alloc_init", "Control[UMFPACK_ALLOC_INIT]", "None", lu->Control[UMFPACK_ALLOC_INIT], &lu->Control[UMFPACK_ALLOC_INIT], NULL));
26326cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_front_alloc_init", "Control[UMFPACK_FRONT_ALLOC_INIT]", "None", lu->Control[UMFPACK_FRONT_ALLOC_INIT], &lu->Control[UMFPACK_ALLOC_INIT], NULL));
26426cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_droptol", "Control[UMFPACK_DROPTOL]", "None", lu->Control[UMFPACK_DROPTOL], &lu->Control[UMFPACK_DROPTOL], NULL));
26526cc229bSBarry Smith 
26626cc229bSBarry Smith   /* Control parameters used by solve */
26726cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_irstep", "Control[UMFPACK_IRSTEP]", "None", lu->Control[UMFPACK_IRSTEP], &lu->Control[UMFPACK_IRSTEP], NULL));
26826cc229bSBarry Smith   PetscOptionsEnd();
26926cc229bSBarry Smith 
2702c7c0729SBarry Smith   if (r) {
2719566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(r, &ra));
2729566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &lu->perm_c));
2737de69702SBarry Smith     /* we cannot simply memcpy on 64-bit archs */
2742d4e2982SHong Zhang     for (i = 0; i < m; i++) lu->perm_c[i] = ra[i];
2759566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(r, &ra));
2761677d0b8SKris Buschelman   }
2771677d0b8SKris Buschelman 
2781677d0b8SKris Buschelman   /* print the control parameters */
2798592901bSBarry Smith   if (lu->Control[UMFPACK_PRL] > 1) umfpack_UMF_report_control(lu->Control);
2801677d0b8SKris Buschelman 
2811677d0b8SKris Buschelman   /* symbolic factorization of A' */
2822c7c0729SBarry Smith   if (r) { /* use Petsc row ordering */
2838592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
284ae26a541SJed Brown     status = umfpack_UMF_qsymbolic(n, m, ai, aj, av, lu->perm_c, &lu->Symbolic, lu->Control, lu->Info);
2852d4e2982SHong Zhang #else
286b0043f70SBarry Smith     status = umfpack_UMF_qsymbolic(n, m, ai, aj, NULL, NULL, lu->perm_c, &lu->Symbolic, lu->Control, lu->Info);
2878592901bSBarry Smith #endif
2882d4e2982SHong Zhang   } else { /* use Umfpack col ordering */
2898592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
290ae26a541SJed Brown     status = umfpack_UMF_symbolic(n, m, ai, aj, av, &lu->Symbolic, lu->Control, lu->Info);
2918592901bSBarry Smith #else
292ae26a541SJed Brown     status = umfpack_UMF_symbolic(n, m, ai, aj, NULL, NULL, &lu->Symbolic, lu->Control, lu->Info);
2938592901bSBarry Smith #endif
2942d4e2982SHong Zhang   }
2952d4e2982SHong Zhang   if (status < 0) {
2968592901bSBarry Smith     umfpack_UMF_report_info(lu->Control, lu->Info);
2978592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
298e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "umfpack_UMF_symbolic failed");
2992d4e2982SHong Zhang   }
3002d4e2982SHong Zhang   /* report sumbolic factorization of A' when Control[PRL] > 3 */
3018592901bSBarry Smith   (void)umfpack_UMF_report_symbolic(lu->Symbolic, lu->Control);
3021677d0b8SKris Buschelman 
3031677d0b8SKris Buschelman   lu->flg            = DIFFERENT_NONZERO_PATTERN;
304e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK = PETSC_TRUE;
3053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3061677d0b8SKris Buschelman }
3071677d0b8SKris Buschelman 
308d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Info_UMFPACK(Mat A, PetscViewer viewer)
309d71ae5a4SJacob Faibussowitsch {
310b9c12af5SBarry Smith   Mat_UMFPACK *lu = (Mat_UMFPACK *)A->data;
311f0c56d0fSKris Buschelman 
3121677d0b8SKris Buschelman   PetscFunctionBegin;
3131677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
3143ba16761SJacob Faibussowitsch   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(PETSC_SUCCESS);
3151677d0b8SKris Buschelman 
3169566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "UMFPACK run parameters:\n"));
3171677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
3189566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_PRL]: %g\n", lu->Control[UMFPACK_PRL]));
3191677d0b8SKris Buschelman 
3201677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
3219566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_STRATEGY]: %g\n", lu->Control[UMFPACK_STRATEGY]));
3229566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_DENSE_COL]: %g\n", lu->Control[UMFPACK_DENSE_COL]));
3239566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_DENSE_ROW]: %g\n", lu->Control[UMFPACK_DENSE_ROW]));
3249566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_AMD_DENSE]: %g\n", lu->Control[UMFPACK_AMD_DENSE]));
3259566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_BLOCK_SIZE]: %g\n", lu->Control[UMFPACK_BLOCK_SIZE]));
3269566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_FIXQ]: %g\n", lu->Control[UMFPACK_FIXQ]));
3279566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_AGGRESSIVE]: %g\n", lu->Control[UMFPACK_AGGRESSIVE]));
3281677d0b8SKris Buschelman 
3291677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
3309566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n", lu->Control[UMFPACK_PIVOT_TOLERANCE]));
3319566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n", lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]));
3329566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_SCALE]: %g\n", lu->Control[UMFPACK_SCALE]));
3339566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_ALLOC_INIT]: %g\n", lu->Control[UMFPACK_ALLOC_INIT]));
3349566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_DROPTOL]: %g\n", lu->Control[UMFPACK_DROPTOL]));
3351677d0b8SKris Buschelman 
3361677d0b8SKris Buschelman   /* Control parameters used by solve */
3379566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_IRSTEP]: %g\n", lu->Control[UMFPACK_IRSTEP]));
3381677d0b8SKris Buschelman 
3391677d0b8SKris Buschelman   /* mat ordering */
34048a46eb9SPierre Jolivet   if (!lu->perm_c) PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_ORDERING]: %s (not using the PETSc ordering)\n", UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]]));
3413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3421677d0b8SKris Buschelman }
3431677d0b8SKris Buschelman 
344d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_UMFPACK(Mat A, PetscViewer viewer)
345d71ae5a4SJacob Faibussowitsch {
346ace3abfcSBarry Smith   PetscBool         iascii;
347d844382dSKris Buschelman   PetscViewerFormat format;
348d844382dSKris Buschelman 
349d844382dSKris Buschelman   PetscFunctionBegin;
3509566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
35132077d6dSBarry Smith   if (iascii) {
3529566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
35348a46eb9SPierre Jolivet     if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_UMFPACK(A, viewer));
354d844382dSKris Buschelman   }
3553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
356d844382dSKris Buschelman }
357d844382dSKris Buschelman 
35866976f2fSJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_seqaij_umfpack(Mat A, MatSolverType *type)
359d71ae5a4SJacob Faibussowitsch {
36035bd34faSBarry Smith   PetscFunctionBegin;
3612692d6eeSBarry Smith   *type = MATSOLVERUMFPACK;
3623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36335bd34faSBarry Smith }
36435bd34faSBarry Smith 
3652f71e704SKris Buschelman /*MC
36611a5261eSBarry Smith   MATSOLVERUMFPACK = "umfpack" - A matrix type providing direct solvers, LU, for sequential matrices
3672f71e704SKris Buschelman   via the external package UMFPACK.
3682f71e704SKris Buschelman 
3692ef1f0ffSBarry Smith   Use `./configure --download-suitesparse` to install PETSc to use UMFPACK
370c2b89b5dSBarry Smith 
3712ef1f0ffSBarry Smith   Use `-pc_type lu` `-pc_factor_mat_solver_type umfpack` to use this direct solver
3722f71e704SKris Buschelman 
3732f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
3742f71e704SKris Buschelman   which correspond to the options database keys below.
3752f71e704SKris Buschelman 
3762f71e704SKris Buschelman   Options Database Keys:
3772ef1f0ffSBarry Smith + -mat_umfpack_ordering                - `CHOLMOD`, `AMD`, `GIVEN`, `METIS`, `BEST`, `NONE`
378ba41fbb6SBarry Smith . -mat_umfpack_prl                     - UMFPACK print level: Control[UMFPACK_PRL]
3792ef1f0ffSBarry Smith . -mat_umfpack_strategy <AUTO>         - (choose one of) `AUTO`, `UNSYMMETRIC`, `SYMMETRIC 2BY2`
3802f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c>     - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
381e08999f5SMatthew G Knepley . -mat_umfpack_dense_row <0.2>         - Control[UMFPACK_DENSE_ROW]
382e08999f5SMatthew G Knepley . -mat_umfpack_amd_dense <10>          - Control[UMFPACK_AMD_DENSE]
3832f71e704SKris Buschelman . -mat_umfpack_block_size <bs>         - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
384e08999f5SMatthew G Knepley . -mat_umfpack_2by2_tolerance <0.01>   - Control[UMFPACK_2BY2_TOLERANCE]
385e08999f5SMatthew G Knepley . -mat_umfpack_fixq <0>                - Control[UMFPACK_FIXQ]
386e08999f5SMatthew G Knepley . -mat_umfpack_aggressive <1>          - Control[UMFPACK_AGGRESSIVE]
3872f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
388e08999f5SMatthew G Knepley . -mat_umfpack_sym_pivot_tolerance <0.001> - Control[UMFPACK_SYM_PIVOT_TOLERANCE]
389e08999f5SMatthew G Knepley . -mat_umfpack_scale <NONE>           - (choose one of) NONE SUM MAX
3902f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta>      - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
391e08999f5SMatthew G Knepley . -mat_umfpack_droptol <0>            - Control[UMFPACK_DROPTOL]
3922f71e704SKris Buschelman - -mat_umfpack_irstep <maxit>          - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
3932f71e704SKris Buschelman 
3942f71e704SKris Buschelman    Level: beginner
395a364b7d2SBarry Smith 
3962ef1f0ffSBarry Smith    Note:
397*1d27aa22SBarry Smith    UMFPACK is part of SuiteSparse <http://faculty.cse.tamu.edu/davis/suitesparse.html>
3982f71e704SKris Buschelman 
3991cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCLU`, `MATSOLVERSUPERLU`, `MATSOLVERMUMPS`, `PCFactorSetMatSolverType()`, `MatSolverType`
4002f71e704SKris Buschelman M*/
401b2573a8aSBarry Smith 
402d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_umfpack(Mat A, MatFactorType ftype, Mat *F)
403d71ae5a4SJacob Faibussowitsch {
4043519fcdcSHong Zhang   Mat          B;
4053519fcdcSHong Zhang   Mat_UMFPACK *lu;
40626cc229bSBarry Smith   PetscInt     m = A->rmap->n, n = A->cmap->n;
4072d4e2982SHong Zhang 
4083519fcdcSHong Zhang   PetscFunctionBegin;
4093519fcdcSHong Zhang   /* Create the factorization matrix F */
4109566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
4119566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n));
4129566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("umfpack", &((PetscObject)B)->type_name));
4139566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
414b9c12af5SBarry Smith 
4154dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&lu));
4162205254eSKarl Rupp 
417b9c12af5SBarry Smith   B->data                  = lu;
418b9c12af5SBarry Smith   B->ops->getinfo          = MatGetInfo_External;
4193519fcdcSHong Zhang   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
4203519fcdcSHong Zhang   B->ops->destroy          = MatDestroy_UMFPACK;
4213519fcdcSHong Zhang   B->ops->view             = MatView_UMFPACK;
4221aef8b4cSStefano Zampini   B->ops->matsolve         = NULL;
4232205254eSKarl Rupp 
4249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_umfpack));
4252205254eSKarl Rupp 
426d5f3da31SBarry Smith   B->factortype   = MAT_FACTOR_LU;
4273519fcdcSHong Zhang   B->assembled    = PETSC_TRUE; /* required by -ksp_view */
4283519fcdcSHong Zhang   B->preallocated = PETSC_TRUE;
4293519fcdcSHong Zhang 
4309566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
4319566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERUMFPACK, &B->solvertype));
432f73b0415SBarry Smith   B->canuseordering = PETSC_TRUE;
4339566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_LU]));
43400c67f3bSHong Zhang 
4353519fcdcSHong Zhang   /* initializations */
4363519fcdcSHong Zhang   /* get the default control parameters */
4378592901bSBarry Smith   umfpack_UMF_defaults(lu->Control);
438aaa8cc7dSPierre Jolivet   lu->perm_c                  = NULL; /* use default UMFPACK col permutation */
4393519fcdcSHong Zhang   lu->Control[UMFPACK_IRSTEP] = 0;    /* max num of iterative refinement steps to attempt */
4403519fcdcSHong Zhang 
4413519fcdcSHong Zhang   *F = B;
4423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4433519fcdcSHong Zhang }
444b2573a8aSBarry Smith 
445db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat, MatFactorType, Mat *);
446db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat, MatFactorType, Mat *);
447db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat, MatFactorType, Mat *);
448a2fc1e05SToby Isaac PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_spqr(Mat, MatFactorType, Mat *);
4492d4e2982SHong Zhang 
450d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_SuiteSparse(void)
451d71ae5a4SJacob Faibussowitsch {
45242c9c57cSBarry Smith   PetscFunctionBegin;
4539566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERUMFPACK, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_seqaij_umfpack));
4549566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERCHOLMOD, MATSEQAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_seqaij_cholmod));
4559566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERCHOLMOD, MATSEQSBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_seqsbaij_cholmod));
4569566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERKLU, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_seqaij_klu));
4579566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSPQR, MATSEQAIJ, MAT_FACTOR_QR, MatGetFactor_seqaij_spqr));
45848a46eb9SPierre Jolivet   if (!PetscDefined(USE_COMPLEX)) PetscCall(MatSolverTypeRegister(MATSOLVERSPQR, MATNORMAL, MAT_FACTOR_QR, MatGetFactor_seqaij_spqr));
4599566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSPQR, MATNORMALHERMITIAN, MAT_FACTOR_QR, MatGetFactor_seqaij_spqr));
4603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46142c9c57cSBarry Smith }
462