11677d0b8SKris Buschelman /*$Id: umfpack.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/ 21677d0b8SKris Buschelman 31677d0b8SKris Buschelman /* 4*59698cb0SHong Zhang Provides an interface to the UMFPACKv4.3 sparse solver 51677d0b8SKris Buschelman */ 61677d0b8SKris Buschelman 71677d0b8SKris Buschelman #include "src/mat/impls/aij/seq/aij.h" 81677d0b8SKris Buschelman 91677d0b8SKris Buschelman EXTERN_C_BEGIN 101677d0b8SKris Buschelman #include "umfpack.h" 111677d0b8SKris Buschelman EXTERN_C_END 121677d0b8SKris Buschelman 131677d0b8SKris Buschelman typedef struct { 141677d0b8SKris Buschelman void *Symbolic, *Numeric; 151677d0b8SKris Buschelman double Info[UMFPACK_INFO], Control[UMFPACK_CONTROL],*W; 161677d0b8SKris Buschelman int *Wi,*ai,*aj,*perm_c; 171677d0b8SKris Buschelman PetscScalar *av; 181677d0b8SKris Buschelman MatStructure flg; 191677d0b8SKris Buschelman PetscTruth PetscMatOdering; 201677d0b8SKris Buschelman 211677d0b8SKris Buschelman /* A few function pointers for inheritance */ 22f0c56d0fSKris Buschelman int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*); 231677d0b8SKris Buschelman int (*MatView)(Mat,PetscViewer); 241677d0b8SKris Buschelman int (*MatAssemblyEnd)(Mat,MatAssemblyType); 25d844382dSKris Buschelman int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*); 261677d0b8SKris Buschelman int (*MatDestroy)(Mat); 271677d0b8SKris Buschelman 281677d0b8SKris Buschelman /* Flag to clean up UMFPACK objects during Destroy */ 291677d0b8SKris Buschelman PetscTruth CleanUpUMFPACK; 30f0c56d0fSKris Buschelman } Mat_UMFPACK; 31f0c56d0fSKris Buschelman 32f0c56d0fSKris Buschelman EXTERN int MatDuplicate_UMFPACK(Mat,MatDuplicateOption,Mat*); 33d844382dSKris Buschelman 34d844382dSKris Buschelman EXTERN_C_BEGIN 35d844382dSKris Buschelman #undef __FUNCT__ 36d844382dSKris Buschelman #define __FUNCT__ "MatConvert_UMFPACK_SeqAIJ" 378e9aea5cSBarry Smith int MatConvert_UMFPACK_SeqAIJ(Mat A,const MatType type,Mat *newmat) { 38d844382dSKris Buschelman /* This routine is only called to convert an unfactored PETSc-UMFPACK matrix */ 39d844382dSKris Buschelman /* to its base PETSc type, so we will ignore 'MatType type'. */ 40d844382dSKris Buschelman int ierr; 41d844382dSKris Buschelman Mat B=*newmat; 42f0c56d0fSKris Buschelman Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr; 43d844382dSKris Buschelman 44d844382dSKris Buschelman PetscFunctionBegin; 45d844382dSKris Buschelman if (B != A) { 46d844382dSKris Buschelman /* This routine was inherited from SeqAIJ. */ 47d844382dSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 48f0c56d0fSKris Buschelman } 49d844382dSKris Buschelman /* Reset the original function pointers */ 50f0c56d0fSKris Buschelman A->ops->duplicate = lu->MatDuplicate; 51d844382dSKris Buschelman A->ops->view = lu->MatView; 52d844382dSKris Buschelman A->ops->assemblyend = lu->MatAssemblyEnd; 53d844382dSKris Buschelman A->ops->lufactorsymbolic = lu->MatLUFactorSymbolic; 54d844382dSKris Buschelman A->ops->destroy = lu->MatDestroy; 55d844382dSKris Buschelman 56d844382dSKris Buschelman ierr = PetscFree(lu);CHKERRQ(ierr); 57d844382dSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr); 58d844382dSKris Buschelman *newmat = B; 59d844382dSKris Buschelman PetscFunctionReturn(0); 60d844382dSKris Buschelman } 61d844382dSKris Buschelman EXTERN_C_END 621677d0b8SKris Buschelman 631677d0b8SKris Buschelman #undef __FUNCT__ 64f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_UMFPACK" 65f0c56d0fSKris Buschelman int MatDestroy_UMFPACK(Mat A) { 66d844382dSKris Buschelman int ierr; 67f0c56d0fSKris Buschelman Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr; 681677d0b8SKris Buschelman 691677d0b8SKris Buschelman PetscFunctionBegin; 70fb731535SKris Buschelman if (lu->CleanUpUMFPACK) { 711677d0b8SKris Buschelman umfpack_di_free_symbolic(&lu->Symbolic) ; 721677d0b8SKris Buschelman umfpack_di_free_numeric(&lu->Numeric) ; 731677d0b8SKris Buschelman ierr = PetscFree(lu->Wi);CHKERRQ(ierr); 741677d0b8SKris Buschelman ierr = PetscFree(lu->W);CHKERRQ(ierr); 751677d0b8SKris Buschelman if (lu->PetscMatOdering) { 761677d0b8SKris Buschelman ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 771677d0b8SKris Buschelman } 781677d0b8SKris Buschelman } 79d844382dSKris Buschelman ierr = MatConvert_UMFPACK_SeqAIJ(A,MATSEQAIJ,&A);CHKERRQ(ierr); 80d844382dSKris Buschelman ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 811677d0b8SKris Buschelman PetscFunctionReturn(0); 821677d0b8SKris Buschelman } 831677d0b8SKris Buschelman 841677d0b8SKris Buschelman #undef __FUNCT__ 85f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_UMFPACK" 86f0c56d0fSKris Buschelman int MatSolve_UMFPACK(Mat A,Vec b,Vec x) { 87f0c56d0fSKris Buschelman Mat_UMFPACK *lu = (Mat_UMFPACK*)A->spptr; 881677d0b8SKris Buschelman PetscScalar *av=lu->av,*ba,*xa; 891677d0b8SKris Buschelman int ierr,*ai=lu->ai,*aj=lu->aj,status; 901677d0b8SKris Buschelman 911677d0b8SKris Buschelman PetscFunctionBegin; 921677d0b8SKris Buschelman /* solve Ax = b by umfpack_di_wsolve */ 931677d0b8SKris Buschelman /* ----------------------------------*/ 941677d0b8SKris Buschelman ierr = VecGetArray(b,&ba); 951677d0b8SKris Buschelman ierr = VecGetArray(x,&xa); 961677d0b8SKris Buschelman 971677d0b8SKris Buschelman status = umfpack_di_wsolve(UMFPACK_At,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W); 981677d0b8SKris Buschelman umfpack_di_report_info(lu->Control, lu->Info); 991677d0b8SKris Buschelman if (status < 0){ 1001677d0b8SKris Buschelman umfpack_di_report_status(lu->Control, status) ; 1011677d0b8SKris Buschelman SETERRQ(1,"umfpack_di_wsolve failed") ; 1021677d0b8SKris Buschelman } 1031677d0b8SKris Buschelman 1041677d0b8SKris Buschelman ierr = VecRestoreArray(b,&ba); 1051677d0b8SKris Buschelman ierr = VecRestoreArray(x,&xa); 1062a325a84SHong Zhang 1071677d0b8SKris Buschelman PetscFunctionReturn(0); 1081677d0b8SKris Buschelman } 1091677d0b8SKris Buschelman 1101677d0b8SKris Buschelman #undef __FUNCT__ 111f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_UMFPACK" 112f0c56d0fSKris Buschelman int MatLUFactorNumeric_UMFPACK(Mat A,Mat *F) { 113f0c56d0fSKris Buschelman Mat_UMFPACK *lu=(Mat_UMFPACK*)(*F)->spptr; 1141677d0b8SKris Buschelman int *ai=lu->ai,*aj=lu->aj,m=A->m,status,ierr; 1151677d0b8SKris Buschelman PetscScalar *av=lu->av; 1161677d0b8SKris Buschelman 1171677d0b8SKris Buschelman PetscFunctionBegin; 1181677d0b8SKris Buschelman /* numeric factorization of A' */ 1191677d0b8SKris Buschelman /* ----------------------------*/ 1202a325a84SHong Zhang if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){ 1212a325a84SHong Zhang umfpack_di_free_numeric(&lu->Numeric) ; 1222a325a84SHong Zhang } 1231677d0b8SKris Buschelman status = umfpack_di_numeric (ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info) ; 1241677d0b8SKris Buschelman if (status < 0) SETERRQ(1,"umfpack_di_numeric failed"); 1251677d0b8SKris Buschelman /* report numeric factorization of A' when Control[PRL] > 3 */ 1261677d0b8SKris Buschelman (void) umfpack_di_report_numeric (lu->Numeric, lu->Control) ; 1271677d0b8SKris Buschelman 1281677d0b8SKris Buschelman if (lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization */ 1291677d0b8SKris Buschelman /* allocate working space to be used by Solve */ 1301677d0b8SKris Buschelman ierr = PetscMalloc(m * sizeof(int), &lu->Wi);CHKERRQ(ierr); 1311677d0b8SKris Buschelman ierr = PetscMalloc(5*m * sizeof(double), &lu->W);CHKERRQ(ierr); 1321677d0b8SKris Buschelman } 1332a325a84SHong Zhang lu->flg = SAME_NONZERO_PATTERN; 1342a325a84SHong Zhang lu->CleanUpUMFPACK = PETSC_TRUE; 1351677d0b8SKris Buschelman PetscFunctionReturn(0); 1361677d0b8SKris Buschelman } 1371677d0b8SKris Buschelman 1381677d0b8SKris Buschelman /* 1391677d0b8SKris Buschelman Note the r permutation is ignored 1401677d0b8SKris Buschelman */ 1411677d0b8SKris Buschelman #undef __FUNCT__ 142f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK" 143f0c56d0fSKris Buschelman int MatLUFactorSymbolic_UMFPACK(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) { 144e1b4c3a1SKris Buschelman Mat B; 1451677d0b8SKris Buschelman Mat_SeqAIJ *mat=(Mat_SeqAIJ*)A->data; 146f0c56d0fSKris Buschelman Mat_UMFPACK *lu; 1471677d0b8SKris Buschelman int ierr,m=A->m,n=A->n,*ai=mat->i,*aj=mat->j,status,*ca; 1481677d0b8SKris Buschelman PetscScalar *av=mat->a; 1491677d0b8SKris Buschelman 1501677d0b8SKris Buschelman PetscFunctionBegin; 1511677d0b8SKris Buschelman /* Create the factorization matrix F */ 152e1b4c3a1SKris Buschelman ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,m,n,&B);CHKERRQ(ierr); 153be5d1d56SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 154e1b4c3a1SKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 1551677d0b8SKris Buschelman 156f0c56d0fSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK; 157f0c56d0fSKris Buschelman B->ops->solve = MatSolve_UMFPACK; 158e1b4c3a1SKris Buschelman B->factor = FACTOR_LU; 15994b7f48cSBarry Smith B->assembled = PETSC_TRUE; /* required by -ksp_view */ 1601677d0b8SKris Buschelman 161f0c56d0fSKris Buschelman lu = (Mat_UMFPACK*)(B->spptr); 1621677d0b8SKris Buschelman 1631677d0b8SKris Buschelman /* initializations */ 1641677d0b8SKris Buschelman /* ------------------------------------------------*/ 1651677d0b8SKris Buschelman /* get the default control parameters */ 1661677d0b8SKris Buschelman umfpack_di_defaults (lu->Control) ; 1671677d0b8SKris Buschelman lu->perm_c = PETSC_NULL; /* use defaul UMFPACK col permutation */ 1681677d0b8SKris Buschelman 1691677d0b8SKris Buschelman ierr = PetscOptionsBegin(A->comm,A->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr); 1701677d0b8SKris Buschelman /* Control parameters used by reporting routiones */ 1711677d0b8SKris Buschelman ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr); 1721677d0b8SKris Buschelman 1731677d0b8SKris Buschelman /* Control parameters for symbolic factorization */ 1741677d0b8SKris Buschelman ierr = PetscOptionsReal("-mat_umfpack_dense_col","Control[UMFPACK_DENSE_COL]","None",lu->Control[UMFPACK_DENSE_COL],&lu->Control[UMFPACK_DENSE_COL],PETSC_NULL);CHKERRQ(ierr); 1751677d0b8SKris Buschelman ierr = PetscOptionsReal("-mat_umfpack_dense_row","Control[UMFPACK_DENSE_ROW]","None",lu->Control[UMFPACK_DENSE_ROW],&lu->Control[UMFPACK_DENSE_ROW],PETSC_NULL);CHKERRQ(ierr); 1761677d0b8SKris Buschelman ierr = PetscOptionsReal("-mat_umfpack_block_size","Control[UMFPACK_BLOCK_SIZE]","None",lu->Control[UMFPACK_BLOCK_SIZE],&lu->Control[UMFPACK_BLOCK_SIZE],PETSC_NULL);CHKERRQ(ierr); 1771677d0b8SKris Buschelman 1781677d0b8SKris Buschelman /* Control parameters used by numeric factorization */ 1791677d0b8SKris Buschelman ierr = PetscOptionsReal("-mat_umfpack_pivot_tolerance","Control[UMFPACK_PIVOT_TOLERANCE]","None",lu->Control[UMFPACK_PIVOT_TOLERANCE],&lu->Control[UMFPACK_PIVOT_TOLERANCE],PETSC_NULL);CHKERRQ(ierr); 1801677d0b8SKris Buschelman ierr = PetscOptionsReal("-mat_umfpack_alloc_init","Control[UMFPACK_ALLOC_INIT]","None",lu->Control[UMFPACK_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],PETSC_NULL);CHKERRQ(ierr); 1811677d0b8SKris Buschelman 1821677d0b8SKris Buschelman /* Control parameters used by solve */ 1831677d0b8SKris Buschelman ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr); 1841677d0b8SKris Buschelman 1851677d0b8SKris Buschelman /* use Petsc mat ordering (notice size is for the transpose) */ 1861677d0b8SKris Buschelman ierr = PetscOptionsHasName(PETSC_NULL,"-pc_lu_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr); 1871677d0b8SKris Buschelman if (lu->PetscMatOdering) { 1881677d0b8SKris Buschelman ierr = ISGetIndices(c,&ca);CHKERRQ(ierr); 1891677d0b8SKris Buschelman ierr = PetscMalloc(A->m*sizeof(int),&lu->perm_c);CHKERRQ(ierr); 1901677d0b8SKris Buschelman ierr = PetscMemcpy(lu->perm_c,ca,A->m*sizeof(int));CHKERRQ(ierr); 1911677d0b8SKris Buschelman ierr = ISRestoreIndices(c,&ca);CHKERRQ(ierr); 1921677d0b8SKris Buschelman } 1931677d0b8SKris Buschelman PetscOptionsEnd(); 1941677d0b8SKris Buschelman 1951677d0b8SKris Buschelman /* print the control parameters */ 1961677d0b8SKris Buschelman if( lu->Control[UMFPACK_PRL] > 1 ) umfpack_di_report_control (lu->Control); 1971677d0b8SKris Buschelman 1981677d0b8SKris Buschelman /* symbolic factorization of A' */ 1991677d0b8SKris Buschelman /* ---------------------------------------------------------------------- */ 2000801d200SKris Buschelman status = umfpack_di_qsymbolic(n,m,ai,aj,PETSC_NULL,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info) ; 2011677d0b8SKris Buschelman if (status < 0){ 2021677d0b8SKris Buschelman umfpack_di_report_info(lu->Control, lu->Info) ; 2031677d0b8SKris Buschelman umfpack_di_report_status(lu->Control, status) ; 2041677d0b8SKris Buschelman SETERRQ(1,"umfpack_di_symbolic failed"); 2051677d0b8SKris Buschelman } 2061677d0b8SKris Buschelman /* report sumbolic factorization of A' when Control[PRL] > 3 */ 2071677d0b8SKris Buschelman (void) umfpack_di_report_symbolic(lu->Symbolic, lu->Control) ; 2081677d0b8SKris Buschelman 2091677d0b8SKris Buschelman lu->flg = DIFFERENT_NONZERO_PATTERN; 2101677d0b8SKris Buschelman lu->ai = ai; 2111677d0b8SKris Buschelman lu->aj = aj; 2121677d0b8SKris Buschelman lu->av = av; 213e1b4c3a1SKris Buschelman lu->CleanUpUMFPACK = PETSC_TRUE; 214e1b4c3a1SKris Buschelman *F = B; 2151677d0b8SKris Buschelman PetscFunctionReturn(0); 2161677d0b8SKris Buschelman } 2171677d0b8SKris Buschelman 2181677d0b8SKris Buschelman #undef __FUNCT__ 219f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_UMFPACK" 220f0c56d0fSKris Buschelman int MatAssemblyEnd_UMFPACK(Mat A,MatAssemblyType mode) { 221d844382dSKris Buschelman int ierr; 222f0c56d0fSKris Buschelman Mat_UMFPACK *lu=(Mat_UMFPACK*)(A->spptr); 223d844382dSKris Buschelman 2241677d0b8SKris Buschelman PetscFunctionBegin; 225d844382dSKris Buschelman ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 226d844382dSKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 227f0c56d0fSKris Buschelman A->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK; 2281677d0b8SKris Buschelman PetscFunctionReturn(0); 2291677d0b8SKris Buschelman } 2301677d0b8SKris Buschelman 23194b7f48cSBarry Smith /* used by -ksp_view */ 2321677d0b8SKris Buschelman #undef __FUNCT__ 233f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_UMFPACK" 234f0c56d0fSKris Buschelman int MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer) { 235f0c56d0fSKris Buschelman Mat_UMFPACK *lu= (Mat_UMFPACK*)A->spptr; 2361677d0b8SKris Buschelman int ierr; 237f0c56d0fSKris Buschelman 2381677d0b8SKris Buschelman PetscFunctionBegin; 2391677d0b8SKris Buschelman /* check if matrix is UMFPACK type */ 240f0c56d0fSKris Buschelman if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0); 2411677d0b8SKris Buschelman 2421677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr); 2431677d0b8SKris Buschelman /* Control parameters used by reporting routiones */ 2441677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr); 2451677d0b8SKris Buschelman 2461677d0b8SKris Buschelman /* Control parameters used by symbolic factorization */ 2471677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr); 2481677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr); 2491677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr); 2501677d0b8SKris Buschelman 2511677d0b8SKris Buschelman /* Control parameters used by numeric factorization */ 2521677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr); 253*59698cb0SHong Zhang 2541677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr); 2551677d0b8SKris Buschelman 2561677d0b8SKris Buschelman /* Control parameters used by solve */ 2571677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr); 2581677d0b8SKris Buschelman 2591677d0b8SKris Buschelman /* mat ordering */ 2601677d0b8SKris Buschelman if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer," UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr); 2611677d0b8SKris Buschelman 2621677d0b8SKris Buschelman PetscFunctionReturn(0); 2631677d0b8SKris Buschelman } 2641677d0b8SKris Buschelman 265d844382dSKris Buschelman #undef __FUNCT__ 266f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_UMFPACK" 267f0c56d0fSKris Buschelman int MatView_UMFPACK(Mat A,PetscViewer viewer) { 268d844382dSKris Buschelman int ierr; 269d844382dSKris Buschelman PetscTruth isascii; 270d844382dSKris Buschelman PetscViewerFormat format; 271f0c56d0fSKris Buschelman Mat_UMFPACK *lu=(Mat_UMFPACK*)(A->spptr); 272d844382dSKris Buschelman 273d844382dSKris Buschelman PetscFunctionBegin; 274d844382dSKris Buschelman ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr); 275d844382dSKris Buschelman 276d844382dSKris Buschelman ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 277d844382dSKris Buschelman if (isascii) { 278d844382dSKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 279d844382dSKris Buschelman if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 280f0c56d0fSKris Buschelman ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr); 281d844382dSKris Buschelman } 282d844382dSKris Buschelman } 283d844382dSKris Buschelman PetscFunctionReturn(0); 284d844382dSKris Buschelman } 285d844382dSKris Buschelman 286d844382dSKris Buschelman EXTERN_C_BEGIN 287d844382dSKris Buschelman #undef __FUNCT__ 288d844382dSKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_UMFPACK" 2898e9aea5cSBarry Smith int MatConvert_SeqAIJ_UMFPACK(Mat A,const MatType type,Mat *newmat) { 290d844382dSKris Buschelman /* This routine is only called to convert to MATUMFPACK */ 291d844382dSKris Buschelman /* from MATSEQAIJ, so we will ignore 'MatType type'. */ 292d844382dSKris Buschelman int ierr; 293d844382dSKris Buschelman Mat B=*newmat; 294f0c56d0fSKris Buschelman Mat_UMFPACK *lu; 295d844382dSKris Buschelman 296d844382dSKris Buschelman PetscFunctionBegin; 297d844382dSKris Buschelman if (B != A) { 298d844382dSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 299d844382dSKris Buschelman } 300d844382dSKris Buschelman 301f0c56d0fSKris Buschelman ierr = PetscNew(Mat_UMFPACK,&lu);CHKERRQ(ierr); 302f0c56d0fSKris Buschelman lu->MatDuplicate = A->ops->duplicate; 303d844382dSKris Buschelman lu->MatView = A->ops->view; 304d844382dSKris Buschelman lu->MatAssemblyEnd = A->ops->assemblyend; 305d844382dSKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 306d844382dSKris Buschelman lu->MatDestroy = A->ops->destroy; 307d844382dSKris Buschelman lu->CleanUpUMFPACK = PETSC_FALSE; 308d844382dSKris Buschelman 309d844382dSKris Buschelman B->spptr = (void*)lu; 310f0c56d0fSKris Buschelman B->ops->duplicate = MatDuplicate_UMFPACK; 311f0c56d0fSKris Buschelman B->ops->view = MatView_UMFPACK; 312f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_UMFPACK; 313f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK; 314f0c56d0fSKris Buschelman B->ops->destroy = MatDestroy_UMFPACK; 315d844382dSKris Buschelman 316d844382dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_umfpack_C", 317d844382dSKris Buschelman "MatConvert_SeqAIJ_UMFPACK",MatConvert_SeqAIJ_UMFPACK);CHKERRQ(ierr); 318d844382dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_umfpack_seqaij_C", 319d844382dSKris Buschelman "MatConvert_UMFPACK_SeqAIJ",MatConvert_UMFPACK_SeqAIJ);CHKERRQ(ierr); 320d844382dSKris Buschelman PetscLogInfo(0,"Using UMFPACK for SeqAIJ LU factorization and solves."); 321d844382dSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATUMFPACK);CHKERRQ(ierr); 322d844382dSKris Buschelman *newmat = B; 323d844382dSKris Buschelman PetscFunctionReturn(0); 324d844382dSKris Buschelman } 325d844382dSKris Buschelman EXTERN_C_END 326d844382dSKris Buschelman 327f0c56d0fSKris Buschelman #undef __FUNCT__ 328f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_UMFPACK" 329f0c56d0fSKris Buschelman int MatDuplicate_UMFPACK(Mat A, MatDuplicateOption op, Mat *M) { 330f0c56d0fSKris Buschelman int ierr; 3318f340917SKris Buschelman Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr; 3328f340917SKris Buschelman 333f0c56d0fSKris Buschelman PetscFunctionBegin; 3348f340917SKris Buschelman ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 3353f953163SKris Buschelman ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_UMFPACK));CHKERRQ(ierr); 336f0c56d0fSKris Buschelman PetscFunctionReturn(0); 337f0c56d0fSKris Buschelman } 338f0c56d0fSKris Buschelman 3392f71e704SKris Buschelman /*MC 340fafad747SKris Buschelman MATUMFPACK - MATUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices 3412f71e704SKris Buschelman via the external package UMFPACK. 3422f71e704SKris Buschelman 3432f71e704SKris Buschelman If UMFPACK is installed (see the manual for 3442f71e704SKris Buschelman instructions on how to declare the existence of external packages), 3452f71e704SKris Buschelman a matrix type can be constructed which invokes UMFPACK solvers. 3462f71e704SKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,UMFPACK). 3472f71e704SKris Buschelman This matrix type is only supported for double precision real. 3482f71e704SKris Buschelman 3492f71e704SKris Buschelman This matrix inherits from MATSEQAIJ. As a result, MatSeqAIJSetPreallocation is 35028b08bd3SKris Buschelman supported for this matrix type. One can also call MatConvert for an inplace conversion to or from 35128b08bd3SKris Buschelman the MATSEQAIJ type without data copy. 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: 3570bad9183SKris Buschelman + -mat_type umfpack - sets the matrix type to "umfpack" during a call to MatSetFromOptions() 3582f71e704SKris Buschelman . -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL] 3592f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL] 3602f71e704SKris Buschelman . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE] 3612f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE] 3622f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT] 3632f71e704SKris Buschelman - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP] 3642f71e704SKris Buschelman 3652f71e704SKris Buschelman Level: beginner 3662f71e704SKris Buschelman 3672f71e704SKris Buschelman .seealso: PCLU 3682f71e704SKris Buschelman M*/ 3692f71e704SKris Buschelman 3701677d0b8SKris Buschelman EXTERN_C_BEGIN 3711677d0b8SKris Buschelman #undef __FUNCT__ 372f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_UMFPACK" 373f0c56d0fSKris Buschelman int MatCreate_UMFPACK(Mat A) { 3741677d0b8SKris Buschelman int ierr; 3751677d0b8SKris Buschelman 3761677d0b8SKris Buschelman PetscFunctionBegin; 3775441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqAIJ and UMFPACK types */ 3785441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATUMFPACK);CHKERRQ(ierr); 3791677d0b8SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 380d844382dSKris Buschelman ierr = MatConvert_SeqAIJ_UMFPACK(A,MATUMFPACK,&A);CHKERRQ(ierr); 3811677d0b8SKris Buschelman PetscFunctionReturn(0); 3821677d0b8SKris Buschelman } 3831677d0b8SKris Buschelman EXTERN_C_END 384