xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision da0e1c470cd0cbbbdc182335ab9ba482b734ddca)
11677d0b8SKris Buschelman 
21677d0b8SKris Buschelman /*
32d4e2982SHong Zhang    Provides an interface to the UMFPACKv5.1 sparse solver
42d4e2982SHong Zhang 
58592901bSBarry Smith    When build with PETSC_USE_64BIT_INDICES this will use UF_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
88592901bSBarry Smith    that UMFPACK UL_Long version MUST be built with 64 bit integers when used.
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 
748592901bSBarry Smith 
758592901bSBarry Smith #define UF_long long long
768592901bSBarry Smith #define UF_long_max LONG_LONG_MAX
778592901bSBarry Smith #define UF_long_id "%lld"
788592901bSBarry Smith 
791677d0b8SKris Buschelman EXTERN_C_BEGIN
80c6db04a5SJed Brown #include <umfpack.h>
811677d0b8SKris Buschelman EXTERN_C_END
821677d0b8SKris Buschelman 
83fcfd50ebSBarry Smith static const char *const UmfpackOrderingTypes[] = {"CHOLMOD","AMD","GIVEN","METIS","BEST","NONE","USER","UmfpackOrderingTypes","UMFPACK_ORDERING_",0};
84fcfd50ebSBarry Smith 
851677d0b8SKris Buschelman typedef struct {
861677d0b8SKris Buschelman   void         *Symbolic, *Numeric;
871677d0b8SKris Buschelman   double       Info[UMFPACK_INFO], Control[UMFPACK_CONTROL],*W;
888592901bSBarry Smith   PetscInt      *Wi,*ai,*aj,*perm_c;
891677d0b8SKris Buschelman   PetscScalar  *av;
901677d0b8SKris Buschelman   MatStructure flg;
91fcfd50ebSBarry Smith   PetscBool    PetscMatOrdering;
921677d0b8SKris Buschelman 
931677d0b8SKris Buschelman   /* Flag to clean up UMFPACK objects during Destroy */
94ace3abfcSBarry Smith   PetscBool  CleanUpUMFPACK;
95f0c56d0fSKris Buschelman } Mat_UMFPACK;
96f0c56d0fSKris Buschelman 
971677d0b8SKris Buschelman #undef __FUNCT__
98f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_UMFPACK"
994b019bd2SJed Brown static PetscErrorCode MatDestroy_UMFPACK(Mat A)
100dfbe8321SBarry Smith {
101dfbe8321SBarry Smith   PetscErrorCode ierr;
102f0c56d0fSKris Buschelman   Mat_UMFPACK    *lu=(Mat_UMFPACK*)A->spptr;
1031677d0b8SKris Buschelman 
1041677d0b8SKris Buschelman   PetscFunctionBegin;
105bf0cc555SLisandro Dalcin   if (lu && lu->CleanUpUMFPACK) {
1068592901bSBarry Smith     umfpack_UMF_free_symbolic(&lu->Symbolic);
1078592901bSBarry Smith     umfpack_UMF_free_numeric(&lu->Numeric);
1081677d0b8SKris Buschelman     ierr = PetscFree(lu->Wi);CHKERRQ(ierr);
1091677d0b8SKris Buschelman     ierr = PetscFree(lu->W);CHKERRQ(ierr);
110fcfd50ebSBarry Smith     if (lu->PetscMatOrdering) {
1111677d0b8SKris Buschelman       ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
1121677d0b8SKris Buschelman     }
1131677d0b8SKris Buschelman   }
114*da0e1c47SShri Abhyankar   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
115b24902e0SBarry Smith   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
1161677d0b8SKris Buschelman   PetscFunctionReturn(0);
1171677d0b8SKris Buschelman }
1181677d0b8SKris Buschelman 
1191677d0b8SKris Buschelman #undef __FUNCT__
1204b019bd2SJed Brown #define __FUNCT__ "MatSolve_UMFPACK_Private"
1214b019bd2SJed Brown static PetscErrorCode MatSolve_UMFPACK_Private(Mat A,Vec b,Vec x,int uflag)
1226849ba73SBarry Smith {
123f0c56d0fSKris Buschelman   Mat_UMFPACK    *lu = (Mat_UMFPACK*)A->spptr;
1241677d0b8SKris Buschelman   PetscScalar    *av=lu->av,*ba,*xa;
125dfbe8321SBarry Smith   PetscErrorCode ierr;
1268592901bSBarry Smith   PetscInt       *ai=lu->ai,*aj=lu->aj,status;
1271677d0b8SKris Buschelman 
1281677d0b8SKris Buschelman   PetscFunctionBegin;
1292d4e2982SHong Zhang   /* solve Ax = b by umfpack_*_wsolve */
1301677d0b8SKris Buschelman   /* ----------------------------------*/
1312d4e2982SHong Zhang 
1321677d0b8SKris Buschelman   ierr = VecGetArray(b,&ba);
1331677d0b8SKris Buschelman   ierr = VecGetArray(x,&xa);
13479c34000SHong Zhang #if defined(PETSC_USE_COMPLEX)
1354b019bd2SJed Brown   status = umfpack_UMF_wsolve(uflag,ai,aj,(PetscReal*)av,NULL,(PetscReal*)xa,NULL,(PetscReal*)ba,NULL,
13679c34000SHong Zhang                               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 
146057dcbc9SJed Brown   ierr = VecRestoreArray(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;
1806849ba73SBarry Smith   PetscErrorCode ierr;
1818592901bSBarry Smith   PetscInt     *ai=lu->ai,*aj=lu->aj,m=A->rmap->n,status;
1821677d0b8SKris Buschelman   PetscScalar *av=lu->av;
1831677d0b8SKris Buschelman 
1841677d0b8SKris Buschelman   PetscFunctionBegin;
1851677d0b8SKris Buschelman   /* numeric factorization of A' */
1861677d0b8SKris Buschelman   /* ----------------------------*/
1872d4e2982SHong Zhang 
1888592901bSBarry Smith   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
1898592901bSBarry Smith     umfpack_UMF_free_numeric(&lu->Numeric);
1908592901bSBarry Smith   }
1912d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1928592901bSBarry Smith   status = umfpack_UMF_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1932d4e2982SHong Zhang #else
1948592901bSBarry Smith   status = umfpack_UMF_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1958592901bSBarry Smith #endif
1969f42a82aSMatthew Knepley   if (status < 0) {
1978592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
198e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_numeric failed");
1999f42a82aSMatthew Knepley   }
2002d4e2982SHong Zhang   /* report numeric factorization of A' when Control[PRL] > 3 */
2018592901bSBarry Smith   (void) umfpack_UMF_report_numeric(lu->Numeric, lu->Control);
2021677d0b8SKris Buschelman 
2031677d0b8SKris Buschelman   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
2041677d0b8SKris Buschelman     /* allocate working space to be used by Solve */
2058592901bSBarry Smith     ierr = PetscMalloc(m * sizeof(PetscInt), &lu->Wi);CHKERRQ(ierr);
2068592901bSBarry Smith     ierr = PetscMalloc(5*m * sizeof(PetscScalar), &lu->W);CHKERRQ(ierr);
2071677d0b8SKris Buschelman   }
2082d4e2982SHong Zhang 
2092a325a84SHong Zhang   lu->flg = SAME_NONZERO_PATTERN;
2102a325a84SHong Zhang   lu->CleanUpUMFPACK = PETSC_TRUE;
2114b019bd2SJed Brown   F->ops->solve          = MatSolve_UMFPACK;
2124b019bd2SJed Brown   F->ops->solvetranspose = MatSolveTranspose_UMFPACK;
2131677d0b8SKris Buschelman   PetscFunctionReturn(0);
2141677d0b8SKris Buschelman }
2151677d0b8SKris Buschelman 
2161677d0b8SKris Buschelman /*
2171677d0b8SKris Buschelman    Note the r permutation is ignored
2181677d0b8SKris Buschelman */
2191677d0b8SKris Buschelman #undef __FUNCT__
220f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK"
2214b019bd2SJed Brown static PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
2226849ba73SBarry Smith {
2231677d0b8SKris Buschelman   Mat_SeqAIJ     *mat=(Mat_SeqAIJ*)A->data;
224719d5645SBarry Smith   Mat_UMFPACK    *lu = (Mat_UMFPACK*)(F->spptr);
225dfbe8321SBarry Smith   PetscErrorCode ierr;
2268592901bSBarry Smith   PetscInt       i,m=A->rmap->n,n=A->cmap->n;
2275d0c19d7SBarry Smith   const PetscInt *ra;
2288592901bSBarry Smith   PetscInt        status;
2291677d0b8SKris Buschelman   PetscScalar    *av=mat->a;
2301677d0b8SKris Buschelman 
2311677d0b8SKris Buschelman   PetscFunctionBegin;
232fcfd50ebSBarry Smith   if (lu->PetscMatOrdering) {
2330f6efc10SHong Zhang     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
2348592901bSBarry Smith     ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr);
2352d4e2982SHong Zhang     /* we cannot simply memcpy on 64 bit archs */
2362d4e2982SHong Zhang     for(i = 0; i < m; i++) lu->perm_c[i] = ra[i];
2370f6efc10SHong Zhang     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
2381677d0b8SKris Buschelman   }
2391677d0b8SKris Buschelman 
2408592901bSBarry Smith   lu->ai = mat->i;
2418592901bSBarry Smith   lu->aj = mat->j;
2422d4e2982SHong Zhang 
2431677d0b8SKris Buschelman   /* print the control parameters */
2448592901bSBarry Smith   if(lu->Control[UMFPACK_PRL] > 1) umfpack_UMF_report_control(lu->Control);
2451677d0b8SKris Buschelman 
2461677d0b8SKris Buschelman   /* symbolic factorization of A' */
2471677d0b8SKris Buschelman   /* ---------------------------------------------------------------------- */
248fcfd50ebSBarry Smith   if (lu->PetscMatOrdering) { /* use Petsc row ordering */
2498592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
2508592901bSBarry Smith     status = umfpack_UMF_qsymbolic(n,m,lu->ai,lu->aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2512d4e2982SHong Zhang #else
25279c34000SHong Zhang     status = umfpack_UMF_qsymbolic(n,m,lu->ai,lu->aj,NULL,NULL,
25379c34000SHong Zhang                                    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)
2578592901bSBarry Smith     status = umfpack_UMF_symbolic(n,m,lu->ai,lu->aj,av,&lu->Symbolic,lu->Control,lu->Info);
2588592901bSBarry Smith #else
25979c34000SHong Zhang     status = umfpack_UMF_symbolic(n,m,lu->ai,lu->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;
2711677d0b8SKris Buschelman   lu->av  = av;
272e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK = PETSC_TRUE;
273719d5645SBarry Smith   (F)->ops->lufactornumeric  = MatLUFactorNumeric_UMFPACK;
2741677d0b8SKris Buschelman   PetscFunctionReturn(0);
2751677d0b8SKris Buschelman }
2761677d0b8SKris Buschelman 
2771677d0b8SKris Buschelman #undef __FUNCT__
278f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_UMFPACK"
2794b019bd2SJed Brown static PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer)
2806849ba73SBarry Smith {
281f0c56d0fSKris Buschelman   Mat_UMFPACK    *lu= (Mat_UMFPACK*)A->spptr;
282dfbe8321SBarry Smith   PetscErrorCode ierr;
283f0c56d0fSKris Buschelman 
2841677d0b8SKris Buschelman   PetscFunctionBegin;
2851677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
286f0c56d0fSKris Buschelman   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
2871677d0b8SKris Buschelman 
2881677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
2891677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
2901677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
2911677d0b8SKris Buschelman 
2921677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
2930f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr);
2941677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
2951677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
2960f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr);
2971677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
2980f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr);
2990f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr);
3001677d0b8SKris Buschelman 
3011677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
3021677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
3030f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr);
3040f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr);
3051677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
3060f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr);
3071677d0b8SKris Buschelman 
3081677d0b8SKris Buschelman   /* Control parameters used by solve */
3091677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
3101677d0b8SKris Buschelman 
3111677d0b8SKris Buschelman   /* mat ordering */
312fcfd50ebSBarry Smith   if (!lu->PetscMatOrdering) {
313fcfd50ebSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ORDERING]: %s (not using the PETSc ordering)\n",UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]]);CHKERRQ(ierr);
314fcfd50ebSBarry Smith   }
3151677d0b8SKris Buschelman 
3161677d0b8SKris Buschelman   PetscFunctionReturn(0);
3171677d0b8SKris Buschelman }
3181677d0b8SKris Buschelman 
319d844382dSKris Buschelman #undef __FUNCT__
320f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_UMFPACK"
3214b019bd2SJed Brown static PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
3226849ba73SBarry Smith {
323dfbe8321SBarry Smith   PetscErrorCode    ierr;
324ace3abfcSBarry Smith   PetscBool         iascii;
325d844382dSKris Buschelman   PetscViewerFormat format;
326d844382dSKris Buschelman 
327d844382dSKris Buschelman   PetscFunctionBegin;
328b24902e0SBarry Smith   ierr = MatView_SeqAIJ(A,viewer);CHKERRQ(ierr);
329d844382dSKris Buschelman 
3302692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
33132077d6dSBarry Smith   if (iascii) {
332d844382dSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
3332f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO) {
334f0c56d0fSKris Buschelman       ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
335d844382dSKris Buschelman     }
336d844382dSKris Buschelman   }
337d844382dSKris Buschelman   PetscFunctionReturn(0);
338d844382dSKris Buschelman }
339d844382dSKris Buschelman 
34035bd34faSBarry Smith EXTERN_C_BEGIN
34135bd34faSBarry Smith #undef __FUNCT__
34235bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_umfpack"
34335bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_seqaij_umfpack(Mat A,const MatSolverPackage *type)
34435bd34faSBarry Smith {
34535bd34faSBarry Smith   PetscFunctionBegin;
3462692d6eeSBarry Smith   *type = MATSOLVERUMFPACK;
34735bd34faSBarry Smith   PetscFunctionReturn(0);
34835bd34faSBarry Smith }
34935bd34faSBarry Smith EXTERN_C_END
35035bd34faSBarry Smith 
35135bd34faSBarry Smith 
3522f71e704SKris Buschelman /*MC
3532692d6eeSBarry Smith   MATSOLVERUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
3542f71e704SKris Buschelman   via the external package UMFPACK.
3552f71e704SKris Buschelman 
356e2e64c6bSBarry Smith   ./configure --download-umfpack to install PETSc to use UMFPACK
3572f71e704SKris Buschelman 
3582f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
3592f71e704SKris Buschelman   which correspond to the options database keys below.
3602f71e704SKris Buschelman 
3612f71e704SKris Buschelman   Options Database Keys:
3625c9eb25fSBarry Smith + -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL]
3633519fcdcSHong Zhang . -mat_umfpack_strategy <AUTO> (choose one of) AUTO UNSYMMETRIC SYMMETRIC 2BY2
3642f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
3653519fcdcSHong Zhang . -mat_umfpack_dense_row <0.2>: Control[UMFPACK_DENSE_ROW]
3663519fcdcSHong Zhang . -mat_umfpack_amd_dense <10>: Control[UMFPACK_AMD_DENSE]
3672f71e704SKris Buschelman . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
3683519fcdcSHong Zhang . -mat_umfpack_2by2_tolerance <0.01>: Control[UMFPACK_2BY2_TOLERANCE]
3693519fcdcSHong Zhang . -mat_umfpack_fixq <0>: Control[UMFPACK_FIXQ]
3703519fcdcSHong Zhang . -mat_umfpack_aggressive <1>: Control[UMFPACK_AGGRESSIVE]
3712f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
3723519fcdcSHong Zhang .  -mat_umfpack_sym_pivot_tolerance <0.001>: Control[UMFPACK_SYM_PIVOT_TOLERANCE]
3733519fcdcSHong Zhang .  -mat_umfpack_scale <NONE> (choose one of) NONE SUM MAX
3742f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
3753519fcdcSHong Zhang .  -mat_umfpack_droptol <0>: Control[UMFPACK_DROPTOL]
3762f71e704SKris Buschelman - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
3772f71e704SKris Buschelman 
3782f71e704SKris Buschelman    Level: beginner
3792f71e704SKris Buschelman 
3802692d6eeSBarry Smith .seealso: PCLU, MATSOLVERSUPERLU, MATSOLVERMUMPS, MATSOLVERSPOOLES, PCFactorSetMatSolverPackage(), MatSolverPackage
3812f71e704SKris Buschelman M*/
3823519fcdcSHong Zhang EXTERN_C_BEGIN
3833519fcdcSHong Zhang #undef __FUNCT__
3843519fcdcSHong Zhang #define __FUNCT__ "MatGetFactor_seqaij_umfpack"
3853519fcdcSHong Zhang PetscErrorCode MatGetFactor_seqaij_umfpack(Mat A,MatFactorType ftype,Mat *F)
3863519fcdcSHong Zhang {
3873519fcdcSHong Zhang   Mat            B;
3883519fcdcSHong Zhang   Mat_UMFPACK    *lu;
3893519fcdcSHong Zhang   PetscErrorCode ierr;
3908592901bSBarry Smith   PetscInt       m=A->rmap->n,n=A->cmap->n,idx;
3912f71e704SKris Buschelman 
392fcfd50ebSBarry Smith   const char     *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC"},
3933519fcdcSHong Zhang                  *scale[]={"NONE","SUM","MAX"};
394ace3abfcSBarry Smith   PetscBool      flg;
3952d4e2982SHong Zhang 
3963519fcdcSHong Zhang   PetscFunctionBegin;
3973519fcdcSHong Zhang   /* Create the factorization matrix F */
3983519fcdcSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
3993519fcdcSHong Zhang   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
4003519fcdcSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
4013519fcdcSHong Zhang   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
4023519fcdcSHong Zhang   ierr = PetscNewLog(B,Mat_UMFPACK,&lu);CHKERRQ(ierr);
4033519fcdcSHong Zhang   B->spptr                 = lu;
4043519fcdcSHong Zhang   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
4053519fcdcSHong Zhang   B->ops->destroy          = MatDestroy_UMFPACK;
4063519fcdcSHong Zhang   B->ops->view             = MatView_UMFPACK;
40735bd34faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_umfpack",MatFactorGetSolverPackage_seqaij_umfpack);CHKERRQ(ierr);
408d5f3da31SBarry Smith   B->factortype            = MAT_FACTOR_LU;
4093519fcdcSHong Zhang   B->assembled             = PETSC_TRUE;  /* required by -ksp_view */
4103519fcdcSHong Zhang   B->preallocated          = PETSC_TRUE;
4113519fcdcSHong Zhang 
4123519fcdcSHong Zhang   /* initializations */
4133519fcdcSHong Zhang   /* ------------------------------------------------*/
4143519fcdcSHong Zhang   /* get the default control parameters */
4158592901bSBarry Smith   umfpack_UMF_defaults(lu->Control);
4163519fcdcSHong Zhang   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
4173519fcdcSHong Zhang   lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */
4183519fcdcSHong Zhang 
4193519fcdcSHong Zhang   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
4203519fcdcSHong Zhang   /* Control parameters used by reporting routiones */
4213519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr);
4223519fcdcSHong Zhang 
4233519fcdcSHong Zhang   /* Control parameters for symbolic factorization */
424fcfd50ebSBarry Smith   ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,3,strategy[0],&idx,&flg);CHKERRQ(ierr);
4253519fcdcSHong Zhang   if (flg) {
4263519fcdcSHong Zhang     switch (idx){
4273519fcdcSHong Zhang     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
4283519fcdcSHong Zhang     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
4293519fcdcSHong Zhang     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
4303519fcdcSHong Zhang     }
4313519fcdcSHong Zhang   }
432fcfd50ebSBarry 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);
433fcfd50ebSBarry Smith   if (flg) lu->Control[UMFPACK_ORDERING] = (int)idx;
4343519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_dense_col","Control[UMFPACK_DENSE_COL]","None",lu->Control[UMFPACK_DENSE_COL],&lu->Control[UMFPACK_DENSE_COL],PETSC_NULL);CHKERRQ(ierr);
4353519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_dense_row","Control[UMFPACK_DENSE_ROW]","None",lu->Control[UMFPACK_DENSE_ROW],&lu->Control[UMFPACK_DENSE_ROW],PETSC_NULL);CHKERRQ(ierr);
4363519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_amd_dense","Control[UMFPACK_AMD_DENSE]","None",lu->Control[UMFPACK_AMD_DENSE],&lu->Control[UMFPACK_AMD_DENSE],PETSC_NULL);CHKERRQ(ierr);
4373519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_block_size","Control[UMFPACK_BLOCK_SIZE]","None",lu->Control[UMFPACK_BLOCK_SIZE],&lu->Control[UMFPACK_BLOCK_SIZE],PETSC_NULL);CHKERRQ(ierr);
4383519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);CHKERRQ(ierr);
4393519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);CHKERRQ(ierr);
4403519fcdcSHong Zhang 
4413519fcdcSHong Zhang   /* Control parameters used by numeric factorization */
4423519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_pivot_tolerance","Control[UMFPACK_PIVOT_TOLERANCE]","None",lu->Control[UMFPACK_PIVOT_TOLERANCE],&lu->Control[UMFPACK_PIVOT_TOLERANCE],PETSC_NULL);CHKERRQ(ierr);
4433519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_sym_pivot_tolerance","Control[UMFPACK_SYM_PIVOT_TOLERANCE]","None",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE],&lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE],PETSC_NULL);CHKERRQ(ierr);
4443519fcdcSHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
4453519fcdcSHong Zhang   if (flg) {
4463519fcdcSHong Zhang     switch (idx){
4473519fcdcSHong Zhang     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
4483519fcdcSHong Zhang     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
4493519fcdcSHong Zhang     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
4503519fcdcSHong Zhang     }
4513519fcdcSHong Zhang   }
4523519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_alloc_init","Control[UMFPACK_ALLOC_INIT]","None",lu->Control[UMFPACK_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],PETSC_NULL);CHKERRQ(ierr);
4533519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_front_alloc_init","Control[UMFPACK_FRONT_ALLOC_INIT]","None",lu->Control[UMFPACK_FRONT_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],PETSC_NULL);CHKERRQ(ierr);
4543519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);CHKERRQ(ierr);
4553519fcdcSHong Zhang 
4563519fcdcSHong Zhang   /* Control parameters used by solve */
4573519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr);
4583519fcdcSHong Zhang 
4593519fcdcSHong Zhang   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */
460fcfd50ebSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-pc_factor_mat_ordering_type",&lu->PetscMatOrdering);CHKERRQ(ierr);
4613519fcdcSHong Zhang   PetscOptionsEnd();
4623519fcdcSHong Zhang   *F = B;
4633519fcdcSHong Zhang   PetscFunctionReturn(0);
4643519fcdcSHong Zhang }
4653519fcdcSHong Zhang EXTERN_C_END
4662d4e2982SHong Zhang 
467