xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 32077d6df18512e94702c86e7d562481ed0cebd0)
11677d0b8SKris Buschelman 
21677d0b8SKris Buschelman /*
359698cb0SHong Zhang         Provides an interface to the UMFPACKv4.3 sparse solver
41677d0b8SKris Buschelman */
51677d0b8SKris Buschelman 
61677d0b8SKris Buschelman #include "src/mat/impls/aij/seq/aij.h"
71677d0b8SKris Buschelman 
81677d0b8SKris Buschelman EXTERN_C_BEGIN
91677d0b8SKris Buschelman #include "umfpack.h"
101677d0b8SKris Buschelman EXTERN_C_END
111677d0b8SKris Buschelman 
121677d0b8SKris Buschelman typedef struct {
131677d0b8SKris Buschelman   void         *Symbolic, *Numeric;
141677d0b8SKris Buschelman   double       Info[UMFPACK_INFO], Control[UMFPACK_CONTROL],*W;
151677d0b8SKris Buschelman   int          *Wi,*ai,*aj,*perm_c;
161677d0b8SKris Buschelman   PetscScalar  *av;
171677d0b8SKris Buschelman   MatStructure flg;
181677d0b8SKris Buschelman   PetscTruth   PetscMatOdering;
191677d0b8SKris Buschelman 
201677d0b8SKris Buschelman   /* A few function pointers for inheritance */
21f0c56d0fSKris Buschelman   int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
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;
29f0c56d0fSKris Buschelman } Mat_UMFPACK;
30f0c56d0fSKris Buschelman 
31f0c56d0fSKris Buschelman EXTERN int MatDuplicate_UMFPACK(Mat,MatDuplicateOption,Mat*);
32d844382dSKris Buschelman 
33d844382dSKris Buschelman EXTERN_C_BEGIN
34d844382dSKris Buschelman #undef __FUNCT__
35d844382dSKris Buschelman #define __FUNCT__ "MatConvert_UMFPACK_SeqAIJ"
368e9aea5cSBarry Smith int MatConvert_UMFPACK_SeqAIJ(Mat A,const MatType type,Mat *newmat) {
37d844382dSKris Buschelman   /* This routine is only called to convert an unfactored PETSc-UMFPACK matrix */
38d844382dSKris Buschelman   /* to its base PETSc type, so we will ignore 'MatType type'. */
39d844382dSKris Buschelman   int         ierr;
40d844382dSKris Buschelman   Mat         B=*newmat;
41f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
42d844382dSKris Buschelman 
43d844382dSKris Buschelman   PetscFunctionBegin;
44d844382dSKris Buschelman   if (B != A) {
45d844382dSKris Buschelman     /* This routine was inherited from SeqAIJ. */
46d844382dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
47f0c56d0fSKris Buschelman   }
48d844382dSKris Buschelman   /* Reset the original function pointers */
49f0c56d0fSKris Buschelman   A->ops->duplicate        = lu->MatDuplicate;
50d844382dSKris Buschelman   A->ops->view             = lu->MatView;
51d844382dSKris Buschelman   A->ops->assemblyend      = lu->MatAssemblyEnd;
52d844382dSKris Buschelman   A->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
53d844382dSKris Buschelman   A->ops->destroy          = lu->MatDestroy;
54d844382dSKris Buschelman 
55d844382dSKris Buschelman   ierr = PetscFree(lu);CHKERRQ(ierr);
56d844382dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
57d844382dSKris Buschelman   *newmat = B;
58d844382dSKris Buschelman   PetscFunctionReturn(0);
59d844382dSKris Buschelman }
60d844382dSKris Buschelman EXTERN_C_END
611677d0b8SKris Buschelman 
621677d0b8SKris Buschelman #undef __FUNCT__
63f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_UMFPACK"
64f0c56d0fSKris Buschelman int MatDestroy_UMFPACK(Mat A) {
65d844382dSKris Buschelman   int         ierr;
66f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
671677d0b8SKris Buschelman 
681677d0b8SKris Buschelman   PetscFunctionBegin;
69fb731535SKris Buschelman   if (lu->CleanUpUMFPACK) {
701677d0b8SKris Buschelman     umfpack_di_free_symbolic(&lu->Symbolic) ;
711677d0b8SKris Buschelman     umfpack_di_free_numeric(&lu->Numeric) ;
721677d0b8SKris Buschelman     ierr = PetscFree(lu->Wi);CHKERRQ(ierr);
731677d0b8SKris Buschelman     ierr = PetscFree(lu->W);CHKERRQ(ierr);
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__
84f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_UMFPACK"
85f0c56d0fSKris Buschelman int MatSolve_UMFPACK(Mat A,Vec b,Vec x) {
86f0c56d0fSKris Buschelman   Mat_UMFPACK *lu = (Mat_UMFPACK*)A->spptr;
871677d0b8SKris Buschelman   PetscScalar *av=lu->av,*ba,*xa;
881677d0b8SKris Buschelman   int         ierr,*ai=lu->ai,*aj=lu->aj,status;
891677d0b8SKris Buschelman 
901677d0b8SKris Buschelman   PetscFunctionBegin;
911677d0b8SKris Buschelman   /* solve Ax = b by umfpack_di_wsolve */
921677d0b8SKris Buschelman   /* ----------------------------------*/
931677d0b8SKris Buschelman   ierr = VecGetArray(b,&ba);
941677d0b8SKris Buschelman   ierr = VecGetArray(x,&xa);
951677d0b8SKris Buschelman 
961677d0b8SKris Buschelman   status = umfpack_di_wsolve(UMFPACK_At,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
971677d0b8SKris Buschelman   umfpack_di_report_info(lu->Control, lu->Info);
981677d0b8SKris Buschelman   if (status < 0){
991677d0b8SKris Buschelman     umfpack_di_report_status(lu->Control, status) ;
1001677d0b8SKris Buschelman     SETERRQ(1,"umfpack_di_wsolve failed") ;
1011677d0b8SKris Buschelman   }
1021677d0b8SKris Buschelman 
1031677d0b8SKris Buschelman   ierr = VecRestoreArray(b,&ba);
1041677d0b8SKris Buschelman   ierr = VecRestoreArray(x,&xa);
1052a325a84SHong Zhang 
1061677d0b8SKris Buschelman   PetscFunctionReturn(0);
1071677d0b8SKris Buschelman }
1081677d0b8SKris Buschelman 
1091677d0b8SKris Buschelman #undef __FUNCT__
110f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_UMFPACK"
111f0c56d0fSKris Buschelman int MatLUFactorNumeric_UMFPACK(Mat A,Mat *F) {
112f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)(*F)->spptr;
1131677d0b8SKris Buschelman   int         *ai=lu->ai,*aj=lu->aj,m=A->m,status,ierr;
1141677d0b8SKris Buschelman   PetscScalar *av=lu->av;
1151677d0b8SKris Buschelman 
1161677d0b8SKris Buschelman   PetscFunctionBegin;
1171677d0b8SKris Buschelman   /* numeric factorization of A' */
1181677d0b8SKris Buschelman   /* ----------------------------*/
1192a325a84SHong Zhang   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
1202a325a84SHong Zhang       umfpack_di_free_numeric(&lu->Numeric) ;
1212a325a84SHong Zhang   }
1221677d0b8SKris Buschelman   status = umfpack_di_numeric (ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info) ;
1231677d0b8SKris Buschelman   if (status < 0) SETERRQ(1,"umfpack_di_numeric failed");
1241677d0b8SKris Buschelman   /* report numeric factorization of A' when Control[PRL] > 3 */
1251677d0b8SKris Buschelman   (void) umfpack_di_report_numeric (lu->Numeric, lu->Control) ;
1261677d0b8SKris Buschelman 
1271677d0b8SKris Buschelman   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
1281677d0b8SKris Buschelman     /* allocate working space to be used by Solve */
1291677d0b8SKris Buschelman     ierr = PetscMalloc(m * sizeof(int), &lu->Wi);CHKERRQ(ierr);
1301677d0b8SKris Buschelman     ierr = PetscMalloc(5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
1311677d0b8SKris Buschelman   }
1322a325a84SHong Zhang   lu->flg = SAME_NONZERO_PATTERN;
1332a325a84SHong Zhang   lu->CleanUpUMFPACK = PETSC_TRUE;
1341677d0b8SKris Buschelman   PetscFunctionReturn(0);
1351677d0b8SKris Buschelman }
1361677d0b8SKris Buschelman 
1371677d0b8SKris Buschelman /*
1381677d0b8SKris Buschelman    Note the r permutation is ignored
1391677d0b8SKris Buschelman */
1401677d0b8SKris Buschelman #undef __FUNCT__
141f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK"
142f0c56d0fSKris Buschelman int MatLUFactorSymbolic_UMFPACK(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
143e1b4c3a1SKris Buschelman   Mat         B;
1441677d0b8SKris Buschelman   Mat_SeqAIJ  *mat=(Mat_SeqAIJ*)A->data;
145f0c56d0fSKris Buschelman   Mat_UMFPACK *lu;
1460f6efc10SHong Zhang   int         ierr,m=A->m,n=A->n,*ai=mat->i,*aj=mat->j,status,*ra,idx;
1471677d0b8SKris Buschelman   PetscScalar *av=mat->a;
1480f6efc10SHong Zhang   const char  *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC","2BY2"},
1490f6efc10SHong Zhang               *scale[]={"NONE","SUM","MAX"};
1500f6efc10SHong Zhang   PetscTruth  flg;
1511677d0b8SKris Buschelman 
1521677d0b8SKris Buschelman   PetscFunctionBegin;
1531677d0b8SKris Buschelman   /* Create the factorization matrix F */
154e1b4c3a1SKris Buschelman   ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,m,n,&B);CHKERRQ(ierr);
155be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
156e1b4c3a1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
1571677d0b8SKris Buschelman 
158f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
159f0c56d0fSKris Buschelman   B->ops->solve           = MatSolve_UMFPACK;
160e1b4c3a1SKris Buschelman   B->factor               = FACTOR_LU;
16194b7f48cSBarry Smith   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */
1621677d0b8SKris Buschelman 
163f0c56d0fSKris Buschelman   lu = (Mat_UMFPACK*)(B->spptr);
1641677d0b8SKris Buschelman 
1651677d0b8SKris Buschelman   /* initializations */
1661677d0b8SKris Buschelman   /* ------------------------------------------------*/
1671677d0b8SKris Buschelman   /* get the default control parameters */
1681677d0b8SKris Buschelman   umfpack_di_defaults (lu->Control) ;
1691677d0b8SKris Buschelman   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
1700f6efc10SHong Zhang   lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */
1711677d0b8SKris Buschelman 
1721677d0b8SKris Buschelman   ierr = PetscOptionsBegin(A->comm,A->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
1731677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
1741677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr);
1751677d0b8SKris Buschelman 
1761677d0b8SKris Buschelman   /* Control parameters for symbolic factorization */
1770f6efc10SHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,4,strategy[0],&idx,&flg);CHKERRQ(ierr);
1780f6efc10SHong Zhang   if (flg) {
1790f6efc10SHong Zhang     switch (idx){
1800f6efc10SHong Zhang     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
1810f6efc10SHong Zhang     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
1820f6efc10SHong Zhang     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
1830f6efc10SHong Zhang     case 3: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_2BY2; break;
1840f6efc10SHong Zhang     }
1850f6efc10SHong Zhang   }
1861677d0b8SKris 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);
1871677d0b8SKris 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);
1880f6efc10SHong 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);
1891677d0b8SKris 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);
1900f6efc10SHong 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);
1910f6efc10SHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);CHKERRQ(ierr);
1920f6efc10SHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);CHKERRQ(ierr);
1931677d0b8SKris Buschelman 
1941677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
1951677d0b8SKris 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);
1960f6efc10SHong 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);
1970f6efc10SHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
1980f6efc10SHong Zhang   if (flg) {
1990f6efc10SHong Zhang     switch (idx){
2000f6efc10SHong Zhang     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
2010f6efc10SHong Zhang     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
2020f6efc10SHong Zhang     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
2030f6efc10SHong Zhang     }
2040f6efc10SHong Zhang   }
2051677d0b8SKris 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);
2060f6efc10SHong 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);
2070f6efc10SHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);CHKERRQ(ierr);
2081677d0b8SKris Buschelman 
2091677d0b8SKris Buschelman   /* Control parameters used by solve */
2101677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr);
2111677d0b8SKris Buschelman 
2120f6efc10SHong Zhang   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */
2131677d0b8SKris Buschelman   ierr = PetscOptionsHasName(PETSC_NULL,"-pc_lu_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr);
2141677d0b8SKris Buschelman   if (lu->PetscMatOdering) {
2150f6efc10SHong Zhang     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
2161677d0b8SKris Buschelman     ierr = PetscMalloc(A->m*sizeof(int),&lu->perm_c);CHKERRQ(ierr);
2170f6efc10SHong Zhang     ierr = PetscMemcpy(lu->perm_c,ra,A->m*sizeof(int));CHKERRQ(ierr);
2180f6efc10SHong Zhang     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
2191677d0b8SKris Buschelman   }
2201677d0b8SKris Buschelman   PetscOptionsEnd();
2211677d0b8SKris Buschelman 
2221677d0b8SKris Buschelman   /* print the control parameters */
2231677d0b8SKris Buschelman   if( lu->Control[UMFPACK_PRL] > 1 ) umfpack_di_report_control (lu->Control);
2241677d0b8SKris Buschelman 
2251677d0b8SKris Buschelman   /* symbolic factorization of A' */
2261677d0b8SKris Buschelman   /* ---------------------------------------------------------------------- */
2270f6efc10SHong Zhang   if (lu->PetscMatOdering) { /* use Petsc row ordering */
2280f6efc10SHong Zhang     status = umfpack_di_qsymbolic(n,m,ai,aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2290f6efc10SHong Zhang   } else { /* use Umfpack col ordering */
2300f6efc10SHong Zhang     status = umfpack_di_symbolic(n,m,ai,aj,av,&lu->Symbolic,lu->Control,lu->Info);
2310f6efc10SHong Zhang   }
2321677d0b8SKris Buschelman   if (status < 0){
2331677d0b8SKris Buschelman     umfpack_di_report_info(lu->Control, lu->Info) ;
2341677d0b8SKris Buschelman     umfpack_di_report_status(lu->Control, status) ;
2351677d0b8SKris Buschelman     SETERRQ(1,"umfpack_di_symbolic failed");
2361677d0b8SKris Buschelman   }
2371677d0b8SKris Buschelman   /* report sumbolic factorization of A' when Control[PRL] > 3 */
2381677d0b8SKris Buschelman   (void) umfpack_di_report_symbolic(lu->Symbolic, lu->Control) ;
2391677d0b8SKris Buschelman 
2401677d0b8SKris Buschelman   lu->flg = DIFFERENT_NONZERO_PATTERN;
2411677d0b8SKris Buschelman   lu->ai  = ai;
2421677d0b8SKris Buschelman   lu->aj  = aj;
2431677d0b8SKris Buschelman   lu->av  = av;
244e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK = PETSC_TRUE;
245e1b4c3a1SKris Buschelman   *F = B;
2461677d0b8SKris Buschelman   PetscFunctionReturn(0);
2471677d0b8SKris Buschelman }
2481677d0b8SKris Buschelman 
2491677d0b8SKris Buschelman #undef __FUNCT__
250f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_UMFPACK"
251f0c56d0fSKris Buschelman int MatAssemblyEnd_UMFPACK(Mat A,MatAssemblyType mode) {
252d844382dSKris Buschelman   int         ierr;
253f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)(A->spptr);
254d844382dSKris Buschelman 
2551677d0b8SKris Buschelman   PetscFunctionBegin;
256d844382dSKris Buschelman   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
257d844382dSKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
258f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
2591677d0b8SKris Buschelman   PetscFunctionReturn(0);
2601677d0b8SKris Buschelman }
2611677d0b8SKris Buschelman 
26294b7f48cSBarry Smith /* used by -ksp_view */
2631677d0b8SKris Buschelman #undef __FUNCT__
264f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_UMFPACK"
265f0c56d0fSKris Buschelman int MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer) {
266f0c56d0fSKris Buschelman   Mat_UMFPACK *lu= (Mat_UMFPACK*)A->spptr;
2671677d0b8SKris Buschelman   int         ierr;
268f0c56d0fSKris Buschelman 
2691677d0b8SKris Buschelman   PetscFunctionBegin;
2701677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
271f0c56d0fSKris Buschelman   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
2721677d0b8SKris Buschelman 
2731677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
2741677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
2751677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
2761677d0b8SKris Buschelman 
2771677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
2780f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr);
2791677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
2801677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
2810f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr);
2821677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
2830f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_2BY2_TOLERANCE]: %g\n",lu->Control[UMFPACK_2BY2_TOLERANCE]);CHKERRQ(ierr);
2840f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr);
2850f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr);
2861677d0b8SKris Buschelman 
2871677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
2881677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
2890f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr);
2900f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr);
2911677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
2920f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr);
2931677d0b8SKris Buschelman 
2941677d0b8SKris Buschelman   /* Control parameters used by solve */
2951677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
2961677d0b8SKris Buschelman 
2971677d0b8SKris Buschelman   /* mat ordering */
2981677d0b8SKris Buschelman   if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer,"  UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr);
2991677d0b8SKris Buschelman 
3001677d0b8SKris Buschelman   PetscFunctionReturn(0);
3011677d0b8SKris Buschelman }
3021677d0b8SKris Buschelman 
303d844382dSKris Buschelman #undef __FUNCT__
304f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_UMFPACK"
305f0c56d0fSKris Buschelman int MatView_UMFPACK(Mat A,PetscViewer viewer) {
306d844382dSKris Buschelman   int               ierr;
307*32077d6dSBarry Smith   PetscTruth        iascii;
308d844382dSKris Buschelman   PetscViewerFormat format;
309f0c56d0fSKris Buschelman   Mat_UMFPACK       *lu=(Mat_UMFPACK*)(A->spptr);
310d844382dSKris Buschelman 
311d844382dSKris Buschelman   PetscFunctionBegin;
312d844382dSKris Buschelman   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
313d844382dSKris Buschelman 
314*32077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
315*32077d6dSBarry Smith   if (iascii) {
316d844382dSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
317d844382dSKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
318f0c56d0fSKris Buschelman       ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
319d844382dSKris Buschelman     }
320d844382dSKris Buschelman   }
321d844382dSKris Buschelman   PetscFunctionReturn(0);
322d844382dSKris Buschelman }
323d844382dSKris Buschelman 
324d844382dSKris Buschelman EXTERN_C_BEGIN
325d844382dSKris Buschelman #undef __FUNCT__
326d844382dSKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_UMFPACK"
3278e9aea5cSBarry Smith int MatConvert_SeqAIJ_UMFPACK(Mat A,const MatType type,Mat *newmat) {
328d844382dSKris Buschelman   /* This routine is only called to convert to MATUMFPACK */
329d844382dSKris Buschelman   /* from MATSEQAIJ, so we will ignore 'MatType type'. */
330d844382dSKris Buschelman   int         ierr;
331d844382dSKris Buschelman   Mat         B=*newmat;
332f0c56d0fSKris Buschelman   Mat_UMFPACK *lu;
333d844382dSKris Buschelman 
334d844382dSKris Buschelman   PetscFunctionBegin;
335d844382dSKris Buschelman   if (B != A) {
336d844382dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
337d844382dSKris Buschelman   }
338d844382dSKris Buschelman 
339f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_UMFPACK,&lu);CHKERRQ(ierr);
340f0c56d0fSKris Buschelman   lu->MatDuplicate         = A->ops->duplicate;
341d844382dSKris Buschelman   lu->MatView              = A->ops->view;
342d844382dSKris Buschelman   lu->MatAssemblyEnd       = A->ops->assemblyend;
343d844382dSKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
344d844382dSKris Buschelman   lu->MatDestroy           = A->ops->destroy;
345d844382dSKris Buschelman   lu->CleanUpUMFPACK       = PETSC_FALSE;
346d844382dSKris Buschelman 
347d844382dSKris Buschelman   B->spptr                 = (void*)lu;
348f0c56d0fSKris Buschelman   B->ops->duplicate        = MatDuplicate_UMFPACK;
349f0c56d0fSKris Buschelman   B->ops->view             = MatView_UMFPACK;
350f0c56d0fSKris Buschelman   B->ops->assemblyend      = MatAssemblyEnd_UMFPACK;
351f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
352f0c56d0fSKris Buschelman   B->ops->destroy          = MatDestroy_UMFPACK;
353d844382dSKris Buschelman 
354d844382dSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_umfpack_C",
355d844382dSKris Buschelman                                            "MatConvert_SeqAIJ_UMFPACK",MatConvert_SeqAIJ_UMFPACK);CHKERRQ(ierr);
356d844382dSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_umfpack_seqaij_C",
357d844382dSKris Buschelman                                            "MatConvert_UMFPACK_SeqAIJ",MatConvert_UMFPACK_SeqAIJ);CHKERRQ(ierr);
358d844382dSKris Buschelman   PetscLogInfo(0,"Using UMFPACK for SeqAIJ LU factorization and solves.");
359d844382dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATUMFPACK);CHKERRQ(ierr);
360d844382dSKris Buschelman   *newmat = B;
361d844382dSKris Buschelman   PetscFunctionReturn(0);
362d844382dSKris Buschelman }
363d844382dSKris Buschelman EXTERN_C_END
364d844382dSKris Buschelman 
365f0c56d0fSKris Buschelman #undef __FUNCT__
366f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_UMFPACK"
367f0c56d0fSKris Buschelman int MatDuplicate_UMFPACK(Mat A, MatDuplicateOption op, Mat *M) {
368f0c56d0fSKris Buschelman   int         ierr;
3698f340917SKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
3708f340917SKris Buschelman 
371f0c56d0fSKris Buschelman   PetscFunctionBegin;
3728f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
3733f953163SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_UMFPACK));CHKERRQ(ierr);
374f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
375f0c56d0fSKris Buschelman }
376f0c56d0fSKris Buschelman 
3772f71e704SKris Buschelman /*MC
378fafad747SKris Buschelman   MATUMFPACK - MATUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
3792f71e704SKris Buschelman   via the external package UMFPACK.
3802f71e704SKris Buschelman 
3812f71e704SKris Buschelman   If UMFPACK is installed (see the manual for
3822f71e704SKris Buschelman   instructions on how to declare the existence of external packages),
3832f71e704SKris Buschelman   a matrix type can be constructed which invokes UMFPACK solvers.
3842f71e704SKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,UMFPACK).
3852f71e704SKris Buschelman   This matrix type is only supported for double precision real.
3862f71e704SKris Buschelman 
3872f71e704SKris Buschelman   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
38828b08bd3SKris Buschelman   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from
38928b08bd3SKris Buschelman   the MATSEQAIJ type without data copy.
3902f71e704SKris Buschelman 
3912f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
3922f71e704SKris Buschelman   which correspond to the options database keys below.
3932f71e704SKris Buschelman 
3942f71e704SKris Buschelman   Options Database Keys:
3950bad9183SKris Buschelman + -mat_type umfpack - sets the matrix type to "umfpack" during a call to MatSetFromOptions()
3962f71e704SKris Buschelman . -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL]
3972f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
3982f71e704SKris Buschelman . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
3992f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
4002f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
4012f71e704SKris Buschelman - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
4022f71e704SKris Buschelman 
4032f71e704SKris Buschelman    Level: beginner
4042f71e704SKris Buschelman 
4052f71e704SKris Buschelman .seealso: PCLU
4062f71e704SKris Buschelman M*/
4072f71e704SKris Buschelman 
4081677d0b8SKris Buschelman EXTERN_C_BEGIN
4091677d0b8SKris Buschelman #undef __FUNCT__
410f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_UMFPACK"
411f0c56d0fSKris Buschelman int MatCreate_UMFPACK(Mat A) {
4121677d0b8SKris Buschelman   int                ierr;
4131677d0b8SKris Buschelman 
4141677d0b8SKris Buschelman   PetscFunctionBegin;
4155441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ and UMFPACK types */
4165441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATUMFPACK);CHKERRQ(ierr);
4171677d0b8SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
418d844382dSKris Buschelman   ierr = MatConvert_SeqAIJ_UMFPACK(A,MATUMFPACK,&A);CHKERRQ(ierr);
4191677d0b8SKris Buschelman   PetscFunctionReturn(0);
4201677d0b8SKris Buschelman }
4211677d0b8SKris Buschelman EXTERN_C_END
422