11677d0b8SKris Buschelman /*$Id: umfpack.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/ 21677d0b8SKris Buschelman 31677d0b8SKris Buschelman /* 459698cb0SHong 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; 147*0f6efc10SHong Zhang int ierr,m=A->m,n=A->n,*ai=mat->i,*aj=mat->j,status,*ra,idx; 1481677d0b8SKris Buschelman PetscScalar *av=mat->a; 149*0f6efc10SHong Zhang const char *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC","2BY2"}, 150*0f6efc10SHong Zhang *scale[]={"NONE","SUM","MAX"}; 151*0f6efc10SHong Zhang PetscTruth flg; 1521677d0b8SKris Buschelman 1531677d0b8SKris Buschelman PetscFunctionBegin; 1541677d0b8SKris Buschelman /* Create the factorization matrix F */ 155e1b4c3a1SKris Buschelman ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,m,n,&B);CHKERRQ(ierr); 156be5d1d56SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 157e1b4c3a1SKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 1581677d0b8SKris Buschelman 159f0c56d0fSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK; 160f0c56d0fSKris Buschelman B->ops->solve = MatSolve_UMFPACK; 161e1b4c3a1SKris Buschelman B->factor = FACTOR_LU; 16294b7f48cSBarry Smith B->assembled = PETSC_TRUE; /* required by -ksp_view */ 1631677d0b8SKris Buschelman 164f0c56d0fSKris Buschelman lu = (Mat_UMFPACK*)(B->spptr); 1651677d0b8SKris Buschelman 1661677d0b8SKris Buschelman /* initializations */ 1671677d0b8SKris Buschelman /* ------------------------------------------------*/ 1681677d0b8SKris Buschelman /* get the default control parameters */ 1691677d0b8SKris Buschelman umfpack_di_defaults (lu->Control) ; 1701677d0b8SKris Buschelman lu->perm_c = PETSC_NULL; /* use defaul UMFPACK col permutation */ 171*0f6efc10SHong Zhang lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */ 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 */ 178*0f6efc10SHong Zhang ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,4,strategy[0],&idx,&flg);CHKERRQ(ierr); 179*0f6efc10SHong Zhang if (flg) { 180*0f6efc10SHong Zhang switch (idx){ 181*0f6efc10SHong Zhang case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break; 182*0f6efc10SHong Zhang case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break; 183*0f6efc10SHong Zhang case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break; 184*0f6efc10SHong Zhang case 3: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_2BY2; break; 185*0f6efc10SHong Zhang } 186*0f6efc10SHong Zhang } 1871677d0b8SKris 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); 1881677d0b8SKris 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); 189*0f6efc10SHong Zhang ierr = PetscOptionsReal("-mat_umfpack_amd_dense","Control[UMFPACK_AMD_DENSE]","None",lu->Control[UMFPACK_AMD_DENSE],&lu->Control[UMFPACK_AMD_DENSE],PETSC_NULL);CHKERRQ(ierr); 1901677d0b8SKris 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); 191*0f6efc10SHong Zhang ierr = PetscOptionsReal("-mat_umfpack_2by2_tolerance","Control[UMFPACK_2BY2_TOLERANCE]","None",lu->Control[UMFPACK_2BY2_TOLERANCE],&lu->Control[UMFPACK_2BY2_TOLERANCE],PETSC_NULL);CHKERRQ(ierr); 192*0f6efc10SHong Zhang ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);CHKERRQ(ierr); 193*0f6efc10SHong Zhang ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);CHKERRQ(ierr); 1941677d0b8SKris Buschelman 1951677d0b8SKris Buschelman /* Control parameters used by numeric factorization */ 1961677d0b8SKris 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); 197*0f6efc10SHong Zhang 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],PETSC_NULL);CHKERRQ(ierr); 198*0f6efc10SHong Zhang ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr); 199*0f6efc10SHong Zhang if (flg) { 200*0f6efc10SHong Zhang switch (idx){ 201*0f6efc10SHong Zhang case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break; 202*0f6efc10SHong Zhang case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break; 203*0f6efc10SHong Zhang case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break; 204*0f6efc10SHong Zhang } 205*0f6efc10SHong Zhang } 2061677d0b8SKris 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); 207*0f6efc10SHong Zhang ierr = PetscOptionsReal("-mat_umfpack_front_alloc_init","Control[UMFPACK_FRONT_ALLOC_INIT]","None",lu->Control[UMFPACK_FRONT_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],PETSC_NULL);CHKERRQ(ierr); 208*0f6efc10SHong Zhang ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);CHKERRQ(ierr); 2091677d0b8SKris Buschelman 2101677d0b8SKris Buschelman /* Control parameters used by solve */ 2111677d0b8SKris Buschelman ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr); 2121677d0b8SKris Buschelman 213*0f6efc10SHong Zhang /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */ 2141677d0b8SKris Buschelman ierr = PetscOptionsHasName(PETSC_NULL,"-pc_lu_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr); 2151677d0b8SKris Buschelman if (lu->PetscMatOdering) { 216*0f6efc10SHong Zhang ierr = ISGetIndices(r,&ra);CHKERRQ(ierr); 2171677d0b8SKris Buschelman ierr = PetscMalloc(A->m*sizeof(int),&lu->perm_c);CHKERRQ(ierr); 218*0f6efc10SHong Zhang ierr = PetscMemcpy(lu->perm_c,ra,A->m*sizeof(int));CHKERRQ(ierr); 219*0f6efc10SHong Zhang ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr); 2201677d0b8SKris Buschelman } 2211677d0b8SKris Buschelman PetscOptionsEnd(); 2221677d0b8SKris Buschelman 2231677d0b8SKris Buschelman /* print the control parameters */ 2241677d0b8SKris Buschelman if( lu->Control[UMFPACK_PRL] > 1 ) umfpack_di_report_control (lu->Control); 2251677d0b8SKris Buschelman 2261677d0b8SKris Buschelman /* symbolic factorization of A' */ 2271677d0b8SKris Buschelman /* ---------------------------------------------------------------------- */ 228*0f6efc10SHong Zhang if (lu->PetscMatOdering) { /* use Petsc row ordering */ 229*0f6efc10SHong Zhang status = umfpack_di_qsymbolic(n,m,ai,aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info); 230*0f6efc10SHong Zhang } else { /* use Umfpack col ordering */ 231*0f6efc10SHong Zhang status = umfpack_di_symbolic(n,m,ai,aj,av,&lu->Symbolic,lu->Control,lu->Info); 232*0f6efc10SHong Zhang } 2331677d0b8SKris Buschelman if (status < 0){ 2341677d0b8SKris Buschelman umfpack_di_report_info(lu->Control, lu->Info) ; 2351677d0b8SKris Buschelman umfpack_di_report_status(lu->Control, status) ; 2361677d0b8SKris Buschelman SETERRQ(1,"umfpack_di_symbolic failed"); 2371677d0b8SKris Buschelman } 2381677d0b8SKris Buschelman /* report sumbolic factorization of A' when Control[PRL] > 3 */ 2391677d0b8SKris Buschelman (void) umfpack_di_report_symbolic(lu->Symbolic, lu->Control) ; 2401677d0b8SKris Buschelman 2411677d0b8SKris Buschelman lu->flg = DIFFERENT_NONZERO_PATTERN; 2421677d0b8SKris Buschelman lu->ai = ai; 2431677d0b8SKris Buschelman lu->aj = aj; 2441677d0b8SKris Buschelman lu->av = av; 245e1b4c3a1SKris Buschelman lu->CleanUpUMFPACK = PETSC_TRUE; 246e1b4c3a1SKris Buschelman *F = B; 2471677d0b8SKris Buschelman PetscFunctionReturn(0); 2481677d0b8SKris Buschelman } 2491677d0b8SKris Buschelman 2501677d0b8SKris Buschelman #undef __FUNCT__ 251f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_UMFPACK" 252f0c56d0fSKris Buschelman int MatAssemblyEnd_UMFPACK(Mat A,MatAssemblyType mode) { 253d844382dSKris Buschelman int ierr; 254f0c56d0fSKris Buschelman Mat_UMFPACK *lu=(Mat_UMFPACK*)(A->spptr); 255d844382dSKris Buschelman 2561677d0b8SKris Buschelman PetscFunctionBegin; 257d844382dSKris Buschelman ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 258d844382dSKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 259f0c56d0fSKris Buschelman A->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK; 2601677d0b8SKris Buschelman PetscFunctionReturn(0); 2611677d0b8SKris Buschelman } 2621677d0b8SKris Buschelman 26394b7f48cSBarry Smith /* used by -ksp_view */ 2641677d0b8SKris Buschelman #undef __FUNCT__ 265f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_UMFPACK" 266f0c56d0fSKris Buschelman int MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer) { 267f0c56d0fSKris Buschelman Mat_UMFPACK *lu= (Mat_UMFPACK*)A->spptr; 2681677d0b8SKris Buschelman int ierr; 269f0c56d0fSKris Buschelman 2701677d0b8SKris Buschelman PetscFunctionBegin; 2711677d0b8SKris Buschelman /* check if matrix is UMFPACK type */ 272f0c56d0fSKris Buschelman if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0); 2731677d0b8SKris Buschelman 2741677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr); 2751677d0b8SKris Buschelman /* Control parameters used by reporting routiones */ 2761677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr); 2771677d0b8SKris Buschelman 2781677d0b8SKris Buschelman /* Control parameters used by symbolic factorization */ 279*0f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr); 2801677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr); 2811677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr); 282*0f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr); 2831677d0b8SKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr); 284*0f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_2BY2_TOLERANCE]: %g\n",lu->Control[UMFPACK_2BY2_TOLERANCE]);CHKERRQ(ierr); 285*0f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr); 286*0f6efc10SHong 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); 290*0f6efc10SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr); 291*0f6efc10SHong 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); 293*0f6efc10SHong 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 */ 2991677d0b8SKris Buschelman if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer," UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr); 3001677d0b8SKris Buschelman 3011677d0b8SKris Buschelman PetscFunctionReturn(0); 3021677d0b8SKris Buschelman } 3031677d0b8SKris Buschelman 304d844382dSKris Buschelman #undef __FUNCT__ 305f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_UMFPACK" 306f0c56d0fSKris Buschelman int MatView_UMFPACK(Mat A,PetscViewer viewer) { 307d844382dSKris Buschelman int ierr; 308d844382dSKris Buschelman PetscTruth isascii; 309d844382dSKris Buschelman PetscViewerFormat format; 310f0c56d0fSKris Buschelman Mat_UMFPACK *lu=(Mat_UMFPACK*)(A->spptr); 311d844382dSKris Buschelman 312d844382dSKris Buschelman PetscFunctionBegin; 313d844382dSKris Buschelman ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr); 314d844382dSKris Buschelman 315d844382dSKris Buschelman ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 316d844382dSKris Buschelman if (isascii) { 317d844382dSKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 318d844382dSKris Buschelman if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 319f0c56d0fSKris Buschelman ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr); 320d844382dSKris Buschelman } 321d844382dSKris Buschelman } 322d844382dSKris Buschelman PetscFunctionReturn(0); 323d844382dSKris Buschelman } 324d844382dSKris Buschelman 325d844382dSKris Buschelman EXTERN_C_BEGIN 326d844382dSKris Buschelman #undef __FUNCT__ 327d844382dSKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_UMFPACK" 3288e9aea5cSBarry Smith int MatConvert_SeqAIJ_UMFPACK(Mat A,const MatType type,Mat *newmat) { 329d844382dSKris Buschelman /* This routine is only called to convert to MATUMFPACK */ 330d844382dSKris Buschelman /* from MATSEQAIJ, so we will ignore 'MatType type'. */ 331d844382dSKris Buschelman int ierr; 332d844382dSKris Buschelman Mat B=*newmat; 333f0c56d0fSKris Buschelman Mat_UMFPACK *lu; 334d844382dSKris Buschelman 335d844382dSKris Buschelman PetscFunctionBegin; 336d844382dSKris Buschelman if (B != A) { 337d844382dSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 338d844382dSKris Buschelman } 339d844382dSKris Buschelman 340f0c56d0fSKris Buschelman ierr = PetscNew(Mat_UMFPACK,&lu);CHKERRQ(ierr); 341f0c56d0fSKris Buschelman lu->MatDuplicate = A->ops->duplicate; 342d844382dSKris Buschelman lu->MatView = A->ops->view; 343d844382dSKris Buschelman lu->MatAssemblyEnd = A->ops->assemblyend; 344d844382dSKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 345d844382dSKris Buschelman lu->MatDestroy = A->ops->destroy; 346d844382dSKris Buschelman lu->CleanUpUMFPACK = PETSC_FALSE; 347d844382dSKris Buschelman 348d844382dSKris Buschelman B->spptr = (void*)lu; 349f0c56d0fSKris Buschelman B->ops->duplicate = MatDuplicate_UMFPACK; 350f0c56d0fSKris Buschelman B->ops->view = MatView_UMFPACK; 351f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_UMFPACK; 352f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK; 353f0c56d0fSKris Buschelman B->ops->destroy = MatDestroy_UMFPACK; 354d844382dSKris Buschelman 355d844382dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_umfpack_C", 356d844382dSKris Buschelman "MatConvert_SeqAIJ_UMFPACK",MatConvert_SeqAIJ_UMFPACK);CHKERRQ(ierr); 357d844382dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_umfpack_seqaij_C", 358d844382dSKris Buschelman "MatConvert_UMFPACK_SeqAIJ",MatConvert_UMFPACK_SeqAIJ);CHKERRQ(ierr); 359d844382dSKris Buschelman PetscLogInfo(0,"Using UMFPACK for SeqAIJ LU factorization and solves."); 360d844382dSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATUMFPACK);CHKERRQ(ierr); 361d844382dSKris Buschelman *newmat = B; 362d844382dSKris Buschelman PetscFunctionReturn(0); 363d844382dSKris Buschelman } 364d844382dSKris Buschelman EXTERN_C_END 365d844382dSKris Buschelman 366f0c56d0fSKris Buschelman #undef __FUNCT__ 367f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_UMFPACK" 368f0c56d0fSKris Buschelman int MatDuplicate_UMFPACK(Mat A, MatDuplicateOption op, Mat *M) { 369f0c56d0fSKris Buschelman int ierr; 3708f340917SKris Buschelman Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr; 3718f340917SKris Buschelman 372f0c56d0fSKris Buschelman PetscFunctionBegin; 3738f340917SKris Buschelman ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 3743f953163SKris Buschelman ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_UMFPACK));CHKERRQ(ierr); 375f0c56d0fSKris Buschelman PetscFunctionReturn(0); 376f0c56d0fSKris Buschelman } 377f0c56d0fSKris Buschelman 3782f71e704SKris Buschelman /*MC 379fafad747SKris Buschelman MATUMFPACK - MATUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices 3802f71e704SKris Buschelman via the external package UMFPACK. 3812f71e704SKris Buschelman 3822f71e704SKris Buschelman If UMFPACK is installed (see the manual for 3832f71e704SKris Buschelman instructions on how to declare the existence of external packages), 3842f71e704SKris Buschelman a matrix type can be constructed which invokes UMFPACK solvers. 3852f71e704SKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,UMFPACK). 3862f71e704SKris Buschelman This matrix type is only supported for double precision real. 3872f71e704SKris Buschelman 3882f71e704SKris Buschelman This matrix inherits from MATSEQAIJ. As a result, MatSeqAIJSetPreallocation is 38928b08bd3SKris Buschelman supported for this matrix type. One can also call MatConvert for an inplace conversion to or from 39028b08bd3SKris Buschelman the MATSEQAIJ type without data copy. 3912f71e704SKris Buschelman 3922f71e704SKris Buschelman Consult UMFPACK documentation for more information about the Control parameters 3932f71e704SKris Buschelman which correspond to the options database keys below. 3942f71e704SKris Buschelman 3952f71e704SKris Buschelman Options Database Keys: 3960bad9183SKris Buschelman + -mat_type umfpack - sets the matrix type to "umfpack" during a call to MatSetFromOptions() 3972f71e704SKris Buschelman . -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL] 3982f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL] 3992f71e704SKris Buschelman . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE] 4002f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE] 4012f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT] 4022f71e704SKris Buschelman - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP] 4032f71e704SKris Buschelman 4042f71e704SKris Buschelman Level: beginner 4052f71e704SKris Buschelman 4062f71e704SKris Buschelman .seealso: PCLU 4072f71e704SKris Buschelman M*/ 4082f71e704SKris Buschelman 4091677d0b8SKris Buschelman EXTERN_C_BEGIN 4101677d0b8SKris Buschelman #undef __FUNCT__ 411f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_UMFPACK" 412f0c56d0fSKris Buschelman int MatCreate_UMFPACK(Mat A) { 4131677d0b8SKris Buschelman int ierr; 4141677d0b8SKris Buschelman 4151677d0b8SKris Buschelman PetscFunctionBegin; 4165441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqAIJ and UMFPACK types */ 4175441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATUMFPACK);CHKERRQ(ierr); 4181677d0b8SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 419d844382dSKris Buschelman ierr = MatConvert_SeqAIJ_UMFPACK(A,MATUMFPACK,&A);CHKERRQ(ierr); 4201677d0b8SKris Buschelman PetscFunctionReturn(0); 4211677d0b8SKris Buschelman } 4221677d0b8SKris Buschelman EXTERN_C_END 423