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 931677d0b8SKris Buschelman #undef __FUNCT__ 94f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_UMFPACK" 954b019bd2SJed Brown static PetscErrorCode MatDestroy_UMFPACK(Mat A) 96dfbe8321SBarry Smith { 97dfbe8321SBarry Smith PetscErrorCode ierr; 98*b9c12af5SBarry Smith Mat_UMFPACK *lu=(Mat_UMFPACK*)A->data; 991677d0b8SKris Buschelman 1001677d0b8SKris Buschelman PetscFunctionBegin; 101*b9c12af5SBarry Smith if (lu->CleanUpUMFPACK) { 1028592901bSBarry Smith umfpack_UMF_free_symbolic(&lu->Symbolic); 1038592901bSBarry Smith umfpack_UMF_free_numeric(&lu->Numeric); 1041677d0b8SKris Buschelman ierr = PetscFree(lu->Wi);CHKERRQ(ierr); 1051677d0b8SKris Buschelman ierr = PetscFree(lu->W);CHKERRQ(ierr); 1061677d0b8SKris Buschelman ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 1071677d0b8SKris Buschelman } 108ae26a541SJed Brown ierr = MatDestroy(&lu->A);CHKERRQ(ierr); 109*b9c12af5SBarry Smith ierr = PetscFree(A->data);CHKERRQ(ierr); 1101677d0b8SKris Buschelman PetscFunctionReturn(0); 1111677d0b8SKris Buschelman } 1121677d0b8SKris Buschelman 1131677d0b8SKris Buschelman #undef __FUNCT__ 1144b019bd2SJed Brown #define __FUNCT__ "MatSolve_UMFPACK_Private" 1154b019bd2SJed Brown static PetscErrorCode MatSolve_UMFPACK_Private(Mat A,Vec b,Vec x,int uflag) 1166849ba73SBarry Smith { 117*b9c12af5SBarry Smith Mat_UMFPACK *lu = (Mat_UMFPACK*)A->data; 118ae26a541SJed Brown Mat_SeqAIJ *a = (Mat_SeqAIJ*)lu->A->data; 119d9ca1df4SBarry Smith PetscScalar *av = a->a,*xa; 120d9ca1df4SBarry Smith const PetscScalar *ba; 121dfbe8321SBarry Smith PetscErrorCode ierr; 122ae26a541SJed Brown PetscInt *ai = a->i,*aj = a->j,status; 123fcf6e9abSJed Brown static PetscBool cite = PETSC_FALSE; 1241677d0b8SKris Buschelman 1251677d0b8SKris Buschelman PetscFunctionBegin; 126fcf6e9abSJed 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); 1272d4e2982SHong Zhang /* solve Ax = b by umfpack_*_wsolve */ 1281677d0b8SKris Buschelman /* ----------------------------------*/ 1292d4e2982SHong Zhang 130ae26a541SJed Brown if (!lu->Wi) { /* first time, allocate working space for wsolve */ 131785e854fSJed Brown ierr = PetscMalloc1(A->rmap->n,&lu->Wi);CHKERRQ(ierr); 132785e854fSJed Brown ierr = PetscMalloc1(5*A->rmap->n,&lu->W);CHKERRQ(ierr); 133ae26a541SJed Brown } 134ae26a541SJed Brown 135d9ca1df4SBarry Smith ierr = VecGetArrayRead(b,&ba); 1361677d0b8SKris Buschelman ierr = VecGetArray(x,&xa); 13779c34000SHong Zhang #if defined(PETSC_USE_COMPLEX) 138b0043f70SBarry 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); 13979c34000SHong Zhang #else 1404b019bd2SJed Brown status = umfpack_UMF_wsolve(uflag,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W); 14179c34000SHong Zhang #endif 1428592901bSBarry Smith umfpack_UMF_report_info(lu->Control, lu->Info); 1431677d0b8SKris Buschelman if (status < 0) { 1448592901bSBarry Smith umfpack_UMF_report_status(lu->Control, status); 145e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_wsolve failed"); 1461677d0b8SKris Buschelman } 1471677d0b8SKris Buschelman 148d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(b,&ba);CHKERRQ(ierr); 149057dcbc9SJed Brown ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr); 1504b019bd2SJed Brown PetscFunctionReturn(0); 1514b019bd2SJed Brown } 1522a325a84SHong Zhang 1534b019bd2SJed Brown #undef __FUNCT__ 1544b019bd2SJed Brown #define __FUNCT__ "MatSolve_UMFPACK" 1554b019bd2SJed Brown static PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x) 1564b019bd2SJed Brown { 1574b019bd2SJed Brown PetscErrorCode ierr; 1584b019bd2SJed Brown 1594b019bd2SJed Brown PetscFunctionBegin; 1604b019bd2SJed Brown /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */ 1614b019bd2SJed Brown ierr = MatSolve_UMFPACK_Private(A,b,x,UMFPACK_Aat);CHKERRQ(ierr); 1624b019bd2SJed Brown PetscFunctionReturn(0); 1634b019bd2SJed Brown } 1644b019bd2SJed Brown 1654b019bd2SJed Brown #undef __FUNCT__ 1664b019bd2SJed Brown #define __FUNCT__ "MatSolveTranspose_UMFPACK" 1674b019bd2SJed Brown static PetscErrorCode MatSolveTranspose_UMFPACK(Mat A,Vec b,Vec x) 1684b019bd2SJed Brown { 1694b019bd2SJed Brown PetscErrorCode ierr; 1704b019bd2SJed Brown 1714b019bd2SJed Brown PetscFunctionBegin; 1724b019bd2SJed Brown /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */ 1734b019bd2SJed Brown ierr = MatSolve_UMFPACK_Private(A,b,x,UMFPACK_A);CHKERRQ(ierr); 1741677d0b8SKris Buschelman PetscFunctionReturn(0); 1751677d0b8SKris Buschelman } 1761677d0b8SKris Buschelman 1771677d0b8SKris Buschelman #undef __FUNCT__ 178f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_UMFPACK" 1794b019bd2SJed Brown static PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat F,Mat A,const MatFactorInfo *info) 1806849ba73SBarry Smith { 181*b9c12af5SBarry Smith Mat_UMFPACK *lu = (Mat_UMFPACK*)(F)->data; 182ae26a541SJed Brown Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 183ae26a541SJed Brown PetscInt *ai = a->i,*aj=a->j,status; 184ae26a541SJed Brown PetscScalar *av = a->a; 1856849ba73SBarry Smith PetscErrorCode ierr; 1861677d0b8SKris Buschelman 1871677d0b8SKris Buschelman PetscFunctionBegin; 1881677d0b8SKris Buschelman /* numeric factorization of A' */ 1891677d0b8SKris Buschelman /* ----------------------------*/ 1902d4e2982SHong Zhang 1918592901bSBarry Smith if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) { 1928592901bSBarry Smith umfpack_UMF_free_numeric(&lu->Numeric); 1938592901bSBarry Smith } 1942d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX) 1958592901bSBarry Smith status = umfpack_UMF_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info); 1962d4e2982SHong Zhang #else 1978592901bSBarry Smith status = umfpack_UMF_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info); 1988592901bSBarry Smith #endif 1999f42a82aSMatthew Knepley if (status < 0) { 2008592901bSBarry Smith umfpack_UMF_report_status(lu->Control, status); 201e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_numeric failed"); 2029f42a82aSMatthew Knepley } 2032d4e2982SHong Zhang /* report numeric factorization of A' when Control[PRL] > 3 */ 2048592901bSBarry Smith (void) umfpack_UMF_report_numeric(lu->Numeric, lu->Control); 2051677d0b8SKris Buschelman 206ae26a541SJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 207ae26a541SJed Brown ierr = MatDestroy(&lu->A);CHKERRQ(ierr); 2082205254eSKarl Rupp 209ae26a541SJed Brown lu->A = A; 2102a325a84SHong Zhang lu->flg = SAME_NONZERO_PATTERN; 2112a325a84SHong Zhang lu->CleanUpUMFPACK = PETSC_TRUE; 2124b019bd2SJed Brown F->ops->solve = MatSolve_UMFPACK; 2134b019bd2SJed Brown F->ops->solvetranspose = MatSolveTranspose_UMFPACK; 2141677d0b8SKris Buschelman PetscFunctionReturn(0); 2151677d0b8SKris Buschelman } 2161677d0b8SKris Buschelman 2171677d0b8SKris Buschelman /* 2181677d0b8SKris Buschelman Note the r permutation is ignored 2191677d0b8SKris Buschelman */ 2201677d0b8SKris Buschelman #undef __FUNCT__ 221f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK" 2224b019bd2SJed Brown static PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 2236849ba73SBarry Smith { 224ae26a541SJed Brown Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 225*b9c12af5SBarry Smith Mat_UMFPACK *lu = (Mat_UMFPACK*)(F->data); 226dfbe8321SBarry Smith PetscErrorCode ierr; 227ae26a541SJed Brown PetscInt i,*ai = a->i,*aj = a->j,m=A->rmap->n,n=A->cmap->n; 228989213a4SSatish Balay #if !defined(PETSC_USE_COMPLEX) 229ae26a541SJed Brown PetscScalar *av = a->a; 230989213a4SSatish Balay #endif 2315d0c19d7SBarry Smith const PetscInt *ra; 2328592901bSBarry Smith PetscInt status; 2331677d0b8SKris Buschelman 2341677d0b8SKris Buschelman PetscFunctionBegin; 235fcfd50ebSBarry Smith if (lu->PetscMatOrdering) { 2360f6efc10SHong Zhang ierr = ISGetIndices(r,&ra);CHKERRQ(ierr); 237785e854fSJed Brown ierr = PetscMalloc1(m,&lu->perm_c);CHKERRQ(ierr); 2382d4e2982SHong Zhang /* we cannot simply memcpy on 64 bit archs */ 2392d4e2982SHong Zhang for (i = 0; i < m; i++) lu->perm_c[i] = ra[i]; 2400f6efc10SHong Zhang ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr); 2411677d0b8SKris Buschelman } 2421677d0b8SKris Buschelman 2431677d0b8SKris Buschelman /* print the control parameters */ 2448592901bSBarry Smith if (lu->Control[UMFPACK_PRL] > 1) umfpack_UMF_report_control(lu->Control); 2451677d0b8SKris Buschelman 2461677d0b8SKris Buschelman /* symbolic factorization of A' */ 2471677d0b8SKris Buschelman /* ---------------------------------------------------------------------- */ 248fcfd50ebSBarry Smith if (lu->PetscMatOrdering) { /* use Petsc row ordering */ 2498592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX) 250ae26a541SJed Brown status = umfpack_UMF_qsymbolic(n,m,ai,aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info); 2512d4e2982SHong Zhang #else 252b0043f70SBarry Smith status = umfpack_UMF_qsymbolic(n,m,ai,aj,NULL,NULL,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info); 2538592901bSBarry Smith #endif 2542d4e2982SHong Zhang } else { /* use Umfpack col ordering */ 2558592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX) 256ae26a541SJed Brown status = umfpack_UMF_symbolic(n,m,ai,aj,av,&lu->Symbolic,lu->Control,lu->Info); 2578592901bSBarry Smith #else 258ae26a541SJed Brown status = umfpack_UMF_symbolic(n,m,ai,aj,NULL,NULL,&lu->Symbolic,lu->Control,lu->Info); 2598592901bSBarry Smith #endif 2602d4e2982SHong Zhang } 2612d4e2982SHong Zhang if (status < 0) { 2628592901bSBarry Smith umfpack_UMF_report_info(lu->Control, lu->Info); 2638592901bSBarry Smith umfpack_UMF_report_status(lu->Control, status); 264e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_symbolic failed"); 2652d4e2982SHong Zhang } 2662d4e2982SHong Zhang /* report sumbolic factorization of A' when Control[PRL] > 3 */ 2678592901bSBarry Smith (void) umfpack_UMF_report_symbolic(lu->Symbolic, lu->Control); 2681677d0b8SKris Buschelman 2691677d0b8SKris Buschelman lu->flg = DIFFERENT_NONZERO_PATTERN; 270e1b4c3a1SKris Buschelman lu->CleanUpUMFPACK = PETSC_TRUE; 271719d5645SBarry Smith (F)->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK; 2721677d0b8SKris Buschelman PetscFunctionReturn(0); 2731677d0b8SKris Buschelman } 2741677d0b8SKris Buschelman 2751677d0b8SKris Buschelman #undef __FUNCT__ 276f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_UMFPACK" 2774b019bd2SJed Brown static PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer) 2786849ba73SBarry Smith { 279*b9c12af5SBarry Smith Mat_UMFPACK *lu= (Mat_UMFPACK*)A->data; 280dfbe8321SBarry Smith PetscErrorCode ierr; 281f0c56d0fSKris Buschelman 2821677d0b8SKris Buschelman PetscFunctionBegin; 2831677d0b8SKris Buschelman /* check if matrix is UMFPACK type */ 284f0c56d0fSKris Buschelman if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0); 2851677d0b8SKris Buschelman 2861677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr); 2871677d0b8SKris Buschelman /* Control parameters used by reporting routiones */ 2881677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr); 2891677d0b8SKris Buschelman 2901677d0b8SKris Buschelman /* Control parameters used by symbolic factorization */ 2910f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr); 2921677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr); 2931677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr); 2940f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr); 2951677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr); 2960f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr); 2970f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr); 2981677d0b8SKris Buschelman 2991677d0b8SKris Buschelman /* Control parameters used by numeric factorization */ 3001677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr); 3010f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr); 3020f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr); 3031677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr); 3040f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr); 3051677d0b8SKris Buschelman 3061677d0b8SKris Buschelman /* Control parameters used by solve */ 3071677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr); 3081677d0b8SKris Buschelman 3091677d0b8SKris Buschelman /* mat ordering */ 310fcfd50ebSBarry Smith if (!lu->PetscMatOrdering) { 311fcfd50ebSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_ORDERING]: %s (not using the PETSc ordering)\n",UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]]);CHKERRQ(ierr); 312fcfd50ebSBarry Smith } 3131677d0b8SKris Buschelman PetscFunctionReturn(0); 3141677d0b8SKris Buschelman } 3151677d0b8SKris Buschelman 316d844382dSKris Buschelman #undef __FUNCT__ 317f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_UMFPACK" 3184b019bd2SJed Brown static PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer) 3196849ba73SBarry Smith { 320dfbe8321SBarry Smith PetscErrorCode ierr; 321ace3abfcSBarry Smith PetscBool iascii; 322d844382dSKris Buschelman PetscViewerFormat format; 323d844382dSKris Buschelman 324d844382dSKris Buschelman PetscFunctionBegin; 325251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 32632077d6dSBarry Smith if (iascii) { 327d844382dSKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 3282f59403fSHong Zhang if (format == PETSC_VIEWER_ASCII_INFO) { 329f0c56d0fSKris Buschelman ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr); 330d844382dSKris Buschelman } 331d844382dSKris Buschelman } 332d844382dSKris Buschelman PetscFunctionReturn(0); 333d844382dSKris Buschelman } 334d844382dSKris Buschelman 33535bd34faSBarry Smith #undef __FUNCT__ 33635bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_umfpack" 33735bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_seqaij_umfpack(Mat A,const MatSolverPackage *type) 33835bd34faSBarry Smith { 33935bd34faSBarry Smith PetscFunctionBegin; 3402692d6eeSBarry Smith *type = MATSOLVERUMFPACK; 34135bd34faSBarry Smith PetscFunctionReturn(0); 34235bd34faSBarry Smith } 34335bd34faSBarry Smith 34435bd34faSBarry Smith 3452f71e704SKris Buschelman /*MC 3462692d6eeSBarry Smith MATSOLVERUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices 3472f71e704SKris Buschelman via the external package UMFPACK. 3482f71e704SKris Buschelman 349c2b89b5dSBarry Smith Use ./configure --download-suitesparse to install PETSc to use UMFPACK 350c2b89b5dSBarry Smith 351c2b89b5dSBarry Smith Use -pc_type lu -pc_factor_mat_solver_package umfpack to us this direct solver 3522f71e704SKris Buschelman 3532f71e704SKris Buschelman Consult UMFPACK documentation for more information about the Control parameters 3542f71e704SKris Buschelman which correspond to the options database keys below. 3552f71e704SKris Buschelman 3562f71e704SKris Buschelman Options Database Keys: 357ba41fbb6SBarry Smith + -mat_umfpack_ordering - CHOLMOD, AMD, GIVEN, METIS, BEST, NONE 358ba41fbb6SBarry Smith . -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL] 359e08999f5SMatthew G Knepley . -mat_umfpack_strategy <AUTO> - (choose one of) AUTO UNSYMMETRIC SYMMETRIC 2BY2 3602f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL] 361e08999f5SMatthew G Knepley . -mat_umfpack_dense_row <0.2> - Control[UMFPACK_DENSE_ROW] 362e08999f5SMatthew G Knepley . -mat_umfpack_amd_dense <10> - Control[UMFPACK_AMD_DENSE] 3632f71e704SKris Buschelman . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE] 364e08999f5SMatthew G Knepley . -mat_umfpack_2by2_tolerance <0.01> - Control[UMFPACK_2BY2_TOLERANCE] 365e08999f5SMatthew G Knepley . -mat_umfpack_fixq <0> - Control[UMFPACK_FIXQ] 366e08999f5SMatthew G Knepley . -mat_umfpack_aggressive <1> - Control[UMFPACK_AGGRESSIVE] 3672f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE] 368e08999f5SMatthew G Knepley . -mat_umfpack_sym_pivot_tolerance <0.001> - Control[UMFPACK_SYM_PIVOT_TOLERANCE] 369e08999f5SMatthew G Knepley . -mat_umfpack_scale <NONE> - (choose one of) NONE SUM MAX 3702f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT] 371e08999f5SMatthew G Knepley . -mat_umfpack_droptol <0> - Control[UMFPACK_DROPTOL] 3722f71e704SKris Buschelman - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP] 3732f71e704SKris Buschelman 3742f71e704SKris Buschelman Level: beginner 375a364b7d2SBarry Smith 376a364b7d2SBarry Smith Note: UMFPACK is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html 3772f71e704SKris Buschelman 378d45987f3SHong Zhang .seealso: PCLU, MATSOLVERSUPERLU, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage 3792f71e704SKris Buschelman M*/ 380b2573a8aSBarry Smith 3813519fcdcSHong Zhang #undef __FUNCT__ 3823519fcdcSHong Zhang #define __FUNCT__ "MatGetFactor_seqaij_umfpack" 3838cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_umfpack(Mat A,MatFactorType ftype,Mat *F) 3843519fcdcSHong Zhang { 3853519fcdcSHong Zhang Mat B; 3863519fcdcSHong Zhang Mat_UMFPACK *lu; 3873519fcdcSHong Zhang PetscErrorCode ierr; 3888592901bSBarry Smith PetscInt m=A->rmap->n,n=A->cmap->n,idx; 3892f71e704SKris Buschelman 3902205254eSKarl Rupp const char *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC"}; 3912205254eSKarl Rupp const char *scale[] ={"NONE","SUM","MAX"}; 392ace3abfcSBarry Smith PetscBool flg; 3932d4e2982SHong Zhang 3943519fcdcSHong Zhang PetscFunctionBegin; 3953519fcdcSHong Zhang /* Create the factorization matrix F */ 396ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 3973519fcdcSHong Zhang ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr); 398*b9c12af5SBarry Smith ierr = PetscStrallocpy("umfpack",&((PetscObject)B)->type_name);CHKERRQ(ierr); 399*b9c12af5SBarry Smith ierr = MatSetUp(B);CHKERRQ(ierr); 400*b9c12af5SBarry Smith 401b00a9115SJed Brown ierr = PetscNewLog(B,&lu);CHKERRQ(ierr); 4022205254eSKarl Rupp 403*b9c12af5SBarry Smith B->data = lu; 404*b9c12af5SBarry Smith B->ops->getinfo = MatGetInfo_External; 4053519fcdcSHong Zhang B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK; 4063519fcdcSHong Zhang B->ops->destroy = MatDestroy_UMFPACK; 4073519fcdcSHong Zhang B->ops->view = MatView_UMFPACK; 4081aef8b4cSStefano Zampini B->ops->matsolve = NULL; 4092205254eSKarl Rupp 410bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_umfpack);CHKERRQ(ierr); 4112205254eSKarl Rupp 412d5f3da31SBarry Smith B->factortype = MAT_FACTOR_LU; 4133519fcdcSHong Zhang B->assembled = PETSC_TRUE; /* required by -ksp_view */ 4143519fcdcSHong Zhang B->preallocated = PETSC_TRUE; 4153519fcdcSHong Zhang 41600c67f3bSHong Zhang ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 41700c67f3bSHong Zhang ierr = PetscStrallocpy(MATSOLVERUMFPACK,&B->solvertype);CHKERRQ(ierr); 41800c67f3bSHong Zhang 4193519fcdcSHong Zhang /* initializations */ 4203519fcdcSHong Zhang /* ------------------------------------------------*/ 4213519fcdcSHong Zhang /* get the default control parameters */ 4228592901bSBarry Smith umfpack_UMF_defaults(lu->Control); 4230298fd71SBarry Smith lu->perm_c = NULL; /* use defaul UMFPACK col permutation */ 4243519fcdcSHong Zhang lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */ 4253519fcdcSHong Zhang 426ce94432eSBarry Smith ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr); 4273519fcdcSHong Zhang /* Control parameters used by reporting routiones */ 4280298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],NULL);CHKERRQ(ierr); 4293519fcdcSHong Zhang 4303519fcdcSHong Zhang /* Control parameters for symbolic factorization */ 431fcfd50ebSBarry Smith ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,3,strategy[0],&idx,&flg);CHKERRQ(ierr); 4323519fcdcSHong Zhang if (flg) { 4333519fcdcSHong Zhang switch (idx) { 4343519fcdcSHong Zhang case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break; 4353519fcdcSHong Zhang case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break; 4363519fcdcSHong Zhang case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break; 4373519fcdcSHong Zhang } 4383519fcdcSHong Zhang } 4398caf3d72SBarry 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); 440fcfd50ebSBarry Smith if (flg) lu->Control[UMFPACK_ORDERING] = (int)idx; 4410298fd71SBarry 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); 4420298fd71SBarry 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); 4430298fd71SBarry 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); 4440298fd71SBarry 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); 4450298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],NULL);CHKERRQ(ierr); 4460298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],NULL);CHKERRQ(ierr); 4473519fcdcSHong Zhang 4483519fcdcSHong Zhang /* Control parameters used by numeric factorization */ 4490298fd71SBarry 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); 4500298fd71SBarry 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); 4513519fcdcSHong Zhang ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr); 4523519fcdcSHong Zhang if (flg) { 4533519fcdcSHong Zhang switch (idx) { 4543519fcdcSHong Zhang case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break; 4553519fcdcSHong Zhang case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break; 4563519fcdcSHong Zhang case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break; 4573519fcdcSHong Zhang } 4583519fcdcSHong Zhang } 4590298fd71SBarry 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); 4600298fd71SBarry 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); 4610298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],NULL);CHKERRQ(ierr); 4623519fcdcSHong Zhang 4633519fcdcSHong Zhang /* Control parameters used by solve */ 4640298fd71SBarry Smith ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],NULL);CHKERRQ(ierr); 4653519fcdcSHong Zhang 4663519fcdcSHong Zhang /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */ 46776a34f28SBarry Smith ierr = PetscOptionsName("-pc_factor_mat_ordering_type","Ordering to do factorization in","MatGetOrdering",&lu->PetscMatOrdering);CHKERRQ(ierr); 4683519fcdcSHong Zhang PetscOptionsEnd(); 4693519fcdcSHong Zhang *F = B; 4703519fcdcSHong Zhang PetscFunctionReturn(0); 4713519fcdcSHong Zhang } 472b2573a8aSBarry Smith 473db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat,MatFactorType,Mat*); 474db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat,MatFactorType,Mat*); 475db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat,MatFactorType,Mat*); 4762d4e2982SHong Zhang 47742c9c57cSBarry Smith #undef __FUNCT__ 47842c9c57cSBarry Smith #define __FUNCT__ "MatSolverPackageRegister_SuiteSparse" 47929b38603SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_SuiteSparse(void) 48042c9c57cSBarry Smith { 48142c9c57cSBarry Smith PetscErrorCode ierr; 48242c9c57cSBarry Smith 48342c9c57cSBarry Smith PetscFunctionBegin; 48442c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERUMFPACK,MATSEQAIJ, MAT_FACTOR_LU,MatGetFactor_seqaij_umfpack);CHKERRQ(ierr); 48542c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERCHOLMOD,MATSEQAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_seqaij_cholmod);CHKERRQ(ierr); 48642c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERCHOLMOD,MATSEQSBAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_seqsbaij_cholmod);CHKERRQ(ierr); 48742c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERKLU,MATSEQAIJ, MAT_FACTOR_LU,MatGetFactor_seqaij_klu);CHKERRQ(ierr); 48842c9c57cSBarry Smith PetscFunctionReturn(0); 48942c9c57cSBarry Smith } 490