11677d0b8SKris Buschelman /*$Id: umfpack.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/ 21677d0b8SKris Buschelman 31677d0b8SKris Buschelman /* 41677d0b8SKris Buschelman Provides an interface to the UMFPACK 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 */ 221677d0b8SKris Buschelman int (*MatView)(Mat,PetscViewer); 231677d0b8SKris Buschelman int (*MatAssemblyEnd)(Mat,MatAssemblyType); 241677d0b8SKris Buschelman int (*MatDestroy)(Mat); 251677d0b8SKris Buschelman 261677d0b8SKris Buschelman /* Flag to clean up UMFPACK objects during Destroy */ 271677d0b8SKris Buschelman PetscTruth CleanUpUMFPACK; 281677d0b8SKris Buschelman } Mat_SeqAIJ_UMFPACK; 291677d0b8SKris Buschelman 301677d0b8SKris Buschelman #undef __FUNCT__ 311677d0b8SKris Buschelman #define __FUNCT__ "MatDestroy_SeqAIJ_UMFPACK" 321677d0b8SKris Buschelman int MatDestroy_SeqAIJ_UMFPACK(Mat A) 331677d0b8SKris Buschelman { 341677d0b8SKris Buschelman Mat_SeqAIJ_UMFPACK *lu = (Mat_SeqAIJ_UMFPACK*)A->spptr; 351677d0b8SKris Buschelman int ierr,(*destroy)(Mat); 361677d0b8SKris Buschelman 371677d0b8SKris Buschelman PetscFunctionBegin; 38*fb731535SKris Buschelman if (lu->CleanUpUMFPACK) { 391677d0b8SKris Buschelman umfpack_di_free_symbolic(&lu->Symbolic) ; 401677d0b8SKris Buschelman umfpack_di_free_numeric(&lu->Numeric) ; 411677d0b8SKris Buschelman ierr = PetscFree(lu->Wi);CHKERRQ(ierr); 421677d0b8SKris Buschelman ierr = PetscFree(lu->W);CHKERRQ(ierr); 431677d0b8SKris Buschelman 441677d0b8SKris Buschelman if (lu->PetscMatOdering) { 451677d0b8SKris Buschelman ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 461677d0b8SKris Buschelman } 471677d0b8SKris Buschelman } 481677d0b8SKris Buschelman 491677d0b8SKris Buschelman destroy = lu->MatDestroy; 501677d0b8SKris Buschelman ierr = PetscFree(lu);CHKERRQ(ierr); 511677d0b8SKris Buschelman ierr = (*destroy)(A);CHKERRQ(ierr); 521677d0b8SKris Buschelman PetscFunctionReturn(0); 531677d0b8SKris Buschelman } 541677d0b8SKris Buschelman 551677d0b8SKris Buschelman #undef __FUNCT__ 561677d0b8SKris Buschelman #define __FUNCT__ "MatView_SeqAIJ_UMFPACK" 571677d0b8SKris Buschelman int MatView_SeqAIJ_UMFPACK(Mat A,PetscViewer viewer) 581677d0b8SKris Buschelman { 591677d0b8SKris Buschelman int ierr; 601677d0b8SKris Buschelman PetscTruth isascii; 611677d0b8SKris Buschelman PetscViewerFormat format; 621677d0b8SKris Buschelman Mat_SeqAIJ_UMFPACK *lu=(Mat_SeqAIJ_UMFPACK*)(A->spptr); 631677d0b8SKris Buschelman 641677d0b8SKris Buschelman PetscFunctionBegin; 651677d0b8SKris Buschelman ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr); 661677d0b8SKris Buschelman 671677d0b8SKris Buschelman ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 681677d0b8SKris Buschelman if (isascii) { 691677d0b8SKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 701677d0b8SKris Buschelman if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 711677d0b8SKris Buschelman ierr = MatSeqAIJFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr); 721677d0b8SKris Buschelman } 731677d0b8SKris Buschelman } 741677d0b8SKris Buschelman PetscFunctionReturn(0); 751677d0b8SKris Buschelman } 761677d0b8SKris Buschelman 771677d0b8SKris Buschelman #undef __FUNCT__ 781677d0b8SKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SeqAIJ_UMFPACK" 791677d0b8SKris Buschelman int MatAssemblyEnd_SeqAIJ_UMFPACK(Mat A,MatAssemblyType mode) { 801677d0b8SKris Buschelman int ierr; 811677d0b8SKris Buschelman Mat_SeqAIJ_UMFPACK *lu=(Mat_SeqAIJ_UMFPACK*)(A->spptr); 821677d0b8SKris Buschelman 831677d0b8SKris Buschelman PetscFunctionBegin; 841677d0b8SKris Buschelman ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 851677d0b8SKris Buschelman ierr = MatUseUMFPACK_SeqAIJ(A);CHKERRQ(ierr); 861677d0b8SKris Buschelman PetscFunctionReturn(0); 871677d0b8SKris Buschelman } 881677d0b8SKris Buschelman 891677d0b8SKris Buschelman #undef __FUNCT__ 901677d0b8SKris Buschelman #define __FUNCT__ "MatSolve_SeqAIJ_UMFPACK" 911677d0b8SKris Buschelman int MatSolve_SeqAIJ_UMFPACK(Mat A,Vec b,Vec x) 921677d0b8SKris Buschelman { 931677d0b8SKris Buschelman Mat_SeqAIJ_UMFPACK *lu = (Mat_SeqAIJ_UMFPACK*)A->spptr; 941677d0b8SKris Buschelman PetscScalar *av=lu->av,*ba,*xa; 951677d0b8SKris Buschelman int ierr,*ai=lu->ai,*aj=lu->aj,status; 961677d0b8SKris Buschelman 971677d0b8SKris Buschelman PetscFunctionBegin; 981677d0b8SKris Buschelman /* solve Ax = b by umfpack_di_wsolve */ 991677d0b8SKris Buschelman /* ----------------------------------*/ 1001677d0b8SKris Buschelman ierr = VecGetArray(b,&ba); 1011677d0b8SKris Buschelman ierr = VecGetArray(x,&xa); 1021677d0b8SKris Buschelman 1031677d0b8SKris Buschelman status = umfpack_di_wsolve(UMFPACK_At,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W); 1041677d0b8SKris Buschelman umfpack_di_report_info(lu->Control, lu->Info); 1051677d0b8SKris Buschelman if (status < 0){ 1061677d0b8SKris Buschelman umfpack_di_report_status(lu->Control, status) ; 1071677d0b8SKris Buschelman SETERRQ(1,"umfpack_di_wsolve failed") ; 1081677d0b8SKris Buschelman } 1091677d0b8SKris Buschelman 1101677d0b8SKris Buschelman ierr = VecRestoreArray(b,&ba); 1111677d0b8SKris Buschelman ierr = VecRestoreArray(x,&xa); 1121677d0b8SKris Buschelman PetscFunctionReturn(0); 1131677d0b8SKris Buschelman } 1141677d0b8SKris Buschelman 1151677d0b8SKris Buschelman #undef __FUNCT__ 1161677d0b8SKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_UMFPACK" 1171677d0b8SKris Buschelman int MatLUFactorNumeric_SeqAIJ_UMFPACK(Mat A,Mat *F) 1181677d0b8SKris Buschelman { 1191677d0b8SKris Buschelman Mat_SeqAIJ_UMFPACK *lu = (Mat_SeqAIJ_UMFPACK*)(*F)->spptr; 1201677d0b8SKris Buschelman int *ai=lu->ai,*aj=lu->aj,m=A->m,status,ierr; 1211677d0b8SKris Buschelman PetscScalar *av=lu->av; 1221677d0b8SKris Buschelman 1231677d0b8SKris Buschelman PetscFunctionBegin; 1241677d0b8SKris Buschelman /* numeric factorization of A' */ 1251677d0b8SKris Buschelman /* ----------------------------*/ 1261677d0b8SKris Buschelman status = umfpack_di_numeric (ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info) ; 1271677d0b8SKris Buschelman if (status < 0) SETERRQ(1,"umfpack_di_numeric failed"); 1281677d0b8SKris Buschelman /* report numeric factorization of A' when Control[PRL] > 3 */ 1291677d0b8SKris Buschelman (void) umfpack_di_report_numeric (lu->Numeric, lu->Control) ; 1301677d0b8SKris Buschelman 1311677d0b8SKris Buschelman if (lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization */ 1321677d0b8SKris Buschelman /* allocate working space to be used by Solve */ 1331677d0b8SKris Buschelman ierr = PetscMalloc(m * sizeof(int), &lu->Wi);CHKERRQ(ierr); 1341677d0b8SKris Buschelman ierr = PetscMalloc(5*m * sizeof(double), &lu->W);CHKERRQ(ierr); 1351677d0b8SKris Buschelman 1361677d0b8SKris Buschelman lu->flg = SAME_NONZERO_PATTERN; 1371677d0b8SKris Buschelman } 1381677d0b8SKris Buschelman 1391677d0b8SKris Buschelman PetscFunctionReturn(0); 1401677d0b8SKris Buschelman } 1411677d0b8SKris Buschelman 1421677d0b8SKris Buschelman /* 1431677d0b8SKris Buschelman Note the r permutation is ignored 1441677d0b8SKris Buschelman */ 1451677d0b8SKris Buschelman #undef __FUNCT__ 1461677d0b8SKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ_UMFPACK" 1471677d0b8SKris Buschelman int MatLUFactorSymbolic_SeqAIJ_UMFPACK(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) 1481677d0b8SKris Buschelman { 1491677d0b8SKris Buschelman Mat_SeqAIJ *mat=(Mat_SeqAIJ*)A->data; 1501677d0b8SKris Buschelman Mat_SeqAIJ_UMFPACK *lu; 1511677d0b8SKris Buschelman int ierr,m=A->m,n=A->n,*ai=mat->i,*aj=mat->j,status,*ca; 1521677d0b8SKris Buschelman PetscScalar *av=mat->a; 1531677d0b8SKris Buschelman 1541677d0b8SKris Buschelman PetscFunctionBegin; 1551677d0b8SKris Buschelman /* Create the factorization matrix F */ 1561677d0b8SKris Buschelman ierr = MatCreateSeqAIJ(A->comm,m,n,PETSC_NULL,PETSC_NULL,F);CHKERRQ(ierr); 1571677d0b8SKris Buschelman 1581677d0b8SKris Buschelman (*F)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_UMFPACK; 1591677d0b8SKris Buschelman (*F)->ops->solve = MatSolve_SeqAIJ_UMFPACK; 1601677d0b8SKris Buschelman (*F)->ops->destroy = MatDestroy_SeqAIJ_UMFPACK; 1611677d0b8SKris Buschelman (*F)->factor = FACTOR_LU; 1621677d0b8SKris Buschelman (*F)->assembled = PETSC_TRUE; /* required by -sles_view */ 1631677d0b8SKris Buschelman 1641677d0b8SKris Buschelman ierr = PetscNew(Mat_SeqAIJ_UMFPACK,&lu);CHKERRQ(ierr); 1651677d0b8SKris Buschelman (*F)->spptr = (void*)lu; 1661677d0b8SKris Buschelman 1671677d0b8SKris Buschelman /* initializations */ 1681677d0b8SKris Buschelman /* ------------------------------------------------*/ 1691677d0b8SKris Buschelman /* get the default control parameters */ 1701677d0b8SKris Buschelman umfpack_di_defaults (lu->Control) ; 1711677d0b8SKris Buschelman lu->perm_c = PETSC_NULL; /* use defaul UMFPACK col permutation */ 1721677d0b8SKris Buschelman 1731677d0b8SKris Buschelman ierr = PetscOptionsBegin(A->comm,A->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr); 1741677d0b8SKris Buschelman /* Control parameters used by reporting routiones */ 1751677d0b8SKris Buschelman ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr); 1761677d0b8SKris Buschelman 1771677d0b8SKris Buschelman /* Control parameters for symbolic factorization */ 1781677d0b8SKris 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); 1791677d0b8SKris 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); 1801677d0b8SKris 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); 1811677d0b8SKris Buschelman 1821677d0b8SKris Buschelman /* Control parameters used by numeric factorization */ 1831677d0b8SKris 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); 184*fb731535SKris Buschelman #if !defined(PETSC_HAVE_UMFPACK_42_OR_NEWER) 1851677d0b8SKris Buschelman ierr = PetscOptionsReal("-mat_umfpack_relaxed_amalgamation","Control[UMFPACK_RELAXED_AMALGAMATION]","None",lu->Control[UMFPACK_RELAXED_AMALGAMATION],&lu->Control[UMFPACK_RELAXED_AMALGAMATION],PETSC_NULL);CHKERRQ(ierr); 1861677d0b8SKris Buschelman ierr = PetscOptionsReal("-mat_umfpack_relaxed2_amalgamation","Control[UMFPACK_RELAXED2_AMALGAMATION]","None",lu->Control[UMFPACK_RELAXED2_AMALGAMATION],&lu->Control[UMFPACK_RELAXED2_AMALGAMATION],PETSC_NULL);CHKERRQ(ierr); 1871677d0b8SKris Buschelman ierr = PetscOptionsReal("-mat_umfpack_relaxed3_amalgamation","Control[UMFPACK_RELAXED3_AMALGAMATION]","None",lu->Control[UMFPACK_RELAXED3_AMALGAMATION],&lu->Control[UMFPACK_RELAXED3_AMALGAMATION],PETSC_NULL);CHKERRQ(ierr); 188*fb731535SKris Buschelman #endif 1891677d0b8SKris 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); 1901677d0b8SKris Buschelman 1911677d0b8SKris Buschelman /* Control parameters used by solve */ 1921677d0b8SKris Buschelman ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr); 1931677d0b8SKris Buschelman 1941677d0b8SKris Buschelman /* use Petsc mat ordering (notice size is for the transpose) */ 1951677d0b8SKris Buschelman ierr = PetscOptionsHasName(PETSC_NULL,"-pc_lu_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr); 1961677d0b8SKris Buschelman if (lu->PetscMatOdering) { 1971677d0b8SKris Buschelman ierr = ISGetIndices(c,&ca);CHKERRQ(ierr); 1981677d0b8SKris Buschelman ierr = PetscMalloc(A->m*sizeof(int),&lu->perm_c);CHKERRQ(ierr); 1991677d0b8SKris Buschelman ierr = PetscMemcpy(lu->perm_c,ca,A->m*sizeof(int));CHKERRQ(ierr); 2001677d0b8SKris Buschelman ierr = ISRestoreIndices(c,&ca);CHKERRQ(ierr); 2011677d0b8SKris Buschelman } 2021677d0b8SKris Buschelman PetscOptionsEnd(); 2031677d0b8SKris Buschelman 2041677d0b8SKris Buschelman /* print the control parameters */ 2051677d0b8SKris Buschelman if( lu->Control[UMFPACK_PRL] > 1 ) umfpack_di_report_control (lu->Control); 2061677d0b8SKris Buschelman 2071677d0b8SKris Buschelman /* symbolic factorization of A' */ 2081677d0b8SKris Buschelman /* ---------------------------------------------------------------------- */ 2091677d0b8SKris Buschelman status = umfpack_di_qsymbolic(n,m,ai,aj,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info) ; 2101677d0b8SKris Buschelman if (status < 0){ 2111677d0b8SKris Buschelman umfpack_di_report_info(lu->Control, lu->Info) ; 2121677d0b8SKris Buschelman umfpack_di_report_status(lu->Control, status) ; 2131677d0b8SKris Buschelman SETERRQ(1,"umfpack_di_symbolic failed"); 2141677d0b8SKris Buschelman } 2151677d0b8SKris Buschelman /* report sumbolic factorization of A' when Control[PRL] > 3 */ 2161677d0b8SKris Buschelman (void) umfpack_di_report_symbolic(lu->Symbolic, lu->Control) ; 2171677d0b8SKris Buschelman 2181677d0b8SKris Buschelman lu->flg = DIFFERENT_NONZERO_PATTERN; 2191677d0b8SKris Buschelman lu->ai = ai; 2201677d0b8SKris Buschelman lu->aj = aj; 2211677d0b8SKris Buschelman lu->av = av; 2221677d0b8SKris Buschelman PetscFunctionReturn(0); 2231677d0b8SKris Buschelman } 2241677d0b8SKris Buschelman 2251677d0b8SKris Buschelman #undef __FUNCT__ 2261677d0b8SKris Buschelman #define __FUNCT__ "MatUseUMFPACK_SeqAIJ" 2271677d0b8SKris Buschelman int MatUseUMFPACK_SeqAIJ(Mat A) 2281677d0b8SKris Buschelman { 2291677d0b8SKris Buschelman PetscFunctionBegin; 2301677d0b8SKris Buschelman A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ_UMFPACK; 2311677d0b8SKris Buschelman PetscFunctionReturn(0); 2321677d0b8SKris Buschelman } 2331677d0b8SKris Buschelman 2341677d0b8SKris Buschelman /* used by -sles_view */ 2351677d0b8SKris Buschelman #undef __FUNCT__ 2361677d0b8SKris Buschelman #define __FUNCT__ "MatSeqAIJFactorInfo_UMFPACK" 2371677d0b8SKris Buschelman int MatSeqAIJFactorInfo_UMFPACK(Mat A,PetscViewer viewer) 2381677d0b8SKris Buschelman { 2391677d0b8SKris Buschelman Mat_SeqAIJ_UMFPACK *lu= (Mat_SeqAIJ_UMFPACK*)A->spptr; 2401677d0b8SKris Buschelman int ierr; 2411677d0b8SKris Buschelman PetscFunctionBegin; 2421677d0b8SKris Buschelman /* check if matrix is UMFPACK type */ 2431677d0b8SKris Buschelman if (A->ops->solve != MatSolve_SeqAIJ_UMFPACK) PetscFunctionReturn(0); 2441677d0b8SKris Buschelman 2451677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr); 2461677d0b8SKris Buschelman /* Control parameters used by reporting routiones */ 2471677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr); 2481677d0b8SKris Buschelman 2491677d0b8SKris Buschelman /* Control parameters used by symbolic factorization */ 2501677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr); 2511677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr); 2521677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr); 2531677d0b8SKris Buschelman 2541677d0b8SKris Buschelman /* Control parameters used by numeric factorization */ 2551677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr); 256*fb731535SKris Buschelman #if !defined(PETSC_HAVE_UMFPACK_42_OR_NEWER) 2571677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_RELAXED_AMALGAMATION]: %g\n",lu->Control[UMFPACK_RELAXED_AMALGAMATION]);CHKERRQ(ierr); 2581677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_RELAXED2_AMALGAMATION]: %g\n",lu->Control[UMFPACK_RELAXED2_AMALGAMATION]);CHKERRQ(ierr); 2591677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_RELAXED3_AMALGAMATION]: %g\n",lu->Control[UMFPACK_RELAXED3_AMALGAMATION]);CHKERRQ(ierr); 260*fb731535SKris Buschelman #endif 2611677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr); 2621677d0b8SKris Buschelman 2631677d0b8SKris Buschelman /* Control parameters used by solve */ 2641677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr); 2651677d0b8SKris Buschelman 2661677d0b8SKris Buschelman /* mat ordering */ 2671677d0b8SKris Buschelman if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer," UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr); 2681677d0b8SKris Buschelman 2691677d0b8SKris Buschelman PetscFunctionReturn(0); 2701677d0b8SKris Buschelman } 2711677d0b8SKris Buschelman 2721677d0b8SKris Buschelman EXTERN_C_BEGIN 2731677d0b8SKris Buschelman #undef __FUNCT__ 2741677d0b8SKris Buschelman #define __FUNCT__ "MatCreate_SeqAIJ_UMFPACK" 2751677d0b8SKris Buschelman int MatCreate_SeqAIJ_UMFPACK(Mat A) { 2761677d0b8SKris Buschelman int ierr; 2771677d0b8SKris Buschelman Mat_SeqAIJ_UMFPACK *lu; 2781677d0b8SKris Buschelman 2791677d0b8SKris Buschelman PetscFunctionBegin; 2801677d0b8SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 2811677d0b8SKris Buschelman ierr = MatUseUMFPACK_SeqAIJ(A);CHKERRQ(ierr); 2821677d0b8SKris Buschelman 2831677d0b8SKris Buschelman ierr = PetscNew(Mat_SeqAIJ_UMFPACK,&lu);CHKERRQ(ierr); 2841677d0b8SKris Buschelman lu->MatView = A->ops->view; 2851677d0b8SKris Buschelman lu->MatAssemblyEnd = A->ops->assemblyend; 2861677d0b8SKris Buschelman lu->MatDestroy = A->ops->destroy; 2871677d0b8SKris Buschelman lu->CleanUpUMFPACK = PETSC_FALSE; 2881677d0b8SKris Buschelman A->spptr = (void*)lu; 2891677d0b8SKris Buschelman A->ops->view = MatView_SeqAIJ_UMFPACK; 2901677d0b8SKris Buschelman A->ops->assemblyend = MatAssemblyEnd_SeqAIJ_UMFPACK; 2911677d0b8SKris Buschelman A->ops->destroy = MatDestroy_SeqAIJ_UMFPACK; 2921677d0b8SKris Buschelman PetscFunctionReturn(0); 2931677d0b8SKris Buschelman } 2941677d0b8SKris Buschelman EXTERN_C_END 295