xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 057dcbc9f699b55ddef3fe6bbdd874f47d40468b)
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"
984b019bd2SJed Brown static 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__
1184b019bd2SJed Brown #define __FUNCT__ "MatSolve_UMFPACK_Private"
1194b019bd2SJed Brown static PetscErrorCode MatSolve_UMFPACK_Private(Mat A,Vec b,Vec x,int uflag)
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 
1301677d0b8SKris Buschelman   ierr = VecGetArray(b,&ba);
1311677d0b8SKris Buschelman   ierr = VecGetArray(x,&xa);
13279c34000SHong Zhang #if defined(PETSC_USE_COMPLEX)
1334b019bd2SJed Brown   status = umfpack_UMF_wsolve(uflag,ai,aj,(PetscReal*)av,NULL,(PetscReal*)xa,NULL,(PetscReal*)ba,NULL,
13479c34000SHong Zhang                               lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
13579c34000SHong Zhang #else
1364b019bd2SJed Brown   status = umfpack_UMF_wsolve(uflag,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
13779c34000SHong Zhang #endif
1388592901bSBarry Smith   umfpack_UMF_report_info(lu->Control, lu->Info);
1391677d0b8SKris Buschelman   if (status < 0){
1408592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
141e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_wsolve failed");
1421677d0b8SKris Buschelman   }
1431677d0b8SKris Buschelman 
144*057dcbc9SJed Brown   ierr = VecRestoreArray(b,&ba);CHKERRQ(ierr);
145*057dcbc9SJed Brown   ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr);
1464b019bd2SJed Brown   PetscFunctionReturn(0);
1474b019bd2SJed Brown }
1482a325a84SHong Zhang 
1494b019bd2SJed Brown #undef __FUNCT__
1504b019bd2SJed Brown #define __FUNCT__ "MatSolve_UMFPACK"
1514b019bd2SJed Brown static PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x)
1524b019bd2SJed Brown {
1534b019bd2SJed Brown   PetscErrorCode ierr;
1544b019bd2SJed Brown 
1554b019bd2SJed Brown   PetscFunctionBegin;
1564b019bd2SJed Brown   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
1574b019bd2SJed Brown   ierr = MatSolve_UMFPACK_Private(A,b,x,UMFPACK_Aat);CHKERRQ(ierr);
1584b019bd2SJed Brown   PetscFunctionReturn(0);
1594b019bd2SJed Brown }
1604b019bd2SJed Brown 
1614b019bd2SJed Brown #undef __FUNCT__
1624b019bd2SJed Brown #define __FUNCT__ "MatSolveTranspose_UMFPACK"
1634b019bd2SJed Brown static PetscErrorCode MatSolveTranspose_UMFPACK(Mat A,Vec b,Vec x)
1644b019bd2SJed Brown {
1654b019bd2SJed Brown   PetscErrorCode ierr;
1664b019bd2SJed Brown 
1674b019bd2SJed Brown   PetscFunctionBegin;
1684b019bd2SJed Brown   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
1694b019bd2SJed Brown   ierr = MatSolve_UMFPACK_Private(A,b,x,UMFPACK_A);CHKERRQ(ierr);
1701677d0b8SKris Buschelman   PetscFunctionReturn(0);
1711677d0b8SKris Buschelman }
1721677d0b8SKris Buschelman 
1731677d0b8SKris Buschelman #undef __FUNCT__
174f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_UMFPACK"
1754b019bd2SJed Brown static PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat F,Mat A,const MatFactorInfo *info)
1766849ba73SBarry Smith {
177719d5645SBarry Smith   Mat_UMFPACK *lu=(Mat_UMFPACK*)(F)->spptr;
1786849ba73SBarry Smith   PetscErrorCode ierr;
1798592901bSBarry Smith   PetscInt     *ai=lu->ai,*aj=lu->aj,m=A->rmap->n,status;
1801677d0b8SKris Buschelman   PetscScalar *av=lu->av;
1811677d0b8SKris Buschelman 
1821677d0b8SKris Buschelman   PetscFunctionBegin;
1831677d0b8SKris Buschelman   /* numeric factorization of A' */
1841677d0b8SKris Buschelman   /* ----------------------------*/
1852d4e2982SHong Zhang 
1868592901bSBarry Smith   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
1878592901bSBarry Smith     umfpack_UMF_free_numeric(&lu->Numeric);
1888592901bSBarry Smith   }
1892d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1908592901bSBarry Smith   status = umfpack_UMF_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1912d4e2982SHong Zhang #else
1928592901bSBarry Smith   status = umfpack_UMF_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1938592901bSBarry Smith #endif
1949f42a82aSMatthew Knepley   if (status < 0) {
1958592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
196e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_numeric failed");
1979f42a82aSMatthew Knepley   }
1982d4e2982SHong Zhang   /* report numeric factorization of A' when Control[PRL] > 3 */
1998592901bSBarry Smith   (void) umfpack_UMF_report_numeric(lu->Numeric, lu->Control);
2001677d0b8SKris Buschelman 
2011677d0b8SKris Buschelman   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
2021677d0b8SKris Buschelman     /* allocate working space to be used by Solve */
2038592901bSBarry Smith     ierr = PetscMalloc(m * sizeof(PetscInt), &lu->Wi);CHKERRQ(ierr);
2048592901bSBarry Smith     ierr = PetscMalloc(5*m * sizeof(PetscScalar), &lu->W);CHKERRQ(ierr);
2051677d0b8SKris Buschelman   }
2062d4e2982SHong Zhang 
2072a325a84SHong Zhang   lu->flg = SAME_NONZERO_PATTERN;
2082a325a84SHong Zhang   lu->CleanUpUMFPACK = PETSC_TRUE;
2094b019bd2SJed Brown   F->ops->solve          = MatSolve_UMFPACK;
2104b019bd2SJed Brown   F->ops->solvetranspose = MatSolveTranspose_UMFPACK;
2111677d0b8SKris Buschelman   PetscFunctionReturn(0);
2121677d0b8SKris Buschelman }
2131677d0b8SKris Buschelman 
2141677d0b8SKris Buschelman /*
2151677d0b8SKris Buschelman    Note the r permutation is ignored
2161677d0b8SKris Buschelman */
2171677d0b8SKris Buschelman #undef __FUNCT__
218f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK"
2194b019bd2SJed Brown static PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
2206849ba73SBarry Smith {
2211677d0b8SKris Buschelman   Mat_SeqAIJ     *mat=(Mat_SeqAIJ*)A->data;
222719d5645SBarry Smith   Mat_UMFPACK    *lu = (Mat_UMFPACK*)(F->spptr);
223dfbe8321SBarry Smith   PetscErrorCode ierr;
2248592901bSBarry Smith   PetscInt       i,m=A->rmap->n,n=A->cmap->n;
2255d0c19d7SBarry Smith   const PetscInt *ra;
2268592901bSBarry Smith   PetscInt        status;
2271677d0b8SKris Buschelman   PetscScalar    *av=mat->a;
2281677d0b8SKris Buschelman 
2291677d0b8SKris Buschelman   PetscFunctionBegin;
2301677d0b8SKris Buschelman   if (lu->PetscMatOdering) {
2310f6efc10SHong Zhang     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
2328592901bSBarry Smith     ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr);
2332d4e2982SHong Zhang     /* we cannot simply memcpy on 64 bit archs */
2342d4e2982SHong Zhang     for(i = 0; i < m; i++) lu->perm_c[i] = ra[i];
2350f6efc10SHong Zhang     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
2361677d0b8SKris Buschelman   }
2371677d0b8SKris Buschelman 
2388592901bSBarry Smith   lu->ai = mat->i;
2398592901bSBarry Smith   lu->aj = mat->j;
2402d4e2982SHong Zhang 
2411677d0b8SKris Buschelman   /* print the control parameters */
2428592901bSBarry Smith   if(lu->Control[UMFPACK_PRL] > 1) umfpack_UMF_report_control(lu->Control);
2431677d0b8SKris Buschelman 
2441677d0b8SKris Buschelman   /* symbolic factorization of A' */
2451677d0b8SKris Buschelman   /* ---------------------------------------------------------------------- */
2460f6efc10SHong Zhang   if (lu->PetscMatOdering) { /* use Petsc row ordering */
2478592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
2488592901bSBarry Smith     status = umfpack_UMF_qsymbolic(n,m,lu->ai,lu->aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2492d4e2982SHong Zhang #else
25079c34000SHong Zhang     status = umfpack_UMF_qsymbolic(n,m,lu->ai,lu->aj,NULL,NULL,
25179c34000SHong Zhang                                    lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2528592901bSBarry Smith #endif
2532d4e2982SHong Zhang   } else { /* use Umfpack col ordering */
2548592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
2558592901bSBarry Smith     status = umfpack_UMF_symbolic(n,m,lu->ai,lu->aj,av,&lu->Symbolic,lu->Control,lu->Info);
2568592901bSBarry Smith #else
25779c34000SHong Zhang     status = umfpack_UMF_symbolic(n,m,lu->ai,lu->aj,NULL,NULL,&lu->Symbolic,lu->Control,lu->Info);
2588592901bSBarry Smith #endif
2592d4e2982SHong Zhang   }
2602d4e2982SHong Zhang   if (status < 0){
2618592901bSBarry Smith     umfpack_UMF_report_info(lu->Control, lu->Info);
2628592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
263e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_symbolic failed");
2642d4e2982SHong Zhang   }
2652d4e2982SHong Zhang   /* report sumbolic factorization of A' when Control[PRL] > 3 */
2668592901bSBarry Smith   (void) umfpack_UMF_report_symbolic(lu->Symbolic, lu->Control);
2671677d0b8SKris Buschelman 
2681677d0b8SKris Buschelman   lu->flg = DIFFERENT_NONZERO_PATTERN;
2691677d0b8SKris Buschelman   lu->av  = av;
270e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK = PETSC_TRUE;
271719d5645SBarry Smith   (F)->ops->lufactornumeric  = MatLUFactorNumeric_UMFPACK;
2721677d0b8SKris Buschelman   PetscFunctionReturn(0);
2731677d0b8SKris Buschelman }
2741677d0b8SKris Buschelman 
2751677d0b8SKris Buschelman #undef __FUNCT__
276f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_UMFPACK"
2774b019bd2SJed Brown static PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer)
2786849ba73SBarry Smith {
279f0c56d0fSKris Buschelman   Mat_UMFPACK    *lu= (Mat_UMFPACK*)A->spptr;
280dfbe8321SBarry Smith   PetscErrorCode ierr;
281f0c56d0fSKris Buschelman 
2821677d0b8SKris Buschelman   PetscFunctionBegin;
2831677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
284f0c56d0fSKris Buschelman   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
2851677d0b8SKris Buschelman 
2861677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
2871677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
2881677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
2891677d0b8SKris Buschelman 
2901677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
2910f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr);
2921677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
2931677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
2940f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr);
2951677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
2960f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_2BY2_TOLERANCE]: %g\n",lu->Control[UMFPACK_2BY2_TOLERANCE]);CHKERRQ(ierr);
2970f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr);
2980f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr);
2991677d0b8SKris Buschelman 
3001677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
3011677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
3020f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr);
3030f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr);
3041677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
3050f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr);
3061677d0b8SKris Buschelman 
3071677d0b8SKris Buschelman   /* Control parameters used by solve */
3081677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
3091677d0b8SKris Buschelman 
3101677d0b8SKris Buschelman   /* mat ordering */
3111677d0b8SKris Buschelman   if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer,"  UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr);
3121677d0b8SKris Buschelman 
3131677d0b8SKris Buschelman   PetscFunctionReturn(0);
3141677d0b8SKris Buschelman }
3151677d0b8SKris Buschelman 
316d844382dSKris Buschelman #undef __FUNCT__
317f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_UMFPACK"
3184b019bd2SJed Brown static PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
3196849ba73SBarry Smith {
320dfbe8321SBarry Smith   PetscErrorCode    ierr;
32132077d6dSBarry Smith   PetscTruth        iascii;
322d844382dSKris Buschelman   PetscViewerFormat format;
323d844382dSKris Buschelman 
324d844382dSKris Buschelman   PetscFunctionBegin;
325b24902e0SBarry Smith   ierr = MatView_SeqAIJ(A,viewer);CHKERRQ(ierr);
326d844382dSKris Buschelman 
3272692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
32832077d6dSBarry Smith   if (iascii) {
329d844382dSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
3302f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO) {
331f0c56d0fSKris Buschelman       ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
332d844382dSKris Buschelman     }
333d844382dSKris Buschelman   }
334d844382dSKris Buschelman   PetscFunctionReturn(0);
335d844382dSKris Buschelman }
336d844382dSKris Buschelman 
33735bd34faSBarry Smith EXTERN_C_BEGIN
33835bd34faSBarry Smith #undef __FUNCT__
33935bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_umfpack"
34035bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_seqaij_umfpack(Mat A,const MatSolverPackage *type)
34135bd34faSBarry Smith {
34235bd34faSBarry Smith   PetscFunctionBegin;
3432692d6eeSBarry Smith   *type = MATSOLVERUMFPACK;
34435bd34faSBarry Smith   PetscFunctionReturn(0);
34535bd34faSBarry Smith }
34635bd34faSBarry Smith EXTERN_C_END
34735bd34faSBarry Smith 
34835bd34faSBarry Smith 
3492f71e704SKris Buschelman /*MC
3502692d6eeSBarry Smith   MATSOLVERUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
3512f71e704SKris Buschelman   via the external package UMFPACK.
3522f71e704SKris Buschelman 
353e2e64c6bSBarry Smith   ./configure --download-umfpack to install PETSc to use UMFPACK
3542f71e704SKris Buschelman 
3552f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
3562f71e704SKris Buschelman   which correspond to the options database keys below.
3572f71e704SKris Buschelman 
3582f71e704SKris Buschelman   Options Database Keys:
3595c9eb25fSBarry Smith + -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL]
3603519fcdcSHong Zhang . -mat_umfpack_strategy <AUTO> (choose one of) AUTO UNSYMMETRIC SYMMETRIC 2BY2
3612f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
3623519fcdcSHong Zhang . -mat_umfpack_dense_row <0.2>: Control[UMFPACK_DENSE_ROW]
3633519fcdcSHong Zhang . -mat_umfpack_amd_dense <10>: Control[UMFPACK_AMD_DENSE]
3642f71e704SKris Buschelman . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
3653519fcdcSHong Zhang . -mat_umfpack_2by2_tolerance <0.01>: Control[UMFPACK_2BY2_TOLERANCE]
3663519fcdcSHong Zhang . -mat_umfpack_fixq <0>: Control[UMFPACK_FIXQ]
3673519fcdcSHong Zhang . -mat_umfpack_aggressive <1>: Control[UMFPACK_AGGRESSIVE]
3682f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
3693519fcdcSHong Zhang .  -mat_umfpack_sym_pivot_tolerance <0.001>: Control[UMFPACK_SYM_PIVOT_TOLERANCE]
3703519fcdcSHong Zhang .  -mat_umfpack_scale <NONE> (choose one of) NONE SUM MAX
3712f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
3723519fcdcSHong Zhang .  -mat_umfpack_droptol <0>: Control[UMFPACK_DROPTOL]
3732f71e704SKris Buschelman - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
3742f71e704SKris Buschelman 
3752f71e704SKris Buschelman    Level: beginner
3762f71e704SKris Buschelman 
3772692d6eeSBarry Smith .seealso: PCLU, MATSOLVERSUPERLU, MATSOLVERMUMPS, MATSOLVERSPOOLES, PCFactorSetMatSolverPackage(), MatSolverPackage
3782f71e704SKris Buschelman M*/
3793519fcdcSHong Zhang EXTERN_C_BEGIN
3803519fcdcSHong Zhang #undef __FUNCT__
3813519fcdcSHong Zhang #define __FUNCT__ "MatGetFactor_seqaij_umfpack"
3823519fcdcSHong Zhang PetscErrorCode MatGetFactor_seqaij_umfpack(Mat A,MatFactorType ftype,Mat *F)
3833519fcdcSHong Zhang {
3843519fcdcSHong Zhang   Mat            B;
3853519fcdcSHong Zhang   Mat_UMFPACK    *lu;
3863519fcdcSHong Zhang   PetscErrorCode ierr;
3878592901bSBarry Smith   PetscInt       m=A->rmap->n,n=A->cmap->n,idx;
3882f71e704SKris Buschelman 
3893519fcdcSHong Zhang   const char     *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC","2BY2"},
3903519fcdcSHong Zhang                  *scale[]={"NONE","SUM","MAX"};
3913519fcdcSHong Zhang   PetscTruth     flg;
3922d4e2982SHong Zhang 
3933519fcdcSHong Zhang   PetscFunctionBegin;
3943519fcdcSHong Zhang   /* Create the factorization matrix F */
3953519fcdcSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
3963519fcdcSHong Zhang   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
3973519fcdcSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
3983519fcdcSHong Zhang   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
3993519fcdcSHong Zhang   ierr = PetscNewLog(B,Mat_UMFPACK,&lu);CHKERRQ(ierr);
4003519fcdcSHong Zhang   B->spptr                 = lu;
4013519fcdcSHong Zhang   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
4023519fcdcSHong Zhang   B->ops->destroy          = MatDestroy_UMFPACK;
4033519fcdcSHong Zhang   B->ops->view             = MatView_UMFPACK;
40435bd34faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_umfpack",MatFactorGetSolverPackage_seqaij_umfpack);CHKERRQ(ierr);
405d5f3da31SBarry Smith   B->factortype            = MAT_FACTOR_LU;
4063519fcdcSHong Zhang   B->assembled             = PETSC_TRUE;  /* required by -ksp_view */
4073519fcdcSHong Zhang   B->preallocated          = PETSC_TRUE;
4083519fcdcSHong Zhang 
4093519fcdcSHong Zhang   /* initializations */
4103519fcdcSHong Zhang   /* ------------------------------------------------*/
4113519fcdcSHong Zhang   /* get the default control parameters */
4128592901bSBarry Smith   umfpack_UMF_defaults(lu->Control);
4133519fcdcSHong Zhang   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
4143519fcdcSHong Zhang   lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */
4153519fcdcSHong Zhang 
4163519fcdcSHong Zhang   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
4173519fcdcSHong Zhang   /* Control parameters used by reporting routiones */
4183519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr);
4193519fcdcSHong Zhang 
4203519fcdcSHong Zhang   /* Control parameters for symbolic factorization */
4213519fcdcSHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,4,strategy[0],&idx,&flg);CHKERRQ(ierr);
4223519fcdcSHong Zhang   if (flg) {
4233519fcdcSHong Zhang     switch (idx){
4243519fcdcSHong Zhang     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
4253519fcdcSHong Zhang     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
4263519fcdcSHong Zhang     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
4273519fcdcSHong Zhang     case 3: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_2BY2; break;
4283519fcdcSHong Zhang     }
4293519fcdcSHong Zhang   }
4303519fcdcSHong 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);
4313519fcdcSHong 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);
4323519fcdcSHong 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);
4333519fcdcSHong 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);
4343519fcdcSHong 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);
4353519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);CHKERRQ(ierr);
4363519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);CHKERRQ(ierr);
4373519fcdcSHong Zhang 
4383519fcdcSHong Zhang   /* Control parameters used by numeric factorization */
4393519fcdcSHong 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);
4403519fcdcSHong 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);
4413519fcdcSHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
4423519fcdcSHong Zhang   if (flg) {
4433519fcdcSHong Zhang     switch (idx){
4443519fcdcSHong Zhang     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
4453519fcdcSHong Zhang     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
4463519fcdcSHong Zhang     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
4473519fcdcSHong Zhang     }
4483519fcdcSHong Zhang   }
4493519fcdcSHong 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);
4503519fcdcSHong 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);
4513519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);CHKERRQ(ierr);
4523519fcdcSHong Zhang 
4533519fcdcSHong Zhang   /* Control parameters used by solve */
4543519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr);
4553519fcdcSHong Zhang 
4563519fcdcSHong Zhang   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */
4573519fcdcSHong Zhang   ierr = PetscOptionsHasName(PETSC_NULL,"-pc_factor_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr);
4583519fcdcSHong Zhang   PetscOptionsEnd();
4593519fcdcSHong Zhang   *F = B;
4603519fcdcSHong Zhang   PetscFunctionReturn(0);
4613519fcdcSHong Zhang }
4623519fcdcSHong Zhang EXTERN_C_END
4632d4e2982SHong Zhang 
464