xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 41c8de11b2c11272bfdf610cc930de5c376c86d6)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
21677d0b8SKris Buschelman 
31677d0b8SKris Buschelman /*
42d4e2982SHong Zhang    Provides an interface to the UMFPACKv5.1 sparse solver
52d4e2982SHong Zhang 
68592901bSBarry Smith    When build with PETSC_USE_64BIT_INDICES this will use UF_Long as the
78592901bSBarry Smith    integer type in UMFPACK, otherwise it will use int. This means
88592901bSBarry Smith    all integers in this file as simply declared as PetscInt. Also it means
98592901bSBarry Smith    that UMFPACK UL_Long version MUST be built with 64 bit integers when used.
102d4e2982SHong Zhang 
111677d0b8SKris Buschelman */
127c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h"
131677d0b8SKris Buschelman 
148592901bSBarry Smith #if defined(PETSC_USE_64BIT_INDICES)
158592901bSBarry Smith #if defined(PETSC_USE_COMPLEX)
168592901bSBarry Smith #define umfpack_UMF_free_symbolic   umfpack_zl_free_symbolic
178592901bSBarry Smith #define umfpack_UMF_free_numeric    umfpack_zl_free_numeric
188592901bSBarry Smith #define umfpack_UMF_wsolve          umfpack_zl_wsolve
198592901bSBarry Smith #define umfpack_UMF_numeric         umfpack_zl_numeric
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
258592901bSBarry Smith #define umfpack_UMF_qsymbolic       umfpack_zl_qsymbolic
268592901bSBarry Smith #define umfpack_UMF_symbolic        umfpack_zl_symbolic
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
328592901bSBarry Smith #define umfpack_UMF_wsolve          umfpack_dl_wsolve
338592901bSBarry Smith #define umfpack_UMF_numeric         umfpack_dl_numeric
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
398592901bSBarry Smith #define umfpack_UMF_qsymbolic       umfpack_dl_qsymbolic
408592901bSBarry Smith #define umfpack_UMF_symbolic        umfpack_dl_symbolic
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 
758592901bSBarry Smith 
768592901bSBarry Smith #define UF_long long long
778592901bSBarry Smith #define UF_long_max LONG_LONG_MAX
788592901bSBarry Smith #define UF_long_id "%lld"
798592901bSBarry Smith 
801677d0b8SKris Buschelman EXTERN_C_BEGIN
811677d0b8SKris Buschelman #include "umfpack.h"
821677d0b8SKris Buschelman EXTERN_C_END
831677d0b8SKris Buschelman 
841677d0b8SKris Buschelman typedef struct {
851677d0b8SKris Buschelman   void         *Symbolic, *Numeric;
861677d0b8SKris Buschelman   double       Info[UMFPACK_INFO], Control[UMFPACK_CONTROL],*W;
878592901bSBarry Smith   PetscInt      *Wi,*ai,*aj,*perm_c;
881677d0b8SKris Buschelman   PetscScalar  *av;
891677d0b8SKris Buschelman   MatStructure flg;
901677d0b8SKris Buschelman   PetscTruth   PetscMatOdering;
911677d0b8SKris Buschelman 
921677d0b8SKris Buschelman   /* Flag to clean up UMFPACK objects during Destroy */
931677d0b8SKris Buschelman   PetscTruth CleanUpUMFPACK;
94f0c56d0fSKris Buschelman } Mat_UMFPACK;
95f0c56d0fSKris Buschelman 
961677d0b8SKris Buschelman #undef __FUNCT__
97f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_UMFPACK"
98dfbe8321SBarry Smith PetscErrorCode MatDestroy_UMFPACK(Mat A)
99dfbe8321SBarry Smith {
100dfbe8321SBarry Smith   PetscErrorCode ierr;
101f0c56d0fSKris Buschelman   Mat_UMFPACK    *lu=(Mat_UMFPACK*)A->spptr;
1021677d0b8SKris Buschelman 
1031677d0b8SKris Buschelman   PetscFunctionBegin;
104fb731535SKris Buschelman   if (lu->CleanUpUMFPACK) {
1058592901bSBarry Smith     umfpack_UMF_free_symbolic(&lu->Symbolic);
1068592901bSBarry Smith     umfpack_UMF_free_numeric(&lu->Numeric);
1071677d0b8SKris Buschelman     ierr = PetscFree(lu->Wi);CHKERRQ(ierr);
1081677d0b8SKris Buschelman     ierr = PetscFree(lu->W);CHKERRQ(ierr);
1091677d0b8SKris Buschelman     if (lu->PetscMatOdering) {
1101677d0b8SKris Buschelman       ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
1111677d0b8SKris Buschelman     }
1121677d0b8SKris Buschelman   }
113b24902e0SBarry Smith   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
1141677d0b8SKris Buschelman   PetscFunctionReturn(0);
1151677d0b8SKris Buschelman }
1161677d0b8SKris Buschelman 
1171677d0b8SKris Buschelman #undef __FUNCT__
118f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_UMFPACK"
1196849ba73SBarry Smith PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x)
1206849ba73SBarry Smith {
121f0c56d0fSKris Buschelman   Mat_UMFPACK    *lu = (Mat_UMFPACK*)A->spptr;
1221677d0b8SKris Buschelman   PetscScalar    *av=lu->av,*ba,*xa;
123dfbe8321SBarry Smith   PetscErrorCode ierr;
1248592901bSBarry Smith   PetscInt       *ai=lu->ai,*aj=lu->aj,status;
1251677d0b8SKris Buschelman 
1261677d0b8SKris Buschelman   PetscFunctionBegin;
1272d4e2982SHong Zhang   /* solve Ax = b by umfpack_*_wsolve */
1281677d0b8SKris Buschelman   /* ----------------------------------*/
1292d4e2982SHong Zhang   ierr = VecConjugate(b);
1302d4e2982SHong Zhang 
1311677d0b8SKris Buschelman   ierr = VecGetArray(b,&ba);
1321677d0b8SKris Buschelman   ierr = VecGetArray(x,&xa);
13379c34000SHong Zhang #if defined(PETSC_USE_COMPLEX)
13479c34000SHong Zhang   status = umfpack_UMF_wsolve(UMFPACK_At,ai,aj,(PetscReal*)av,NULL,(PetscReal*)xa,NULL,(PetscReal*)ba,NULL,
13579c34000SHong Zhang                               lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
13679c34000SHong Zhang #else
1378592901bSBarry Smith   status = umfpack_UMF_wsolve(UMFPACK_At,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
13879c34000SHong Zhang #endif
1398592901bSBarry Smith   umfpack_UMF_report_info(lu->Control, lu->Info);
1401677d0b8SKris Buschelman   if (status < 0){
1418592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
1428592901bSBarry Smith     SETERRQ(PETSC_ERR_LIB,"umfpack_UMF_wsolve failed");
1431677d0b8SKris Buschelman   }
1441677d0b8SKris Buschelman 
1451677d0b8SKris Buschelman   ierr = VecRestoreArray(b,&ba);
1461677d0b8SKris Buschelman   ierr = VecRestoreArray(x,&xa);
1472a325a84SHong Zhang 
1482d4e2982SHong Zhang   ierr = VecConjugate(b);
1495b2a9f60SHong Zhang   ierr = VecConjugate(x);
1501677d0b8SKris Buschelman   PetscFunctionReturn(0);
1511677d0b8SKris Buschelman }
1521677d0b8SKris Buschelman 
1531677d0b8SKris Buschelman #undef __FUNCT__
154f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_UMFPACK"
1550481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat F,Mat A,const MatFactorInfo *info)
1566849ba73SBarry Smith {
157719d5645SBarry Smith   Mat_UMFPACK *lu=(Mat_UMFPACK*)(F)->spptr;
1586849ba73SBarry Smith   PetscErrorCode ierr;
1598592901bSBarry Smith   PetscInt     *ai=lu->ai,*aj=lu->aj,m=A->rmap->n,status;
1601677d0b8SKris Buschelman   PetscScalar *av=lu->av;
1611677d0b8SKris Buschelman 
1621677d0b8SKris Buschelman   PetscFunctionBegin;
1631677d0b8SKris Buschelman   /* numeric factorization of A' */
1641677d0b8SKris Buschelman   /* ----------------------------*/
1652d4e2982SHong Zhang 
1668592901bSBarry Smith   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
1678592901bSBarry Smith     umfpack_UMF_free_numeric(&lu->Numeric);
1688592901bSBarry Smith   }
1692d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1708592901bSBarry Smith   status = umfpack_UMF_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1712d4e2982SHong Zhang #else
1728592901bSBarry Smith   status = umfpack_UMF_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1738592901bSBarry Smith #endif
1749f42a82aSMatthew Knepley   if (status < 0) {
1758592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
1768592901bSBarry Smith     SETERRQ(PETSC_ERR_LIB,"umfpack_UMF_numeric failed");
1779f42a82aSMatthew Knepley   }
1782d4e2982SHong Zhang   /* report numeric factorization of A' when Control[PRL] > 3 */
1798592901bSBarry Smith   (void) umfpack_UMF_report_numeric(lu->Numeric, lu->Control);
1801677d0b8SKris Buschelman 
1811677d0b8SKris Buschelman   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
1821677d0b8SKris Buschelman     /* allocate working space to be used by Solve */
1838592901bSBarry Smith     ierr = PetscMalloc(m * sizeof(PetscInt), &lu->Wi);CHKERRQ(ierr);
1848592901bSBarry Smith     ierr = PetscMalloc(5*m * sizeof(PetscScalar), &lu->W);CHKERRQ(ierr);
1851677d0b8SKris Buschelman   }
1862d4e2982SHong Zhang 
1872a325a84SHong Zhang   lu->flg = SAME_NONZERO_PATTERN;
1882a325a84SHong Zhang   lu->CleanUpUMFPACK = PETSC_TRUE;
189719d5645SBarry Smith   (F)->ops->solve            = MatSolve_UMFPACK;
1901677d0b8SKris Buschelman   PetscFunctionReturn(0);
1911677d0b8SKris Buschelman }
1921677d0b8SKris Buschelman 
1931677d0b8SKris Buschelman /*
1941677d0b8SKris Buschelman    Note the r permutation is ignored
1951677d0b8SKris Buschelman */
1961677d0b8SKris Buschelman #undef __FUNCT__
197f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK"
1980481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1996849ba73SBarry Smith {
2001677d0b8SKris Buschelman   Mat_SeqAIJ     *mat=(Mat_SeqAIJ*)A->data;
201719d5645SBarry Smith   Mat_UMFPACK    *lu = (Mat_UMFPACK*)(F->spptr);
202dfbe8321SBarry Smith   PetscErrorCode ierr;
2038592901bSBarry Smith   PetscInt       i,m=A->rmap->n,n=A->cmap->n;
2045d0c19d7SBarry Smith   const PetscInt *ra;
2058592901bSBarry Smith   PetscInt        status;
2061677d0b8SKris Buschelman   PetscScalar    *av=mat->a;
2071677d0b8SKris Buschelman 
2081677d0b8SKris Buschelman   PetscFunctionBegin;
2091677d0b8SKris Buschelman   if (lu->PetscMatOdering) {
2100f6efc10SHong Zhang     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
2118592901bSBarry Smith     ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr);
2122d4e2982SHong Zhang     /* we cannot simply memcpy on 64 bit archs */
2132d4e2982SHong Zhang     for(i = 0; i < m; i++) lu->perm_c[i] = ra[i];
2140f6efc10SHong Zhang     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
2151677d0b8SKris Buschelman   }
2161677d0b8SKris Buschelman 
2178592901bSBarry Smith   lu->ai = mat->i;
2188592901bSBarry Smith   lu->aj = mat->j;
2192d4e2982SHong Zhang 
2201677d0b8SKris Buschelman   /* print the control parameters */
2218592901bSBarry Smith   if(lu->Control[UMFPACK_PRL] > 1) umfpack_UMF_report_control(lu->Control);
2221677d0b8SKris Buschelman 
2231677d0b8SKris Buschelman   /* symbolic factorization of A' */
2241677d0b8SKris Buschelman   /* ---------------------------------------------------------------------- */
2250f6efc10SHong Zhang   if (lu->PetscMatOdering) { /* use Petsc row ordering */
2268592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
2278592901bSBarry Smith     status = umfpack_UMF_qsymbolic(n,m,lu->ai,lu->aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2282d4e2982SHong Zhang #else
22979c34000SHong Zhang     status = umfpack_UMF_qsymbolic(n,m,lu->ai,lu->aj,NULL,NULL,
23079c34000SHong Zhang                                    lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2318592901bSBarry Smith #endif
2322d4e2982SHong Zhang   } else { /* use Umfpack col ordering */
2338592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
2348592901bSBarry Smith     status = umfpack_UMF_symbolic(n,m,lu->ai,lu->aj,av,&lu->Symbolic,lu->Control,lu->Info);
2358592901bSBarry Smith #else
23679c34000SHong Zhang     status = umfpack_UMF_symbolic(n,m,lu->ai,lu->aj,NULL,NULL,&lu->Symbolic,lu->Control,lu->Info);
2378592901bSBarry Smith #endif
2382d4e2982SHong Zhang   }
2392d4e2982SHong Zhang   if (status < 0){
2408592901bSBarry Smith     umfpack_UMF_report_info(lu->Control, lu->Info);
2418592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
2428592901bSBarry Smith     SETERRQ(PETSC_ERR_LIB,"umfpack_UMF_symbolic failed");
2432d4e2982SHong Zhang   }
2442d4e2982SHong Zhang   /* report sumbolic factorization of A' when Control[PRL] > 3 */
2458592901bSBarry Smith   (void) umfpack_UMF_report_symbolic(lu->Symbolic, lu->Control);
2461677d0b8SKris Buschelman 
2471677d0b8SKris Buschelman   lu->flg = DIFFERENT_NONZERO_PATTERN;
2481677d0b8SKris Buschelman   lu->av  = av;
249e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK = PETSC_TRUE;
250719d5645SBarry Smith   (F)->ops->lufactornumeric  = MatLUFactorNumeric_UMFPACK;
2511677d0b8SKris Buschelman   PetscFunctionReturn(0);
2521677d0b8SKris Buschelman }
2531677d0b8SKris Buschelman 
2541677d0b8SKris Buschelman #undef __FUNCT__
255f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_UMFPACK"
2566849ba73SBarry Smith PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer)
2576849ba73SBarry Smith {
258f0c56d0fSKris Buschelman   Mat_UMFPACK    *lu= (Mat_UMFPACK*)A->spptr;
259dfbe8321SBarry Smith   PetscErrorCode ierr;
260f0c56d0fSKris Buschelman 
2611677d0b8SKris Buschelman   PetscFunctionBegin;
2621677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
263f0c56d0fSKris Buschelman   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
2641677d0b8SKris Buschelman 
2651677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
2661677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
2671677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
2681677d0b8SKris Buschelman 
2691677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
2700f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr);
2711677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
2721677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
2730f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr);
2741677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
2750f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_2BY2_TOLERANCE]: %g\n",lu->Control[UMFPACK_2BY2_TOLERANCE]);CHKERRQ(ierr);
2760f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr);
2770f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr);
2781677d0b8SKris Buschelman 
2791677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
2801677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
2810f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr);
2820f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr);
2831677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
2840f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr);
2851677d0b8SKris Buschelman 
2861677d0b8SKris Buschelman   /* Control parameters used by solve */
2871677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
2881677d0b8SKris Buschelman 
2891677d0b8SKris Buschelman   /* mat ordering */
2901677d0b8SKris Buschelman   if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer,"  UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr);
2911677d0b8SKris Buschelman 
2921677d0b8SKris Buschelman   PetscFunctionReturn(0);
2931677d0b8SKris Buschelman }
2941677d0b8SKris Buschelman 
295d844382dSKris Buschelman #undef __FUNCT__
296f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_UMFPACK"
2976849ba73SBarry Smith PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
2986849ba73SBarry Smith {
299dfbe8321SBarry Smith   PetscErrorCode    ierr;
30032077d6dSBarry Smith   PetscTruth        iascii;
301d844382dSKris Buschelman   PetscViewerFormat format;
302d844382dSKris Buschelman 
303d844382dSKris Buschelman   PetscFunctionBegin;
304b24902e0SBarry Smith   ierr = MatView_SeqAIJ(A,viewer);CHKERRQ(ierr);
305d844382dSKris Buschelman 
30632077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
30732077d6dSBarry Smith   if (iascii) {
308d844382dSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
3092f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO) {
310f0c56d0fSKris Buschelman       ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
311d844382dSKris Buschelman     }
312d844382dSKris Buschelman   }
313d844382dSKris Buschelman   PetscFunctionReturn(0);
314d844382dSKris Buschelman }
315d844382dSKris Buschelman 
31635bd34faSBarry Smith EXTERN_C_BEGIN
31735bd34faSBarry Smith #undef __FUNCT__
31835bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_umfpack"
31935bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_seqaij_umfpack(Mat A,const MatSolverPackage *type)
32035bd34faSBarry Smith {
32135bd34faSBarry Smith   PetscFunctionBegin;
32235bd34faSBarry Smith   *type = MAT_SOLVER_UMFPACK;
32335bd34faSBarry Smith   PetscFunctionReturn(0);
32435bd34faSBarry Smith }
32535bd34faSBarry Smith EXTERN_C_END
32635bd34faSBarry Smith 
32735bd34faSBarry Smith 
3282f71e704SKris Buschelman /*MC
3295c9eb25fSBarry Smith   MAT_SOLVER_UMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
3302f71e704SKris Buschelman   via the external package UMFPACK.
3312f71e704SKris Buschelman 
3325c9eb25fSBarry Smith   config/configure.py --download-umfpack to install PETSc to use UMFPACK
3332f71e704SKris Buschelman 
3342f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
3352f71e704SKris Buschelman   which correspond to the options database keys below.
3362f71e704SKris Buschelman 
3372f71e704SKris Buschelman   Options Database Keys:
3385c9eb25fSBarry Smith + -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL]
3393519fcdcSHong Zhang . -mat_umfpack_strategy <AUTO> (choose one of) AUTO UNSYMMETRIC SYMMETRIC 2BY2
3402f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
3413519fcdcSHong Zhang . -mat_umfpack_dense_row <0.2>: Control[UMFPACK_DENSE_ROW]
3423519fcdcSHong Zhang . -mat_umfpack_amd_dense <10>: Control[UMFPACK_AMD_DENSE]
3432f71e704SKris Buschelman . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
3443519fcdcSHong Zhang . -mat_umfpack_2by2_tolerance <0.01>: Control[UMFPACK_2BY2_TOLERANCE]
3453519fcdcSHong Zhang . -mat_umfpack_fixq <0>: Control[UMFPACK_FIXQ]
3463519fcdcSHong Zhang . -mat_umfpack_aggressive <1>: Control[UMFPACK_AGGRESSIVE]
3472f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
3483519fcdcSHong Zhang .  -mat_umfpack_sym_pivot_tolerance <0.001>: Control[UMFPACK_SYM_PIVOT_TOLERANCE]
3493519fcdcSHong Zhang .  -mat_umfpack_scale <NONE> (choose one of) NONE SUM MAX
3502f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
3513519fcdcSHong Zhang .  -mat_umfpack_droptol <0>: Control[UMFPACK_DROPTOL]
3522f71e704SKris Buschelman - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
3532f71e704SKris Buschelman 
3542f71e704SKris Buschelman    Level: beginner
3552f71e704SKris Buschelman 
356*41c8de11SBarry Smith .seealso: PCLU, MAT_SOLVER_SUPERLU, MAT_SOLVER_MUMPS, MAT_SOLVER_SPOOLES, PCFactorSetMatSolverPackage(), MatSolverPackage
3572f71e704SKris Buschelman M*/
3583519fcdcSHong Zhang EXTERN_C_BEGIN
3593519fcdcSHong Zhang #undef __FUNCT__
3603519fcdcSHong Zhang #define __FUNCT__ "MatGetFactor_seqaij_umfpack"
3613519fcdcSHong Zhang PetscErrorCode MatGetFactor_seqaij_umfpack(Mat A,MatFactorType ftype,Mat *F)
3623519fcdcSHong Zhang {
3633519fcdcSHong Zhang   Mat            B;
3643519fcdcSHong Zhang   Mat_UMFPACK    *lu;
3653519fcdcSHong Zhang   PetscErrorCode ierr;
3668592901bSBarry Smith   PetscInt       m=A->rmap->n,n=A->cmap->n,idx;
3672f71e704SKris Buschelman 
3683519fcdcSHong Zhang   const char     *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC","2BY2"},
3693519fcdcSHong Zhang                  *scale[]={"NONE","SUM","MAX"};
3703519fcdcSHong Zhang   PetscTruth     flg;
3712d4e2982SHong Zhang 
3723519fcdcSHong Zhang   PetscFunctionBegin;
3733519fcdcSHong Zhang   /* Create the factorization matrix F */
3743519fcdcSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
3753519fcdcSHong Zhang   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
3763519fcdcSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
3773519fcdcSHong Zhang   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
3783519fcdcSHong Zhang   ierr = PetscNewLog(B,Mat_UMFPACK,&lu);CHKERRQ(ierr);
3793519fcdcSHong Zhang   B->spptr                 = lu;
3803519fcdcSHong Zhang   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
3813519fcdcSHong Zhang   B->ops->destroy          = MatDestroy_UMFPACK;
3823519fcdcSHong Zhang   B->ops->view             = MatView_UMFPACK;
38335bd34faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_umfpack",MatFactorGetSolverPackage_seqaij_umfpack);CHKERRQ(ierr);
3843519fcdcSHong Zhang   B->factor                = MAT_FACTOR_LU;
3853519fcdcSHong Zhang   B->assembled             = PETSC_TRUE;  /* required by -ksp_view */
3863519fcdcSHong Zhang   B->preallocated          = PETSC_TRUE;
3873519fcdcSHong Zhang 
3883519fcdcSHong Zhang   /* initializations */
3893519fcdcSHong Zhang   /* ------------------------------------------------*/
3903519fcdcSHong Zhang   /* get the default control parameters */
3918592901bSBarry Smith   umfpack_UMF_defaults(lu->Control);
3923519fcdcSHong Zhang   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
3933519fcdcSHong Zhang   lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */
3943519fcdcSHong Zhang 
3953519fcdcSHong Zhang   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
3963519fcdcSHong Zhang   /* Control parameters used by reporting routiones */
3973519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr);
3983519fcdcSHong Zhang 
3993519fcdcSHong Zhang   /* Control parameters for symbolic factorization */
4003519fcdcSHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,4,strategy[0],&idx,&flg);CHKERRQ(ierr);
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     case 3: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_2BY2; break;
4073519fcdcSHong Zhang     }
4083519fcdcSHong Zhang   }
4093519fcdcSHong 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);
4103519fcdcSHong 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);
4113519fcdcSHong 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);
4123519fcdcSHong 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);
4133519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_2by2_tolerance","Control[UMFPACK_2BY2_TOLERANCE]","None",lu->Control[UMFPACK_2BY2_TOLERANCE],&lu->Control[UMFPACK_2BY2_TOLERANCE],PETSC_NULL);CHKERRQ(ierr);
4143519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);CHKERRQ(ierr);
4153519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);CHKERRQ(ierr);
4163519fcdcSHong Zhang 
4173519fcdcSHong Zhang   /* Control parameters used by numeric factorization */
4183519fcdcSHong 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);
4193519fcdcSHong 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);
4203519fcdcSHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
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   }
4283519fcdcSHong 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);
4293519fcdcSHong 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);
4303519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);CHKERRQ(ierr);
4313519fcdcSHong Zhang 
4323519fcdcSHong Zhang   /* Control parameters used by solve */
4333519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr);
4343519fcdcSHong Zhang 
4353519fcdcSHong Zhang   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */
4363519fcdcSHong Zhang   ierr = PetscOptionsHasName(PETSC_NULL,"-pc_factor_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr);
4373519fcdcSHong Zhang   PetscOptionsEnd();
4383519fcdcSHong Zhang   *F = B;
4393519fcdcSHong Zhang   PetscFunctionReturn(0);
4403519fcdcSHong Zhang }
4413519fcdcSHong Zhang EXTERN_C_END
4422d4e2982SHong Zhang 
443