xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision dd39110ba674af879e763a5939b18b2b35b4c87f)
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;
871677d0b8SKris Buschelman 
881677d0b8SKris Buschelman   /* Flag to clean up UMFPACK objects during Destroy */
89ace3abfcSBarry Smith   PetscBool CleanUpUMFPACK;
90f0c56d0fSKris Buschelman } Mat_UMFPACK;
91f0c56d0fSKris Buschelman 
924b019bd2SJed Brown static PetscErrorCode MatDestroy_UMFPACK(Mat A)
93dfbe8321SBarry Smith {
94b9c12af5SBarry Smith   Mat_UMFPACK    *lu=(Mat_UMFPACK*)A->data;
951677d0b8SKris Buschelman 
961677d0b8SKris Buschelman   PetscFunctionBegin;
97b9c12af5SBarry Smith   if (lu->CleanUpUMFPACK) {
988592901bSBarry Smith     umfpack_UMF_free_symbolic(&lu->Symbolic);
998592901bSBarry Smith     umfpack_UMF_free_numeric(&lu->Numeric);
1009566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->Wi));
1019566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->W));
1029566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->perm_c));
1031677d0b8SKris Buschelman   }
1049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lu->A));
1059566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
1061677d0b8SKris Buschelman   PetscFunctionReturn(0);
1071677d0b8SKris Buschelman }
1081677d0b8SKris Buschelman 
1094b019bd2SJed Brown static PetscErrorCode MatSolve_UMFPACK_Private(Mat A,Vec b,Vec x,int uflag)
1106849ba73SBarry Smith {
111b9c12af5SBarry Smith   Mat_UMFPACK       *lu = (Mat_UMFPACK*)A->data;
112ae26a541SJed Brown   Mat_SeqAIJ        *a  = (Mat_SeqAIJ*)lu->A->data;
113d9ca1df4SBarry Smith   PetscScalar       *av = a->a,*xa;
114d9ca1df4SBarry Smith   const PetscScalar *ba;
115ae26a541SJed Brown   PetscInt          *ai = a->i,*aj = a->j,status;
116fcf6e9abSJed Brown   static PetscBool  cite = PETSC_FALSE;
1171677d0b8SKris Buschelman 
1181677d0b8SKris Buschelman   PetscFunctionBegin;
119a643c07aSBarry Smith   if (!A->rmap->n) PetscFunctionReturn(0);
1209566063dSJacob Faibussowitsch   PetscCall(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));
1212d4e2982SHong Zhang   /* solve Ax = b by umfpack_*_wsolve */
1221677d0b8SKris Buschelman   /* ----------------------------------*/
1232d4e2982SHong Zhang 
124ae26a541SJed Brown   if (!lu->Wi) {  /* first time, allocate working space for wsolve */
1259566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(A->rmap->n,&lu->Wi));
1269566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5*A->rmap->n,&lu->W));
127ae26a541SJed Brown   }
128ae26a541SJed Brown 
1299566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(b,&ba));
1309566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x,&xa));
13179c34000SHong Zhang #if defined(PETSC_USE_COMPLEX)
132b0043f70SBarry 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);
13379c34000SHong Zhang #else
1344b019bd2SJed Brown   status = umfpack_UMF_wsolve(uflag,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
13579c34000SHong Zhang #endif
1368592901bSBarry Smith   umfpack_UMF_report_info(lu->Control, lu->Info);
1371677d0b8SKris Buschelman   if (status < 0) {
1388592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
139e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_wsolve failed");
1401677d0b8SKris Buschelman   }
1411677d0b8SKris Buschelman 
1429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(b,&ba));
1439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x,&xa));
1444b019bd2SJed Brown   PetscFunctionReturn(0);
1454b019bd2SJed Brown }
1462a325a84SHong Zhang 
1474b019bd2SJed Brown static PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x)
1484b019bd2SJed Brown {
1494b019bd2SJed Brown   PetscFunctionBegin;
1504b019bd2SJed Brown   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
1519566063dSJacob Faibussowitsch   PetscCall(MatSolve_UMFPACK_Private(A,b,x,UMFPACK_Aat));
1524b019bd2SJed Brown   PetscFunctionReturn(0);
1534b019bd2SJed Brown }
1544b019bd2SJed Brown 
1554b019bd2SJed Brown static PetscErrorCode MatSolveTranspose_UMFPACK(Mat A,Vec b,Vec x)
1564b019bd2SJed Brown {
1574b019bd2SJed Brown   PetscFunctionBegin;
1584b019bd2SJed Brown   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
1599566063dSJacob Faibussowitsch   PetscCall(MatSolve_UMFPACK_Private(A,b,x,UMFPACK_A));
1601677d0b8SKris Buschelman   PetscFunctionReturn(0);
1611677d0b8SKris Buschelman }
1621677d0b8SKris Buschelman 
1634b019bd2SJed Brown static PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat F,Mat A,const MatFactorInfo *info)
1646849ba73SBarry Smith {
165b9c12af5SBarry Smith   Mat_UMFPACK    *lu = (Mat_UMFPACK*)(F)->data;
166ae26a541SJed Brown   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
167ae26a541SJed Brown   PetscInt       *ai = a->i,*aj=a->j,status;
168ae26a541SJed Brown   PetscScalar    *av = a->a;
1691677d0b8SKris Buschelman 
1701677d0b8SKris Buschelman   PetscFunctionBegin;
171a643c07aSBarry Smith   if (!A->rmap->n) PetscFunctionReturn(0);
1721677d0b8SKris Buschelman   /* numeric factorization of A' */
1731677d0b8SKris Buschelman   /* ----------------------------*/
1742d4e2982SHong Zhang 
1758592901bSBarry Smith   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) {
1768592901bSBarry Smith     umfpack_UMF_free_numeric(&lu->Numeric);
1778592901bSBarry Smith   }
1782d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1798592901bSBarry Smith   status = umfpack_UMF_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1802d4e2982SHong Zhang #else
1818592901bSBarry Smith   status = umfpack_UMF_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1828592901bSBarry Smith #endif
1839f42a82aSMatthew Knepley   if (status < 0) {
1848592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
185e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_numeric failed");
1869f42a82aSMatthew Knepley   }
1872d4e2982SHong Zhang   /* report numeric factorization of A' when Control[PRL] > 3 */
1888592901bSBarry Smith   (void) umfpack_UMF_report_numeric(lu->Numeric, lu->Control);
1891677d0b8SKris Buschelman 
1909566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)A));
1919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lu->A));
1922205254eSKarl Rupp 
193ae26a541SJed Brown   lu->A                  = A;
1942a325a84SHong Zhang   lu->flg                = SAME_NONZERO_PATTERN;
1952a325a84SHong Zhang   lu->CleanUpUMFPACK     = PETSC_TRUE;
1964b019bd2SJed Brown   F->ops->solve          = MatSolve_UMFPACK;
1974b019bd2SJed Brown   F->ops->solvetranspose = MatSolveTranspose_UMFPACK;
1981677d0b8SKris Buschelman   PetscFunctionReturn(0);
1991677d0b8SKris Buschelman }
2001677d0b8SKris Buschelman 
2014b019bd2SJed Brown static PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
2026849ba73SBarry Smith {
203ae26a541SJed Brown   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
204b9c12af5SBarry Smith   Mat_UMFPACK    *lu = (Mat_UMFPACK*)(F->data);
205ae26a541SJed Brown   PetscInt       i,*ai = a->i,*aj = a->j,m=A->rmap->n,n=A->cmap->n;
206989213a4SSatish Balay #if !defined(PETSC_USE_COMPLEX)
207ae26a541SJed Brown   PetscScalar    *av = a->a;
208989213a4SSatish Balay #endif
2095d0c19d7SBarry Smith   const PetscInt *ra;
2108592901bSBarry Smith   PetscInt       status;
2111677d0b8SKris Buschelman 
2121677d0b8SKris Buschelman   PetscFunctionBegin;
213a643c07aSBarry Smith   (F)->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
214a643c07aSBarry Smith   if (!n) PetscFunctionReturn(0);
2152c7c0729SBarry Smith   if (r) {
2169566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(r,&ra));
2179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m,&lu->perm_c));
2182d4e2982SHong Zhang     /* we cannot simply memcpy on 64 bit archs */
2192d4e2982SHong Zhang     for (i = 0; i < m; i++) lu->perm_c[i] = ra[i];
2209566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(r,&ra));
2211677d0b8SKris Buschelman   }
2221677d0b8SKris Buschelman 
2231677d0b8SKris Buschelman   /* print the control parameters */
2248592901bSBarry Smith   if (lu->Control[UMFPACK_PRL] > 1) umfpack_UMF_report_control(lu->Control);
2251677d0b8SKris Buschelman 
2261677d0b8SKris Buschelman   /* symbolic factorization of A' */
2271677d0b8SKris Buschelman   /* ---------------------------------------------------------------------- */
2282c7c0729SBarry Smith   if (r) { /* use Petsc row ordering */
2298592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
230ae26a541SJed Brown     status = umfpack_UMF_qsymbolic(n,m,ai,aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2312d4e2982SHong Zhang #else
232b0043f70SBarry Smith     status = umfpack_UMF_qsymbolic(n,m,ai,aj,NULL,NULL,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2338592901bSBarry Smith #endif
2342d4e2982SHong Zhang   } else { /* use Umfpack col ordering */
2358592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
236ae26a541SJed Brown     status = umfpack_UMF_symbolic(n,m,ai,aj,av,&lu->Symbolic,lu->Control,lu->Info);
2378592901bSBarry Smith #else
238ae26a541SJed Brown     status = umfpack_UMF_symbolic(n,m,ai,aj,NULL,NULL,&lu->Symbolic,lu->Control,lu->Info);
2398592901bSBarry Smith #endif
2402d4e2982SHong Zhang   }
2412d4e2982SHong Zhang   if (status < 0) {
2428592901bSBarry Smith     umfpack_UMF_report_info(lu->Control, lu->Info);
2438592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
244e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_symbolic failed");
2452d4e2982SHong Zhang   }
2462d4e2982SHong Zhang   /* report sumbolic factorization of A' when Control[PRL] > 3 */
2478592901bSBarry Smith   (void) umfpack_UMF_report_symbolic(lu->Symbolic, lu->Control);
2481677d0b8SKris Buschelman 
2491677d0b8SKris Buschelman   lu->flg                   = DIFFERENT_NONZERO_PATTERN;
250e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK        = PETSC_TRUE;
2511677d0b8SKris Buschelman   PetscFunctionReturn(0);
2521677d0b8SKris Buschelman }
2531677d0b8SKris Buschelman 
254860c79edSBarry Smith static PetscErrorCode MatView_Info_UMFPACK(Mat A,PetscViewer viewer)
2556849ba73SBarry Smith {
256b9c12af5SBarry Smith   Mat_UMFPACK    *lu= (Mat_UMFPACK*)A->data;
257f0c56d0fSKris Buschelman 
2581677d0b8SKris Buschelman   PetscFunctionBegin;
2591677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
260f0c56d0fSKris Buschelman   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
2611677d0b8SKris Buschelman 
2629566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n"));
2631677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
2649566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]));
2651677d0b8SKris Buschelman 
2661677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
2679566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]));
2689566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]));
2699566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]));
2709566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]));
2719566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]));
2729566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]));
2739566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]));
2741677d0b8SKris Buschelman 
2751677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
2769566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]));
2779566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]));
2789566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]));
2799566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]));
2809566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]));
2811677d0b8SKris Buschelman 
2821677d0b8SKris Buschelman   /* Control parameters used by solve */
2839566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]));
2841677d0b8SKris Buschelman 
2851677d0b8SKris Buschelman   /* mat ordering */
2862c7c0729SBarry Smith   if (!lu->perm_c) {
2879566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ORDERING]: %s (not using the PETSc ordering)\n",UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]]));
288fcfd50ebSBarry Smith   }
2891677d0b8SKris Buschelman   PetscFunctionReturn(0);
2901677d0b8SKris Buschelman }
2911677d0b8SKris Buschelman 
2924b019bd2SJed Brown static PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
2936849ba73SBarry Smith {
294ace3abfcSBarry Smith   PetscBool         iascii;
295d844382dSKris Buschelman   PetscViewerFormat format;
296d844382dSKris Buschelman 
297d844382dSKris Buschelman   PetscFunctionBegin;
2989566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
29932077d6dSBarry Smith   if (iascii) {
3009566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer,&format));
3012f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO) {
3029566063dSJacob Faibussowitsch       PetscCall(MatView_Info_UMFPACK(A,viewer));
303d844382dSKris Buschelman     }
304d844382dSKris Buschelman   }
305d844382dSKris Buschelman   PetscFunctionReturn(0);
306d844382dSKris Buschelman }
307d844382dSKris Buschelman 
308ea799195SBarry Smith PetscErrorCode MatFactorGetSolverType_seqaij_umfpack(Mat A,MatSolverType *type)
30935bd34faSBarry Smith {
31035bd34faSBarry Smith   PetscFunctionBegin;
3112692d6eeSBarry Smith   *type = MATSOLVERUMFPACK;
31235bd34faSBarry Smith   PetscFunctionReturn(0);
31335bd34faSBarry Smith }
31435bd34faSBarry Smith 
3152f71e704SKris Buschelman /*MC
3162692d6eeSBarry Smith   MATSOLVERUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
3172f71e704SKris Buschelman   via the external package UMFPACK.
3182f71e704SKris Buschelman 
319c2b89b5dSBarry Smith   Use ./configure --download-suitesparse to install PETSc to use UMFPACK
320c2b89b5dSBarry Smith 
3213ca39a21SBarry Smith   Use -pc_type lu -pc_factor_mat_solver_type umfpack to use this direct solver
3222f71e704SKris Buschelman 
3232f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
3242f71e704SKris Buschelman   which correspond to the options database keys below.
3252f71e704SKris Buschelman 
3262f71e704SKris Buschelman   Options Database Keys:
327ba41fbb6SBarry Smith + -mat_umfpack_ordering                - CHOLMOD, AMD, GIVEN, METIS, BEST, NONE
328ba41fbb6SBarry Smith . -mat_umfpack_prl                     - UMFPACK print level: Control[UMFPACK_PRL]
329e08999f5SMatthew G Knepley . -mat_umfpack_strategy <AUTO>         - (choose one of) AUTO UNSYMMETRIC SYMMETRIC 2BY2
3302f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c>     - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
331e08999f5SMatthew G Knepley . -mat_umfpack_dense_row <0.2>         - Control[UMFPACK_DENSE_ROW]
332e08999f5SMatthew G Knepley . -mat_umfpack_amd_dense <10>          - Control[UMFPACK_AMD_DENSE]
3332f71e704SKris Buschelman . -mat_umfpack_block_size <bs>         - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
334e08999f5SMatthew G Knepley . -mat_umfpack_2by2_tolerance <0.01>   - Control[UMFPACK_2BY2_TOLERANCE]
335e08999f5SMatthew G Knepley . -mat_umfpack_fixq <0>                - Control[UMFPACK_FIXQ]
336e08999f5SMatthew G Knepley . -mat_umfpack_aggressive <1>          - Control[UMFPACK_AGGRESSIVE]
3372f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
338e08999f5SMatthew G Knepley . -mat_umfpack_sym_pivot_tolerance <0.001> - Control[UMFPACK_SYM_PIVOT_TOLERANCE]
339e08999f5SMatthew G Knepley . -mat_umfpack_scale <NONE>           - (choose one of) NONE SUM MAX
3402f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta>      - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
341e08999f5SMatthew G Knepley . -mat_umfpack_droptol <0>            - Control[UMFPACK_DROPTOL]
3422f71e704SKris Buschelman - -mat_umfpack_irstep <maxit>          - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
3432f71e704SKris Buschelman 
3442f71e704SKris Buschelman    Level: beginner
345a364b7d2SBarry Smith 
346a364b7d2SBarry Smith    Note: UMFPACK is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html
3472f71e704SKris Buschelman 
3483ca39a21SBarry Smith .seealso: PCLU, MATSOLVERSUPERLU, MATSOLVERMUMPS, PCFactorSetMatSolverType(), MatSolverType
3492f71e704SKris Buschelman M*/
350b2573a8aSBarry Smith 
3518cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_umfpack(Mat A,MatFactorType ftype,Mat *F)
3523519fcdcSHong Zhang {
3533519fcdcSHong Zhang   Mat            B;
3543519fcdcSHong Zhang   Mat_UMFPACK    *lu;
3553519fcdcSHong Zhang   PetscErrorCode ierr;
3568592901bSBarry Smith   PetscInt       m=A->rmap->n,n=A->cmap->n,idx;
3572205254eSKarl Rupp   const char     *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC"};
3582205254eSKarl Rupp   const char     *scale[]   ={"NONE","SUM","MAX"};
359ace3abfcSBarry Smith   PetscBool      flg;
3602d4e2982SHong Zhang 
3613519fcdcSHong Zhang   PetscFunctionBegin;
3623519fcdcSHong Zhang   /* Create the factorization matrix F */
3639566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B));
3649566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n));
3659566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("umfpack",&((PetscObject)B)->type_name));
3669566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
367b9c12af5SBarry Smith 
3689566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(B,&lu));
3692205254eSKarl Rupp 
370b9c12af5SBarry Smith   B->data                   = lu;
371b9c12af5SBarry Smith   B->ops->getinfo          = MatGetInfo_External;
3723519fcdcSHong Zhang   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
3733519fcdcSHong Zhang   B->ops->destroy          = MatDestroy_UMFPACK;
3743519fcdcSHong Zhang   B->ops->view             = MatView_UMFPACK;
3751aef8b4cSStefano Zampini   B->ops->matsolve         = NULL;
3762205254eSKarl Rupp 
3779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_umfpack));
3782205254eSKarl Rupp 
379d5f3da31SBarry Smith   B->factortype   = MAT_FACTOR_LU;
3803519fcdcSHong Zhang   B->assembled    = PETSC_TRUE;           /* required by -ksp_view */
3813519fcdcSHong Zhang   B->preallocated = PETSC_TRUE;
3823519fcdcSHong Zhang 
3839566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
3849566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERUMFPACK,&B->solvertype));
385f73b0415SBarry Smith   B->canuseordering = PETSC_TRUE;
3869566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&B->preferredordering[MAT_FACTOR_LU]));
38700c67f3bSHong Zhang 
3883519fcdcSHong Zhang   /* initializations */
3893519fcdcSHong Zhang   /* ------------------------------------------------*/
3903519fcdcSHong Zhang   /* get the default control parameters */
3918592901bSBarry Smith   umfpack_UMF_defaults(lu->Control);
3920298fd71SBarry Smith   lu->perm_c                  = NULL; /* use defaul UMFPACK col permutation */
3933519fcdcSHong Zhang   lu->Control[UMFPACK_IRSTEP] = 0;          /* max num of iterative refinement steps to attempt */
3943519fcdcSHong Zhang 
3959566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"UMFPACK Options","Mat");PetscCall(ierr);
3963519fcdcSHong Zhang   /* Control parameters used by reporting routiones */
3979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],NULL));
3983519fcdcSHong Zhang 
3993519fcdcSHong Zhang   /* Control parameters for symbolic factorization */
4009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,3,strategy[0],&idx,&flg));
4013519fcdcSHong Zhang   if (flg) {
4023519fcdcSHong Zhang     switch (idx) {
4033519fcdcSHong Zhang     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
4043519fcdcSHong Zhang     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
4053519fcdcSHong Zhang     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
4063519fcdcSHong Zhang     }
4073519fcdcSHong Zhang   }
408*dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-mat_umfpack_ordering","Internal ordering method","None",UmfpackOrderingTypes,PETSC_STATIC_ARRAY_LENGTH(UmfpackOrderingTypes),UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]],&idx,&flg));
409fcfd50ebSBarry Smith   if (flg) lu->Control[UMFPACK_ORDERING] = (int)idx;
4109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_umfpack_dense_col","Control[UMFPACK_DENSE_COL]","None",lu->Control[UMFPACK_DENSE_COL],&lu->Control[UMFPACK_DENSE_COL],NULL));
4119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_umfpack_dense_row","Control[UMFPACK_DENSE_ROW]","None",lu->Control[UMFPACK_DENSE_ROW],&lu->Control[UMFPACK_DENSE_ROW],NULL));
4129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_umfpack_amd_dense","Control[UMFPACK_AMD_DENSE]","None",lu->Control[UMFPACK_AMD_DENSE],&lu->Control[UMFPACK_AMD_DENSE],NULL));
4139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_umfpack_block_size","Control[UMFPACK_BLOCK_SIZE]","None",lu->Control[UMFPACK_BLOCK_SIZE],&lu->Control[UMFPACK_BLOCK_SIZE],NULL));
4149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],NULL));
4159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],NULL));
4163519fcdcSHong Zhang 
4173519fcdcSHong Zhang   /* Control parameters used by numeric factorization */
4189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_umfpack_pivot_tolerance","Control[UMFPACK_PIVOT_TOLERANCE]","None",lu->Control[UMFPACK_PIVOT_TOLERANCE],&lu->Control[UMFPACK_PIVOT_TOLERANCE],NULL));
4199566063dSJacob Faibussowitsch   PetscCall(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));
4209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg));
4213519fcdcSHong Zhang   if (flg) {
4223519fcdcSHong Zhang     switch (idx) {
4233519fcdcSHong Zhang     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
4243519fcdcSHong Zhang     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
4253519fcdcSHong Zhang     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
4263519fcdcSHong Zhang     }
4273519fcdcSHong Zhang   }
4289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_umfpack_alloc_init","Control[UMFPACK_ALLOC_INIT]","None",lu->Control[UMFPACK_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],NULL));
4299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_umfpack_front_alloc_init","Control[UMFPACK_FRONT_ALLOC_INIT]","None",lu->Control[UMFPACK_FRONT_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],NULL));
4309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],NULL));
4313519fcdcSHong Zhang 
4323519fcdcSHong Zhang   /* Control parameters used by solve */
4339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],NULL));
4349566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
4353519fcdcSHong Zhang   *F = B;
4363519fcdcSHong Zhang   PetscFunctionReturn(0);
4373519fcdcSHong Zhang }
438b2573a8aSBarry Smith 
439db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat,MatFactorType,Mat*);
440db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat,MatFactorType,Mat*);
441db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat,MatFactorType,Mat*);
442a2fc1e05SToby Isaac PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_spqr(Mat,MatFactorType,Mat*);
4432d4e2982SHong Zhang 
4443ca39a21SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_SuiteSparse(void)
44542c9c57cSBarry Smith {
44642c9c57cSBarry Smith   PetscFunctionBegin;
4479566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERUMFPACK,MATSEQAIJ,  MAT_FACTOR_LU,MatGetFactor_seqaij_umfpack));
4489566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERCHOLMOD,MATSEQAIJ,  MAT_FACTOR_CHOLESKY,MatGetFactor_seqaij_cholmod));
4499566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERCHOLMOD,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_seqsbaij_cholmod));
4509566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERKLU,MATSEQAIJ,      MAT_FACTOR_LU,MatGetFactor_seqaij_klu));
4519566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSPQR,MATSEQAIJ,     MAT_FACTOR_QR,MatGetFactor_seqaij_spqr));
45265f45395SPierre Jolivet   if (!PetscDefined(USE_COMPLEX)) {
4539566063dSJacob Faibussowitsch     PetscCall(MatSolverTypeRegister(MATSOLVERSPQR,MATNORMAL,   MAT_FACTOR_QR,MatGetFactor_seqaij_spqr));
45465f45395SPierre Jolivet   }
4559566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSPQR,MATNORMALHERMITIAN, MAT_FACTOR_QR,MatGetFactor_seqaij_spqr));
45642c9c57cSBarry Smith   PetscFunctionReturn(0);
45742c9c57cSBarry Smith }
458