xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 2d4e298286c851afd491cc748c20829c18377a93)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
21677d0b8SKris Buschelman 
31677d0b8SKris Buschelman /*
4*2d4e2982SHong Zhang    Provides an interface to the UMFPACKv5.1 sparse solver
5*2d4e2982SHong Zhang 
6*2d4e2982SHong Zhang    This interface uses the "UF_long version" of the UMFPACK API
7*2d4e2982SHong Zhang    (*_dl_* and *_zl_* routines, instead of *_di_* and *_zi_* routines)
8*2d4e2982SHong Zhang    so that UMFPACK can address more than 2Gb of memory on 64 bit
9*2d4e2982SHong Zhang    machines.
10*2d4e2982SHong Zhang 
11*2d4e2982SHong Zhang    If sizeof(UF_long) == 32 the interface only allocates the memory
12*2d4e2982SHong Zhang    necessary for UMFPACK's working arrays (W, Wi and possibly
13*2d4e2982SHong Zhang    perm_c). If sizeof(UF_long) == 64, in addition to allocating the
14*2d4e2982SHong Zhang    working arrays, the interface also has to re-allocate the matrix
15*2d4e2982SHong Zhang    index arrays (ai and aj, which must be stored as UF_long).
16*2d4e2982SHong Zhang 
17*2d4e2982SHong Zhang    The interface is implemented for both real and complex
18*2d4e2982SHong Zhang    arithmetic. Complex numbers are assumed to be in packed format,
19*2d4e2982SHong Zhang    which requires UMFPACK >= 4.4.
20*2d4e2982SHong Zhang 
21*2d4e2982SHong Zhang    We thank Christophe Geuzaine <geuzaine@acm.caltech.edu> for upgrading this interface to the UMFPACKv5.1
221677d0b8SKris Buschelman */
231677d0b8SKris Buschelman 
241677d0b8SKris Buschelman #include "src/mat/impls/aij/seq/aij.h"
251677d0b8SKris Buschelman 
261677d0b8SKris Buschelman EXTERN_C_BEGIN
271677d0b8SKris Buschelman #include "umfpack.h"
281677d0b8SKris Buschelman EXTERN_C_END
291677d0b8SKris Buschelman 
301677d0b8SKris Buschelman typedef struct {
311677d0b8SKris Buschelman   void         *Symbolic, *Numeric;
321677d0b8SKris Buschelman   double       Info[UMFPACK_INFO], Control[UMFPACK_CONTROL],*W;
33*2d4e2982SHong Zhang   UF_long      *Wi,*ai,*aj,*perm_c;
341677d0b8SKris Buschelman   PetscScalar  *av;
351677d0b8SKris Buschelman   MatStructure flg;
361677d0b8SKris Buschelman   PetscTruth   PetscMatOdering;
371677d0b8SKris Buschelman 
381677d0b8SKris Buschelman   /* A few function pointers for inheritance */
396849ba73SBarry Smith   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
406849ba73SBarry Smith   PetscErrorCode (*MatView)(Mat,PetscViewer);
416849ba73SBarry Smith   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
426849ba73SBarry Smith   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
436849ba73SBarry Smith   PetscErrorCode (*MatDestroy)(Mat);
441677d0b8SKris Buschelman 
451677d0b8SKris Buschelman   /* Flag to clean up UMFPACK objects during Destroy */
461677d0b8SKris Buschelman   PetscTruth CleanUpUMFPACK;
47f0c56d0fSKris Buschelman } Mat_UMFPACK;
48f0c56d0fSKris Buschelman 
49dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_UMFPACK(Mat,MatDuplicateOption,Mat*);
50d844382dSKris Buschelman 
51d844382dSKris Buschelman EXTERN_C_BEGIN
52d844382dSKris Buschelman #undef __FUNCT__
53d844382dSKris Buschelman #define __FUNCT__ "MatConvert_UMFPACK_SeqAIJ"
5475179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_UMFPACK_SeqAIJ(Mat A,MatType type,MatReuse reuse,Mat *newmat)
556849ba73SBarry Smith {
56dfbe8321SBarry Smith   PetscErrorCode ierr;
57d844382dSKris Buschelman   Mat         B=*newmat;
58f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
59d844382dSKris Buschelman 
60d844382dSKris Buschelman   PetscFunctionBegin;
61ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
62d844382dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
63f0c56d0fSKris Buschelman   }
64d844382dSKris Buschelman   /* Reset the original function pointers */
65901853e0SKris Buschelman   B->ops->duplicate        = lu->MatDuplicate;
66901853e0SKris Buschelman   B->ops->view             = lu->MatView;
67901853e0SKris Buschelman   B->ops->assemblyend      = lu->MatAssemblyEnd;
68901853e0SKris Buschelman   B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
69901853e0SKris Buschelman   B->ops->destroy          = lu->MatDestroy;
70901853e0SKris Buschelman 
71901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_umfpack_C","",PETSC_NULL);CHKERRQ(ierr);
72901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_umfpack_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
73901853e0SKris Buschelman 
74d844382dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
75d844382dSKris Buschelman   *newmat = B;
76d844382dSKris Buschelman   PetscFunctionReturn(0);
77d844382dSKris Buschelman }
78d844382dSKris Buschelman EXTERN_C_END
791677d0b8SKris Buschelman 
801677d0b8SKris Buschelman #undef __FUNCT__
81f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_UMFPACK"
82dfbe8321SBarry Smith PetscErrorCode MatDestroy_UMFPACK(Mat A)
83dfbe8321SBarry Smith {
84dfbe8321SBarry Smith   PetscErrorCode ierr;
85f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
861677d0b8SKris Buschelman 
871677d0b8SKris Buschelman   PetscFunctionBegin;
88fb731535SKris Buschelman   if (lu->CleanUpUMFPACK) {
89*2d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
90*2d4e2982SHong Zhang     umfpack_zl_free_symbolic(&lu->Symbolic);
91*2d4e2982SHong Zhang     umfpack_zl_free_numeric(&lu->Numeric);
92*2d4e2982SHong Zhang #else
93*2d4e2982SHong Zhang     umfpack_dl_free_symbolic(&lu->Symbolic);
94*2d4e2982SHong Zhang     umfpack_dl_free_numeric(&lu->Numeric);
95*2d4e2982SHong Zhang #endif
961677d0b8SKris Buschelman     ierr = PetscFree(lu->Wi);CHKERRQ(ierr);
971677d0b8SKris Buschelman     ierr = PetscFree(lu->W);CHKERRQ(ierr);
98*2d4e2982SHong Zhang     if(sizeof(UF_long) != sizeof(int)){
99*2d4e2982SHong Zhang       ierr = PetscFree(lu->ai);CHKERRQ(ierr);
100*2d4e2982SHong Zhang       ierr = PetscFree(lu->aj);CHKERRQ(ierr);
101*2d4e2982SHong Zhang     }
1021677d0b8SKris Buschelman     if (lu->PetscMatOdering) {
1031677d0b8SKris Buschelman       ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
1041677d0b8SKris Buschelman     }
1051677d0b8SKris Buschelman   }
106ceb03754SKris Buschelman   ierr = MatConvert_UMFPACK_SeqAIJ(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
107d844382dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
1081677d0b8SKris Buschelman   PetscFunctionReturn(0);
1091677d0b8SKris Buschelman }
1101677d0b8SKris Buschelman 
1111677d0b8SKris Buschelman #undef __FUNCT__
112f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_UMFPACK"
1136849ba73SBarry Smith PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x)
1146849ba73SBarry Smith {
115f0c56d0fSKris Buschelman   Mat_UMFPACK *lu = (Mat_UMFPACK*)A->spptr;
1161677d0b8SKris Buschelman   PetscScalar *av=lu->av,*ba,*xa;
117dfbe8321SBarry Smith   PetscErrorCode ierr;
118*2d4e2982SHong Zhang   UF_long     *ai=lu->ai,*aj=lu->aj,status;
1191677d0b8SKris Buschelman 
1201677d0b8SKris Buschelman   PetscFunctionBegin;
121*2d4e2982SHong Zhang   /* solve Ax = b by umfpack_*_wsolve */
1221677d0b8SKris Buschelman   /* ----------------------------------*/
123*2d4e2982SHong Zhang 
124*2d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
125*2d4e2982SHong Zhang   ierr = VecConjugate(b);
126*2d4e2982SHong Zhang #endif
127*2d4e2982SHong Zhang 
1281677d0b8SKris Buschelman   ierr = VecGetArray(b,&ba);
1291677d0b8SKris Buschelman   ierr = VecGetArray(x,&xa);
1301677d0b8SKris Buschelman 
131*2d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
132*2d4e2982SHong Zhang   status = umfpack_zl_wsolve(UMFPACK_At,ai,aj,(double*)av,NULL,(double*)xa,NULL,(double*)ba,NULL,
133*2d4e2982SHong Zhang 			     lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
134*2d4e2982SHong Zhang   umfpack_zl_report_info(lu->Control, lu->Info);
1351677d0b8SKris Buschelman   if (status < 0){
136*2d4e2982SHong Zhang     umfpack_zl_report_status(lu->Control, status);
137*2d4e2982SHong Zhang     SETERRQ(PETSC_ERR_LIB,"umfpack_zl_wsolve failed");
1381677d0b8SKris Buschelman   }
139*2d4e2982SHong Zhang #else
140*2d4e2982SHong Zhang   status = umfpack_dl_wsolve(UMFPACK_At,ai,aj,av,xa,ba,
141*2d4e2982SHong Zhang 			     lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
142*2d4e2982SHong Zhang   umfpack_dl_report_info(lu->Control, lu->Info);
143*2d4e2982SHong Zhang   if (status < 0){
144*2d4e2982SHong Zhang     umfpack_dl_report_status(lu->Control, status);
145*2d4e2982SHong Zhang     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_wsolve failed");
146*2d4e2982SHong Zhang   }
147*2d4e2982SHong Zhang #endif
1481677d0b8SKris Buschelman 
1491677d0b8SKris Buschelman   ierr = VecRestoreArray(b,&ba);
1501677d0b8SKris Buschelman   ierr = VecRestoreArray(x,&xa);
1512a325a84SHong Zhang 
152*2d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
153*2d4e2982SHong Zhang   ierr = VecConjugate(b);
154*2d4e2982SHong Zhang #endif
155*2d4e2982SHong Zhang 
1561677d0b8SKris Buschelman   PetscFunctionReturn(0);
1571677d0b8SKris Buschelman }
1581677d0b8SKris Buschelman 
1591677d0b8SKris Buschelman #undef __FUNCT__
160f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_UMFPACK"
161af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat A,MatFactorInfo *info,Mat *F)
1626849ba73SBarry Smith {
163f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)(*F)->spptr;
1646849ba73SBarry Smith   PetscErrorCode ierr;
165*2d4e2982SHong Zhang   UF_long     *ai=lu->ai,*aj=lu->aj,m=A->rmap.n,status;
1661677d0b8SKris Buschelman   PetscScalar *av=lu->av;
1671677d0b8SKris Buschelman 
1681677d0b8SKris Buschelman   PetscFunctionBegin;
1691677d0b8SKris Buschelman   /* numeric factorization of A' */
1701677d0b8SKris Buschelman   /* ----------------------------*/
171*2d4e2982SHong Zhang 
172*2d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1732a325a84SHong Zhang   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
174*2d4e2982SHong Zhang     umfpack_zl_free_numeric(&lu->Numeric);
1752a325a84SHong Zhang   }
176*2d4e2982SHong Zhang   status = umfpack_zl_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
177*2d4e2982SHong Zhang   if (status < 0) SETERRQ(PETSC_ERR_LIB,"umfpack_zl_numeric failed");
1781677d0b8SKris Buschelman   /* report numeric factorization of A' when Control[PRL] > 3 */
179*2d4e2982SHong Zhang   (void) umfpack_zl_report_numeric(lu->Numeric, lu->Control);
180*2d4e2982SHong Zhang #else
181*2d4e2982SHong Zhang   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
182*2d4e2982SHong Zhang     umfpack_dl_free_numeric(&lu->Numeric);
183*2d4e2982SHong Zhang   }
184*2d4e2982SHong Zhang   status = umfpack_dl_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
185*2d4e2982SHong Zhang   if (status < 0) SETERRQ(PETSC_ERR_LIB,"umfpack_dl_numeric failed");
186*2d4e2982SHong Zhang   /* report numeric factorization of A' when Control[PRL] > 3 */
187*2d4e2982SHong Zhang   (void) umfpack_dl_report_numeric(lu->Numeric, lu->Control);
188*2d4e2982SHong Zhang #endif
1891677d0b8SKris Buschelman 
1901677d0b8SKris Buschelman   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
1911677d0b8SKris Buschelman     /* allocate working space to be used by Solve */
192*2d4e2982SHong Zhang     ierr = PetscMalloc(m * sizeof(UF_long), &lu->Wi);CHKERRQ(ierr);
193*2d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
194*2d4e2982SHong Zhang     ierr = PetscMalloc(2*5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
195*2d4e2982SHong Zhang #else
1961677d0b8SKris Buschelman     ierr = PetscMalloc(5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
197*2d4e2982SHong Zhang #endif
1981677d0b8SKris Buschelman   }
199*2d4e2982SHong Zhang 
2002a325a84SHong Zhang   lu->flg = SAME_NONZERO_PATTERN;
2012a325a84SHong Zhang   lu->CleanUpUMFPACK = PETSC_TRUE;
2021677d0b8SKris Buschelman   PetscFunctionReturn(0);
2031677d0b8SKris Buschelman }
2041677d0b8SKris Buschelman 
2051677d0b8SKris Buschelman /*
2061677d0b8SKris Buschelman    Note the r permutation is ignored
2071677d0b8SKris Buschelman */
2081677d0b8SKris Buschelman #undef __FUNCT__
209f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK"
2106849ba73SBarry Smith PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
2116849ba73SBarry Smith {
212e1b4c3a1SKris Buschelman   Mat         B;
2131677d0b8SKris Buschelman   Mat_SeqAIJ  *mat=(Mat_SeqAIJ*)A->data;
214f0c56d0fSKris Buschelman   Mat_UMFPACK *lu;
215dfbe8321SBarry Smith   PetscErrorCode ierr;
216*2d4e2982SHong Zhang   int         i,m=A->rmap.n,n=A->cmap.n,*ra,idx;
217*2d4e2982SHong Zhang   UF_long     status;
218*2d4e2982SHong Zhang 
2191677d0b8SKris Buschelman   PetscScalar *av=mat->a;
2200f6efc10SHong Zhang   const char  *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC","2BY2"},
2210f6efc10SHong Zhang               *scale[]={"NONE","SUM","MAX"};
2220f6efc10SHong Zhang   PetscTruth  flg;
2231677d0b8SKris Buschelman 
2241677d0b8SKris Buschelman   PetscFunctionBegin;
2251677d0b8SKris Buschelman   /* Create the factorization matrix F */
226f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
227f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
228be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
229e1b4c3a1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
2301677d0b8SKris Buschelman 
231f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
232f0c56d0fSKris Buschelman   B->ops->solve           = MatSolve_UMFPACK;
233e1b4c3a1SKris Buschelman   B->factor               = FACTOR_LU;
23494b7f48cSBarry Smith   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */
2351677d0b8SKris Buschelman 
236f0c56d0fSKris Buschelman   lu = (Mat_UMFPACK*)(B->spptr);
2371677d0b8SKris Buschelman 
2381677d0b8SKris Buschelman   /* initializations */
2391677d0b8SKris Buschelman   /* ------------------------------------------------*/
2401677d0b8SKris Buschelman   /* get the default control parameters */
241*2d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
242*2d4e2982SHong Zhang   umfpack_zl_defaults(lu->Control);
243*2d4e2982SHong Zhang #else
244*2d4e2982SHong Zhang   umfpack_dl_defaults(lu->Control);
245*2d4e2982SHong Zhang #endif
2461677d0b8SKris Buschelman   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
2470f6efc10SHong Zhang   lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */
2481677d0b8SKris Buschelman 
2491677d0b8SKris Buschelman   ierr = PetscOptionsBegin(A->comm,A->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
2501677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
2511677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr);
2521677d0b8SKris Buschelman 
2531677d0b8SKris Buschelman   /* Control parameters for symbolic factorization */
2540f6efc10SHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,4,strategy[0],&idx,&flg);CHKERRQ(ierr);
2550f6efc10SHong Zhang   if (flg) {
2560f6efc10SHong Zhang     switch (idx){
2570f6efc10SHong Zhang     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
2580f6efc10SHong Zhang     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
2590f6efc10SHong Zhang     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
2600f6efc10SHong Zhang     case 3: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_2BY2; break;
2610f6efc10SHong Zhang     }
2620f6efc10SHong Zhang   }
2631677d0b8SKris 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);
2641677d0b8SKris 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);
2650f6efc10SHong 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);
2661677d0b8SKris 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);
2670f6efc10SHong 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);
2680f6efc10SHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);CHKERRQ(ierr);
2690f6efc10SHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);CHKERRQ(ierr);
2701677d0b8SKris Buschelman 
2711677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
2721677d0b8SKris 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);
2730f6efc10SHong 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);
2740f6efc10SHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
2750f6efc10SHong Zhang   if (flg) {
2760f6efc10SHong Zhang     switch (idx){
2770f6efc10SHong Zhang     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
2780f6efc10SHong Zhang     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
2790f6efc10SHong Zhang     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
2800f6efc10SHong Zhang     }
2810f6efc10SHong Zhang   }
2821677d0b8SKris 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);
2830f6efc10SHong 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);
2840f6efc10SHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);CHKERRQ(ierr);
2851677d0b8SKris Buschelman 
2861677d0b8SKris Buschelman   /* Control parameters used by solve */
2871677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr);
2881677d0b8SKris Buschelman 
2890f6efc10SHong Zhang   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */
2902401956bSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-pc_factor_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr);
2911677d0b8SKris Buschelman   if (lu->PetscMatOdering) {
2920f6efc10SHong Zhang     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
293*2d4e2982SHong Zhang     ierr = PetscMalloc(m*sizeof(UF_long),&lu->perm_c);CHKERRQ(ierr);
294*2d4e2982SHong Zhang     /* we cannot simply memcpy on 64 bit archs */
295*2d4e2982SHong Zhang     for(i = 0; i < m; i++) lu->perm_c[i] = ra[i];
2960f6efc10SHong Zhang     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
2971677d0b8SKris Buschelman   }
2981677d0b8SKris Buschelman   PetscOptionsEnd();
2991677d0b8SKris Buschelman 
300*2d4e2982SHong Zhang   if(sizeof(UF_long) != sizeof(int)){
301*2d4e2982SHong Zhang     /* we cannot directly use mat->i and mat->j on 64 bit archs */
302*2d4e2982SHong Zhang     ierr = PetscMalloc((m+1)*sizeof(UF_long),&lu->ai);CHKERRQ(ierr);
303*2d4e2982SHong Zhang     ierr = PetscMalloc(mat->nz*sizeof(UF_long),&lu->aj);CHKERRQ(ierr);
304*2d4e2982SHong Zhang     for(i = 0; i < m + 1; i++) lu->ai[i] = mat->i[i];
305*2d4e2982SHong Zhang     for(i = 0; i < mat->nz; i++) lu->aj[i] = mat->j[i];
306*2d4e2982SHong Zhang   }
307*2d4e2982SHong Zhang   else{
308*2d4e2982SHong Zhang     lu->ai = (UF_long*)mat->i;
309*2d4e2982SHong Zhang     lu->aj = (UF_long*)mat->j;
310*2d4e2982SHong Zhang   }
311*2d4e2982SHong Zhang 
3121677d0b8SKris Buschelman   /* print the control parameters */
313*2d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
314*2d4e2982SHong Zhang   if(lu->Control[UMFPACK_PRL] > 1) umfpack_zl_report_control(lu->Control);
315*2d4e2982SHong Zhang #else
316*2d4e2982SHong Zhang   if(lu->Control[UMFPACK_PRL] > 1) umfpack_dl_report_control(lu->Control);
317*2d4e2982SHong Zhang #endif
3181677d0b8SKris Buschelman 
3191677d0b8SKris Buschelman   /* symbolic factorization of A' */
3201677d0b8SKris Buschelman   /* ---------------------------------------------------------------------- */
321*2d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
3220f6efc10SHong Zhang   if (lu->PetscMatOdering) { /* use Petsc row ordering */
323*2d4e2982SHong Zhang     status = umfpack_zl_qsymbolic(n,m,lu->ai,lu->aj,(double*)av,NULL,
324*2d4e2982SHong Zhang 				  lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
3250f6efc10SHong Zhang   } else { /* use Umfpack col ordering */
326*2d4e2982SHong Zhang     status = umfpack_zl_symbolic(n,m,lu->ai,lu->aj,(double*)av,NULL,
327*2d4e2982SHong Zhang 				 &lu->Symbolic,lu->Control,lu->Info);
3280f6efc10SHong Zhang   }
3291677d0b8SKris Buschelman   if (status < 0){
330*2d4e2982SHong Zhang     umfpack_zl_report_info(lu->Control, lu->Info);
331*2d4e2982SHong Zhang     umfpack_zl_report_status(lu->Control, status);
332*2d4e2982SHong Zhang     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_symbolic failed");
3331677d0b8SKris Buschelman   }
3341677d0b8SKris Buschelman   /* report sumbolic factorization of A' when Control[PRL] > 3 */
335*2d4e2982SHong Zhang   (void) umfpack_zl_report_symbolic(lu->Symbolic, lu->Control);
336*2d4e2982SHong Zhang #else
337*2d4e2982SHong Zhang   if (lu->PetscMatOdering) { /* use Petsc row ordering */
338*2d4e2982SHong Zhang     status = umfpack_dl_qsymbolic(n,m,lu->ai,lu->aj,av,
339*2d4e2982SHong Zhang 				  lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
340*2d4e2982SHong Zhang   } else { /* use Umfpack col ordering */
341*2d4e2982SHong Zhang     status = umfpack_dl_symbolic(n,m,lu->ai,lu->aj,av,
342*2d4e2982SHong Zhang 				 &lu->Symbolic,lu->Control,lu->Info);
343*2d4e2982SHong Zhang   }
344*2d4e2982SHong Zhang   if (status < 0){
345*2d4e2982SHong Zhang     umfpack_dl_report_info(lu->Control, lu->Info);
346*2d4e2982SHong Zhang     umfpack_dl_report_status(lu->Control, status);
347*2d4e2982SHong Zhang     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_symbolic failed");
348*2d4e2982SHong Zhang   }
349*2d4e2982SHong Zhang   /* report sumbolic factorization of A' when Control[PRL] > 3 */
350*2d4e2982SHong Zhang   (void) umfpack_dl_report_symbolic(lu->Symbolic, lu->Control);
351*2d4e2982SHong Zhang #endif
3521677d0b8SKris Buschelman 
3531677d0b8SKris Buschelman   lu->flg = DIFFERENT_NONZERO_PATTERN;
3541677d0b8SKris Buschelman   lu->av  = av;
355e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK = PETSC_TRUE;
356e1b4c3a1SKris Buschelman   *F = B;
3571677d0b8SKris Buschelman   PetscFunctionReturn(0);
3581677d0b8SKris Buschelman }
3591677d0b8SKris Buschelman 
3601677d0b8SKris Buschelman #undef __FUNCT__
361f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_UMFPACK"
3626849ba73SBarry Smith PetscErrorCode MatAssemblyEnd_UMFPACK(Mat A,MatAssemblyType mode)
3636849ba73SBarry Smith {
364dfbe8321SBarry Smith   PetscErrorCode ierr;
365f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)(A->spptr);
366d844382dSKris Buschelman 
3671677d0b8SKris Buschelman   PetscFunctionBegin;
368d844382dSKris Buschelman   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
369d844382dSKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
370f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
3711677d0b8SKris Buschelman   PetscFunctionReturn(0);
3721677d0b8SKris Buschelman }
3731677d0b8SKris Buschelman 
37494b7f48cSBarry Smith /* used by -ksp_view */
3751677d0b8SKris Buschelman #undef __FUNCT__
376f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_UMFPACK"
3776849ba73SBarry Smith PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer)
3786849ba73SBarry Smith {
379f0c56d0fSKris Buschelman   Mat_UMFPACK *lu= (Mat_UMFPACK*)A->spptr;
380dfbe8321SBarry Smith   PetscErrorCode ierr;
381f0c56d0fSKris Buschelman 
3821677d0b8SKris Buschelman   PetscFunctionBegin;
3831677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
384f0c56d0fSKris Buschelman   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
3851677d0b8SKris Buschelman 
3861677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
3871677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
3881677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
3891677d0b8SKris Buschelman 
3901677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
3910f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr);
3921677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
3931677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
3940f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr);
3951677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
3960f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_2BY2_TOLERANCE]: %g\n",lu->Control[UMFPACK_2BY2_TOLERANCE]);CHKERRQ(ierr);
3970f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr);
3980f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr);
3991677d0b8SKris Buschelman 
4001677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
4011677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
4020f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr);
4030f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr);
4041677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
4050f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr);
4061677d0b8SKris Buschelman 
4071677d0b8SKris Buschelman   /* Control parameters used by solve */
4081677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
4091677d0b8SKris Buschelman 
4101677d0b8SKris Buschelman   /* mat ordering */
4111677d0b8SKris Buschelman   if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer,"  UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr);
4121677d0b8SKris Buschelman 
4131677d0b8SKris Buschelman   PetscFunctionReturn(0);
4141677d0b8SKris Buschelman }
4151677d0b8SKris Buschelman 
416d844382dSKris Buschelman #undef __FUNCT__
417f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_UMFPACK"
4186849ba73SBarry Smith PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
4196849ba73SBarry Smith {
420dfbe8321SBarry Smith   PetscErrorCode ierr;
42132077d6dSBarry Smith   PetscTruth        iascii;
422d844382dSKris Buschelman   PetscViewerFormat format;
423f0c56d0fSKris Buschelman   Mat_UMFPACK       *lu=(Mat_UMFPACK*)(A->spptr);
424d844382dSKris Buschelman 
425d844382dSKris Buschelman   PetscFunctionBegin;
426d844382dSKris Buschelman   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
427d844382dSKris Buschelman 
42832077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
42932077d6dSBarry Smith   if (iascii) {
430d844382dSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
4312f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO) {
432f0c56d0fSKris Buschelman       ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
433d844382dSKris Buschelman     }
434d844382dSKris Buschelman   }
435d844382dSKris Buschelman   PetscFunctionReturn(0);
436d844382dSKris Buschelman }
437d844382dSKris Buschelman 
438d844382dSKris Buschelman EXTERN_C_BEGIN
439d844382dSKris Buschelman #undef __FUNCT__
440d844382dSKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_UMFPACK"
44175179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_UMFPACK(Mat A,MatType type,MatReuse reuse,Mat *newmat)
4426849ba73SBarry Smith {
443dfbe8321SBarry Smith   PetscErrorCode ierr;
444d844382dSKris Buschelman   Mat            B=*newmat;
445f0c56d0fSKris Buschelman   Mat_UMFPACK    *lu;
446d844382dSKris Buschelman 
447d844382dSKris Buschelman   PetscFunctionBegin;
448ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
449d844382dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
450d844382dSKris Buschelman   }
451d844382dSKris Buschelman 
452f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_UMFPACK,&lu);CHKERRQ(ierr);
453f0c56d0fSKris Buschelman   lu->MatDuplicate         = A->ops->duplicate;
454d844382dSKris Buschelman   lu->MatView              = A->ops->view;
455d844382dSKris Buschelman   lu->MatAssemblyEnd       = A->ops->assemblyend;
456d844382dSKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
457d844382dSKris Buschelman   lu->MatDestroy           = A->ops->destroy;
458d844382dSKris Buschelman   lu->CleanUpUMFPACK       = PETSC_FALSE;
459d844382dSKris Buschelman 
460d844382dSKris Buschelman   B->spptr                 = (void*)lu;
461f0c56d0fSKris Buschelman   B->ops->duplicate        = MatDuplicate_UMFPACK;
462f0c56d0fSKris Buschelman   B->ops->view             = MatView_UMFPACK;
463f0c56d0fSKris Buschelman   B->ops->assemblyend      = MatAssemblyEnd_UMFPACK;
464f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
465f0c56d0fSKris Buschelman   B->ops->destroy          = MatDestroy_UMFPACK;
466d844382dSKris Buschelman 
467d844382dSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_umfpack_C",
468d844382dSKris Buschelman                                            "MatConvert_SeqAIJ_UMFPACK",MatConvert_SeqAIJ_UMFPACK);CHKERRQ(ierr);
469d844382dSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_umfpack_seqaij_C",
470d844382dSKris Buschelman                                            "MatConvert_UMFPACK_SeqAIJ",MatConvert_UMFPACK_SeqAIJ);CHKERRQ(ierr);
471ae15b995SBarry Smith   ierr = PetscInfo(0,"Using UMFPACK for SeqAIJ LU factorization and solves.\n");CHKERRQ(ierr);
472d844382dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATUMFPACK);CHKERRQ(ierr);
473d844382dSKris Buschelman   *newmat = B;
474d844382dSKris Buschelman   PetscFunctionReturn(0);
475d844382dSKris Buschelman }
476d844382dSKris Buschelman EXTERN_C_END
477d844382dSKris Buschelman 
478f0c56d0fSKris Buschelman #undef __FUNCT__
479f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_UMFPACK"
4806849ba73SBarry Smith PetscErrorCode MatDuplicate_UMFPACK(Mat A, MatDuplicateOption op, Mat *M)
4816849ba73SBarry Smith {
482dfbe8321SBarry Smith   PetscErrorCode ierr;
4838f340917SKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
4848f340917SKris Buschelman 
485f0c56d0fSKris Buschelman   PetscFunctionBegin;
4868f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
4873f953163SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_UMFPACK));CHKERRQ(ierr);
488f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
489f0c56d0fSKris Buschelman }
490f0c56d0fSKris Buschelman 
4912f71e704SKris Buschelman /*MC
492fafad747SKris Buschelman   MATUMFPACK - MATUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
4932f71e704SKris Buschelman   via the external package UMFPACK.
4942f71e704SKris Buschelman 
4952f71e704SKris Buschelman   If UMFPACK is installed (see the manual for
4962f71e704SKris Buschelman   instructions on how to declare the existence of external packages),
4972f71e704SKris Buschelman   a matrix type can be constructed which invokes UMFPACK solvers.
4982f71e704SKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,UMFPACK).
4992f71e704SKris Buschelman 
5002f71e704SKris Buschelman   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
50128b08bd3SKris Buschelman   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from
50228b08bd3SKris Buschelman   the MATSEQAIJ type without data copy.
5032f71e704SKris Buschelman 
5042f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
5052f71e704SKris Buschelman   which correspond to the options database keys below.
5062f71e704SKris Buschelman 
5072f71e704SKris Buschelman   Options Database Keys:
5080bad9183SKris Buschelman + -mat_type umfpack - sets the matrix type to "umfpack" during a call to MatSetFromOptions()
5092f71e704SKris Buschelman . -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL]
5102f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
5112f71e704SKris Buschelman . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
5122f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
5132f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
5142f71e704SKris Buschelman - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
5152f71e704SKris Buschelman 
5162f71e704SKris Buschelman    Level: beginner
5172f71e704SKris Buschelman 
5182f71e704SKris Buschelman .seealso: PCLU
5192f71e704SKris Buschelman M*/
5202f71e704SKris Buschelman 
5211677d0b8SKris Buschelman EXTERN_C_BEGIN
5221677d0b8SKris Buschelman #undef __FUNCT__
523f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_UMFPACK"
524be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_UMFPACK(Mat A)
525dfbe8321SBarry Smith {
526dfbe8321SBarry Smith   PetscErrorCode ierr;
5271677d0b8SKris Buschelman 
5281677d0b8SKris Buschelman   PetscFunctionBegin;
5291677d0b8SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
530ceb03754SKris Buschelman   ierr = MatConvert_SeqAIJ_UMFPACK(A,MATUMFPACK,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
5311677d0b8SKris Buschelman   PetscFunctionReturn(0);
5321677d0b8SKris Buschelman }
5331677d0b8SKris Buschelman EXTERN_C_END
534*2d4e2982SHong Zhang 
535*2d4e2982SHong Zhang 
536*2d4e2982SHong Zhang 
537