xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
11677d0b8SKris Buschelman 
21677d0b8SKris Buschelman /*
3d515b9b4SShri Abhyankar    Provides an interface to the UMFPACK sparse solver available through SuiteSparse version 4.2.1
42d4e2982SHong Zhang 
59e475b0dSSatish Balay    When build with PETSC_USE_64BIT_INDICES this will use Suitesparse_long as the
68592901bSBarry Smith    integer type in UMFPACK, otherwise it will use int. This means
78592901bSBarry Smith    all integers in this file as simply declared as PetscInt. Also it means
87de69702SBarry Smith    that one cannot use 64BIT_INDICES on 32-bit pointer systems [as Suitesparse_long is 32-bit only]
92d4e2982SHong Zhang 
101677d0b8SKris Buschelman */
11c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
121677d0b8SKris Buschelman 
138592901bSBarry Smith #if defined(PETSC_USE_64BIT_INDICES)
148592901bSBarry Smith   #if defined(PETSC_USE_COMPLEX)
158592901bSBarry Smith     #define umfpack_UMF_free_symbolic umfpack_zl_free_symbolic
168592901bSBarry Smith     #define umfpack_UMF_free_numeric  umfpack_zl_free_numeric
17d755ee67SBarry Smith     /* the type casts are needed because PetscInt is long long while SuiteSparse_long is long and compilers warn even when they are identical */
18d755ee67SBarry 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)
19d755ee67SBarry 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)
208592901bSBarry Smith     #define umfpack_UMF_report_numeric                                   umfpack_zl_report_numeric
218592901bSBarry Smith     #define umfpack_UMF_report_control                                   umfpack_zl_report_control
228592901bSBarry Smith     #define umfpack_UMF_report_status                                    umfpack_zl_report_status
238592901bSBarry Smith     #define umfpack_UMF_report_info                                      umfpack_zl_report_info
248592901bSBarry Smith     #define umfpack_UMF_report_symbolic                                  umfpack_zl_report_symbolic
25d755ee67SBarry 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)
26d755ee67SBarry 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)
278592901bSBarry Smith     #define umfpack_UMF_defaults                                         umfpack_zl_defaults
288592901bSBarry Smith 
298592901bSBarry Smith   #else
308592901bSBarry Smith     #define umfpack_UMF_free_symbolic                           umfpack_dl_free_symbolic
318592901bSBarry Smith     #define umfpack_UMF_free_numeric                            umfpack_dl_free_numeric
32d755ee67SBarry 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)
33d755ee67SBarry 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)
348592901bSBarry Smith     #define umfpack_UMF_report_numeric                          umfpack_dl_report_numeric
358592901bSBarry Smith     #define umfpack_UMF_report_control                          umfpack_dl_report_control
368592901bSBarry Smith     #define umfpack_UMF_report_status                           umfpack_dl_report_status
378592901bSBarry Smith     #define umfpack_UMF_report_info                             umfpack_dl_report_info
388592901bSBarry Smith     #define umfpack_UMF_report_symbolic                         umfpack_dl_report_symbolic
39d755ee67SBarry 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)
40d755ee67SBarry 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)
418592901bSBarry Smith     #define umfpack_UMF_defaults                                umfpack_dl_defaults
428592901bSBarry Smith   #endif
438592901bSBarry Smith 
448592901bSBarry Smith #else
458592901bSBarry Smith   #if defined(PETSC_USE_COMPLEX)
468592901bSBarry Smith     #define umfpack_UMF_free_symbolic   umfpack_zi_free_symbolic
478592901bSBarry Smith     #define umfpack_UMF_free_numeric    umfpack_zi_free_numeric
488592901bSBarry Smith     #define umfpack_UMF_wsolve          umfpack_zi_wsolve
498592901bSBarry Smith     #define umfpack_UMF_numeric         umfpack_zi_numeric
508592901bSBarry Smith     #define umfpack_UMF_report_numeric  umfpack_zi_report_numeric
518592901bSBarry Smith     #define umfpack_UMF_report_control  umfpack_zi_report_control
528592901bSBarry Smith     #define umfpack_UMF_report_status   umfpack_zi_report_status
538592901bSBarry Smith     #define umfpack_UMF_report_info     umfpack_zi_report_info
548592901bSBarry Smith     #define umfpack_UMF_report_symbolic umfpack_zi_report_symbolic
558592901bSBarry Smith     #define umfpack_UMF_qsymbolic       umfpack_zi_qsymbolic
568592901bSBarry Smith     #define umfpack_UMF_symbolic        umfpack_zi_symbolic
578592901bSBarry Smith     #define umfpack_UMF_defaults        umfpack_zi_defaults
588592901bSBarry Smith 
598592901bSBarry Smith   #else
608592901bSBarry Smith     #define umfpack_UMF_free_symbolic   umfpack_di_free_symbolic
618592901bSBarry Smith     #define umfpack_UMF_free_numeric    umfpack_di_free_numeric
628592901bSBarry Smith     #define umfpack_UMF_wsolve          umfpack_di_wsolve
638592901bSBarry Smith     #define umfpack_UMF_numeric         umfpack_di_numeric
648592901bSBarry Smith     #define umfpack_UMF_report_numeric  umfpack_di_report_numeric
658592901bSBarry Smith     #define umfpack_UMF_report_control  umfpack_di_report_control
668592901bSBarry Smith     #define umfpack_UMF_report_status   umfpack_di_report_status
678592901bSBarry Smith     #define umfpack_UMF_report_info     umfpack_di_report_info
688592901bSBarry Smith     #define umfpack_UMF_report_symbolic umfpack_di_report_symbolic
698592901bSBarry Smith     #define umfpack_UMF_qsymbolic       umfpack_di_qsymbolic
708592901bSBarry Smith     #define umfpack_UMF_symbolic        umfpack_di_symbolic
718592901bSBarry Smith     #define umfpack_UMF_defaults        umfpack_di_defaults
728592901bSBarry Smith   #endif
738592901bSBarry Smith #endif
748592901bSBarry Smith 
751677d0b8SKris Buschelman EXTERN_C_BEGIN
76c6db04a5SJed Brown #include <umfpack.h>
771677d0b8SKris Buschelman EXTERN_C_END
781677d0b8SKris Buschelman 
79fcfd50ebSBarry Smith static const char *const UmfpackOrderingTypes[] = {"CHOLMOD", "AMD", "GIVEN", "METIS", "BEST", "NONE", "USER", "UmfpackOrderingTypes", "UMFPACK_ORDERING_", 0};
80fcfd50ebSBarry Smith 
811677d0b8SKris Buschelman typedef struct {
821677d0b8SKris Buschelman   void        *Symbolic, *Numeric;
831677d0b8SKris Buschelman   double       Info[UMFPACK_INFO], Control[UMFPACK_CONTROL], *W;
84ae26a541SJed Brown   PetscInt    *Wi, *perm_c;
85ae26a541SJed Brown   Mat          A; /* Matrix used for factorization */
861677d0b8SKris Buschelman   MatStructure flg;
871677d0b8SKris Buschelman 
881677d0b8SKris Buschelman   /* Flag to clean up UMFPACK objects during Destroy */
89ace3abfcSBarry Smith   PetscBool CleanUpUMFPACK;
90f0c56d0fSKris Buschelman } Mat_UMFPACK;
91f0c56d0fSKris Buschelman 
92d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_UMFPACK(Mat A)
93d71ae5a4SJacob Faibussowitsch {
94b9c12af5SBarry Smith   Mat_UMFPACK *lu = (Mat_UMFPACK *)A->data;
951677d0b8SKris Buschelman 
961677d0b8SKris Buschelman   PetscFunctionBegin;
97b9c12af5SBarry Smith   if (lu->CleanUpUMFPACK) {
988592901bSBarry Smith     umfpack_UMF_free_symbolic(&lu->Symbolic);
998592901bSBarry Smith     umfpack_UMF_free_numeric(&lu->Numeric);
1009566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->Wi));
1019566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->W));
1029566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->perm_c));
1031677d0b8SKris Buschelman   }
1049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lu->A));
1052e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL));
1069566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
1073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1081677d0b8SKris Buschelman }
1091677d0b8SKris Buschelman 
110d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_UMFPACK_Private(Mat A, Vec b, Vec x, int uflag)
111d71ae5a4SJacob Faibussowitsch {
112b9c12af5SBarry Smith   Mat_UMFPACK       *lu = (Mat_UMFPACK *)A->data;
113ae26a541SJed Brown   Mat_SeqAIJ        *a  = (Mat_SeqAIJ *)lu->A->data;
114d9ca1df4SBarry Smith   PetscScalar       *av = a->a, *xa;
115d9ca1df4SBarry Smith   const PetscScalar *ba;
116ae26a541SJed Brown   PetscInt          *ai = a->i, *aj = a->j, status;
117fcf6e9abSJed Brown   static PetscBool   cite = PETSC_FALSE;
1181677d0b8SKris Buschelman 
1191677d0b8SKris Buschelman   PetscFunctionBegin;
1203ba16761SJacob Faibussowitsch   if (!A->rmap->n) PetscFunctionReturn(PETSC_SUCCESS);
1219371c9d4SSatish 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  "
1229371c9d4SSatish Balay                                    "volume={30},\n  number={2},\n  pages={196--199},\n  year={2004},\n  publisher={ACM}\n}\n",
1239371c9d4SSatish Balay                                    &cite));
1242d4e2982SHong Zhang   /* solve Ax = b by umfpack_*_wsolve */
1252d4e2982SHong Zhang 
126ae26a541SJed Brown   if (!lu->Wi) { /* first time, allocate working space for wsolve */
1279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(A->rmap->n, &lu->Wi));
1289566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5 * A->rmap->n, &lu->W));
129ae26a541SJed Brown   }
130ae26a541SJed Brown 
1319566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(b, &ba));
1329566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &xa));
13379c34000SHong Zhang #if defined(PETSC_USE_COMPLEX)
134b0043f70SBarry 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);
13579c34000SHong Zhang #else
1364b019bd2SJed Brown   status = umfpack_UMF_wsolve(uflag, ai, aj, av, xa, ba, lu->Numeric, lu->Control, lu->Info, lu->Wi, lu->W);
13779c34000SHong Zhang #endif
1388592901bSBarry Smith   umfpack_UMF_report_info(lu->Control, lu->Info);
1391677d0b8SKris Buschelman   if (status < 0) {
1408592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
141e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "umfpack_UMF_wsolve failed");
1421677d0b8SKris Buschelman   }
1431677d0b8SKris Buschelman 
1449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(b, &ba));
1459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &xa));
1463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1474b019bd2SJed Brown }
1482a325a84SHong Zhang 
149d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_UMFPACK(Mat A, Vec b, Vec x)
150d71ae5a4SJacob Faibussowitsch {
1514b019bd2SJed Brown   PetscFunctionBegin;
1524b019bd2SJed Brown   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
1539566063dSJacob Faibussowitsch   PetscCall(MatSolve_UMFPACK_Private(A, b, x, UMFPACK_Aat));
1543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1554b019bd2SJed Brown }
1564b019bd2SJed Brown 
157d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolveTranspose_UMFPACK(Mat A, Vec b, Vec x)
158d71ae5a4SJacob Faibussowitsch {
1594b019bd2SJed Brown   PetscFunctionBegin;
1604b019bd2SJed Brown   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
1619566063dSJacob Faibussowitsch   PetscCall(MatSolve_UMFPACK_Private(A, b, x, UMFPACK_A));
1623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1631677d0b8SKris Buschelman }
1641677d0b8SKris Buschelman 
165d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat F, Mat A, const MatFactorInfo *info)
166d71ae5a4SJacob Faibussowitsch {
167b9c12af5SBarry Smith   Mat_UMFPACK *lu = (Mat_UMFPACK *)(F)->data;
168ae26a541SJed Brown   Mat_SeqAIJ  *a  = (Mat_SeqAIJ *)A->data;
169ae26a541SJed Brown   PetscInt    *ai = a->i, *aj = a->j, status;
170ae26a541SJed Brown   PetscScalar *av = a->a;
1711677d0b8SKris Buschelman 
1721677d0b8SKris Buschelman   PetscFunctionBegin;
1733ba16761SJacob Faibussowitsch   if (!A->rmap->n) PetscFunctionReturn(PETSC_SUCCESS);
1741677d0b8SKris Buschelman   /* numeric factorization of A' */
1752d4e2982SHong Zhang 
176ad540459SPierre Jolivet   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) umfpack_UMF_free_numeric(&lu->Numeric);
1772d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1788592901bSBarry Smith   status = umfpack_UMF_numeric(ai, aj, (double *)av, NULL, lu->Symbolic, &lu->Numeric, lu->Control, lu->Info);
1792d4e2982SHong Zhang #else
1808592901bSBarry Smith   status = umfpack_UMF_numeric(ai, aj, av, lu->Symbolic, &lu->Numeric, lu->Control, lu->Info);
1818592901bSBarry Smith #endif
1829f42a82aSMatthew Knepley   if (status < 0) {
1838592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
184e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "umfpack_UMF_numeric failed");
1859f42a82aSMatthew Knepley   }
1862d4e2982SHong Zhang   /* report numeric factorization of A' when Control[PRL] > 3 */
1878592901bSBarry Smith   (void)umfpack_UMF_report_numeric(lu->Numeric, lu->Control);
1881677d0b8SKris Buschelman 
1899566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)A));
1909566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lu->A));
1912205254eSKarl Rupp 
192ae26a541SJed Brown   lu->A                  = A;
1932a325a84SHong Zhang   lu->flg                = SAME_NONZERO_PATTERN;
1942a325a84SHong Zhang   lu->CleanUpUMFPACK     = PETSC_TRUE;
1954b019bd2SJed Brown   F->ops->solve          = MatSolve_UMFPACK;
1964b019bd2SJed Brown   F->ops->solvetranspose = MatSolveTranspose_UMFPACK;
1973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1981677d0b8SKris Buschelman }
1991677d0b8SKris Buschelman 
200d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info)
201d71ae5a4SJacob Faibussowitsch {
202ae26a541SJed Brown   Mat_SeqAIJ  *a  = (Mat_SeqAIJ *)A->data;
203b9c12af5SBarry Smith   Mat_UMFPACK *lu = (Mat_UMFPACK *)(F->data);
20426cc229bSBarry Smith   PetscInt     i, *ai = a->i, *aj = a->j, m = A->rmap->n, n = A->cmap->n, status, idx;
205989213a4SSatish Balay #if !defined(PETSC_USE_COMPLEX)
206ae26a541SJed Brown   PetscScalar *av = a->a;
207989213a4SSatish Balay #endif
2085d0c19d7SBarry Smith   const PetscInt *ra;
20926cc229bSBarry Smith   const char     *strategy[] = {"AUTO", "UNSYMMETRIC", "SYMMETRIC"};
21026cc229bSBarry Smith   const char     *scale[]    = {"NONE", "SUM", "MAX"};
21126cc229bSBarry Smith   PetscBool       flg;
2121677d0b8SKris Buschelman 
2131677d0b8SKris Buschelman   PetscFunctionBegin;
214a643c07aSBarry Smith   (F)->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
2153ba16761SJacob Faibussowitsch   if (!n) PetscFunctionReturn(PETSC_SUCCESS);
21626cc229bSBarry Smith 
21726cc229bSBarry Smith   /* Set options to F */
21826cc229bSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "UMFPACK Options", "Mat");
21926cc229bSBarry Smith   /* Control parameters used by reporting routiones */
22026cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_prl", "Control[UMFPACK_PRL]", "None", lu->Control[UMFPACK_PRL], &lu->Control[UMFPACK_PRL], NULL));
22126cc229bSBarry Smith 
22226cc229bSBarry Smith   /* Control parameters for symbolic factorization */
22326cc229bSBarry Smith   PetscCall(PetscOptionsEList("-mat_umfpack_strategy", "ordering and pivoting strategy", "None", strategy, 3, strategy[0], &idx, &flg));
22426cc229bSBarry Smith   if (flg) {
22526cc229bSBarry Smith     switch (idx) {
226d71ae5a4SJacob Faibussowitsch     case 0:
227d71ae5a4SJacob Faibussowitsch       lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO;
228d71ae5a4SJacob Faibussowitsch       break;
229d71ae5a4SJacob Faibussowitsch     case 1:
230d71ae5a4SJacob Faibussowitsch       lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC;
231d71ae5a4SJacob Faibussowitsch       break;
232d71ae5a4SJacob Faibussowitsch     case 2:
233d71ae5a4SJacob Faibussowitsch       lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC;
234d71ae5a4SJacob Faibussowitsch       break;
23526cc229bSBarry Smith     }
23626cc229bSBarry Smith   }
23726cc229bSBarry Smith   PetscCall(PetscOptionsEList("-mat_umfpack_ordering", "Internal ordering method", "None", UmfpackOrderingTypes, PETSC_STATIC_ARRAY_LENGTH(UmfpackOrderingTypes), UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]], &idx, &flg));
23826cc229bSBarry Smith   if (flg) lu->Control[UMFPACK_ORDERING] = (int)idx;
23926cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_dense_col", "Control[UMFPACK_DENSE_COL]", "None", lu->Control[UMFPACK_DENSE_COL], &lu->Control[UMFPACK_DENSE_COL], NULL));
24026cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_dense_row", "Control[UMFPACK_DENSE_ROW]", "None", lu->Control[UMFPACK_DENSE_ROW], &lu->Control[UMFPACK_DENSE_ROW], NULL));
24126cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_amd_dense", "Control[UMFPACK_AMD_DENSE]", "None", lu->Control[UMFPACK_AMD_DENSE], &lu->Control[UMFPACK_AMD_DENSE], NULL));
24226cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_block_size", "Control[UMFPACK_BLOCK_SIZE]", "None", lu->Control[UMFPACK_BLOCK_SIZE], &lu->Control[UMFPACK_BLOCK_SIZE], NULL));
24326cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_fixq", "Control[UMFPACK_FIXQ]", "None", lu->Control[UMFPACK_FIXQ], &lu->Control[UMFPACK_FIXQ], NULL));
24426cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_aggressive", "Control[UMFPACK_AGGRESSIVE]", "None", lu->Control[UMFPACK_AGGRESSIVE], &lu->Control[UMFPACK_AGGRESSIVE], NULL));
24526cc229bSBarry Smith 
24626cc229bSBarry Smith   /* Control parameters used by numeric factorization */
24726cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_pivot_tolerance", "Control[UMFPACK_PIVOT_TOLERANCE]", "None", lu->Control[UMFPACK_PIVOT_TOLERANCE], &lu->Control[UMFPACK_PIVOT_TOLERANCE], NULL));
24826cc229bSBarry 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));
24926cc229bSBarry Smith   PetscCall(PetscOptionsEList("-mat_umfpack_scale", "Control[UMFPACK_SCALE]", "None", scale, 3, scale[0], &idx, &flg));
25026cc229bSBarry Smith   if (flg) {
25126cc229bSBarry Smith     switch (idx) {
252d71ae5a4SJacob Faibussowitsch     case 0:
253d71ae5a4SJacob Faibussowitsch       lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE;
254d71ae5a4SJacob Faibussowitsch       break;
255d71ae5a4SJacob Faibussowitsch     case 1:
256d71ae5a4SJacob Faibussowitsch       lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM;
257d71ae5a4SJacob Faibussowitsch       break;
258d71ae5a4SJacob Faibussowitsch     case 2:
259d71ae5a4SJacob Faibussowitsch       lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX;
260d71ae5a4SJacob Faibussowitsch       break;
26126cc229bSBarry Smith     }
26226cc229bSBarry Smith   }
26326cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_alloc_init", "Control[UMFPACK_ALLOC_INIT]", "None", lu->Control[UMFPACK_ALLOC_INIT], &lu->Control[UMFPACK_ALLOC_INIT], NULL));
26426cc229bSBarry 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));
26526cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_droptol", "Control[UMFPACK_DROPTOL]", "None", lu->Control[UMFPACK_DROPTOL], &lu->Control[UMFPACK_DROPTOL], NULL));
26626cc229bSBarry Smith 
26726cc229bSBarry Smith   /* Control parameters used by solve */
26826cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_irstep", "Control[UMFPACK_IRSTEP]", "None", lu->Control[UMFPACK_IRSTEP], &lu->Control[UMFPACK_IRSTEP], NULL));
26926cc229bSBarry Smith   PetscOptionsEnd();
27026cc229bSBarry Smith 
2712c7c0729SBarry Smith   if (r) {
2729566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(r, &ra));
2739566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &lu->perm_c));
2747de69702SBarry Smith     /* we cannot simply memcpy on 64-bit archs */
2752d4e2982SHong Zhang     for (i = 0; i < m; i++) lu->perm_c[i] = ra[i];
2769566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(r, &ra));
2771677d0b8SKris Buschelman   }
2781677d0b8SKris Buschelman 
2791677d0b8SKris Buschelman   /* print the control parameters */
2808592901bSBarry Smith   if (lu->Control[UMFPACK_PRL] > 1) umfpack_UMF_report_control(lu->Control);
2811677d0b8SKris Buschelman 
2821677d0b8SKris Buschelman   /* symbolic factorization of A' */
2832c7c0729SBarry Smith   if (r) { /* use Petsc row ordering */
2848592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
285ae26a541SJed Brown     status = umfpack_UMF_qsymbolic(n, m, ai, aj, av, lu->perm_c, &lu->Symbolic, lu->Control, lu->Info);
2862d4e2982SHong Zhang #else
287b0043f70SBarry Smith     status = umfpack_UMF_qsymbolic(n, m, ai, aj, NULL, NULL, lu->perm_c, &lu->Symbolic, lu->Control, lu->Info);
2888592901bSBarry Smith #endif
2892d4e2982SHong Zhang   } else { /* use Umfpack col ordering */
2908592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
291ae26a541SJed Brown     status = umfpack_UMF_symbolic(n, m, ai, aj, av, &lu->Symbolic, lu->Control, lu->Info);
2928592901bSBarry Smith #else
293ae26a541SJed Brown     status = umfpack_UMF_symbolic(n, m, ai, aj, NULL, NULL, &lu->Symbolic, lu->Control, lu->Info);
2948592901bSBarry Smith #endif
2952d4e2982SHong Zhang   }
2962d4e2982SHong Zhang   if (status < 0) {
2978592901bSBarry Smith     umfpack_UMF_report_info(lu->Control, lu->Info);
2988592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
299e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "umfpack_UMF_symbolic failed");
3002d4e2982SHong Zhang   }
3012d4e2982SHong Zhang   /* report sumbolic factorization of A' when Control[PRL] > 3 */
3028592901bSBarry Smith   (void)umfpack_UMF_report_symbolic(lu->Symbolic, lu->Control);
3031677d0b8SKris Buschelman 
3041677d0b8SKris Buschelman   lu->flg            = DIFFERENT_NONZERO_PATTERN;
305e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK = PETSC_TRUE;
3063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3071677d0b8SKris Buschelman }
3081677d0b8SKris Buschelman 
309d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Info_UMFPACK(Mat A, PetscViewer viewer)
310d71ae5a4SJacob Faibussowitsch {
311b9c12af5SBarry Smith   Mat_UMFPACK *lu = (Mat_UMFPACK *)A->data;
312f0c56d0fSKris Buschelman 
3131677d0b8SKris Buschelman   PetscFunctionBegin;
3141677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
3153ba16761SJacob Faibussowitsch   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(PETSC_SUCCESS);
3161677d0b8SKris Buschelman 
3179566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "UMFPACK run parameters:\n"));
3181677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
3199566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_PRL]: %g\n", lu->Control[UMFPACK_PRL]));
3201677d0b8SKris Buschelman 
3211677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
3229566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_STRATEGY]: %g\n", lu->Control[UMFPACK_STRATEGY]));
3239566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_DENSE_COL]: %g\n", lu->Control[UMFPACK_DENSE_COL]));
3249566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_DENSE_ROW]: %g\n", lu->Control[UMFPACK_DENSE_ROW]));
3259566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_AMD_DENSE]: %g\n", lu->Control[UMFPACK_AMD_DENSE]));
3269566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_BLOCK_SIZE]: %g\n", lu->Control[UMFPACK_BLOCK_SIZE]));
3279566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_FIXQ]: %g\n", lu->Control[UMFPACK_FIXQ]));
3289566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_AGGRESSIVE]: %g\n", lu->Control[UMFPACK_AGGRESSIVE]));
3291677d0b8SKris Buschelman 
3301677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
3319566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n", lu->Control[UMFPACK_PIVOT_TOLERANCE]));
3329566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n", lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]));
3339566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_SCALE]: %g\n", lu->Control[UMFPACK_SCALE]));
3349566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_ALLOC_INIT]: %g\n", lu->Control[UMFPACK_ALLOC_INIT]));
3359566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_DROPTOL]: %g\n", lu->Control[UMFPACK_DROPTOL]));
3361677d0b8SKris Buschelman 
3371677d0b8SKris Buschelman   /* Control parameters used by solve */
3389566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_IRSTEP]: %g\n", lu->Control[UMFPACK_IRSTEP]));
3391677d0b8SKris Buschelman 
3401677d0b8SKris Buschelman   /* mat ordering */
34148a46eb9SPierre Jolivet   if (!lu->perm_c) PetscCall(PetscViewerASCIIPrintf(viewer, "  Control[UMFPACK_ORDERING]: %s (not using the PETSc ordering)\n", UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]]));
3423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3431677d0b8SKris Buschelman }
3441677d0b8SKris Buschelman 
345d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_UMFPACK(Mat A, PetscViewer viewer)
346d71ae5a4SJacob Faibussowitsch {
347ace3abfcSBarry Smith   PetscBool         iascii;
348d844382dSKris Buschelman   PetscViewerFormat format;
349d844382dSKris Buschelman 
350d844382dSKris Buschelman   PetscFunctionBegin;
3519566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
35232077d6dSBarry Smith   if (iascii) {
3539566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
35448a46eb9SPierre Jolivet     if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_UMFPACK(A, viewer));
355d844382dSKris Buschelman   }
3563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
357d844382dSKris Buschelman }
358d844382dSKris Buschelman 
359d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFactorGetSolverType_seqaij_umfpack(Mat A, MatSolverType *type)
360d71ae5a4SJacob Faibussowitsch {
36135bd34faSBarry Smith   PetscFunctionBegin;
3622692d6eeSBarry Smith   *type = MATSOLVERUMFPACK;
3633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36435bd34faSBarry Smith }
36535bd34faSBarry Smith 
3662f71e704SKris Buschelman /*MC
36711a5261eSBarry Smith   MATSOLVERUMFPACK = "umfpack" - A matrix type providing direct solvers, LU, for sequential matrices
3682f71e704SKris Buschelman   via the external package UMFPACK.
3692f71e704SKris Buschelman 
3702ef1f0ffSBarry Smith   Use `./configure --download-suitesparse` to install PETSc to use UMFPACK
371c2b89b5dSBarry Smith 
3722ef1f0ffSBarry Smith   Use `-pc_type lu` `-pc_factor_mat_solver_type umfpack` to use this direct solver
3732f71e704SKris Buschelman 
3742f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
3752f71e704SKris Buschelman   which correspond to the options database keys below.
3762f71e704SKris Buschelman 
3772f71e704SKris Buschelman   Options Database Keys:
3782ef1f0ffSBarry Smith + -mat_umfpack_ordering                - `CHOLMOD`, `AMD`, `GIVEN`, `METIS`, `BEST`, `NONE`
379ba41fbb6SBarry Smith . -mat_umfpack_prl                     - UMFPACK print level: Control[UMFPACK_PRL]
3802ef1f0ffSBarry Smith . -mat_umfpack_strategy <AUTO>         - (choose one of) `AUTO`, `UNSYMMETRIC`, `SYMMETRIC 2BY2`
3812f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c>     - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
382e08999f5SMatthew G Knepley . -mat_umfpack_dense_row <0.2>         - Control[UMFPACK_DENSE_ROW]
383e08999f5SMatthew G Knepley . -mat_umfpack_amd_dense <10>          - Control[UMFPACK_AMD_DENSE]
3842f71e704SKris Buschelman . -mat_umfpack_block_size <bs>         - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
385e08999f5SMatthew G Knepley . -mat_umfpack_2by2_tolerance <0.01>   - Control[UMFPACK_2BY2_TOLERANCE]
386e08999f5SMatthew G Knepley . -mat_umfpack_fixq <0>                - Control[UMFPACK_FIXQ]
387e08999f5SMatthew G Knepley . -mat_umfpack_aggressive <1>          - Control[UMFPACK_AGGRESSIVE]
3882f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
389e08999f5SMatthew G Knepley . -mat_umfpack_sym_pivot_tolerance <0.001> - Control[UMFPACK_SYM_PIVOT_TOLERANCE]
390e08999f5SMatthew G Knepley . -mat_umfpack_scale <NONE>           - (choose one of) NONE SUM MAX
3912f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta>      - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
392e08999f5SMatthew G Knepley . -mat_umfpack_droptol <0>            - Control[UMFPACK_DROPTOL]
3932f71e704SKris Buschelman - -mat_umfpack_irstep <maxit>          - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
3942f71e704SKris Buschelman 
3952f71e704SKris Buschelman    Level: beginner
396a364b7d2SBarry Smith 
3972ef1f0ffSBarry Smith    Note:
3982ef1f0ffSBarry Smith    UMFPACK is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html
3992f71e704SKris Buschelman 
400*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCLU`, `MATSOLVERSUPERLU`, `MATSOLVERMUMPS`, `PCFactorSetMatSolverType()`, `MatSolverType`
4012f71e704SKris Buschelman M*/
402b2573a8aSBarry Smith 
403d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_umfpack(Mat A, MatFactorType ftype, Mat *F)
404d71ae5a4SJacob Faibussowitsch {
4053519fcdcSHong Zhang   Mat          B;
4063519fcdcSHong Zhang   Mat_UMFPACK *lu;
40726cc229bSBarry Smith   PetscInt     m = A->rmap->n, n = A->cmap->n;
4082d4e2982SHong Zhang 
4093519fcdcSHong Zhang   PetscFunctionBegin;
4103519fcdcSHong Zhang   /* Create the factorization matrix F */
4119566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
4129566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n));
4139566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("umfpack", &((PetscObject)B)->type_name));
4149566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
415b9c12af5SBarry Smith 
4164dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&lu));
4172205254eSKarl Rupp 
418b9c12af5SBarry Smith   B->data                  = lu;
419b9c12af5SBarry Smith   B->ops->getinfo          = MatGetInfo_External;
4203519fcdcSHong Zhang   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
4213519fcdcSHong Zhang   B->ops->destroy          = MatDestroy_UMFPACK;
4223519fcdcSHong Zhang   B->ops->view             = MatView_UMFPACK;
4231aef8b4cSStefano Zampini   B->ops->matsolve         = NULL;
4242205254eSKarl Rupp 
4259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_umfpack));
4262205254eSKarl Rupp 
427d5f3da31SBarry Smith   B->factortype   = MAT_FACTOR_LU;
4283519fcdcSHong Zhang   B->assembled    = PETSC_TRUE; /* required by -ksp_view */
4293519fcdcSHong Zhang   B->preallocated = PETSC_TRUE;
4303519fcdcSHong Zhang 
4319566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
4329566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERUMFPACK, &B->solvertype));
433f73b0415SBarry Smith   B->canuseordering = PETSC_TRUE;
4349566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_LU]));
43500c67f3bSHong Zhang 
4363519fcdcSHong Zhang   /* initializations */
4373519fcdcSHong Zhang   /* get the default control parameters */
4388592901bSBarry Smith   umfpack_UMF_defaults(lu->Control);
439aaa8cc7dSPierre Jolivet   lu->perm_c                  = NULL; /* use default UMFPACK col permutation */
4403519fcdcSHong Zhang   lu->Control[UMFPACK_IRSTEP] = 0;    /* max num of iterative refinement steps to attempt */
4413519fcdcSHong Zhang 
4423519fcdcSHong Zhang   *F = B;
4433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4443519fcdcSHong Zhang }
445b2573a8aSBarry Smith 
446db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat, MatFactorType, Mat *);
447db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat, MatFactorType, Mat *);
448db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat, MatFactorType, Mat *);
449a2fc1e05SToby Isaac PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_spqr(Mat, MatFactorType, Mat *);
4502d4e2982SHong Zhang 
451d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_SuiteSparse(void)
452d71ae5a4SJacob Faibussowitsch {
45342c9c57cSBarry Smith   PetscFunctionBegin;
4549566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERUMFPACK, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_seqaij_umfpack));
4559566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERCHOLMOD, MATSEQAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_seqaij_cholmod));
4569566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERCHOLMOD, MATSEQSBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_seqsbaij_cholmod));
4579566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERKLU, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_seqaij_klu));
4589566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSPQR, MATSEQAIJ, MAT_FACTOR_QR, MatGetFactor_seqaij_spqr));
45948a46eb9SPierre Jolivet   if (!PetscDefined(USE_COMPLEX)) PetscCall(MatSolverTypeRegister(MATSOLVERSPQR, MATNORMAL, MAT_FACTOR_QR, MatGetFactor_seqaij_spqr));
4609566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSPQR, MATNORMALHERMITIAN, MAT_FACTOR_QR, MatGetFactor_seqaij_spqr));
4613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46242c9c57cSBarry Smith }
463