xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
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 */
21*6849ba73SBarry Smith   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
22*6849ba73SBarry Smith   PetscErrorCode (*MatView)(Mat,PetscViewer);
23*6849ba73SBarry Smith   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
24*6849ba73SBarry Smith   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
25*6849ba73SBarry Smith   PetscErrorCode (*MatDestroy)(Mat);
261677d0b8SKris Buschelman 
271677d0b8SKris Buschelman   /* Flag to clean up UMFPACK objects during Destroy */
281677d0b8SKris Buschelman   PetscTruth CleanUpUMFPACK;
29f0c56d0fSKris Buschelman } Mat_UMFPACK;
30f0c56d0fSKris Buschelman 
31dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_UMFPACK(Mat,MatDuplicateOption,Mat*);
32d844382dSKris Buschelman 
33d844382dSKris Buschelman EXTERN_C_BEGIN
34d844382dSKris Buschelman #undef __FUNCT__
35d844382dSKris Buschelman #define __FUNCT__ "MatConvert_UMFPACK_SeqAIJ"
36*6849ba73SBarry Smith PetscErrorCode MatConvert_UMFPACK_SeqAIJ(Mat A,const MatType type,Mat *newmat)
37*6849ba73SBarry Smith {
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'. */
40dfbe8321SBarry Smith   PetscErrorCode 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"
65dfbe8321SBarry Smith PetscErrorCode MatDestroy_UMFPACK(Mat A)
66dfbe8321SBarry Smith {
67dfbe8321SBarry Smith   PetscErrorCode ierr;
68f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
691677d0b8SKris Buschelman 
701677d0b8SKris Buschelman   PetscFunctionBegin;
71fb731535SKris Buschelman   if (lu->CleanUpUMFPACK) {
721677d0b8SKris Buschelman     umfpack_di_free_symbolic(&lu->Symbolic);
731677d0b8SKris Buschelman     umfpack_di_free_numeric(&lu->Numeric);
741677d0b8SKris Buschelman     ierr = PetscFree(lu->Wi);CHKERRQ(ierr);
751677d0b8SKris Buschelman     ierr = PetscFree(lu->W);CHKERRQ(ierr);
761677d0b8SKris Buschelman     if (lu->PetscMatOdering) {
771677d0b8SKris Buschelman       ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
781677d0b8SKris Buschelman     }
791677d0b8SKris Buschelman   }
80d844382dSKris Buschelman   ierr = MatConvert_UMFPACK_SeqAIJ(A,MATSEQAIJ,&A);CHKERRQ(ierr);
81d844382dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
821677d0b8SKris Buschelman   PetscFunctionReturn(0);
831677d0b8SKris Buschelman }
841677d0b8SKris Buschelman 
851677d0b8SKris Buschelman #undef __FUNCT__
86f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_UMFPACK"
87*6849ba73SBarry Smith PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x)
88*6849ba73SBarry Smith {
89f0c56d0fSKris Buschelman   Mat_UMFPACK *lu = (Mat_UMFPACK*)A->spptr;
901677d0b8SKris Buschelman   PetscScalar *av=lu->av,*ba,*xa;
91dfbe8321SBarry Smith   PetscErrorCode ierr;
92dfbe8321SBarry Smith   int          *ai=lu->ai,*aj=lu->aj,status;
931677d0b8SKris Buschelman 
941677d0b8SKris Buschelman   PetscFunctionBegin;
951677d0b8SKris Buschelman   /* solve Ax = b by umfpack_di_wsolve */
961677d0b8SKris Buschelman   /* ----------------------------------*/
971677d0b8SKris Buschelman   ierr = VecGetArray(b,&ba);
981677d0b8SKris Buschelman   ierr = VecGetArray(x,&xa);
991677d0b8SKris Buschelman 
1001677d0b8SKris Buschelman   status = umfpack_di_wsolve(UMFPACK_At,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
1011677d0b8SKris Buschelman   umfpack_di_report_info(lu->Control, lu->Info);
1021677d0b8SKris Buschelman   if (status < 0){
1031677d0b8SKris Buschelman     umfpack_di_report_status(lu->Control, status);
1041677d0b8SKris Buschelman     SETERRQ(1,"umfpack_di_wsolve failed");
1051677d0b8SKris Buschelman   }
1061677d0b8SKris Buschelman 
1071677d0b8SKris Buschelman   ierr = VecRestoreArray(b,&ba);
1081677d0b8SKris Buschelman   ierr = VecRestoreArray(x,&xa);
1092a325a84SHong Zhang 
1101677d0b8SKris Buschelman   PetscFunctionReturn(0);
1111677d0b8SKris Buschelman }
1121677d0b8SKris Buschelman 
1131677d0b8SKris Buschelman #undef __FUNCT__
114f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_UMFPACK"
115*6849ba73SBarry Smith PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat A,Mat *F)
116*6849ba73SBarry Smith {
117f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)(*F)->spptr;
118*6849ba73SBarry Smith   PetscErrorCode ierr;
119*6849ba73SBarry Smith   int         *ai=lu->ai,*aj=lu->aj,m=A->m,status;
1201677d0b8SKris Buschelman   PetscScalar *av=lu->av;
1211677d0b8SKris Buschelman 
1221677d0b8SKris Buschelman   PetscFunctionBegin;
1231677d0b8SKris Buschelman   /* numeric factorization of A' */
1241677d0b8SKris Buschelman   /* ----------------------------*/
1252a325a84SHong Zhang   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
1262a325a84SHong Zhang       umfpack_di_free_numeric(&lu->Numeric);
1272a325a84SHong Zhang   }
1281677d0b8SKris Buschelman   status = umfpack_di_numeric (ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1291677d0b8SKris Buschelman   if (status < 0) SETERRQ(1,"umfpack_di_numeric failed");
1301677d0b8SKris Buschelman   /* report numeric factorization of A' when Control[PRL] > 3 */
1311677d0b8SKris Buschelman   (void) umfpack_di_report_numeric (lu->Numeric, lu->Control);
1321677d0b8SKris Buschelman 
1331677d0b8SKris Buschelman   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
1341677d0b8SKris Buschelman     /* allocate working space to be used by Solve */
1351677d0b8SKris Buschelman     ierr = PetscMalloc(m * sizeof(int), &lu->Wi);CHKERRQ(ierr);
1361677d0b8SKris Buschelman     ierr = PetscMalloc(5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
1371677d0b8SKris Buschelman   }
1382a325a84SHong Zhang   lu->flg = SAME_NONZERO_PATTERN;
1392a325a84SHong Zhang   lu->CleanUpUMFPACK = PETSC_TRUE;
1401677d0b8SKris Buschelman   PetscFunctionReturn(0);
1411677d0b8SKris Buschelman }
1421677d0b8SKris Buschelman 
1431677d0b8SKris Buschelman /*
1441677d0b8SKris Buschelman    Note the r permutation is ignored
1451677d0b8SKris Buschelman */
1461677d0b8SKris Buschelman #undef __FUNCT__
147f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK"
148*6849ba73SBarry Smith PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
149*6849ba73SBarry Smith {
150e1b4c3a1SKris Buschelman   Mat         B;
1511677d0b8SKris Buschelman   Mat_SeqAIJ  *mat=(Mat_SeqAIJ*)A->data;
152f0c56d0fSKris Buschelman   Mat_UMFPACK *lu;
153dfbe8321SBarry Smith   PetscErrorCode ierr;
154dfbe8321SBarry Smith   int          m=A->m,n=A->n,*ai=mat->i,*aj=mat->j,status,*ra,idx;
1551677d0b8SKris Buschelman   PetscScalar *av=mat->a;
1560f6efc10SHong Zhang   const char  *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC","2BY2"},
1570f6efc10SHong Zhang               *scale[]={"NONE","SUM","MAX"};
1580f6efc10SHong Zhang   PetscTruth  flg;
1591677d0b8SKris Buschelman 
1601677d0b8SKris Buschelman   PetscFunctionBegin;
1611677d0b8SKris Buschelman   /* Create the factorization matrix F */
162e1b4c3a1SKris Buschelman   ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,m,n,&B);CHKERRQ(ierr);
163be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
164e1b4c3a1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
1651677d0b8SKris Buschelman 
166f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
167f0c56d0fSKris Buschelman   B->ops->solve           = MatSolve_UMFPACK;
168e1b4c3a1SKris Buschelman   B->factor               = FACTOR_LU;
16994b7f48cSBarry Smith   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */
1701677d0b8SKris Buschelman 
171f0c56d0fSKris Buschelman   lu = (Mat_UMFPACK*)(B->spptr);
1721677d0b8SKris Buschelman 
1731677d0b8SKris Buschelman   /* initializations */
1741677d0b8SKris Buschelman   /* ------------------------------------------------*/
1751677d0b8SKris Buschelman   /* get the default control parameters */
1761677d0b8SKris Buschelman   umfpack_di_defaults (lu->Control);
1771677d0b8SKris Buschelman   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
1780f6efc10SHong Zhang   lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */
1791677d0b8SKris Buschelman 
1801677d0b8SKris Buschelman   ierr = PetscOptionsBegin(A->comm,A->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
1811677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
1821677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr);
1831677d0b8SKris Buschelman 
1841677d0b8SKris Buschelman   /* Control parameters for symbolic factorization */
1850f6efc10SHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,4,strategy[0],&idx,&flg);CHKERRQ(ierr);
1860f6efc10SHong Zhang   if (flg) {
1870f6efc10SHong Zhang     switch (idx){
1880f6efc10SHong Zhang     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
1890f6efc10SHong Zhang     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
1900f6efc10SHong Zhang     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
1910f6efc10SHong Zhang     case 3: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_2BY2; break;
1920f6efc10SHong Zhang     }
1930f6efc10SHong Zhang   }
1941677d0b8SKris 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);
1951677d0b8SKris 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);
1960f6efc10SHong 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);
1971677d0b8SKris 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);
1980f6efc10SHong 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);
1990f6efc10SHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);CHKERRQ(ierr);
2000f6efc10SHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);CHKERRQ(ierr);
2011677d0b8SKris Buschelman 
2021677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
2031677d0b8SKris 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);
2040f6efc10SHong 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);
2050f6efc10SHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
2060f6efc10SHong Zhang   if (flg) {
2070f6efc10SHong Zhang     switch (idx){
2080f6efc10SHong Zhang     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
2090f6efc10SHong Zhang     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
2100f6efc10SHong Zhang     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
2110f6efc10SHong Zhang     }
2120f6efc10SHong Zhang   }
2131677d0b8SKris 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);
2140f6efc10SHong 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);
2150f6efc10SHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);CHKERRQ(ierr);
2161677d0b8SKris Buschelman 
2171677d0b8SKris Buschelman   /* Control parameters used by solve */
2181677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr);
2191677d0b8SKris Buschelman 
2200f6efc10SHong Zhang   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */
2211677d0b8SKris Buschelman   ierr = PetscOptionsHasName(PETSC_NULL,"-pc_lu_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr);
2221677d0b8SKris Buschelman   if (lu->PetscMatOdering) {
2230f6efc10SHong Zhang     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
2241677d0b8SKris Buschelman     ierr = PetscMalloc(A->m*sizeof(int),&lu->perm_c);CHKERRQ(ierr);
2250f6efc10SHong Zhang     ierr = PetscMemcpy(lu->perm_c,ra,A->m*sizeof(int));CHKERRQ(ierr);
2260f6efc10SHong Zhang     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
2271677d0b8SKris Buschelman   }
2281677d0b8SKris Buschelman   PetscOptionsEnd();
2291677d0b8SKris Buschelman 
2301677d0b8SKris Buschelman   /* print the control parameters */
2311677d0b8SKris Buschelman   if( lu->Control[UMFPACK_PRL] > 1 ) umfpack_di_report_control (lu->Control);
2321677d0b8SKris Buschelman 
2331677d0b8SKris Buschelman   /* symbolic factorization of A' */
2341677d0b8SKris Buschelman   /* ---------------------------------------------------------------------- */
2350f6efc10SHong Zhang   if (lu->PetscMatOdering) { /* use Petsc row ordering */
2360f6efc10SHong Zhang     status = umfpack_di_qsymbolic(n,m,ai,aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2370f6efc10SHong Zhang   } else { /* use Umfpack col ordering */
2380f6efc10SHong Zhang     status = umfpack_di_symbolic(n,m,ai,aj,av,&lu->Symbolic,lu->Control,lu->Info);
2390f6efc10SHong Zhang   }
2401677d0b8SKris Buschelman   if (status < 0){
2411677d0b8SKris Buschelman     umfpack_di_report_info(lu->Control, lu->Info);
2421677d0b8SKris Buschelman     umfpack_di_report_status(lu->Control, status);
2431677d0b8SKris Buschelman     SETERRQ(1,"umfpack_di_symbolic failed");
2441677d0b8SKris Buschelman   }
2451677d0b8SKris Buschelman   /* report sumbolic factorization of A' when Control[PRL] > 3 */
2461677d0b8SKris Buschelman   (void) umfpack_di_report_symbolic(lu->Symbolic, lu->Control);
2471677d0b8SKris Buschelman 
2481677d0b8SKris Buschelman   lu->flg = DIFFERENT_NONZERO_PATTERN;
2491677d0b8SKris Buschelman   lu->ai  = ai;
2501677d0b8SKris Buschelman   lu->aj  = aj;
2511677d0b8SKris Buschelman   lu->av  = av;
252e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK = PETSC_TRUE;
253e1b4c3a1SKris Buschelman   *F = B;
2541677d0b8SKris Buschelman   PetscFunctionReturn(0);
2551677d0b8SKris Buschelman }
2561677d0b8SKris Buschelman 
2571677d0b8SKris Buschelman #undef __FUNCT__
258f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_UMFPACK"
259*6849ba73SBarry Smith PetscErrorCode MatAssemblyEnd_UMFPACK(Mat A,MatAssemblyType mode)
260*6849ba73SBarry Smith {
261dfbe8321SBarry Smith   PetscErrorCode ierr;
262f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)(A->spptr);
263d844382dSKris Buschelman 
2641677d0b8SKris Buschelman   PetscFunctionBegin;
265d844382dSKris Buschelman   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
266d844382dSKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
267f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
2681677d0b8SKris Buschelman   PetscFunctionReturn(0);
2691677d0b8SKris Buschelman }
2701677d0b8SKris Buschelman 
27194b7f48cSBarry Smith /* used by -ksp_view */
2721677d0b8SKris Buschelman #undef __FUNCT__
273f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_UMFPACK"
274*6849ba73SBarry Smith PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer)
275*6849ba73SBarry Smith {
276f0c56d0fSKris Buschelman   Mat_UMFPACK *lu= (Mat_UMFPACK*)A->spptr;
277dfbe8321SBarry Smith   PetscErrorCode ierr;
278f0c56d0fSKris Buschelman 
2791677d0b8SKris Buschelman   PetscFunctionBegin;
2801677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
281f0c56d0fSKris Buschelman   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
2821677d0b8SKris Buschelman 
2831677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
2841677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
2851677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
2861677d0b8SKris Buschelman 
2871677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
2880f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr);
2891677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
2901677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
2910f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr);
2921677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
2930f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_2BY2_TOLERANCE]: %g\n",lu->Control[UMFPACK_2BY2_TOLERANCE]);CHKERRQ(ierr);
2940f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr);
2950f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr);
2961677d0b8SKris Buschelman 
2971677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
2981677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
2990f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr);
3000f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr);
3011677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
3020f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr);
3031677d0b8SKris Buschelman 
3041677d0b8SKris Buschelman   /* Control parameters used by solve */
3051677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
3061677d0b8SKris Buschelman 
3071677d0b8SKris Buschelman   /* mat ordering */
3081677d0b8SKris Buschelman   if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer,"  UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr);
3091677d0b8SKris Buschelman 
3101677d0b8SKris Buschelman   PetscFunctionReturn(0);
3111677d0b8SKris Buschelman }
3121677d0b8SKris Buschelman 
313d844382dSKris Buschelman #undef __FUNCT__
314f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_UMFPACK"
315*6849ba73SBarry Smith PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
316*6849ba73SBarry Smith {
317dfbe8321SBarry Smith   PetscErrorCode ierr;
31832077d6dSBarry Smith   PetscTruth        iascii;
319d844382dSKris Buschelman   PetscViewerFormat format;
320f0c56d0fSKris Buschelman   Mat_UMFPACK       *lu=(Mat_UMFPACK*)(A->spptr);
321d844382dSKris Buschelman 
322d844382dSKris Buschelman   PetscFunctionBegin;
323d844382dSKris Buschelman   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
324d844382dSKris Buschelman 
32532077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
32632077d6dSBarry Smith   if (iascii) {
327d844382dSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
328d844382dSKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
329f0c56d0fSKris Buschelman       ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
330d844382dSKris Buschelman     }
331d844382dSKris Buschelman   }
332d844382dSKris Buschelman   PetscFunctionReturn(0);
333d844382dSKris Buschelman }
334d844382dSKris Buschelman 
335d844382dSKris Buschelman EXTERN_C_BEGIN
336d844382dSKris Buschelman #undef __FUNCT__
337d844382dSKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_UMFPACK"
338*6849ba73SBarry Smith PetscErrorCode MatConvert_SeqAIJ_UMFPACK(Mat A,const MatType type,Mat *newmat)
339*6849ba73SBarry Smith {
340d844382dSKris Buschelman   /* This routine is only called to convert to MATUMFPACK */
341d844382dSKris Buschelman   /* from MATSEQAIJ, so we will ignore 'MatType type'. */
342dfbe8321SBarry Smith   PetscErrorCode ierr;
343d844382dSKris Buschelman   Mat         B=*newmat;
344f0c56d0fSKris Buschelman   Mat_UMFPACK *lu;
345d844382dSKris Buschelman 
346d844382dSKris Buschelman   PetscFunctionBegin;
347d844382dSKris Buschelman   if (B != A) {
348d844382dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
349d844382dSKris Buschelman   }
350d844382dSKris Buschelman 
351f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_UMFPACK,&lu);CHKERRQ(ierr);
352f0c56d0fSKris Buschelman   lu->MatDuplicate         = A->ops->duplicate;
353d844382dSKris Buschelman   lu->MatView              = A->ops->view;
354d844382dSKris Buschelman   lu->MatAssemblyEnd       = A->ops->assemblyend;
355d844382dSKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
356d844382dSKris Buschelman   lu->MatDestroy           = A->ops->destroy;
357d844382dSKris Buschelman   lu->CleanUpUMFPACK       = PETSC_FALSE;
358d844382dSKris Buschelman 
359d844382dSKris Buschelman   B->spptr                 = (void*)lu;
360f0c56d0fSKris Buschelman   B->ops->duplicate        = MatDuplicate_UMFPACK;
361f0c56d0fSKris Buschelman   B->ops->view             = MatView_UMFPACK;
362f0c56d0fSKris Buschelman   B->ops->assemblyend      = MatAssemblyEnd_UMFPACK;
363f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
364f0c56d0fSKris Buschelman   B->ops->destroy          = MatDestroy_UMFPACK;
365d844382dSKris Buschelman 
366d844382dSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_umfpack_C",
367d844382dSKris Buschelman                                            "MatConvert_SeqAIJ_UMFPACK",MatConvert_SeqAIJ_UMFPACK);CHKERRQ(ierr);
368d844382dSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_umfpack_seqaij_C",
369d844382dSKris Buschelman                                            "MatConvert_UMFPACK_SeqAIJ",MatConvert_UMFPACK_SeqAIJ);CHKERRQ(ierr);
370d844382dSKris Buschelman   PetscLogInfo(0,"Using UMFPACK for SeqAIJ LU factorization and solves.");
371d844382dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATUMFPACK);CHKERRQ(ierr);
372d844382dSKris Buschelman   *newmat = B;
373d844382dSKris Buschelman   PetscFunctionReturn(0);
374d844382dSKris Buschelman }
375d844382dSKris Buschelman EXTERN_C_END
376d844382dSKris Buschelman 
377f0c56d0fSKris Buschelman #undef __FUNCT__
378f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_UMFPACK"
379*6849ba73SBarry Smith PetscErrorCode MatDuplicate_UMFPACK(Mat A, MatDuplicateOption op, Mat *M)
380*6849ba73SBarry Smith {
381dfbe8321SBarry Smith   PetscErrorCode ierr;
3828f340917SKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
3838f340917SKris Buschelman 
384f0c56d0fSKris Buschelman   PetscFunctionBegin;
3858f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
3863f953163SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_UMFPACK));CHKERRQ(ierr);
387f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
388f0c56d0fSKris Buschelman }
389f0c56d0fSKris Buschelman 
3902f71e704SKris Buschelman /*MC
391fafad747SKris Buschelman   MATUMFPACK - MATUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
3922f71e704SKris Buschelman   via the external package UMFPACK.
3932f71e704SKris Buschelman 
3942f71e704SKris Buschelman   If UMFPACK is installed (see the manual for
3952f71e704SKris Buschelman   instructions on how to declare the existence of external packages),
3962f71e704SKris Buschelman   a matrix type can be constructed which invokes UMFPACK solvers.
3972f71e704SKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,UMFPACK).
3982f71e704SKris Buschelman   This matrix type is only supported for double precision real.
3992f71e704SKris Buschelman 
4002f71e704SKris Buschelman   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
40128b08bd3SKris Buschelman   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from
40228b08bd3SKris Buschelman   the MATSEQAIJ type without data copy.
4032f71e704SKris Buschelman 
4042f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
4052f71e704SKris Buschelman   which correspond to the options database keys below.
4062f71e704SKris Buschelman 
4072f71e704SKris Buschelman   Options Database Keys:
4080bad9183SKris Buschelman + -mat_type umfpack - sets the matrix type to "umfpack" during a call to MatSetFromOptions()
4092f71e704SKris Buschelman . -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL]
4102f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
4112f71e704SKris Buschelman . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
4122f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
4132f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
4142f71e704SKris Buschelman - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
4152f71e704SKris Buschelman 
4162f71e704SKris Buschelman    Level: beginner
4172f71e704SKris Buschelman 
4182f71e704SKris Buschelman .seealso: PCLU
4192f71e704SKris Buschelman M*/
4202f71e704SKris Buschelman 
4211677d0b8SKris Buschelman EXTERN_C_BEGIN
4221677d0b8SKris Buschelman #undef __FUNCT__
423f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_UMFPACK"
424dfbe8321SBarry Smith PetscErrorCode MatCreate_UMFPACK(Mat A)
425dfbe8321SBarry Smith {
426dfbe8321SBarry Smith   PetscErrorCode ierr;
4271677d0b8SKris Buschelman 
4281677d0b8SKris Buschelman   PetscFunctionBegin;
4295441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ and UMFPACK types */
4305441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATUMFPACK);CHKERRQ(ierr);
4311677d0b8SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
432d844382dSKris Buschelman   ierr = MatConvert_SeqAIJ_UMFPACK(A,MATUMFPACK,&A);CHKERRQ(ierr);
4331677d0b8SKris Buschelman   PetscFunctionReturn(0);
4341677d0b8SKris Buschelman }
4351677d0b8SKris Buschelman EXTERN_C_END
436