Actual source code: klu.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  2: /*
  3:    Provides an interface to the KLUv1.2 sparse solver

  5:    When build with PETSC_USE_64BIT_INDICES this will use SuiteSparse_long as the
  6:    integer type in KLU, otherwise it will use int. This means
  7:    all integers in this file are simply declared as PetscInt. Also it means
  8:    that KLU SuiteSparse_long version MUST be built with 64 bit integers when used.

 10: */
 11: #include <../src/mat/impls/aij/seq/aij.h>

 13: #if defined(PETSC_USE_64BIT_INDICES)
 14: #define klu_K_defaults                klu_l_defaults
 15: #define klu_K_analyze                 klu_l_analyze
 16: #define klu_K_analyze_given           klu_l_analyze_given
 17: #define klu_K_free_symbolic           klu_l_free_symbolic
 18: #define klu_K_free_numeric            klu_l_free_numeric
 19: #define klu_K_common                  klu_l_common
 20: #define klu_K_symbolic                klu_l_symbolic
 21: #define klu_K_numeric                 klu_l_numeric
 22: #if defined(PETSC_USE_COMPLEX)
 23: #define klu_K_factor                  klu_zl_factor
 24: #define klu_K_solve                   klu_zl_solve
 25: #define klu_K_tsolve                  klu_zl_tsolve
 26: #define klu_K_refactor                klu_zl_refactor
 27: #define klu_K_sort                    klu_zl_sort
 28: #define klu_K_flops                   klu_zl_flops
 29: #define klu_K_rgrowth                 klu_zl_rgrowth
 30: #define klu_K_condest                 klu_zl_condest
 31: #define klu_K_rcond                   klu_zl_rcond
 32: #define klu_K_scale                   klu_zl_scale
 33: #else
 34: #define klu_K_factor                  klu_l_factor
 35: #define klu_K_solve                   klu_l_solve
 36: #define klu_K_tsolve                  klu_l_tsolve
 37: #define klu_K_refactor                klu_l_refactor
 38: #define klu_K_sort                    klu_l_sort
 39: #define klu_K_flops                   klu_l_flops
 40: #define klu_K_rgrowth                 klu_l_rgrowth
 41: #define klu_K_condest                 klu_l_condest
 42: #define klu_K_rcond                   klu_l_rcond
 43: #define klu_K_scale                   klu_l_scale
 44: #endif
 45: #else
 46: #define klu_K_defaults                klu_defaults
 47: #define klu_K_analyze                 klu_analyze
 48: #define klu_K_analyze_given           klu_analyze_given
 49: #define klu_K_free_symbolic           klu_free_symbolic
 50: #define klu_K_free_numeric            klu_free_numeric
 51: #define klu_K_common                  klu_common
 52: #define klu_K_symbolic                klu_symbolic
 53: #define klu_K_numeric                 klu_numeric
 54: #if defined(PETSC_USE_COMPLEX)
 55: #define klu_K_factor                  klu_z_factor
 56: #define klu_K_solve                   klu_z_solve
 57: #define klu_K_tsolve                  klu_z_tsolve
 58: #define klu_K_refactor                klu_z_refactor
 59: #define klu_K_sort                    klu_z_sort
 60: #define klu_K_flops                   klu_z_flops
 61: #define klu_K_rgrowth                 klu_z_rgrowth
 62: #define klu_K_condest                 klu_z_condest
 63: #define klu_K_rcond                   klu_z_rcond
 64: #define klu_K_scale                   klu_z_scale
 65: #else
 66: #define klu_K_factor                  klu_factor
 67: #define klu_K_solve                   klu_solve
 68: #define klu_K_tsolve                  klu_tsolve
 69: #define klu_K_refactor                klu_refactor
 70: #define klu_K_sort                    klu_sort
 71: #define klu_K_flops                   klu_flops
 72: #define klu_K_rgrowth                 klu_rgrowth
 73: #define klu_K_condest                 klu_condest
 74: #define klu_K_rcond                   klu_rcond
 75: #define klu_K_scale                   klu_scale
 76: #endif
 77: #endif


 80: #define SuiteSparse_long long long
 81: #define SuiteSparse_long_max LONG_LONG_MAX
 82: #define SuiteSparse_long_id "%lld"

 84: EXTERN_C_BEGIN
 85: #include <klu.h>
 86: EXTERN_C_END

 88: static const char *KluOrderingTypes[] = {"AMD","COLAMD","PETSC"};
 89: static const char *scale[] ={"NONE","SUM","MAX"};

 91: typedef struct {
 92:   klu_K_common   Common;
 93:   klu_K_symbolic *Symbolic;
 94:   klu_K_numeric  *Numeric;
 95:   PetscInt     *perm_c,*perm_r;
 96:   MatStructure flg;
 97:   PetscBool    PetscMatOrdering;

 99:   /* Flag to clean up KLU objects during Destroy */
100:   PetscBool CleanUpKLU;
101: } Mat_KLU;

105: static PetscErrorCode MatDestroy_KLU(Mat A)
106: {
108:   Mat_KLU    *lu=(Mat_KLU*)A->spptr;

111:   if (lu && lu->CleanUpKLU) {
112:     klu_K_free_symbolic(&lu->Symbolic,&lu->Common);
113:     klu_K_free_numeric(&lu->Numeric,&lu->Common);
114:     PetscFree2(lu->perm_r,lu->perm_c);
115:   }
116:   PetscFree(A->spptr);
117:   MatDestroy_SeqAIJ(A);
118:   return(0);
119: }

123: static PetscErrorCode MatSolveTranspose_KLU(Mat A,Vec b,Vec x)
124: {
125:   Mat_KLU       *lu = (Mat_KLU*)A->spptr;
126:   PetscScalar    *xa;
128:   PetscInt       status;

131:   /* KLU uses a column major format, solve Ax = b by klu_*_solve */
132:   /* ----------------------------------*/
133:   VecCopy(b,x); /* klu_solve stores the solution in rhs */
134:   VecGetArray(x,&xa);
135:   status = klu_K_solve(lu->Symbolic,lu->Numeric,A->rmap->n,1,(PetscReal*)xa,&lu->Common);
136:   if (status != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Solve failed");
137:   VecRestoreArray(x,&xa);
138:   return(0);
139: }

143: static PetscErrorCode MatSolve_KLU(Mat A,Vec b,Vec x)
144: {
145:   Mat_KLU       *lu = (Mat_KLU*)A->spptr;
146:   PetscScalar    *xa;
148:   PetscInt       status;

151:   /* KLU uses a column major format, solve A^Tx = b by klu_*_tsolve */
152:   /* ----------------------------------*/
153:   VecCopy(b,x); /* klu_solve stores the solution in rhs */
154:   VecGetArray(x,&xa);
155: #if defined(PETSC_USE_COMPLEX)
156:   PetscInt conj_solve=1;
157:   status = klu_K_tsolve(lu->Symbolic,lu->Numeric,A->rmap->n,1,(PetscReal*)xa,conj_solve,&lu->Common); /* conjugate solve */
158: #else
159:   status = klu_K_tsolve(lu->Symbolic,lu->Numeric,A->rmap->n,1,xa,&lu->Common);
160: #endif  
161:   if (status != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Solve failed");
162:   VecRestoreArray(x,&xa);
163:   return(0);
164: }

168: static PetscErrorCode MatLUFactorNumeric_KLU(Mat F,Mat A,const MatFactorInfo *info)
169: {
170:   Mat_KLU        *lu = (Mat_KLU*)(F)->spptr;
171:   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
172:   PetscInt       *ai = a->i,*aj=a->j;
173:   PetscScalar    *av = a->a;

176:   /* numeric factorization of A' */
177:   /* ----------------------------*/

179:   if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) {
180:     klu_K_free_numeric(&lu->Numeric,&lu->Common);
181:   }
182:   lu->Numeric = klu_K_factor(ai,aj,(PetscReal*)av,lu->Symbolic,&lu->Common);
183:   if(!lu->Numeric) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Numeric factorization failed");

185:   lu->flg                = SAME_NONZERO_PATTERN;
186:   lu->CleanUpKLU         = PETSC_TRUE;
187:   F->ops->solve          = MatSolve_KLU;
188:   F->ops->solvetranspose = MatSolveTranspose_KLU;
189:   return(0);
190: }

194: static PetscErrorCode MatLUFactorSymbolic_KLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
195: {
196:   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
197:   Mat_KLU       *lu = (Mat_KLU*)(F->spptr);
199:   PetscInt       i,*ai = a->i,*aj = a->j,m=A->rmap->n,n=A->cmap->n;
200:   const PetscInt *ra,*ca;

203:   if (lu->PetscMatOrdering) {
204:     ISGetIndices(r,&ra);
205:     ISGetIndices(c,&ca);
206:     PetscMalloc2(m,&lu->perm_r,n,&lu->perm_c);
207:     /* we cannot simply memcpy on 64 bit archs */
208:     for (i = 0; i < m; i++) lu->perm_r[i] = ra[i];
209:     for (i = 0; i < n; i++) lu->perm_c[i] = ca[i];
210:     ISRestoreIndices(r,&ra);
211:     ISRestoreIndices(c,&ca);
212:   }

214:   /* symbolic factorization of A' */
215:   /* ---------------------------------------------------------------------- */
216:   if (lu->PetscMatOrdering) { /* use Petsc ordering */
217:     lu->Symbolic = klu_K_analyze_given(n,ai,aj,lu->perm_c,lu->perm_r,&lu->Common);
218:   } else { /* use klu internal ordering */
219:     lu->Symbolic = klu_K_analyze(n,ai,aj,&lu->Common);
220:   }
221:   if (!lu->Symbolic) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Symbolic Factorization failed");

223:   lu->flg                   = DIFFERENT_NONZERO_PATTERN;
224:   lu->CleanUpKLU            = PETSC_TRUE;
225:   (F)->ops->lufactornumeric = MatLUFactorNumeric_KLU;
226:   return(0);
227: }

231: static PetscErrorCode MatFactorInfo_KLU(Mat A,PetscViewer viewer)
232: {
233:   Mat_KLU       *lu= (Mat_KLU*)A->spptr;
234:   klu_K_numeric *Numeric=(klu_K_numeric*)lu->Numeric;

238:   /* check if matrix is KLU type */
239:   if (A->ops->solve != MatSolve_KLU) return(0);

241:   PetscViewerASCIIPrintf(viewer,"KLU stats:\n");
242:   PetscViewerASCIIPrintf(viewer,"  Number of diagonal blocks: %d\n",Numeric->nblocks);
243:   PetscViewerASCIIPrintf(viewer,"  Total nonzeros=%d\n",Numeric->lnz+Numeric->unz);

245:   PetscViewerASCIIPrintf(viewer,"KLU runtime parameters:\n");

247:   /* Control parameters used by numeric factorization */
248:   PetscViewerASCIIPrintf(viewer,"  Partial pivoting tolerance: %g\n",lu->Common.tol);
249:   /* BTF preordering */
250:   PetscViewerASCIIPrintf(viewer,"  BTF preordering enabled: %d\n",lu->Common.btf);
251:   /* mat ordering */
252:   if (!lu->PetscMatOrdering) {
253:     PetscViewerASCIIPrintf(viewer,"  Ordering: %s (not using the PETSc ordering)\n",KluOrderingTypes[(int)lu->Common.ordering]);
254:   } else {
255:     PetscViewerASCIIPrintf(viewer,"  Using PETSc ordering\n");
256:   }
257:   /* matrix row scaling */
258:   PetscViewerASCIIPrintf(viewer, "  Matrix row scaling: %s\n",scale[(int)lu->Common.scale]);
259:   return(0);
260: }

264: static PetscErrorCode MatView_KLU(Mat A,PetscViewer viewer)
265: {
266:   PetscErrorCode    ierr;
267:   PetscBool         iascii;
268:   PetscViewerFormat format;

271:   MatView_SeqAIJ(A,viewer);

273:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
274:   if (iascii) {
275:     PetscViewerGetFormat(viewer,&format);
276:     if (format == PETSC_VIEWER_ASCII_INFO) {
277:       MatFactorInfo_KLU(A,viewer);
278:     }
279:   }
280:   return(0);
281: }

285: PetscErrorCode MatFactorGetSolverPackage_seqaij_klu(Mat A,const MatSolverPackage *type)
286: {
288:   *type = MATSOLVERKLU;
289:   return(0);
290: }


293: /*MC
294:   MATSOLVERKLU = "klu" - A matrix type providing direct solvers (LU) for sequential matrices
295:   via the external package KLU.

297:   ./configure --download-suitesparse to install PETSc to use KLU

299:   Use -pc_type lu -pc_factor_mat_solver_package klu to us this direct solver

301:   Consult KLU documentation for more information on the options database keys below.

303:   Options Database Keys:
304: + -mat_klu_pivot_tol <0.001>                  - Partial pivoting tolerance
305: . -mat_klu_use_btf <1>                        - Use BTF preordering
306: . -mat_klu_ordering <AMD>                     - KLU reordering scheme to reduce fill-in (choose one of) AMD COLAMD PETSC
307: - -mat_klu_row_scale <NONE>                   - Matrix row scaling (choose one of) NONE SUM MAX 

309:    Note: KLU is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html

311:    Level: beginner

313: .seealso: PCLU, MATSOLVERUMFPACK, MATSOLVERCHOLMOD, PCFactorSetMatSolverPackage(), MatSolverPackage
314: M*/

318: PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat A,MatFactorType ftype,Mat *F)
319: {
320:   Mat            B;
321:   Mat_KLU       *lu;
323:   PetscInt       m=A->rmap->n,n=A->cmap->n,idx,status;
324:   PetscBool      flg;

327:   /* Create the factorization matrix F */
328:   MatCreate(PetscObjectComm((PetscObject)A),&B);
329:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);
330:   MatSetType(B,((PetscObject)A)->type_name);
331:   MatSeqAIJSetPreallocation(B,0,NULL);
332:   PetscNewLog(B,&lu);

334:   B->spptr                 = lu;
335:   B->ops->lufactorsymbolic = MatLUFactorSymbolic_KLU;
336:   B->ops->destroy          = MatDestroy_KLU;
337:   B->ops->view             = MatView_KLU;

339:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_klu);

341:   B->factortype   = MAT_FACTOR_LU;
342:   B->assembled    = PETSC_TRUE;           /* required by -ksp_view */
343:   B->preallocated = PETSC_TRUE;

345:   PetscFree(B->solvertype);
346:   PetscStrallocpy(MATSOLVERKLU,&B->solvertype);

348:   /* initializations */
349:   /* ------------------------------------------------*/
350:   /* get the default control parameters */
351:   status = klu_K_defaults(&lu->Common);
352:   if(status <= 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Initialization failed");

354:   lu->Common.scale = 0; /* No row scaling */

356:   PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"KLU Options","Mat");
357:   /* Partial pivoting tolerance */
358:   PetscOptionsReal("-mat_klu_pivot_tol","Partial pivoting tolerance","None",lu->Common.tol,&lu->Common.tol,NULL);
359:   /* BTF pre-ordering */
360:   PetscOptionsInt("-mat_klu_use_btf","Enable BTF preordering","None",lu->Common.btf,&lu->Common.btf,NULL);
361:   /* Matrix reordering */
362:   PetscOptionsEList("-mat_klu_ordering","Internal ordering method","None",KluOrderingTypes,sizeof(KluOrderingTypes)/sizeof(KluOrderingTypes[0]),KluOrderingTypes[0],&idx,&flg);
363:   if (flg) {
364:     if ((int)idx == 2) lu->PetscMatOrdering = PETSC_TRUE;   /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Klu perm_c) */
365:     else lu->Common.ordering = (int)idx;
366:   }
367:   /* Matrix row scaling */
368:   PetscOptionsEList("-mat_klu_row_scale","Matrix row scaling","None",scale,3,scale[0],&idx,&flg);
369:   PetscOptionsEnd();
370:   *F = B;
371:   return(0);
372: }