xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 5b2a9f606ee6c957b065fbf8ce7f96f1b2adcab9)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
21677d0b8SKris Buschelman 
31677d0b8SKris Buschelman /*
42d4e2982SHong Zhang    Provides an interface to the UMFPACKv5.1 sparse solver
52d4e2982SHong Zhang 
62d4e2982SHong Zhang    This interface uses the "UF_long version" of the UMFPACK API
72d4e2982SHong Zhang    (*_dl_* and *_zl_* routines, instead of *_di_* and *_zi_* routines)
82d4e2982SHong Zhang    so that UMFPACK can address more than 2Gb of memory on 64 bit
92d4e2982SHong Zhang    machines.
102d4e2982SHong Zhang 
112d4e2982SHong Zhang    If sizeof(UF_long) == 32 the interface only allocates the memory
122d4e2982SHong Zhang    necessary for UMFPACK's working arrays (W, Wi and possibly
132d4e2982SHong Zhang    perm_c). If sizeof(UF_long) == 64, in addition to allocating the
142d4e2982SHong Zhang    working arrays, the interface also has to re-allocate the matrix
152d4e2982SHong Zhang    index arrays (ai and aj, which must be stored as UF_long).
162d4e2982SHong Zhang 
172d4e2982SHong Zhang    The interface is implemented for both real and complex
182d4e2982SHong Zhang    arithmetic. Complex numbers are assumed to be in packed format,
192d4e2982SHong Zhang    which requires UMFPACK >= 4.4.
202d4e2982SHong Zhang 
212d4e2982SHong 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;
332d4e2982SHong 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) {
892d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
902d4e2982SHong Zhang     umfpack_zl_free_symbolic(&lu->Symbolic);
912d4e2982SHong Zhang     umfpack_zl_free_numeric(&lu->Numeric);
922d4e2982SHong Zhang #else
932d4e2982SHong Zhang     umfpack_dl_free_symbolic(&lu->Symbolic);
942d4e2982SHong Zhang     umfpack_dl_free_numeric(&lu->Numeric);
952d4e2982SHong Zhang #endif
961677d0b8SKris Buschelman     ierr = PetscFree(lu->Wi);CHKERRQ(ierr);
971677d0b8SKris Buschelman     ierr = PetscFree(lu->W);CHKERRQ(ierr);
982d4e2982SHong Zhang     if(sizeof(UF_long) != sizeof(int)){
992d4e2982SHong Zhang       ierr = PetscFree(lu->ai);CHKERRQ(ierr);
1002d4e2982SHong Zhang       ierr = PetscFree(lu->aj);CHKERRQ(ierr);
1012d4e2982SHong 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;
1182d4e2982SHong Zhang   UF_long     *ai=lu->ai,*aj=lu->aj,status;
1191677d0b8SKris Buschelman 
1201677d0b8SKris Buschelman   PetscFunctionBegin;
1212d4e2982SHong Zhang   /* solve Ax = b by umfpack_*_wsolve */
1221677d0b8SKris Buschelman   /* ----------------------------------*/
1232d4e2982SHong Zhang 
1242d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1252d4e2982SHong Zhang   ierr = VecConjugate(b);
1262d4e2982SHong Zhang #endif
1272d4e2982SHong Zhang 
1281677d0b8SKris Buschelman   ierr = VecGetArray(b,&ba);
1291677d0b8SKris Buschelman   ierr = VecGetArray(x,&xa);
1301677d0b8SKris Buschelman 
1312d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1322d4e2982SHong Zhang   status = umfpack_zl_wsolve(UMFPACK_At,ai,aj,(double*)av,NULL,(double*)xa,NULL,(double*)ba,NULL,
1332d4e2982SHong Zhang 			     lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
1342d4e2982SHong Zhang   umfpack_zl_report_info(lu->Control, lu->Info);
1351677d0b8SKris Buschelman   if (status < 0){
1362d4e2982SHong Zhang     umfpack_zl_report_status(lu->Control, status);
1372d4e2982SHong Zhang     SETERRQ(PETSC_ERR_LIB,"umfpack_zl_wsolve failed");
1381677d0b8SKris Buschelman   }
1392d4e2982SHong Zhang #else
1402d4e2982SHong Zhang   status = umfpack_dl_wsolve(UMFPACK_At,ai,aj,av,xa,ba,
1412d4e2982SHong Zhang 			     lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
1422d4e2982SHong Zhang   umfpack_dl_report_info(lu->Control, lu->Info);
1432d4e2982SHong Zhang   if (status < 0){
1442d4e2982SHong Zhang     umfpack_dl_report_status(lu->Control, status);
1452d4e2982SHong Zhang     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_wsolve failed");
1462d4e2982SHong Zhang   }
1472d4e2982SHong Zhang #endif
1481677d0b8SKris Buschelman 
1491677d0b8SKris Buschelman   ierr = VecRestoreArray(b,&ba);
1501677d0b8SKris Buschelman   ierr = VecRestoreArray(x,&xa);
1512a325a84SHong Zhang 
1522d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1532d4e2982SHong Zhang   ierr = VecConjugate(b);
154*5b2a9f60SHong Zhang   ierr = VecConjugate(x);
1552d4e2982SHong Zhang #endif
1562d4e2982SHong Zhang 
1571677d0b8SKris Buschelman   PetscFunctionReturn(0);
1581677d0b8SKris Buschelman }
1591677d0b8SKris Buschelman 
1601677d0b8SKris Buschelman #undef __FUNCT__
161f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_UMFPACK"
162af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat A,MatFactorInfo *info,Mat *F)
1636849ba73SBarry Smith {
164f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)(*F)->spptr;
1656849ba73SBarry Smith   PetscErrorCode ierr;
1662d4e2982SHong Zhang   UF_long     *ai=lu->ai,*aj=lu->aj,m=A->rmap.n,status;
1671677d0b8SKris Buschelman   PetscScalar *av=lu->av;
1681677d0b8SKris Buschelman 
1691677d0b8SKris Buschelman   PetscFunctionBegin;
1701677d0b8SKris Buschelman   /* numeric factorization of A' */
1711677d0b8SKris Buschelman   /* ----------------------------*/
1722d4e2982SHong Zhang 
1732d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1742a325a84SHong Zhang   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
1752d4e2982SHong Zhang     umfpack_zl_free_numeric(&lu->Numeric);
1762a325a84SHong Zhang   }
1772d4e2982SHong Zhang   status = umfpack_zl_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1782d4e2982SHong Zhang   if (status < 0) SETERRQ(PETSC_ERR_LIB,"umfpack_zl_numeric failed");
1791677d0b8SKris Buschelman   /* report numeric factorization of A' when Control[PRL] > 3 */
1802d4e2982SHong Zhang   (void) umfpack_zl_report_numeric(lu->Numeric, lu->Control);
1812d4e2982SHong Zhang #else
1822d4e2982SHong Zhang   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
1832d4e2982SHong Zhang     umfpack_dl_free_numeric(&lu->Numeric);
1842d4e2982SHong Zhang   }
1852d4e2982SHong Zhang   status = umfpack_dl_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1862d4e2982SHong Zhang   if (status < 0) SETERRQ(PETSC_ERR_LIB,"umfpack_dl_numeric failed");
1872d4e2982SHong Zhang   /* report numeric factorization of A' when Control[PRL] > 3 */
1882d4e2982SHong Zhang   (void) umfpack_dl_report_numeric(lu->Numeric, lu->Control);
1892d4e2982SHong Zhang #endif
1901677d0b8SKris Buschelman 
1911677d0b8SKris Buschelman   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
1921677d0b8SKris Buschelman     /* allocate working space to be used by Solve */
1932d4e2982SHong Zhang     ierr = PetscMalloc(m * sizeof(UF_long), &lu->Wi);CHKERRQ(ierr);
1942d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1952d4e2982SHong Zhang     ierr = PetscMalloc(2*5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
1962d4e2982SHong Zhang #else
1971677d0b8SKris Buschelman     ierr = PetscMalloc(5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
1982d4e2982SHong Zhang #endif
1991677d0b8SKris Buschelman   }
2002d4e2982SHong Zhang 
2012a325a84SHong Zhang   lu->flg = SAME_NONZERO_PATTERN;
2022a325a84SHong Zhang   lu->CleanUpUMFPACK = PETSC_TRUE;
2031677d0b8SKris Buschelman   PetscFunctionReturn(0);
2041677d0b8SKris Buschelman }
2051677d0b8SKris Buschelman 
2061677d0b8SKris Buschelman /*
2071677d0b8SKris Buschelman    Note the r permutation is ignored
2081677d0b8SKris Buschelman */
2091677d0b8SKris Buschelman #undef __FUNCT__
210f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK"
2116849ba73SBarry Smith PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
2126849ba73SBarry Smith {
213e1b4c3a1SKris Buschelman   Mat         B;
2141677d0b8SKris Buschelman   Mat_SeqAIJ  *mat=(Mat_SeqAIJ*)A->data;
215f0c56d0fSKris Buschelman   Mat_UMFPACK *lu;
216dfbe8321SBarry Smith   PetscErrorCode ierr;
2172d4e2982SHong Zhang   int         i,m=A->rmap.n,n=A->cmap.n,*ra,idx;
2182d4e2982SHong Zhang   UF_long     status;
2192d4e2982SHong Zhang 
2201677d0b8SKris Buschelman   PetscScalar *av=mat->a;
2210f6efc10SHong Zhang   const char  *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC","2BY2"},
2220f6efc10SHong Zhang               *scale[]={"NONE","SUM","MAX"};
2230f6efc10SHong Zhang   PetscTruth  flg;
2241677d0b8SKris Buschelman 
2251677d0b8SKris Buschelman   PetscFunctionBegin;
2261677d0b8SKris Buschelman   /* Create the factorization matrix F */
227f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
228f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
229be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
230e1b4c3a1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
2311677d0b8SKris Buschelman 
232f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
233f0c56d0fSKris Buschelman   B->ops->solve           = MatSolve_UMFPACK;
234e1b4c3a1SKris Buschelman   B->factor               = FACTOR_LU;
23594b7f48cSBarry Smith   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */
2361677d0b8SKris Buschelman 
237f0c56d0fSKris Buschelman   lu = (Mat_UMFPACK*)(B->spptr);
2381677d0b8SKris Buschelman 
2391677d0b8SKris Buschelman   /* initializations */
2401677d0b8SKris Buschelman   /* ------------------------------------------------*/
2411677d0b8SKris Buschelman   /* get the default control parameters */
2422d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
2432d4e2982SHong Zhang   umfpack_zl_defaults(lu->Control);
2442d4e2982SHong Zhang #else
2452d4e2982SHong Zhang   umfpack_dl_defaults(lu->Control);
2462d4e2982SHong Zhang #endif
2471677d0b8SKris Buschelman   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
2480f6efc10SHong Zhang   lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */
2491677d0b8SKris Buschelman 
2501677d0b8SKris Buschelman   ierr = PetscOptionsBegin(A->comm,A->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
2511677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
2521677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr);
2531677d0b8SKris Buschelman 
2541677d0b8SKris Buschelman   /* Control parameters for symbolic factorization */
2550f6efc10SHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,4,strategy[0],&idx,&flg);CHKERRQ(ierr);
2560f6efc10SHong Zhang   if (flg) {
2570f6efc10SHong Zhang     switch (idx){
2580f6efc10SHong Zhang     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
2590f6efc10SHong Zhang     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
2600f6efc10SHong Zhang     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
2610f6efc10SHong Zhang     case 3: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_2BY2; break;
2620f6efc10SHong Zhang     }
2630f6efc10SHong Zhang   }
2641677d0b8SKris 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);
2651677d0b8SKris 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);
2660f6efc10SHong 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);
2671677d0b8SKris 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);
2680f6efc10SHong 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);
2690f6efc10SHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);CHKERRQ(ierr);
2700f6efc10SHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);CHKERRQ(ierr);
2711677d0b8SKris Buschelman 
2721677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
2731677d0b8SKris 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);
2740f6efc10SHong 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);
2750f6efc10SHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
2760f6efc10SHong Zhang   if (flg) {
2770f6efc10SHong Zhang     switch (idx){
2780f6efc10SHong Zhang     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
2790f6efc10SHong Zhang     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
2800f6efc10SHong Zhang     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
2810f6efc10SHong Zhang     }
2820f6efc10SHong Zhang   }
2831677d0b8SKris 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);
2840f6efc10SHong 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);
2850f6efc10SHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);CHKERRQ(ierr);
2861677d0b8SKris Buschelman 
2871677d0b8SKris Buschelman   /* Control parameters used by solve */
2881677d0b8SKris Buschelman   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr);
2891677d0b8SKris Buschelman 
2900f6efc10SHong Zhang   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */
2912401956bSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-pc_factor_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr);
2921677d0b8SKris Buschelman   if (lu->PetscMatOdering) {
2930f6efc10SHong Zhang     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
2942d4e2982SHong Zhang     ierr = PetscMalloc(m*sizeof(UF_long),&lu->perm_c);CHKERRQ(ierr);
2952d4e2982SHong Zhang     /* we cannot simply memcpy on 64 bit archs */
2962d4e2982SHong Zhang     for(i = 0; i < m; i++) lu->perm_c[i] = ra[i];
2970f6efc10SHong Zhang     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
2981677d0b8SKris Buschelman   }
2991677d0b8SKris Buschelman   PetscOptionsEnd();
3001677d0b8SKris Buschelman 
3012d4e2982SHong Zhang   if(sizeof(UF_long) != sizeof(int)){
3022d4e2982SHong Zhang     /* we cannot directly use mat->i and mat->j on 64 bit archs */
3032d4e2982SHong Zhang     ierr = PetscMalloc((m+1)*sizeof(UF_long),&lu->ai);CHKERRQ(ierr);
3042d4e2982SHong Zhang     ierr = PetscMalloc(mat->nz*sizeof(UF_long),&lu->aj);CHKERRQ(ierr);
3052d4e2982SHong Zhang     for(i = 0; i < m + 1; i++) lu->ai[i] = mat->i[i];
3062d4e2982SHong Zhang     for(i = 0; i < mat->nz; i++) lu->aj[i] = mat->j[i];
3072d4e2982SHong Zhang   }
3082d4e2982SHong Zhang   else{
3092d4e2982SHong Zhang     lu->ai = (UF_long*)mat->i;
3102d4e2982SHong Zhang     lu->aj = (UF_long*)mat->j;
3112d4e2982SHong Zhang   }
3122d4e2982SHong Zhang 
3131677d0b8SKris Buschelman   /* print the control parameters */
3142d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
3152d4e2982SHong Zhang   if(lu->Control[UMFPACK_PRL] > 1) umfpack_zl_report_control(lu->Control);
3162d4e2982SHong Zhang #else
3172d4e2982SHong Zhang   if(lu->Control[UMFPACK_PRL] > 1) umfpack_dl_report_control(lu->Control);
3182d4e2982SHong Zhang #endif
3191677d0b8SKris Buschelman 
3201677d0b8SKris Buschelman   /* symbolic factorization of A' */
3211677d0b8SKris Buschelman   /* ---------------------------------------------------------------------- */
3222d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
3230f6efc10SHong Zhang   if (lu->PetscMatOdering) { /* use Petsc row ordering */
3242d4e2982SHong Zhang     status = umfpack_zl_qsymbolic(n,m,lu->ai,lu->aj,(double*)av,NULL,
3252d4e2982SHong Zhang 				  lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
3260f6efc10SHong Zhang   } else { /* use Umfpack col ordering */
3272d4e2982SHong Zhang     status = umfpack_zl_symbolic(n,m,lu->ai,lu->aj,(double*)av,NULL,
3282d4e2982SHong Zhang 				 &lu->Symbolic,lu->Control,lu->Info);
3290f6efc10SHong Zhang   }
3301677d0b8SKris Buschelman   if (status < 0){
3312d4e2982SHong Zhang     umfpack_zl_report_info(lu->Control, lu->Info);
3322d4e2982SHong Zhang     umfpack_zl_report_status(lu->Control, status);
3332d4e2982SHong Zhang     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_symbolic failed");
3341677d0b8SKris Buschelman   }
3351677d0b8SKris Buschelman   /* report sumbolic factorization of A' when Control[PRL] > 3 */
3362d4e2982SHong Zhang   (void) umfpack_zl_report_symbolic(lu->Symbolic, lu->Control);
3372d4e2982SHong Zhang #else
3382d4e2982SHong Zhang   if (lu->PetscMatOdering) { /* use Petsc row ordering */
3392d4e2982SHong Zhang     status = umfpack_dl_qsymbolic(n,m,lu->ai,lu->aj,av,
3402d4e2982SHong Zhang 				  lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
3412d4e2982SHong Zhang   } else { /* use Umfpack col ordering */
3422d4e2982SHong Zhang     status = umfpack_dl_symbolic(n,m,lu->ai,lu->aj,av,
3432d4e2982SHong Zhang 				 &lu->Symbolic,lu->Control,lu->Info);
3442d4e2982SHong Zhang   }
3452d4e2982SHong Zhang   if (status < 0){
3462d4e2982SHong Zhang     umfpack_dl_report_info(lu->Control, lu->Info);
3472d4e2982SHong Zhang     umfpack_dl_report_status(lu->Control, status);
3482d4e2982SHong Zhang     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_symbolic failed");
3492d4e2982SHong Zhang   }
3502d4e2982SHong Zhang   /* report sumbolic factorization of A' when Control[PRL] > 3 */
3512d4e2982SHong Zhang   (void) umfpack_dl_report_symbolic(lu->Symbolic, lu->Control);
3522d4e2982SHong Zhang #endif
3531677d0b8SKris Buschelman 
3541677d0b8SKris Buschelman   lu->flg = DIFFERENT_NONZERO_PATTERN;
3551677d0b8SKris Buschelman   lu->av  = av;
356e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK = PETSC_TRUE;
357e1b4c3a1SKris Buschelman   *F = B;
3581677d0b8SKris Buschelman   PetscFunctionReturn(0);
3591677d0b8SKris Buschelman }
3601677d0b8SKris Buschelman 
3611677d0b8SKris Buschelman #undef __FUNCT__
362f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_UMFPACK"
3636849ba73SBarry Smith PetscErrorCode MatAssemblyEnd_UMFPACK(Mat A,MatAssemblyType mode)
3646849ba73SBarry Smith {
365dfbe8321SBarry Smith   PetscErrorCode ierr;
366f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)(A->spptr);
367d844382dSKris Buschelman 
3681677d0b8SKris Buschelman   PetscFunctionBegin;
369d844382dSKris Buschelman   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
370d844382dSKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
371f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
3721677d0b8SKris Buschelman   PetscFunctionReturn(0);
3731677d0b8SKris Buschelman }
3741677d0b8SKris Buschelman 
37594b7f48cSBarry Smith /* used by -ksp_view */
3761677d0b8SKris Buschelman #undef __FUNCT__
377f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_UMFPACK"
3786849ba73SBarry Smith PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer)
3796849ba73SBarry Smith {
380f0c56d0fSKris Buschelman   Mat_UMFPACK *lu= (Mat_UMFPACK*)A->spptr;
381dfbe8321SBarry Smith   PetscErrorCode ierr;
382f0c56d0fSKris Buschelman 
3831677d0b8SKris Buschelman   PetscFunctionBegin;
3841677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
385f0c56d0fSKris Buschelman   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
3861677d0b8SKris Buschelman 
3871677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
3881677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
3891677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
3901677d0b8SKris Buschelman 
3911677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
3920f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr);
3931677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
3941677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
3950f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr);
3961677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
3970f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_2BY2_TOLERANCE]: %g\n",lu->Control[UMFPACK_2BY2_TOLERANCE]);CHKERRQ(ierr);
3980f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr);
3990f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr);
4001677d0b8SKris Buschelman 
4011677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
4021677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
4030f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr);
4040f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr);
4051677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
4060f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr);
4071677d0b8SKris Buschelman 
4081677d0b8SKris Buschelman   /* Control parameters used by solve */
4091677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
4101677d0b8SKris Buschelman 
4111677d0b8SKris Buschelman   /* mat ordering */
4121677d0b8SKris Buschelman   if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer,"  UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr);
4131677d0b8SKris Buschelman 
4141677d0b8SKris Buschelman   PetscFunctionReturn(0);
4151677d0b8SKris Buschelman }
4161677d0b8SKris Buschelman 
417d844382dSKris Buschelman #undef __FUNCT__
418f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_UMFPACK"
4196849ba73SBarry Smith PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
4206849ba73SBarry Smith {
421dfbe8321SBarry Smith   PetscErrorCode ierr;
42232077d6dSBarry Smith   PetscTruth        iascii;
423d844382dSKris Buschelman   PetscViewerFormat format;
424f0c56d0fSKris Buschelman   Mat_UMFPACK       *lu=(Mat_UMFPACK*)(A->spptr);
425d844382dSKris Buschelman 
426d844382dSKris Buschelman   PetscFunctionBegin;
427d844382dSKris Buschelman   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
428d844382dSKris Buschelman 
42932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
43032077d6dSBarry Smith   if (iascii) {
431d844382dSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
4322f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO) {
433f0c56d0fSKris Buschelman       ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
434d844382dSKris Buschelman     }
435d844382dSKris Buschelman   }
436d844382dSKris Buschelman   PetscFunctionReturn(0);
437d844382dSKris Buschelman }
438d844382dSKris Buschelman 
439d844382dSKris Buschelman EXTERN_C_BEGIN
440d844382dSKris Buschelman #undef __FUNCT__
441d844382dSKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_UMFPACK"
44275179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_UMFPACK(Mat A,MatType type,MatReuse reuse,Mat *newmat)
4436849ba73SBarry Smith {
444dfbe8321SBarry Smith   PetscErrorCode ierr;
445d844382dSKris Buschelman   Mat            B=*newmat;
446f0c56d0fSKris Buschelman   Mat_UMFPACK    *lu;
447d844382dSKris Buschelman 
448d844382dSKris Buschelman   PetscFunctionBegin;
449ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
450d844382dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
451d844382dSKris Buschelman   }
452d844382dSKris Buschelman 
453f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_UMFPACK,&lu);CHKERRQ(ierr);
454f0c56d0fSKris Buschelman   lu->MatDuplicate         = A->ops->duplicate;
455d844382dSKris Buschelman   lu->MatView              = A->ops->view;
456d844382dSKris Buschelman   lu->MatAssemblyEnd       = A->ops->assemblyend;
457d844382dSKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
458d844382dSKris Buschelman   lu->MatDestroy           = A->ops->destroy;
459d844382dSKris Buschelman   lu->CleanUpUMFPACK       = PETSC_FALSE;
460d844382dSKris Buschelman 
461d844382dSKris Buschelman   B->spptr                 = (void*)lu;
462f0c56d0fSKris Buschelman   B->ops->duplicate        = MatDuplicate_UMFPACK;
463f0c56d0fSKris Buschelman   B->ops->view             = MatView_UMFPACK;
464f0c56d0fSKris Buschelman   B->ops->assemblyend      = MatAssemblyEnd_UMFPACK;
465f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
466f0c56d0fSKris Buschelman   B->ops->destroy          = MatDestroy_UMFPACK;
467d844382dSKris Buschelman 
468d844382dSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_umfpack_C",
469d844382dSKris Buschelman                                            "MatConvert_SeqAIJ_UMFPACK",MatConvert_SeqAIJ_UMFPACK);CHKERRQ(ierr);
470d844382dSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_umfpack_seqaij_C",
471d844382dSKris Buschelman                                            "MatConvert_UMFPACK_SeqAIJ",MatConvert_UMFPACK_SeqAIJ);CHKERRQ(ierr);
472ae15b995SBarry Smith   ierr = PetscInfo(0,"Using UMFPACK for SeqAIJ LU factorization and solves.\n");CHKERRQ(ierr);
473d844382dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATUMFPACK);CHKERRQ(ierr);
474d844382dSKris Buschelman   *newmat = B;
475d844382dSKris Buschelman   PetscFunctionReturn(0);
476d844382dSKris Buschelman }
477d844382dSKris Buschelman EXTERN_C_END
478d844382dSKris Buschelman 
479f0c56d0fSKris Buschelman #undef __FUNCT__
480f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_UMFPACK"
4816849ba73SBarry Smith PetscErrorCode MatDuplicate_UMFPACK(Mat A, MatDuplicateOption op, Mat *M)
4826849ba73SBarry Smith {
483dfbe8321SBarry Smith   PetscErrorCode ierr;
4848f340917SKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)A->spptr;
4858f340917SKris Buschelman 
486f0c56d0fSKris Buschelman   PetscFunctionBegin;
4878f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
4883f953163SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_UMFPACK));CHKERRQ(ierr);
489f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
490f0c56d0fSKris Buschelman }
491f0c56d0fSKris Buschelman 
4922f71e704SKris Buschelman /*MC
493fafad747SKris Buschelman   MATUMFPACK - MATUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
4942f71e704SKris Buschelman   via the external package UMFPACK.
4952f71e704SKris Buschelman 
4962f71e704SKris Buschelman   If UMFPACK is installed (see the manual for
4972f71e704SKris Buschelman   instructions on how to declare the existence of external packages),
4982f71e704SKris Buschelman   a matrix type can be constructed which invokes UMFPACK solvers.
4992f71e704SKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,UMFPACK).
5002f71e704SKris Buschelman 
5012f71e704SKris Buschelman   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
50228b08bd3SKris Buschelman   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from
50328b08bd3SKris Buschelman   the MATSEQAIJ type without data copy.
5042f71e704SKris Buschelman 
5052f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
5062f71e704SKris Buschelman   which correspond to the options database keys below.
5072f71e704SKris Buschelman 
5082f71e704SKris Buschelman   Options Database Keys:
5090bad9183SKris Buschelman + -mat_type umfpack - sets the matrix type to "umfpack" during a call to MatSetFromOptions()
5102f71e704SKris Buschelman . -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL]
5112f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
5122f71e704SKris Buschelman . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
5132f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
5142f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
5152f71e704SKris Buschelman - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
5162f71e704SKris Buschelman 
5172f71e704SKris Buschelman    Level: beginner
5182f71e704SKris Buschelman 
5192f71e704SKris Buschelman .seealso: PCLU
5202f71e704SKris Buschelman M*/
5212f71e704SKris Buschelman 
5221677d0b8SKris Buschelman EXTERN_C_BEGIN
5231677d0b8SKris Buschelman #undef __FUNCT__
524f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_UMFPACK"
525be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_UMFPACK(Mat A)
526dfbe8321SBarry Smith {
527dfbe8321SBarry Smith   PetscErrorCode ierr;
5281677d0b8SKris Buschelman 
5291677d0b8SKris Buschelman   PetscFunctionBegin;
5301677d0b8SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
531ceb03754SKris Buschelman   ierr = MatConvert_SeqAIJ_UMFPACK(A,MATUMFPACK,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
5321677d0b8SKris Buschelman   PetscFunctionReturn(0);
5331677d0b8SKris Buschelman }
5341677d0b8SKris Buschelman EXTERN_C_END
5352d4e2982SHong Zhang 
5362d4e2982SHong Zhang 
5372d4e2982SHong Zhang 
538