xref: /petsc/src/mat/impls/aij/seq/klu/klu.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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
869e15a41SShri Abhyankar    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 
9669e15a41SShri Abhyankar static PetscErrorCode MatDestroy_KLU(Mat A)
9769e15a41SShri Abhyankar {
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);
104*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree2(lu->perm_r,lu->perm_c));
10569e15a41SShri Abhyankar   }
106*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(A->data));
10769e15a41SShri Abhyankar   PetscFunctionReturn(0);
10869e15a41SShri Abhyankar }
10969e15a41SShri Abhyankar 
11069e15a41SShri Abhyankar static PetscErrorCode MatSolveTranspose_KLU(Mat A,Vec b,Vec x)
11169e15a41SShri Abhyankar {
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 */
11869e15a41SShri Abhyankar   /* ----------------------------------*/
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(b,x)); /* klu_solve stores the solution in rhs */
120*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(x,&xa));
12169e15a41SShri Abhyankar   status = klu_K_solve(lu->Symbolic,lu->Numeric,A->rmap->n,1,(PetscReal*)xa,&lu->Common);
122*5f80ce2aSJacob Faibussowitsch   PetscCheck(status == 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Solve failed");
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(x,&xa));
12469e15a41SShri Abhyankar   PetscFunctionReturn(0);
12569e15a41SShri Abhyankar }
12669e15a41SShri Abhyankar 
12769e15a41SShri Abhyankar static PetscErrorCode MatSolve_KLU(Mat A,Vec b,Vec x)
12869e15a41SShri Abhyankar {
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 */
13569e15a41SShri Abhyankar   /* ----------------------------------*/
136*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(b,x)); /* klu_solve stores the solution in rhs */
137*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(x,&xa));
13869e15a41SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
13969e15a41SShri Abhyankar   PetscInt conj_solve=1;
14069e15a41SShri Abhyankar   status = klu_K_tsolve(lu->Symbolic,lu->Numeric,A->rmap->n,1,(PetscReal*)xa,conj_solve,&lu->Common); /* conjugate solve */
14169e15a41SShri Abhyankar #else
14269e15a41SShri Abhyankar   status = klu_K_tsolve(lu->Symbolic,lu->Numeric,A->rmap->n,1,xa,&lu->Common);
14369e15a41SShri Abhyankar #endif
144*5f80ce2aSJacob Faibussowitsch   PetscCheck(status == 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Solve failed");
145*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(x,&xa));
14669e15a41SShri Abhyankar   PetscFunctionReturn(0);
14769e15a41SShri Abhyankar }
14869e15a41SShri Abhyankar 
14969e15a41SShri Abhyankar static PetscErrorCode MatLUFactorNumeric_KLU(Mat F,Mat A,const MatFactorInfo *info)
15069e15a41SShri Abhyankar {
151569c4379SBarry Smith   Mat_KLU        *lu = (Mat_KLU*)(F)->data;
15269e15a41SShri Abhyankar   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
15369e15a41SShri Abhyankar   PetscInt       *ai = a->i,*aj=a->j;
15469e15a41SShri Abhyankar   PetscScalar    *av = a->a;
15569e15a41SShri Abhyankar 
15669e15a41SShri Abhyankar   PetscFunctionBegin;
15769e15a41SShri Abhyankar   /* numeric factorization of A' */
15869e15a41SShri Abhyankar   /* ----------------------------*/
15969e15a41SShri Abhyankar 
16069e15a41SShri Abhyankar   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) {
16169e15a41SShri Abhyankar     klu_K_free_numeric(&lu->Numeric,&lu->Common);
16269e15a41SShri Abhyankar   }
16369e15a41SShri Abhyankar   lu->Numeric = klu_K_factor(ai,aj,(PetscReal*)av,lu->Symbolic,&lu->Common);
164*5f80ce2aSJacob Faibussowitsch   PetscCheck(lu->Numeric,PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Numeric factorization failed");
16569e15a41SShri Abhyankar 
16669e15a41SShri Abhyankar   lu->flg                = SAME_NONZERO_PATTERN;
16769e15a41SShri Abhyankar   lu->CleanUpKLU         = PETSC_TRUE;
16869e15a41SShri Abhyankar   F->ops->solve          = MatSolve_KLU;
16969e15a41SShri Abhyankar   F->ops->solvetranspose = MatSolveTranspose_KLU;
17069e15a41SShri Abhyankar   PetscFunctionReturn(0);
17169e15a41SShri Abhyankar }
17269e15a41SShri Abhyankar 
17369e15a41SShri Abhyankar static PetscErrorCode MatLUFactorSymbolic_KLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
17469e15a41SShri Abhyankar {
17569e15a41SShri Abhyankar   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
176569c4379SBarry Smith   Mat_KLU        *lu = (Mat_KLU*)(F->data);
17769e15a41SShri Abhyankar   PetscInt       i,*ai = a->i,*aj = a->j,m=A->rmap->n,n=A->cmap->n;
17869e15a41SShri Abhyankar   const PetscInt *ra,*ca;
17969e15a41SShri Abhyankar 
18069e15a41SShri Abhyankar   PetscFunctionBegin;
18169e15a41SShri Abhyankar   if (lu->PetscMatOrdering) {
182*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(r,&ra));
183*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(c,&ca));
184*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc2(m,&lu->perm_r,n,&lu->perm_c));
18569e15a41SShri Abhyankar     /* we cannot simply memcpy on 64 bit archs */
18669e15a41SShri Abhyankar     for (i = 0; i < m; i++) lu->perm_r[i] = ra[i];
18769e15a41SShri Abhyankar     for (i = 0; i < n; i++) lu->perm_c[i] = ca[i];
188*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(r,&ra));
189*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(c,&ca));
19069e15a41SShri Abhyankar   }
19169e15a41SShri Abhyankar 
19269e15a41SShri Abhyankar   /* symbolic factorization of A' */
19369e15a41SShri Abhyankar   /* ---------------------------------------------------------------------- */
1944ac6704cSBarry Smith   if (r) {
1954ac6704cSBarry Smith     lu->PetscMatOrdering = PETSC_TRUE;
19669e15a41SShri Abhyankar     lu->Symbolic = klu_K_analyze_given(n,ai,aj,lu->perm_c,lu->perm_r,&lu->Common);
19769e15a41SShri Abhyankar   } else { /* use klu internal ordering */
19869e15a41SShri Abhyankar     lu->Symbolic = klu_K_analyze(n,ai,aj,&lu->Common);
19969e15a41SShri Abhyankar   }
200*5f80ce2aSJacob Faibussowitsch   PetscCheck(lu->Symbolic,PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Symbolic Factorization failed");
20169e15a41SShri Abhyankar 
20269e15a41SShri Abhyankar   lu->flg                   = DIFFERENT_NONZERO_PATTERN;
20369e15a41SShri Abhyankar   lu->CleanUpKLU            = PETSC_TRUE;
20469e15a41SShri Abhyankar   (F)->ops->lufactornumeric = MatLUFactorNumeric_KLU;
20569e15a41SShri Abhyankar   PetscFunctionReturn(0);
20669e15a41SShri Abhyankar }
20769e15a41SShri Abhyankar 
208860c79edSBarry Smith static PetscErrorCode MatView_Info_KLU(Mat A,PetscViewer viewer)
20969e15a41SShri Abhyankar {
210569c4379SBarry Smith   Mat_KLU       *lu= (Mat_KLU*)A->data;
21169e15a41SShri Abhyankar   klu_K_numeric *Numeric=(klu_K_numeric*)lu->Numeric;
21269e15a41SShri Abhyankar 
21369e15a41SShri Abhyankar   PetscFunctionBegin;
214*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"KLU stats:\n"));
215*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"  Number of diagonal blocks: %" PetscInt_FMT "\n",(PetscInt)(Numeric->nblocks)));
216*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"  Total nonzeros=%" PetscInt_FMT "\n",(PetscInt)(Numeric->lnz+Numeric->unz)));
217*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"KLU runtime parameters:\n"));
21869e15a41SShri Abhyankar   /* Control parameters used by numeric factorization */
219*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"  Partial pivoting tolerance: %g\n",lu->Common.tol));
22069e15a41SShri Abhyankar   /* BTF preordering */
221*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"  BTF preordering enabled: %" PetscInt_FMT "\n",(PetscInt)(lu->Common.btf)));
22269e15a41SShri Abhyankar   /* mat ordering */
22369e15a41SShri Abhyankar   if (!lu->PetscMatOrdering) {
224*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"  Ordering: %s (not using the PETSc ordering)\n",KluOrderingTypes[(int)lu->Common.ordering]));
22569e15a41SShri Abhyankar   }
22669e15a41SShri Abhyankar   /* matrix row scaling */
227*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer, "  Matrix row scaling: %s\n",scale[(int)lu->Common.scale]));
22869e15a41SShri Abhyankar   PetscFunctionReturn(0);
22969e15a41SShri Abhyankar }
23069e15a41SShri Abhyankar 
23169e15a41SShri Abhyankar static PetscErrorCode MatView_KLU(Mat A,PetscViewer viewer)
23269e15a41SShri Abhyankar {
23369e15a41SShri Abhyankar   PetscBool         iascii;
23469e15a41SShri Abhyankar   PetscViewerFormat format;
23569e15a41SShri Abhyankar 
23669e15a41SShri Abhyankar   PetscFunctionBegin;
237*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
23869e15a41SShri Abhyankar   if (iascii) {
239*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerGetFormat(viewer,&format));
24069e15a41SShri Abhyankar     if (format == PETSC_VIEWER_ASCII_INFO) {
241*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatView_Info_KLU(A,viewer));
24269e15a41SShri Abhyankar     }
24369e15a41SShri Abhyankar   }
24469e15a41SShri Abhyankar   PetscFunctionReturn(0);
24569e15a41SShri Abhyankar }
24669e15a41SShri Abhyankar 
247ea799195SBarry Smith PetscErrorCode MatFactorGetSolverType_seqaij_klu(Mat A,MatSolverType *type)
24869e15a41SShri Abhyankar {
24969e15a41SShri Abhyankar   PetscFunctionBegin;
25069e15a41SShri Abhyankar   *type = MATSOLVERKLU;
25169e15a41SShri Abhyankar   PetscFunctionReturn(0);
25269e15a41SShri Abhyankar }
25369e15a41SShri Abhyankar 
25469e15a41SShri Abhyankar /*MC
25569e15a41SShri Abhyankar   MATSOLVERKLU = "klu" - A matrix type providing direct solvers (LU) for sequential matrices
25669e15a41SShri Abhyankar   via the external package KLU.
25769e15a41SShri Abhyankar 
258057b0833SShri Abhyankar   ./configure --download-suitesparse to install PETSc to use KLU
25969e15a41SShri Abhyankar 
2603ca39a21SBarry Smith   Use -pc_type lu -pc_factor_mat_solver_type klu to use this direct solver
261c2b89b5dSBarry Smith 
26269e15a41SShri Abhyankar   Consult KLU documentation for more information on the options database keys below.
26369e15a41SShri Abhyankar 
26469e15a41SShri Abhyankar   Options Database Keys:
26569e15a41SShri Abhyankar + -mat_klu_pivot_tol <0.001>                  - Partial pivoting tolerance
26669e15a41SShri Abhyankar . -mat_klu_use_btf <1>                        - Use BTF preordering
26769e15a41SShri Abhyankar . -mat_klu_ordering <AMD>                     - KLU reordering scheme to reduce fill-in (choose one of) AMD COLAMD PETSC
26869e15a41SShri Abhyankar - -mat_klu_row_scale <NONE>                   - Matrix row scaling (choose one of) NONE SUM MAX
26969e15a41SShri Abhyankar 
270a364b7d2SBarry Smith    Note: KLU is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html
271a364b7d2SBarry Smith 
27269e15a41SShri Abhyankar    Level: beginner
27369e15a41SShri Abhyankar 
2743ca39a21SBarry Smith .seealso: PCLU, MATSOLVERUMFPACK, MATSOLVERCHOLMOD, PCFactorSetMatSolverType(), MatSolverType
27569e15a41SShri Abhyankar M*/
27669e15a41SShri Abhyankar 
277db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat A,MatFactorType ftype,Mat *F)
27869e15a41SShri Abhyankar {
27969e15a41SShri Abhyankar   Mat            B;
28069e15a41SShri Abhyankar   Mat_KLU       *lu;
28169e15a41SShri Abhyankar   PetscErrorCode ierr;
2824ac6704cSBarry Smith   PetscInt       m=A->rmap->n,n=A->cmap->n,idx = 0,status;
28369e15a41SShri Abhyankar   PetscBool      flg;
28469e15a41SShri Abhyankar 
28569e15a41SShri Abhyankar   PetscFunctionBegin;
28669e15a41SShri Abhyankar   /* Create the factorization matrix F */
287*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PetscObjectComm((PetscObject)A),&B));
288*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n));
289*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrallocpy("klu",&((PetscObject)B)->type_name));
290*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(B));
291569c4379SBarry Smith 
292*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNewLog(B,&lu));
29369e15a41SShri Abhyankar 
294569c4379SBarry Smith   B->data                  = lu;
295569c4379SBarry Smith   B->ops->getinfo          = MatGetInfo_External;
29669e15a41SShri Abhyankar   B->ops->lufactorsymbolic = MatLUFactorSymbolic_KLU;
29769e15a41SShri Abhyankar   B->ops->destroy          = MatDestroy_KLU;
29869e15a41SShri Abhyankar   B->ops->view             = MatView_KLU;
29969e15a41SShri Abhyankar 
300*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_klu));
30169e15a41SShri Abhyankar 
30269e15a41SShri Abhyankar   B->factortype   = MAT_FACTOR_LU;
30369e15a41SShri Abhyankar   B->assembled    = PETSC_TRUE;           /* required by -ksp_view */
30469e15a41SShri Abhyankar   B->preallocated = PETSC_TRUE;
30569e15a41SShri Abhyankar 
306*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(B->solvertype));
307*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrallocpy(MATSOLVERKLU,&B->solvertype));
308f73b0415SBarry Smith   B->canuseordering = PETSC_TRUE;
309*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&B->preferredordering[MAT_FACTOR_LU]));
31000c67f3bSHong Zhang 
31169e15a41SShri Abhyankar   /* initializations */
31269e15a41SShri Abhyankar   /* ------------------------------------------------*/
31369e15a41SShri Abhyankar   /* get the default control parameters */
31469e15a41SShri Abhyankar   status = klu_K_defaults(&lu->Common);
315*5f80ce2aSJacob Faibussowitsch   PetscCheck(status > 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Initialization failed");
3166c4ed002SBarry Smith 
31769e15a41SShri Abhyankar   lu->Common.scale = 0; /* No row scaling */
31869e15a41SShri Abhyankar 
31969e15a41SShri Abhyankar   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"KLU Options","Mat");CHKERRQ(ierr);
32069e15a41SShri Abhyankar   /* Partial pivoting tolerance */
321*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-mat_klu_pivot_tol","Partial pivoting tolerance","None",lu->Common.tol,&lu->Common.tol,NULL));
32269e15a41SShri Abhyankar   /* BTF pre-ordering */
323*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-mat_klu_use_btf","Enable BTF preordering","None",(PetscInt)lu->Common.btf,(PetscInt*)&lu->Common.btf,NULL));
32469e15a41SShri Abhyankar   /* Matrix reordering */
325*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsEList("-mat_klu_ordering","Internal ordering method","None",KluOrderingTypes,sizeof(KluOrderingTypes)/sizeof(KluOrderingTypes[0]),KluOrderingTypes[0],&idx,&flg));
3264ac6704cSBarry Smith   lu->Common.ordering = (int)idx;
32769e15a41SShri Abhyankar   /* Matrix row scaling */
328*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsEList("-mat_klu_row_scale","Matrix row scaling","None",scale,3,scale[0],&idx,&flg));
329*5f80ce2aSJacob Faibussowitsch   ierr = PetscOptionsEnd();CHKERRQ(ierr);
33069e15a41SShri Abhyankar   *F = B;
33169e15a41SShri Abhyankar   PetscFunctionReturn(0);
33269e15a41SShri Abhyankar }
333