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 89e475b0dSSatish Balay that one cannot use 64BIT_INDICES on 32bit machines [as Suitesparse_long is 32bit 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; 87fcfd50ebSBarry Smith PetscBool PetscMatOrdering; 881677d0b8SKris Buschelman 891677d0b8SKris Buschelman /* Flag to clean up UMFPACK objects during Destroy */ 90ace3abfcSBarry Smith PetscBool CleanUpUMFPACK; 91f0c56d0fSKris Buschelman } Mat_UMFPACK; 92f0c56d0fSKris Buschelman 934b019bd2SJed Brown static PetscErrorCode MatDestroy_UMFPACK(Mat A) 94dfbe8321SBarry Smith { 95dfbe8321SBarry Smith PetscErrorCode ierr; 96b9c12af5SBarry Smith Mat_UMFPACK *lu=(Mat_UMFPACK*)A->data; 971677d0b8SKris Buschelman 981677d0b8SKris Buschelman PetscFunctionBegin; 99b9c12af5SBarry Smith if (lu->CleanUpUMFPACK) { 1008592901bSBarry Smith umfpack_UMF_free_symbolic(&lu->Symbolic); 1018592901bSBarry Smith umfpack_UMF_free_numeric(&lu->Numeric); 1021677d0b8SKris Buschelman ierr = PetscFree(lu->Wi);CHKERRQ(ierr); 1031677d0b8SKris Buschelman ierr = PetscFree(lu->W);CHKERRQ(ierr); 1041677d0b8SKris Buschelman ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 1051677d0b8SKris Buschelman } 106ae26a541SJed Brown ierr = MatDestroy(&lu->A);CHKERRQ(ierr); 107b9c12af5SBarry Smith ierr = PetscFree(A->data);CHKERRQ(ierr); 1081677d0b8SKris Buschelman PetscFunctionReturn(0); 1091677d0b8SKris Buschelman } 1101677d0b8SKris Buschelman 1114b019bd2SJed Brown static PetscErrorCode MatSolve_UMFPACK_Private(Mat A,Vec b,Vec x,int uflag) 1126849ba73SBarry Smith { 113b9c12af5SBarry Smith Mat_UMFPACK *lu = (Mat_UMFPACK*)A->data; 114ae26a541SJed Brown Mat_SeqAIJ *a = (Mat_SeqAIJ*)lu->A->data; 115d9ca1df4SBarry Smith PetscScalar *av = a->a,*xa; 116d9ca1df4SBarry Smith const PetscScalar *ba; 117dfbe8321SBarry Smith PetscErrorCode ierr; 118ae26a541SJed Brown PetscInt *ai = a->i,*aj = a->j,status; 119fcf6e9abSJed Brown static PetscBool cite = PETSC_FALSE; 1201677d0b8SKris Buschelman 1211677d0b8SKris Buschelman PetscFunctionBegin; 122a643c07aSBarry Smith if (!A->rmap->n) PetscFunctionReturn(0); 123fcf6e9abSJed Brown ierr = 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 volume={30},\n number={2},\n pages={196--199},\n year={2004},\n publisher={ACM}\n}\n",&cite);CHKERRQ(ierr); 1242d4e2982SHong Zhang /* solve Ax = b by umfpack_*_wsolve */ 1251677d0b8SKris Buschelman /* ----------------------------------*/ 1262d4e2982SHong Zhang 127ae26a541SJed Brown if (!lu->Wi) { /* first time, allocate working space for wsolve */ 128785e854fSJed Brown ierr = PetscMalloc1(A->rmap->n,&lu->Wi);CHKERRQ(ierr); 129785e854fSJed Brown ierr = PetscMalloc1(5*A->rmap->n,&lu->W);CHKERRQ(ierr); 130ae26a541SJed Brown } 131ae26a541SJed Brown 132d9ca1df4SBarry Smith ierr = VecGetArrayRead(b,&ba); 1331677d0b8SKris Buschelman ierr = VecGetArray(x,&xa); 13479c34000SHong Zhang #if defined(PETSC_USE_COMPLEX) 135b0043f70SBarry 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); 13679c34000SHong Zhang #else 1374b019bd2SJed Brown status = umfpack_UMF_wsolve(uflag,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W); 13879c34000SHong Zhang #endif 1398592901bSBarry Smith umfpack_UMF_report_info(lu->Control, lu->Info); 1401677d0b8SKris Buschelman if (status < 0) { 1418592901bSBarry Smith umfpack_UMF_report_status(lu->Control, status); 142e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_wsolve failed"); 1431677d0b8SKris Buschelman } 1441677d0b8SKris Buschelman 145d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(b,&ba);CHKERRQ(ierr); 146057dcbc9SJed Brown ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr); 1474b019bd2SJed Brown PetscFunctionReturn(0); 1484b019bd2SJed Brown } 1492a325a84SHong Zhang 1504b019bd2SJed Brown static PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x) 1514b019bd2SJed Brown { 1524b019bd2SJed Brown PetscErrorCode ierr; 1534b019bd2SJed Brown 1544b019bd2SJed Brown PetscFunctionBegin; 1554b019bd2SJed Brown /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */ 1564b019bd2SJed Brown ierr = MatSolve_UMFPACK_Private(A,b,x,UMFPACK_Aat);CHKERRQ(ierr); 1574b019bd2SJed Brown PetscFunctionReturn(0); 1584b019bd2SJed Brown } 1594b019bd2SJed Brown 1604b019bd2SJed Brown static PetscErrorCode MatSolveTranspose_UMFPACK(Mat A,Vec b,Vec x) 1614b019bd2SJed Brown { 1624b019bd2SJed Brown PetscErrorCode ierr; 1634b019bd2SJed Brown 1644b019bd2SJed Brown PetscFunctionBegin; 1654b019bd2SJed Brown /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */ 1664b019bd2SJed Brown ierr = MatSolve_UMFPACK_Private(A,b,x,UMFPACK_A);CHKERRQ(ierr); 1671677d0b8SKris Buschelman PetscFunctionReturn(0); 1681677d0b8SKris Buschelman } 1691677d0b8SKris Buschelman 1704b019bd2SJed Brown static PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat F,Mat A,const MatFactorInfo *info) 1716849ba73SBarry Smith { 172b9c12af5SBarry Smith Mat_UMFPACK *lu = (Mat_UMFPACK*)(F)->data; 173ae26a541SJed Brown Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 174ae26a541SJed Brown PetscInt *ai = a->i,*aj=a->j,status; 175ae26a541SJed Brown PetscScalar *av = a->a; 1766849ba73SBarry Smith PetscErrorCode ierr; 1771677d0b8SKris Buschelman 1781677d0b8SKris Buschelman PetscFunctionBegin; 179a643c07aSBarry Smith if (!A->rmap->n) PetscFunctionReturn(0); 1801677d0b8SKris Buschelman /* numeric factorization of A' */ 1811677d0b8SKris Buschelman /* ----------------------------*/ 1822d4e2982SHong Zhang 1838592901bSBarry Smith if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) { 1848592901bSBarry Smith umfpack_UMF_free_numeric(&lu->Numeric); 1858592901bSBarry Smith } 1862d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX) 1878592901bSBarry Smith status = umfpack_UMF_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info); 1882d4e2982SHong Zhang #else 1898592901bSBarry Smith status = umfpack_UMF_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info); 1908592901bSBarry Smith #endif 1919f42a82aSMatthew Knepley if (status < 0) { 1928592901bSBarry Smith umfpack_UMF_report_status(lu->Control, status); 193e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_numeric failed"); 1949f42a82aSMatthew Knepley } 1952d4e2982SHong Zhang /* report numeric factorization of A' when Control[PRL] > 3 */ 1968592901bSBarry Smith (void) umfpack_UMF_report_numeric(lu->Numeric, lu->Control); 1971677d0b8SKris Buschelman 198ae26a541SJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 199ae26a541SJed Brown ierr = MatDestroy(&lu->A);CHKERRQ(ierr); 2002205254eSKarl Rupp 201ae26a541SJed Brown lu->A = A; 2022a325a84SHong Zhang lu->flg = SAME_NONZERO_PATTERN; 2032a325a84SHong Zhang lu->CleanUpUMFPACK = PETSC_TRUE; 2044b019bd2SJed Brown F->ops->solve = MatSolve_UMFPACK; 2054b019bd2SJed Brown F->ops->solvetranspose = MatSolveTranspose_UMFPACK; 2061677d0b8SKris Buschelman PetscFunctionReturn(0); 2071677d0b8SKris Buschelman } 2081677d0b8SKris Buschelman 2091677d0b8SKris Buschelman /* 2101677d0b8SKris Buschelman Note the r permutation is ignored 2111677d0b8SKris Buschelman */ 2124b019bd2SJed Brown static PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 2136849ba73SBarry Smith { 214ae26a541SJed Brown Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 215b9c12af5SBarry Smith Mat_UMFPACK *lu = (Mat_UMFPACK*)(F->data); 216dfbe8321SBarry Smith PetscErrorCode ierr; 217ae26a541SJed Brown PetscInt i,*ai = a->i,*aj = a->j,m=A->rmap->n,n=A->cmap->n; 218989213a4SSatish Balay #if !defined(PETSC_USE_COMPLEX) 219ae26a541SJed Brown PetscScalar *av = a->a; 220989213a4SSatish Balay #endif 2215d0c19d7SBarry Smith const PetscInt *ra; 2228592901bSBarry Smith PetscInt status; 2231677d0b8SKris Buschelman 2241677d0b8SKris Buschelman PetscFunctionBegin; 225a643c07aSBarry Smith (F)->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK; 226a643c07aSBarry Smith if (!n) PetscFunctionReturn(0); 227fcfd50ebSBarry Smith if (lu->PetscMatOrdering) { 2280f6efc10SHong Zhang ierr = ISGetIndices(r,&ra);CHKERRQ(ierr); 229785e854fSJed Brown ierr = PetscMalloc1(m,&lu->perm_c);CHKERRQ(ierr); 2302d4e2982SHong Zhang /* we cannot simply memcpy on 64 bit archs */ 2312d4e2982SHong Zhang for (i = 0; i < m; i++) lu->perm_c[i] = ra[i]; 2320f6efc10SHong Zhang ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr); 2331677d0b8SKris Buschelman } 2341677d0b8SKris Buschelman 2351677d0b8SKris Buschelman /* print the control parameters */ 2368592901bSBarry Smith if (lu->Control[UMFPACK_PRL] > 1) umfpack_UMF_report_control(lu->Control); 2371677d0b8SKris Buschelman 2381677d0b8SKris Buschelman /* symbolic factorization of A' */ 2391677d0b8SKris Buschelman /* ---------------------------------------------------------------------- */ 240fcfd50ebSBarry Smith if (lu->PetscMatOrdering) { /* use Petsc row ordering */ 2418592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX) 242ae26a541SJed Brown status = umfpack_UMF_qsymbolic(n,m,ai,aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info); 2432d4e2982SHong Zhang #else 244b0043f70SBarry Smith status = umfpack_UMF_qsymbolic(n,m,ai,aj,NULL,NULL,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info); 2458592901bSBarry Smith #endif 2462d4e2982SHong Zhang } else { /* use Umfpack col ordering */ 2478592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX) 248ae26a541SJed Brown status = umfpack_UMF_symbolic(n,m,ai,aj,av,&lu->Symbolic,lu->Control,lu->Info); 2498592901bSBarry Smith #else 250ae26a541SJed Brown status = umfpack_UMF_symbolic(n,m,ai,aj,NULL,NULL,&lu->Symbolic,lu->Control,lu->Info); 2518592901bSBarry Smith #endif 2522d4e2982SHong Zhang } 2532d4e2982SHong Zhang if (status < 0) { 2548592901bSBarry Smith umfpack_UMF_report_info(lu->Control, lu->Info); 2558592901bSBarry Smith umfpack_UMF_report_status(lu->Control, status); 256e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_symbolic failed"); 2572d4e2982SHong Zhang } 2582d4e2982SHong Zhang /* report sumbolic factorization of A' when Control[PRL] > 3 */ 2598592901bSBarry Smith (void) umfpack_UMF_report_symbolic(lu->Symbolic, lu->Control); 2601677d0b8SKris Buschelman 2611677d0b8SKris Buschelman lu->flg = DIFFERENT_NONZERO_PATTERN; 262e1b4c3a1SKris Buschelman lu->CleanUpUMFPACK = PETSC_TRUE; 2631677d0b8SKris Buschelman PetscFunctionReturn(0); 2641677d0b8SKris Buschelman } 2651677d0b8SKris Buschelman 266*860c79edSBarry Smith static PetscErrorCode MatView_Info_UMFPACK(Mat A,PetscViewer viewer) 2676849ba73SBarry Smith { 268b9c12af5SBarry Smith Mat_UMFPACK *lu= (Mat_UMFPACK*)A->data; 269dfbe8321SBarry Smith PetscErrorCode ierr; 270f0c56d0fSKris Buschelman 2711677d0b8SKris Buschelman PetscFunctionBegin; 2721677d0b8SKris Buschelman /* check if matrix is UMFPACK type */ 273f0c56d0fSKris Buschelman if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0); 2741677d0b8SKris Buschelman 2751677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr); 2761677d0b8SKris Buschelman /* Control parameters used by reporting routiones */ 2771677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr); 2781677d0b8SKris Buschelman 2791677d0b8SKris Buschelman /* Control parameters used by symbolic factorization */ 2800f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr); 2811677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr); 2821677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr); 2830f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr); 2841677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr); 2850f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr); 2860f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr); 2871677d0b8SKris Buschelman 2881677d0b8SKris Buschelman /* Control parameters used by numeric factorization */ 2891677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr); 2900f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr); 2910f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr); 2921677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr); 2930f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr); 2941677d0b8SKris Buschelman 2951677d0b8SKris Buschelman /* Control parameters used by solve */ 2961677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr); 2971677d0b8SKris Buschelman 2981677d0b8SKris Buschelman /* mat ordering */ 299fcfd50ebSBarry Smith if (!lu->PetscMatOrdering) { 300fcfd50ebSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_ORDERING]: %s (not using the PETSc ordering)\n",UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]]);CHKERRQ(ierr); 301fcfd50ebSBarry Smith } 3021677d0b8SKris Buschelman PetscFunctionReturn(0); 3031677d0b8SKris Buschelman } 3041677d0b8SKris Buschelman 3054b019bd2SJed Brown static PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer) 3066849ba73SBarry Smith { 307dfbe8321SBarry Smith PetscErrorCode ierr; 308ace3abfcSBarry Smith PetscBool iascii; 309d844382dSKris Buschelman PetscViewerFormat format; 310d844382dSKris Buschelman 311d844382dSKris Buschelman PetscFunctionBegin; 312251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 31332077d6dSBarry Smith if (iascii) { 314d844382dSKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 3152f59403fSHong Zhang if (format == PETSC_VIEWER_ASCII_INFO) { 316*860c79edSBarry Smith ierr = MatView_Info_UMFPACK(A,viewer);CHKERRQ(ierr); 317d844382dSKris Buschelman } 318d844382dSKris Buschelman } 319d844382dSKris Buschelman PetscFunctionReturn(0); 320d844382dSKris Buschelman } 321d844382dSKris Buschelman 32235bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_seqaij_umfpack(Mat A,const MatSolverPackage *type) 32335bd34faSBarry Smith { 32435bd34faSBarry Smith PetscFunctionBegin; 3252692d6eeSBarry Smith *type = MATSOLVERUMFPACK; 32635bd34faSBarry Smith PetscFunctionReturn(0); 32735bd34faSBarry Smith } 32835bd34faSBarry Smith 32935bd34faSBarry Smith 3302f71e704SKris Buschelman /*MC 3312692d6eeSBarry Smith MATSOLVERUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices 3322f71e704SKris Buschelman via the external package UMFPACK. 3332f71e704SKris Buschelman 334c2b89b5dSBarry Smith Use ./configure --download-suitesparse to install PETSc to use UMFPACK 335c2b89b5dSBarry Smith 336c2b89b5dSBarry Smith Use -pc_type lu -pc_factor_mat_solver_package umfpack to us this direct solver 3372f71e704SKris Buschelman 3382f71e704SKris Buschelman Consult UMFPACK documentation for more information about the Control parameters 3392f71e704SKris Buschelman which correspond to the options database keys below. 3402f71e704SKris Buschelman 3412f71e704SKris Buschelman Options Database Keys: 342ba41fbb6SBarry Smith + -mat_umfpack_ordering - CHOLMOD, AMD, GIVEN, METIS, BEST, NONE 343ba41fbb6SBarry Smith . -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL] 344e08999f5SMatthew G Knepley . -mat_umfpack_strategy <AUTO> - (choose one of) AUTO UNSYMMETRIC SYMMETRIC 2BY2 3452f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL] 346e08999f5SMatthew G Knepley . -mat_umfpack_dense_row <0.2> - Control[UMFPACK_DENSE_ROW] 347e08999f5SMatthew G Knepley . -mat_umfpack_amd_dense <10> - Control[UMFPACK_AMD_DENSE] 3482f71e704SKris Buschelman . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE] 349e08999f5SMatthew G Knepley . -mat_umfpack_2by2_tolerance <0.01> - Control[UMFPACK_2BY2_TOLERANCE] 350e08999f5SMatthew G Knepley . -mat_umfpack_fixq <0> - Control[UMFPACK_FIXQ] 351e08999f5SMatthew G Knepley . -mat_umfpack_aggressive <1> - Control[UMFPACK_AGGRESSIVE] 3522f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE] 353e08999f5SMatthew G Knepley . -mat_umfpack_sym_pivot_tolerance <0.001> - Control[UMFPACK_SYM_PIVOT_TOLERANCE] 354e08999f5SMatthew G Knepley . -mat_umfpack_scale <NONE> - (choose one of) NONE SUM MAX 3552f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT] 356e08999f5SMatthew G Knepley . -mat_umfpack_droptol <0> - Control[UMFPACK_DROPTOL] 3572f71e704SKris Buschelman - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP] 3582f71e704SKris Buschelman 3592f71e704SKris Buschelman Level: beginner 360a364b7d2SBarry Smith 361a364b7d2SBarry Smith Note: UMFPACK is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html 3622f71e704SKris Buschelman 363d45987f3SHong Zhang .seealso: PCLU, MATSOLVERSUPERLU, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage 3642f71e704SKris Buschelman M*/ 365b2573a8aSBarry Smith 3668cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_umfpack(Mat A,MatFactorType ftype,Mat *F) 3673519fcdcSHong Zhang { 3683519fcdcSHong Zhang Mat B; 3693519fcdcSHong Zhang Mat_UMFPACK *lu; 3703519fcdcSHong Zhang PetscErrorCode ierr; 3718592901bSBarry Smith PetscInt m=A->rmap->n,n=A->cmap->n,idx; 3722f71e704SKris Buschelman 3732205254eSKarl Rupp const char *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC"}; 3742205254eSKarl Rupp const char *scale[] ={"NONE","SUM","MAX"}; 375ace3abfcSBarry Smith PetscBool flg; 3762d4e2982SHong Zhang 3773519fcdcSHong Zhang PetscFunctionBegin; 3783519fcdcSHong Zhang /* Create the factorization matrix F */ 379ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 3803519fcdcSHong Zhang ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr); 381b9c12af5SBarry Smith ierr = PetscStrallocpy("umfpack",&((PetscObject)B)->type_name);CHKERRQ(ierr); 382b9c12af5SBarry Smith ierr = MatSetUp(B);CHKERRQ(ierr); 383b9c12af5SBarry Smith 384b00a9115SJed Brown ierr = PetscNewLog(B,&lu);CHKERRQ(ierr); 3852205254eSKarl Rupp 386b9c12af5SBarry Smith B->data = lu; 387b9c12af5SBarry Smith B->ops->getinfo = MatGetInfo_External; 3883519fcdcSHong Zhang B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK; 3893519fcdcSHong Zhang B->ops->destroy = MatDestroy_UMFPACK; 3903519fcdcSHong Zhang B->ops->view = MatView_UMFPACK; 3911aef8b4cSStefano Zampini B->ops->matsolve = NULL; 3922205254eSKarl Rupp 393bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_umfpack);CHKERRQ(ierr); 3942205254eSKarl Rupp 395d5f3da31SBarry Smith B->factortype = MAT_FACTOR_LU; 3963519fcdcSHong Zhang B->assembled = PETSC_TRUE; /* required by -ksp_view */ 3973519fcdcSHong Zhang B->preallocated = PETSC_TRUE; 3983519fcdcSHong Zhang 39900c67f3bSHong Zhang ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 40000c67f3bSHong Zhang ierr = PetscStrallocpy(MATSOLVERUMFPACK,&B->solvertype);CHKERRQ(ierr); 40100c67f3bSHong Zhang 4023519fcdcSHong Zhang /* initializations */ 4033519fcdcSHong Zhang /* ------------------------------------------------*/ 4043519fcdcSHong Zhang /* get the default control parameters */ 4058592901bSBarry Smith umfpack_UMF_defaults(lu->Control); 4060298fd71SBarry Smith lu->perm_c = NULL; /* use defaul UMFPACK col permutation */ 4073519fcdcSHong Zhang lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */ 4083519fcdcSHong Zhang 409ce94432eSBarry Smith ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr); 4103519fcdcSHong Zhang /* Control parameters used by reporting routiones */ 4110298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],NULL);CHKERRQ(ierr); 4123519fcdcSHong Zhang 4133519fcdcSHong Zhang /* Control parameters for symbolic factorization */ 414fcfd50ebSBarry Smith ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,3,strategy[0],&idx,&flg);CHKERRQ(ierr); 4153519fcdcSHong Zhang if (flg) { 4163519fcdcSHong Zhang switch (idx) { 4173519fcdcSHong Zhang case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break; 4183519fcdcSHong Zhang case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break; 4193519fcdcSHong Zhang case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break; 4203519fcdcSHong Zhang } 4213519fcdcSHong Zhang } 4228caf3d72SBarry Smith ierr = PetscOptionsEList("-mat_umfpack_ordering","Internal ordering method","None",UmfpackOrderingTypes,sizeof(UmfpackOrderingTypes)/sizeof(UmfpackOrderingTypes[0]),UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]],&idx,&flg);CHKERRQ(ierr); 423fcfd50ebSBarry Smith if (flg) lu->Control[UMFPACK_ORDERING] = (int)idx; 4240298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_dense_col","Control[UMFPACK_DENSE_COL]","None",lu->Control[UMFPACK_DENSE_COL],&lu->Control[UMFPACK_DENSE_COL],NULL);CHKERRQ(ierr); 4250298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_dense_row","Control[UMFPACK_DENSE_ROW]","None",lu->Control[UMFPACK_DENSE_ROW],&lu->Control[UMFPACK_DENSE_ROW],NULL);CHKERRQ(ierr); 4260298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_amd_dense","Control[UMFPACK_AMD_DENSE]","None",lu->Control[UMFPACK_AMD_DENSE],&lu->Control[UMFPACK_AMD_DENSE],NULL);CHKERRQ(ierr); 4270298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_block_size","Control[UMFPACK_BLOCK_SIZE]","None",lu->Control[UMFPACK_BLOCK_SIZE],&lu->Control[UMFPACK_BLOCK_SIZE],NULL);CHKERRQ(ierr); 4280298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],NULL);CHKERRQ(ierr); 4290298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],NULL);CHKERRQ(ierr); 4303519fcdcSHong Zhang 4313519fcdcSHong Zhang /* Control parameters used by numeric factorization */ 4320298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_pivot_tolerance","Control[UMFPACK_PIVOT_TOLERANCE]","None",lu->Control[UMFPACK_PIVOT_TOLERANCE],&lu->Control[UMFPACK_PIVOT_TOLERANCE],NULL);CHKERRQ(ierr); 4330298fd71SBarry Smith ierr = 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);CHKERRQ(ierr); 4343519fcdcSHong Zhang ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr); 4353519fcdcSHong Zhang if (flg) { 4363519fcdcSHong Zhang switch (idx) { 4373519fcdcSHong Zhang case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break; 4383519fcdcSHong Zhang case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break; 4393519fcdcSHong Zhang case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break; 4403519fcdcSHong Zhang } 4413519fcdcSHong Zhang } 4420298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_alloc_init","Control[UMFPACK_ALLOC_INIT]","None",lu->Control[UMFPACK_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],NULL);CHKERRQ(ierr); 4430298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_front_alloc_init","Control[UMFPACK_FRONT_ALLOC_INIT]","None",lu->Control[UMFPACK_FRONT_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],NULL);CHKERRQ(ierr); 4440298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],NULL);CHKERRQ(ierr); 4453519fcdcSHong Zhang 4463519fcdcSHong Zhang /* Control parameters used by solve */ 4470298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],NULL);CHKERRQ(ierr); 4483519fcdcSHong Zhang 4493519fcdcSHong Zhang /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */ 45076a34f28SBarry Smith ierr = PetscOptionsName("-pc_factor_mat_ordering_type","Ordering to do factorization in","MatGetOrdering",&lu->PetscMatOrdering);CHKERRQ(ierr); 4513519fcdcSHong Zhang PetscOptionsEnd(); 4523519fcdcSHong Zhang *F = B; 4533519fcdcSHong Zhang PetscFunctionReturn(0); 4543519fcdcSHong Zhang } 455b2573a8aSBarry Smith 456db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat,MatFactorType,Mat*); 457db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat,MatFactorType,Mat*); 458db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat,MatFactorType,Mat*); 4592d4e2982SHong Zhang 46029b38603SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_SuiteSparse(void) 46142c9c57cSBarry Smith { 46242c9c57cSBarry Smith PetscErrorCode ierr; 46342c9c57cSBarry Smith 46442c9c57cSBarry Smith PetscFunctionBegin; 46542c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERUMFPACK,MATSEQAIJ, MAT_FACTOR_LU,MatGetFactor_seqaij_umfpack);CHKERRQ(ierr); 46642c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERCHOLMOD,MATSEQAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_seqaij_cholmod);CHKERRQ(ierr); 46742c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERCHOLMOD,MATSEQSBAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_seqsbaij_cholmod);CHKERRQ(ierr); 46842c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERKLU,MATSEQAIJ, MAT_FACTOR_LU,MatGetFactor_seqaij_klu);CHKERRQ(ierr); 46942c9c57cSBarry Smith PetscFunctionReturn(0); 47042c9c57cSBarry Smith } 471