169e15a41SShri Abhyankar 269e15a41SShri Abhyankar /* 369e15a41SShri Abhyankar Provides an interface to the KLUv1.2 sparse solver 469e15a41SShri Abhyankar 569e15a41SShri Abhyankar When build with PETSC_USE_64BIT_INDICES this will use SuiteSparse_long as the 669e15a41SShri Abhyankar integer type in KLU, otherwise it will use int. This means 769e15a41SShri Abhyankar all integers in this file are simply declared as PetscInt. Also it means 87de69702SBarry Smith that KLU SuiteSparse_long version MUST be built with 64-bit integers when used. 969e15a41SShri Abhyankar 1069e15a41SShri Abhyankar */ 1169e15a41SShri Abhyankar #include <../src/mat/impls/aij/seq/aij.h> 1269e15a41SShri Abhyankar 1369e15a41SShri Abhyankar #if defined(PETSC_USE_64BIT_INDICES) 1469e15a41SShri Abhyankar #define klu_K_defaults klu_l_defaults 1530704e1fSBarry Smith #define klu_K_analyze(a, b, c, d) klu_l_analyze((SuiteSparse_long)a, (SuiteSparse_long *)b, (SuiteSparse_long *)c, d) 1630704e1fSBarry Smith #define klu_K_analyze_given(a, b, c, d, e, f) klu_l_analyze_given((SuiteSparse_long)a, (SuiteSparse_long *)b, (SuiteSparse_long *)c, (SuiteSparse_long *)d, (SuiteSparse_long *)e, f) 1769e15a41SShri Abhyankar #define klu_K_free_symbolic klu_l_free_symbolic 1869e15a41SShri Abhyankar #define klu_K_free_numeric klu_l_free_numeric 1969e15a41SShri Abhyankar #define klu_K_common klu_l_common 2069e15a41SShri Abhyankar #define klu_K_symbolic klu_l_symbolic 2169e15a41SShri Abhyankar #define klu_K_numeric klu_l_numeric 2269e15a41SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 2330704e1fSBarry Smith #define klu_K_factor(a, b, c, d, e) klu_zl_factor((SuiteSparse_long *)a, (SuiteSparse_long *)b, c, d, e); 2469e15a41SShri Abhyankar #define klu_K_solve klu_zl_solve 2569e15a41SShri Abhyankar #define klu_K_tsolve klu_zl_tsolve 2669e15a41SShri Abhyankar #define klu_K_refactor klu_zl_refactor 2769e15a41SShri Abhyankar #define klu_K_sort klu_zl_sort 2869e15a41SShri Abhyankar #define klu_K_flops klu_zl_flops 2969e15a41SShri Abhyankar #define klu_K_rgrowth klu_zl_rgrowth 3069e15a41SShri Abhyankar #define klu_K_condest klu_zl_condest 3169e15a41SShri Abhyankar #define klu_K_rcond klu_zl_rcond 3269e15a41SShri Abhyankar #define klu_K_scale klu_zl_scale 3369e15a41SShri Abhyankar #else 3430704e1fSBarry Smith #define klu_K_factor(a, b, c, d, e) klu_l_factor((SuiteSparse_long *)a, (SuiteSparse_long *)b, c, d, e); 3569e15a41SShri Abhyankar #define klu_K_solve klu_l_solve 3669e15a41SShri Abhyankar #define klu_K_tsolve klu_l_tsolve 3769e15a41SShri Abhyankar #define klu_K_refactor klu_l_refactor 3869e15a41SShri Abhyankar #define klu_K_sort klu_l_sort 3969e15a41SShri Abhyankar #define klu_K_flops klu_l_flops 4069e15a41SShri Abhyankar #define klu_K_rgrowth klu_l_rgrowth 4169e15a41SShri Abhyankar #define klu_K_condest klu_l_condest 4269e15a41SShri Abhyankar #define klu_K_rcond klu_l_rcond 4369e15a41SShri Abhyankar #define klu_K_scale klu_l_scale 4469e15a41SShri Abhyankar #endif 4569e15a41SShri Abhyankar #else 4669e15a41SShri Abhyankar #define klu_K_defaults klu_defaults 4769e15a41SShri Abhyankar #define klu_K_analyze klu_analyze 4869e15a41SShri Abhyankar #define klu_K_analyze_given klu_analyze_given 4969e15a41SShri Abhyankar #define klu_K_free_symbolic klu_free_symbolic 5069e15a41SShri Abhyankar #define klu_K_free_numeric klu_free_numeric 5169e15a41SShri Abhyankar #define klu_K_common klu_common 5269e15a41SShri Abhyankar #define klu_K_symbolic klu_symbolic 5369e15a41SShri Abhyankar #define klu_K_numeric klu_numeric 5469e15a41SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 5569e15a41SShri Abhyankar #define klu_K_factor klu_z_factor 5669e15a41SShri Abhyankar #define klu_K_solve klu_z_solve 5769e15a41SShri Abhyankar #define klu_K_tsolve klu_z_tsolve 5869e15a41SShri Abhyankar #define klu_K_refactor klu_z_refactor 5969e15a41SShri Abhyankar #define klu_K_sort klu_z_sort 6069e15a41SShri Abhyankar #define klu_K_flops klu_z_flops 6169e15a41SShri Abhyankar #define klu_K_rgrowth klu_z_rgrowth 6269e15a41SShri Abhyankar #define klu_K_condest klu_z_condest 6369e15a41SShri Abhyankar #define klu_K_rcond klu_z_rcond 6469e15a41SShri Abhyankar #define klu_K_scale klu_z_scale 6569e15a41SShri Abhyankar #else 6669e15a41SShri Abhyankar #define klu_K_factor klu_factor 6769e15a41SShri Abhyankar #define klu_K_solve klu_solve 6869e15a41SShri Abhyankar #define klu_K_tsolve klu_tsolve 6969e15a41SShri Abhyankar #define klu_K_refactor klu_refactor 7069e15a41SShri Abhyankar #define klu_K_sort klu_sort 7169e15a41SShri Abhyankar #define klu_K_flops klu_flops 7269e15a41SShri Abhyankar #define klu_K_rgrowth klu_rgrowth 7369e15a41SShri Abhyankar #define klu_K_condest klu_condest 7469e15a41SShri Abhyankar #define klu_K_rcond klu_rcond 7569e15a41SShri Abhyankar #define klu_K_scale klu_scale 7669e15a41SShri Abhyankar #endif 7769e15a41SShri Abhyankar #endif 7869e15a41SShri Abhyankar 7969e15a41SShri Abhyankar EXTERN_C_BEGIN 8069e15a41SShri Abhyankar #include <klu.h> 8169e15a41SShri Abhyankar EXTERN_C_END 8269e15a41SShri Abhyankar 834ac6704cSBarry Smith static const char *KluOrderingTypes[] = {"AMD", "COLAMD"}; 8469e15a41SShri Abhyankar static const char *scale[] = {"NONE", "SUM", "MAX"}; 8569e15a41SShri Abhyankar 8669e15a41SShri Abhyankar typedef struct { 8769e15a41SShri Abhyankar klu_K_common Common; 8869e15a41SShri Abhyankar klu_K_symbolic *Symbolic; 8969e15a41SShri Abhyankar klu_K_numeric *Numeric; 9069e15a41SShri Abhyankar PetscInt *perm_c, *perm_r; 9169e15a41SShri Abhyankar MatStructure flg; 9269e15a41SShri Abhyankar PetscBool PetscMatOrdering; 9369e15a41SShri Abhyankar PetscBool CleanUpKLU; 9469e15a41SShri Abhyankar } Mat_KLU; 9569e15a41SShri Abhyankar 96d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_KLU(Mat A) 97d71ae5a4SJacob Faibussowitsch { 98569c4379SBarry Smith Mat_KLU *lu = (Mat_KLU *)A->data; 9969e15a41SShri Abhyankar 10069e15a41SShri Abhyankar PetscFunctionBegin; 101569c4379SBarry Smith if (lu->CleanUpKLU) { 10269e15a41SShri Abhyankar klu_K_free_symbolic(&lu->Symbolic, &lu->Common); 10369e15a41SShri Abhyankar klu_K_free_numeric(&lu->Numeric, &lu->Common); 1049566063dSJacob Faibussowitsch PetscCall(PetscFree2(lu->perm_r, lu->perm_c)); 10569e15a41SShri Abhyankar } 1062e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL)); 1079566063dSJacob Faibussowitsch PetscCall(PetscFree(A->data)); 1083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10969e15a41SShri Abhyankar } 11069e15a41SShri Abhyankar 111d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolveTranspose_KLU(Mat A, Vec b, Vec x) 112d71ae5a4SJacob Faibussowitsch { 113569c4379SBarry Smith Mat_KLU *lu = (Mat_KLU *)A->data; 11469e15a41SShri Abhyankar PetscScalar *xa; 11569e15a41SShri Abhyankar PetscInt status; 11669e15a41SShri Abhyankar 11769e15a41SShri Abhyankar PetscFunctionBegin; 11869e15a41SShri Abhyankar /* KLU uses a column major format, solve Ax = b by klu_*_solve */ 1199566063dSJacob Faibussowitsch PetscCall(VecCopy(b, x)); /* klu_solve stores the solution in rhs */ 1209566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &xa)); 12169e15a41SShri Abhyankar status = klu_K_solve(lu->Symbolic, lu->Numeric, A->rmap->n, 1, (PetscReal *)xa, &lu->Common); 1225f80ce2aSJacob Faibussowitsch PetscCheck(status == 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "KLU Solve failed"); 1239566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &xa)); 1243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12569e15a41SShri Abhyankar } 12669e15a41SShri Abhyankar 127d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_KLU(Mat A, Vec b, Vec x) 128d71ae5a4SJacob Faibussowitsch { 129569c4379SBarry Smith Mat_KLU *lu = (Mat_KLU *)A->data; 13069e15a41SShri Abhyankar PetscScalar *xa; 13169e15a41SShri Abhyankar PetscInt status; 13269e15a41SShri Abhyankar 13369e15a41SShri Abhyankar PetscFunctionBegin; 13469e15a41SShri Abhyankar /* KLU uses a column major format, solve A^Tx = b by klu_*_tsolve */ 1359566063dSJacob Faibussowitsch PetscCall(VecCopy(b, x)); /* klu_solve stores the solution in rhs */ 1369566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &xa)); 13769e15a41SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 13869e15a41SShri Abhyankar PetscInt conj_solve = 1; 13969e15a41SShri Abhyankar status = klu_K_tsolve(lu->Symbolic, lu->Numeric, A->rmap->n, 1, (PetscReal *)xa, conj_solve, &lu->Common); /* conjugate solve */ 14069e15a41SShri Abhyankar #else 14169e15a41SShri Abhyankar status = klu_K_tsolve(lu->Symbolic, lu->Numeric, A->rmap->n, 1, xa, &lu->Common); 14269e15a41SShri Abhyankar #endif 1435f80ce2aSJacob Faibussowitsch PetscCheck(status == 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "KLU Solve failed"); 1449566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &xa)); 1453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14669e15a41SShri Abhyankar } 14769e15a41SShri Abhyankar 148d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorNumeric_KLU(Mat F, Mat A, const MatFactorInfo *info) 149d71ae5a4SJacob Faibussowitsch { 150569c4379SBarry Smith Mat_KLU *lu = (Mat_KLU *)(F)->data; 15169e15a41SShri Abhyankar Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 15269e15a41SShri Abhyankar PetscInt *ai = a->i, *aj = a->j; 15369e15a41SShri Abhyankar PetscScalar *av = a->a; 15469e15a41SShri Abhyankar 15569e15a41SShri Abhyankar PetscFunctionBegin; 15669e15a41SShri Abhyankar /* numeric factorization of A' */ 15769e15a41SShri Abhyankar 158ad540459SPierre Jolivet if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) klu_K_free_numeric(&lu->Numeric, &lu->Common); 15969e15a41SShri Abhyankar lu->Numeric = klu_K_factor(ai, aj, (PetscReal *)av, lu->Symbolic, &lu->Common); 1605f80ce2aSJacob Faibussowitsch PetscCheck(lu->Numeric, PETSC_COMM_SELF, PETSC_ERR_LIB, "KLU Numeric factorization failed"); 16169e15a41SShri Abhyankar 16269e15a41SShri Abhyankar lu->flg = SAME_NONZERO_PATTERN; 16369e15a41SShri Abhyankar lu->CleanUpKLU = PETSC_TRUE; 16469e15a41SShri Abhyankar F->ops->solve = MatSolve_KLU; 16569e15a41SShri Abhyankar F->ops->solvetranspose = MatSolveTranspose_KLU; 1663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16769e15a41SShri Abhyankar } 16869e15a41SShri Abhyankar 169d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorSymbolic_KLU(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info) 170d71ae5a4SJacob Faibussowitsch { 17169e15a41SShri Abhyankar Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 172569c4379SBarry Smith Mat_KLU *lu = (Mat_KLU *)(F->data); 17369e15a41SShri Abhyankar PetscInt i, *ai = a->i, *aj = a->j, m = A->rmap->n, n = A->cmap->n; 17469e15a41SShri Abhyankar const PetscInt *ra, *ca; 17569e15a41SShri Abhyankar 17669e15a41SShri Abhyankar PetscFunctionBegin; 17769e15a41SShri Abhyankar if (lu->PetscMatOrdering) { 1789566063dSJacob Faibussowitsch PetscCall(ISGetIndices(r, &ra)); 1799566063dSJacob Faibussowitsch PetscCall(ISGetIndices(c, &ca)); 1809566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(m, &lu->perm_r, n, &lu->perm_c)); 1817de69702SBarry Smith /* we cannot simply memcpy on 64-bit archs */ 18269e15a41SShri Abhyankar for (i = 0; i < m; i++) lu->perm_r[i] = ra[i]; 18369e15a41SShri Abhyankar for (i = 0; i < n; i++) lu->perm_c[i] = ca[i]; 1849566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(r, &ra)); 1859566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(c, &ca)); 18669e15a41SShri Abhyankar } 18769e15a41SShri Abhyankar 18869e15a41SShri Abhyankar /* symbolic factorization of A' */ 1894ac6704cSBarry Smith if (r) { 1904ac6704cSBarry Smith lu->PetscMatOrdering = PETSC_TRUE; 19169e15a41SShri Abhyankar lu->Symbolic = klu_K_analyze_given(n, ai, aj, lu->perm_c, lu->perm_r, &lu->Common); 19269e15a41SShri Abhyankar } else { /* use klu internal ordering */ 19369e15a41SShri Abhyankar lu->Symbolic = klu_K_analyze(n, ai, aj, &lu->Common); 19469e15a41SShri Abhyankar } 1955f80ce2aSJacob Faibussowitsch PetscCheck(lu->Symbolic, PETSC_COMM_SELF, PETSC_ERR_LIB, "KLU Symbolic Factorization failed"); 19669e15a41SShri Abhyankar 19769e15a41SShri Abhyankar lu->flg = DIFFERENT_NONZERO_PATTERN; 19869e15a41SShri Abhyankar lu->CleanUpKLU = PETSC_TRUE; 19969e15a41SShri Abhyankar (F)->ops->lufactornumeric = MatLUFactorNumeric_KLU; 2003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20169e15a41SShri Abhyankar } 20269e15a41SShri Abhyankar 203d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Info_KLU(Mat A, PetscViewer viewer) 204d71ae5a4SJacob Faibussowitsch { 205569c4379SBarry Smith Mat_KLU *lu = (Mat_KLU *)A->data; 20669e15a41SShri Abhyankar klu_K_numeric *Numeric = (klu_K_numeric *)lu->Numeric; 20769e15a41SShri Abhyankar 20869e15a41SShri Abhyankar PetscFunctionBegin; 2099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "KLU stats:\n")); 2109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Number of diagonal blocks: %" PetscInt_FMT "\n", (PetscInt)(Numeric->nblocks))); 2119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Total nonzeros=%" PetscInt_FMT "\n", (PetscInt)(Numeric->lnz + Numeric->unz))); 2129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "KLU runtime parameters:\n")); 21369e15a41SShri Abhyankar /* Control parameters used by numeric factorization */ 2149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Partial pivoting tolerance: %g\n", lu->Common.tol)); 21569e15a41SShri Abhyankar /* BTF preordering */ 2169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " BTF preordering enabled: %" PetscInt_FMT "\n", (PetscInt)(lu->Common.btf))); 21769e15a41SShri Abhyankar /* mat ordering */ 21848a46eb9SPierre Jolivet if (!lu->PetscMatOrdering) PetscCall(PetscViewerASCIIPrintf(viewer, " Ordering: %s (not using the PETSc ordering)\n", KluOrderingTypes[(int)lu->Common.ordering])); 21969e15a41SShri Abhyankar /* matrix row scaling */ 2209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Matrix row scaling: %s\n", scale[(int)lu->Common.scale])); 2213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22269e15a41SShri Abhyankar } 22369e15a41SShri Abhyankar 224d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_KLU(Mat A, PetscViewer viewer) 225d71ae5a4SJacob Faibussowitsch { 22669e15a41SShri Abhyankar PetscBool iascii; 22769e15a41SShri Abhyankar PetscViewerFormat format; 22869e15a41SShri Abhyankar 22969e15a41SShri Abhyankar PetscFunctionBegin; 2309566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 23169e15a41SShri Abhyankar if (iascii) { 2329566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 23348a46eb9SPierre Jolivet if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_KLU(A, viewer)); 23469e15a41SShri Abhyankar } 2353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 23669e15a41SShri Abhyankar } 23769e15a41SShri Abhyankar 238d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFactorGetSolverType_seqaij_klu(Mat A, MatSolverType *type) 239d71ae5a4SJacob Faibussowitsch { 24069e15a41SShri Abhyankar PetscFunctionBegin; 24169e15a41SShri Abhyankar *type = MATSOLVERKLU; 2423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 24369e15a41SShri Abhyankar } 24469e15a41SShri Abhyankar 24569e15a41SShri Abhyankar /*MC 24611a5261eSBarry Smith MATSOLVERKLU = "klu" - A matrix type providing direct solvers, LU, for sequential matrices 24769e15a41SShri Abhyankar via the external package KLU. 24869e15a41SShri Abhyankar 2492ef1f0ffSBarry Smith `./configure --download-suitesparse` to install PETSc to use KLU 25069e15a41SShri Abhyankar 2512ef1f0ffSBarry Smith Use `-pc_type lu` `-pc_factor_mat_solver_type klu` to use this direct solver 252c2b89b5dSBarry Smith 25369e15a41SShri Abhyankar Consult KLU documentation for more information on the options database keys below. 25469e15a41SShri Abhyankar 25569e15a41SShri Abhyankar Options Database Keys: 25669e15a41SShri Abhyankar + -mat_klu_pivot_tol <0.001> - Partial pivoting tolerance 25769e15a41SShri Abhyankar . -mat_klu_use_btf <1> - Use BTF preordering 2582ef1f0ffSBarry Smith . -mat_klu_ordering <AMD> - KLU reordering scheme to reduce fill-in (choose one of) `AMD`, `COLAMD`, `PETSC` 2592ef1f0ffSBarry Smith - -mat_klu_row_scale <NONE> - Matrix row scaling (choose one of) `NONE`, `SUM`, `MAX` 260a364b7d2SBarry Smith 26169e15a41SShri Abhyankar Level: beginner 26269e15a41SShri Abhyankar 2632ef1f0ffSBarry Smith Note: 2642ef1f0ffSBarry Smith KLU is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html 2652ef1f0ffSBarry Smith 266*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCLU`, `MATSOLVERUMFPACK`, `MATSOLVERCHOLMOD`, `PCFactorSetMatSolverType()`, `MatSolverType` 26769e15a41SShri Abhyankar M*/ 26869e15a41SShri Abhyankar 269d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat A, MatFactorType ftype, Mat *F) 270d71ae5a4SJacob Faibussowitsch { 27169e15a41SShri Abhyankar Mat B; 27269e15a41SShri Abhyankar Mat_KLU *lu; 2734ac6704cSBarry Smith PetscInt m = A->rmap->n, n = A->cmap->n, idx = 0, status; 27469e15a41SShri Abhyankar PetscBool flg; 27569e15a41SShri Abhyankar 27669e15a41SShri Abhyankar PetscFunctionBegin; 27769e15a41SShri Abhyankar /* Create the factorization matrix F */ 2789566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 2799566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n)); 2809566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("klu", &((PetscObject)B)->type_name)); 2819566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 282569c4379SBarry Smith 2834dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&lu)); 28469e15a41SShri Abhyankar 285569c4379SBarry Smith B->data = lu; 286569c4379SBarry Smith B->ops->getinfo = MatGetInfo_External; 28769e15a41SShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_KLU; 28869e15a41SShri Abhyankar B->ops->destroy = MatDestroy_KLU; 28969e15a41SShri Abhyankar B->ops->view = MatView_KLU; 29069e15a41SShri Abhyankar 2919566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_klu)); 29269e15a41SShri Abhyankar 29369e15a41SShri Abhyankar B->factortype = MAT_FACTOR_LU; 29469e15a41SShri Abhyankar B->assembled = PETSC_TRUE; /* required by -ksp_view */ 29569e15a41SShri Abhyankar B->preallocated = PETSC_TRUE; 29669e15a41SShri Abhyankar 2979566063dSJacob Faibussowitsch PetscCall(PetscFree(B->solvertype)); 2989566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERKLU, &B->solvertype)); 299f73b0415SBarry Smith B->canuseordering = PETSC_TRUE; 3009566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_LU])); 30100c67f3bSHong Zhang 30269e15a41SShri Abhyankar /* initializations */ 30369e15a41SShri Abhyankar /* get the default control parameters */ 30469e15a41SShri Abhyankar status = klu_K_defaults(&lu->Common); 3055f80ce2aSJacob Faibussowitsch PetscCheck(status > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "KLU Initialization failed"); 3066c4ed002SBarry Smith 30769e15a41SShri Abhyankar lu->Common.scale = 0; /* No row scaling */ 30869e15a41SShri Abhyankar 30926cc229bSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)B), ((PetscObject)B)->prefix, "KLU Options", "Mat"); 31069e15a41SShri Abhyankar /* Partial pivoting tolerance */ 3119566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_klu_pivot_tol", "Partial pivoting tolerance", "None", lu->Common.tol, &lu->Common.tol, NULL)); 31269e15a41SShri Abhyankar /* BTF pre-ordering */ 3139566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_klu_use_btf", "Enable BTF preordering", "None", (PetscInt)lu->Common.btf, (PetscInt *)&lu->Common.btf, NULL)); 31469e15a41SShri Abhyankar /* Matrix reordering */ 315dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-mat_klu_ordering", "Internal ordering method", "None", KluOrderingTypes, PETSC_STATIC_ARRAY_LENGTH(KluOrderingTypes), KluOrderingTypes[0], &idx, &flg)); 3164ac6704cSBarry Smith lu->Common.ordering = (int)idx; 31769e15a41SShri Abhyankar /* Matrix row scaling */ 3189566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-mat_klu_row_scale", "Matrix row scaling", "None", scale, 3, scale[0], &idx, &flg)); 319d0609cedSBarry Smith PetscOptionsEnd(); 32069e15a41SShri Abhyankar *F = B; 3213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 32269e15a41SShri Abhyankar } 323