xref: /petsc/src/mat/impls/aij/seq/umfpack/umfpack.c (revision 26cc229b24afffee8e9a203ddc8754b40d660fc7)
11677d0b8SKris Buschelman 
21677d0b8SKris Buschelman /*
3d515b9b4SShri Abhyankar    Provides an interface to the UMFPACK sparse solver available through SuiteSparse version 4.2.1
42d4e2982SHong Zhang 
59e475b0dSSatish Balay    When build with PETSC_USE_64BIT_INDICES this will use Suitesparse_long as the
68592901bSBarry Smith    integer type in UMFPACK, otherwise it will use int. This means
78592901bSBarry Smith    all integers in this file as simply declared as PetscInt. Also it means
89e475b0dSSatish Balay    that one cannot use 64BIT_INDICES on 32bit machines [as Suitesparse_long is 32bit only]
92d4e2982SHong Zhang 
101677d0b8SKris Buschelman */
11c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
121677d0b8SKris Buschelman 
138592901bSBarry Smith #if defined(PETSC_USE_64BIT_INDICES)
148592901bSBarry Smith #if defined(PETSC_USE_COMPLEX)
158592901bSBarry Smith #define umfpack_UMF_free_symbolic                      umfpack_zl_free_symbolic
168592901bSBarry Smith #define umfpack_UMF_free_numeric                       umfpack_zl_free_numeric
17d755ee67SBarry Smith /* the type casts are needed because PetscInt is long long while SuiteSparse_long is long and compilers warn even when they are identical */
18d755ee67SBarry Smith #define umfpack_UMF_wsolve(a,b,c,d,e,f,g,h,i,j,k,l,m,n) umfpack_zl_wsolve(a,(SuiteSparse_long*)b,(SuiteSparse_long*)c,d,e,f,g,h,i,(SuiteSparse_long*)j,k,l,(SuiteSparse_long*)m,n)
19d755ee67SBarry Smith #define umfpack_UMF_numeric(a,b,c,d,e,f,g,h)          umfpack_zl_numeric((SuiteSparse_long*)a,(SuiteSparse_long*)b,c,d,e,f,g,h)
208592901bSBarry Smith #define umfpack_UMF_report_numeric                    umfpack_zl_report_numeric
218592901bSBarry Smith #define umfpack_UMF_report_control                    umfpack_zl_report_control
228592901bSBarry Smith #define umfpack_UMF_report_status                     umfpack_zl_report_status
238592901bSBarry Smith #define umfpack_UMF_report_info                       umfpack_zl_report_info
248592901bSBarry Smith #define umfpack_UMF_report_symbolic                   umfpack_zl_report_symbolic
25d755ee67SBarry Smith #define umfpack_UMF_qsymbolic(a,b,c,d,e,f,g,h,i,j)    umfpack_zl_qsymbolic(a,b,(SuiteSparse_long*)c,(SuiteSparse_long*)d,e,f,(SuiteSparse_long*)g,h,i,j)
26d755ee67SBarry Smith #define umfpack_UMF_symbolic(a,b,c,d,e,f,g,h,i)       umfpack_zl_symbolic(a,b,(SuiteSparse_long*)c,(SuiteSparse_long*)d,e,f,g,h,i)
278592901bSBarry Smith #define umfpack_UMF_defaults                          umfpack_zl_defaults
288592901bSBarry Smith 
298592901bSBarry Smith #else
308592901bSBarry Smith #define umfpack_UMF_free_symbolic                  umfpack_dl_free_symbolic
318592901bSBarry Smith #define umfpack_UMF_free_numeric                   umfpack_dl_free_numeric
32d755ee67SBarry Smith #define umfpack_UMF_wsolve(a,b,c,d,e,f,g,h,i,j,k)  umfpack_dl_wsolve(a,(SuiteSparse_long*)b,(SuiteSparse_long*)c,d,e,f,g,h,i,(SuiteSparse_long*)j,k)
33d755ee67SBarry Smith #define umfpack_UMF_numeric(a,b,c,d,e,f,g)         umfpack_dl_numeric((SuiteSparse_long*)a,(SuiteSparse_long*)b,c,d,e,f,g)
348592901bSBarry Smith #define umfpack_UMF_report_numeric                 umfpack_dl_report_numeric
358592901bSBarry Smith #define umfpack_UMF_report_control                 umfpack_dl_report_control
368592901bSBarry Smith #define umfpack_UMF_report_status                  umfpack_dl_report_status
378592901bSBarry Smith #define umfpack_UMF_report_info                    umfpack_dl_report_info
388592901bSBarry Smith #define umfpack_UMF_report_symbolic                umfpack_dl_report_symbolic
39d755ee67SBarry Smith #define umfpack_UMF_qsymbolic(a,b,c,d,e,f,g,h,i)   umfpack_dl_qsymbolic(a,b,(SuiteSparse_long*)c,(SuiteSparse_long*)d,e,(SuiteSparse_long*)f,g,h,i)
40d755ee67SBarry Smith #define umfpack_UMF_symbolic(a,b,c,d,e,f,g,h)      umfpack_dl_symbolic(a,b,(SuiteSparse_long*)c,(SuiteSparse_long*)d,e,f,g,h)
418592901bSBarry Smith #define umfpack_UMF_defaults                       umfpack_dl_defaults
428592901bSBarry Smith #endif
438592901bSBarry Smith 
448592901bSBarry Smith #else
458592901bSBarry Smith #if defined(PETSC_USE_COMPLEX)
468592901bSBarry Smith #define umfpack_UMF_free_symbolic   umfpack_zi_free_symbolic
478592901bSBarry Smith #define umfpack_UMF_free_numeric    umfpack_zi_free_numeric
488592901bSBarry Smith #define umfpack_UMF_wsolve          umfpack_zi_wsolve
498592901bSBarry Smith #define umfpack_UMF_numeric         umfpack_zi_numeric
508592901bSBarry Smith #define umfpack_UMF_report_numeric  umfpack_zi_report_numeric
518592901bSBarry Smith #define umfpack_UMF_report_control  umfpack_zi_report_control
528592901bSBarry Smith #define umfpack_UMF_report_status   umfpack_zi_report_status
538592901bSBarry Smith #define umfpack_UMF_report_info     umfpack_zi_report_info
548592901bSBarry Smith #define umfpack_UMF_report_symbolic umfpack_zi_report_symbolic
558592901bSBarry Smith #define umfpack_UMF_qsymbolic       umfpack_zi_qsymbolic
568592901bSBarry Smith #define umfpack_UMF_symbolic        umfpack_zi_symbolic
578592901bSBarry Smith #define umfpack_UMF_defaults        umfpack_zi_defaults
588592901bSBarry Smith 
598592901bSBarry Smith #else
608592901bSBarry Smith #define umfpack_UMF_free_symbolic   umfpack_di_free_symbolic
618592901bSBarry Smith #define umfpack_UMF_free_numeric    umfpack_di_free_numeric
628592901bSBarry Smith #define umfpack_UMF_wsolve          umfpack_di_wsolve
638592901bSBarry Smith #define umfpack_UMF_numeric         umfpack_di_numeric
648592901bSBarry Smith #define umfpack_UMF_report_numeric  umfpack_di_report_numeric
658592901bSBarry Smith #define umfpack_UMF_report_control  umfpack_di_report_control
668592901bSBarry Smith #define umfpack_UMF_report_status   umfpack_di_report_status
678592901bSBarry Smith #define umfpack_UMF_report_info     umfpack_di_report_info
688592901bSBarry Smith #define umfpack_UMF_report_symbolic umfpack_di_report_symbolic
698592901bSBarry Smith #define umfpack_UMF_qsymbolic       umfpack_di_qsymbolic
708592901bSBarry Smith #define umfpack_UMF_symbolic        umfpack_di_symbolic
718592901bSBarry Smith #define umfpack_UMF_defaults        umfpack_di_defaults
728592901bSBarry Smith #endif
738592901bSBarry Smith #endif
748592901bSBarry Smith 
751677d0b8SKris Buschelman EXTERN_C_BEGIN
76c6db04a5SJed Brown #include <umfpack.h>
771677d0b8SKris Buschelman EXTERN_C_END
781677d0b8SKris Buschelman 
79fcfd50ebSBarry Smith static const char *const UmfpackOrderingTypes[] = {"CHOLMOD","AMD","GIVEN","METIS","BEST","NONE","USER","UmfpackOrderingTypes","UMFPACK_ORDERING_",0};
80fcfd50ebSBarry Smith 
811677d0b8SKris Buschelman typedef struct {
821677d0b8SKris Buschelman   void         *Symbolic, *Numeric;
831677d0b8SKris Buschelman   double       Info[UMFPACK_INFO], Control[UMFPACK_CONTROL],*W;
84ae26a541SJed Brown   PetscInt     *Wi,*perm_c;
85ae26a541SJed Brown   Mat          A;               /* Matrix used for factorization */
861677d0b8SKris Buschelman   MatStructure flg;
871677d0b8SKris Buschelman 
881677d0b8SKris Buschelman   /* Flag to clean up UMFPACK objects during Destroy */
89ace3abfcSBarry Smith   PetscBool CleanUpUMFPACK;
90f0c56d0fSKris Buschelman } Mat_UMFPACK;
91f0c56d0fSKris Buschelman 
924b019bd2SJed Brown static PetscErrorCode MatDestroy_UMFPACK(Mat A)
93dfbe8321SBarry Smith {
94b9c12af5SBarry Smith   Mat_UMFPACK    *lu=(Mat_UMFPACK*)A->data;
951677d0b8SKris Buschelman 
961677d0b8SKris Buschelman   PetscFunctionBegin;
97b9c12af5SBarry Smith   if (lu->CleanUpUMFPACK) {
988592901bSBarry Smith     umfpack_UMF_free_symbolic(&lu->Symbolic);
998592901bSBarry Smith     umfpack_UMF_free_numeric(&lu->Numeric);
1009566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->Wi));
1019566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->W));
1029566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->perm_c));
1031677d0b8SKris Buschelman   }
1049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lu->A));
1052e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL));
1069566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
1071677d0b8SKris Buschelman   PetscFunctionReturn(0);
1081677d0b8SKris Buschelman }
1091677d0b8SKris Buschelman 
1104b019bd2SJed Brown static PetscErrorCode MatSolve_UMFPACK_Private(Mat A,Vec b,Vec x,int uflag)
1116849ba73SBarry Smith {
112b9c12af5SBarry Smith   Mat_UMFPACK       *lu = (Mat_UMFPACK*)A->data;
113ae26a541SJed Brown   Mat_SeqAIJ        *a  = (Mat_SeqAIJ*)lu->A->data;
114d9ca1df4SBarry Smith   PetscScalar       *av = a->a,*xa;
115d9ca1df4SBarry Smith   const PetscScalar *ba;
116ae26a541SJed Brown   PetscInt          *ai = a->i,*aj = a->j,status;
117fcf6e9abSJed Brown   static PetscBool  cite = PETSC_FALSE;
1181677d0b8SKris Buschelman 
1191677d0b8SKris Buschelman   PetscFunctionBegin;
120a643c07aSBarry Smith   if (!A->rmap->n) PetscFunctionReturn(0);
1219566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister("@article{davis2004algorithm,\n  title={Algorithm 832: {UMFPACK} V4.3---An Unsymmetric-Pattern Multifrontal Method},\n  author={Davis, Timothy A},\n  journal={ACM Transactions on Mathematical Software (TOMS)},\n  volume={30},\n  number={2},\n  pages={196--199},\n  year={2004},\n  publisher={ACM}\n}\n",&cite));
1222d4e2982SHong Zhang   /* solve Ax = b by umfpack_*_wsolve */
1231677d0b8SKris Buschelman   /* ----------------------------------*/
1242d4e2982SHong Zhang 
125ae26a541SJed Brown   if (!lu->Wi) {  /* first time, allocate working space for wsolve */
1269566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(A->rmap->n,&lu->Wi));
1279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5*A->rmap->n,&lu->W));
128ae26a541SJed Brown   }
129ae26a541SJed Brown 
1309566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(b,&ba));
1319566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x,&xa));
13279c34000SHong Zhang #if defined(PETSC_USE_COMPLEX)
133b0043f70SBarry Smith   status = umfpack_UMF_wsolve(uflag,ai,aj,(PetscReal*)av,NULL,(PetscReal*)xa,NULL,(PetscReal*)ba,NULL,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
13479c34000SHong Zhang #else
1354b019bd2SJed Brown   status = umfpack_UMF_wsolve(uflag,ai,aj,av,xa,ba,lu->Numeric,lu->Control,lu->Info,lu->Wi,lu->W);
13679c34000SHong Zhang #endif
1378592901bSBarry Smith   umfpack_UMF_report_info(lu->Control, lu->Info);
1381677d0b8SKris Buschelman   if (status < 0) {
1398592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
140e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_wsolve failed");
1411677d0b8SKris Buschelman   }
1421677d0b8SKris Buschelman 
1439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(b,&ba));
1449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x,&xa));
1454b019bd2SJed Brown   PetscFunctionReturn(0);
1464b019bd2SJed Brown }
1472a325a84SHong Zhang 
1484b019bd2SJed Brown static PetscErrorCode MatSolve_UMFPACK(Mat A,Vec b,Vec x)
1494b019bd2SJed Brown {
1504b019bd2SJed Brown   PetscFunctionBegin;
1514b019bd2SJed Brown   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
1529566063dSJacob Faibussowitsch   PetscCall(MatSolve_UMFPACK_Private(A,b,x,UMFPACK_Aat));
1534b019bd2SJed Brown   PetscFunctionReturn(0);
1544b019bd2SJed Brown }
1554b019bd2SJed Brown 
1564b019bd2SJed Brown static PetscErrorCode MatSolveTranspose_UMFPACK(Mat A,Vec b,Vec x)
1574b019bd2SJed Brown {
1584b019bd2SJed Brown   PetscFunctionBegin;
1594b019bd2SJed Brown   /* We gave UMFPACK the algebraic transpose (because it assumes column alignment) */
1609566063dSJacob Faibussowitsch   PetscCall(MatSolve_UMFPACK_Private(A,b,x,UMFPACK_A));
1611677d0b8SKris Buschelman   PetscFunctionReturn(0);
1621677d0b8SKris Buschelman }
1631677d0b8SKris Buschelman 
1644b019bd2SJed Brown static PetscErrorCode MatLUFactorNumeric_UMFPACK(Mat F,Mat A,const MatFactorInfo *info)
1656849ba73SBarry Smith {
166b9c12af5SBarry Smith   Mat_UMFPACK    *lu = (Mat_UMFPACK*)(F)->data;
167ae26a541SJed Brown   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
168ae26a541SJed Brown   PetscInt       *ai = a->i,*aj=a->j,status;
169ae26a541SJed Brown   PetscScalar    *av = a->a;
1701677d0b8SKris Buschelman 
1711677d0b8SKris Buschelman   PetscFunctionBegin;
172a643c07aSBarry Smith   if (!A->rmap->n) PetscFunctionReturn(0);
1731677d0b8SKris Buschelman   /* numeric factorization of A' */
1741677d0b8SKris Buschelman   /* ----------------------------*/
1752d4e2982SHong Zhang 
1768592901bSBarry Smith   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) {
1778592901bSBarry Smith     umfpack_UMF_free_numeric(&lu->Numeric);
1788592901bSBarry Smith   }
1792d4e2982SHong Zhang #if defined(PETSC_USE_COMPLEX)
1808592901bSBarry Smith   status = umfpack_UMF_numeric(ai,aj,(double*)av,NULL,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1812d4e2982SHong Zhang #else
1828592901bSBarry Smith   status = umfpack_UMF_numeric(ai,aj,av,lu->Symbolic,&lu->Numeric,lu->Control,lu->Info);
1838592901bSBarry Smith #endif
1849f42a82aSMatthew Knepley   if (status < 0) {
1858592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
186e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_numeric failed");
1879f42a82aSMatthew Knepley   }
1882d4e2982SHong Zhang   /* report numeric factorization of A' when Control[PRL] > 3 */
1898592901bSBarry Smith   (void) umfpack_UMF_report_numeric(lu->Numeric, lu->Control);
1901677d0b8SKris Buschelman 
1919566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)A));
1929566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&lu->A));
1932205254eSKarl Rupp 
194ae26a541SJed Brown   lu->A                  = A;
1952a325a84SHong Zhang   lu->flg                = SAME_NONZERO_PATTERN;
1962a325a84SHong Zhang   lu->CleanUpUMFPACK     = PETSC_TRUE;
1974b019bd2SJed Brown   F->ops->solve          = MatSolve_UMFPACK;
1984b019bd2SJed Brown   F->ops->solvetranspose = MatSolveTranspose_UMFPACK;
1991677d0b8SKris Buschelman   PetscFunctionReturn(0);
2001677d0b8SKris Buschelman }
2011677d0b8SKris Buschelman 
2024b019bd2SJed Brown static PetscErrorCode MatLUFactorSymbolic_UMFPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
2036849ba73SBarry Smith {
204ae26a541SJed Brown   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
205b9c12af5SBarry Smith   Mat_UMFPACK    *lu = (Mat_UMFPACK*)(F->data);
206*26cc229bSBarry Smith   PetscInt       i,*ai = a->i,*aj = a->j,m=A->rmap->n,n=A->cmap->n,status,idx;
207989213a4SSatish Balay #if !defined(PETSC_USE_COMPLEX)
208ae26a541SJed Brown   PetscScalar    *av = a->a;
209989213a4SSatish Balay #endif
2105d0c19d7SBarry Smith   const PetscInt *ra;
211*26cc229bSBarry Smith   const char     *strategy[]={"AUTO","UNSYMMETRIC","SYMMETRIC"};
212*26cc229bSBarry Smith   const char     *scale[]   ={"NONE","SUM","MAX"};
213*26cc229bSBarry Smith   PetscBool      flg;
2141677d0b8SKris Buschelman 
2151677d0b8SKris Buschelman   PetscFunctionBegin;
216a643c07aSBarry Smith   (F)->ops->lufactornumeric = MatLUFactorNumeric_UMFPACK;
217a643c07aSBarry Smith   if (!n) PetscFunctionReturn(0);
218*26cc229bSBarry Smith 
219*26cc229bSBarry Smith   /* Set options to F */
220*26cc229bSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)F),((PetscObject)F)->prefix,"UMFPACK Options","Mat");
221*26cc229bSBarry Smith   /* Control parameters used by reporting routiones */
222*26cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_prl","Control[UMFPACK_PRL]","None",lu->Control[UMFPACK_PRL],&lu->Control[UMFPACK_PRL],NULL));
223*26cc229bSBarry Smith 
224*26cc229bSBarry Smith   /* Control parameters for symbolic factorization */
225*26cc229bSBarry Smith   PetscCall(PetscOptionsEList("-mat_umfpack_strategy","ordering and pivoting strategy","None",strategy,3,strategy[0],&idx,&flg));
226*26cc229bSBarry Smith   if (flg) {
227*26cc229bSBarry Smith     switch (idx) {
228*26cc229bSBarry Smith     case 0: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_AUTO; break;
229*26cc229bSBarry Smith     case 1: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_UNSYMMETRIC; break;
230*26cc229bSBarry Smith     case 2: lu->Control[UMFPACK_STRATEGY] = UMFPACK_STRATEGY_SYMMETRIC; break;
231*26cc229bSBarry Smith     }
232*26cc229bSBarry Smith   }
233*26cc229bSBarry Smith   PetscCall(PetscOptionsEList("-mat_umfpack_ordering","Internal ordering method","None",UmfpackOrderingTypes,PETSC_STATIC_ARRAY_LENGTH(UmfpackOrderingTypes),UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]],&idx,&flg));
234*26cc229bSBarry Smith   if (flg) lu->Control[UMFPACK_ORDERING] = (int)idx;
235*26cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_dense_col","Control[UMFPACK_DENSE_COL]","None",lu->Control[UMFPACK_DENSE_COL],&lu->Control[UMFPACK_DENSE_COL],NULL));
236*26cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_dense_row","Control[UMFPACK_DENSE_ROW]","None",lu->Control[UMFPACK_DENSE_ROW],&lu->Control[UMFPACK_DENSE_ROW],NULL));
237*26cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_amd_dense","Control[UMFPACK_AMD_DENSE]","None",lu->Control[UMFPACK_AMD_DENSE],&lu->Control[UMFPACK_AMD_DENSE],NULL));
238*26cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_block_size","Control[UMFPACK_BLOCK_SIZE]","None",lu->Control[UMFPACK_BLOCK_SIZE],&lu->Control[UMFPACK_BLOCK_SIZE],NULL));
239*26cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_fixq","Control[UMFPACK_FIXQ]","None",lu->Control[UMFPACK_FIXQ],&lu->Control[UMFPACK_FIXQ],NULL));
240*26cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_aggressive","Control[UMFPACK_AGGRESSIVE]","None",lu->Control[UMFPACK_AGGRESSIVE],&lu->Control[UMFPACK_AGGRESSIVE],NULL));
241*26cc229bSBarry Smith 
242*26cc229bSBarry Smith   /* Control parameters used by numeric factorization */
243*26cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_pivot_tolerance","Control[UMFPACK_PIVOT_TOLERANCE]","None",lu->Control[UMFPACK_PIVOT_TOLERANCE],&lu->Control[UMFPACK_PIVOT_TOLERANCE],NULL));
244*26cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_sym_pivot_tolerance","Control[UMFPACK_SYM_PIVOT_TOLERANCE]","None",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE],&lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE],NULL));
245*26cc229bSBarry Smith   PetscCall(PetscOptionsEList("-mat_umfpack_scale","Control[UMFPACK_SCALE]","None",scale,3,scale[0],&idx,&flg));
246*26cc229bSBarry Smith   if (flg) {
247*26cc229bSBarry Smith     switch (idx) {
248*26cc229bSBarry Smith     case 0: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_NONE; break;
249*26cc229bSBarry Smith     case 1: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_SUM; break;
250*26cc229bSBarry Smith     case 2: lu->Control[UMFPACK_SCALE] = UMFPACK_SCALE_MAX; break;
251*26cc229bSBarry Smith     }
252*26cc229bSBarry Smith   }
253*26cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_alloc_init","Control[UMFPACK_ALLOC_INIT]","None",lu->Control[UMFPACK_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],NULL));
254*26cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_front_alloc_init","Control[UMFPACK_FRONT_ALLOC_INIT]","None",lu->Control[UMFPACK_FRONT_ALLOC_INIT],&lu->Control[UMFPACK_ALLOC_INIT],NULL));
255*26cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_droptol","Control[UMFPACK_DROPTOL]","None",lu->Control[UMFPACK_DROPTOL],&lu->Control[UMFPACK_DROPTOL],NULL));
256*26cc229bSBarry Smith 
257*26cc229bSBarry Smith   /* Control parameters used by solve */
258*26cc229bSBarry Smith   PetscCall(PetscOptionsReal("-mat_umfpack_irstep","Control[UMFPACK_IRSTEP]","None",lu->Control[UMFPACK_IRSTEP],&lu->Control[UMFPACK_IRSTEP],NULL));
259*26cc229bSBarry Smith   PetscOptionsEnd();
260*26cc229bSBarry Smith 
2612c7c0729SBarry Smith   if (r) {
2629566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(r,&ra));
2639566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m,&lu->perm_c));
2642d4e2982SHong Zhang     /* we cannot simply memcpy on 64 bit archs */
2652d4e2982SHong Zhang     for (i = 0; i < m; i++) lu->perm_c[i] = ra[i];
2669566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(r,&ra));
2671677d0b8SKris Buschelman   }
2681677d0b8SKris Buschelman 
2691677d0b8SKris Buschelman   /* print the control parameters */
2708592901bSBarry Smith   if (lu->Control[UMFPACK_PRL] > 1) umfpack_UMF_report_control(lu->Control);
2711677d0b8SKris Buschelman 
2721677d0b8SKris Buschelman   /* symbolic factorization of A' */
2731677d0b8SKris Buschelman   /* ---------------------------------------------------------------------- */
2742c7c0729SBarry Smith   if (r) { /* use Petsc row ordering */
2758592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
276ae26a541SJed Brown     status = umfpack_UMF_qsymbolic(n,m,ai,aj,av,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2772d4e2982SHong Zhang #else
278b0043f70SBarry Smith     status = umfpack_UMF_qsymbolic(n,m,ai,aj,NULL,NULL,lu->perm_c,&lu->Symbolic,lu->Control,lu->Info);
2798592901bSBarry Smith #endif
2802d4e2982SHong Zhang   } else { /* use Umfpack col ordering */
2818592901bSBarry Smith #if !defined(PETSC_USE_COMPLEX)
282ae26a541SJed Brown     status = umfpack_UMF_symbolic(n,m,ai,aj,av,&lu->Symbolic,lu->Control,lu->Info);
2838592901bSBarry Smith #else
284ae26a541SJed Brown     status = umfpack_UMF_symbolic(n,m,ai,aj,NULL,NULL,&lu->Symbolic,lu->Control,lu->Info);
2858592901bSBarry Smith #endif
2862d4e2982SHong Zhang   }
2872d4e2982SHong Zhang   if (status < 0) {
2888592901bSBarry Smith     umfpack_UMF_report_info(lu->Control, lu->Info);
2898592901bSBarry Smith     umfpack_UMF_report_status(lu->Control, status);
290e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"umfpack_UMF_symbolic failed");
2912d4e2982SHong Zhang   }
2922d4e2982SHong Zhang   /* report sumbolic factorization of A' when Control[PRL] > 3 */
2938592901bSBarry Smith   (void) umfpack_UMF_report_symbolic(lu->Symbolic, lu->Control);
2941677d0b8SKris Buschelman 
2951677d0b8SKris Buschelman   lu->flg                   = DIFFERENT_NONZERO_PATTERN;
296e1b4c3a1SKris Buschelman   lu->CleanUpUMFPACK        = PETSC_TRUE;
2971677d0b8SKris Buschelman   PetscFunctionReturn(0);
2981677d0b8SKris Buschelman }
2991677d0b8SKris Buschelman 
300860c79edSBarry Smith static PetscErrorCode MatView_Info_UMFPACK(Mat A,PetscViewer viewer)
3016849ba73SBarry Smith {
302b9c12af5SBarry Smith   Mat_UMFPACK    *lu= (Mat_UMFPACK*)A->data;
303f0c56d0fSKris Buschelman 
3041677d0b8SKris Buschelman   PetscFunctionBegin;
3051677d0b8SKris Buschelman   /* check if matrix is UMFPACK type */
306f0c56d0fSKris Buschelman   if (A->ops->solve != MatSolve_UMFPACK) PetscFunctionReturn(0);
3071677d0b8SKris Buschelman 
3089566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"UMFPACK run parameters:\n"));
3091677d0b8SKris Buschelman   /* Control parameters used by reporting routiones */
3109566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PRL]: %g\n",lu->Control[UMFPACK_PRL]));
3111677d0b8SKris Buschelman 
3121677d0b8SKris Buschelman   /* Control parameters used by symbolic factorization */
3139566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_STRATEGY]: %g\n",lu->Control[UMFPACK_STRATEGY]));
3149566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_COL]: %g\n",lu->Control[UMFPACK_DENSE_COL]));
3159566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DENSE_ROW]: %g\n",lu->Control[UMFPACK_DENSE_ROW]));
3169566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AMD_DENSE]: %g\n",lu->Control[UMFPACK_AMD_DENSE]));
3179566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_BLOCK_SIZE]: %g\n",lu->Control[UMFPACK_BLOCK_SIZE]));
3189566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_FIXQ]: %g\n",lu->Control[UMFPACK_FIXQ]));
3199566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_AGGRESSIVE]: %g\n",lu->Control[UMFPACK_AGGRESSIVE]));
3201677d0b8SKris Buschelman 
3211677d0b8SKris Buschelman   /* Control parameters used by numeric factorization */
3229566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_PIVOT_TOLERANCE]));
3239566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SYM_PIVOT_TOLERANCE]: %g\n",lu->Control[UMFPACK_SYM_PIVOT_TOLERANCE]));
3249566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_SCALE]: %g\n",lu->Control[UMFPACK_SCALE]));
3259566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ALLOC_INIT]: %g\n",lu->Control[UMFPACK_ALLOC_INIT]));
3269566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_DROPTOL]: %g\n",lu->Control[UMFPACK_DROPTOL]));
3271677d0b8SKris Buschelman 
3281677d0b8SKris Buschelman   /* Control parameters used by solve */
3299566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_IRSTEP]: %g\n",lu->Control[UMFPACK_IRSTEP]));
3301677d0b8SKris Buschelman 
3311677d0b8SKris Buschelman   /* mat ordering */
3322c7c0729SBarry Smith   if (!lu->perm_c) {
3339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Control[UMFPACK_ORDERING]: %s (not using the PETSc ordering)\n",UmfpackOrderingTypes[(int)lu->Control[UMFPACK_ORDERING]]));
334fcfd50ebSBarry Smith   }
3351677d0b8SKris Buschelman   PetscFunctionReturn(0);
3361677d0b8SKris Buschelman }
3371677d0b8SKris Buschelman 
3384b019bd2SJed Brown static PetscErrorCode MatView_UMFPACK(Mat A,PetscViewer viewer)
3396849ba73SBarry Smith {
340ace3abfcSBarry Smith   PetscBool         iascii;
341d844382dSKris Buschelman   PetscViewerFormat format;
342d844382dSKris Buschelman 
343d844382dSKris Buschelman   PetscFunctionBegin;
3449566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
34532077d6dSBarry Smith   if (iascii) {
3469566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer,&format));
3472f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO) {
3489566063dSJacob Faibussowitsch       PetscCall(MatView_Info_UMFPACK(A,viewer));
349d844382dSKris Buschelman     }
350d844382dSKris Buschelman   }
351d844382dSKris Buschelman   PetscFunctionReturn(0);
352d844382dSKris Buschelman }
353d844382dSKris Buschelman 
354ea799195SBarry Smith PetscErrorCode MatFactorGetSolverType_seqaij_umfpack(Mat A,MatSolverType *type)
35535bd34faSBarry Smith {
35635bd34faSBarry Smith   PetscFunctionBegin;
3572692d6eeSBarry Smith   *type = MATSOLVERUMFPACK;
35835bd34faSBarry Smith   PetscFunctionReturn(0);
35935bd34faSBarry Smith }
36035bd34faSBarry Smith 
3612f71e704SKris Buschelman /*MC
3622692d6eeSBarry Smith   MATSOLVERUMFPACK = "umfpack" - A matrix type providing direct solvers (LU) for sequential matrices
3632f71e704SKris Buschelman   via the external package UMFPACK.
3642f71e704SKris Buschelman 
365c2b89b5dSBarry Smith   Use ./configure --download-suitesparse to install PETSc to use UMFPACK
366c2b89b5dSBarry Smith 
3673ca39a21SBarry Smith   Use -pc_type lu -pc_factor_mat_solver_type umfpack to use this direct solver
3682f71e704SKris Buschelman 
3692f71e704SKris Buschelman   Consult UMFPACK documentation for more information about the Control parameters
3702f71e704SKris Buschelman   which correspond to the options database keys below.
3712f71e704SKris Buschelman 
3722f71e704SKris Buschelman   Options Database Keys:
373ba41fbb6SBarry Smith + -mat_umfpack_ordering                - CHOLMOD, AMD, GIVEN, METIS, BEST, NONE
374ba41fbb6SBarry Smith . -mat_umfpack_prl                     - UMFPACK print level: Control[UMFPACK_PRL]
375e08999f5SMatthew G Knepley . -mat_umfpack_strategy <AUTO>         - (choose one of) AUTO UNSYMMETRIC SYMMETRIC 2BY2
3762f71e704SKris Buschelman . -mat_umfpack_dense_col <alpha_c>     - UMFPACK dense column threshold: Control[UMFPACK_DENSE_COL]
377e08999f5SMatthew G Knepley . -mat_umfpack_dense_row <0.2>         - Control[UMFPACK_DENSE_ROW]
378e08999f5SMatthew G Knepley . -mat_umfpack_amd_dense <10>          - Control[UMFPACK_AMD_DENSE]
3792f71e704SKris Buschelman . -mat_umfpack_block_size <bs>         - UMFPACK block size for BLAS-Level 3 calls: Control[UMFPACK_BLOCK_SIZE]
380e08999f5SMatthew G Knepley . -mat_umfpack_2by2_tolerance <0.01>   - Control[UMFPACK_2BY2_TOLERANCE]
381e08999f5SMatthew G Knepley . -mat_umfpack_fixq <0>                - Control[UMFPACK_FIXQ]
382e08999f5SMatthew G Knepley . -mat_umfpack_aggressive <1>          - Control[UMFPACK_AGGRESSIVE]
3832f71e704SKris Buschelman . -mat_umfpack_pivot_tolerance <delta> - UMFPACK partial pivot tolerance: Control[UMFPACK_PIVOT_TOLERANCE]
384e08999f5SMatthew G Knepley . -mat_umfpack_sym_pivot_tolerance <0.001> - Control[UMFPACK_SYM_PIVOT_TOLERANCE]
385e08999f5SMatthew G Knepley . -mat_umfpack_scale <NONE>           - (choose one of) NONE SUM MAX
3862f71e704SKris Buschelman . -mat_umfpack_alloc_init <delta>      - UMFPACK factorized matrix allocation modifier: Control[UMFPACK_ALLOC_INIT]
387e08999f5SMatthew G Knepley . -mat_umfpack_droptol <0>            - Control[UMFPACK_DROPTOL]
3882f71e704SKris Buschelman - -mat_umfpack_irstep <maxit>          - UMFPACK maximum number of iterative refinement steps: Control[UMFPACK_IRSTEP]
3892f71e704SKris Buschelman 
3902f71e704SKris Buschelman    Level: beginner
391a364b7d2SBarry Smith 
392a364b7d2SBarry Smith    Note: UMFPACK is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html
3932f71e704SKris Buschelman 
394db781477SPatrick Sanan .seealso: `PCLU`, `MATSOLVERSUPERLU`, `MATSOLVERMUMPS`, `PCFactorSetMatSolverType()`, `MatSolverType`
3952f71e704SKris Buschelman M*/
396b2573a8aSBarry Smith 
3978cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_umfpack(Mat A,MatFactorType ftype,Mat *F)
3983519fcdcSHong Zhang {
3993519fcdcSHong Zhang   Mat            B;
4003519fcdcSHong Zhang   Mat_UMFPACK    *lu;
401*26cc229bSBarry Smith   PetscInt       m=A->rmap->n,n=A->cmap->n;
4022d4e2982SHong Zhang 
4033519fcdcSHong Zhang   PetscFunctionBegin;
4043519fcdcSHong Zhang   /* Create the factorization matrix F */
4059566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B));
4069566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n));
4079566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("umfpack",&((PetscObject)B)->type_name));
4089566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
409b9c12af5SBarry Smith 
4109566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(B,&lu));
4112205254eSKarl Rupp 
412b9c12af5SBarry Smith   B->data                   = lu;
413b9c12af5SBarry Smith   B->ops->getinfo          = MatGetInfo_External;
4143519fcdcSHong Zhang   B->ops->lufactorsymbolic = MatLUFactorSymbolic_UMFPACK;
4153519fcdcSHong Zhang   B->ops->destroy          = MatDestroy_UMFPACK;
4163519fcdcSHong Zhang   B->ops->view             = MatView_UMFPACK;
4171aef8b4cSStefano Zampini   B->ops->matsolve         = NULL;
4182205254eSKarl Rupp 
4199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_umfpack));
4202205254eSKarl Rupp 
421d5f3da31SBarry Smith   B->factortype   = MAT_FACTOR_LU;
4223519fcdcSHong Zhang   B->assembled    = PETSC_TRUE;           /* required by -ksp_view */
4233519fcdcSHong Zhang   B->preallocated = PETSC_TRUE;
4243519fcdcSHong Zhang 
4259566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
4269566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERUMFPACK,&B->solvertype));
427f73b0415SBarry Smith   B->canuseordering = PETSC_TRUE;
4289566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&B->preferredordering[MAT_FACTOR_LU]));
42900c67f3bSHong Zhang 
4303519fcdcSHong Zhang   /* initializations */
4313519fcdcSHong Zhang   /* ------------------------------------------------*/
4323519fcdcSHong Zhang   /* get the default control parameters */
4338592901bSBarry Smith   umfpack_UMF_defaults(lu->Control);
4340298fd71SBarry Smith   lu->perm_c                  = NULL; /* use defaul UMFPACK col permutation */
4353519fcdcSHong Zhang   lu->Control[UMFPACK_IRSTEP] = 0;          /* max num of iterative refinement steps to attempt */
4363519fcdcSHong Zhang 
4373519fcdcSHong Zhang   *F = B;
4383519fcdcSHong Zhang   PetscFunctionReturn(0);
4393519fcdcSHong Zhang }
440b2573a8aSBarry Smith 
441db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat,MatFactorType,Mat*);
442db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat,MatFactorType,Mat*);
443db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat,MatFactorType,Mat*);
444a2fc1e05SToby Isaac PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_spqr(Mat,MatFactorType,Mat*);
4452d4e2982SHong Zhang 
4463ca39a21SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_SuiteSparse(void)
44742c9c57cSBarry Smith {
44842c9c57cSBarry Smith   PetscFunctionBegin;
4499566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERUMFPACK,MATSEQAIJ,  MAT_FACTOR_LU,MatGetFactor_seqaij_umfpack));
4509566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERCHOLMOD,MATSEQAIJ,  MAT_FACTOR_CHOLESKY,MatGetFactor_seqaij_cholmod));
4519566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERCHOLMOD,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_seqsbaij_cholmod));
4529566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERKLU,MATSEQAIJ,      MAT_FACTOR_LU,MatGetFactor_seqaij_klu));
4539566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSPQR,MATSEQAIJ,     MAT_FACTOR_QR,MatGetFactor_seqaij_spqr));
45465f45395SPierre Jolivet   if (!PetscDefined(USE_COMPLEX)) {
4559566063dSJacob Faibussowitsch     PetscCall(MatSolverTypeRegister(MATSOLVERSPQR,MATNORMAL,   MAT_FACTOR_QR,MatGetFactor_seqaij_spqr));
45665f45395SPierre Jolivet   }
4579566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERSPQR,MATNORMALHERMITIAN, MAT_FACTOR_QR,MatGetFactor_seqaij_spqr));
45842c9c57cSBarry Smith   PetscFunctionReturn(0);
45942c9c57cSBarry Smith }
460