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