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; 115cf27e480SPierre Jolivet PetscInt *ai = a->i, *aj = a->j; 116cf27e480SPierre Jolivet int 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 { 16757508eceSPierre Jolivet Mat_UMFPACK *lu = (Mat_UMFPACK *)F->data; 168ae26a541SJed Brown Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 169cf27e480SPierre Jolivet PetscInt *ai = a->i, *aj = a->j; 170cf27e480SPierre Jolivet int status; 171ae26a541SJed Brown PetscScalar *av = a->a; 1721677d0b8SKris Buschelman 1731677d0b8SKris Buschelman PetscFunctionBegin; 1743ba16761SJacob Faibussowitsch if (!A->rmap->n) PetscFunctionReturn(PETSC_SUCCESS); 1751677d0b8SKris Buschelman /* numeric factorization of A' */ 1762d4e2982SHong Zhang 177ad540459SPierre Jolivet if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) umfpack_UMF_free_numeric(&lu->Numeric); 1782d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX) 1798592901bSBarry Smith status = umfpack_UMF_numeric(ai, aj, (double *)av, NULL, lu->Symbolic, &lu->Numeric, lu->Control, lu->Info); 1802d4e2982SHong Zhang #else 1818592901bSBarry Smith status = umfpack_UMF_numeric(ai, aj, av, lu->Symbolic, &lu->Numeric, lu->Control, lu->Info); 1828592901bSBarry Smith #endif 1839f42a82aSMatthew Knepley if (status < 0) { 1848592901bSBarry Smith umfpack_UMF_report_status(lu->Control, status); 185e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "umfpack_UMF_numeric failed"); 1869f42a82aSMatthew Knepley } 1872d4e2982SHong Zhang /* report numeric factorization of A' when Control[PRL] > 3 */ 1888592901bSBarry Smith (void)umfpack_UMF_report_numeric(lu->Numeric, lu->Control); 1891677d0b8SKris Buschelman 1909566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 1919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&lu->A)); 1922205254eSKarl Rupp 193ae26a541SJed Brown lu->A = A; 1942a325a84SHong Zhang lu->flg = SAME_NONZERO_PATTERN; 1952a325a84SHong Zhang lu->CleanUpUMFPACK = PETSC_TRUE; 1964b019bd2SJed Brown F->ops->solve = MatSolve_UMFPACK; 1974b019bd2SJed Brown F->ops->solvetranspose = MatSolveTranspose_UMFPACK; 1983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1991677d0b8SKris Buschelman } 2001677d0b8SKris Buschelman 201d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info) 202d71ae5a4SJacob Faibussowitsch { 203ae26a541SJed Brown Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 204f4f49eeaSPierre Jolivet Mat_UMFPACK *lu = (Mat_UMFPACK *)F->data; 205cf27e480SPierre Jolivet PetscInt i, *ai = a->i, *aj = a->j, m = A->rmap->n, n = A->cmap->n, idx; 206cf27e480SPierre Jolivet int status; 207989213a4SSatish Balay #if !defined(PETSC_USE_COMPLEX) 208ae26a541SJed Brown PetscScalar *av = a->a; 209989213a4SSatish Balay #endif 2105d0c19d7SBarry Smith const PetscInt *ra; 21126cc229bSBarry Smith const char *strategy[] = {"AUTO", "UNSYMMETRIC", "SYMMETRIC"}; 21226cc229bSBarry Smith const char *scale[] = {"NONE", "SUM", "MAX"}; 21326cc229bSBarry Smith PetscBool flg; 2141677d0b8SKris Buschelman 2151677d0b8SKris Buschelman PetscFunctionBegin; 21657508eceSPierre Jolivet F->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK; 2173ba16761SJacob Faibussowitsch if (!n) PetscFunctionReturn(PETSC_SUCCESS); 21826cc229bSBarry Smith 21926cc229bSBarry Smith /* Set options to F */ 22026cc229bSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "UMFPACK Options", "Mat"); 22126cc229bSBarry Smith /* Control parameters used by reporting routiones */ 22226cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_umfpack_prl", "Control[UMFPACK_PRL]", "None", lu->Control[UMFPACK_PRL], &lu->Control[UMFPACK_PRL], NULL)); 22326cc229bSBarry Smith 22426cc229bSBarry Smith /* Control parameters for symbolic factorization */ 22526cc229bSBarry Smith PetscCall(PetscOptionsEList("-mat_umfpack_strategy", "ordering and pivoting strategy", "None", strategy, 3, strategy[0], &idx, &flg)); 22626cc229bSBarry Smith if (flg) { 22726cc229bSBarry Smith switch (idx) { 228d71ae5a4SJacob Faibussowitsch case 0: 229d71ae5a4SJacob Faibussowitsch lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; 230d71ae5a4SJacob Faibussowitsch break; 231d71ae5a4SJacob Faibussowitsch case 1: 232d71ae5a4SJacob Faibussowitsch lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; 233d71ae5a4SJacob Faibussowitsch break; 234d71ae5a4SJacob Faibussowitsch case 2: 235d71ae5a4SJacob Faibussowitsch lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; 236d71ae5a4SJacob Faibussowitsch break; 23726cc229bSBarry Smith } 23826cc229bSBarry Smith } 23926cc229bSBarry Smith PetscCall(PetscOptionsEList("-mat_umfpack_ordering", "Internal ordering method", "None", UmfpackOrderingTypes, PETSC_STATIC_ARRAY_LENGTH(UmfpackOrderingTypes), UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]], &idx, &flg)); 24026cc229bSBarry Smith if (flg) lu->Control[UMFPACK_ORDERING] = (int)idx; 24126cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_umfpack_dense_col", "Control[UMFPACK_DENSE_COL]", "None", lu->Control[UMFPACK_DENSE_COL], &lu->Control[UMFPACK_DENSE_COL], NULL)); 24226cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_umfpack_dense_row", "Control[UMFPACK_DENSE_ROW]", "None", lu->Control[UMFPACK_DENSE_ROW], &lu->Control[UMFPACK_DENSE_ROW], NULL)); 24326cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_umfpack_amd_dense", "Control[UMFPACK_AMD_DENSE]", "None", lu->Control[UMFPACK_AMD_DENSE], &lu->Control[UMFPACK_AMD_DENSE], NULL)); 24426cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_umfpack_block_size", "Control[UMFPACK_BLOCK_SIZE]", "None", lu->Control[UMFPACK_BLOCK_SIZE], &lu->Control[UMFPACK_BLOCK_SIZE], NULL)); 24526cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_umfpack_fixq", "Control[UMFPACK_FIXQ]", "None", lu->Control[UMFPACK_FIXQ], &lu->Control[UMFPACK_FIXQ], NULL)); 24626cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_umfpack_aggressive", "Control[UMFPACK_AGGRESSIVE]", "None", lu->Control[UMFPACK_AGGRESSIVE], &lu->Control[UMFPACK_AGGRESSIVE], NULL)); 24726cc229bSBarry Smith 24826cc229bSBarry Smith /* Control parameters used by numeric factorization */ 24926cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_umfpack_pivot_tolerance", "Control[UMFPACK_PIVOT_TOLERANCE]", "None", lu->Control[UMFPACK_PIVOT_TOLERANCE], &lu->Control[UMFPACK_PIVOT_TOLERANCE], NULL)); 25026cc229bSBarry 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)); 25126cc229bSBarry Smith PetscCall(PetscOptionsEList("-mat_umfpack_scale", "Control[UMFPACK_SCALE]", "None", scale, 3, scale[0], &idx, &flg)); 25226cc229bSBarry Smith if (flg) { 25326cc229bSBarry Smith switch (idx) { 254d71ae5a4SJacob Faibussowitsch case 0: 255d71ae5a4SJacob Faibussowitsch lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; 256d71ae5a4SJacob Faibussowitsch break; 257d71ae5a4SJacob Faibussowitsch case 1: 258d71ae5a4SJacob Faibussowitsch lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; 259d71ae5a4SJacob Faibussowitsch break; 260d71ae5a4SJacob Faibussowitsch case 2: 261d71ae5a4SJacob Faibussowitsch lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; 262d71ae5a4SJacob Faibussowitsch break; 26326cc229bSBarry Smith } 26426cc229bSBarry Smith } 26526cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_umfpack_alloc_init", "Control[UMFPACK_ALLOC_INIT]", "None", lu->Control[UMFPACK_ALLOC_INIT], &lu->Control[UMFPACK_ALLOC_INIT], NULL)); 26626cc229bSBarry 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)); 26726cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_umfpack_droptol", "Control[UMFPACK_DROPTOL]", "None", lu->Control[UMFPACK_DROPTOL], &lu->Control[UMFPACK_DROPTOL], NULL)); 26826cc229bSBarry Smith 26926cc229bSBarry Smith /* Control parameters used by solve */ 27026cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_umfpack_irstep", "Control[UMFPACK_IRSTEP]", "None", lu->Control[UMFPACK_IRSTEP], &lu->Control[UMFPACK_IRSTEP], NULL)); 27126cc229bSBarry Smith PetscOptionsEnd(); 27226cc229bSBarry Smith 2732c7c0729SBarry Smith if (r) { 2749566063dSJacob Faibussowitsch PetscCall(ISGetIndices(r, &ra)); 2759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &lu->perm_c)); 2767de69702SBarry Smith /* we cannot simply memcpy on 64-bit archs */ 2772d4e2982SHong Zhang for (i = 0; i < m; i++) lu->perm_c[i] = ra[i]; 2789566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(r, &ra)); 2791677d0b8SKris Buschelman } 2801677d0b8SKris Buschelman 2811677d0b8SKris Buschelman /* print the control parameters */ 2828592901bSBarry Smith if (lu->Control[UMFPACK_PRL] > 1) umfpack_UMF_report_control(lu->Control); 2831677d0b8SKris Buschelman 2841677d0b8SKris Buschelman /* symbolic factorization of A' */ 285f0b74427SPierre Jolivet if (r) { /* use PETSc row ordering */ 2868592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX) 287ae26a541SJed Brown status = umfpack_UMF_qsymbolic(n, m, ai, aj, av, lu->perm_c, &lu->Symbolic, lu->Control, lu->Info); 2882d4e2982SHong Zhang #else 289b0043f70SBarry Smith status = umfpack_UMF_qsymbolic(n, m, ai, aj, NULL, NULL, lu->perm_c, &lu->Symbolic, lu->Control, lu->Info); 2908592901bSBarry Smith #endif 2912d4e2982SHong Zhang } else { /* use Umfpack col ordering */ 2928592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX) 293ae26a541SJed Brown status = umfpack_UMF_symbolic(n, m, ai, aj, av, &lu->Symbolic, lu->Control, lu->Info); 2948592901bSBarry Smith #else 295ae26a541SJed Brown status = umfpack_UMF_symbolic(n, m, ai, aj, NULL, NULL, &lu->Symbolic, lu->Control, lu->Info); 2968592901bSBarry Smith #endif 2972d4e2982SHong Zhang } 2982d4e2982SHong Zhang if (status < 0) { 2998592901bSBarry Smith umfpack_UMF_report_info(lu->Control, lu->Info); 3008592901bSBarry Smith umfpack_UMF_report_status(lu->Control, status); 301e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "umfpack_UMF_symbolic failed"); 3022d4e2982SHong Zhang } 3032d4e2982SHong Zhang /* report sumbolic factorization of A' when Control[PRL] > 3 */ 3048592901bSBarry Smith (void)umfpack_UMF_report_symbolic(lu->Symbolic, lu->Control); 3051677d0b8SKris Buschelman 3061677d0b8SKris Buschelman lu->flg = DIFFERENT_NONZERO_PATTERN; 307e1b4c3a1SKris Buschelman lu->CleanUpUMFPACK = PETSC_TRUE; 3083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3091677d0b8SKris Buschelman } 3101677d0b8SKris Buschelman 311d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Info_UMFPACK(Mat A, PetscViewer viewer) 312d71ae5a4SJacob Faibussowitsch { 313b9c12af5SBarry Smith Mat_UMFPACK *lu = (Mat_UMFPACK *)A->data; 314f0c56d0fSKris Buschelman 3151677d0b8SKris Buschelman PetscFunctionBegin; 3161677d0b8SKris Buschelman /* check if matrix is UMFPACK type */ 3173ba16761SJacob Faibussowitsch if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(PETSC_SUCCESS); 3181677d0b8SKris Buschelman 3199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "UMFPACK run parameters:\n")); 3201677d0b8SKris Buschelman /* Control parameters used by reporting routiones */ 3219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Control[UMFPACK_PRL]: %g\n", lu->Control[UMFPACK_PRL])); 3221677d0b8SKris Buschelman 3231677d0b8SKris Buschelman /* Control parameters used by symbolic factorization */ 3249566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Control[UMFPACK_STRATEGY]: %g\n", lu->Control[UMFPACK_STRATEGY])); 3259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Control[UMFPACK_DENSE_COL]: %g\n", lu->Control[UMFPACK_DENSE_COL])); 3269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Control[UMFPACK_DENSE_ROW]: %g\n", lu->Control[UMFPACK_DENSE_ROW])); 3279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Control[UMFPACK_AMD_DENSE]: %g\n", lu->Control[UMFPACK_AMD_DENSE])); 3289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Control[UMFPACK_BLOCK_SIZE]: %g\n", lu->Control[UMFPACK_BLOCK_SIZE])); 3299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Control[UMFPACK_FIXQ]: %g\n", lu->Control[UMFPACK_FIXQ])); 3309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Control[UMFPACK_AGGRESSIVE]: %g\n", lu->Control[UMFPACK_AGGRESSIVE])); 3311677d0b8SKris Buschelman 3321677d0b8SKris Buschelman /* Control parameters used by numeric factorization */ 3339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Control[UMFPACK_PIVOT_TOLERANCE]: %g\n", lu->Control[UMFPACK_PIVOT_TOLERANCE])); 3349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n", lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE])); 3359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Control[UMFPACK_SCALE]: %g\n", lu->Control[UMFPACK_SCALE])); 3369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Control[UMFPACK_ALLOC_INIT]: %g\n", lu->Control[UMFPACK_ALLOC_INIT])); 3379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Control[UMFPACK_DROPTOL]: %g\n", lu->Control[UMFPACK_DROPTOL])); 3381677d0b8SKris Buschelman 3391677d0b8SKris Buschelman /* Control parameters used by solve */ 3409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Control[UMFPACK_IRSTEP]: %g\n", lu->Control[UMFPACK_IRSTEP])); 3411677d0b8SKris Buschelman 3421677d0b8SKris Buschelman /* mat ordering */ 34348a46eb9SPierre Jolivet if (!lu->perm_c) PetscCall(PetscViewerASCIIPrintf(viewer, " Control[UMFPACK_ORDERING]: %s (not using the PETSc ordering)\n", UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]])); 3443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3451677d0b8SKris Buschelman } 3461677d0b8SKris Buschelman 347d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_UMFPACK(Mat A, PetscViewer viewer) 348d71ae5a4SJacob Faibussowitsch { 349*9f196a02SMartin Diehl PetscBool isascii; 350d844382dSKris Buschelman PetscViewerFormat format; 351d844382dSKris Buschelman 352d844382dSKris Buschelman PetscFunctionBegin; 353*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 354*9f196a02SMartin Diehl if (isascii) { 3559566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 35648a46eb9SPierre Jolivet if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_UMFPACK(A, viewer)); 357d844382dSKris Buschelman } 3583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 359d844382dSKris Buschelman } 360d844382dSKris Buschelman 36166976f2fSJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_seqaij_umfpack(Mat A, MatSolverType *type) 362d71ae5a4SJacob Faibussowitsch { 36335bd34faSBarry Smith PetscFunctionBegin; 3642692d6eeSBarry Smith *type = MATSOLVERUMFPACK; 3653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 36635bd34faSBarry Smith } 36735bd34faSBarry Smith 3682f71e704SKris Buschelman /*MC 36911a5261eSBarry Smith MATSOLVERUMFPACK = "umfpack" - A matrix type providing direct solvers, LU, for sequential matrices 3702f71e704SKris Buschelman via the external package UMFPACK. 3712f71e704SKris Buschelman 3722ef1f0ffSBarry Smith Use `./configure --download-suitesparse` to install PETSc to use UMFPACK 373c2b89b5dSBarry Smith 3742ef1f0ffSBarry Smith Use `-pc_type lu` `-pc_factor_mat_solver_type umfpack` to use this direct solver 3752f71e704SKris Buschelman 3762f71e704SKris Buschelman Consult UMFPACK documentation for more information about the Control parameters 3772f71e704SKris Buschelman which correspond to the options database keys below. 3782f71e704SKris Buschelman 3792f71e704SKris Buschelman Options Database Keys: 3802ef1f0ffSBarry Smith + -mat_umfpack_ordering - `CHOLMOD`, `AMD`, `GIVEN`, `METIS`, `BEST`, `NONE` 381ba41fbb6SBarry Smith . -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL] 3822ef1f0ffSBarry Smith . -mat_umfpack_strategy <AUTO> - (choose one of) `AUTO`, `UNSYMMETRIC`, `SYMMETRIC 2BY2` 3832f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL] 384e08999f5SMatthew G Knepley . -mat_umfpack_dense_row <0.2> - Control[UMFPACK_DENSE_ROW] 385e08999f5SMatthew G Knepley . -mat_umfpack_amd_dense <10> - Control[UMFPACK_AMD_DENSE] 3862f71e704SKris Buschelman . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE] 387e08999f5SMatthew G Knepley . -mat_umfpack_2by2_tolerance <0.01> - Control[UMFPACK_2BY2_TOLERANCE] 388e08999f5SMatthew G Knepley . -mat_umfpack_fixq <0> - Control[UMFPACK_FIXQ] 389e08999f5SMatthew G Knepley . -mat_umfpack_aggressive <1> - Control[UMFPACK_AGGRESSIVE] 3902f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE] 391e08999f5SMatthew G Knepley . -mat_umfpack_sym_pivot_tolerance <0.001> - Control[UMFPACK_SYM_PIVOT_TOLERANCE] 392e08999f5SMatthew G Knepley . -mat_umfpack_scale <NONE> - (choose one of) NONE SUM MAX 3932f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT] 394e08999f5SMatthew G Knepley . -mat_umfpack_droptol <0> - Control[UMFPACK_DROPTOL] 3952f71e704SKris Buschelman - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP] 3962f71e704SKris Buschelman 3972f71e704SKris Buschelman Level: beginner 398a364b7d2SBarry Smith 3992ef1f0ffSBarry Smith Note: 4001d27aa22SBarry Smith UMFPACK is part of SuiteSparse <http://faculty.cse.tamu.edu/davis/suitesparse.html> 4012f71e704SKris Buschelman 4021cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCLU`, `MATSOLVERSUPERLU`, `MATSOLVERMUMPS`, `PCFactorSetMatSolverType()`, `MatSolverType` 4032f71e704SKris Buschelman M*/ 404b2573a8aSBarry Smith 405d1f0640dSPierre Jolivet static PetscErrorCode MatGetFactor_seqaij_umfpack(Mat A, MatFactorType ftype, Mat *F) 406d71ae5a4SJacob Faibussowitsch { 4073519fcdcSHong Zhang Mat B; 4083519fcdcSHong Zhang Mat_UMFPACK *lu; 40926cc229bSBarry Smith PetscInt m = A->rmap->n, n = A->cmap->n; 4102d4e2982SHong Zhang 4113519fcdcSHong Zhang PetscFunctionBegin; 4123519fcdcSHong Zhang /* Create the factorization matrix F */ 4139566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 4149566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n)); 4159566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("umfpack", &((PetscObject)B)->type_name)); 4169566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 417b9c12af5SBarry Smith 4184dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&lu)); 4192205254eSKarl Rupp 420b9c12af5SBarry Smith B->data = lu; 421b9c12af5SBarry Smith B->ops->getinfo = MatGetInfo_External; 4223519fcdcSHong Zhang B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK; 4233519fcdcSHong Zhang B->ops->destroy = MatDestroy_UMFPACK; 4243519fcdcSHong Zhang B->ops->view = MatView_UMFPACK; 4251aef8b4cSStefano Zampini B->ops->matsolve = NULL; 4262205254eSKarl Rupp 4279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_umfpack)); 4282205254eSKarl Rupp 429d5f3da31SBarry Smith B->factortype = MAT_FACTOR_LU; 4303519fcdcSHong Zhang B->assembled = PETSC_TRUE; /* required by -ksp_view */ 4313519fcdcSHong Zhang B->preallocated = PETSC_TRUE; 4323519fcdcSHong Zhang 4339566063dSJacob Faibussowitsch PetscCall(PetscFree(B->solvertype)); 4349566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERUMFPACK, &B->solvertype)); 435f73b0415SBarry Smith B->canuseordering = PETSC_TRUE; 4369566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_LU])); 43700c67f3bSHong Zhang 4383519fcdcSHong Zhang /* initializations */ 4393519fcdcSHong Zhang /* get the default control parameters */ 4408592901bSBarry Smith umfpack_UMF_defaults(lu->Control); 441aaa8cc7dSPierre Jolivet lu->perm_c = NULL; /* use default UMFPACK col permutation */ 4423519fcdcSHong Zhang lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */ 4433519fcdcSHong Zhang 4443519fcdcSHong Zhang *F = B; 4453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4463519fcdcSHong Zhang } 447b2573a8aSBarry Smith 448db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat, MatFactorType, Mat *); 449db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat, MatFactorType, Mat *); 450db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat, MatFactorType, Mat *); 451a2fc1e05SToby Isaac PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_spqr(Mat, MatFactorType, Mat *); 4522d4e2982SHong Zhang 453d1f0640dSPierre Jolivet PETSC_INTERN PetscErrorCode MatSolverTypeRegister_SuiteSparse(void) 454d71ae5a4SJacob Faibussowitsch { 45542c9c57cSBarry Smith PetscFunctionBegin; 4569566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERUMFPACK, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_seqaij_umfpack)); 4579566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERCHOLMOD, MATSEQAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_seqaij_cholmod)); 4589566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERCHOLMOD, MATSEQSBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_seqsbaij_cholmod)); 4599566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERKLU, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_seqaij_klu)); 4609566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSPQR, MATSEQAIJ, MAT_FACTOR_QR, MatGetFactor_seqaij_spqr)); 46148a46eb9SPierre Jolivet if (!PetscDefined(USE_COMPLEX)) PetscCall(MatSolverTypeRegister(MATSOLVERSPQR, MATNORMAL, MAT_FACTOR_QR, MatGetFactor_seqaij_spqr)); 4629566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSPQR, MATNORMALHERMITIAN, MAT_FACTOR_QR, MatGetFactor_seqaij_spqr)); 4633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46442c9c57cSBarry Smith } 465