xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 3519fcdc3ccd010619df22db39613e3a1abe34d6)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
21677d0b8SKris Buschelman 
31677d0b8SKris Buschelman /*
42d4e2982SHong Zhang    Provides an interface to the UMFPACKv5.1 sparse solver
52d4e2982SHong Zhang 
62d4e2982SHong Zhang    This interface uses the "UF_long version" of the UMFPACK API
72d4e2982SHong Zhang    (*_dl_* and *_zl_* routines, instead of *_di_* and *_zi_* routines)
82d4e2982SHong Zhang    so that UMFPACK can address more than 2Gb of memory on 64 bit
92d4e2982SHong Zhang    machines.
102d4e2982SHong Zhang 
112d4e2982SHong Zhang    If sizeof(UF_long) == 32 the interface only allocates the memory
122d4e2982SHong Zhang    necessary for UMFPACK's working arrays (W, Wi and possibly
132d4e2982SHong Zhang    perm_c). If sizeof(UF_long) == 64, in addition to allocating the
142d4e2982SHong Zhang    working arrays, the interface also has to re-allocate the matrix
152d4e2982SHong Zhang    index arrays (ai and aj, which must be stored as UF_long).
162d4e2982SHong Zhang 
172d4e2982SHong Zhang    The interface is implemented for both real and complex
182d4e2982SHong Zhang    arithmetic. Complex numbers are assumed to be in packed format,
192d4e2982SHong Zhang    which requires UMFPACK >= 4.4.
202d4e2982SHong Zhang 
212d4e2982SHong Zhang    We thank Christophe Geuzaine <geuzaine@acm.caltech.edu> for upgrading this interface to the UMFPACKv5.1
221677d0b8SKris Buschelman */
231677d0b8SKris Buschelman 
241677d0b8SKris Buschelman #include "src/mat/impls/aij/seq/aij.h"
251677d0b8SKris Buschelman 
261677d0b8SKris Buschelman EXTERN_C_BEGIN
271677d0b8SKris Buschelman #include "umfpack.h"
281677d0b8SKris Buschelman EXTERN_C_END
291677d0b8SKris Buschelman 
301677d0b8SKris Buschelman typedef struct {
311677d0b8SKris Buschelman   void         *Symbolic, *Numeric;
321677d0b8SKris Buschelman   double       Info[UMFPACK_INFO], Control[UMFPACK_CONTROL],*W;
332d4e2982SHong Zhang   UF_long      *Wi,*ai,*aj,*perm_c;
341677d0b8SKris Buschelman   PetscScalar  *av;
351677d0b8SKris Buschelman   MatStructure flg;
361677d0b8SKris Buschelman   PetscTruth   PetscMatOdering;
371677d0b8SKris Buschelman 
381677d0b8SKris Buschelman   /* Flag to clean up UMFPACK objects during Destroy */
391677d0b8SKris Buschelman   PetscTruth CleanUpUMFPACK;
40f0c56d0fSKris Buschelman } Mat_UMFPACK;
41f0c56d0fSKris Buschelman 
421677d0b8SKris Buschelman #undef __FUNCT__
43f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_UMFPACK"
44dfbe8321SBarry Smith PetscErrorCode MatDestroy_UMFPACK(Mat A)
45dfbe8321SBarry Smith {
46dfbe8321SBarry Smith   PetscErrorCode ierr;
47f0c56d0fSKris Buschelman   Mat_UMFPACK    *lu=(Mat_UMFPACK*)A->spptr;
481677d0b8SKris Buschelman 
491677d0b8SKris Buschelman   PetscFunctionBegin;
50fb731535SKris Buschelman   if (lu->CleanUpUMFPACK) {
512d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
522d4e2982SHong Zhang     umfpack_zl_free_symbolic(&lu->Symbolic);
532d4e2982SHong Zhang     umfpack_zl_free_numeric(&lu->Numeric);
542d4e2982SHong Zhang #else
552d4e2982SHong Zhang     umfpack_dl_free_symbolic(&lu->Symbolic);
562d4e2982SHong Zhang     umfpack_dl_free_numeric(&lu->Numeric);
572d4e2982SHong Zhang #endif
581677d0b8SKris Buschelman     ierr = PetscFree(lu->Wi);CHKERRQ(ierr);
591677d0b8SKris Buschelman     ierr = PetscFree(lu->W);CHKERRQ(ierr);
602d4e2982SHong Zhang     if(sizeof(UF_long) != sizeof(int)){
612d4e2982SHong Zhang       ierr = PetscFree(lu->ai);CHKERRQ(ierr);
622d4e2982SHong Zhang       ierr = PetscFree(lu->aj);CHKERRQ(ierr);
632d4e2982SHong Zhang     }
641677d0b8SKris Buschelman     if (lu->PetscMatOdering) {
651677d0b8SKris Buschelman       ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
661677d0b8SKris Buschelman     }
671677d0b8SKris Buschelman   }
68b24902e0SBarry Smith   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
691677d0b8SKris Buschelman   PetscFunctionReturn(0);
701677d0b8SKris Buschelman }
711677d0b8SKris Buschelman 
721677d0b8SKris Buschelman #undef __FUNCT__
73f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_UMFPACK"
746849ba73SBarry Smith PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x)
756849ba73SBarry Smith {
76f0c56d0fSKris Buschelman   Mat_UMFPACK *lu = (Mat_UMFPACK*)A->spptr;
771677d0b8SKris Buschelman   PetscScalar *av=lu->av,*ba,*xa;
78dfbe8321SBarry Smith   PetscErrorCode ierr;
792d4e2982SHong Zhang   UF_long     *ai=lu->ai,*aj=lu->aj,status;
801677d0b8SKris Buschelman 
811677d0b8SKris Buschelman   PetscFunctionBegin;
822d4e2982SHong Zhang   /* solve Ax = b by umfpack_*_wsolve */
831677d0b8SKris Buschelman   /* ----------------------------------*/
842d4e2982SHong Zhang 
852d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
862d4e2982SHong Zhang   ierr = VecConjugate(b);
872d4e2982SHong Zhang #endif
882d4e2982SHong Zhang 
891677d0b8SKris Buschelman   ierr = VecGetArray(b,&ba);
901677d0b8SKris Buschelman   ierr = VecGetArray(x,&xa);
911677d0b8SKris Buschelman 
922d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
932d4e2982SHong Zhang   status = umfpack_zl_wsolve(UMFPACK_At,ai,aj,(double*)av,NULL,(double*)xa,NULL,(double*)ba,NULL,
942d4e2982SHong Zhang 			     lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
952d4e2982SHong Zhang   umfpack_zl_report_info(lu->Control, lu->Info);
961677d0b8SKris Buschelman   if (status < 0){
972d4e2982SHong Zhang     umfpack_zl_report_status(lu->Control, status);
982d4e2982SHong Zhang     SETERRQ(PETSC_ERR_LIB,"umfpack_zl_wsolve failed");
991677d0b8SKris Buschelman   }
1002d4e2982SHong Zhang #else
1012d4e2982SHong Zhang   status = umfpack_dl_wsolve(UMFPACK_At,ai,aj,av,xa,ba,
1022d4e2982SHong Zhang 			     lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
1032d4e2982SHong Zhang   umfpack_dl_report_info(lu->Control, lu->Info);
1042d4e2982SHong Zhang   if (status < 0){
1052d4e2982SHong Zhang     umfpack_dl_report_status(lu->Control, status);
1062d4e2982SHong Zhang     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_wsolve failed");
1072d4e2982SHong Zhang   }
1082d4e2982SHong Zhang #endif
1091677d0b8SKris Buschelman 
1101677d0b8SKris Buschelman   ierr = VecRestoreArray(b,&ba);
1111677d0b8SKris Buschelman   ierr = VecRestoreArray(x,&xa);
1122a325a84SHong Zhang 
1132d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1142d4e2982SHong Zhang   ierr = VecConjugate(b);
1155b2a9f60SHong Zhang   ierr = VecConjugate(x);
1162d4e2982SHong Zhang #endif
1172d4e2982SHong Zhang 
1181677d0b8SKris Buschelman   PetscFunctionReturn(0);
1191677d0b8SKris Buschelman }
1201677d0b8SKris Buschelman 
1211677d0b8SKris Buschelman #undef __FUNCT__
122f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_UMFPACK"
123af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat A,MatFactorInfo *info,Mat *F)
1246849ba73SBarry Smith {
125f0c56d0fSKris Buschelman   Mat_UMFPACK *lu=(Mat_UMFPACK*)(*F)->spptr;
1266849ba73SBarry Smith   PetscErrorCode ierr;
127d0f46423SBarry Smith   UF_long     *ai=lu->ai,*aj=lu->aj,m=A->rmap->n,status;
1281677d0b8SKris Buschelman   PetscScalar *av=lu->av;
1291677d0b8SKris Buschelman 
1301677d0b8SKris Buschelman   PetscFunctionBegin;
1311677d0b8SKris Buschelman   /* numeric factorization of A' */
1321677d0b8SKris Buschelman   /* ----------------------------*/
1332d4e2982SHong Zhang 
1342d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1352a325a84SHong Zhang   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
1362d4e2982SHong Zhang     umfpack_zl_free_numeric(&lu->Numeric);
1372a325a84SHong Zhang   }
1382d4e2982SHong Zhang   status = umfpack_zl_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1399f42a82aSMatthew Knepley   if (status < 0) {
1409f42a82aSMatthew Knepley     umfpack_zl_report_status(lu->Control, status);
1419f42a82aSMatthew Knepley     SETERRQ(PETSC_ERR_LIB,"umfpack_zl_numeric failed");
1429f42a82aSMatthew Knepley   }
1431677d0b8SKris Buschelman   /* report numeric factorization of A' when Control[PRL] > 3 */
1442d4e2982SHong Zhang   (void) umfpack_zl_report_numeric(lu->Numeric, lu->Control);
1452d4e2982SHong Zhang #else
1462d4e2982SHong Zhang   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric){
1472d4e2982SHong Zhang     umfpack_dl_free_numeric(&lu->Numeric);
1482d4e2982SHong Zhang   }
1492d4e2982SHong Zhang   status = umfpack_dl_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1509f42a82aSMatthew Knepley   if (status < 0) {
1519f42a82aSMatthew Knepley     umfpack_zl_report_status(lu->Control, status);
1529f42a82aSMatthew Knepley     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_numeric failed");
1539f42a82aSMatthew Knepley   }
1542d4e2982SHong Zhang   /* report numeric factorization of A' when Control[PRL] > 3 */
1552d4e2982SHong Zhang   (void) umfpack_dl_report_numeric(lu->Numeric, lu->Control);
1562d4e2982SHong Zhang #endif
1571677d0b8SKris Buschelman 
1581677d0b8SKris Buschelman   if (lu->flg == DIFFERENT_NONZERO_PATTERN){  /* first numeric factorization */
1591677d0b8SKris Buschelman     /* allocate working space to be used by Solve */
1602d4e2982SHong Zhang     ierr = PetscMalloc(m * sizeof(UF_long), &lu->Wi);CHKERRQ(ierr);
1612d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1622d4e2982SHong Zhang     ierr = PetscMalloc(2*5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
1632d4e2982SHong Zhang #else
1641677d0b8SKris Buschelman     ierr = PetscMalloc(5*m * sizeof(double), &lu->W);CHKERRQ(ierr);
1652d4e2982SHong Zhang #endif
1661677d0b8SKris Buschelman   }
1672d4e2982SHong Zhang 
1682a325a84SHong Zhang   lu->flg = SAME_NONZERO_PATTERN;
1692a325a84SHong Zhang   lu->CleanUpUMFPACK = PETSC_TRUE;
1701677d0b8SKris Buschelman   PetscFunctionReturn(0);
1711677d0b8SKris Buschelman }
1721677d0b8SKris Buschelman 
1731677d0b8SKris Buschelman /*
1741677d0b8SKris Buschelman    Note the r permutation is ignored
1751677d0b8SKris Buschelman */
1761677d0b8SKris Buschelman #undef __FUNCT__
177f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_UMFPACK"
1786849ba73SBarry Smith PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
1796849ba73SBarry Smith {
1801677d0b8SKris Buschelman   Mat_SeqAIJ     *mat=(Mat_SeqAIJ*)A->data;
181b24902e0SBarry Smith   Mat_UMFPACK    *lu = (Mat_UMFPACK*)((*F)->spptr);
182dfbe8321SBarry Smith   PetscErrorCode ierr;
183d0f46423SBarry Smith   int            i,m=A->rmap->n,n=A->cmap->n,*ra;
1842d4e2982SHong Zhang   UF_long        status;
1851677d0b8SKris Buschelman   PetscScalar    *av=mat->a;
1861677d0b8SKris Buschelman 
1871677d0b8SKris Buschelman   PetscFunctionBegin;
1881677d0b8SKris Buschelman   if (lu->PetscMatOdering) {
1890f6efc10SHong Zhang     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
1902d4e2982SHong Zhang     ierr = PetscMalloc(m*sizeof(UF_long),&lu->perm_c);CHKERRQ(ierr);
1912d4e2982SHong Zhang     /* we cannot simply memcpy on 64 bit archs */
1922d4e2982SHong Zhang     for(i = 0; i < m; i++) lu->perm_c[i] = ra[i];
1930f6efc10SHong Zhang     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
1941677d0b8SKris Buschelman   }
1951677d0b8SKris Buschelman 
1962d4e2982SHong Zhang   if(sizeof(UF_long) != sizeof(int)){
1972d4e2982SHong Zhang     /* we cannot directly use mat->i and mat->j on 64 bit archs */
1982d4e2982SHong Zhang     ierr = PetscMalloc((m+1)*sizeof(UF_long),&lu->ai);CHKERRQ(ierr);
1992d4e2982SHong Zhang     ierr = PetscMalloc(mat->nz*sizeof(UF_long),&lu->aj);CHKERRQ(ierr);
2002d4e2982SHong Zhang     for(i = 0; i < m + 1; i++) lu->ai[i] = mat->i[i];
2012d4e2982SHong Zhang     for(i = 0; i < mat->nz; i++) lu->aj[i] = mat->j[i];
2022d4e2982SHong Zhang   }
2032d4e2982SHong Zhang   else{
2042d4e2982SHong Zhang     lu->ai = (UF_long*)mat->i;
2052d4e2982SHong Zhang     lu->aj = (UF_long*)mat->j;
2062d4e2982SHong Zhang   }
2072d4e2982SHong Zhang 
2081677d0b8SKris Buschelman   /* print the control parameters */
2092d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
2102d4e2982SHong Zhang   if(lu->Control[UMFPACK_PRL] > 1) umfpack_zl_report_control(lu->Control);
2112d4e2982SHong Zhang #else
2122d4e2982SHong Zhang   if(lu->Control[UMFPACK_PRL] > 1) umfpack_dl_report_control(lu->Control);
2132d4e2982SHong Zhang #endif
2141677d0b8SKris Buschelman 
2151677d0b8SKris Buschelman   /* symbolic factorization of A' */
2161677d0b8SKris Buschelman   /* ---------------------------------------------------------------------- */
2172d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
2180f6efc10SHong Zhang   if (lu->PetscMatOdering) { /* use Petsc row ordering */
2192d4e2982SHong Zhang     status = umfpack_zl_qsymbolic(n,m,lu->ai,lu->aj,(double*)av,NULL,
2202d4e2982SHong Zhang 				  lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2210f6efc10SHong Zhang   } else { /* use Umfpack col ordering */
2222d4e2982SHong Zhang     status = umfpack_zl_symbolic(n,m,lu->ai,lu->aj,(double*)av,NULL,
2232d4e2982SHong Zhang 				 &lu->Symbolic,lu->Control,lu->Info);
2240f6efc10SHong Zhang   }
2251677d0b8SKris Buschelman   if (status < 0){
2262d4e2982SHong Zhang     umfpack_zl_report_info(lu->Control, lu->Info);
2272d4e2982SHong Zhang     umfpack_zl_report_status(lu->Control, status);
2282d4e2982SHong Zhang     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_symbolic failed");
2291677d0b8SKris Buschelman   }
2301677d0b8SKris Buschelman   /* report sumbolic factorization of A' when Control[PRL] > 3 */
2312d4e2982SHong Zhang   (void) umfpack_zl_report_symbolic(lu->Symbolic, lu->Control);
2322d4e2982SHong Zhang #else
2332d4e2982SHong Zhang   if (lu->PetscMatOdering) { /* use Petsc row ordering */
2342d4e2982SHong Zhang     status = umfpack_dl_qsymbolic(n,m,lu->ai,lu->aj,av,
2352d4e2982SHong Zhang 				  lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2362d4e2982SHong Zhang   } else { /* use Umfpack col ordering */
2372d4e2982SHong Zhang     status = umfpack_dl_symbolic(n,m,lu->ai,lu->aj,av,
2382d4e2982SHong Zhang 				 &lu->Symbolic,lu->Control,lu->Info);
2392d4e2982SHong Zhang   }
2402d4e2982SHong Zhang   if (status < 0){
2412d4e2982SHong Zhang     umfpack_dl_report_info(lu->Control, lu->Info);
2422d4e2982SHong Zhang     umfpack_dl_report_status(lu->Control, status);
2432d4e2982SHong Zhang     SETERRQ(PETSC_ERR_LIB,"umfpack_dl_symbolic failed");
2442d4e2982SHong Zhang   }
2452d4e2982SHong Zhang   /* report sumbolic factorization of A' when Control[PRL] > 3 */
2462d4e2982SHong Zhang   (void) umfpack_dl_report_symbolic(lu->Symbolic, lu->Control);
2472d4e2982SHong Zhang #endif
2481677d0b8SKris Buschelman 
2491677d0b8SKris Buschelman   lu->flg = DIFFERENT_NONZERO_PATTERN;
2501677d0b8SKris Buschelman   lu->av  = av;
251e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK = PETSC_TRUE;
252db4efbfdSBarry Smith   (*F)->ops->lufactornumeric  = MatLUFactorNumeric_UMFPACK;
253db4efbfdSBarry Smith   (*F)->ops->solve            = MatSolve_UMFPACK;
2541677d0b8SKris Buschelman   PetscFunctionReturn(0);
2551677d0b8SKris Buschelman }
2561677d0b8SKris Buschelman 
2571677d0b8SKris Buschelman #undef __FUNCT__
258f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_UMFPACK"
2596849ba73SBarry Smith PetscErrorCode MatFactorInfo_UMFPACK(Mat A,PetscViewer viewer)
2606849ba73SBarry Smith {
261f0c56d0fSKris Buschelman   Mat_UMFPACK    *lu= (Mat_UMFPACK*)A->spptr;
262dfbe8321SBarry Smith   PetscErrorCode ierr;
263f0c56d0fSKris Buschelman 
2641677d0b8SKris Buschelman   PetscFunctionBegin;
2651677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
266f0c56d0fSKris Buschelman   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
2671677d0b8SKris Buschelman 
2681677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n");CHKERRQ(ierr);
2691677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
2701677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]);CHKERRQ(ierr);
2711677d0b8SKris Buschelman 
2721677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
2730f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]);CHKERRQ(ierr);
2741677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]);CHKERRQ(ierr);
2751677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]);CHKERRQ(ierr);
2760f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]);CHKERRQ(ierr);
2771677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]);CHKERRQ(ierr);
2780f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_2BY2_TOLERANCE]: %g\n",lu->Control[UMFPACK_2BY2_TOLERANCE]);CHKERRQ(ierr);
2790f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]);CHKERRQ(ierr);
2800f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]);CHKERRQ(ierr);
2811677d0b8SKris Buschelman 
2821677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
2831677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]);CHKERRQ(ierr);
2840f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]);CHKERRQ(ierr);
2850f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]);CHKERRQ(ierr);
2861677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]);CHKERRQ(ierr);
2870f6efc10SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]);CHKERRQ(ierr);
2881677d0b8SKris Buschelman 
2891677d0b8SKris Buschelman   /* Control parameters used by solve */
2901677d0b8SKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]);CHKERRQ(ierr);
2911677d0b8SKris Buschelman 
2921677d0b8SKris Buschelman   /* mat ordering */
2931677d0b8SKris Buschelman   if(!lu->PetscMatOdering) ierr = PetscViewerASCIIPrintf(viewer,"  UMFPACK default matrix ordering is used (not the PETSc matrix ordering) \n");CHKERRQ(ierr);
2941677d0b8SKris Buschelman 
2951677d0b8SKris Buschelman   PetscFunctionReturn(0);
2961677d0b8SKris Buschelman }
2971677d0b8SKris Buschelman 
298d844382dSKris Buschelman #undef __FUNCT__
299f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_UMFPACK"
3006849ba73SBarry Smith PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
3016849ba73SBarry Smith {
302dfbe8321SBarry Smith   PetscErrorCode    ierr;
30332077d6dSBarry Smith   PetscTruth        iascii;
304d844382dSKris Buschelman   PetscViewerFormat format;
305d844382dSKris Buschelman 
306d844382dSKris Buschelman   PetscFunctionBegin;
307b24902e0SBarry Smith   ierr = MatView_SeqAIJ(A,viewer);CHKERRQ(ierr);
308d844382dSKris Buschelman 
30932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
31032077d6dSBarry Smith   if (iascii) {
311d844382dSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
3122f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO) {
313f0c56d0fSKris Buschelman       ierr = MatFactorInfo_UMFPACK(A,viewer);CHKERRQ(ierr);
314d844382dSKris Buschelman     }
315d844382dSKris Buschelman   }
316d844382dSKris Buschelman   PetscFunctionReturn(0);
317d844382dSKris Buschelman }
318d844382dSKris Buschelman 
3192f71e704SKris Buschelman /*MC
3205c9eb25fSBarry Smith   MAT_SOLVER_UMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
3212f71e704SKris Buschelman   via the external package UMFPACK.
3222f71e704SKris Buschelman 
3235c9eb25fSBarry Smith   config/configure.py --download-umfpack to install PETSc to use UMFPACK
3242f71e704SKris Buschelman 
3252f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
3262f71e704SKris Buschelman   which correspond to the options database keys below.
3272f71e704SKris Buschelman 
3282f71e704SKris Buschelman   Options Database Keys:
3295c9eb25fSBarry Smith + -mat_umfpack_prl - UMFPACK print level: Control[UMFPACK_PRL]
330*3519fcdcSHong Zhang . -mat_umfpack_strategy <AUTO> (choose one of) AUTO UNSYMMETRIC SYMMETRIC 2BY2
3312f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c> - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
332*3519fcdcSHong Zhang . -mat_umfpack_dense_row <0.2>: Control[UMFPACK_DENSE_ROW]
333*3519fcdcSHong Zhang . -mat_umfpack_amd_dense <10>: Control[UMFPACK_AMD_DENSE]
3342f71e704SKris Buschelman . -mat_umfpack_block_size <bs> - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
335*3519fcdcSHong Zhang . -mat_umfpack_2by2_tolerance <0.01>: Control[UMFPACK_2BY2_TOLERANCE]
336*3519fcdcSHong Zhang . -mat_umfpack_fixq <0>: Control[UMFPACK_FIXQ]
337*3519fcdcSHong Zhang . -mat_umfpack_aggressive <1>: Control[UMFPACK_AGGRESSIVE]
3382f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
339*3519fcdcSHong Zhang .  -mat_umfpack_sym_pivot_tolerance <0.001>: Control[UMFPACK_SYM_PIVOT_TOLERANCE]
340*3519fcdcSHong Zhang .  -mat_umfpack_scale <NONE> (choose one of) NONE SUM MAX
3412f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta> - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
342*3519fcdcSHong Zhang .  -mat_umfpack_droptol <0>: Control[UMFPACK_DROPTOL]
3432f71e704SKris Buschelman - -mat_umfpack_irstep <maxit> - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
3442f71e704SKris Buschelman 
3452f71e704SKris Buschelman    Level: beginner
3462f71e704SKris Buschelman 
347*3519fcdcSHong Zhang .seealso: PCLU, MAT_SOLVER_SUPERLU, MAT_SOLVER_MUMPS, MAT_SOLVER_SPOOLES
3482f71e704SKris Buschelman M*/
349*3519fcdcSHong Zhang EXTERN_C_BEGIN
350*3519fcdcSHong Zhang #undef __FUNCT__
351*3519fcdcSHong Zhang #define __FUNCT__ "MatGetFactor_seqaij_umfpack"
352*3519fcdcSHong Zhang PetscErrorCode MatGetFactor_seqaij_umfpack(Mat A,MatFactorType ftype,Mat *F)
353*3519fcdcSHong Zhang {
354*3519fcdcSHong Zhang   Mat            B;
355*3519fcdcSHong Zhang   Mat_UMFPACK    *lu;
356*3519fcdcSHong Zhang   PetscErrorCode ierr;
357*3519fcdcSHong Zhang   int            m=A->rmap->n,n=A->cmap->n,idx;
3582f71e704SKris Buschelman 
359*3519fcdcSHong Zhang   const char     *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC","2BY2"},
360*3519fcdcSHong Zhang                  *scale[]={"NONE","SUM","MAX"};
361*3519fcdcSHong Zhang   PetscTruth     flg;
3622d4e2982SHong Zhang 
363*3519fcdcSHong Zhang   PetscFunctionBegin;
364*3519fcdcSHong Zhang   /* Create the factorization matrix F */
365*3519fcdcSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
366*3519fcdcSHong Zhang   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
367*3519fcdcSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
368*3519fcdcSHong Zhang   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
369*3519fcdcSHong Zhang   ierr = PetscNewLog(B,Mat_UMFPACK,&lu);CHKERRQ(ierr);
370*3519fcdcSHong Zhang   B->spptr                 = lu;
371*3519fcdcSHong Zhang   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
372*3519fcdcSHong Zhang   B->ops->destroy          = MatDestroy_UMFPACK;
373*3519fcdcSHong Zhang   B->ops->view             = MatView_UMFPACK;
374*3519fcdcSHong Zhang   B->factor                = MAT_FACTOR_LU;
375*3519fcdcSHong Zhang   B->assembled             = PETSC_TRUE;  /* required by -ksp_view */
376*3519fcdcSHong Zhang   B->preallocated          = PETSC_TRUE;
377*3519fcdcSHong Zhang 
378*3519fcdcSHong Zhang   /* initializations */
379*3519fcdcSHong Zhang   /* ------------------------------------------------*/
380*3519fcdcSHong Zhang   /* get the default control parameters */
381*3519fcdcSHong Zhang #if defined(PETSC_USE_COMPLEX)
382*3519fcdcSHong Zhang   umfpack_zl_defaults(lu->Control);
383*3519fcdcSHong Zhang #else
384*3519fcdcSHong Zhang   umfpack_dl_defaults(lu->Control);
385*3519fcdcSHong Zhang #endif
386*3519fcdcSHong Zhang   lu->perm_c = PETSC_NULL;  /* use defaul UMFPACK col permutation */
387*3519fcdcSHong Zhang   lu->Control[UMFPACK_IRSTEP] = 0; /* max num of iterative refinement steps to attempt */
388*3519fcdcSHong Zhang 
389*3519fcdcSHong Zhang   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"UMFPACK Options","Mat");CHKERRQ(ierr);
390*3519fcdcSHong Zhang   /* Control parameters used by reporting routiones */
391*3519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],PETSC_NULL);CHKERRQ(ierr);
392*3519fcdcSHong Zhang 
393*3519fcdcSHong Zhang   /* Control parameters for symbolic factorization */
394*3519fcdcSHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,4,strategy[0],&idx,&flg);CHKERRQ(ierr);
395*3519fcdcSHong Zhang   if (flg) {
396*3519fcdcSHong Zhang     switch (idx){
397*3519fcdcSHong Zhang     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
398*3519fcdcSHong Zhang     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
399*3519fcdcSHong Zhang     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
400*3519fcdcSHong Zhang     case 3: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_2BY2; break;
401*3519fcdcSHong Zhang     }
402*3519fcdcSHong Zhang   }
403*3519fcdcSHong 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);
404*3519fcdcSHong 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);
405*3519fcdcSHong 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);
406*3519fcdcSHong 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);
407*3519fcdcSHong 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);
408*3519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],PETSC_NULL);CHKERRQ(ierr);
409*3519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],PETSC_NULL);CHKERRQ(ierr);
410*3519fcdcSHong Zhang 
411*3519fcdcSHong Zhang   /* Control parameters used by numeric factorization */
412*3519fcdcSHong 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);
413*3519fcdcSHong 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);
414*3519fcdcSHong Zhang   ierr = PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
415*3519fcdcSHong Zhang   if (flg) {
416*3519fcdcSHong Zhang     switch (idx){
417*3519fcdcSHong Zhang     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
418*3519fcdcSHong Zhang     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
419*3519fcdcSHong Zhang     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
420*3519fcdcSHong Zhang     }
421*3519fcdcSHong Zhang   }
422*3519fcdcSHong 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);
423*3519fcdcSHong 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);
424*3519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],PETSC_NULL);CHKERRQ(ierr);
425*3519fcdcSHong Zhang 
426*3519fcdcSHong Zhang   /* Control parameters used by solve */
427*3519fcdcSHong Zhang   ierr = PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],PETSC_NULL);CHKERRQ(ierr);
428*3519fcdcSHong Zhang 
429*3519fcdcSHong Zhang   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Umfpack perm_c) */
430*3519fcdcSHong Zhang   ierr = PetscOptionsHasName(PETSC_NULL,"-pc_factor_mat_ordering_type",&lu->PetscMatOdering);CHKERRQ(ierr);
431*3519fcdcSHong Zhang   PetscOptionsEnd();
432*3519fcdcSHong Zhang   *F = B;
433*3519fcdcSHong Zhang   PetscFunctionReturn(0);
434*3519fcdcSHong Zhang }
435*3519fcdcSHong Zhang EXTERN_C_END
4362d4e2982SHong Zhang 
437