xref: /petsc/src/mat/impls/aij/seq/klu/klu.c (revision 941063a0c6742edd6ba3de0434f3e015738e856a)
169e15a41SShri Abhyankar /*
269e15a41SShri Abhyankar    Provides an interface to the KLUv1.2 sparse solver
369e15a41SShri Abhyankar 
469e15a41SShri Abhyankar    When build with PETSC_USE_64BIT_INDICES this will use SuiteSparse_long as the
569e15a41SShri Abhyankar    integer type in KLU, otherwise it will use int. This means
669e15a41SShri Abhyankar    all integers in this file are simply declared as PetscInt. Also it means
77de69702SBarry Smith    that KLU SuiteSparse_long version MUST be built with 64-bit integers when used.
869e15a41SShri Abhyankar 
969e15a41SShri Abhyankar */
1069e15a41SShri Abhyankar #include <../src/mat/impls/aij/seq/aij.h>
1169e15a41SShri Abhyankar 
1269e15a41SShri Abhyankar #if defined(PETSC_USE_64BIT_INDICES)
1369e15a41SShri Abhyankar   #define klu_K_defaults                        klu_l_defaults
1430704e1fSBarry Smith   #define klu_K_analyze(a, b, c, d)             klu_l_analyze((SuiteSparse_long)a, (SuiteSparse_long *)b, (SuiteSparse_long *)c, d)
1530704e1fSBarry 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)
1669e15a41SShri Abhyankar   #define klu_K_free_symbolic                   klu_l_free_symbolic
1769e15a41SShri Abhyankar   #define klu_K_free_numeric                    klu_l_free_numeric
1869e15a41SShri Abhyankar   #define klu_K_common                          klu_l_common
1969e15a41SShri Abhyankar   #define klu_K_symbolic                        klu_l_symbolic
2069e15a41SShri Abhyankar   #define klu_K_numeric                         klu_l_numeric
2169e15a41SShri Abhyankar   #if defined(PETSC_USE_COMPLEX)
2230704e1fSBarry Smith     #define klu_K_factor(a, b, c, d, e)       klu_zl_factor((SuiteSparse_long *)a, (SuiteSparse_long *)b, c, d, e);
2369e15a41SShri Abhyankar     #define klu_K_solve                       klu_zl_solve
24*407b358cSPierre Jolivet     #define klu_K_tsolve(a, b, c, d, e, f, g) klu_zl_tsolve((a), (b), (c), (d), (PetscReal *)(e), (f), (g))
2569e15a41SShri Abhyankar     #define klu_K_refactor                    klu_zl_refactor
2669e15a41SShri Abhyankar     #define klu_K_sort                        klu_zl_sort
2769e15a41SShri Abhyankar     #define klu_K_flops                       klu_zl_flops
2869e15a41SShri Abhyankar     #define klu_K_rgrowth                     klu_zl_rgrowth
2969e15a41SShri Abhyankar     #define klu_K_condest                     klu_zl_condest
3069e15a41SShri Abhyankar     #define klu_K_rcond                       klu_zl_rcond
3169e15a41SShri Abhyankar     #define klu_K_scale                       klu_zl_scale
3269e15a41SShri Abhyankar   #else
3330704e1fSBarry Smith     #define klu_K_factor(a, b, c, d, e)       klu_l_factor((SuiteSparse_long *)a, (SuiteSparse_long *)b, c, d, e);
3469e15a41SShri Abhyankar     #define klu_K_solve                       klu_l_solve
35*407b358cSPierre Jolivet     #define klu_K_tsolve(a, b, c, d, e, f, g) klu_l_tsolve((a), (b), (c), (d), (e), (g))
3669e15a41SShri Abhyankar     #define klu_K_refactor                    klu_l_refactor
3769e15a41SShri Abhyankar     #define klu_K_sort                        klu_l_sort
3869e15a41SShri Abhyankar     #define klu_K_flops                       klu_l_flops
3969e15a41SShri Abhyankar     #define klu_K_rgrowth                     klu_l_rgrowth
4069e15a41SShri Abhyankar     #define klu_K_condest                     klu_l_condest
4169e15a41SShri Abhyankar     #define klu_K_rcond                       klu_l_rcond
4269e15a41SShri Abhyankar     #define klu_K_scale                       klu_l_scale
4369e15a41SShri Abhyankar   #endif
4469e15a41SShri Abhyankar #else
4569e15a41SShri Abhyankar   #define klu_K_defaults      klu_defaults
4669e15a41SShri Abhyankar   #define klu_K_analyze       klu_analyze
4769e15a41SShri Abhyankar   #define klu_K_analyze_given klu_analyze_given
4869e15a41SShri Abhyankar   #define klu_K_free_symbolic klu_free_symbolic
4969e15a41SShri Abhyankar   #define klu_K_free_numeric  klu_free_numeric
5069e15a41SShri Abhyankar   #define klu_K_common        klu_common
5169e15a41SShri Abhyankar   #define klu_K_symbolic      klu_symbolic
5269e15a41SShri Abhyankar   #define klu_K_numeric       klu_numeric
5369e15a41SShri Abhyankar   #if defined(PETSC_USE_COMPLEX)
5469e15a41SShri Abhyankar     #define klu_K_factor                      klu_z_factor
5569e15a41SShri Abhyankar     #define klu_K_solve                       klu_z_solve
56*407b358cSPierre Jolivet     #define klu_K_tsolve(a, b, c, d, e, f, g) klu_z_tsolve((a), (b), (c), (d), (PetscReal *)(e), (f), (g))
5769e15a41SShri Abhyankar     #define klu_K_refactor                    klu_z_refactor
5869e15a41SShri Abhyankar     #define klu_K_sort                        klu_z_sort
5969e15a41SShri Abhyankar     #define klu_K_flops                       klu_z_flops
6069e15a41SShri Abhyankar     #define klu_K_rgrowth                     klu_z_rgrowth
6169e15a41SShri Abhyankar     #define klu_K_condest                     klu_z_condest
6269e15a41SShri Abhyankar     #define klu_K_rcond                       klu_z_rcond
6369e15a41SShri Abhyankar     #define klu_K_scale                       klu_z_scale
6469e15a41SShri Abhyankar   #else
6569e15a41SShri Abhyankar     #define klu_K_factor                      klu_factor
6669e15a41SShri Abhyankar     #define klu_K_solve                       klu_solve
67*407b358cSPierre Jolivet     #define klu_K_tsolve(a, b, c, d, e, f, g) klu_tsolve((a), (b), (c), (d), (e), (g))
6869e15a41SShri Abhyankar     #define klu_K_refactor                    klu_refactor
6969e15a41SShri Abhyankar     #define klu_K_sort                        klu_sort
7069e15a41SShri Abhyankar     #define klu_K_flops                       klu_flops
7169e15a41SShri Abhyankar     #define klu_K_rgrowth                     klu_rgrowth
7269e15a41SShri Abhyankar     #define klu_K_condest                     klu_condest
7369e15a41SShri Abhyankar     #define klu_K_rcond                       klu_rcond
7469e15a41SShri Abhyankar     #define klu_K_scale                       klu_scale
7569e15a41SShri Abhyankar   #endif
7669e15a41SShri Abhyankar #endif
7769e15a41SShri Abhyankar 
7869e15a41SShri Abhyankar EXTERN_C_BEGIN
7969e15a41SShri Abhyankar #include <klu.h>
8069e15a41SShri Abhyankar EXTERN_C_END
8169e15a41SShri Abhyankar 
824ac6704cSBarry Smith static const char *KluOrderingTypes[] = {"AMD", "COLAMD"};
8369e15a41SShri Abhyankar static const char *scale[]            = {"NONE", "SUM", "MAX"};
8469e15a41SShri Abhyankar 
8569e15a41SShri Abhyankar typedef struct {
8669e15a41SShri Abhyankar   klu_K_common    Common;
8769e15a41SShri Abhyankar   klu_K_symbolic *Symbolic;
8869e15a41SShri Abhyankar   klu_K_numeric  *Numeric;
8969e15a41SShri Abhyankar   PetscInt       *perm_c, *perm_r;
9069e15a41SShri Abhyankar   MatStructure    flg;
9169e15a41SShri Abhyankar   PetscBool       PetscMatOrdering;
9269e15a41SShri Abhyankar   PetscBool       CleanUpKLU;
9369e15a41SShri Abhyankar } Mat_KLU;
9469e15a41SShri Abhyankar 
MatDestroy_KLU(Mat A)95d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_KLU(Mat A)
96d71ae5a4SJacob Faibussowitsch {
97569c4379SBarry Smith   Mat_KLU *lu = (Mat_KLU *)A->data;
9869e15a41SShri Abhyankar 
9969e15a41SShri Abhyankar   PetscFunctionBegin;
100569c4379SBarry Smith   if (lu->CleanUpKLU) {
10169e15a41SShri Abhyankar     klu_K_free_symbolic(&lu->Symbolic, &lu->Common);
10269e15a41SShri Abhyankar     klu_K_free_numeric(&lu->Numeric, &lu->Common);
1039566063dSJacob Faibussowitsch     PetscCall(PetscFree2(lu->perm_r, lu->perm_c));
10469e15a41SShri Abhyankar   }
1052e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL));
1069566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
1073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10869e15a41SShri Abhyankar }
10969e15a41SShri Abhyankar 
MatSolveTranspose_KLU(Mat A,Vec b,Vec x)110d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolveTranspose_KLU(Mat A, Vec b, Vec x)
111d71ae5a4SJacob Faibussowitsch {
112569c4379SBarry Smith   Mat_KLU     *lu = (Mat_KLU *)A->data;
11369e15a41SShri Abhyankar   PetscScalar *xa;
11469e15a41SShri Abhyankar   PetscInt     status;
11569e15a41SShri Abhyankar 
11669e15a41SShri Abhyankar   PetscFunctionBegin;
11769e15a41SShri Abhyankar   /* KLU uses a column major format, solve Ax = b by klu_*_solve */
1189566063dSJacob Faibussowitsch   PetscCall(VecCopy(b, x)); /* klu_solve stores the solution in rhs */
1199566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &xa));
12069e15a41SShri Abhyankar   status = klu_K_solve(lu->Symbolic, lu->Numeric, A->rmap->n, 1, (PetscReal *)xa, &lu->Common);
1215f80ce2aSJacob Faibussowitsch   PetscCheck(status == 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "KLU Solve failed");
1229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &xa));
1233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12469e15a41SShri Abhyankar }
12569e15a41SShri Abhyankar 
MatSolve_KLU(Mat A,Vec b,Vec x)126d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_KLU(Mat A, Vec b, Vec x)
127d71ae5a4SJacob Faibussowitsch {
128569c4379SBarry Smith   Mat_KLU     *lu = (Mat_KLU *)A->data;
12969e15a41SShri Abhyankar   PetscScalar *xa;
13069e15a41SShri Abhyankar 
13169e15a41SShri Abhyankar   PetscFunctionBegin;
13269e15a41SShri Abhyankar   /* KLU uses a column major format, solve A^Tx = b by klu_*_tsolve */
1339566063dSJacob Faibussowitsch   PetscCall(VecCopy(b, x)); /* klu_solve stores the solution in rhs */
1349566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &xa));
135*407b358cSPierre Jolivet   PetscCheck(klu_K_tsolve(lu->Symbolic, lu->Numeric, A->rmap->n, 1, xa, 1, &lu->Common), PETSC_COMM_SELF, PETSC_ERR_LIB, "KLU Solve failed"); /* conjugate solve */
1369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &xa));
1373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13869e15a41SShri Abhyankar }
13969e15a41SShri Abhyankar 
MatLUFactorNumeric_KLU(Mat F,Mat A,const MatFactorInfo * info)140d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorNumeric_KLU(Mat F, Mat A, const MatFactorInfo *info)
141d71ae5a4SJacob Faibussowitsch {
14257508eceSPierre Jolivet   Mat_KLU     *lu = (Mat_KLU *)F->data;
14369e15a41SShri Abhyankar   Mat_SeqAIJ  *a  = (Mat_SeqAIJ *)A->data;
14469e15a41SShri Abhyankar   PetscInt    *ai = a->i, *aj = a->j;
14569e15a41SShri Abhyankar   PetscScalar *av = a->a;
14669e15a41SShri Abhyankar 
14769e15a41SShri Abhyankar   PetscFunctionBegin;
14869e15a41SShri Abhyankar   /* numeric factorization of A' */
14969e15a41SShri Abhyankar 
150ad540459SPierre Jolivet   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) klu_K_free_numeric(&lu->Numeric, &lu->Common);
15169e15a41SShri Abhyankar   lu->Numeric = klu_K_factor(ai, aj, (PetscReal *)av, lu->Symbolic, &lu->Common);
1525f80ce2aSJacob Faibussowitsch   PetscCheck(lu->Numeric, PETSC_COMM_SELF, PETSC_ERR_LIB, "KLU Numeric factorization failed");
15369e15a41SShri Abhyankar 
15469e15a41SShri Abhyankar   lu->flg                = SAME_NONZERO_PATTERN;
15569e15a41SShri Abhyankar   lu->CleanUpKLU         = PETSC_TRUE;
15669e15a41SShri Abhyankar   F->ops->solve          = MatSolve_KLU;
15769e15a41SShri Abhyankar   F->ops->solvetranspose = MatSolveTranspose_KLU;
1583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15969e15a41SShri Abhyankar }
16069e15a41SShri Abhyankar 
MatLUFactorSymbolic_KLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo * info)161d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorSymbolic_KLU(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info)
162d71ae5a4SJacob Faibussowitsch {
16369e15a41SShri Abhyankar   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *)A->data;
164f4f49eeaSPierre Jolivet   Mat_KLU        *lu = (Mat_KLU *)F->data;
16569e15a41SShri Abhyankar   PetscInt        i, *ai = a->i, *aj = a->j, m = A->rmap->n, n = A->cmap->n;
16669e15a41SShri Abhyankar   const PetscInt *ra, *ca;
16769e15a41SShri Abhyankar 
16869e15a41SShri Abhyankar   PetscFunctionBegin;
16969e15a41SShri Abhyankar   if (lu->PetscMatOrdering) {
1709566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(r, &ra));
1719566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(c, &ca));
1729566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(m, &lu->perm_r, n, &lu->perm_c));
1737de69702SBarry Smith     /* we cannot simply memcpy on 64-bit archs */
17469e15a41SShri Abhyankar     for (i = 0; i < m; i++) lu->perm_r[i] = ra[i];
17569e15a41SShri Abhyankar     for (i = 0; i < n; i++) lu->perm_c[i] = ca[i];
1769566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(r, &ra));
1779566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(c, &ca));
17869e15a41SShri Abhyankar   }
17969e15a41SShri Abhyankar 
18069e15a41SShri Abhyankar   /* symbolic factorization of A' */
1814ac6704cSBarry Smith   if (r) {
1824ac6704cSBarry Smith     lu->PetscMatOrdering = PETSC_TRUE;
18369e15a41SShri Abhyankar     lu->Symbolic         = klu_K_analyze_given(n, ai, aj, lu->perm_c, lu->perm_r, &lu->Common);
18469e15a41SShri Abhyankar   } else { /* use klu internal ordering */
18569e15a41SShri Abhyankar     lu->Symbolic = klu_K_analyze(n, ai, aj, &lu->Common);
18669e15a41SShri Abhyankar   }
1875f80ce2aSJacob Faibussowitsch   PetscCheck(lu->Symbolic, PETSC_COMM_SELF, PETSC_ERR_LIB, "KLU Symbolic Factorization failed");
18869e15a41SShri Abhyankar 
18969e15a41SShri Abhyankar   lu->flg                 = DIFFERENT_NONZERO_PATTERN;
19069e15a41SShri Abhyankar   lu->CleanUpKLU          = PETSC_TRUE;
19157508eceSPierre Jolivet   F->ops->lufactornumeric = MatLUFactorNumeric_KLU;
1923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19369e15a41SShri Abhyankar }
19469e15a41SShri Abhyankar 
MatView_Info_KLU(Mat A,PetscViewer viewer)195d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Info_KLU(Mat A, PetscViewer viewer)
196d71ae5a4SJacob Faibussowitsch {
197569c4379SBarry Smith   Mat_KLU       *lu      = (Mat_KLU *)A->data;
19869e15a41SShri Abhyankar   klu_K_numeric *Numeric = (klu_K_numeric *)lu->Numeric;
19969e15a41SShri Abhyankar 
20069e15a41SShri Abhyankar   PetscFunctionBegin;
2019566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "KLU stats:\n"));
202f4f49eeaSPierre Jolivet   PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of diagonal blocks: %" PetscInt_FMT "\n", (PetscInt)Numeric->nblocks));
2039566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Total nonzeros=%" PetscInt_FMT "\n", (PetscInt)(Numeric->lnz + Numeric->unz)));
2049566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "KLU runtime parameters:\n"));
20569e15a41SShri Abhyankar   /* Control parameters used by numeric factorization */
2069566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Partial pivoting tolerance: %g\n", lu->Common.tol));
20769e15a41SShri Abhyankar   /* BTF preordering */
208f4f49eeaSPierre Jolivet   PetscCall(PetscViewerASCIIPrintf(viewer, "  BTF preordering enabled: %" PetscInt_FMT "\n", (PetscInt)lu->Common.btf));
20969e15a41SShri Abhyankar   /* mat ordering */
21048a46eb9SPierre Jolivet   if (!lu->PetscMatOrdering) PetscCall(PetscViewerASCIIPrintf(viewer, "  Ordering: %s (not using the PETSc ordering)\n", KluOrderingTypes[(int)lu->Common.ordering]));
21169e15a41SShri Abhyankar   /* matrix row scaling */
2129566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Matrix row scaling: %s\n", scale[(int)lu->Common.scale]));
2133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21469e15a41SShri Abhyankar }
21569e15a41SShri Abhyankar 
MatView_KLU(Mat A,PetscViewer viewer)216d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_KLU(Mat A, PetscViewer viewer)
217d71ae5a4SJacob Faibussowitsch {
2189f196a02SMartin Diehl   PetscBool         isascii;
21969e15a41SShri Abhyankar   PetscViewerFormat format;
22069e15a41SShri Abhyankar 
22169e15a41SShri Abhyankar   PetscFunctionBegin;
2229f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2239f196a02SMartin Diehl   if (isascii) {
2249566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
22548a46eb9SPierre Jolivet     if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_KLU(A, viewer));
22669e15a41SShri Abhyankar   }
2273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22869e15a41SShri Abhyankar }
22969e15a41SShri Abhyankar 
MatFactorGetSolverType_seqaij_klu(Mat A,MatSolverType * type)23066976f2fSJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_seqaij_klu(Mat A, MatSolverType *type)
231d71ae5a4SJacob Faibussowitsch {
23269e15a41SShri Abhyankar   PetscFunctionBegin;
23369e15a41SShri Abhyankar   *type = MATSOLVERKLU;
2343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23569e15a41SShri Abhyankar }
23669e15a41SShri Abhyankar 
23769e15a41SShri Abhyankar /*MC
23811a5261eSBarry Smith   MATSOLVERKLU = "klu" - A matrix type providing direct solvers, LU, for sequential matrices
23969e15a41SShri Abhyankar   via the external package KLU.
24069e15a41SShri Abhyankar 
2412ef1f0ffSBarry Smith   `./configure --download-suitesparse` to install PETSc to use KLU
24269e15a41SShri Abhyankar 
2432ef1f0ffSBarry Smith   Use `-pc_type lu` `-pc_factor_mat_solver_type klu` to use this direct solver
244c2b89b5dSBarry Smith 
24569e15a41SShri Abhyankar   Consult KLU documentation for more information on the options database keys below.
24669e15a41SShri Abhyankar 
24769e15a41SShri Abhyankar   Options Database Keys:
24869e15a41SShri Abhyankar + -mat_klu_pivot_tol <0.001>                  - Partial pivoting tolerance
24969e15a41SShri Abhyankar . -mat_klu_use_btf <1>                        - Use BTF preordering
2502ef1f0ffSBarry Smith . -mat_klu_ordering <AMD>                     - KLU reordering scheme to reduce fill-in (choose one of) `AMD`, `COLAMD`, `PETSC`
2512ef1f0ffSBarry Smith - -mat_klu_row_scale <NONE>                   - Matrix row scaling (choose one of) `NONE`, `SUM`, `MAX`
252a364b7d2SBarry Smith 
25369e15a41SShri Abhyankar    Level: beginner
25469e15a41SShri Abhyankar 
2552ef1f0ffSBarry Smith    Note:
2561d27aa22SBarry Smith    KLU is part of SuiteSparse <http://faculty.cse.tamu.edu/davis/suitesparse.html>
2572ef1f0ffSBarry Smith 
2581cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCLU`, `MATSOLVERUMFPACK`, `MATSOLVERCHOLMOD`, `PCFactorSetMatSolverType()`, `MatSolverType`
25969e15a41SShri Abhyankar M*/
26069e15a41SShri Abhyankar 
MatGetFactor_seqaij_klu(Mat A,MatFactorType ftype,Mat * F)261d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat A, MatFactorType ftype, Mat *F)
262d71ae5a4SJacob Faibussowitsch {
26369e15a41SShri Abhyankar   Mat       B;
26469e15a41SShri Abhyankar   Mat_KLU  *lu;
2654ac6704cSBarry Smith   PetscInt  m = A->rmap->n, n = A->cmap->n, idx = 0, status;
26669e15a41SShri Abhyankar   PetscBool flg;
26769e15a41SShri Abhyankar 
26869e15a41SShri Abhyankar   PetscFunctionBegin;
26969e15a41SShri Abhyankar   /* Create the factorization matrix F */
2709566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
2719566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n));
2729566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("klu", &((PetscObject)B)->type_name));
2739566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
274569c4379SBarry Smith 
2754dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&lu));
27669e15a41SShri Abhyankar 
277569c4379SBarry Smith   B->data                  = lu;
278569c4379SBarry Smith   B->ops->getinfo          = MatGetInfo_External;
27969e15a41SShri Abhyankar   B->ops->lufactorsymbolic = MatLUFactorSymbolic_KLU;
28069e15a41SShri Abhyankar   B->ops->destroy          = MatDestroy_KLU;
28169e15a41SShri Abhyankar   B->ops->view             = MatView_KLU;
28269e15a41SShri Abhyankar 
2839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_klu));
28469e15a41SShri Abhyankar 
28569e15a41SShri Abhyankar   B->factortype   = MAT_FACTOR_LU;
28669e15a41SShri Abhyankar   B->assembled    = PETSC_TRUE; /* required by -ksp_view */
28769e15a41SShri Abhyankar   B->preallocated = PETSC_TRUE;
28869e15a41SShri Abhyankar 
2899566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
2909566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERKLU, &B->solvertype));
291f73b0415SBarry Smith   B->canuseordering = PETSC_TRUE;
2929566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_LU]));
29300c67f3bSHong Zhang 
29469e15a41SShri Abhyankar   /* initializations */
29569e15a41SShri Abhyankar   /* get the default control parameters */
29669e15a41SShri Abhyankar   status = klu_K_defaults(&lu->Common);
2975f80ce2aSJacob Faibussowitsch   PetscCheck(status > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "KLU Initialization failed");
2986c4ed002SBarry Smith 
29969e15a41SShri Abhyankar   lu->Common.scale = 0; /* No row scaling */
30069e15a41SShri Abhyankar 
30126cc229bSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)B), ((PetscObject)B)->prefix, "KLU Options", "Mat");
30269e15a41SShri Abhyankar   /* Partial pivoting tolerance */
3039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_klu_pivot_tol", "Partial pivoting tolerance", "None", lu->Common.tol, &lu->Common.tol, NULL));
30469e15a41SShri Abhyankar   /* BTF pre-ordering */
3059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_klu_use_btf", "Enable BTF preordering", "None", (PetscInt)lu->Common.btf, (PetscInt *)&lu->Common.btf, NULL));
30669e15a41SShri Abhyankar   /* Matrix reordering */
307dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-mat_klu_ordering", "Internal ordering method", "None", KluOrderingTypes, PETSC_STATIC_ARRAY_LENGTH(KluOrderingTypes), KluOrderingTypes[0], &idx, &flg));
3084ac6704cSBarry Smith   lu->Common.ordering = (int)idx;
30969e15a41SShri Abhyankar   /* Matrix row scaling */
3109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-mat_klu_row_scale", "Matrix row scaling", "None", scale, 3, scale[0], &idx, &flg));
311d0609cedSBarry Smith   PetscOptionsEnd();
31269e15a41SShri Abhyankar   *F = B;
3133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31469e15a41SShri Abhyankar }
315