xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 2a325a84c8111bef5aae3ce3ac98ba23198eb906)
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 */
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;
90*2a325a84SHong Zhang   int         m;
911677d0b8SKris Buschelman 
921677d0b8SKris Buschelman   PetscFunctionBegin;
931677d0b8SKris Buschelman   /* solve Ax = b by umfpack_di_wsolve */
941677d0b8SKris Buschelman   /* ----------------------------------*/
951677d0b8SKris Buschelman   ierr = VecGetArray(b,&ba);
961677d0b8SKris Buschelman   ierr = VecGetArray(x,&xa);
971677d0b8SKris Buschelman 
981677d0b8SKris Buschelman   status = umfpack_di_wsolve(UMFPACK_At,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
991677d0b8SKris Buschelman   umfpack_di_report_info(lu->Control, lu->Info);
1001677d0b8SKris Buschelman   if (status < 0){
1011677d0b8SKris Buschelman     umfpack_di_report_status(lu->Control, status) ;
1021677d0b8SKris Buschelman     SETERRQ(1,"umfpack_di_wsolve failed") ;
1031677d0b8SKris Buschelman   }
1041677d0b8SKris Buschelman 
1051677d0b8SKris Buschelman   ierr = VecRestoreArray(b,&ba);
1061677d0b8SKris Buschelman   ierr = VecRestoreArray(x,&xa);
107*2a325a84SHong Zhang 
1081677d0b8SKris Buschelman   PetscFunctionReturn(0);
1091677d0b8SKris Buschelman }
1101677d0b8SKris Buschelman 
1111677d0b8SKris Buschelman #undef __FUNCT__
112f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_UMFPACK"
113f0c56d0fSKris Buschelman int MatLUFactorNumeric_UMFPACK(Mat A,Mat *F) {
114f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)(*F)->spptr;
1151677d0b8SKris Buschelman   int         *ai=lu->ai,*aj=lu->aj,m=A->m,status,ierr;
1161677d0b8SKris Buschelman   PetscScalar *av=lu->av;
1171677d0b8SKris Buschelman 
1181677d0b8SKris Buschelman   PetscFunctionBegin;
1191677d0b8SKris Buschelman   /* numeric factorization of A' */
1201677d0b8SKris Buschelman   /* ----------------------------*/
121*2a325a84SHong Zhang   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
122*2a325a84SHong Zhang       umfpack_di_free_numeric(&lu->Numeric) ;
123*2a325a84SHong Zhang   }
1241677d0b8SKris Buschelman   status = umfpack_di_numeric (ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info) ;
1251677d0b8SKris Buschelman   if (status < 0) SETERRQ(1,"umfpack_di_numeric failed");
1261677d0b8SKris Buschelman   /* report numeric factorization of A' when Control[PRL] > 3 */
1271677d0b8SKris Buschelman   (void) umfpack_di_report_numeric (lu->Numeric, lu->Control) ;
1281677d0b8SKris Buschelman 
1291677d0b8SKris Buschelman   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
1301677d0b8SKris Buschelman     /* allocate working space to be used by Solve */
1311677d0b8SKris Buschelman     ierr = PetscMalloc(m * sizeof(int), &lu->Wi);CHKERRQ(ierr);
1321677d0b8SKris Buschelman     ierr = PetscMalloc(5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
1331677d0b8SKris Buschelman   }
134*2a325a84SHong Zhang   lu->flg = SAME_NONZERO_PATTERN;
135*2a325a84SHong Zhang   lu->CleanUpUMFPACK = PETSC_TRUE;
1361677d0b8SKris Buschelman   PetscFunctionReturn(0);
1371677d0b8SKris Buschelman }
1381677d0b8SKris Buschelman 
1391677d0b8SKris Buschelman /*
1401677d0b8SKris Buschelman    Note the r permutation is ignored
1411677d0b8SKris Buschelman */
1421677d0b8SKris Buschelman #undef __FUNCT__
143f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK"
144f0c56d0fSKris Buschelman int MatLUFactorSymbolic_UMFPACK(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
145e1b4c3a1SKris Buschelman   Mat         B;
1461677d0b8SKris Buschelman   Mat_SeqAIJ  *mat=(Mat_SeqAIJ*)A->data;
147f0c56d0fSKris Buschelman   Mat_UMFPACK *lu;
1481677d0b8SKris Buschelman   int         ierr,m=A->m,n=A->n,*ai=mat->i,*aj=mat->j,status,*ca;
1491677d0b8SKris Buschelman   PetscScalar *av=mat->a;
1501677d0b8SKris Buschelman 
1511677d0b8SKris Buschelman   PetscFunctionBegin;
1521677d0b8SKris Buschelman   /* Create the factorization matrix F */
153e1b4c3a1SKris Buschelman   ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,m,n,&B);CHKERRQ(ierr);
154be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
155e1b4c3a1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
1561677d0b8SKris Buschelman 
157f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
158f0c56d0fSKris Buschelman   B->ops->solve           = MatSolve_UMFPACK;
159e1b4c3a1SKris Buschelman   B->factor               = FACTOR_LU;
16094b7f48cSBarry Smith   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */
1611677d0b8SKris Buschelman 
162f0c56d0fSKris Buschelman   lu = (Mat_UMFPACK*)(B->spptr);
1631677d0b8SKris Buschelman 
1641677d0b8SKris Buschelman   /* initializations */
1651677d0b8SKris Buschelman   /* ------------------------------------------------*/
1661677d0b8SKris Buschelman   /* get the default control parameters */
1671677d0b8SKris Buschelman   umfpack_di_defaults (lu->Control) ;
1681677d0b8SKris Buschelman   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
1691677d0b8SKris Buschelman 
1701677d0b8SKris Buschelman   ierr = PetscOptionsBegin(A->comm,A->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
1711677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
1721677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr);
1731677d0b8SKris Buschelman 
1741677d0b8SKris Buschelman   /* Control parameters for symbolic factorization */
1751677d0b8SKris 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);
1761677d0b8SKris 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);
1771677d0b8SKris 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);
1781677d0b8SKris Buschelman 
1791677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
1801677d0b8SKris 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);
1810801d200SKris Buschelman #if !defined(PETSC_HAVE_UMFPACK_41_OR_NEWER)
1821677d0b8SKris 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);
1831677d0b8SKris 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);
1841677d0b8SKris 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);
185fb731535SKris Buschelman #endif
1861677d0b8SKris 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);
1871677d0b8SKris Buschelman 
1881677d0b8SKris Buschelman   /* Control parameters used by solve */
1891677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr);
1901677d0b8SKris Buschelman 
1911677d0b8SKris Buschelman   /* use Petsc mat ordering (notice size is for the transpose) */
1921677d0b8SKris Buschelman   ierr = PetscOptionsHasName(PETSC_NULL,"-pc_lu_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr);
1931677d0b8SKris Buschelman   if (lu->PetscMatOdering) {
1941677d0b8SKris Buschelman     ierr = ISGetIndices(c,&ca);CHKERRQ(ierr);
1951677d0b8SKris Buschelman     ierr = PetscMalloc(A->m*sizeof(int),&lu->perm_c);CHKERRQ(ierr);
1961677d0b8SKris Buschelman     ierr = PetscMemcpy(lu->perm_c,ca,A->m*sizeof(int));CHKERRQ(ierr);
1971677d0b8SKris Buschelman     ierr = ISRestoreIndices(c,&ca);CHKERRQ(ierr);
1981677d0b8SKris Buschelman   }
1991677d0b8SKris Buschelman   PetscOptionsEnd();
2001677d0b8SKris Buschelman 
2011677d0b8SKris Buschelman   /* print the control parameters */
2021677d0b8SKris Buschelman   if( lu->Control[UMFPACK_PRL] > 1 ) umfpack_di_report_control (lu->Control);
2031677d0b8SKris Buschelman 
2041677d0b8SKris Buschelman   /* symbolic factorization of A' */
2051677d0b8SKris Buschelman   /* ---------------------------------------------------------------------- */
2060801d200SKris Buschelman #if defined(PETSC_HAVE_UMFPACK_41_OR_NEWER)
2070801d200SKris Buschelman   status = umfpack_di_qsymbolic(n,m,ai,aj,PETSC_NULL,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info) ;
2080801d200SKris Buschelman #else
2091677d0b8SKris Buschelman   status = umfpack_di_qsymbolic(n,m,ai,aj,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info) ;
2100801d200SKris Buschelman #endif
2111677d0b8SKris Buschelman   if (status < 0){
2121677d0b8SKris Buschelman     umfpack_di_report_info(lu->Control, lu->Info) ;
2131677d0b8SKris Buschelman     umfpack_di_report_status(lu->Control, status) ;
2141677d0b8SKris Buschelman     SETERRQ(1,"umfpack_di_symbolic failed");
2151677d0b8SKris Buschelman   }
2161677d0b8SKris Buschelman   /* report sumbolic factorization of A' when Control[PRL] > 3 */
2171677d0b8SKris Buschelman   (void) umfpack_di_report_symbolic(lu->Symbolic, lu->Control) ;
2181677d0b8SKris Buschelman 
2191677d0b8SKris Buschelman   lu->flg = DIFFERENT_NONZERO_PATTERN;
2201677d0b8SKris Buschelman   lu->ai  = ai;
2211677d0b8SKris Buschelman   lu->aj  = aj;
2221677d0b8SKris Buschelman   lu->av  = av;
223e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK = PETSC_TRUE;
224e1b4c3a1SKris Buschelman   *F = B;
2251677d0b8SKris Buschelman   PetscFunctionReturn(0);
2261677d0b8SKris Buschelman }
2271677d0b8SKris Buschelman 
2281677d0b8SKris Buschelman #undef __FUNCT__
229f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_UMFPACK"
230f0c56d0fSKris Buschelman int MatAssemblyEnd_UMFPACK(Mat A,MatAssemblyType mode) {
231d844382dSKris Buschelman   int         ierr;
232f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)(A->spptr);
233d844382dSKris Buschelman 
2341677d0b8SKris Buschelman   PetscFunctionBegin;
235d844382dSKris Buschelman   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
236d844382dSKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
237f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
2381677d0b8SKris Buschelman   PetscFunctionReturn(0);
2391677d0b8SKris Buschelman }
2401677d0b8SKris Buschelman 
24194b7f48cSBarry Smith /* used by -ksp_view */
2421677d0b8SKris Buschelman #undef __FUNCT__
243f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_UMFPACK"
244f0c56d0fSKris Buschelman int MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer) {
245f0c56d0fSKris Buschelman   Mat_UMFPACK *lu= (Mat_UMFPACK*)A->spptr;
2461677d0b8SKris Buschelman   int         ierr;
247f0c56d0fSKris Buschelman 
2481677d0b8SKris Buschelman   PetscFunctionBegin;
2491677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
250f0c56d0fSKris Buschelman   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
2511677d0b8SKris Buschelman 
2521677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
2531677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
2541677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
2551677d0b8SKris Buschelman 
2561677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
2571677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
2581677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
2591677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
2601677d0b8SKris Buschelman 
2611677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
2621677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
26354c494e3SKris Buschelman #if !defined(PETSC_HAVE_UMFPACK_41_OR_NEWER)
2641677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_RELAXED_AMALGAMATION]: %g\n",lu->Control[UMFPACK_RELAXED_AMALGAMATION]);CHKERRQ(ierr);
2651677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_RELAXED2_AMALGAMATION]: %g\n",lu->Control[UMFPACK_RELAXED2_AMALGAMATION]);CHKERRQ(ierr);
2661677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_RELAXED3_AMALGAMATION]: %g\n",lu->Control[UMFPACK_RELAXED3_AMALGAMATION]);CHKERRQ(ierr);
267fb731535SKris Buschelman #endif
2681677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
2691677d0b8SKris Buschelman 
2701677d0b8SKris Buschelman   /* Control parameters used by solve */
2711677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
2721677d0b8SKris Buschelman 
2731677d0b8SKris Buschelman   /* mat ordering */
2741677d0b8SKris Buschelman   if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer,"  UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr);
2751677d0b8SKris Buschelman 
2761677d0b8SKris Buschelman   PetscFunctionReturn(0);
2771677d0b8SKris Buschelman }
2781677d0b8SKris Buschelman 
279d844382dSKris Buschelman #undef __FUNCT__
280f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_UMFPACK"
281f0c56d0fSKris Buschelman int MatView_UMFPACK(Mat A,PetscViewer viewer) {
282d844382dSKris Buschelman   int               ierr;
283d844382dSKris Buschelman   PetscTruth        isascii;
284d844382dSKris Buschelman   PetscViewerFormat format;
285f0c56d0fSKris Buschelman   Mat_UMFPACK       *lu=(Mat_UMFPACK*)(A->spptr);
286d844382dSKris Buschelman 
287d844382dSKris Buschelman   PetscFunctionBegin;
288d844382dSKris Buschelman   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
289d844382dSKris Buschelman 
290d844382dSKris Buschelman   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
291d844382dSKris Buschelman   if (isascii) {
292d844382dSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
293d844382dSKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
294f0c56d0fSKris Buschelman       ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
295d844382dSKris Buschelman     }
296d844382dSKris Buschelman   }
297d844382dSKris Buschelman   PetscFunctionReturn(0);
298d844382dSKris Buschelman }
299d844382dSKris Buschelman 
300d844382dSKris Buschelman EXTERN_C_BEGIN
301d844382dSKris Buschelman #undef __FUNCT__
302d844382dSKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_UMFPACK"
3038e9aea5cSBarry Smith int MatConvert_SeqAIJ_UMFPACK(Mat A,const MatType type,Mat *newmat) {
304d844382dSKris Buschelman   /* This routine is only called to convert to MATUMFPACK */
305d844382dSKris Buschelman   /* from MATSEQAIJ, so we will ignore 'MatType type'. */
306d844382dSKris Buschelman   int         ierr;
307d844382dSKris Buschelman   Mat         B=*newmat;
308f0c56d0fSKris Buschelman   Mat_UMFPACK *lu;
309d844382dSKris Buschelman 
310d844382dSKris Buschelman   PetscFunctionBegin;
311d844382dSKris Buschelman   if (B != A) {
312d844382dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
313d844382dSKris Buschelman   }
314d844382dSKris Buschelman 
315f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_UMFPACK,&lu);CHKERRQ(ierr);
316f0c56d0fSKris Buschelman   lu->MatDuplicate         = A->ops->duplicate;
317d844382dSKris Buschelman   lu->MatView              = A->ops->view;
318d844382dSKris Buschelman   lu->MatAssemblyEnd       = A->ops->assemblyend;
319d844382dSKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
320d844382dSKris Buschelman   lu->MatDestroy           = A->ops->destroy;
321d844382dSKris Buschelman   lu->CleanUpUMFPACK       = PETSC_FALSE;
322d844382dSKris Buschelman 
323d844382dSKris Buschelman   B->spptr                 = (void*)lu;
324f0c56d0fSKris Buschelman   B->ops->duplicate        = MatDuplicate_UMFPACK;
325f0c56d0fSKris Buschelman   B->ops->view             = MatView_UMFPACK;
326f0c56d0fSKris Buschelman   B->ops->assemblyend      = MatAssemblyEnd_UMFPACK;
327f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
328f0c56d0fSKris Buschelman   B->ops->destroy          = MatDestroy_UMFPACK;
329d844382dSKris Buschelman 
330d844382dSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_umfpack_C",
331d844382dSKris Buschelman                                            "MatConvert_SeqAIJ_UMFPACK",MatConvert_SeqAIJ_UMFPACK);CHKERRQ(ierr);
332d844382dSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_umfpack_seqaij_C",
333d844382dSKris Buschelman                                            "MatConvert_UMFPACK_SeqAIJ",MatConvert_UMFPACK_SeqAIJ);CHKERRQ(ierr);
334d844382dSKris Buschelman   PetscLogInfo(0,"Using UMFPACK for SeqAIJ LU factorization and solves.");
335d844382dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATUMFPACK);CHKERRQ(ierr);
336d844382dSKris Buschelman   *newmat = B;
337d844382dSKris Buschelman   PetscFunctionReturn(0);
338d844382dSKris Buschelman }
339d844382dSKris Buschelman EXTERN_C_END
340d844382dSKris Buschelman 
341f0c56d0fSKris Buschelman #undef __FUNCT__
342f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_UMFPACK"
343f0c56d0fSKris Buschelman int MatDuplicate_UMFPACK(Mat A, MatDuplicateOption op, Mat *M) {
344f0c56d0fSKris Buschelman   int         ierr;
3458f340917SKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
3468f340917SKris Buschelman 
347f0c56d0fSKris Buschelman   PetscFunctionBegin;
3488f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
3493f953163SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_UMFPACK));CHKERRQ(ierr);
350f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
351f0c56d0fSKris Buschelman }
352f0c56d0fSKris Buschelman 
3532f71e704SKris Buschelman /*MC
354fafad747SKris Buschelman   MATUMFPACK - MATUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
3552f71e704SKris Buschelman   via the external package UMFPACK.
3562f71e704SKris Buschelman 
3572f71e704SKris Buschelman   If UMFPACK is installed (see the manual for
3582f71e704SKris Buschelman   instructions on how to declare the existence of external packages),
3592f71e704SKris Buschelman   a matrix type can be constructed which invokes UMFPACK solvers.
3602f71e704SKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,UMFPACK).
3612f71e704SKris Buschelman   This matrix type is only supported for double precision real.
3622f71e704SKris Buschelman 
3632f71e704SKris Buschelman   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
36428b08bd3SKris Buschelman   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from
36528b08bd3SKris Buschelman   the MATSEQAIJ type without data copy.
3662f71e704SKris Buschelman 
3672f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
3682f71e704SKris Buschelman   which correspond to the options database keys below.
3692f71e704SKris Buschelman 
3702f71e704SKris Buschelman   Options Database Keys:
3710bad9183SKris Buschelman + -mat_type umfpack - sets the matrix type to "umfpack" during a call to MatSetFromOptions()
3722f71e704SKris Buschelman . -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL]
3732f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
3742f71e704SKris Buschelman . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
3752f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
3762f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
3772f71e704SKris Buschelman - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
3782f71e704SKris Buschelman 
3792f71e704SKris Buschelman    Level: beginner
3802f71e704SKris Buschelman 
3812f71e704SKris Buschelman .seealso: PCLU
3822f71e704SKris Buschelman M*/
3832f71e704SKris Buschelman 
3841677d0b8SKris Buschelman EXTERN_C_BEGIN
3851677d0b8SKris Buschelman #undef __FUNCT__
386f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_UMFPACK"
387f0c56d0fSKris Buschelman int MatCreate_UMFPACK(Mat A) {
3881677d0b8SKris Buschelman   int                ierr;
3891677d0b8SKris Buschelman 
3901677d0b8SKris Buschelman   PetscFunctionBegin;
3915441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ and UMFPACK types */
3925441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATUMFPACK);CHKERRQ(ierr);
3931677d0b8SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
394d844382dSKris Buschelman   ierr = MatConvert_SeqAIJ_UMFPACK(A,MATUMFPACK,&A);CHKERRQ(ierr);
3951677d0b8SKris Buschelman   PetscFunctionReturn(0);
3961677d0b8SKris Buschelman }
3971677d0b8SKris Buschelman EXTERN_C_END
398