xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision fcf6e9ab63b11f47efe2e1e248ac0bcd297e1e64)
11677d0b8SKris Buschelman 
21677d0b8SKris Buschelman /*
3d515b9b4SShri Abhyankar    Provides an interface to the UMFPACK sparse solver available through SuiteSparse version 4.2.1
42d4e2982SHong Zhang 
59e475b0dSSatish Balay    When build with PETSC_USE_64BIT_INDICES this will use Suitesparse_long as the
68592901bSBarry Smith    integer type in UMFPACK, otherwise it will use int. This means
78592901bSBarry Smith    all integers in this file as simply declared as PetscInt. Also it means
89e475b0dSSatish Balay    that one cannot use 64BIT_INDICES on 32bit machines [as Suitesparse_long is 32bit only]
92d4e2982SHong Zhang 
101677d0b8SKris Buschelman */
11c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
121677d0b8SKris Buschelman 
138592901bSBarry Smith #if defined(PETSC_USE_64BIT_INDICES)
148592901bSBarry Smith #if defined(PETSC_USE_COMPLEX)
158592901bSBarry Smith #define umfpack_UMF_free_symbolic                      umfpack_zl_free_symbolic
168592901bSBarry Smith #define umfpack_UMF_free_numeric                       umfpack_zl_free_numeric
17d755ee67SBarry Smith /* the type casts are needed because PetscInt is long long while SuiteSparse_long is long and compilers warn even when they are identical */
18d755ee67SBarry Smith #define umfpack_UMF_wsolve(a,b,c,d,e,f,g,h,i,j,k,l,m,n) umfpack_zl_wsolve(a,(SuiteSparse_long*)b,(SuiteSparse_long*)c,d,e,f,g,h,i,(SuiteSparse_long*)j,k,l,(SuiteSparse_long*)m,n)
19d755ee67SBarry Smith #define umfpack_UMF_numeric(a,b,c,d,e,f,g,h)          umfpack_zl_numeric((SuiteSparse_long*)a,(SuiteSparse_long*)b,c,d,e,f,g,h)
208592901bSBarry Smith #define umfpack_UMF_report_numeric                    umfpack_zl_report_numeric
218592901bSBarry Smith #define umfpack_UMF_report_control                    umfpack_zl_report_control
228592901bSBarry Smith #define umfpack_UMF_report_status                     umfpack_zl_report_status
238592901bSBarry Smith #define umfpack_UMF_report_info                       umfpack_zl_report_info
248592901bSBarry Smith #define umfpack_UMF_report_symbolic                   umfpack_zl_report_symbolic
25d755ee67SBarry Smith #define umfpack_UMF_qsymbolic(a,b,c,d,e,f,g,h,i,j)    umfpack_zl_qsymbolic(a,b,(SuiteSparse_long*)c,(SuiteSparse_long*)d,e,f,(SuiteSparse_long*)g,h,i,j)
26d755ee67SBarry Smith #define umfpack_UMF_symbolic(a,b,c,d,e,f,g,h,i)       umfpack_zl_symbolic(a,b,(SuiteSparse_long*)c,(SuiteSparse_long*)d,e,f,g,h,i)
278592901bSBarry Smith #define umfpack_UMF_defaults                          umfpack_zl_defaults
288592901bSBarry Smith 
298592901bSBarry Smith #else
308592901bSBarry Smith #define umfpack_UMF_free_symbolic                  umfpack_dl_free_symbolic
318592901bSBarry Smith #define umfpack_UMF_free_numeric                   umfpack_dl_free_numeric
32d755ee67SBarry Smith #define umfpack_UMF_wsolve(a,b,c,d,e,f,g,h,i,j,k)  umfpack_dl_wsolve(a,(SuiteSparse_long*)b,(SuiteSparse_long*)c,d,e,f,g,h,i,(SuiteSparse_long*)j,k)
33d755ee67SBarry Smith #define umfpack_UMF_numeric(a,b,c,d,e,f,g)         umfpack_dl_numeric((SuiteSparse_long*)a,(SuiteSparse_long*)b,c,d,e,f,g)
348592901bSBarry Smith #define umfpack_UMF_report_numeric                 umfpack_dl_report_numeric
358592901bSBarry Smith #define umfpack_UMF_report_control                 umfpack_dl_report_control
368592901bSBarry Smith #define umfpack_UMF_report_status                  umfpack_dl_report_status
378592901bSBarry Smith #define umfpack_UMF_report_info                    umfpack_dl_report_info
388592901bSBarry Smith #define umfpack_UMF_report_symbolic                umfpack_dl_report_symbolic
39d755ee67SBarry Smith #define umfpack_UMF_qsymbolic(a,b,c,d,e,f,g,h,i)   umfpack_dl_qsymbolic(a,b,(SuiteSparse_long*)c,(SuiteSparse_long*)d,e,(SuiteSparse_long*)f,g,h,i)
40d755ee67SBarry Smith #define umfpack_UMF_symbolic(a,b,c,d,e,f,g,h)      umfpack_dl_symbolic(a,b,(SuiteSparse_long*)c,(SuiteSparse_long*)d,e,f,g,h)
418592901bSBarry Smith #define umfpack_UMF_defaults                       umfpack_dl_defaults
428592901bSBarry Smith #endif
438592901bSBarry Smith 
448592901bSBarry Smith #else
458592901bSBarry Smith #if defined(PETSC_USE_COMPLEX)
468592901bSBarry Smith #define umfpack_UMF_free_symbolic   umfpack_zi_free_symbolic
478592901bSBarry Smith #define umfpack_UMF_free_numeric    umfpack_zi_free_numeric
488592901bSBarry Smith #define umfpack_UMF_wsolve          umfpack_zi_wsolve
498592901bSBarry Smith #define umfpack_UMF_numeric         umfpack_zi_numeric
508592901bSBarry Smith #define umfpack_UMF_report_numeric  umfpack_zi_report_numeric
518592901bSBarry Smith #define umfpack_UMF_report_control  umfpack_zi_report_control
528592901bSBarry Smith #define umfpack_UMF_report_status   umfpack_zi_report_status
538592901bSBarry Smith #define umfpack_UMF_report_info     umfpack_zi_report_info
548592901bSBarry Smith #define umfpack_UMF_report_symbolic umfpack_zi_report_symbolic
558592901bSBarry Smith #define umfpack_UMF_qsymbolic       umfpack_zi_qsymbolic
568592901bSBarry Smith #define umfpack_UMF_symbolic        umfpack_zi_symbolic
578592901bSBarry Smith #define umfpack_UMF_defaults        umfpack_zi_defaults
588592901bSBarry Smith 
598592901bSBarry Smith #else
608592901bSBarry Smith #define umfpack_UMF_free_symbolic   umfpack_di_free_symbolic
618592901bSBarry Smith #define umfpack_UMF_free_numeric    umfpack_di_free_numeric
628592901bSBarry Smith #define umfpack_UMF_wsolve          umfpack_di_wsolve
638592901bSBarry Smith #define umfpack_UMF_numeric         umfpack_di_numeric
648592901bSBarry Smith #define umfpack_UMF_report_numeric  umfpack_di_report_numeric
658592901bSBarry Smith #define umfpack_UMF_report_control  umfpack_di_report_control
668592901bSBarry Smith #define umfpack_UMF_report_status   umfpack_di_report_status
678592901bSBarry Smith #define umfpack_UMF_report_info     umfpack_di_report_info
688592901bSBarry Smith #define umfpack_UMF_report_symbolic umfpack_di_report_symbolic
698592901bSBarry Smith #define umfpack_UMF_qsymbolic       umfpack_di_qsymbolic
708592901bSBarry Smith #define umfpack_UMF_symbolic        umfpack_di_symbolic
718592901bSBarry Smith #define umfpack_UMF_defaults        umfpack_di_defaults
728592901bSBarry Smith #endif
738592901bSBarry Smith #endif
748592901bSBarry Smith 
751677d0b8SKris Buschelman EXTERN_C_BEGIN
76c6db04a5SJed Brown #include <umfpack.h>
771677d0b8SKris Buschelman EXTERN_C_END
781677d0b8SKris Buschelman 
79fcfd50ebSBarry Smith static const char *const UmfpackOrderingTypes[] = {"CHOLMOD","AMD","GIVEN","METIS","BEST","NONE","USER","UmfpackOrderingTypes","UMFPACK_ORDERING_",0};
80fcfd50ebSBarry Smith 
811677d0b8SKris Buschelman typedef struct {
821677d0b8SKris Buschelman   void         *Symbolic, *Numeric;
831677d0b8SKris Buschelman   double       Info[UMFPACK_INFO], Control[UMFPACK_CONTROL],*W;
84ae26a541SJed Brown   PetscInt     *Wi,*perm_c;
85ae26a541SJed Brown   Mat          A;               /* Matrix used for factorization */
861677d0b8SKris Buschelman   MatStructure flg;
87fcfd50ebSBarry Smith   PetscBool    PetscMatOrdering;
881677d0b8SKris Buschelman 
891677d0b8SKris Buschelman   /* Flag to clean up UMFPACK objects during Destroy */
90ace3abfcSBarry Smith   PetscBool CleanUpUMFPACK;
91f0c56d0fSKris Buschelman } Mat_UMFPACK;
92f0c56d0fSKris Buschelman 
931677d0b8SKris Buschelman #undef __FUNCT__
94f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_UMFPACK"
954b019bd2SJed Brown static PetscErrorCode MatDestroy_UMFPACK(Mat A)
96dfbe8321SBarry Smith {
97dfbe8321SBarry Smith   PetscErrorCode ierr;
98f0c56d0fSKris Buschelman   Mat_UMFPACK    *lu=(Mat_UMFPACK*)A->spptr;
991677d0b8SKris Buschelman 
1001677d0b8SKris Buschelman   PetscFunctionBegin;
101bf0cc555SLisandro Dalcin   if (lu && lu->CleanUpUMFPACK) {
1028592901bSBarry Smith     umfpack_UMF_free_symbolic(&lu->Symbolic);
1038592901bSBarry Smith     umfpack_UMF_free_numeric(&lu->Numeric);
1041677d0b8SKris Buschelman     ierr = PetscFree(lu->Wi);CHKERRQ(ierr);
1051677d0b8SKris Buschelman     ierr = PetscFree(lu->W);CHKERRQ(ierr);
1061677d0b8SKris Buschelman     ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
1071677d0b8SKris Buschelman   }
108ae26a541SJed Brown   ierr = MatDestroy(&lu->A);CHKERRQ(ierr);
109da0e1c47SShri Abhyankar   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
110b24902e0SBarry Smith   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
1111677d0b8SKris Buschelman   PetscFunctionReturn(0);
1121677d0b8SKris Buschelman }
1131677d0b8SKris Buschelman 
1141677d0b8SKris Buschelman #undef __FUNCT__
1154b019bd2SJed Brown #define __FUNCT__ "MatSolve_UMFPACK_Private"
1164b019bd2SJed Brown static PetscErrorCode MatSolve_UMFPACK_Private(Mat A,Vec b,Vec x,int uflag)
1176849ba73SBarry Smith {
118f0c56d0fSKris Buschelman   Mat_UMFPACK       *lu = (Mat_UMFPACK*)A->spptr;
119ae26a541SJed Brown   Mat_SeqAIJ        *a  = (Mat_SeqAIJ*)lu->A->data;
120d9ca1df4SBarry Smith   PetscScalar       *av = a->a,*xa;
121d9ca1df4SBarry Smith   const PetscScalar *ba;
122dfbe8321SBarry Smith   PetscErrorCode    ierr;
123ae26a541SJed Brown   PetscInt          *ai = a->i,*aj = a->j,status;
124*fcf6e9abSJed Brown   static PetscBool  cite = PETSC_FALSE;
1251677d0b8SKris Buschelman 
1261677d0b8SKris Buschelman   PetscFunctionBegin;
127*fcf6e9abSJed Brown   ierr = PetscCitationsRegister("@article{davis2004algorithm,\n  title={Algorithm 832: {UMFPACK} V4.3---An Unsymmetric-Pattern Multifrontal Method},\n  author={Davis, Timothy A},\n  journal={ACM Transactions on Mathematical Software (TOMS)},\n  volume={30},\n  number={2},\n  pages={196--199},\n  year={2004},\n  publisher={ACM}\n}\n",&cite);CHKERRQ(ierr);
1282d4e2982SHong Zhang   /* solve Ax = b by umfpack_*_wsolve */
1291677d0b8SKris Buschelman   /* ----------------------------------*/
1302d4e2982SHong Zhang 
131ae26a541SJed Brown   if (!lu->Wi) {  /* first time, allocate working space for wsolve */
132785e854fSJed Brown     ierr = PetscMalloc1(A->rmap->n,&lu->Wi);CHKERRQ(ierr);
133785e854fSJed Brown     ierr = PetscMalloc1(5*A->rmap->n,&lu->W);CHKERRQ(ierr);
134ae26a541SJed Brown   }
135ae26a541SJed Brown 
136d9ca1df4SBarry Smith   ierr = VecGetArrayRead(b,&ba);
1371677d0b8SKris Buschelman   ierr = VecGetArray(x,&xa);
13879c34000SHong Zhang #if defined(PETSC_USE_COMPLEX)
139b0043f70SBarry Smith   status = umfpack_UMF_wsolve(uflag,ai,aj,(PetscReal*)av,NULL,(PetscReal*)xa,NULL,(PetscReal*)ba,NULL,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
14079c34000SHong Zhang #else
1414b019bd2SJed Brown   status = umfpack_UMF_wsolve(uflag,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
14279c34000SHong Zhang #endif
1438592901bSBarry Smith   umfpack_UMF_report_info(lu->Control, lu->Info);
1441677d0b8SKris Buschelman   if (status < 0) {
1458592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
146e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_wsolve failed");
1471677d0b8SKris Buschelman   }
1481677d0b8SKris Buschelman 
149d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(b,&ba);CHKERRQ(ierr);
150057dcbc9SJed Brown   ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr);
1514b019bd2SJed Brown   PetscFunctionReturn(0);
1524b019bd2SJed Brown }
1532a325a84SHong Zhang 
1544b019bd2SJed Brown #undef __FUNCT__
1554b019bd2SJed Brown #define __FUNCT__ "MatSolve_UMFPACK"
1564b019bd2SJed Brown static PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x)
1574b019bd2SJed Brown {
1584b019bd2SJed Brown   PetscErrorCode ierr;
1594b019bd2SJed Brown 
1604b019bd2SJed Brown   PetscFunctionBegin;
1614b019bd2SJed Brown   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
1624b019bd2SJed Brown   ierr = MatSolve_UMFPACK_Private(A,b,x,UMFPACK_Aat);CHKERRQ(ierr);
1634b019bd2SJed Brown   PetscFunctionReturn(0);
1644b019bd2SJed Brown }
1654b019bd2SJed Brown 
1664b019bd2SJed Brown #undef __FUNCT__
1674b019bd2SJed Brown #define __FUNCT__ "MatSolveTranspose_UMFPACK"
1684b019bd2SJed Brown static PetscErrorCode MatSolveTranspose_UMFPACK(Mat A,Vec b,Vec x)
1694b019bd2SJed Brown {
1704b019bd2SJed Brown   PetscErrorCode ierr;
1714b019bd2SJed Brown 
1724b019bd2SJed Brown   PetscFunctionBegin;
1734b019bd2SJed Brown   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
1744b019bd2SJed Brown   ierr = MatSolve_UMFPACK_Private(A,b,x,UMFPACK_A);CHKERRQ(ierr);
1751677d0b8SKris Buschelman   PetscFunctionReturn(0);
1761677d0b8SKris Buschelman }
1771677d0b8SKris Buschelman 
1781677d0b8SKris Buschelman #undef __FUNCT__
179f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_UMFPACK"
1804b019bd2SJed Brown static PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat F,Mat A,const MatFactorInfo *info)
1816849ba73SBarry Smith {
182719d5645SBarry Smith   Mat_UMFPACK    *lu = (Mat_UMFPACK*)(F)->spptr;
183ae26a541SJed Brown   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
184ae26a541SJed Brown   PetscInt       *ai = a->i,*aj=a->j,status;
185ae26a541SJed Brown   PetscScalar    *av = a->a;
1866849ba73SBarry Smith   PetscErrorCode ierr;
1871677d0b8SKris Buschelman 
1881677d0b8SKris Buschelman   PetscFunctionBegin;
1891677d0b8SKris Buschelman   /* numeric factorization of A' */
1901677d0b8SKris Buschelman   /* ----------------------------*/
1912d4e2982SHong Zhang 
1928592901bSBarry Smith   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) {
1938592901bSBarry Smith     umfpack_UMF_free_numeric(&lu->Numeric);
1948592901bSBarry Smith   }
1952d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1968592901bSBarry Smith   status = umfpack_UMF_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1972d4e2982SHong Zhang #else
1988592901bSBarry Smith   status = umfpack_UMF_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1998592901bSBarry Smith #endif
2009f42a82aSMatthew Knepley   if (status < 0) {
2018592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
202e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_numeric failed");
2039f42a82aSMatthew Knepley   }
2042d4e2982SHong Zhang   /* report numeric factorization of A' when Control[PRL] > 3 */
2058592901bSBarry Smith   (void) umfpack_UMF_report_numeric(lu->Numeric, lu->Control);
2061677d0b8SKris Buschelman 
207ae26a541SJed Brown   ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
208ae26a541SJed Brown   ierr = MatDestroy(&lu->A);CHKERRQ(ierr);
2092205254eSKarl Rupp 
210ae26a541SJed Brown   lu->A                  = A;
2112a325a84SHong Zhang   lu->flg                = SAME_NONZERO_PATTERN;
2122a325a84SHong Zhang   lu->CleanUpUMFPACK     = PETSC_TRUE;
2134b019bd2SJed Brown   F->ops->solve          = MatSolve_UMFPACK;
2144b019bd2SJed Brown   F->ops->solvetranspose = MatSolveTranspose_UMFPACK;
2151677d0b8SKris Buschelman   PetscFunctionReturn(0);
2161677d0b8SKris Buschelman }
2171677d0b8SKris Buschelman 
2181677d0b8SKris Buschelman /*
2191677d0b8SKris Buschelman    Note the r permutation is ignored
2201677d0b8SKris Buschelman */
2211677d0b8SKris Buschelman #undef __FUNCT__
222f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK"
2234b019bd2SJed Brown static PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
2246849ba73SBarry Smith {
225ae26a541SJed Brown   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
226719d5645SBarry Smith   Mat_UMFPACK    *lu = (Mat_UMFPACK*)(F->spptr);
227dfbe8321SBarry Smith   PetscErrorCode ierr;
228ae26a541SJed Brown   PetscInt       i,*ai = a->i,*aj = a->j,m=A->rmap->n,n=A->cmap->n;
229989213a4SSatish Balay #if !defined(PETSC_USE_COMPLEX)
230ae26a541SJed Brown   PetscScalar    *av = a->a;
231989213a4SSatish Balay #endif
2325d0c19d7SBarry Smith   const PetscInt *ra;
2338592901bSBarry Smith   PetscInt       status;
2341677d0b8SKris Buschelman 
2351677d0b8SKris Buschelman   PetscFunctionBegin;
236fcfd50ebSBarry Smith   if (lu->PetscMatOrdering) {
2370f6efc10SHong Zhang     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
238785e854fSJed Brown     ierr = PetscMalloc1(m,&lu->perm_c);CHKERRQ(ierr);
2392d4e2982SHong Zhang     /* we cannot simply memcpy on 64 bit archs */
2402d4e2982SHong Zhang     for (i = 0; i < m; i++) lu->perm_c[i] = ra[i];
2410f6efc10SHong Zhang     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
2421677d0b8SKris Buschelman   }
2431677d0b8SKris Buschelman 
2441677d0b8SKris Buschelman   /* print the control parameters */
2458592901bSBarry Smith   if (lu->Control[UMFPACK_PRL] > 1) umfpack_UMF_report_control(lu->Control);
2461677d0b8SKris Buschelman 
2471677d0b8SKris Buschelman   /* symbolic factorization of A' */
2481677d0b8SKris Buschelman   /* ---------------------------------------------------------------------- */
249fcfd50ebSBarry Smith   if (lu->PetscMatOrdering) { /* use Petsc row ordering */
2508592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
251ae26a541SJed Brown     status = umfpack_UMF_qsymbolic(n,m,ai,aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2522d4e2982SHong Zhang #else
253b0043f70SBarry Smith     status = umfpack_UMF_qsymbolic(n,m,ai,aj,NULL,NULL,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2548592901bSBarry Smith #endif
2552d4e2982SHong Zhang   } else { /* use Umfpack col ordering */
2568592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
257ae26a541SJed Brown     status = umfpack_UMF_symbolic(n,m,ai,aj,av,&lu->Symbolic,lu->Control,lu->Info);
2588592901bSBarry Smith #else
259ae26a541SJed Brown     status = umfpack_UMF_symbolic(n,m,ai,aj,NULL,NULL,&lu->Symbolic,lu->Control,lu->Info);
2608592901bSBarry Smith #endif
2612d4e2982SHong Zhang   }
2622d4e2982SHong Zhang   if (status < 0) {
2638592901bSBarry Smith     umfpack_UMF_report_info(lu->Control, lu->Info);
2648592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
265e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_symbolic failed");
2662d4e2982SHong Zhang   }
2672d4e2982SHong Zhang   /* report sumbolic factorization of A' when Control[PRL] > 3 */
2688592901bSBarry Smith   (void) umfpack_UMF_report_symbolic(lu->Symbolic, lu->Control);
2691677d0b8SKris Buschelman 
2701677d0b8SKris Buschelman   lu->flg                   = DIFFERENT_NONZERO_PATTERN;
271e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK        = PETSC_TRUE;
272719d5645SBarry Smith   (F)->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
2731677d0b8SKris Buschelman   PetscFunctionReturn(0);
2741677d0b8SKris Buschelman }
2751677d0b8SKris Buschelman 
2761677d0b8SKris Buschelman #undef __FUNCT__
277f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_UMFPACK"
2784b019bd2SJed Brown static PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer)
2796849ba73SBarry Smith {
280f0c56d0fSKris Buschelman   Mat_UMFPACK    *lu= (Mat_UMFPACK*)A->spptr;
281dfbe8321SBarry Smith   PetscErrorCode ierr;
282f0c56d0fSKris Buschelman 
2831677d0b8SKris Buschelman   PetscFunctionBegin;
2841677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
285f0c56d0fSKris Buschelman   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
2861677d0b8SKris Buschelman 
2871677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
2881677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
2891677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
2901677d0b8SKris Buschelman 
2911677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
2920f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr);
2931677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
2941677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
2950f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr);
2961677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
2970f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr);
2980f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr);
2991677d0b8SKris Buschelman 
3001677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
3011677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
3020f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr);
3030f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr);
3041677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
3050f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr);
3061677d0b8SKris Buschelman 
3071677d0b8SKris Buschelman   /* Control parameters used by solve */
3081677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
3091677d0b8SKris Buschelman 
3101677d0b8SKris Buschelman   /* mat ordering */
311fcfd50ebSBarry Smith   if (!lu->PetscMatOrdering) {
312fcfd50ebSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ORDERING]: %s (not using the PETSc ordering)\n",UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]]);CHKERRQ(ierr);
313fcfd50ebSBarry Smith   }
3141677d0b8SKris Buschelman   PetscFunctionReturn(0);
3151677d0b8SKris Buschelman }
3161677d0b8SKris Buschelman 
317d844382dSKris Buschelman #undef __FUNCT__
318f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_UMFPACK"
3194b019bd2SJed Brown static PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
3206849ba73SBarry Smith {
321dfbe8321SBarry Smith   PetscErrorCode    ierr;
322ace3abfcSBarry Smith   PetscBool         iascii;
323d844382dSKris Buschelman   PetscViewerFormat format;
324d844382dSKris Buschelman 
325d844382dSKris Buschelman   PetscFunctionBegin;
326b24902e0SBarry Smith   ierr = MatView_SeqAIJ(A,viewer);CHKERRQ(ierr);
327d844382dSKris Buschelman 
328251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
32932077d6dSBarry Smith   if (iascii) {
330d844382dSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
3312f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO) {
332f0c56d0fSKris Buschelman       ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
333d844382dSKris Buschelman     }
334d844382dSKris Buschelman   }
335d844382dSKris Buschelman   PetscFunctionReturn(0);
336d844382dSKris Buschelman }
337d844382dSKris Buschelman 
33835bd34faSBarry Smith #undef __FUNCT__
33935bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_umfpack"
34035bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_seqaij_umfpack(Mat A,const MatSolverPackage *type)
34135bd34faSBarry Smith {
34235bd34faSBarry Smith   PetscFunctionBegin;
3432692d6eeSBarry Smith   *type = MATSOLVERUMFPACK;
34435bd34faSBarry Smith   PetscFunctionReturn(0);
34535bd34faSBarry Smith }
34635bd34faSBarry Smith 
34735bd34faSBarry Smith 
3482f71e704SKris Buschelman /*MC
3492692d6eeSBarry Smith   MATSOLVERUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
3502f71e704SKris Buschelman   via the external package UMFPACK.
3512f71e704SKris Buschelman 
352c2b89b5dSBarry Smith   Use ./configure --download-suitesparse to install PETSc to use UMFPACK
353c2b89b5dSBarry Smith 
354c2b89b5dSBarry Smith   Use -pc_type lu -pc_factor_mat_solver_package umfpack to us this direct solver
3552f71e704SKris Buschelman 
3562f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
3572f71e704SKris Buschelman   which correspond to the options database keys below.
3582f71e704SKris Buschelman 
3592f71e704SKris Buschelman   Options Database Keys:
360ba41fbb6SBarry Smith + -mat_umfpack_ordering                - CHOLMOD, AMD, GIVEN, METIS, BEST, NONE
361ba41fbb6SBarry Smith . -mat_umfpack_prl                     - UMFPACK print level: Control[UMFPACK_PRL]
362e08999f5SMatthew G Knepley . -mat_umfpack_strategy <AUTO>         - (choose one of) AUTO UNSYMMETRIC SYMMETRIC 2BY2
3632f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c>     - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
364e08999f5SMatthew G Knepley . -mat_umfpack_dense_row <0.2>         - Control[UMFPACK_DENSE_ROW]
365e08999f5SMatthew G Knepley . -mat_umfpack_amd_dense <10>          - Control[UMFPACK_AMD_DENSE]
3662f71e704SKris Buschelman . -mat_umfpack_block_size <bs>         - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
367e08999f5SMatthew G Knepley . -mat_umfpack_2by2_tolerance <0.01>   - Control[UMFPACK_2BY2_TOLERANCE]
368e08999f5SMatthew G Knepley . -mat_umfpack_fixq <0>                - Control[UMFPACK_FIXQ]
369e08999f5SMatthew G Knepley . -mat_umfpack_aggressive <1>          - Control[UMFPACK_AGGRESSIVE]
3702f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
371e08999f5SMatthew G Knepley . -mat_umfpack_sym_pivot_tolerance <0.001> - Control[UMFPACK_SYM_PIVOT_TOLERANCE]
372e08999f5SMatthew G Knepley . -mat_umfpack_scale <NONE>           - (choose one of) NONE SUM MAX
3732f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta>      - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
374e08999f5SMatthew G Knepley . -mat_umfpack_droptol <0>            - Control[UMFPACK_DROPTOL]
3752f71e704SKris Buschelman - -mat_umfpack_irstep <maxit>          - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
3762f71e704SKris Buschelman 
3772f71e704SKris Buschelman    Level: beginner
378a364b7d2SBarry Smith 
379a364b7d2SBarry Smith    Note: UMFPACK is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html
3802f71e704SKris Buschelman 
381d45987f3SHong Zhang .seealso: PCLU, MATSOLVERSUPERLU, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage
3822f71e704SKris Buschelman M*/
383b2573a8aSBarry Smith 
3843519fcdcSHong Zhang #undef __FUNCT__
3853519fcdcSHong Zhang #define __FUNCT__ "MatGetFactor_seqaij_umfpack"
3868cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_umfpack(Mat A,MatFactorType ftype,Mat *F)
3873519fcdcSHong Zhang {
3883519fcdcSHong Zhang   Mat            B;
3893519fcdcSHong Zhang   Mat_UMFPACK    *lu;
3903519fcdcSHong Zhang   PetscErrorCode ierr;
3918592901bSBarry Smith   PetscInt       m=A->rmap->n,n=A->cmap->n,idx;
3922f71e704SKris Buschelman 
3932205254eSKarl Rupp   const char *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC"};
3942205254eSKarl Rupp   const char *scale[]   ={"NONE","SUM","MAX"};
395ace3abfcSBarry Smith   PetscBool  flg;
3962d4e2982SHong Zhang 
3973519fcdcSHong Zhang   PetscFunctionBegin;
3983519fcdcSHong Zhang   /* Create the factorization matrix F */
399ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
4003519fcdcSHong Zhang   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
4013519fcdcSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
4020298fd71SBarry Smith   ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
403b00a9115SJed Brown   ierr = PetscNewLog(B,&lu);CHKERRQ(ierr);
4042205254eSKarl Rupp 
4053519fcdcSHong Zhang   B->spptr                 = lu;
4063519fcdcSHong Zhang   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
4073519fcdcSHong Zhang   B->ops->destroy          = MatDestroy_UMFPACK;
4083519fcdcSHong Zhang   B->ops->view             = MatView_UMFPACK;
4091aef8b4cSStefano Zampini   B->ops->matsolve         = NULL;
4102205254eSKarl Rupp 
411bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_umfpack);CHKERRQ(ierr);
4122205254eSKarl Rupp 
413d5f3da31SBarry Smith   B->factortype   = MAT_FACTOR_LU;
4143519fcdcSHong Zhang   B->assembled    = PETSC_TRUE;           /* required by -ksp_view */
4153519fcdcSHong Zhang   B->preallocated = PETSC_TRUE;
4163519fcdcSHong Zhang 
4173519fcdcSHong Zhang   /* initializations */
4183519fcdcSHong Zhang   /* ------------------------------------------------*/
4193519fcdcSHong Zhang   /* get the default control parameters */
4208592901bSBarry Smith   umfpack_UMF_defaults(lu->Control);
4210298fd71SBarry Smith   lu->perm_c                  = NULL; /* use defaul UMFPACK col permutation */
4223519fcdcSHong Zhang   lu->Control[UMFPACK_IRSTEP] = 0;          /* max num of iterative refinement steps to attempt */
4233519fcdcSHong Zhang 
424ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
4253519fcdcSHong Zhang   /* Control parameters used by reporting routiones */
4260298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],NULL);CHKERRQ(ierr);
4273519fcdcSHong Zhang 
4283519fcdcSHong Zhang   /* Control parameters for symbolic factorization */
429fcfd50ebSBarry Smith   ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,3,strategy[0],&idx,&flg);CHKERRQ(ierr);
4303519fcdcSHong Zhang   if (flg) {
4313519fcdcSHong Zhang     switch (idx) {
4323519fcdcSHong Zhang     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
4333519fcdcSHong Zhang     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
4343519fcdcSHong Zhang     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
4353519fcdcSHong Zhang     }
4363519fcdcSHong Zhang   }
4378caf3d72SBarry Smith   ierr = PetscOptionsEList("-mat_umfpack_ordering","Internal ordering method","None",UmfpackOrderingTypes,sizeof(UmfpackOrderingTypes)/sizeof(UmfpackOrderingTypes[0]),UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]],&idx,&flg);CHKERRQ(ierr);
438fcfd50ebSBarry Smith   if (flg) lu->Control[UMFPACK_ORDERING] = (int)idx;
4390298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_dense_col","Control[UMFPACK_DENSE_COL]","None",lu->Control[UMFPACK_DENSE_COL],&lu->Control[UMFPACK_DENSE_COL],NULL);CHKERRQ(ierr);
4400298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_dense_row","Control[UMFPACK_DENSE_ROW]","None",lu->Control[UMFPACK_DENSE_ROW],&lu->Control[UMFPACK_DENSE_ROW],NULL);CHKERRQ(ierr);
4410298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_amd_dense","Control[UMFPACK_AMD_DENSE]","None",lu->Control[UMFPACK_AMD_DENSE],&lu->Control[UMFPACK_AMD_DENSE],NULL);CHKERRQ(ierr);
4420298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_block_size","Control[UMFPACK_BLOCK_SIZE]","None",lu->Control[UMFPACK_BLOCK_SIZE],&lu->Control[UMFPACK_BLOCK_SIZE],NULL);CHKERRQ(ierr);
4430298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],NULL);CHKERRQ(ierr);
4440298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],NULL);CHKERRQ(ierr);
4453519fcdcSHong Zhang 
4463519fcdcSHong Zhang   /* Control parameters used by numeric factorization */
4470298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_pivot_tolerance","Control[UMFPACK_PIVOT_TOLERANCE]","None",lu->Control[UMFPACK_PIVOT_TOLERANCE],&lu->Control[UMFPACK_PIVOT_TOLERANCE],NULL);CHKERRQ(ierr);
4480298fd71SBarry Smith   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],NULL);CHKERRQ(ierr);
4493519fcdcSHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
4503519fcdcSHong Zhang   if (flg) {
4513519fcdcSHong Zhang     switch (idx) {
4523519fcdcSHong Zhang     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
4533519fcdcSHong Zhang     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
4543519fcdcSHong Zhang     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
4553519fcdcSHong Zhang     }
4563519fcdcSHong Zhang   }
4570298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_alloc_init","Control[UMFPACK_ALLOC_INIT]","None",lu->Control[UMFPACK_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],NULL);CHKERRQ(ierr);
4580298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_front_alloc_init","Control[UMFPACK_FRONT_ALLOC_INIT]","None",lu->Control[UMFPACK_FRONT_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],NULL);CHKERRQ(ierr);
4590298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],NULL);CHKERRQ(ierr);
4603519fcdcSHong Zhang 
4613519fcdcSHong Zhang   /* Control parameters used by solve */
4620298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],NULL);CHKERRQ(ierr);
4633519fcdcSHong Zhang 
4643519fcdcSHong Zhang   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */
46576a34f28SBarry Smith   ierr = PetscOptionsName("-pc_factor_mat_ordering_type","Ordering to do factorization in","MatGetOrdering",&lu->PetscMatOrdering);CHKERRQ(ierr);
4663519fcdcSHong Zhang   PetscOptionsEnd();
4673519fcdcSHong Zhang   *F = B;
4683519fcdcSHong Zhang   PetscFunctionReturn(0);
4693519fcdcSHong Zhang }
470b2573a8aSBarry Smith 
47142c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat,MatFactorType,Mat*);
47242c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat,MatFactorType,Mat*);
47342c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat,MatFactorType,Mat*);
4742d4e2982SHong Zhang 
47542c9c57cSBarry Smith #undef __FUNCT__
47642c9c57cSBarry Smith #define __FUNCT__ "MatSolverPackageRegister_SuiteSparse"
47729b38603SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_SuiteSparse(void)
47842c9c57cSBarry Smith {
47942c9c57cSBarry Smith   PetscErrorCode ierr;
48042c9c57cSBarry Smith 
48142c9c57cSBarry Smith   PetscFunctionBegin;
48242c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERUMFPACK,MATSEQAIJ,      MAT_FACTOR_LU,MatGetFactor_seqaij_umfpack);CHKERRQ(ierr);
48342c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERCHOLMOD,MATSEQAIJ,      MAT_FACTOR_CHOLESKY,MatGetFactor_seqaij_cholmod);CHKERRQ(ierr);
48442c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERCHOLMOD,MATSEQSBAIJ,      MAT_FACTOR_CHOLESKY,MatGetFactor_seqsbaij_cholmod);CHKERRQ(ierr);
48542c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERKLU,MATSEQAIJ,          MAT_FACTOR_LU,MatGetFactor_seqaij_klu);CHKERRQ(ierr);
48642c9c57cSBarry Smith   PetscFunctionReturn(0);
48742c9c57cSBarry Smith }
488