xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 2f71e70437ff1f4f201a034ed19298a7e2c6d424)
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);
24d844382dSKris Buschelman   int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
251677d0b8SKris Buschelman   int (*MatDestroy)(Mat);
261677d0b8SKris Buschelman 
271677d0b8SKris Buschelman   /* Flag to clean up UMFPACK objects during Destroy */
281677d0b8SKris Buschelman   PetscTruth CleanUpUMFPACK;
291677d0b8SKris Buschelman } Mat_SeqAIJ_UMFPACK;
30d844382dSKris Buschelman 
31d844382dSKris Buschelman EXTERN_C_BEGIN
32d844382dSKris Buschelman #undef __FUNCT__
33d844382dSKris Buschelman #define __FUNCT__ "MatConvert_UMFPACK_SeqAIJ"
34d844382dSKris Buschelman int MatConvert_UMFPACK_SeqAIJ(Mat A,MatType type,Mat *newmat) {
35d844382dSKris Buschelman   /* This routine is only called to convert an unfactored PETSc-UMFPACK matrix */
36d844382dSKris Buschelman   /* to its base PETSc type, so we will ignore 'MatType type'. */
37d844382dSKris Buschelman   int                ierr;
38d844382dSKris Buschelman   Mat                B=*newmat;
39d844382dSKris Buschelman   Mat_SeqAIJ_UMFPACK *lu=(Mat_SeqAIJ_UMFPACK*)A->spptr;
40d844382dSKris Buschelman 
41d844382dSKris Buschelman   PetscFunctionBegin;
42d844382dSKris Buschelman   if (B != A) {
43d844382dSKris Buschelman     /* This routine was inherited from SeqAIJ. */
44d844382dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
45d844382dSKris Buschelman   } else {
46d844382dSKris Buschelman     /* Reset the original function pointers */
47d844382dSKris Buschelman     A->ops->view             = lu->MatView;
48d844382dSKris Buschelman     A->ops->assemblyend      = lu->MatAssemblyEnd;
49d844382dSKris Buschelman     A->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
50d844382dSKris Buschelman     A->ops->destroy          = lu->MatDestroy;
51d844382dSKris Buschelman 
52d844382dSKris Buschelman     ierr = PetscFree(lu);CHKERRQ(ierr);
53d844382dSKris Buschelman     ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
54d844382dSKris Buschelman   }
55d844382dSKris Buschelman   *newmat = B;
56d844382dSKris Buschelman   PetscFunctionReturn(0);
57d844382dSKris Buschelman }
58d844382dSKris Buschelman EXTERN_C_END
591677d0b8SKris Buschelman 
601677d0b8SKris Buschelman #undef __FUNCT__
611677d0b8SKris Buschelman #define __FUNCT__ "MatDestroy_SeqAIJ_UMFPACK"
621677d0b8SKris Buschelman int MatDestroy_SeqAIJ_UMFPACK(Mat A)
631677d0b8SKris Buschelman {
64d844382dSKris Buschelman   int                ierr;
651677d0b8SKris Buschelman   Mat_SeqAIJ_UMFPACK *lu = (Mat_SeqAIJ_UMFPACK*)A->spptr;
661677d0b8SKris Buschelman 
671677d0b8SKris Buschelman   PetscFunctionBegin;
68fb731535SKris Buschelman   if (lu->CleanUpUMFPACK) {
691677d0b8SKris Buschelman     umfpack_di_free_symbolic(&lu->Symbolic) ;
701677d0b8SKris Buschelman     umfpack_di_free_numeric(&lu->Numeric) ;
711677d0b8SKris Buschelman     ierr = PetscFree(lu->Wi);CHKERRQ(ierr);
721677d0b8SKris Buschelman     ierr = PetscFree(lu->W);CHKERRQ(ierr);
731677d0b8SKris Buschelman 
741677d0b8SKris Buschelman     if (lu->PetscMatOdering) {
751677d0b8SKris Buschelman       ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
761677d0b8SKris Buschelman     }
771677d0b8SKris Buschelman   }
78d844382dSKris Buschelman   ierr = MatConvert_UMFPACK_SeqAIJ(A,MATSEQAIJ,&A);CHKERRQ(ierr);
79d844382dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
801677d0b8SKris Buschelman   PetscFunctionReturn(0);
811677d0b8SKris Buschelman }
821677d0b8SKris Buschelman 
831677d0b8SKris Buschelman #undef __FUNCT__
841677d0b8SKris Buschelman #define __FUNCT__ "MatSolve_SeqAIJ_UMFPACK"
851677d0b8SKris Buschelman int MatSolve_SeqAIJ_UMFPACK(Mat A,Vec b,Vec x)
861677d0b8SKris Buschelman {
871677d0b8SKris Buschelman   Mat_SeqAIJ_UMFPACK *lu = (Mat_SeqAIJ_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);
1061677d0b8SKris Buschelman   PetscFunctionReturn(0);
1071677d0b8SKris Buschelman }
1081677d0b8SKris Buschelman 
1091677d0b8SKris Buschelman #undef __FUNCT__
1101677d0b8SKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_UMFPACK"
1111677d0b8SKris Buschelman int MatLUFactorNumeric_SeqAIJ_UMFPACK(Mat A,Mat *F)
1121677d0b8SKris Buschelman {
1131677d0b8SKris Buschelman   Mat_SeqAIJ_UMFPACK *lu = (Mat_SeqAIJ_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   /* ----------------------------*/
1201677d0b8SKris Buschelman   status = umfpack_di_numeric (ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info) ;
1211677d0b8SKris Buschelman   if (status < 0) SETERRQ(1,"umfpack_di_numeric failed");
1221677d0b8SKris Buschelman   /* report numeric factorization of A' when Control[PRL] > 3 */
1231677d0b8SKris Buschelman   (void) umfpack_di_report_numeric (lu->Numeric, lu->Control) ;
1241677d0b8SKris Buschelman 
1251677d0b8SKris Buschelman   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
1261677d0b8SKris Buschelman     /* allocate working space to be used by Solve */
1271677d0b8SKris Buschelman     ierr = PetscMalloc(m * sizeof(int), &lu->Wi);CHKERRQ(ierr);
1281677d0b8SKris Buschelman     ierr = PetscMalloc(5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
1291677d0b8SKris Buschelman 
1301677d0b8SKris Buschelman     lu->flg = SAME_NONZERO_PATTERN;
1311677d0b8SKris Buschelman   }
1321677d0b8SKris Buschelman 
1331677d0b8SKris Buschelman   PetscFunctionReturn(0);
1341677d0b8SKris Buschelman }
1351677d0b8SKris Buschelman 
1361677d0b8SKris Buschelman /*
1371677d0b8SKris Buschelman    Note the r permutation is ignored
1381677d0b8SKris Buschelman */
1391677d0b8SKris Buschelman #undef __FUNCT__
1401677d0b8SKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ_UMFPACK"
1411677d0b8SKris Buschelman int MatLUFactorSymbolic_SeqAIJ_UMFPACK(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
1421677d0b8SKris Buschelman {
143e1b4c3a1SKris Buschelman   Mat                 B;
1441677d0b8SKris Buschelman   Mat_SeqAIJ          *mat=(Mat_SeqAIJ*)A->data;
1451677d0b8SKris Buschelman   Mat_SeqAIJ_UMFPACK  *lu;
1461677d0b8SKris Buschelman   int                 ierr,m=A->m,n=A->n,*ai=mat->i,*aj=mat->j,status,*ca;
1471677d0b8SKris Buschelman   PetscScalar         *av=mat->a;
1481677d0b8SKris Buschelman 
1491677d0b8SKris Buschelman   PetscFunctionBegin;
1501677d0b8SKris Buschelman   /* Create the factorization matrix F */
151e1b4c3a1SKris Buschelman   ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,m,n,&B);CHKERRQ(ierr);
152e1b4c3a1SKris Buschelman   ierr = MatSetType(B,MATUMFPACK);CHKERRQ(ierr);
153e1b4c3a1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
1541677d0b8SKris Buschelman 
155e1b4c3a1SKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_UMFPACK;
156e1b4c3a1SKris Buschelman   B->ops->solve           = MatSolve_SeqAIJ_UMFPACK;
157e1b4c3a1SKris Buschelman   B->factor               = FACTOR_LU;
158e1b4c3a1SKris Buschelman   B->assembled            = PETSC_TRUE;  /* required by -sles_view */
1591677d0b8SKris Buschelman 
160e1b4c3a1SKris Buschelman   lu = (Mat_SeqAIJ_UMFPACK*)(B->spptr);
1611677d0b8SKris Buschelman 
1621677d0b8SKris Buschelman   /* initializations */
1631677d0b8SKris Buschelman   /* ------------------------------------------------*/
1641677d0b8SKris Buschelman   /* get the default control parameters */
1651677d0b8SKris Buschelman   umfpack_di_defaults (lu->Control) ;
1661677d0b8SKris Buschelman   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
1671677d0b8SKris Buschelman 
1681677d0b8SKris Buschelman   ierr = PetscOptionsBegin(A->comm,A->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
1691677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
1701677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr);
1711677d0b8SKris Buschelman 
1721677d0b8SKris Buschelman   /* Control parameters for symbolic factorization */
1731677d0b8SKris 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);
1741677d0b8SKris 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);
1751677d0b8SKris 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);
1761677d0b8SKris Buschelman 
1771677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
1781677d0b8SKris 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);
1790801d200SKris Buschelman #if !defined(PETSC_HAVE_UMFPACK_41_OR_NEWER)
1801677d0b8SKris 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);
1811677d0b8SKris 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);
1821677d0b8SKris 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);
183fb731535SKris Buschelman #endif
1841677d0b8SKris 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);
1851677d0b8SKris Buschelman 
1861677d0b8SKris Buschelman   /* Control parameters used by solve */
1871677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr);
1881677d0b8SKris Buschelman 
1891677d0b8SKris Buschelman   /* use Petsc mat ordering (notice size is for the transpose) */
1901677d0b8SKris Buschelman   ierr = PetscOptionsHasName(PETSC_NULL,"-pc_lu_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr);
1911677d0b8SKris Buschelman   if (lu->PetscMatOdering) {
1921677d0b8SKris Buschelman     ierr = ISGetIndices(c,&ca);CHKERRQ(ierr);
1931677d0b8SKris Buschelman     ierr = PetscMalloc(A->m*sizeof(int),&lu->perm_c);CHKERRQ(ierr);
1941677d0b8SKris Buschelman     ierr = PetscMemcpy(lu->perm_c,ca,A->m*sizeof(int));CHKERRQ(ierr);
1951677d0b8SKris Buschelman     ierr = ISRestoreIndices(c,&ca);CHKERRQ(ierr);
1961677d0b8SKris Buschelman   }
1971677d0b8SKris Buschelman   PetscOptionsEnd();
1981677d0b8SKris Buschelman 
1991677d0b8SKris Buschelman   /* print the control parameters */
2001677d0b8SKris Buschelman   if( lu->Control[UMFPACK_PRL] > 1 ) umfpack_di_report_control (lu->Control);
2011677d0b8SKris Buschelman 
2021677d0b8SKris Buschelman   /* symbolic factorization of A' */
2031677d0b8SKris Buschelman   /* ---------------------------------------------------------------------- */
2040801d200SKris Buschelman #if defined(PETSC_HAVE_UMFPACK_41_OR_NEWER)
2050801d200SKris Buschelman   status = umfpack_di_qsymbolic(n,m,ai,aj,PETSC_NULL,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info) ;
2060801d200SKris Buschelman #else
2071677d0b8SKris Buschelman   status = umfpack_di_qsymbolic(n,m,ai,aj,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info) ;
2080801d200SKris Buschelman #endif
2091677d0b8SKris Buschelman   if (status < 0){
2101677d0b8SKris Buschelman     umfpack_di_report_info(lu->Control, lu->Info) ;
2111677d0b8SKris Buschelman     umfpack_di_report_status(lu->Control, status) ;
2121677d0b8SKris Buschelman     SETERRQ(1,"umfpack_di_symbolic failed");
2131677d0b8SKris Buschelman   }
2141677d0b8SKris Buschelman   /* report sumbolic factorization of A' when Control[PRL] > 3 */
2151677d0b8SKris Buschelman   (void) umfpack_di_report_symbolic(lu->Symbolic, lu->Control) ;
2161677d0b8SKris Buschelman 
2171677d0b8SKris Buschelman   lu->flg = DIFFERENT_NONZERO_PATTERN;
2181677d0b8SKris Buschelman   lu->ai  = ai;
2191677d0b8SKris Buschelman   lu->aj  = aj;
2201677d0b8SKris Buschelman   lu->av  = av;
221e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK = PETSC_TRUE;
222e1b4c3a1SKris Buschelman   *F = B;
2231677d0b8SKris Buschelman   PetscFunctionReturn(0);
2241677d0b8SKris Buschelman }
2251677d0b8SKris Buschelman 
2261677d0b8SKris Buschelman #undef __FUNCT__
227d844382dSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SeqAIJ_UMFPACK"
228d844382dSKris Buschelman int MatAssemblyEnd_SeqAIJ_UMFPACK(Mat A,MatAssemblyType mode) {
229d844382dSKris Buschelman   int                ierr;
230d844382dSKris Buschelman   Mat_SeqAIJ_UMFPACK *lu=(Mat_SeqAIJ_UMFPACK*)(A->spptr);
231d844382dSKris Buschelman 
2321677d0b8SKris Buschelman   PetscFunctionBegin;
233d844382dSKris Buschelman   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
234d844382dSKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
2351677d0b8SKris Buschelman   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ_UMFPACK;
2361677d0b8SKris Buschelman   PetscFunctionReturn(0);
2371677d0b8SKris Buschelman }
2381677d0b8SKris Buschelman 
2391677d0b8SKris Buschelman /* used by -sles_view */
2401677d0b8SKris Buschelman #undef __FUNCT__
2411677d0b8SKris Buschelman #define __FUNCT__ "MatSeqAIJFactorInfo_UMFPACK"
2421677d0b8SKris Buschelman int MatSeqAIJFactorInfo_UMFPACK(Mat A,PetscViewer viewer)
2431677d0b8SKris Buschelman {
2441677d0b8SKris Buschelman   Mat_SeqAIJ_UMFPACK      *lu= (Mat_SeqAIJ_UMFPACK*)A->spptr;
2451677d0b8SKris Buschelman   int                     ierr;
2461677d0b8SKris Buschelman   PetscFunctionBegin;
2471677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
2481677d0b8SKris Buschelman   if (A->ops->solve != MatSolve_SeqAIJ_UMFPACK) PetscFunctionReturn(0);
2491677d0b8SKris Buschelman 
2501677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
2511677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
2521677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
2531677d0b8SKris Buschelman 
2541677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
2551677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
2561677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
2571677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
2581677d0b8SKris Buschelman 
2591677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
2601677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
26154c494e3SKris Buschelman #if !defined(PETSC_HAVE_UMFPACK_41_OR_NEWER)
2621677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_RELAXED_AMALGAMATION]: %g\n",lu->Control[UMFPACK_RELAXED_AMALGAMATION]);CHKERRQ(ierr);
2631677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_RELAXED2_AMALGAMATION]: %g\n",lu->Control[UMFPACK_RELAXED2_AMALGAMATION]);CHKERRQ(ierr);
2641677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_RELAXED3_AMALGAMATION]: %g\n",lu->Control[UMFPACK_RELAXED3_AMALGAMATION]);CHKERRQ(ierr);
265fb731535SKris Buschelman #endif
2661677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
2671677d0b8SKris Buschelman 
2681677d0b8SKris Buschelman   /* Control parameters used by solve */
2691677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
2701677d0b8SKris Buschelman 
2711677d0b8SKris Buschelman   /* mat ordering */
2721677d0b8SKris Buschelman   if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer,"  UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr);
2731677d0b8SKris Buschelman 
2741677d0b8SKris Buschelman   PetscFunctionReturn(0);
2751677d0b8SKris Buschelman }
2761677d0b8SKris Buschelman 
277d844382dSKris Buschelman #undef __FUNCT__
278d844382dSKris Buschelman #define __FUNCT__ "MatView_SeqAIJ_UMFPACK"
279d844382dSKris Buschelman int MatView_SeqAIJ_UMFPACK(Mat A,PetscViewer viewer)
280d844382dSKris Buschelman {
281d844382dSKris Buschelman   int                 ierr;
282d844382dSKris Buschelman   PetscTruth          isascii;
283d844382dSKris Buschelman   PetscViewerFormat   format;
284d844382dSKris Buschelman   Mat_SeqAIJ_UMFPACK  *lu=(Mat_SeqAIJ_UMFPACK*)(A->spptr);
285d844382dSKris Buschelman 
286d844382dSKris Buschelman   PetscFunctionBegin;
287d844382dSKris Buschelman   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
288d844382dSKris Buschelman 
289d844382dSKris Buschelman   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
290d844382dSKris Buschelman   if (isascii) {
291d844382dSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
292d844382dSKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
293d844382dSKris Buschelman       ierr = MatSeqAIJFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
294d844382dSKris Buschelman     }
295d844382dSKris Buschelman   }
296d844382dSKris Buschelman   PetscFunctionReturn(0);
297d844382dSKris Buschelman }
298d844382dSKris Buschelman 
299d844382dSKris Buschelman EXTERN_C_BEGIN
300d844382dSKris Buschelman #undef __FUNCT__
301d844382dSKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_UMFPACK"
302d844382dSKris Buschelman int MatConvert_SeqAIJ_UMFPACK(Mat A,MatType type,Mat *newmat) {
303d844382dSKris Buschelman   /* This routine is only called to convert to MATUMFPACK */
304d844382dSKris Buschelman   /* from MATSEQAIJ, so we will ignore 'MatType type'. */
305d844382dSKris Buschelman   int                ierr;
306d844382dSKris Buschelman   Mat                B=*newmat;
307d844382dSKris Buschelman   Mat_SeqAIJ_UMFPACK *lu;
308d844382dSKris Buschelman 
309d844382dSKris Buschelman   PetscFunctionBegin;
310d844382dSKris Buschelman   if (B != A) {
311d844382dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
312d844382dSKris Buschelman   }
313d844382dSKris Buschelman 
314d844382dSKris Buschelman   ierr = PetscNew(Mat_SeqAIJ_UMFPACK,&lu);CHKERRQ(ierr);
315d844382dSKris Buschelman   lu->MatView              = A->ops->view;
316d844382dSKris Buschelman   lu->MatAssemblyEnd       = A->ops->assemblyend;
317d844382dSKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
318d844382dSKris Buschelman   lu->MatDestroy           = A->ops->destroy;
319d844382dSKris Buschelman   lu->CleanUpUMFPACK       = PETSC_FALSE;
320d844382dSKris Buschelman 
321d844382dSKris Buschelman   B->spptr                 = (void*)lu;
322d844382dSKris Buschelman   B->ops->view             = MatView_SeqAIJ_UMFPACK;
323d844382dSKris Buschelman   B->ops->assemblyend      = MatAssemblyEnd_SeqAIJ_UMFPACK;
324d844382dSKris Buschelman   B->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ_UMFPACK;
325d844382dSKris Buschelman   B->ops->destroy          = MatDestroy_SeqAIJ_UMFPACK;
326d844382dSKris Buschelman 
327d844382dSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_umfpack_C",
328d844382dSKris Buschelman                                            "MatConvert_SeqAIJ_UMFPACK",MatConvert_SeqAIJ_UMFPACK);CHKERRQ(ierr);
329d844382dSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_umfpack_seqaij_C",
330d844382dSKris Buschelman                                            "MatConvert_UMFPACK_SeqAIJ",MatConvert_UMFPACK_SeqAIJ);CHKERRQ(ierr);
331d844382dSKris Buschelman   PetscLogInfo(0,"Using UMFPACK for SeqAIJ LU factorization and solves.");
332d844382dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATUMFPACK);CHKERRQ(ierr);
333d844382dSKris Buschelman   *newmat = B;
334d844382dSKris Buschelman   PetscFunctionReturn(0);
335d844382dSKris Buschelman }
336d844382dSKris Buschelman EXTERN_C_END
337d844382dSKris Buschelman 
338*2f71e704SKris Buschelman /*MC
339*2f71e704SKris Buschelman   MATUMFPACK - a matrix type providing direct solvers (LU) for sequential matrices
340*2f71e704SKris Buschelman   via the external package UMFPACK.
341*2f71e704SKris Buschelman 
342*2f71e704SKris Buschelman   If UMFPACK is installed (see the manual for
343*2f71e704SKris Buschelman   instructions on how to declare the existence of external packages),
344*2f71e704SKris Buschelman   a matrix type can be constructed which invokes UMFPACK solvers.
345*2f71e704SKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,UMFPACK).
346*2f71e704SKris Buschelman   This matrix type is only supported for double precision real.
347*2f71e704SKris Buschelman 
348*2f71e704SKris Buschelman   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
349*2f71e704SKris Buschelman   supported for this matrix type.
350*2f71e704SKris Buschelman 
351*2f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
352*2f71e704SKris Buschelman   which correspond to the options database keys below.
353*2f71e704SKris Buschelman 
354*2f71e704SKris Buschelman   Options Database Keys:
355*2f71e704SKris Buschelman + -mat_type umfpack - sets the matrix type to umfpack during a call to MatSetFromOptions()
356*2f71e704SKris Buschelman . -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL]
357*2f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
358*2f71e704SKris Buschelman . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
359*2f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
360*2f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
361*2f71e704SKris Buschelman - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
362*2f71e704SKris Buschelman 
363*2f71e704SKris Buschelman    Level: beginner
364*2f71e704SKris Buschelman 
365*2f71e704SKris Buschelman .seealso: PCLU
366*2f71e704SKris Buschelman M*/
367*2f71e704SKris Buschelman 
3681677d0b8SKris Buschelman EXTERN_C_BEGIN
3691677d0b8SKris Buschelman #undef __FUNCT__
3701677d0b8SKris Buschelman #define __FUNCT__ "MatCreate_SeqAIJ_UMFPACK"
3711677d0b8SKris Buschelman int MatCreate_SeqAIJ_UMFPACK(Mat A) {
3721677d0b8SKris Buschelman   int                ierr;
3731677d0b8SKris Buschelman 
3741677d0b8SKris Buschelman   PetscFunctionBegin;
3751677d0b8SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
376d844382dSKris Buschelman   ierr = MatConvert_SeqAIJ_UMFPACK(A,MATUMFPACK,&A);CHKERRQ(ierr);
3771677d0b8SKris Buschelman   PetscFunctionReturn(0);
3781677d0b8SKris Buschelman }
3791677d0b8SKris Buschelman EXTERN_C_END
380