1 2 /* 3 Provides an interface to the KLUv1.2 sparse solver 4 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. 9 10 */ 11 #include <../src/mat/impls/aij/seq/aij.h> 12 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 78 79 #define SuiteSparse_long long long 80 #define SuiteSparse_long_max LONG_LONG_MAX 81 #define SuiteSparse_long_id "%lld" 82 83 EXTERN_C_BEGIN 84 #include <klu.h> 85 EXTERN_C_END 86 87 static const char *KluOrderingTypes[] = {"AMD","COLAMD","PETSC"}; 88 static const char *scale[] ={"NONE","SUM","MAX"}; 89 90 typedef struct { 91 klu_K_common Common; 92 klu_K_symbolic *Symbolic; 93 klu_K_numeric *Numeric; 94 PetscInt *perm_c,*perm_r; 95 MatStructure flg; 96 PetscBool PetscMatOrdering; 97 PetscBool CleanUpKLU; 98 } Mat_KLU; 99 100 static PetscErrorCode MatDestroy_KLU(Mat A) 101 { 102 PetscErrorCode ierr; 103 Mat_KLU *lu=(Mat_KLU*)A->data; 104 105 PetscFunctionBegin; 106 if (lu->CleanUpKLU) { 107 klu_K_free_symbolic(&lu->Symbolic,&lu->Common); 108 klu_K_free_numeric(&lu->Numeric,&lu->Common); 109 ierr = PetscFree2(lu->perm_r,lu->perm_c);CHKERRQ(ierr); 110 } 111 ierr = PetscFree(A->data);CHKERRQ(ierr); 112 PetscFunctionReturn(0); 113 } 114 115 static PetscErrorCode MatSolveTranspose_KLU(Mat A,Vec b,Vec x) 116 { 117 Mat_KLU *lu = (Mat_KLU*)A->data; 118 PetscScalar *xa; 119 PetscErrorCode ierr; 120 PetscInt status; 121 122 PetscFunctionBegin; 123 /* KLU uses a column major format, solve Ax = b by klu_*_solve */ 124 /* ----------------------------------*/ 125 ierr = VecCopy(b,x); /* klu_solve stores the solution in rhs */ 126 ierr = VecGetArray(x,&xa); 127 status = klu_K_solve(lu->Symbolic,lu->Numeric,A->rmap->n,1,(PetscReal*)xa,&lu->Common); 128 if (status != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Solve failed"); 129 ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr); 130 PetscFunctionReturn(0); 131 } 132 133 static PetscErrorCode MatSolve_KLU(Mat A,Vec b,Vec x) 134 { 135 Mat_KLU *lu = (Mat_KLU*)A->data; 136 PetscScalar *xa; 137 PetscErrorCode ierr; 138 PetscInt status; 139 140 PetscFunctionBegin; 141 /* KLU uses a column major format, solve A^Tx = b by klu_*_tsolve */ 142 /* ----------------------------------*/ 143 ierr = VecCopy(b,x); /* klu_solve stores the solution in rhs */ 144 ierr = VecGetArray(x,&xa); 145 #if defined(PETSC_USE_COMPLEX) 146 PetscInt conj_solve=1; 147 status = klu_K_tsolve(lu->Symbolic,lu->Numeric,A->rmap->n,1,(PetscReal*)xa,conj_solve,&lu->Common); /* conjugate solve */ 148 #else 149 status = klu_K_tsolve(lu->Symbolic,lu->Numeric,A->rmap->n,1,xa,&lu->Common); 150 #endif 151 if (status != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Solve failed"); 152 ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr); 153 PetscFunctionReturn(0); 154 } 155 156 static PetscErrorCode MatLUFactorNumeric_KLU(Mat F,Mat A,const MatFactorInfo *info) 157 { 158 Mat_KLU *lu = (Mat_KLU*)(F)->data; 159 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 160 PetscInt *ai = a->i,*aj=a->j; 161 PetscScalar *av = a->a; 162 163 PetscFunctionBegin; 164 /* numeric factorization of A' */ 165 /* ----------------------------*/ 166 167 if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) { 168 klu_K_free_numeric(&lu->Numeric,&lu->Common); 169 } 170 lu->Numeric = klu_K_factor(ai,aj,(PetscReal*)av,lu->Symbolic,&lu->Common); 171 if(!lu->Numeric) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Numeric factorization failed"); 172 173 lu->flg = SAME_NONZERO_PATTERN; 174 lu->CleanUpKLU = PETSC_TRUE; 175 F->ops->solve = MatSolve_KLU; 176 F->ops->solvetranspose = MatSolveTranspose_KLU; 177 PetscFunctionReturn(0); 178 } 179 180 static PetscErrorCode MatLUFactorSymbolic_KLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 181 { 182 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 183 Mat_KLU *lu = (Mat_KLU*)(F->data); 184 PetscErrorCode ierr; 185 PetscInt i,*ai = a->i,*aj = a->j,m=A->rmap->n,n=A->cmap->n; 186 const PetscInt *ra,*ca; 187 188 PetscFunctionBegin; 189 if (lu->PetscMatOrdering) { 190 ierr = ISGetIndices(r,&ra);CHKERRQ(ierr); 191 ierr = ISGetIndices(c,&ca);CHKERRQ(ierr); 192 ierr = PetscMalloc2(m,&lu->perm_r,n,&lu->perm_c);CHKERRQ(ierr); 193 /* we cannot simply memcpy on 64 bit archs */ 194 for (i = 0; i < m; i++) lu->perm_r[i] = ra[i]; 195 for (i = 0; i < n; i++) lu->perm_c[i] = ca[i]; 196 ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr); 197 ierr = ISRestoreIndices(c,&ca);CHKERRQ(ierr); 198 } 199 200 /* symbolic factorization of A' */ 201 /* ---------------------------------------------------------------------- */ 202 if (lu->PetscMatOrdering) { /* use Petsc ordering */ 203 lu->Symbolic = klu_K_analyze_given(n,ai,aj,lu->perm_c,lu->perm_r,&lu->Common); 204 } else { /* use klu internal ordering */ 205 lu->Symbolic = klu_K_analyze(n,ai,aj,&lu->Common); 206 } 207 if (!lu->Symbolic) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Symbolic Factorization failed"); 208 209 lu->flg = DIFFERENT_NONZERO_PATTERN; 210 lu->CleanUpKLU = PETSC_TRUE; 211 (F)->ops->lufactornumeric = MatLUFactorNumeric_KLU; 212 PetscFunctionReturn(0); 213 } 214 215 static PetscErrorCode MatFactorInfo_KLU(Mat A,PetscViewer viewer) 216 { 217 Mat_KLU *lu= (Mat_KLU*)A->data; 218 klu_K_numeric *Numeric=(klu_K_numeric*)lu->Numeric; 219 PetscErrorCode ierr; 220 221 PetscFunctionBegin; 222 /* check if matrix is KLU type */ 223 if (A->ops->solve != MatSolve_KLU) PetscFunctionReturn(0); 224 225 ierr = PetscViewerASCIIPrintf(viewer,"KLU stats:\n");CHKERRQ(ierr); 226 ierr = PetscViewerASCIIPrintf(viewer," Number of diagonal blocks: %d\n",Numeric->nblocks); 227 ierr = PetscViewerASCIIPrintf(viewer," Total nonzeros=%d\n",Numeric->lnz+Numeric->unz);CHKERRQ(ierr); 228 229 ierr = PetscViewerASCIIPrintf(viewer,"KLU runtime parameters:\n");CHKERRQ(ierr); 230 231 /* Control parameters used by numeric factorization */ 232 ierr = PetscViewerASCIIPrintf(viewer," Partial pivoting tolerance: %g\n",lu->Common.tol);CHKERRQ(ierr); 233 /* BTF preordering */ 234 ierr = PetscViewerASCIIPrintf(viewer," BTF preordering enabled: %d\n",lu->Common.btf);CHKERRQ(ierr); 235 /* mat ordering */ 236 if (!lu->PetscMatOrdering) { 237 ierr = PetscViewerASCIIPrintf(viewer," Ordering: %s (not using the PETSc ordering)\n",KluOrderingTypes[(int)lu->Common.ordering]);CHKERRQ(ierr); 238 } else { 239 ierr = PetscViewerASCIIPrintf(viewer," Using PETSc ordering\n");CHKERRQ(ierr); 240 } 241 /* matrix row scaling */ 242 ierr = PetscViewerASCIIPrintf(viewer, " Matrix row scaling: %s\n",scale[(int)lu->Common.scale]);CHKERRQ(ierr); 243 PetscFunctionReturn(0); 244 } 245 246 static PetscErrorCode MatView_KLU(Mat A,PetscViewer viewer) 247 { 248 PetscErrorCode ierr; 249 PetscBool iascii; 250 PetscViewerFormat format; 251 252 PetscFunctionBegin; 253 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 254 if (iascii) { 255 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 256 if (format == PETSC_VIEWER_ASCII_INFO) { 257 ierr = MatFactorInfo_KLU(A,viewer);CHKERRQ(ierr); 258 } 259 } 260 PetscFunctionReturn(0); 261 } 262 263 PetscErrorCode MatFactorGetSolverPackage_seqaij_klu(Mat A,const MatSolverPackage *type) 264 { 265 PetscFunctionBegin; 266 *type = MATSOLVERKLU; 267 PetscFunctionReturn(0); 268 } 269 270 271 /*MC 272 MATSOLVERKLU = "klu" - A matrix type providing direct solvers (LU) for sequential matrices 273 via the external package KLU. 274 275 ./configure --download-suitesparse to install PETSc to use KLU 276 277 Use -pc_type lu -pc_factor_mat_solver_package klu to us this direct solver 278 279 Consult KLU documentation for more information on the options database keys below. 280 281 Options Database Keys: 282 + -mat_klu_pivot_tol <0.001> - Partial pivoting tolerance 283 . -mat_klu_use_btf <1> - Use BTF preordering 284 . -mat_klu_ordering <AMD> - KLU reordering scheme to reduce fill-in (choose one of) AMD COLAMD PETSC 285 - -mat_klu_row_scale <NONE> - Matrix row scaling (choose one of) NONE SUM MAX 286 287 Note: KLU is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html 288 289 Level: beginner 290 291 .seealso: PCLU, MATSOLVERUMFPACK, MATSOLVERCHOLMOD, PCFactorSetMatSolverPackage(), MatSolverPackage 292 M*/ 293 294 PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat A,MatFactorType ftype,Mat *F) 295 { 296 Mat B; 297 Mat_KLU *lu; 298 PetscErrorCode ierr; 299 PetscInt m=A->rmap->n,n=A->cmap->n,idx,status; 300 PetscBool flg; 301 302 PetscFunctionBegin; 303 /* Create the factorization matrix F */ 304 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 305 ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr); 306 ierr = PetscStrallocpy("klu",&((PetscObject)B)->type_name);CHKERRQ(ierr); 307 ierr = MatSetUp(B);CHKERRQ(ierr); 308 309 ierr = PetscNewLog(B,&lu);CHKERRQ(ierr); 310 311 B->data = lu; 312 B->ops->getinfo = MatGetInfo_External; 313 B->ops->lufactorsymbolic = MatLUFactorSymbolic_KLU; 314 B->ops->destroy = MatDestroy_KLU; 315 B->ops->view = MatView_KLU; 316 317 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_klu);CHKERRQ(ierr); 318 319 B->factortype = MAT_FACTOR_LU; 320 B->assembled = PETSC_TRUE; /* required by -ksp_view */ 321 B->preallocated = PETSC_TRUE; 322 323 ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 324 ierr = PetscStrallocpy(MATSOLVERKLU,&B->solvertype);CHKERRQ(ierr); 325 326 /* initializations */ 327 /* ------------------------------------------------*/ 328 /* get the default control parameters */ 329 status = klu_K_defaults(&lu->Common); 330 if(status <= 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Initialization failed"); 331 332 lu->Common.scale = 0; /* No row scaling */ 333 334 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"KLU Options","Mat");CHKERRQ(ierr); 335 /* Partial pivoting tolerance */ 336 ierr = PetscOptionsReal("-mat_klu_pivot_tol","Partial pivoting tolerance","None",lu->Common.tol,&lu->Common.tol,NULL);CHKERRQ(ierr); 337 /* BTF pre-ordering */ 338 ierr = PetscOptionsInt("-mat_klu_use_btf","Enable BTF preordering","None",lu->Common.btf,&lu->Common.btf,NULL);CHKERRQ(ierr); 339 /* Matrix reordering */ 340 ierr = PetscOptionsEList("-mat_klu_ordering","Internal ordering method","None",KluOrderingTypes,sizeof(KluOrderingTypes)/sizeof(KluOrderingTypes[0]),KluOrderingTypes[0],&idx,&flg);CHKERRQ(ierr); 341 if (flg) { 342 if ((int)idx == 2) lu->PetscMatOrdering = PETSC_TRUE; /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Klu perm_c) */ 343 else lu->Common.ordering = (int)idx; 344 } 345 /* Matrix row scaling */ 346 ierr = PetscOptionsEList("-mat_klu_row_scale","Matrix row scaling","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr); 347 PetscOptionsEnd(); 348 *F = B; 349 PetscFunctionReturn(0); 350 } 351