xref: /petsc/src/mat/impls/aij/seq/klu/klu.c (revision 69e15a41f829c371e4f5c4ef2075fdc7e71d4503)
1*69e15a41SShri Abhyankar 
2*69e15a41SShri Abhyankar /*
3*69e15a41SShri Abhyankar    Provides an interface to the KLUv1.2 sparse solver
4*69e15a41SShri Abhyankar 
5*69e15a41SShri Abhyankar    When build with PETSC_USE_64BIT_INDICES this will use SuiteSparse_long as the
6*69e15a41SShri Abhyankar    integer type in KLU, otherwise it will use int. This means
7*69e15a41SShri Abhyankar    all integers in this file are simply declared as PetscInt. Also it means
8*69e15a41SShri Abhyankar    that KLU SuiteSparse_long version MUST be built with 64 bit integers when used.
9*69e15a41SShri Abhyankar 
10*69e15a41SShri Abhyankar */
11*69e15a41SShri Abhyankar #include <../src/mat/impls/aij/seq/aij.h>
12*69e15a41SShri Abhyankar 
13*69e15a41SShri Abhyankar #if defined(PETSC_USE_64BIT_INDICES)
14*69e15a41SShri Abhyankar #define klu_K_defaults                klu_l_defaults
15*69e15a41SShri Abhyankar #define klu_K_analyze                 klu_l_analyze
16*69e15a41SShri Abhyankar #define klu_K_analyze_given           klu_l_analyze_given
17*69e15a41SShri Abhyankar #define klu_K_free_symbolic           klu_l_free_symbolic
18*69e15a41SShri Abhyankar #define klu_K_free_numeric            klu_l_free_numeric
19*69e15a41SShri Abhyankar #define klu_K_common                  klu_l_common
20*69e15a41SShri Abhyankar #define klu_K_symbolic                klu_l_symbolic
21*69e15a41SShri Abhyankar #define klu_K_numeric                 klu_l_numeric
22*69e15a41SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
23*69e15a41SShri Abhyankar #define klu_K_factor                  klu_zl_factor
24*69e15a41SShri Abhyankar #define klu_K_solve                   klu_zl_solve
25*69e15a41SShri Abhyankar #define klu_K_tsolve                  klu_zl_tsolve
26*69e15a41SShri Abhyankar #define klu_K_refactor                klu_zl_refactor
27*69e15a41SShri Abhyankar #define klu_K_sort                    klu_zl_sort
28*69e15a41SShri Abhyankar #define klu_K_flops                   klu_zl_flops
29*69e15a41SShri Abhyankar #define klu_K_rgrowth                 klu_zl_rgrowth
30*69e15a41SShri Abhyankar #define klu_K_condest                 klu_zl_condest
31*69e15a41SShri Abhyankar #define klu_K_rcond                   klu_zl_rcond
32*69e15a41SShri Abhyankar #define klu_K_scale                   klu_zl_scale
33*69e15a41SShri Abhyankar #else
34*69e15a41SShri Abhyankar #define klu_K_factor                  klu_l_factor
35*69e15a41SShri Abhyankar #define klu_K_solve                   klu_l_solve
36*69e15a41SShri Abhyankar #define klu_K_tsolve                  klu_l_tsolve
37*69e15a41SShri Abhyankar #define klu_K_refactor                klu_l_refactor
38*69e15a41SShri Abhyankar #define klu_K_sort                    klu_l_sort
39*69e15a41SShri Abhyankar #define klu_K_flops                   klu_l_flops
40*69e15a41SShri Abhyankar #define klu_K_rgrowth                 klu_l_rgrowth
41*69e15a41SShri Abhyankar #define klu_K_condest                 klu_l_condest
42*69e15a41SShri Abhyankar #define klu_K_rcond                   klu_l_rcond
43*69e15a41SShri Abhyankar #define klu_K_scale                   klu_l_scale
44*69e15a41SShri Abhyankar #endif
45*69e15a41SShri Abhyankar #else
46*69e15a41SShri Abhyankar #define klu_K_defaults                klu_defaults
47*69e15a41SShri Abhyankar #define klu_K_analyze                 klu_analyze
48*69e15a41SShri Abhyankar #define klu_K_analyze_given           klu_analyze_given
49*69e15a41SShri Abhyankar #define klu_K_free_symbolic           klu_free_symbolic
50*69e15a41SShri Abhyankar #define klu_K_free_numeric            klu_free_numeric
51*69e15a41SShri Abhyankar #define klu_K_common                  klu_common
52*69e15a41SShri Abhyankar #define klu_K_symbolic                klu_symbolic
53*69e15a41SShri Abhyankar #define klu_K_numeric                 klu_numeric
54*69e15a41SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
55*69e15a41SShri Abhyankar #define klu_K_factor                  klu_z_factor
56*69e15a41SShri Abhyankar #define klu_K_solve                   klu_z_solve
57*69e15a41SShri Abhyankar #define klu_K_tsolve                  klu_z_tsolve
58*69e15a41SShri Abhyankar #define klu_K_refactor                klu_z_refactor
59*69e15a41SShri Abhyankar #define klu_K_sort                    klu_z_sort
60*69e15a41SShri Abhyankar #define klu_K_flops                   klu_z_flops
61*69e15a41SShri Abhyankar #define klu_K_rgrowth                 klu_z_rgrowth
62*69e15a41SShri Abhyankar #define klu_K_condest                 klu_z_condest
63*69e15a41SShri Abhyankar #define klu_K_rcond                   klu_z_rcond
64*69e15a41SShri Abhyankar #define klu_K_scale                   klu_z_scale
65*69e15a41SShri Abhyankar #else
66*69e15a41SShri Abhyankar #define klu_K_factor                  klu_factor
67*69e15a41SShri Abhyankar #define klu_K_solve                   klu_solve
68*69e15a41SShri Abhyankar #define klu_K_tsolve                  klu_tsolve
69*69e15a41SShri Abhyankar #define klu_K_refactor                klu_refactor
70*69e15a41SShri Abhyankar #define klu_K_sort                    klu_sort
71*69e15a41SShri Abhyankar #define klu_K_flops                   klu_flops
72*69e15a41SShri Abhyankar #define klu_K_rgrowth                 klu_rgrowth
73*69e15a41SShri Abhyankar #define klu_K_condest                 klu_condest
74*69e15a41SShri Abhyankar #define klu_K_rcond                   klu_rcond
75*69e15a41SShri Abhyankar #define klu_K_scale                   klu_scale
76*69e15a41SShri Abhyankar #endif
77*69e15a41SShri Abhyankar #endif
78*69e15a41SShri Abhyankar 
79*69e15a41SShri Abhyankar 
80*69e15a41SShri Abhyankar #define SuiteSparse_long long long
81*69e15a41SShri Abhyankar #define SuiteSparse_long_max LONG_LONG_MAX
82*69e15a41SShri Abhyankar #define SuiteSparse_long_id "%lld"
83*69e15a41SShri Abhyankar 
84*69e15a41SShri Abhyankar EXTERN_C_BEGIN
85*69e15a41SShri Abhyankar #include <klu.h>
86*69e15a41SShri Abhyankar EXTERN_C_END
87*69e15a41SShri Abhyankar 
88*69e15a41SShri Abhyankar static const char *KluOrderingTypes[] = {"AMD","COLAMD","PETSC"};
89*69e15a41SShri Abhyankar static const char *scale[] ={"NONE","SUM","MAX"};
90*69e15a41SShri Abhyankar 
91*69e15a41SShri Abhyankar typedef struct {
92*69e15a41SShri Abhyankar   klu_K_common   Common;
93*69e15a41SShri Abhyankar   klu_K_symbolic *Symbolic;
94*69e15a41SShri Abhyankar   klu_K_numeric  *Numeric;
95*69e15a41SShri Abhyankar   PetscInt     *perm_c,*perm_r;
96*69e15a41SShri Abhyankar   MatStructure flg;
97*69e15a41SShri Abhyankar   PetscBool    PetscMatOrdering;
98*69e15a41SShri Abhyankar 
99*69e15a41SShri Abhyankar   /* Flag to clean up KLU objects during Destroy */
100*69e15a41SShri Abhyankar   PetscBool CleanUpKLU;
101*69e15a41SShri Abhyankar } Mat_KLU;
102*69e15a41SShri Abhyankar 
103*69e15a41SShri Abhyankar #undef __FUNCT__
104*69e15a41SShri Abhyankar #define __FUNCT__ "MatDestroy_KLU"
105*69e15a41SShri Abhyankar static PetscErrorCode MatDestroy_KLU(Mat A)
106*69e15a41SShri Abhyankar {
107*69e15a41SShri Abhyankar   PetscErrorCode ierr;
108*69e15a41SShri Abhyankar   Mat_KLU    *lu=(Mat_KLU*)A->spptr;
109*69e15a41SShri Abhyankar 
110*69e15a41SShri Abhyankar   PetscFunctionBegin;
111*69e15a41SShri Abhyankar   if (lu && lu->CleanUpKLU) {
112*69e15a41SShri Abhyankar     klu_K_free_symbolic(&lu->Symbolic,&lu->Common);
113*69e15a41SShri Abhyankar     klu_K_free_numeric(&lu->Numeric,&lu->Common);
114*69e15a41SShri Abhyankar     ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
115*69e15a41SShri Abhyankar     ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
116*69e15a41SShri Abhyankar   }
117*69e15a41SShri Abhyankar   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
118*69e15a41SShri Abhyankar   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
119*69e15a41SShri Abhyankar   PetscFunctionReturn(0);
120*69e15a41SShri Abhyankar }
121*69e15a41SShri Abhyankar 
122*69e15a41SShri Abhyankar #undef __FUNCT__
123*69e15a41SShri Abhyankar #define __FUNCT__ "MatSolveTranspose_KLU"
124*69e15a41SShri Abhyankar static PetscErrorCode MatSolveTranspose_KLU(Mat A,Vec b,Vec x)
125*69e15a41SShri Abhyankar {
126*69e15a41SShri Abhyankar   Mat_KLU       *lu = (Mat_KLU*)A->spptr;
127*69e15a41SShri Abhyankar   PetscScalar    *xa;
128*69e15a41SShri Abhyankar   PetscErrorCode ierr;
129*69e15a41SShri Abhyankar   PetscInt       status;
130*69e15a41SShri Abhyankar 
131*69e15a41SShri Abhyankar   PetscFunctionBegin;
132*69e15a41SShri Abhyankar   /* KLU uses a column major format, solve Ax = b by klu_*_solve */
133*69e15a41SShri Abhyankar   /* ----------------------------------*/
134*69e15a41SShri Abhyankar   ierr = VecCopy(b,x); /* klu_solve stores the solution in rhs */
135*69e15a41SShri Abhyankar   ierr = VecGetArray(x,&xa);
136*69e15a41SShri Abhyankar   status = klu_K_solve(lu->Symbolic,lu->Numeric,A->rmap->n,1,(PetscReal*)xa,&lu->Common);
137*69e15a41SShri Abhyankar   if (status != 1) {
138*69e15a41SShri Abhyankar     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Solve failed");
139*69e15a41SShri Abhyankar   }
140*69e15a41SShri Abhyankar   ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr);
141*69e15a41SShri Abhyankar   PetscFunctionReturn(0);
142*69e15a41SShri Abhyankar }
143*69e15a41SShri Abhyankar 
144*69e15a41SShri Abhyankar #undef __FUNCT__
145*69e15a41SShri Abhyankar #define __FUNCT__ "MatSolve_KLU"
146*69e15a41SShri Abhyankar static PetscErrorCode MatSolve_KLU(Mat A,Vec b,Vec x)
147*69e15a41SShri Abhyankar {
148*69e15a41SShri Abhyankar   Mat_KLU       *lu = (Mat_KLU*)A->spptr;
149*69e15a41SShri Abhyankar   PetscScalar    *xa;
150*69e15a41SShri Abhyankar   PetscErrorCode ierr;
151*69e15a41SShri Abhyankar   PetscInt       status;
152*69e15a41SShri Abhyankar 
153*69e15a41SShri Abhyankar   PetscFunctionBegin;
154*69e15a41SShri Abhyankar   /* KLU uses a column major format, solve A^Tx = b by klu_*_tsolve */
155*69e15a41SShri Abhyankar   /* ----------------------------------*/
156*69e15a41SShri Abhyankar   ierr = VecCopy(b,x); /* klu_solve stores the solution in rhs */
157*69e15a41SShri Abhyankar   ierr = VecGetArray(x,&xa);
158*69e15a41SShri Abhyankar #if defined(PETSC_USE_COMPLEX)
159*69e15a41SShri Abhyankar   PetscInt conj_solve=1;
160*69e15a41SShri Abhyankar   status = klu_K_tsolve(lu->Symbolic,lu->Numeric,A->rmap->n,1,(PetscReal*)xa,conj_solve,&lu->Common); /* conjugate solve */
161*69e15a41SShri Abhyankar #else
162*69e15a41SShri Abhyankar   status = klu_K_tsolve(lu->Symbolic,lu->Numeric,A->rmap->n,1,xa,&lu->Common);
163*69e15a41SShri Abhyankar #endif
164*69e15a41SShri Abhyankar   if (status != 1) {
165*69e15a41SShri Abhyankar     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Solve failed");
166*69e15a41SShri Abhyankar   }
167*69e15a41SShri Abhyankar   ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr);
168*69e15a41SShri Abhyankar   PetscFunctionReturn(0);
169*69e15a41SShri Abhyankar }
170*69e15a41SShri Abhyankar 
171*69e15a41SShri Abhyankar #undef __FUNCT__
172*69e15a41SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_KLU"
173*69e15a41SShri Abhyankar static PetscErrorCode MatLUFactorNumeric_KLU(Mat F,Mat A,const MatFactorInfo *info)
174*69e15a41SShri Abhyankar {
175*69e15a41SShri Abhyankar   Mat_KLU        *lu = (Mat_KLU*)(F)->spptr;
176*69e15a41SShri Abhyankar   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
177*69e15a41SShri Abhyankar   PetscInt       *ai = a->i,*aj=a->j;
178*69e15a41SShri Abhyankar   PetscScalar    *av = a->a;
179*69e15a41SShri Abhyankar 
180*69e15a41SShri Abhyankar   PetscFunctionBegin;
181*69e15a41SShri Abhyankar   /* numeric factorization of A' */
182*69e15a41SShri Abhyankar   /* ----------------------------*/
183*69e15a41SShri Abhyankar 
184*69e15a41SShri Abhyankar   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) {
185*69e15a41SShri Abhyankar     klu_K_free_numeric(&lu->Numeric,&lu->Common);
186*69e15a41SShri Abhyankar   }
187*69e15a41SShri Abhyankar   lu->Numeric = klu_K_factor(ai,aj,(PetscReal*)av,lu->Symbolic,&lu->Common);
188*69e15a41SShri Abhyankar   if(!lu->Numeric) {
189*69e15a41SShri Abhyankar     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Numeric factorization failed");
190*69e15a41SShri Abhyankar   }
191*69e15a41SShri Abhyankar 
192*69e15a41SShri Abhyankar   lu->flg                = SAME_NONZERO_PATTERN;
193*69e15a41SShri Abhyankar   lu->CleanUpKLU         = PETSC_TRUE;
194*69e15a41SShri Abhyankar   F->ops->solve          = MatSolve_KLU;
195*69e15a41SShri Abhyankar   F->ops->solvetranspose = MatSolveTranspose_KLU;
196*69e15a41SShri Abhyankar   PetscFunctionReturn(0);
197*69e15a41SShri Abhyankar }
198*69e15a41SShri Abhyankar 
199*69e15a41SShri Abhyankar #undef __FUNCT__
200*69e15a41SShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_KLU"
201*69e15a41SShri Abhyankar static PetscErrorCode MatLUFactorSymbolic_KLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
202*69e15a41SShri Abhyankar {
203*69e15a41SShri Abhyankar   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
204*69e15a41SShri Abhyankar   Mat_KLU       *lu = (Mat_KLU*)(F->spptr);
205*69e15a41SShri Abhyankar   PetscErrorCode ierr;
206*69e15a41SShri Abhyankar   PetscInt       i,*ai = a->i,*aj = a->j,m=A->rmap->n,n=A->cmap->n;
207*69e15a41SShri Abhyankar   const PetscInt *ra,*ca;
208*69e15a41SShri Abhyankar 
209*69e15a41SShri Abhyankar   PetscFunctionBegin;
210*69e15a41SShri Abhyankar   if (lu->PetscMatOrdering) {
211*69e15a41SShri Abhyankar     ierr = ISGetIndices(r,&ra);CHKERRQ(ierr);
212*69e15a41SShri Abhyankar     ierr = ISGetIndices(c,&ca);CHKERRQ(ierr);
213*69e15a41SShri Abhyankar     ierr = PetscMalloc2(m,&lu->perm_c,n,&lu->perm_r);CHKERRQ(ierr);
214*69e15a41SShri Abhyankar     /* we cannot simply memcpy on 64 bit archs */
215*69e15a41SShri Abhyankar     for (i = 0; i < m; i++) lu->perm_r[i] = ra[i];
216*69e15a41SShri Abhyankar     for (i = 0; i < n; i++) lu->perm_c[i] = ca[i];
217*69e15a41SShri Abhyankar     ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr);
218*69e15a41SShri Abhyankar     ierr = ISRestoreIndices(c,&ca);CHKERRQ(ierr);
219*69e15a41SShri Abhyankar   }
220*69e15a41SShri Abhyankar 
221*69e15a41SShri Abhyankar   /* symbolic factorization of A' */
222*69e15a41SShri Abhyankar   /* ---------------------------------------------------------------------- */
223*69e15a41SShri Abhyankar   if (lu->PetscMatOrdering) { /* use Petsc ordering */
224*69e15a41SShri Abhyankar     lu->Symbolic = klu_K_analyze_given(n,ai,aj,lu->perm_c,lu->perm_r,&lu->Common);
225*69e15a41SShri Abhyankar   } else { /* use klu internal ordering */
226*69e15a41SShri Abhyankar     lu->Symbolic = klu_K_analyze(n,ai,aj,&lu->Common);
227*69e15a41SShri Abhyankar   }
228*69e15a41SShri Abhyankar   if (!lu->Symbolic) {
229*69e15a41SShri Abhyankar     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Symbolic Factorization failed");
230*69e15a41SShri Abhyankar   }
231*69e15a41SShri Abhyankar 
232*69e15a41SShri Abhyankar   lu->flg                   = DIFFERENT_NONZERO_PATTERN;
233*69e15a41SShri Abhyankar   lu->CleanUpKLU            = PETSC_TRUE;
234*69e15a41SShri Abhyankar   (F)->ops->lufactornumeric = MatLUFactorNumeric_KLU;
235*69e15a41SShri Abhyankar   PetscFunctionReturn(0);
236*69e15a41SShri Abhyankar }
237*69e15a41SShri Abhyankar 
238*69e15a41SShri Abhyankar #undef __FUNCT__
239*69e15a41SShri Abhyankar #define __FUNCT__ "MatFactorInfo_KLU"
240*69e15a41SShri Abhyankar static PetscErrorCode MatFactorInfo_KLU(Mat A,PetscViewer viewer)
241*69e15a41SShri Abhyankar {
242*69e15a41SShri Abhyankar   Mat_KLU       *lu= (Mat_KLU*)A->spptr;
243*69e15a41SShri Abhyankar   klu_K_numeric *Numeric=(klu_K_numeric*)lu->Numeric;
244*69e15a41SShri Abhyankar   PetscErrorCode ierr;
245*69e15a41SShri Abhyankar 
246*69e15a41SShri Abhyankar   PetscFunctionBegin;
247*69e15a41SShri Abhyankar   /* check if matrix is KLU type */
248*69e15a41SShri Abhyankar   if (A->ops->solve != MatSolve_KLU) PetscFunctionReturn(0);
249*69e15a41SShri Abhyankar 
250*69e15a41SShri Abhyankar   ierr = PetscViewerASCIIPrintf(viewer,"KLU stats:\n");CHKERRQ(ierr);
251*69e15a41SShri Abhyankar   ierr = PetscViewerASCIIPrintf(viewer,"  Number of diagonal blocks: %d\n",Numeric->nblocks);
252*69e15a41SShri Abhyankar   ierr = PetscViewerASCIIPrintf(viewer,"  Total nonzeros=%d\n",Numeric->lnz+Numeric->unz);CHKERRQ(ierr);
253*69e15a41SShri Abhyankar 
254*69e15a41SShri Abhyankar   ierr = PetscViewerASCIIPrintf(viewer,"KLU runtime parameters:\n");CHKERRQ(ierr);
255*69e15a41SShri Abhyankar 
256*69e15a41SShri Abhyankar   /* Control parameters used by numeric factorization */
257*69e15a41SShri Abhyankar   ierr = PetscViewerASCIIPrintf(viewer,"  Partial pivoting tolerance: %g\n",lu->Common.tol);CHKERRQ(ierr);
258*69e15a41SShri Abhyankar   /* BTF preordering */
259*69e15a41SShri Abhyankar   ierr = PetscViewerASCIIPrintf(viewer,"  BTF preordering enabled: %d\n",lu->Common.btf);CHKERRQ(ierr);
260*69e15a41SShri Abhyankar   /* mat ordering */
261*69e15a41SShri Abhyankar   if (!lu->PetscMatOrdering) {
262*69e15a41SShri Abhyankar     ierr = PetscViewerASCIIPrintf(viewer,"  Ordering: %s (not using the PETSc ordering)\n",KluOrderingTypes[(int)lu->Common.ordering]);CHKERRQ(ierr);
263*69e15a41SShri Abhyankar   } else {
264*69e15a41SShri Abhyankar     ierr = PetscViewerASCIIPrintf(viewer,"  Using PETSc ordering\n");CHKERRQ(ierr);
265*69e15a41SShri Abhyankar   }
266*69e15a41SShri Abhyankar   /* matrix row scaling */
267*69e15a41SShri Abhyankar   ierr = PetscViewerASCIIPrintf(viewer, "  Matrix row scaling: %s\n",scale[(int)lu->Common.scale]);CHKERRQ(ierr);
268*69e15a41SShri Abhyankar   PetscFunctionReturn(0);
269*69e15a41SShri Abhyankar }
270*69e15a41SShri Abhyankar 
271*69e15a41SShri Abhyankar #undef __FUNCT__
272*69e15a41SShri Abhyankar #define __FUNCT__ "MatView_KLU"
273*69e15a41SShri Abhyankar static PetscErrorCode MatView_KLU(Mat A,PetscViewer viewer)
274*69e15a41SShri Abhyankar {
275*69e15a41SShri Abhyankar   PetscErrorCode    ierr;
276*69e15a41SShri Abhyankar   PetscBool         iascii;
277*69e15a41SShri Abhyankar   PetscViewerFormat format;
278*69e15a41SShri Abhyankar 
279*69e15a41SShri Abhyankar   PetscFunctionBegin;
280*69e15a41SShri Abhyankar   ierr = MatView_SeqAIJ(A,viewer);CHKERRQ(ierr);
281*69e15a41SShri Abhyankar 
282*69e15a41SShri Abhyankar   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
283*69e15a41SShri Abhyankar   if (iascii) {
284*69e15a41SShri Abhyankar     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
285*69e15a41SShri Abhyankar     if (format == PETSC_VIEWER_ASCII_INFO) {
286*69e15a41SShri Abhyankar       ierr = MatFactorInfo_KLU(A,viewer);CHKERRQ(ierr);
287*69e15a41SShri Abhyankar     }
288*69e15a41SShri Abhyankar   }
289*69e15a41SShri Abhyankar   PetscFunctionReturn(0);
290*69e15a41SShri Abhyankar }
291*69e15a41SShri Abhyankar 
292*69e15a41SShri Abhyankar #undef __FUNCT__
293*69e15a41SShri Abhyankar #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_klu"
294*69e15a41SShri Abhyankar PetscErrorCode MatFactorGetSolverPackage_seqaij_klu(Mat A,const MatSolverPackage *type)
295*69e15a41SShri Abhyankar {
296*69e15a41SShri Abhyankar   PetscFunctionBegin;
297*69e15a41SShri Abhyankar   *type = MATSOLVERKLU;
298*69e15a41SShri Abhyankar   PetscFunctionReturn(0);
299*69e15a41SShri Abhyankar }
300*69e15a41SShri Abhyankar 
301*69e15a41SShri Abhyankar 
302*69e15a41SShri Abhyankar /*MC
303*69e15a41SShri Abhyankar   MATSOLVERKLU = "klu" - A matrix type providing direct solvers (LU) for sequential matrices
304*69e15a41SShri Abhyankar   via the external package KLU.
305*69e15a41SShri Abhyankar 
306*69e15a41SShri Abhyankar   ./configure --download-SuiteSparse --SuiteSparse-enable-klu to install PETSc to use KLU
307*69e15a41SShri Abhyankar 
308*69e15a41SShri Abhyankar   Consult KLU documentation for more information on the options database keys below.
309*69e15a41SShri Abhyankar 
310*69e15a41SShri Abhyankar   Options Database Keys:
311*69e15a41SShri Abhyankar + -mat_klu_pivot_tol <0.001>                  - Partial pivoting tolerance
312*69e15a41SShri Abhyankar . -mat_klu_use_btf <1>                        - Use BTF preordering
313*69e15a41SShri Abhyankar . -mat_klu_ordering <AMD>                     - KLU reordering scheme to reduce fill-in (choose one of) AMD COLAMD PETSC
314*69e15a41SShri Abhyankar - -mat_klu_row_scale <NONE>                   - Matrix row scaling (choose one of) NONE SUM MAX
315*69e15a41SShri Abhyankar 
316*69e15a41SShri Abhyankar    Level: beginner
317*69e15a41SShri Abhyankar 
318*69e15a41SShri Abhyankar .seealso: PCLU, MATSOLVERUMFPACK, MATSOLVERCHOLMOD, PCFactorSetMatSolverPackage(), MatSolverPackage
319*69e15a41SShri Abhyankar M*/
320*69e15a41SShri Abhyankar 
321*69e15a41SShri Abhyankar #undef __FUNCT__
322*69e15a41SShri Abhyankar #define __FUNCT__ "MatGetFactor_seqaij_klu"
323*69e15a41SShri Abhyankar PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat A,MatFactorType ftype,Mat *F)
324*69e15a41SShri Abhyankar {
325*69e15a41SShri Abhyankar   Mat            B;
326*69e15a41SShri Abhyankar   Mat_KLU       *lu;
327*69e15a41SShri Abhyankar   PetscErrorCode ierr;
328*69e15a41SShri Abhyankar   PetscInt       m=A->rmap->n,n=A->cmap->n,idx,status;
329*69e15a41SShri Abhyankar   PetscBool      flg;
330*69e15a41SShri Abhyankar 
331*69e15a41SShri Abhyankar   PetscFunctionBegin;
332*69e15a41SShri Abhyankar   /* Create the factorization matrix F */
333*69e15a41SShri Abhyankar   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
334*69e15a41SShri Abhyankar   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
335*69e15a41SShri Abhyankar   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
336*69e15a41SShri Abhyankar   ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
337*69e15a41SShri Abhyankar   ierr = PetscNewLog(B,&lu);CHKERRQ(ierr);
338*69e15a41SShri Abhyankar 
339*69e15a41SShri Abhyankar   B->spptr                 = lu;
340*69e15a41SShri Abhyankar   B->ops->lufactorsymbolic = MatLUFactorSymbolic_KLU;
341*69e15a41SShri Abhyankar   B->ops->destroy          = MatDestroy_KLU;
342*69e15a41SShri Abhyankar   B->ops->view             = MatView_KLU;
343*69e15a41SShri Abhyankar 
344*69e15a41SShri Abhyankar   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_klu);CHKERRQ(ierr);
345*69e15a41SShri Abhyankar 
346*69e15a41SShri Abhyankar   B->factortype   = MAT_FACTOR_LU;
347*69e15a41SShri Abhyankar   B->assembled    = PETSC_TRUE;           /* required by -ksp_view */
348*69e15a41SShri Abhyankar   B->preallocated = PETSC_TRUE;
349*69e15a41SShri Abhyankar 
350*69e15a41SShri Abhyankar   /* initializations */
351*69e15a41SShri Abhyankar   /* ------------------------------------------------*/
352*69e15a41SShri Abhyankar   /* get the default control parameters */
353*69e15a41SShri Abhyankar   status = klu_K_defaults(&lu->Common);
354*69e15a41SShri Abhyankar   if(status <= 0) {
355*69e15a41SShri Abhyankar     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Initialization failed");
356*69e15a41SShri Abhyankar   }
357*69e15a41SShri Abhyankar   lu->Common.scale = 0; /* No row scaling */
358*69e15a41SShri Abhyankar 
359*69e15a41SShri Abhyankar   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"KLU Options","Mat");CHKERRQ(ierr);
360*69e15a41SShri Abhyankar   /* Partial pivoting tolerance */
361*69e15a41SShri Abhyankar   ierr = PetscOptionsReal("-mat_klu_pivot_tol","Partial pivoting tolerance","None",lu->Common.tol,&lu->Common.tol,NULL);CHKERRQ(ierr);
362*69e15a41SShri Abhyankar   /* BTF pre-ordering */
363*69e15a41SShri Abhyankar   ierr = PetscOptionsInt("-mat_klu_use_btf","Enable BTF preordering","None",lu->Common.btf,&lu->Common.btf,NULL);CHKERRQ(ierr);
364*69e15a41SShri Abhyankar   /* Matrix reordering */
365*69e15a41SShri Abhyankar   ierr = PetscOptionsEList("-mat_klu_ordering","Internal ordering method","None",KluOrderingTypes,sizeof(KluOrderingTypes)/sizeof(KluOrderingTypes[0]),KluOrderingTypes[0],&idx,&flg);CHKERRQ(ierr);
366*69e15a41SShri Abhyankar   if (flg) {
367*69e15a41SShri Abhyankar     if ((int)idx == 2) lu->PetscMatOrdering = PETSC_TRUE;   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Klu perm_c) */
368*69e15a41SShri Abhyankar     else lu->Common.ordering = (int)idx;
369*69e15a41SShri Abhyankar   }
370*69e15a41SShri Abhyankar   /* Matrix row scaling */
371*69e15a41SShri Abhyankar   ierr = PetscOptionsEList("-mat_klu_row_scale","Matrix row scaling","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr);
372*69e15a41SShri Abhyankar   PetscOptionsEnd();
373*69e15a41SShri Abhyankar   *F = B;
374*69e15a41SShri Abhyankar   PetscFunctionReturn(0);
375*69e15a41SShri Abhyankar }
376