1 /* 2 Defines a Cholesky factorization preconditioner for any Mat implementation. 3 Presently only provided for MPIRowbs format (i.e. BlockSolve). 4 */ 5 6 #include "src/ksp/pc/impls/factor/icc/icc.h" /*I "petscpc.h" I*/ 7 8 EXTERN_C_BEGIN 9 #undef __FUNCT__ 10 #define __FUNCT__ "PCICCSetMatOrdering_ICC" 11 PetscErrorCode PCICCSetMatOrdering_ICC(PC pc,MatOrderingType ordering) 12 { 13 PC_ICC *dir = (PC_ICC*)pc->data; 14 PetscErrorCode ierr; 15 16 PetscFunctionBegin; 17 ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 18 ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 19 PetscFunctionReturn(0); 20 } 21 EXTERN_C_END 22 23 EXTERN_C_BEGIN 24 #undef __FUNCT__ 25 #define __FUNCT__ "PCICCSetDamping_ICC" 26 PetscErrorCode PCICCSetDamping_ICC(PC pc,PetscReal damping) 27 { 28 PC_ICC *dir; 29 30 PetscFunctionBegin; 31 dir = (PC_ICC*)pc->data; 32 if (damping == (PetscReal) PETSC_DECIDE) { 33 dir->info.damping = 1.e-12; 34 } else { 35 dir->info.damping = damping; 36 } 37 PetscFunctionReturn(0); 38 } 39 EXTERN_C_END 40 41 EXTERN_C_BEGIN 42 #undef __FUNCT__ 43 #define __FUNCT__ "PCICCSetShift_ICC" 44 PetscErrorCode PCICCSetShift_ICC(PC pc,PetscTruth shift) 45 { 46 PC_ICC *dir; 47 48 PetscFunctionBegin; 49 dir = (PC_ICC*)pc->data; 50 dir->info.shift = shift; 51 if (shift) dir->info.shift_fraction = 0.0; 52 PetscFunctionReturn(0); 53 } 54 EXTERN_C_END 55 56 EXTERN_C_BEGIN 57 #undef __FUNCT__ 58 #define __FUNCT__ "PCICCSetSetZeroPivot_ICC" 59 PetscErrorCode PCICCSetZeroPivot_ICC(PC pc,PetscReal z) 60 { 61 PC_ICC *lu; 62 63 PetscFunctionBegin; 64 lu = (PC_ICC*)pc->data; 65 lu->info.zeropivot = z; 66 PetscFunctionReturn(0); 67 } 68 EXTERN_C_END 69 70 EXTERN_C_BEGIN 71 #undef __FUNCT__ 72 #define __FUNCT__ "PCICCSetFill_ICC" 73 PetscErrorCode PCICCSetFill_ICC(PC pc,PetscReal fill) 74 { 75 PC_ICC *dir; 76 77 PetscFunctionBegin; 78 dir = (PC_ICC*)pc->data; 79 dir->info.fill = fill; 80 PetscFunctionReturn(0); 81 } 82 EXTERN_C_END 83 84 EXTERN_C_BEGIN 85 #undef __FUNCT__ 86 #define __FUNCT__ "PCICCSetLevels_ICC" 87 PetscErrorCode PCICCSetLevels_ICC(PC pc,PetscInt levels) 88 { 89 PC_ICC *icc; 90 91 PetscFunctionBegin; 92 icc = (PC_ICC*)pc->data; 93 icc->info.levels = levels; 94 PetscFunctionReturn(0); 95 } 96 EXTERN_C_END 97 98 #undef __FUNCT__ 99 #define __FUNCT__ "PCICCSetMatOrdering" 100 /*@ 101 PCICCSetMatOrdering - Sets the ordering routine (to reduce fill) to 102 be used it the ICC factorization. 103 104 Collective on PC 105 106 Input Parameters: 107 + pc - the preconditioner context 108 - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 109 110 Options Database Key: 111 . -pc_icc_mat_ordering_type <nd,rcm,...> - Sets ordering routine 112 113 Level: intermediate 114 115 .seealso: PCLUSetMatOrdering() 116 117 .keywords: PC, ICC, set, matrix, reordering 118 119 @*/ 120 PetscErrorCode PCICCSetMatOrdering(PC pc,MatOrderingType ordering) 121 { 122 PetscErrorCode ierr,(*f)(PC,MatOrderingType); 123 124 PetscFunctionBegin; 125 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 126 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 127 if (f) { 128 ierr = (*f)(pc,ordering);CHKERRQ(ierr); 129 } 130 PetscFunctionReturn(0); 131 } 132 133 #undef __FUNCT__ 134 #define __FUNCT__ "PCICCSetLevels" 135 /*@ 136 PCICCSetLevels - Sets the number of levels of fill to use. 137 138 Collective on PC 139 140 Input Parameters: 141 + pc - the preconditioner context 142 - levels - number of levels of fill 143 144 Options Database Key: 145 . -pc_icc_levels <levels> - Sets fill level 146 147 Level: intermediate 148 149 Concepts: ICC^setting levels of fill 150 151 @*/ 152 PetscErrorCode PCICCSetLevels(PC pc,PetscInt levels) 153 { 154 PetscErrorCode ierr,(*f)(PC,PetscInt); 155 156 PetscFunctionBegin; 157 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 158 if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 159 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr); 160 if (f) { 161 ierr = (*f)(pc,levels);CHKERRQ(ierr); 162 } 163 PetscFunctionReturn(0); 164 } 165 166 #undef __FUNCT__ 167 #define __FUNCT__ "PCICCSetFill" 168 /*@ 169 PCICCSetFill - Indicate the amount of fill you expect in the factored matrix, 170 where fill = number nonzeros in factor/number nonzeros in original matrix. 171 172 Collective on PC 173 174 Input Parameters: 175 + pc - the preconditioner context 176 - fill - amount of expected fill 177 178 Options Database Key: 179 $ -pc_icc_fill <fill> 180 181 Note: 182 For sparse matrix factorizations it is difficult to predict how much 183 fill to expect. By running with the option -log_info PETSc will print the 184 actual amount of fill used; allowing you to set the value accurately for 185 future runs. But default PETSc uses a value of 1.0 186 187 Level: intermediate 188 189 .keywords: PC, set, factorization, direct, fill 190 191 .seealso: PCLUSetFill() 192 @*/ 193 PetscErrorCode PCICCSetFill(PC pc,PetscReal fill) 194 { 195 PetscErrorCode ierr,(*f)(PC,PetscReal); 196 197 PetscFunctionBegin; 198 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 199 if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less than 1.0"); 200 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 201 if (f) { 202 ierr = (*f)(pc,fill);CHKERRQ(ierr); 203 } 204 PetscFunctionReturn(0); 205 } 206 207 #undef __FUNCT__ 208 #define __FUNCT__ "PCICCSetDamping" 209 /*@ 210 PCICCSetDamping - adds this quantity to the diagonal of the matrix during the 211 ICC numerical factorization 212 213 Collective on PC 214 215 Input Parameters: 216 + pc - the preconditioner context 217 - damping - amount of damping 218 219 Options Database Key: 220 . -pc_icc_damping <damping> - Sets damping amount or PETSC_DECIDE for the default 221 222 Note: If 0.0 is given, then no damping is used. If a diagonal element is classified as a zero 223 pivot, then the damping is doubled until this is alleviated. 224 225 Level: intermediate 226 227 .keywords: PC, set, factorization, direct, fill 228 229 .seealso: PCICCSetFill(), PCLUSetDamping() 230 @*/ 231 PetscErrorCode PCICCSetDamping(PC pc,PetscReal damping) 232 { 233 PetscErrorCode ierr,(*f)(PC,PetscReal); 234 235 PetscFunctionBegin; 236 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 237 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 238 if (f) { 239 ierr = (*f)(pc,damping);CHKERRQ(ierr); 240 } 241 PetscFunctionReturn(0); 242 } 243 244 #undef __FUNCT__ 245 #define __FUNCT__ "PCICCSetShift" 246 /*@ 247 PCICCSetShift - specify whether to use Manteuffel shifting of ICC. 248 If an ICC factorisation breaks down because of nonpositive pivots, 249 adding sufficient identity to the diagonal will remedy this. 250 251 Manteuffel shifting for ICC uses a different algorithm than the ILU case. 252 Here we base the shift on the lack of diagonal dominance when a negative 253 pivot occurs. 254 255 Input parameters: 256 + pc - the preconditioner context 257 - shifting - PETSC_TRUE to set shift else PETSC_FALSE 258 259 Options Database Key: 260 . -pc_icc_shift - Activate PCICCSetShift() 261 262 Level: intermediate 263 264 .keywords: PC, indefinite, factorization, incomplete, ICC 265 266 .seealso: PCILUSetShift() 267 @*/ 268 PetscErrorCode PCICCSetShift(PC pc,PetscTruth shift) 269 { 270 PetscErrorCode ierr,(*f)(PC,PetscTruth); 271 272 PetscFunctionBegin; 273 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 274 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetShift_C",(void (**)(void))&f);CHKERRQ(ierr); 275 if (f) { 276 ierr = (*f)(pc,shift);CHKERRQ(ierr); 277 } 278 PetscFunctionReturn(0); 279 } 280 281 #undef __FUNCT__ 282 #define __FUNCT__ "PCICCSetZeroPivot" 283 /*@ 284 PCICCSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 285 286 Collective on PC 287 288 Input Parameters: 289 + pc - the preconditioner context 290 - zero - all pivots smaller than this will be considered zero 291 292 Options Database Key: 293 . -pc_ilu_zeropivot <zero> - Sets the zero pivot size 294 295 Level: intermediate 296 297 .keywords: PC, set, factorization, direct, fill 298 299 .seealso: PCICCSetFill(), PCLUSetDamp(), PCLUSetZeroPivot() 300 @*/ 301 PetscErrorCode PCICCSetZeroPivot(PC pc,PetscReal zero) 302 { 303 PetscErrorCode ierr,(*f)(PC,PetscReal); 304 305 PetscFunctionBegin; 306 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 307 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetZeroPivot_C",(void (**)(void))&f);CHKERRQ(ierr); 308 if (f) { 309 ierr = (*f)(pc,zero);CHKERRQ(ierr); 310 } 311 PetscFunctionReturn(0); 312 } 313 314 #undef __FUNCT__ 315 #define __FUNCT__ "PCSetup_ICC" 316 static PetscErrorCode PCSetup_ICC(PC pc) 317 { 318 PC_ICC *icc = (PC_ICC*)pc->data; 319 IS perm,cperm; 320 PetscErrorCode ierr; 321 322 PetscFunctionBegin; 323 ierr = MatGetOrdering(pc->pmat,icc->ordering,&perm,&cperm);CHKERRQ(ierr); 324 325 if (!pc->setupcalled) { 326 ierr = MatICCFactorSymbolic(pc->pmat,perm,&icc->info,&icc->fact);CHKERRQ(ierr); 327 } else if (pc->flag != SAME_NONZERO_PATTERN) { 328 ierr = MatDestroy(icc->fact);CHKERRQ(ierr); 329 ierr = MatICCFactorSymbolic(pc->pmat,perm,&icc->info,&icc->fact);CHKERRQ(ierr); 330 } 331 ierr = ISDestroy(cperm);CHKERRQ(ierr); 332 ierr = ISDestroy(perm);CHKERRQ(ierr); 333 ierr = MatCholeskyFactorNumeric(pc->pmat,&icc->info,&icc->fact);CHKERRQ(ierr); 334 PetscFunctionReturn(0); 335 } 336 337 #undef __FUNCT__ 338 #define __FUNCT__ "PCDestroy_ICC" 339 static PetscErrorCode PCDestroy_ICC(PC pc) 340 { 341 PC_ICC *icc = (PC_ICC*)pc->data; 342 PetscErrorCode ierr; 343 344 PetscFunctionBegin; 345 if (icc->fact) {ierr = MatDestroy(icc->fact);CHKERRQ(ierr);} 346 ierr = PetscStrfree(icc->ordering);CHKERRQ(ierr); 347 ierr = PetscFree(icc);CHKERRQ(ierr); 348 PetscFunctionReturn(0); 349 } 350 351 #undef __FUNCT__ 352 #define __FUNCT__ "PCApply_ICC" 353 static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y) 354 { 355 PC_ICC *icc = (PC_ICC*)pc->data; 356 PetscErrorCode ierr; 357 358 PetscFunctionBegin; 359 ierr = MatSolve(icc->fact,x,y);CHKERRQ(ierr); 360 PetscFunctionReturn(0); 361 } 362 363 #undef __FUNCT__ 364 #define __FUNCT__ "PCApplySymmetricLeft_ICC" 365 static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y) 366 { 367 PetscErrorCode ierr; 368 PC_ICC *icc = (PC_ICC*)pc->data; 369 370 PetscFunctionBegin; 371 ierr = MatForwardSolve(icc->fact,x,y);CHKERRQ(ierr); 372 PetscFunctionReturn(0); 373 } 374 375 #undef __FUNCT__ 376 #define __FUNCT__ "PCApplySymmetricRight_ICC" 377 static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y) 378 { 379 PetscErrorCode ierr; 380 PC_ICC *icc = (PC_ICC*)pc->data; 381 382 PetscFunctionBegin; 383 ierr = MatBackwardSolve(icc->fact,x,y);CHKERRQ(ierr); 384 PetscFunctionReturn(0); 385 } 386 387 #undef __FUNCT__ 388 #define __FUNCT__ "PCGetFactoredMatrix_ICC" 389 static PetscErrorCode PCGetFactoredMatrix_ICC(PC pc,Mat *mat) 390 { 391 PC_ICC *icc = (PC_ICC*)pc->data; 392 393 PetscFunctionBegin; 394 *mat = icc->fact; 395 PetscFunctionReturn(0); 396 } 397 398 #undef __FUNCT__ 399 #define __FUNCT__ "PCSetFromOptions_ICC" 400 static PetscErrorCode PCSetFromOptions_ICC(PC pc) 401 { 402 PC_ICC *icc = (PC_ICC*)pc->data; 403 char tname[256]; 404 PetscTruth flg; 405 PetscErrorCode ierr; 406 PetscFList ordlist; 407 408 PetscFunctionBegin; 409 ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 410 ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr); 411 ierr = PetscOptionsReal("-pc_icc_levels","levels of fill","PCICCSetLevels",icc->info.levels,&icc->info.levels,&flg);CHKERRQ(ierr); 412 ierr = PetscOptionsReal("-pc_icc_fill","Expected fill in factorization","PCICCSetFill",icc->info.fill,&icc->info.fill,&flg);CHKERRQ(ierr); 413 ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 414 ierr = PetscOptionsList("-pc_icc_mat_ordering_type","Reorder to reduce nonzeros in ICC","PCICCSetMatOrdering",ordlist,icc->ordering,tname,256,&flg);CHKERRQ(ierr); 415 if (flg) { 416 ierr = PCICCSetMatOrdering(pc,tname);CHKERRQ(ierr); 417 } 418 ierr = PetscOptionsName("-pc_icc_damping","Damping added to diagonal","PCICCSetDamping",&flg);CHKERRQ(ierr); 419 if (flg) { 420 ierr = PCICCSetDamping(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr); 421 } 422 ierr = PetscOptionsReal("-pc_icc_damping","Damping added to diagonal","PCICCSetDamping",icc->info.damping,&icc->info.damping,0);CHKERRQ(ierr); 423 ierr = PetscOptionsName("-pc_icc_shift","Manteuffel shift applied to diagonal","PCICCSetShift",&flg);CHKERRQ(ierr); 424 if (flg) { 425 ierr = PCICCSetShift(pc,PETSC_TRUE);CHKERRQ(ierr); 426 } else { 427 ierr = PCICCSetShift(pc,PETSC_FALSE);CHKERRQ(ierr); 428 } 429 ierr = PetscOptionsReal("-pc_icc_zeropivot","Pivot is considered zero if less than","PCICCSetZeroPivot",icc->info.zeropivot,&icc->info.zeropivot,0);CHKERRQ(ierr); 430 431 ierr = PetscOptionsTail();CHKERRQ(ierr); 432 PetscFunctionReturn(0); 433 } 434 435 #undef __FUNCT__ 436 #define __FUNCT__ "PCView_ICC" 437 static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer) 438 { 439 PC_ICC *icc = (PC_ICC*)pc->data; 440 PetscErrorCode ierr; 441 PetscTruth isstring,iascii; 442 443 PetscFunctionBegin; 444 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 445 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 446 if (iascii) { 447 if (icc->info.levels == 1) { 448 ierr = PetscViewerASCIIPrintf(viewer," ICC: %D level of fill\n",(PetscInt)icc->info.levels);CHKERRQ(ierr); 449 } else { 450 ierr = PetscViewerASCIIPrintf(viewer," ICC: %D levels of fill\n",(PetscInt)icc->info.levels);CHKERRQ(ierr); 451 } 452 ierr = PetscViewerASCIIPrintf(viewer," ICC: max fill ratio allocated %g\n",icc->info.fill);CHKERRQ(ierr); 453 if (icc->info.shift) {ierr = PetscViewerASCIIPrintf(viewer," ICC: using Manteuffel shift\n");CHKERRQ(ierr);} 454 } else if (isstring) { 455 ierr = PetscViewerStringSPrintf(viewer," lvls=%D",(PetscInt)icc->info.levels);CHKERRQ(ierr);CHKERRQ(ierr); 456 } else { 457 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCICC",((PetscObject)viewer)->type_name); 458 } 459 PetscFunctionReturn(0); 460 } 461 462 /*MC 463 PCICC - Incomplete Cholesky factorization preconditioners. 464 465 Options Database Keys: 466 + -pc_icc_levels <k> - number of levels of fill for ICC(k) 467 . -pc_icc_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for 468 its factorization (overwrites original matrix) 469 . -pc_icc_damping - add damping to diagonal to prevent zero (or very small) pivots 470 . -pc_icc_shift - apply Manteuffel shift to diagonal to force positive definite preconditioner 471 . -pc_icc_zeropivot <tol> - set tolerance for what is considered a zero pivot 472 . -pc_icc_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 473 - -pc_icc_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 474 475 Level: beginner 476 477 Concepts: incomplete Cholesky factorization 478 479 Notes: Only implemented for some matrix formats. Not implemented in parallel 480 481 For BAIJ matrices this implements a point block ICC. 482 483 The Manteuffel shift is only implemented for matrices with block size 1 484 485 By default, the Manteuffel is applied (for matrices with block size 1). Call PCICCSetShift(pc,PETSC_FALSE); 486 to turn off the shift. 487 488 489 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 490 PCICCSetSetZeroPivot(), PCICCSetDamping(), PCICCSetShift(), 491 PCICCSetFill(), PCICCSetMatOrdering(), PCICCSetReuseOrdering(), 492 PCICCSetLevels() 493 494 M*/ 495 496 EXTERN_C_BEGIN 497 #undef __FUNCT__ 498 #define __FUNCT__ "PCCreate_ICC" 499 PetscErrorCode PCCreate_ICC(PC pc) 500 { 501 PetscErrorCode ierr; 502 PC_ICC *icc; 503 504 PetscFunctionBegin; 505 ierr = PetscNew(PC_ICC,&icc);CHKERRQ(ierr); 506 PetscLogObjectMemory(pc,sizeof(PC_ICC)); 507 508 icc->fact = 0; 509 ierr = PetscStrallocpy(MATORDERING_NATURAL,&icc->ordering);CHKERRQ(ierr); 510 ierr = MatFactorInfoInitialize(&icc->info);CHKERRQ(ierr); 511 icc->info.levels = 0; 512 icc->info.fill = 1.0; 513 icc->implctx = 0; 514 515 icc->info.dtcol = PETSC_DEFAULT; 516 icc->info.damping = 0.0; 517 icc->info.shift = PETSC_TRUE; 518 icc->info.shift_fraction = 0.0; 519 icc->info.zeropivot = 1.e-12; 520 pc->data = (void*)icc; 521 522 pc->ops->apply = PCApply_ICC; 523 pc->ops->setup = PCSetup_ICC; 524 pc->ops->destroy = PCDestroy_ICC; 525 pc->ops->setfromoptions = PCSetFromOptions_ICC; 526 pc->ops->view = PCView_ICC; 527 pc->ops->getfactoredmatrix = PCGetFactoredMatrix_ICC; 528 pc->ops->applysymmetricleft = PCApplySymmetricLeft_ICC; 529 pc->ops->applysymmetricright = PCApplySymmetricRight_ICC; 530 531 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetLevels_C","PCICCSetLevels_ICC", 532 PCICCSetLevels_ICC);CHKERRQ(ierr); 533 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetFill_C","PCICCSetFill_ICC", 534 PCICCSetFill_ICC);CHKERRQ(ierr); 535 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetDamping_C","PCICCSetDamping_ICC", 536 PCICCSetDamping_ICC);CHKERRQ(ierr); 537 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetShift_C","PCICCSetShift_ICC", 538 PCICCSetShift_ICC);CHKERRQ(ierr); 539 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetMatOrdering_C","PCICCSetMatOrdering_ICC", 540 PCICCSetMatOrdering_ICC);CHKERRQ(ierr); 541 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetZeroPivot_C","PCICCSetZeroPivot_ICC", 542 PCICCSetZeroPivot_ICC);CHKERRQ(ierr); 543 PetscFunctionReturn(0); 544 } 545 EXTERN_C_END 546 547 548