xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 989213a406d9b6f86c7753f25e76e52d89245077)
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
178592901bSBarry Smith #define umfpack_UMF_wsolve          umfpack_zl_wsolve
188592901bSBarry Smith #define umfpack_UMF_numeric         umfpack_zl_numeric
198592901bSBarry Smith #define umfpack_UMF_report_numeric  umfpack_zl_report_numeric
208592901bSBarry Smith #define umfpack_UMF_report_control  umfpack_zl_report_control
218592901bSBarry Smith #define umfpack_UMF_report_status   umfpack_zl_report_status
228592901bSBarry Smith #define umfpack_UMF_report_info     umfpack_zl_report_info
238592901bSBarry Smith #define umfpack_UMF_report_symbolic umfpack_zl_report_symbolic
248592901bSBarry Smith #define umfpack_UMF_qsymbolic       umfpack_zl_qsymbolic
258592901bSBarry Smith #define umfpack_UMF_symbolic        umfpack_zl_symbolic
268592901bSBarry Smith #define umfpack_UMF_defaults        umfpack_zl_defaults
278592901bSBarry Smith 
288592901bSBarry Smith #else
298592901bSBarry Smith #define umfpack_UMF_free_symbolic   umfpack_dl_free_symbolic
308592901bSBarry Smith #define umfpack_UMF_free_numeric    umfpack_dl_free_numeric
318592901bSBarry Smith #define umfpack_UMF_wsolve          umfpack_dl_wsolve
328592901bSBarry Smith #define umfpack_UMF_numeric         umfpack_dl_numeric
338592901bSBarry Smith #define umfpack_UMF_report_numeric  umfpack_dl_report_numeric
348592901bSBarry Smith #define umfpack_UMF_report_control  umfpack_dl_report_control
358592901bSBarry Smith #define umfpack_UMF_report_status   umfpack_dl_report_status
368592901bSBarry Smith #define umfpack_UMF_report_info     umfpack_dl_report_info
378592901bSBarry Smith #define umfpack_UMF_report_symbolic umfpack_dl_report_symbolic
388592901bSBarry Smith #define umfpack_UMF_qsymbolic       umfpack_dl_qsymbolic
398592901bSBarry Smith #define umfpack_UMF_symbolic        umfpack_dl_symbolic
408592901bSBarry Smith #define umfpack_UMF_defaults        umfpack_dl_defaults
418592901bSBarry Smith #endif
428592901bSBarry Smith 
438592901bSBarry Smith #else
448592901bSBarry Smith #if defined(PETSC_USE_COMPLEX)
458592901bSBarry Smith #define umfpack_UMF_free_symbolic   umfpack_zi_free_symbolic
468592901bSBarry Smith #define umfpack_UMF_free_numeric    umfpack_zi_free_numeric
478592901bSBarry Smith #define umfpack_UMF_wsolve          umfpack_zi_wsolve
488592901bSBarry Smith #define umfpack_UMF_numeric         umfpack_zi_numeric
498592901bSBarry Smith #define umfpack_UMF_report_numeric  umfpack_zi_report_numeric
508592901bSBarry Smith #define umfpack_UMF_report_control  umfpack_zi_report_control
518592901bSBarry Smith #define umfpack_UMF_report_status   umfpack_zi_report_status
528592901bSBarry Smith #define umfpack_UMF_report_info     umfpack_zi_report_info
538592901bSBarry Smith #define umfpack_UMF_report_symbolic umfpack_zi_report_symbolic
548592901bSBarry Smith #define umfpack_UMF_qsymbolic       umfpack_zi_qsymbolic
558592901bSBarry Smith #define umfpack_UMF_symbolic        umfpack_zi_symbolic
568592901bSBarry Smith #define umfpack_UMF_defaults        umfpack_zi_defaults
578592901bSBarry Smith 
588592901bSBarry Smith #else
598592901bSBarry Smith #define umfpack_UMF_free_symbolic   umfpack_di_free_symbolic
608592901bSBarry Smith #define umfpack_UMF_free_numeric    umfpack_di_free_numeric
618592901bSBarry Smith #define umfpack_UMF_wsolve          umfpack_di_wsolve
628592901bSBarry Smith #define umfpack_UMF_numeric         umfpack_di_numeric
638592901bSBarry Smith #define umfpack_UMF_report_numeric  umfpack_di_report_numeric
648592901bSBarry Smith #define umfpack_UMF_report_control  umfpack_di_report_control
658592901bSBarry Smith #define umfpack_UMF_report_status   umfpack_di_report_status
668592901bSBarry Smith #define umfpack_UMF_report_info     umfpack_di_report_info
678592901bSBarry Smith #define umfpack_UMF_report_symbolic umfpack_di_report_symbolic
688592901bSBarry Smith #define umfpack_UMF_qsymbolic       umfpack_di_qsymbolic
698592901bSBarry Smith #define umfpack_UMF_symbolic        umfpack_di_symbolic
708592901bSBarry Smith #define umfpack_UMF_defaults        umfpack_di_defaults
718592901bSBarry Smith #endif
728592901bSBarry Smith #endif
738592901bSBarry Smith 
741677d0b8SKris Buschelman EXTERN_C_BEGIN
75c6db04a5SJed Brown #include <umfpack.h>
761677d0b8SKris Buschelman EXTERN_C_END
771677d0b8SKris Buschelman 
78fcfd50ebSBarry Smith static const char *const UmfpackOrderingTypes[] = {"CHOLMOD","AMD","GIVEN","METIS","BEST","NONE","USER","UmfpackOrderingTypes","UMFPACK_ORDERING_",0};
79fcfd50ebSBarry Smith 
801677d0b8SKris Buschelman typedef struct {
811677d0b8SKris Buschelman   void         *Symbolic, *Numeric;
821677d0b8SKris Buschelman   double       Info[UMFPACK_INFO], Control[UMFPACK_CONTROL],*W;
83ae26a541SJed Brown   PetscInt     *Wi,*perm_c;
84ae26a541SJed Brown   Mat          A;               /* Matrix used for factorization */
851677d0b8SKris Buschelman   MatStructure flg;
86fcfd50ebSBarry Smith   PetscBool    PetscMatOrdering;
871677d0b8SKris Buschelman 
881677d0b8SKris Buschelman   /* Flag to clean up UMFPACK objects during Destroy */
89ace3abfcSBarry Smith   PetscBool CleanUpUMFPACK;
90f0c56d0fSKris Buschelman } Mat_UMFPACK;
91f0c56d0fSKris Buschelman 
921677d0b8SKris Buschelman #undef __FUNCT__
93f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_UMFPACK"
944b019bd2SJed Brown static PetscErrorCode MatDestroy_UMFPACK(Mat A)
95dfbe8321SBarry Smith {
96dfbe8321SBarry Smith   PetscErrorCode ierr;
97f0c56d0fSKris Buschelman   Mat_UMFPACK    *lu=(Mat_UMFPACK*)A->spptr;
981677d0b8SKris Buschelman 
991677d0b8SKris Buschelman   PetscFunctionBegin;
100bf0cc555SLisandro Dalcin   if (lu && lu->CleanUpUMFPACK) {
1018592901bSBarry Smith     umfpack_UMF_free_symbolic(&lu->Symbolic);
1028592901bSBarry Smith     umfpack_UMF_free_numeric(&lu->Numeric);
1031677d0b8SKris Buschelman     ierr = PetscFree(lu->Wi);CHKERRQ(ierr);
1041677d0b8SKris Buschelman     ierr = PetscFree(lu->W);CHKERRQ(ierr);
1051677d0b8SKris Buschelman     ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
1061677d0b8SKris Buschelman   }
107ae26a541SJed Brown   ierr = MatDestroy(&lu->A);CHKERRQ(ierr);
108da0e1c47SShri Abhyankar   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
109b24902e0SBarry Smith   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
1101677d0b8SKris Buschelman   PetscFunctionReturn(0);
1111677d0b8SKris Buschelman }
1121677d0b8SKris Buschelman 
1131677d0b8SKris Buschelman #undef __FUNCT__
1144b019bd2SJed Brown #define __FUNCT__ "MatSolve_UMFPACK_Private"
1154b019bd2SJed Brown static PetscErrorCode MatSolve_UMFPACK_Private(Mat A,Vec b,Vec x,int uflag)
1166849ba73SBarry Smith {
117f0c56d0fSKris Buschelman   Mat_UMFPACK       *lu = (Mat_UMFPACK*)A->spptr;
118ae26a541SJed Brown   Mat_SeqAIJ        *a  = (Mat_SeqAIJ*)lu->A->data;
119d9ca1df4SBarry Smith   PetscScalar       *av = a->a,*xa;
120d9ca1df4SBarry Smith   const PetscScalar *ba;
121dfbe8321SBarry Smith   PetscErrorCode    ierr;
122ae26a541SJed Brown   PetscInt          *ai = a->i,*aj = a->j,status;
1231677d0b8SKris Buschelman 
1241677d0b8SKris Buschelman   PetscFunctionBegin;
1252d4e2982SHong Zhang   /* solve Ax = b by umfpack_*_wsolve */
1261677d0b8SKris Buschelman   /* ----------------------------------*/
1272d4e2982SHong Zhang 
128ae26a541SJed Brown   if (!lu->Wi) {  /* first time, allocate working space for wsolve */
129785e854fSJed Brown     ierr = PetscMalloc1(A->rmap->n,&lu->Wi);CHKERRQ(ierr);
130785e854fSJed Brown     ierr = PetscMalloc1(5*A->rmap->n,&lu->W);CHKERRQ(ierr);
131ae26a541SJed Brown   }
132ae26a541SJed Brown 
133d9ca1df4SBarry Smith   ierr = VecGetArrayRead(b,&ba);
1341677d0b8SKris Buschelman   ierr = VecGetArray(x,&xa);
13579c34000SHong Zhang #if defined(PETSC_USE_COMPLEX)
136b0043f70SBarry 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);
13779c34000SHong Zhang #else
1384b019bd2SJed Brown   status = umfpack_UMF_wsolve(uflag,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
13979c34000SHong Zhang #endif
1408592901bSBarry Smith   umfpack_UMF_report_info(lu->Control, lu->Info);
1411677d0b8SKris Buschelman   if (status < 0) {
1428592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
143e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_wsolve failed");
1441677d0b8SKris Buschelman   }
1451677d0b8SKris Buschelman 
146d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(b,&ba);CHKERRQ(ierr);
147057dcbc9SJed Brown   ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr);
1484b019bd2SJed Brown   PetscFunctionReturn(0);
1494b019bd2SJed Brown }
1502a325a84SHong Zhang 
1514b019bd2SJed Brown #undef __FUNCT__
1524b019bd2SJed Brown #define __FUNCT__ "MatSolve_UMFPACK"
1534b019bd2SJed Brown static PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x)
1544b019bd2SJed Brown {
1554b019bd2SJed Brown   PetscErrorCode ierr;
1564b019bd2SJed Brown 
1574b019bd2SJed Brown   PetscFunctionBegin;
1584b019bd2SJed Brown   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
1594b019bd2SJed Brown   ierr = MatSolve_UMFPACK_Private(A,b,x,UMFPACK_Aat);CHKERRQ(ierr);
1604b019bd2SJed Brown   PetscFunctionReturn(0);
1614b019bd2SJed Brown }
1624b019bd2SJed Brown 
1634b019bd2SJed Brown #undef __FUNCT__
1644b019bd2SJed Brown #define __FUNCT__ "MatSolveTranspose_UMFPACK"
1654b019bd2SJed Brown static PetscErrorCode MatSolveTranspose_UMFPACK(Mat A,Vec b,Vec x)
1664b019bd2SJed Brown {
1674b019bd2SJed Brown   PetscErrorCode ierr;
1684b019bd2SJed Brown 
1694b019bd2SJed Brown   PetscFunctionBegin;
1704b019bd2SJed Brown   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
1714b019bd2SJed Brown   ierr = MatSolve_UMFPACK_Private(A,b,x,UMFPACK_A);CHKERRQ(ierr);
1721677d0b8SKris Buschelman   PetscFunctionReturn(0);
1731677d0b8SKris Buschelman }
1741677d0b8SKris Buschelman 
1751677d0b8SKris Buschelman #undef __FUNCT__
176f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_UMFPACK"
1774b019bd2SJed Brown static PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat F,Mat A,const MatFactorInfo *info)
1786849ba73SBarry Smith {
179719d5645SBarry Smith   Mat_UMFPACK    *lu = (Mat_UMFPACK*)(F)->spptr;
180ae26a541SJed Brown   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
181ae26a541SJed Brown   PetscInt       *ai = a->i,*aj=a->j,status;
182ae26a541SJed Brown   PetscScalar    *av = a->a;
1836849ba73SBarry Smith   PetscErrorCode ierr;
1841677d0b8SKris Buschelman 
1851677d0b8SKris Buschelman   PetscFunctionBegin;
1861677d0b8SKris Buschelman   /* numeric factorization of A' */
1871677d0b8SKris Buschelman   /* ----------------------------*/
1882d4e2982SHong Zhang 
1898592901bSBarry Smith   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) {
1908592901bSBarry Smith     umfpack_UMF_free_numeric(&lu->Numeric);
1918592901bSBarry Smith   }
1922d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1938592901bSBarry Smith   status = umfpack_UMF_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1942d4e2982SHong Zhang #else
1958592901bSBarry Smith   status = umfpack_UMF_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1968592901bSBarry Smith #endif
1979f42a82aSMatthew Knepley   if (status < 0) {
1988592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
199e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_numeric failed");
2009f42a82aSMatthew Knepley   }
2012d4e2982SHong Zhang   /* report numeric factorization of A' when Control[PRL] > 3 */
2028592901bSBarry Smith   (void) umfpack_UMF_report_numeric(lu->Numeric, lu->Control);
2031677d0b8SKris Buschelman 
204ae26a541SJed Brown   ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
205ae26a541SJed Brown   ierr = MatDestroy(&lu->A);CHKERRQ(ierr);
2062205254eSKarl Rupp 
207ae26a541SJed Brown   lu->A                  = A;
2082a325a84SHong Zhang   lu->flg                = SAME_NONZERO_PATTERN;
2092a325a84SHong Zhang   lu->CleanUpUMFPACK     = PETSC_TRUE;
2104b019bd2SJed Brown   F->ops->solve          = MatSolve_UMFPACK;
2114b019bd2SJed Brown   F->ops->solvetranspose = MatSolveTranspose_UMFPACK;
2121677d0b8SKris Buschelman   PetscFunctionReturn(0);
2131677d0b8SKris Buschelman }
2141677d0b8SKris Buschelman 
2151677d0b8SKris Buschelman /*
2161677d0b8SKris Buschelman    Note the r permutation is ignored
2171677d0b8SKris Buschelman */
2181677d0b8SKris Buschelman #undef __FUNCT__
219f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK"
2204b019bd2SJed Brown static PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
2216849ba73SBarry Smith {
222ae26a541SJed Brown   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
223719d5645SBarry Smith   Mat_UMFPACK    *lu = (Mat_UMFPACK*)(F->spptr);
224dfbe8321SBarry Smith   PetscErrorCode ierr;
225ae26a541SJed Brown   PetscInt       i,*ai = a->i,*aj = a->j,m=A->rmap->n,n=A->cmap->n;
226*989213a4SSatish Balay #if !defined(PETSC_USE_COMPLEX)
227ae26a541SJed Brown   PetscScalar    *av = a->a;
228*989213a4SSatish Balay #endif
2295d0c19d7SBarry Smith   const PetscInt *ra;
2308592901bSBarry Smith   PetscInt       status;
2311677d0b8SKris Buschelman 
2321677d0b8SKris Buschelman   PetscFunctionBegin;
233fcfd50ebSBarry Smith   if (lu->PetscMatOrdering) {
2340f6efc10SHong Zhang     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
235785e854fSJed Brown     ierr = PetscMalloc1(m,&lu->perm_c);CHKERRQ(ierr);
2362d4e2982SHong Zhang     /* we cannot simply memcpy on 64 bit archs */
2372d4e2982SHong Zhang     for (i = 0; i < m; i++) lu->perm_c[i] = ra[i];
2380f6efc10SHong Zhang     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
2391677d0b8SKris Buschelman   }
2401677d0b8SKris Buschelman 
2411677d0b8SKris Buschelman   /* print the control parameters */
2428592901bSBarry Smith   if (lu->Control[UMFPACK_PRL] > 1) umfpack_UMF_report_control(lu->Control);
2431677d0b8SKris Buschelman 
2441677d0b8SKris Buschelman   /* symbolic factorization of A' */
2451677d0b8SKris Buschelman   /* ---------------------------------------------------------------------- */
246fcfd50ebSBarry Smith   if (lu->PetscMatOrdering) { /* use Petsc row ordering */
2478592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
248ae26a541SJed Brown     status = umfpack_UMF_qsymbolic(n,m,ai,aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2492d4e2982SHong Zhang #else
250b0043f70SBarry Smith     status = umfpack_UMF_qsymbolic(n,m,ai,aj,NULL,NULL,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2518592901bSBarry Smith #endif
2522d4e2982SHong Zhang   } else { /* use Umfpack col ordering */
2538592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
254ae26a541SJed Brown     status = umfpack_UMF_symbolic(n,m,ai,aj,av,&lu->Symbolic,lu->Control,lu->Info);
2558592901bSBarry Smith #else
256ae26a541SJed Brown     status = umfpack_UMF_symbolic(n,m,ai,aj,NULL,NULL,&lu->Symbolic,lu->Control,lu->Info);
2578592901bSBarry Smith #endif
2582d4e2982SHong Zhang   }
2592d4e2982SHong Zhang   if (status < 0) {
2608592901bSBarry Smith     umfpack_UMF_report_info(lu->Control, lu->Info);
2618592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
262e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_symbolic failed");
2632d4e2982SHong Zhang   }
2642d4e2982SHong Zhang   /* report sumbolic factorization of A' when Control[PRL] > 3 */
2658592901bSBarry Smith   (void) umfpack_UMF_report_symbolic(lu->Symbolic, lu->Control);
2661677d0b8SKris Buschelman 
2671677d0b8SKris Buschelman   lu->flg                   = DIFFERENT_NONZERO_PATTERN;
268e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK        = PETSC_TRUE;
269719d5645SBarry Smith   (F)->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
2701677d0b8SKris Buschelman   PetscFunctionReturn(0);
2711677d0b8SKris Buschelman }
2721677d0b8SKris Buschelman 
2731677d0b8SKris Buschelman #undef __FUNCT__
274f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_UMFPACK"
2754b019bd2SJed Brown static PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer)
2766849ba73SBarry Smith {
277f0c56d0fSKris Buschelman   Mat_UMFPACK    *lu= (Mat_UMFPACK*)A->spptr;
278dfbe8321SBarry Smith   PetscErrorCode ierr;
279f0c56d0fSKris Buschelman 
2801677d0b8SKris Buschelman   PetscFunctionBegin;
2811677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
282f0c56d0fSKris Buschelman   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
2831677d0b8SKris Buschelman 
2841677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
2851677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
2861677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
2871677d0b8SKris Buschelman 
2881677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
2890f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr);
2901677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
2911677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
2920f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr);
2931677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
2940f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr);
2950f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr);
2961677d0b8SKris Buschelman 
2971677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
2981677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
2990f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr);
3000f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr);
3011677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
3020f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr);
3031677d0b8SKris Buschelman 
3041677d0b8SKris Buschelman   /* Control parameters used by solve */
3051677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
3061677d0b8SKris Buschelman 
3071677d0b8SKris Buschelman   /* mat ordering */
308fcfd50ebSBarry Smith   if (!lu->PetscMatOrdering) {
309fcfd50ebSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ORDERING]: %s (not using the PETSc ordering)\n",UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]]);CHKERRQ(ierr);
310fcfd50ebSBarry Smith   }
3111677d0b8SKris Buschelman   PetscFunctionReturn(0);
3121677d0b8SKris Buschelman }
3131677d0b8SKris Buschelman 
314d844382dSKris Buschelman #undef __FUNCT__
315f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_UMFPACK"
3164b019bd2SJed Brown static PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
3176849ba73SBarry Smith {
318dfbe8321SBarry Smith   PetscErrorCode    ierr;
319ace3abfcSBarry Smith   PetscBool         iascii;
320d844382dSKris Buschelman   PetscViewerFormat format;
321d844382dSKris Buschelman 
322d844382dSKris Buschelman   PetscFunctionBegin;
323b24902e0SBarry Smith   ierr = MatView_SeqAIJ(A,viewer);CHKERRQ(ierr);
324d844382dSKris Buschelman 
325251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
32632077d6dSBarry Smith   if (iascii) {
327d844382dSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
3282f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO) {
329f0c56d0fSKris Buschelman       ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
330d844382dSKris Buschelman     }
331d844382dSKris Buschelman   }
332d844382dSKris Buschelman   PetscFunctionReturn(0);
333d844382dSKris Buschelman }
334d844382dSKris Buschelman 
33535bd34faSBarry Smith #undef __FUNCT__
33635bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_umfpack"
33735bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_seqaij_umfpack(Mat A,const MatSolverPackage *type)
33835bd34faSBarry Smith {
33935bd34faSBarry Smith   PetscFunctionBegin;
3402692d6eeSBarry Smith   *type = MATSOLVERUMFPACK;
34135bd34faSBarry Smith   PetscFunctionReturn(0);
34235bd34faSBarry Smith }
34335bd34faSBarry Smith 
34435bd34faSBarry Smith 
3452f71e704SKris Buschelman /*MC
3462692d6eeSBarry Smith   MATSOLVERUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
3472f71e704SKris Buschelman   via the external package UMFPACK.
3482f71e704SKris Buschelman 
349f02a4714SShri   ./configure --download-suitesparse to install PETSc to use UMFPACK
3502f71e704SKris Buschelman 
3512f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
3522f71e704SKris Buschelman   which correspond to the options database keys below.
3532f71e704SKris Buschelman 
3542f71e704SKris Buschelman   Options Database Keys:
355ba41fbb6SBarry Smith + -mat_umfpack_ordering                - CHOLMOD, AMD, GIVEN, METIS, BEST, NONE
356ba41fbb6SBarry Smith . -mat_umfpack_prl                     - UMFPACK print level: Control[UMFPACK_PRL]
357e08999f5SMatthew G Knepley . -mat_umfpack_strategy <AUTO>         - (choose one of) AUTO UNSYMMETRIC SYMMETRIC 2BY2
3582f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c>     - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
359e08999f5SMatthew G Knepley . -mat_umfpack_dense_row <0.2>         - Control[UMFPACK_DENSE_ROW]
360e08999f5SMatthew G Knepley . -mat_umfpack_amd_dense <10>          - Control[UMFPACK_AMD_DENSE]
3612f71e704SKris Buschelman . -mat_umfpack_block_size <bs>         - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
362e08999f5SMatthew G Knepley . -mat_umfpack_2by2_tolerance <0.01>   - Control[UMFPACK_2BY2_TOLERANCE]
363e08999f5SMatthew G Knepley . -mat_umfpack_fixq <0>                - Control[UMFPACK_FIXQ]
364e08999f5SMatthew G Knepley . -mat_umfpack_aggressive <1>          - Control[UMFPACK_AGGRESSIVE]
3652f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
366e08999f5SMatthew G Knepley . -mat_umfpack_sym_pivot_tolerance <0.001> - Control[UMFPACK_SYM_PIVOT_TOLERANCE]
367e08999f5SMatthew G Knepley . -mat_umfpack_scale <NONE>           - (choose one of) NONE SUM MAX
3682f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta>      - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
369e08999f5SMatthew G Knepley . -mat_umfpack_droptol <0>            - Control[UMFPACK_DROPTOL]
3702f71e704SKris Buschelman - -mat_umfpack_irstep <maxit>          - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
3712f71e704SKris Buschelman 
3722f71e704SKris Buschelman    Level: beginner
373a364b7d2SBarry Smith 
374a364b7d2SBarry Smith    Note: UMFPACK is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html
3752f71e704SKris Buschelman 
376d45987f3SHong Zhang .seealso: PCLU, MATSOLVERSUPERLU, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage
3772f71e704SKris Buschelman M*/
378b2573a8aSBarry Smith 
3793519fcdcSHong Zhang #undef __FUNCT__
3803519fcdcSHong Zhang #define __FUNCT__ "MatGetFactor_seqaij_umfpack"
3818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_umfpack(Mat A,MatFactorType ftype,Mat *F)
3823519fcdcSHong Zhang {
3833519fcdcSHong Zhang   Mat            B;
3843519fcdcSHong Zhang   Mat_UMFPACK    *lu;
3853519fcdcSHong Zhang   PetscErrorCode ierr;
3868592901bSBarry Smith   PetscInt       m=A->rmap->n,n=A->cmap->n,idx;
3872f71e704SKris Buschelman 
3882205254eSKarl Rupp   const char *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC"};
3892205254eSKarl Rupp   const char *scale[]   ={"NONE","SUM","MAX"};
390ace3abfcSBarry Smith   PetscBool  flg;
3912d4e2982SHong Zhang 
3923519fcdcSHong Zhang   PetscFunctionBegin;
3933519fcdcSHong Zhang   /* Create the factorization matrix F */
394ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
3953519fcdcSHong Zhang   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
3963519fcdcSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
3970298fd71SBarry Smith   ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
398b00a9115SJed Brown   ierr = PetscNewLog(B,&lu);CHKERRQ(ierr);
3992205254eSKarl Rupp 
4003519fcdcSHong Zhang   B->spptr                 = lu;
4013519fcdcSHong Zhang   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
4023519fcdcSHong Zhang   B->ops->destroy          = MatDestroy_UMFPACK;
4033519fcdcSHong Zhang   B->ops->view             = MatView_UMFPACK;
4042205254eSKarl Rupp 
405bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_umfpack);CHKERRQ(ierr);
4062205254eSKarl Rupp 
407d5f3da31SBarry Smith   B->factortype   = MAT_FACTOR_LU;
4083519fcdcSHong Zhang   B->assembled    = PETSC_TRUE;           /* required by -ksp_view */
4093519fcdcSHong Zhang   B->preallocated = PETSC_TRUE;
4103519fcdcSHong Zhang 
4113519fcdcSHong Zhang   /* initializations */
4123519fcdcSHong Zhang   /* ------------------------------------------------*/
4133519fcdcSHong Zhang   /* get the default control parameters */
4148592901bSBarry Smith   umfpack_UMF_defaults(lu->Control);
4150298fd71SBarry Smith   lu->perm_c                  = NULL; /* use defaul UMFPACK col permutation */
4163519fcdcSHong Zhang   lu->Control[UMFPACK_IRSTEP] = 0;          /* max num of iterative refinement steps to attempt */
4173519fcdcSHong Zhang 
418ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
4193519fcdcSHong Zhang   /* Control parameters used by reporting routiones */
4200298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],NULL);CHKERRQ(ierr);
4213519fcdcSHong Zhang 
4223519fcdcSHong Zhang   /* Control parameters for symbolic factorization */
423fcfd50ebSBarry Smith   ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,3,strategy[0],&idx,&flg);CHKERRQ(ierr);
4243519fcdcSHong Zhang   if (flg) {
4253519fcdcSHong Zhang     switch (idx) {
4263519fcdcSHong Zhang     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
4273519fcdcSHong Zhang     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
4283519fcdcSHong Zhang     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
4293519fcdcSHong Zhang     }
4303519fcdcSHong Zhang   }
4318caf3d72SBarry 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);
432fcfd50ebSBarry Smith   if (flg) lu->Control[UMFPACK_ORDERING] = (int)idx;
4330298fd71SBarry 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);
4340298fd71SBarry 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);
4350298fd71SBarry 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);
4360298fd71SBarry 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);
4370298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],NULL);CHKERRQ(ierr);
4380298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],NULL);CHKERRQ(ierr);
4393519fcdcSHong Zhang 
4403519fcdcSHong Zhang   /* Control parameters used by numeric factorization */
4410298fd71SBarry 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);
4420298fd71SBarry 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);
4433519fcdcSHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
4443519fcdcSHong Zhang   if (flg) {
4453519fcdcSHong Zhang     switch (idx) {
4463519fcdcSHong Zhang     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
4473519fcdcSHong Zhang     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
4483519fcdcSHong Zhang     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
4493519fcdcSHong Zhang     }
4503519fcdcSHong Zhang   }
4510298fd71SBarry 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);
4520298fd71SBarry 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);
4530298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],NULL);CHKERRQ(ierr);
4543519fcdcSHong Zhang 
4553519fcdcSHong Zhang   /* Control parameters used by solve */
4560298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],NULL);CHKERRQ(ierr);
4573519fcdcSHong Zhang 
4583519fcdcSHong Zhang   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */
4590298fd71SBarry Smith   ierr = PetscOptionsHasName(NULL,"-pc_factor_mat_ordering_type",&lu->PetscMatOrdering);CHKERRQ(ierr);
4603519fcdcSHong Zhang   PetscOptionsEnd();
4613519fcdcSHong Zhang   *F = B;
4623519fcdcSHong Zhang   PetscFunctionReturn(0);
4633519fcdcSHong Zhang }
464b2573a8aSBarry Smith 
46542c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat,MatFactorType,Mat*);
46642c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat,MatFactorType,Mat*);
46742c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat,MatFactorType,Mat*);
4682d4e2982SHong Zhang 
46942c9c57cSBarry Smith #undef __FUNCT__
47042c9c57cSBarry Smith #define __FUNCT__ "MatSolverPackageRegister_SuiteSparse"
47129b38603SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_SuiteSparse(void)
47242c9c57cSBarry Smith {
47342c9c57cSBarry Smith   PetscErrorCode ierr;
47442c9c57cSBarry Smith 
47542c9c57cSBarry Smith   PetscFunctionBegin;
47642c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERUMFPACK,MATSEQAIJ,      MAT_FACTOR_LU,MatGetFactor_seqaij_umfpack);CHKERRQ(ierr);
47742c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERCHOLMOD,MATSEQAIJ,      MAT_FACTOR_CHOLESKY,MatGetFactor_seqaij_cholmod);CHKERRQ(ierr);
47842c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERCHOLMOD,MATSEQSBAIJ,      MAT_FACTOR_CHOLESKY,MatGetFactor_seqsbaij_cholmod);CHKERRQ(ierr);
47942c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERKLU,MATSEQAIJ,          MAT_FACTOR_LU,MatGetFactor_seqaij_klu);CHKERRQ(ierr);
48042c9c57cSBarry Smith   PetscFunctionReturn(0);
48142c9c57cSBarry Smith }
482