xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision ae15b995b5732fffd2de5a75cf61ef7190c6fef1)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
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 */
226849ba73SBarry Smith   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
236849ba73SBarry Smith   PetscErrorCode (*MatView)(Mat,PetscViewer);
246849ba73SBarry Smith   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
256849ba73SBarry Smith   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
266849ba73SBarry Smith   PetscErrorCode (*MatDestroy)(Mat);
271677d0b8SKris Buschelman 
281677d0b8SKris Buschelman   /* Flag to clean up UMFPACK objects during Destroy */
291677d0b8SKris Buschelman   PetscTruth CleanUpUMFPACK;
30f0c56d0fSKris Buschelman } Mat_UMFPACK;
31f0c56d0fSKris Buschelman 
32dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_UMFPACK(Mat,MatDuplicateOption,Mat*);
33d844382dSKris Buschelman 
34d844382dSKris Buschelman EXTERN_C_BEGIN
35d844382dSKris Buschelman #undef __FUNCT__
36d844382dSKris Buschelman #define __FUNCT__ "MatConvert_UMFPACK_SeqAIJ"
3775179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_UMFPACK_SeqAIJ(Mat A,MatType type,MatReuse reuse,Mat *newmat)
386849ba73SBarry Smith {
39d844382dSKris Buschelman   /* This routine is only called to convert an unfactored PETSc-UMFPACK matrix */
40d844382dSKris Buschelman   /* to its base PETSc type, so we will ignore 'MatType type'. */
41dfbe8321SBarry Smith   PetscErrorCode ierr;
42d844382dSKris Buschelman   Mat         B=*newmat;
43f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
44d844382dSKris Buschelman 
45d844382dSKris Buschelman   PetscFunctionBegin;
46ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
47d844382dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
48f0c56d0fSKris Buschelman   }
49d844382dSKris Buschelman   /* Reset the original function pointers */
50901853e0SKris Buschelman   B->ops->duplicate        = lu->MatDuplicate;
51901853e0SKris Buschelman   B->ops->view             = lu->MatView;
52901853e0SKris Buschelman   B->ops->assemblyend      = lu->MatAssemblyEnd;
53901853e0SKris Buschelman   B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
54901853e0SKris Buschelman   B->ops->destroy          = lu->MatDestroy;
55d844382dSKris Buschelman 
56d844382dSKris Buschelman   ierr = PetscFree(lu);CHKERRQ(ierr);
57901853e0SKris Buschelman 
58901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_umfpack_C","",PETSC_NULL);CHKERRQ(ierr);
59901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_umfpack_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
60901853e0SKris Buschelman 
61d844382dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
62d844382dSKris Buschelman   *newmat = B;
63d844382dSKris Buschelman   PetscFunctionReturn(0);
64d844382dSKris Buschelman }
65d844382dSKris Buschelman EXTERN_C_END
661677d0b8SKris Buschelman 
671677d0b8SKris Buschelman #undef __FUNCT__
68f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_UMFPACK"
69dfbe8321SBarry Smith PetscErrorCode MatDestroy_UMFPACK(Mat A)
70dfbe8321SBarry Smith {
71dfbe8321SBarry Smith   PetscErrorCode ierr;
72f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
731677d0b8SKris Buschelman 
741677d0b8SKris Buschelman   PetscFunctionBegin;
75fb731535SKris Buschelman   if (lu->CleanUpUMFPACK) {
761677d0b8SKris Buschelman     umfpack_di_free_symbolic(&lu->Symbolic);
771677d0b8SKris Buschelman     umfpack_di_free_numeric(&lu->Numeric);
781677d0b8SKris Buschelman     ierr = PetscFree(lu->Wi);CHKERRQ(ierr);
791677d0b8SKris Buschelman     ierr = PetscFree(lu->W);CHKERRQ(ierr);
801677d0b8SKris Buschelman     if (lu->PetscMatOdering) {
811677d0b8SKris Buschelman       ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
821677d0b8SKris Buschelman     }
831677d0b8SKris Buschelman   }
84ceb03754SKris Buschelman   ierr = MatConvert_UMFPACK_SeqAIJ(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
85d844382dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
861677d0b8SKris Buschelman   PetscFunctionReturn(0);
871677d0b8SKris Buschelman }
881677d0b8SKris Buschelman 
891677d0b8SKris Buschelman #undef __FUNCT__
90f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_UMFPACK"
916849ba73SBarry Smith PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x)
926849ba73SBarry Smith {
93f0c56d0fSKris Buschelman   Mat_UMFPACK *lu = (Mat_UMFPACK*)A->spptr;
941677d0b8SKris Buschelman   PetscScalar *av=lu->av,*ba,*xa;
95dfbe8321SBarry Smith   PetscErrorCode ierr;
96dfbe8321SBarry Smith   int          *ai=lu->ai,*aj=lu->aj,status;
971677d0b8SKris Buschelman 
981677d0b8SKris Buschelman   PetscFunctionBegin;
991677d0b8SKris Buschelman   /* solve Ax = b by umfpack_di_wsolve */
1001677d0b8SKris Buschelman   /* ----------------------------------*/
1011677d0b8SKris Buschelman   ierr = VecGetArray(b,&ba);
1021677d0b8SKris Buschelman   ierr = VecGetArray(x,&xa);
1031677d0b8SKris Buschelman 
1041677d0b8SKris Buschelman   status = umfpack_di_wsolve(UMFPACK_At,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
1051677d0b8SKris Buschelman   umfpack_di_report_info(lu->Control, lu->Info);
1061677d0b8SKris Buschelman   if (status < 0){
1071677d0b8SKris Buschelman     umfpack_di_report_status(lu->Control, status);
108e005ede5SBarry Smith     SETERRQ(PETSC_ERR_LIB,"umfpack_di_wsolve failed");
1091677d0b8SKris Buschelman   }
1101677d0b8SKris Buschelman 
1111677d0b8SKris Buschelman   ierr = VecRestoreArray(b,&ba);
1121677d0b8SKris Buschelman   ierr = VecRestoreArray(x,&xa);
1132a325a84SHong Zhang 
1141677d0b8SKris Buschelman   PetscFunctionReturn(0);
1151677d0b8SKris Buschelman }
1161677d0b8SKris Buschelman 
1171677d0b8SKris Buschelman #undef __FUNCT__
118f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_UMFPACK"
119af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat A,MatFactorInfo *info,Mat *F)
1206849ba73SBarry Smith {
121f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)(*F)->spptr;
1226849ba73SBarry Smith   PetscErrorCode ierr;
1236849ba73SBarry Smith   int         *ai=lu->ai,*aj=lu->aj,m=A->m,status;
1241677d0b8SKris Buschelman   PetscScalar *av=lu->av;
1251677d0b8SKris Buschelman 
1261677d0b8SKris Buschelman   PetscFunctionBegin;
1271677d0b8SKris Buschelman   /* numeric factorization of A' */
1281677d0b8SKris Buschelman   /* ----------------------------*/
1292a325a84SHong Zhang   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
1302a325a84SHong Zhang       umfpack_di_free_numeric(&lu->Numeric);
1312a325a84SHong Zhang   }
1321677d0b8SKris Buschelman   status = umfpack_di_numeric (ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
133e005ede5SBarry Smith   if (status < 0) SETERRQ(PETSC_ERR_LIB,"umfpack_di_numeric failed");
1341677d0b8SKris Buschelman   /* report numeric factorization of A' when Control[PRL] > 3 */
1351677d0b8SKris Buschelman   (void) umfpack_di_report_numeric (lu->Numeric, lu->Control);
1361677d0b8SKris Buschelman 
1371677d0b8SKris Buschelman   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
1381677d0b8SKris Buschelman     /* allocate working space to be used by Solve */
1391677d0b8SKris Buschelman     ierr = PetscMalloc(m * sizeof(int), &lu->Wi);CHKERRQ(ierr);
1401677d0b8SKris Buschelman     ierr = PetscMalloc(5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
1411677d0b8SKris Buschelman   }
1422a325a84SHong Zhang   lu->flg = SAME_NONZERO_PATTERN;
1432a325a84SHong Zhang   lu->CleanUpUMFPACK = PETSC_TRUE;
1441677d0b8SKris Buschelman   PetscFunctionReturn(0);
1451677d0b8SKris Buschelman }
1461677d0b8SKris Buschelman 
1471677d0b8SKris Buschelman /*
1481677d0b8SKris Buschelman    Note the r permutation is ignored
1491677d0b8SKris Buschelman */
1501677d0b8SKris Buschelman #undef __FUNCT__
151f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK"
1526849ba73SBarry Smith PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
1536849ba73SBarry Smith {
154e1b4c3a1SKris Buschelman   Mat         B;
1551677d0b8SKris Buschelman   Mat_SeqAIJ  *mat=(Mat_SeqAIJ*)A->data;
156f0c56d0fSKris Buschelman   Mat_UMFPACK *lu;
157dfbe8321SBarry Smith   PetscErrorCode ierr;
158dfbe8321SBarry Smith   int          m=A->m,n=A->n,*ai=mat->i,*aj=mat->j,status,*ra,idx;
1591677d0b8SKris Buschelman   PetscScalar *av=mat->a;
1600f6efc10SHong Zhang   const char  *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC","2BY2"},
1610f6efc10SHong Zhang               *scale[]={"NONE","SUM","MAX"};
1620f6efc10SHong Zhang   PetscTruth  flg;
1631677d0b8SKris Buschelman 
1641677d0b8SKris Buschelman   PetscFunctionBegin;
1651677d0b8SKris Buschelman   /* Create the factorization matrix F */
166f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
167f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
168be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
169e1b4c3a1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
1701677d0b8SKris Buschelman 
171f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
172f0c56d0fSKris Buschelman   B->ops->solve           = MatSolve_UMFPACK;
173e1b4c3a1SKris Buschelman   B->factor               = FACTOR_LU;
17494b7f48cSBarry Smith   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */
1751677d0b8SKris Buschelman 
176f0c56d0fSKris Buschelman   lu = (Mat_UMFPACK*)(B->spptr);
1771677d0b8SKris Buschelman 
1781677d0b8SKris Buschelman   /* initializations */
1791677d0b8SKris Buschelman   /* ------------------------------------------------*/
1801677d0b8SKris Buschelman   /* get the default control parameters */
1811677d0b8SKris Buschelman   umfpack_di_defaults (lu->Control);
1821677d0b8SKris Buschelman   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
1830f6efc10SHong Zhang   lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */
1841677d0b8SKris Buschelman 
1851677d0b8SKris Buschelman   ierr = PetscOptionsBegin(A->comm,A->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
1861677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
1871677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr);
1881677d0b8SKris Buschelman 
1891677d0b8SKris Buschelman   /* Control parameters for symbolic factorization */
1900f6efc10SHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,4,strategy[0],&idx,&flg);CHKERRQ(ierr);
1910f6efc10SHong Zhang   if (flg) {
1920f6efc10SHong Zhang     switch (idx){
1930f6efc10SHong Zhang     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
1940f6efc10SHong Zhang     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
1950f6efc10SHong Zhang     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
1960f6efc10SHong Zhang     case 3: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_2BY2; break;
1970f6efc10SHong Zhang     }
1980f6efc10SHong Zhang   }
1991677d0b8SKris 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);
2001677d0b8SKris 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);
2010f6efc10SHong 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);
2021677d0b8SKris 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);
2030f6efc10SHong 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);
2040f6efc10SHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);CHKERRQ(ierr);
2050f6efc10SHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);CHKERRQ(ierr);
2061677d0b8SKris Buschelman 
2071677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
2081677d0b8SKris 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);
2090f6efc10SHong 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);
2100f6efc10SHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
2110f6efc10SHong Zhang   if (flg) {
2120f6efc10SHong Zhang     switch (idx){
2130f6efc10SHong Zhang     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
2140f6efc10SHong Zhang     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
2150f6efc10SHong Zhang     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
2160f6efc10SHong Zhang     }
2170f6efc10SHong Zhang   }
2181677d0b8SKris 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);
2190f6efc10SHong 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);
2200f6efc10SHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);CHKERRQ(ierr);
2211677d0b8SKris Buschelman 
2221677d0b8SKris Buschelman   /* Control parameters used by solve */
2231677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr);
2241677d0b8SKris Buschelman 
2250f6efc10SHong Zhang   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */
2262401956bSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-pc_factor_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr);
2271677d0b8SKris Buschelman   if (lu->PetscMatOdering) {
2280f6efc10SHong Zhang     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
2291677d0b8SKris Buschelman     ierr = PetscMalloc(A->m*sizeof(int),&lu->perm_c);CHKERRQ(ierr);
2300f6efc10SHong Zhang     ierr = PetscMemcpy(lu->perm_c,ra,A->m*sizeof(int));CHKERRQ(ierr);
2310f6efc10SHong Zhang     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
2321677d0b8SKris Buschelman   }
2331677d0b8SKris Buschelman   PetscOptionsEnd();
2341677d0b8SKris Buschelman 
2351677d0b8SKris Buschelman   /* print the control parameters */
2361677d0b8SKris Buschelman   if( lu->Control[UMFPACK_PRL] > 1 ) umfpack_di_report_control (lu->Control);
2371677d0b8SKris Buschelman 
2381677d0b8SKris Buschelman   /* symbolic factorization of A' */
2391677d0b8SKris Buschelman   /* ---------------------------------------------------------------------- */
2400f6efc10SHong Zhang   if (lu->PetscMatOdering) { /* use Petsc row ordering */
2410f6efc10SHong Zhang     status = umfpack_di_qsymbolic(n,m,ai,aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2420f6efc10SHong Zhang   } else { /* use Umfpack col ordering */
2430f6efc10SHong Zhang     status = umfpack_di_symbolic(n,m,ai,aj,av,&lu->Symbolic,lu->Control,lu->Info);
2440f6efc10SHong Zhang   }
2451677d0b8SKris Buschelman   if (status < 0){
2461677d0b8SKris Buschelman     umfpack_di_report_info(lu->Control, lu->Info);
2471677d0b8SKris Buschelman     umfpack_di_report_status(lu->Control, status);
248e005ede5SBarry Smith     SETERRQ(PETSC_ERR_LIB,"umfpack_di_symbolic failed");
2491677d0b8SKris Buschelman   }
2501677d0b8SKris Buschelman   /* report sumbolic factorization of A' when Control[PRL] > 3 */
2511677d0b8SKris Buschelman   (void) umfpack_di_report_symbolic(lu->Symbolic, lu->Control);
2521677d0b8SKris Buschelman 
2531677d0b8SKris Buschelman   lu->flg = DIFFERENT_NONZERO_PATTERN;
2541677d0b8SKris Buschelman   lu->ai  = ai;
2551677d0b8SKris Buschelman   lu->aj  = aj;
2561677d0b8SKris Buschelman   lu->av  = av;
257e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK = PETSC_TRUE;
258e1b4c3a1SKris Buschelman   *F = B;
2591677d0b8SKris Buschelman   PetscFunctionReturn(0);
2601677d0b8SKris Buschelman }
2611677d0b8SKris Buschelman 
2621677d0b8SKris Buschelman #undef __FUNCT__
263f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_UMFPACK"
2646849ba73SBarry Smith PetscErrorCode MatAssemblyEnd_UMFPACK(Mat A,MatAssemblyType mode)
2656849ba73SBarry Smith {
266dfbe8321SBarry Smith   PetscErrorCode ierr;
267f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)(A->spptr);
268d844382dSKris Buschelman 
2691677d0b8SKris Buschelman   PetscFunctionBegin;
270d844382dSKris Buschelman   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
271d844382dSKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
272f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
2731677d0b8SKris Buschelman   PetscFunctionReturn(0);
2741677d0b8SKris Buschelman }
2751677d0b8SKris Buschelman 
27694b7f48cSBarry Smith /* used by -ksp_view */
2771677d0b8SKris Buschelman #undef __FUNCT__
278f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_UMFPACK"
2796849ba73SBarry Smith PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer)
2806849ba73SBarry Smith {
281f0c56d0fSKris Buschelman   Mat_UMFPACK *lu= (Mat_UMFPACK*)A->spptr;
282dfbe8321SBarry Smith   PetscErrorCode ierr;
283f0c56d0fSKris Buschelman 
2841677d0b8SKris Buschelman   PetscFunctionBegin;
2851677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
286f0c56d0fSKris Buschelman   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
2871677d0b8SKris Buschelman 
2881677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
2891677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
2901677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
2911677d0b8SKris Buschelman 
2921677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
2930f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr);
2941677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
2951677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
2960f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr);
2971677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
2980f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_2BY2_TOLERANCE]: %g\n",lu->Control[UMFPACK_2BY2_TOLERANCE]);CHKERRQ(ierr);
2990f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr);
3000f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr);
3011677d0b8SKris Buschelman 
3021677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
3031677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
3040f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr);
3050f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr);
3061677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
3070f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr);
3081677d0b8SKris Buschelman 
3091677d0b8SKris Buschelman   /* Control parameters used by solve */
3101677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
3111677d0b8SKris Buschelman 
3121677d0b8SKris Buschelman   /* mat ordering */
3131677d0b8SKris Buschelman   if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer,"  UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr);
3141677d0b8SKris Buschelman 
3151677d0b8SKris Buschelman   PetscFunctionReturn(0);
3161677d0b8SKris Buschelman }
3171677d0b8SKris Buschelman 
318d844382dSKris Buschelman #undef __FUNCT__
319f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_UMFPACK"
3206849ba73SBarry Smith PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
3216849ba73SBarry Smith {
322dfbe8321SBarry Smith   PetscErrorCode ierr;
32332077d6dSBarry Smith   PetscTruth        iascii;
324d844382dSKris Buschelman   PetscViewerFormat format;
325f0c56d0fSKris Buschelman   Mat_UMFPACK       *lu=(Mat_UMFPACK*)(A->spptr);
326d844382dSKris Buschelman 
327d844382dSKris Buschelman   PetscFunctionBegin;
328d844382dSKris Buschelman   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
329d844382dSKris Buschelman 
33032077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
33132077d6dSBarry Smith   if (iascii) {
332d844382dSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
333d844382dSKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
334f0c56d0fSKris Buschelman       ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
335d844382dSKris Buschelman     }
336d844382dSKris Buschelman   }
337d844382dSKris Buschelman   PetscFunctionReturn(0);
338d844382dSKris Buschelman }
339d844382dSKris Buschelman 
340d844382dSKris Buschelman EXTERN_C_BEGIN
341d844382dSKris Buschelman #undef __FUNCT__
342d844382dSKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_UMFPACK"
34375179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_UMFPACK(Mat A,MatType type,MatReuse reuse,Mat *newmat)
3446849ba73SBarry Smith {
345d844382dSKris Buschelman   /* This routine is only called to convert to MATUMFPACK */
346d844382dSKris Buschelman   /* from MATSEQAIJ, so we will ignore 'MatType type'. */
347dfbe8321SBarry Smith   PetscErrorCode ierr;
348d844382dSKris Buschelman   Mat         B=*newmat;
349f0c56d0fSKris Buschelman   Mat_UMFPACK *lu;
350d844382dSKris Buschelman 
351d844382dSKris Buschelman   PetscFunctionBegin;
352ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
353d844382dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
354d844382dSKris Buschelman   }
355d844382dSKris Buschelman 
356f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_UMFPACK,&lu);CHKERRQ(ierr);
357f0c56d0fSKris Buschelman   lu->MatDuplicate         = A->ops->duplicate;
358d844382dSKris Buschelman   lu->MatView              = A->ops->view;
359d844382dSKris Buschelman   lu->MatAssemblyEnd       = A->ops->assemblyend;
360d844382dSKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
361d844382dSKris Buschelman   lu->MatDestroy           = A->ops->destroy;
362d844382dSKris Buschelman   lu->CleanUpUMFPACK       = PETSC_FALSE;
363d844382dSKris Buschelman 
364d844382dSKris Buschelman   B->spptr                 = (void*)lu;
365f0c56d0fSKris Buschelman   B->ops->duplicate        = MatDuplicate_UMFPACK;
366f0c56d0fSKris Buschelman   B->ops->view             = MatView_UMFPACK;
367f0c56d0fSKris Buschelman   B->ops->assemblyend      = MatAssemblyEnd_UMFPACK;
368f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
369f0c56d0fSKris Buschelman   B->ops->destroy          = MatDestroy_UMFPACK;
370d844382dSKris Buschelman 
371d844382dSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_umfpack_C",
372d844382dSKris Buschelman                                            "MatConvert_SeqAIJ_UMFPACK",MatConvert_SeqAIJ_UMFPACK);CHKERRQ(ierr);
373d844382dSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_umfpack_seqaij_C",
374d844382dSKris Buschelman                                            "MatConvert_UMFPACK_SeqAIJ",MatConvert_UMFPACK_SeqAIJ);CHKERRQ(ierr);
375*ae15b995SBarry Smith   ierr = PetscInfo(0,"Using UMFPACK for SeqAIJ LU factorization and solves.\n");CHKERRQ(ierr);
376d844382dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATUMFPACK);CHKERRQ(ierr);
377d844382dSKris Buschelman   *newmat = B;
378d844382dSKris Buschelman   PetscFunctionReturn(0);
379d844382dSKris Buschelman }
380d844382dSKris Buschelman EXTERN_C_END
381d844382dSKris Buschelman 
382f0c56d0fSKris Buschelman #undef __FUNCT__
383f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_UMFPACK"
3846849ba73SBarry Smith PetscErrorCode MatDuplicate_UMFPACK(Mat A, MatDuplicateOption op, Mat *M)
3856849ba73SBarry Smith {
386dfbe8321SBarry Smith   PetscErrorCode ierr;
3878f340917SKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
3888f340917SKris Buschelman 
389f0c56d0fSKris Buschelman   PetscFunctionBegin;
3908f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
3913f953163SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_UMFPACK));CHKERRQ(ierr);
392f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
393f0c56d0fSKris Buschelman }
394f0c56d0fSKris Buschelman 
3952f71e704SKris Buschelman /*MC
396fafad747SKris Buschelman   MATUMFPACK - MATUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
3972f71e704SKris Buschelman   via the external package UMFPACK.
3982f71e704SKris Buschelman 
3992f71e704SKris Buschelman   If UMFPACK is installed (see the manual for
4002f71e704SKris Buschelman   instructions on how to declare the existence of external packages),
4012f71e704SKris Buschelman   a matrix type can be constructed which invokes UMFPACK solvers.
4022f71e704SKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,UMFPACK).
4032f71e704SKris Buschelman   This matrix type is only supported for double precision real.
4042f71e704SKris Buschelman 
4052f71e704SKris Buschelman   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
40628b08bd3SKris Buschelman   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from
40728b08bd3SKris Buschelman   the MATSEQAIJ type without data copy.
4082f71e704SKris Buschelman 
4092f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
4102f71e704SKris Buschelman   which correspond to the options database keys below.
4112f71e704SKris Buschelman 
4122f71e704SKris Buschelman   Options Database Keys:
4130bad9183SKris Buschelman + -mat_type umfpack - sets the matrix type to "umfpack" during a call to MatSetFromOptions()
4142f71e704SKris Buschelman . -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL]
4152f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
4162f71e704SKris Buschelman . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
4172f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
4182f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
4192f71e704SKris Buschelman - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
4202f71e704SKris Buschelman 
4212f71e704SKris Buschelman    Level: beginner
4222f71e704SKris Buschelman 
4232f71e704SKris Buschelman .seealso: PCLU
4242f71e704SKris Buschelman M*/
4252f71e704SKris Buschelman 
4261677d0b8SKris Buschelman EXTERN_C_BEGIN
4271677d0b8SKris Buschelman #undef __FUNCT__
428f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_UMFPACK"
429be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_UMFPACK(Mat A)
430dfbe8321SBarry Smith {
431dfbe8321SBarry Smith   PetscErrorCode ierr;
4321677d0b8SKris Buschelman 
4331677d0b8SKris Buschelman   PetscFunctionBegin;
4345441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ and UMFPACK types */
4355441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATUMFPACK);CHKERRQ(ierr);
4361677d0b8SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
437ceb03754SKris Buschelman   ierr = MatConvert_SeqAIJ_UMFPACK(A,MATUMFPACK,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
4381677d0b8SKris Buschelman   PetscFunctionReturn(0);
4391677d0b8SKris Buschelman }
4401677d0b8SKris Buschelman EXTERN_C_END
441