1 /* 2 Defines a ILU factorization preconditioner for any Mat implementation 3 */ 4 #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 5 #include "src/ksp/pc/impls/factor/ilu/ilu.h" 6 #include "src/mat/matimpl.h" 7 8 /* ------------------------------------------------------------------------------------------*/ 9 10 EXTERN_C_BEGIN 11 #undef __FUNCT__ 12 #define __FUNCT__ "PCILUReorderForNonzeroDiagonal_ILU" 13 PetscErrorCode PCILUReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z) 14 { 15 PC_ILU *ilu = (PC_ILU*)pc->data; 16 17 PetscFunctionBegin; 18 ilu->nonzerosalongdiagonal = PETSC_TRUE; 19 if (z == PETSC_DECIDE) { 20 ilu->nonzerosalongdiagonaltol = 1.e-10; 21 } else { 22 ilu->nonzerosalongdiagonaltol = z; 23 } 24 PetscFunctionReturn(0); 25 } 26 EXTERN_C_END 27 28 #undef __FUNCT__ 29 #define __FUNCT__ "PCDestroy_ILU_Internal" 30 PetscErrorCode PCDestroy_ILU_Internal(PC pc) 31 { 32 PC_ILU *ilu = (PC_ILU*)pc->data; 33 PetscErrorCode ierr; 34 35 PetscFunctionBegin; 36 if (!ilu->inplace && ilu->fact) {ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);} 37 if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 38 if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 39 PetscFunctionReturn(0); 40 } 41 42 EXTERN_C_BEGIN 43 #undef __FUNCT__ 44 #define __FUNCT__ "PCILUSetUseDropTolerance_ILU" 45 PetscErrorCode PCILUSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 46 { 47 PC_ILU *ilu; 48 PetscErrorCode ierr; 49 50 PetscFunctionBegin; 51 ilu = (PC_ILU*)pc->data; 52 if (pc->setupcalled && (!ilu->usedt || ilu->info.dt != dt || ilu->info.dtcol != dtcol || ilu->info.dtcount != dtcount)) { 53 pc->setupcalled = 0; 54 ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 55 } 56 ilu->usedt = PETSC_TRUE; 57 ilu->info.dt = dt; 58 ilu->info.dtcol = dtcol; 59 ilu->info.dtcount = dtcount; 60 ilu->info.fill = PETSC_DEFAULT; 61 PetscFunctionReturn(0); 62 } 63 EXTERN_C_END 64 65 EXTERN_C_BEGIN 66 #undef __FUNCT__ 67 #define __FUNCT__ "PCILUSetFill_ILU" 68 PetscErrorCode PCILUSetFill_ILU(PC pc,PetscReal fill) 69 { 70 PC_ILU *dir; 71 72 PetscFunctionBegin; 73 dir = (PC_ILU*)pc->data; 74 dir->info.fill = fill; 75 PetscFunctionReturn(0); 76 } 77 EXTERN_C_END 78 79 EXTERN_C_BEGIN 80 #undef __FUNCT__ 81 #define __FUNCT__ "PCILUSetMatOrdering_ILU" 82 PetscErrorCode PCILUSetMatOrdering_ILU(PC pc,MatOrderingType ordering) 83 { 84 PC_ILU *dir = (PC_ILU*)pc->data; 85 PetscErrorCode ierr; 86 PetscTruth flg; 87 88 PetscFunctionBegin; 89 if (!pc->setupcalled) { 90 ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 91 ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 92 } else { 93 ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr); 94 if (!flg) { 95 pc->setupcalled = 0; 96 ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 97 ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 98 /* free the data structures, then create them again */ 99 ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 100 } 101 } 102 PetscFunctionReturn(0); 103 } 104 EXTERN_C_END 105 106 EXTERN_C_BEGIN 107 #undef __FUNCT__ 108 #define __FUNCT__ "PCILUSetReuseOrdering_ILU" 109 PetscErrorCode PCILUSetReuseOrdering_ILU(PC pc,PetscTruth flag) 110 { 111 PC_ILU *ilu; 112 113 PetscFunctionBegin; 114 ilu = (PC_ILU*)pc->data; 115 ilu->reuseordering = flag; 116 PetscFunctionReturn(0); 117 } 118 EXTERN_C_END 119 120 EXTERN_C_BEGIN 121 #undef __FUNCT__ 122 #define __FUNCT__ "PCILUDTSetReuseFill_ILUDT" 123 PetscErrorCode PCILUDTSetReuseFill_ILUDT(PC pc,PetscTruth flag) 124 { 125 PC_ILU *ilu; 126 127 PetscFunctionBegin; 128 ilu = (PC_ILU*)pc->data; 129 ilu->reusefill = flag; 130 if (flag) SETERRQ(PETSC_ERR_SUP,"Not yet supported"); 131 PetscFunctionReturn(0); 132 } 133 EXTERN_C_END 134 135 EXTERN_C_BEGIN 136 #undef __FUNCT__ 137 #define __FUNCT__ "PCILUSetLevels_ILU" 138 PetscErrorCode PCILUSetLevels_ILU(PC pc,PetscInt levels) 139 { 140 PC_ILU *ilu; 141 PetscErrorCode ierr; 142 143 PetscFunctionBegin; 144 ilu = (PC_ILU*)pc->data; 145 146 if (!pc->setupcalled) { 147 ilu->info.levels = levels; 148 } else if (ilu->usedt || ilu->info.levels != levels) { 149 ilu->info.levels = levels; 150 pc->setupcalled = 0; 151 ilu->usedt = PETSC_FALSE; 152 ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 153 } 154 PetscFunctionReturn(0); 155 } 156 EXTERN_C_END 157 158 EXTERN_C_BEGIN 159 #undef __FUNCT__ 160 #define __FUNCT__ "PCILUSetUseInPlace_ILU" 161 PetscErrorCode PCILUSetUseInPlace_ILU(PC pc) 162 { 163 PC_ILU *dir; 164 165 PetscFunctionBegin; 166 dir = (PC_ILU*)pc->data; 167 dir->inplace = PETSC_TRUE; 168 PetscFunctionReturn(0); 169 } 170 EXTERN_C_END 171 172 EXTERN_C_BEGIN 173 #undef __FUNCT__ 174 #define __FUNCT__ "PCILUSetAllowDiagonalFill" 175 PetscErrorCode PCILUSetAllowDiagonalFill_ILU(PC pc) 176 { 177 PC_ILU *dir; 178 179 PetscFunctionBegin; 180 dir = (PC_ILU*)pc->data; 181 dir->info.diagonal_fill = 1; 182 PetscFunctionReturn(0); 183 } 184 EXTERN_C_END 185 186 #undef __FUNCT__ 187 #define __FUNCT__ "PCILUSetUseDropTolerance" 188 /*@ 189 PCILUSetUseDropTolerance - The preconditioner will use an ILU 190 based on a drop tolerance. 191 192 Collective on PC 193 194 Input Parameters: 195 + pc - the preconditioner context 196 . dt - the drop tolerance, try from 1.e-10 to .1 197 . dtcol - tolerance for column pivot, good values [0.1 to 0.01] 198 - maxrowcount - the max number of nonzeros allowed in a row, best value 199 depends on the number of nonzeros in row of original matrix 200 201 Options Database Key: 202 . -pc_ilu_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance 203 204 Level: intermediate 205 206 Notes: 207 This uses the iludt() code of Saad's SPARSKIT package 208 209 .keywords: PC, levels, reordering, factorization, incomplete, ILU 210 @*/ 211 PetscErrorCode PCILUSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount) 212 { 213 PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt); 214 215 PetscFunctionBegin; 216 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 217 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 218 if (f) { 219 ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr); 220 } 221 PetscFunctionReturn(0); 222 } 223 224 #undef __FUNCT__ 225 #define __FUNCT__ "PCILUSetFill" 226 /*@ 227 PCILUSetFill - Indicate the amount of fill you expect in the factored matrix, 228 where fill = number nonzeros in factor/number nonzeros in original matrix. 229 230 Collective on PC 231 232 Input Parameters: 233 + pc - the preconditioner context 234 - fill - amount of expected fill 235 236 Options Database Key: 237 $ -pc_ilu_fill <fill> 238 239 Note: 240 For sparse matrix factorizations it is difficult to predict how much 241 fill to expect. By running with the option -log_info PETSc will print the 242 actual amount of fill used; allowing you to set the value accurately for 243 future runs. But default PETSc uses a value of 1.0 244 245 Level: intermediate 246 247 .keywords: PC, set, factorization, direct, fill 248 249 .seealso: PCLUSetFill() 250 @*/ 251 PetscErrorCode PCILUSetFill(PC pc,PetscReal fill) 252 { 253 PetscErrorCode ierr,(*f)(PC,PetscReal); 254 255 PetscFunctionBegin; 256 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 257 if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less than 1.0"); 258 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 259 if (f) { 260 ierr = (*f)(pc,fill);CHKERRQ(ierr); 261 } 262 PetscFunctionReturn(0); 263 } 264 265 #undef __FUNCT__ 266 #define __FUNCT__ "PCILUSetMatOrdering" 267 /*@C 268 PCILUSetMatOrdering - Sets the ordering routine (to reduce fill) to 269 be used in the ILU factorization. 270 271 Collective on PC 272 273 Input Parameters: 274 + pc - the preconditioner context 275 - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 276 277 Options Database Key: 278 . -pc_ilu_mat_ordering_type <nd,rcm,...> - Sets ordering routine 279 280 Level: intermediate 281 282 Notes: natural ordering is used by default 283 284 .seealso: PCLUSetMatOrdering() 285 286 .keywords: PC, ILU, set, matrix, reordering 287 288 @*/ 289 PetscErrorCode PCILUSetMatOrdering(PC pc,MatOrderingType ordering) 290 { 291 PetscErrorCode ierr,(*f)(PC,MatOrderingType); 292 293 PetscFunctionBegin; 294 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 295 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 296 if (f) { 297 ierr = (*f)(pc,ordering);CHKERRQ(ierr); 298 } 299 PetscFunctionReturn(0); 300 } 301 302 #undef __FUNCT__ 303 #define __FUNCT__ "PCILUSetReuseOrdering" 304 /*@ 305 PCILUSetReuseOrdering - When similar matrices are factored, this 306 causes the ordering computed in the first factor to be used for all 307 following factors; applies to both fill and drop tolerance ILUs. 308 309 Collective on PC 310 311 Input Parameters: 312 + pc - the preconditioner context 313 - flag - PETSC_TRUE to reuse else PETSC_FALSE 314 315 Options Database Key: 316 . -pc_ilu_reuse_ordering - Activate PCILUSetReuseOrdering() 317 318 Level: intermediate 319 320 .keywords: PC, levels, reordering, factorization, incomplete, ILU 321 322 .seealso: PCILUDTSetReuseFill(), PCLUSetReuseOrdering(), PCLUSetReuseFill() 323 @*/ 324 PetscErrorCode PCILUSetReuseOrdering(PC pc,PetscTruth flag) 325 { 326 PetscErrorCode ierr,(*f)(PC,PetscTruth); 327 328 PetscFunctionBegin; 329 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 330 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 331 if (f) { 332 ierr = (*f)(pc,flag);CHKERRQ(ierr); 333 } 334 PetscFunctionReturn(0); 335 } 336 337 #undef __FUNCT__ 338 #define __FUNCT__ "PCILUDTSetReuseFill" 339 /*@ 340 PCILUDTSetReuseFill - When matrices with same nonzero structure are ILUDT factored, 341 this causes later ones to use the fill computed in the initial factorization. 342 343 Collective on PC 344 345 Input Parameters: 346 + pc - the preconditioner context 347 - flag - PETSC_TRUE to reuse else PETSC_FALSE 348 349 Options Database Key: 350 . -pc_iludt_reuse_fill - Activates PCILUDTSetReuseFill() 351 352 Level: intermediate 353 354 .keywords: PC, levels, reordering, factorization, incomplete, ILU 355 356 .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCLUSetReuseFill() 357 @*/ 358 PetscErrorCode PCILUDTSetReuseFill(PC pc,PetscTruth flag) 359 { 360 PetscErrorCode ierr,(*f)(PC,PetscTruth); 361 362 PetscFunctionBegin; 363 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 364 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUDTSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr); 365 if (f) { 366 ierr = (*f)(pc,flag);CHKERRQ(ierr); 367 } 368 PetscFunctionReturn(0); 369 } 370 371 #undef __FUNCT__ 372 #define __FUNCT__ "PCILUSetLevels" 373 /*@ 374 PCILUSetLevels - Sets the number of levels of fill to use. 375 376 Collective on PC 377 378 Input Parameters: 379 + pc - the preconditioner context 380 - levels - number of levels of fill 381 382 Options Database Key: 383 . -pc_ilu_levels <levels> - Sets fill level 384 385 Level: intermediate 386 387 .keywords: PC, levels, fill, factorization, incomplete, ILU 388 @*/ 389 PetscErrorCode PCILUSetLevels(PC pc,PetscInt levels) 390 { 391 PetscErrorCode ierr,(*f)(PC,PetscInt); 392 393 PetscFunctionBegin; 394 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 395 if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 396 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr); 397 if (f) { 398 ierr = (*f)(pc,levels);CHKERRQ(ierr); 399 } 400 PetscFunctionReturn(0); 401 } 402 403 #undef __FUNCT__ 404 #define __FUNCT__ "PCILUSetAllowDiagonalFill" 405 /*@ 406 PCILUSetAllowDiagonalFill - Causes all diagonal matrix entries to be 407 treated as level 0 fill even if there is no non-zero location. 408 409 Collective on PC 410 411 Input Parameters: 412 + pc - the preconditioner context 413 414 Options Database Key: 415 . -pc_ilu_diagonal_fill 416 417 Notes: 418 Does not apply with 0 fill. 419 420 Level: intermediate 421 422 .keywords: PC, levels, fill, factorization, incomplete, ILU 423 @*/ 424 PetscErrorCode PCILUSetAllowDiagonalFill(PC pc) 425 { 426 PetscErrorCode ierr,(*f)(PC); 427 428 PetscFunctionBegin; 429 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 430 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr); 431 if (f) { 432 ierr = (*f)(pc);CHKERRQ(ierr); 433 } 434 PetscFunctionReturn(0); 435 } 436 437 #undef __FUNCT__ 438 #define __FUNCT__ "PCILUSetUseInPlace" 439 /*@ 440 PCILUSetUseInPlace - Tells the system to do an in-place incomplete factorization. 441 Collective on PC 442 443 Input Parameters: 444 . pc - the preconditioner context 445 446 Options Database Key: 447 . -pc_ilu_in_place - Activates in-place factorization 448 449 Notes: 450 PCILUSetUseInPlace() is intended for use with matrix-free variants of 451 Krylov methods, or when a different matrices are employed for the linear 452 system and preconditioner, or with ASM preconditioning. Do NOT use 453 this option if the linear system 454 matrix also serves as the preconditioning matrix, since the factored 455 matrix would then overwrite the original matrix. 456 457 Only works well with ILU(0). 458 459 Level: intermediate 460 461 .keywords: PC, set, factorization, inplace, in-place, ILU 462 463 .seealso: PCLUSetUseInPlace() 464 @*/ 465 PetscErrorCode PCILUSetUseInPlace(PC pc) 466 { 467 PetscErrorCode ierr,(*f)(PC); 468 469 PetscFunctionBegin; 470 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 471 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr); 472 if (f) { 473 ierr = (*f)(pc);CHKERRQ(ierr); 474 } 475 PetscFunctionReturn(0); 476 } 477 478 #undef __FUNCT__ 479 #define __FUNCT__ "PCILUReorderForNonzeroDiagonal" 480 /*@ 481 PCILUReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 482 483 Collective on PC 484 485 Input Parameters: 486 + pc - the preconditioner context 487 - tol - diagonal entries smaller than this in absolute value are considered zero 488 489 Options Database Key: 490 . -pc_lu_nonzeros_along_diagonal 491 492 Level: intermediate 493 494 .keywords: PC, set, factorization, direct, fill 495 496 .seealso: PCILUSetFill(), MatReorderForNonzeroDiagonal() 497 @*/ 498 PetscErrorCode PCILUReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 499 { 500 PetscErrorCode ierr,(*f)(PC,PetscReal); 501 502 PetscFunctionBegin; 503 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 504 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr); 505 if (f) { 506 ierr = (*f)(pc,rtol);CHKERRQ(ierr); 507 } 508 PetscFunctionReturn(0); 509 } 510 511 #undef __FUNCT__ 512 #define __FUNCT__ "PCILUSetPivotInBlocks" 513 /*@ 514 PCILUSetPivotInBlocks - Determines if pivoting is done while factoring each block 515 with BAIJ or SBAIJ matrices 516 517 Collective on PC 518 519 Input Parameters: 520 + pc - the preconditioner context 521 - pivot - PETSC_TRUE or PETSC_FALSE 522 523 Options Database Key: 524 . -pc_ilu_pivot_in_blocks <true,false> 525 526 Level: intermediate 527 528 .seealso: PCIILUSetMatOrdering(), PCILUSetPivoting() 529 @*/ 530 PetscErrorCode PCILUSetPivotInBlocks(PC pc,PetscTruth pivot) 531 { 532 PetscErrorCode ierr,(*f)(PC,PetscTruth); 533 534 PetscFunctionBegin; 535 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 536 if (f) { 537 ierr = (*f)(pc,pivot);CHKERRQ(ierr); 538 } 539 PetscFunctionReturn(0); 540 } 541 542 /* ------------------------------------------------------------------------------------------*/ 543 544 EXTERN_C_BEGIN 545 #undef __FUNCT__ 546 #define __FUNCT__ "PCILUSetPivotInBlocks_ILU" 547 PetscErrorCode PCILUSetPivotInBlocks_ILU(PC pc,PetscTruth pivot) 548 { 549 PC_ILU *dir = (PC_ILU*)pc->data; 550 551 PetscFunctionBegin; 552 dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 553 PetscFunctionReturn(0); 554 } 555 EXTERN_C_END 556 557 #undef __FUNCT__ 558 #define __FUNCT__ "PCSetFromOptions_ILU" 559 static PetscErrorCode PCSetFromOptions_ILU(PC pc) 560 { 561 PetscErrorCode ierr; 562 PetscInt dtmax = 3,itmp; 563 PetscTruth flg,set; 564 PetscReal dt[3]; 565 char tname[256]; 566 PC_ILU *ilu = (PC_ILU*)pc->data; 567 PetscFList ordlist; 568 PetscReal tol; 569 570 PetscFunctionBegin; 571 ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 572 ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr); 573 ierr = PetscOptionsInt("-pc_ilu_levels","levels of fill","PCILUSetLevels",(PetscInt)ilu->info.levels,&itmp,&flg);CHKERRQ(ierr); 574 if (flg) ilu->info.levels = itmp; 575 ierr = PetscOptionsName("-pc_ilu_in_place","do factorization in place","PCILUSetUseInPlace",&ilu->inplace);CHKERRQ(ierr); 576 ierr = PetscOptionsName("-pc_ilu_diagonal_fill","Allow fill into empty diagonal entry","PCILUSetAllowDiagonalFill",&flg);CHKERRQ(ierr); 577 ilu->info.diagonal_fill = (double) flg; 578 ierr = PetscOptionsName("-pc_iludt_reuse_fill","Reuse fill from previous ILUdt","PCILUDTSetReuseFill",&ilu->reusefill);CHKERRQ(ierr); 579 ierr = PetscOptionsName("-pc_ilu_reuse_ordering","Reuse previous reordering","PCILUSetReuseOrdering",&ilu->reuseordering);CHKERRQ(ierr); 580 581 ierr = PetscOptionsName("-pc_factor_shiftnonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 582 if (flg) { 583 ierr = PCFactorSetShiftNonzero((PetscReal) PETSC_DECIDE,&ilu->info);CHKERRQ(ierr); 584 } 585 ierr = PetscOptionsReal("-pc_factor_shiftnonzero","Shift added to diagonal","PCFactorSetShiftNonzero",ilu->info.shiftnz,&ilu->info.shiftnz,0);CHKERRQ(ierr); 586 587 ierr = PetscOptionsName("-pc_factor_shiftpd","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr); 588 if (flg) { 589 ierr = PetscOptionsInt("-pc_factor_shiftpd","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",(PetscInt)ilu->info.shiftpd,&itmp,&flg); CHKERRQ(ierr); 590 if (flg && !itmp) { 591 ierr = PCFactorSetShiftPd(PETSC_FALSE,&ilu->info);CHKERRQ(ierr); 592 } else { 593 ierr = PCFactorSetShiftPd(PETSC_TRUE,&ilu->info);CHKERRQ(ierr); 594 } 595 } 596 ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",ilu->info.zeropivot,&ilu->info.zeropivot,0);CHKERRQ(ierr); 597 598 dt[0] = ilu->info.dt; 599 dt[1] = ilu->info.dtcol; 600 dt[2] = ilu->info.dtcount; 601 ierr = PetscOptionsRealArray("-pc_ilu_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCILUSetUseDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr); 602 if (flg) { 603 ierr = PCILUSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr); 604 } 605 ierr = PetscOptionsReal("-pc_ilu_fill","Expected fill in factorization","PCILUSetFill",ilu->info.fill,&ilu->info.fill,&flg);CHKERRQ(ierr); 606 ierr = PetscOptionsName("-pc_ilu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCILUReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 607 if (flg) { 608 tol = PETSC_DECIDE; 609 ierr = PetscOptionsReal("-pc_ilu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCILUReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 610 ierr = PCILUReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 611 } 612 613 ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 614 ierr = PetscOptionsList("-pc_ilu_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCILUSetMatOrdering",ordlist,ilu->ordering,tname,256,&flg);CHKERRQ(ierr); 615 if (flg) { 616 ierr = PCILUSetMatOrdering(pc,tname);CHKERRQ(ierr); 617 } 618 flg = ilu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 619 ierr = PetscOptionsLogical("-pc_ilu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCILUSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 620 if (set) { 621 ierr = PCILUSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 622 } 623 ierr = PetscOptionsTail();CHKERRQ(ierr); 624 PetscFunctionReturn(0); 625 } 626 627 #undef __FUNCT__ 628 #define __FUNCT__ "PCView_ILU" 629 static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer) 630 { 631 PC_ILU *ilu = (PC_ILU*)pc->data; 632 PetscErrorCode ierr; 633 PetscTruth isstring,iascii; 634 635 PetscFunctionBegin; 636 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 637 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 638 if (iascii) { 639 if (ilu->usedt) { 640 ierr = PetscViewerASCIIPrintf(viewer," ILU: drop tolerance %g\n",ilu->info.dt);CHKERRQ(ierr); 641 ierr = PetscViewerASCIIPrintf(viewer," ILU: max nonzeros per row %D\n",(PetscInt)ilu->info.dtcount);CHKERRQ(ierr); 642 ierr = PetscViewerASCIIPrintf(viewer," ILU: column permutation tolerance %g\n",ilu->info.dtcol);CHKERRQ(ierr); 643 } else if (ilu->info.levels == 1) { 644 ierr = PetscViewerASCIIPrintf(viewer," ILU: %D level of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr); 645 } else { 646 ierr = PetscViewerASCIIPrintf(viewer," ILU: %D levels of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr); 647 } 648 ierr = PetscViewerASCIIPrintf(viewer," ILU: max fill ratio allocated %g\n",ilu->info.fill);CHKERRQ(ierr); 649 ierr = PetscViewerASCIIPrintf(viewer," ILU: tolerance for zero pivot %g\n",ilu->info.zeropivot);CHKERRQ(ierr); 650 if (ilu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," ILU: using Manteuffel shift\n");CHKERRQ(ierr);} 651 if (ilu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," in-place factorization\n");CHKERRQ(ierr);} 652 else {ierr = PetscViewerASCIIPrintf(viewer," out-of-place factorization\n");CHKERRQ(ierr);} 653 ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",ilu->ordering);CHKERRQ(ierr); 654 if (ilu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 655 if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 656 if (ilu->fact) { 657 ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 658 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 659 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 660 ierr = MatView(ilu->fact,viewer);CHKERRQ(ierr); 661 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 662 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 663 } 664 } else if (isstring) { 665 ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)ilu->info.levels,ilu->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 666 } else { 667 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name); 668 } 669 PetscFunctionReturn(0); 670 } 671 672 #undef __FUNCT__ 673 #define __FUNCT__ "PCSetUp_ILU" 674 static PetscErrorCode PCSetUp_ILU(PC pc) 675 { 676 PetscErrorCode ierr; 677 PC_ILU *ilu = (PC_ILU*)pc->data; 678 679 PetscFunctionBegin; 680 if (ilu->inplace) { 681 if (!pc->setupcalled) { 682 683 /* In-place factorization only makes sense with the natural ordering, 684 so we only need to get the ordering once, even if nonzero structure changes */ 685 ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 686 if (ilu->row) PetscLogObjectParent(pc,ilu->row); 687 if (ilu->col) PetscLogObjectParent(pc,ilu->col); 688 } 689 690 /* In place ILU only makes sense with fill factor of 1.0 because 691 cannot have levels of fill */ 692 ilu->info.fill = 1.0; 693 ilu->info.diagonal_fill = 0; 694 ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&ilu->info);CHKERRQ(ierr); 695 ilu->fact = pc->pmat; 696 } else if (ilu->usedt) { 697 if (!pc->setupcalled) { 698 ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 699 if (ilu->row) PetscLogObjectParent(pc,ilu->row); 700 if (ilu->col) PetscLogObjectParent(pc,ilu->col); 701 ierr = MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);CHKERRQ(ierr); 702 PetscLogObjectParent(pc,ilu->fact); 703 } else if (pc->flag != SAME_NONZERO_PATTERN) { 704 ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 705 if (!ilu->reuseordering) { 706 if (ilu->row) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 707 if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 708 ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 709 if (ilu->row) PetscLogObjectParent(pc,ilu->row); 710 if (ilu->col) PetscLogObjectParent(pc,ilu->col); 711 } 712 ierr = MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);CHKERRQ(ierr); 713 PetscLogObjectParent(pc,ilu->fact); 714 } else if (!ilu->reusefill) { 715 ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 716 ierr = MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);CHKERRQ(ierr); 717 PetscLogObjectParent(pc,ilu->fact); 718 } else { 719 ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr); 720 } 721 } else { 722 if (!pc->setupcalled) { 723 /* first time in so compute reordering and symbolic factorization */ 724 ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 725 if (ilu->row) PetscLogObjectParent(pc,ilu->row); 726 if (ilu->col) PetscLogObjectParent(pc,ilu->col); 727 /* Remove zeros along diagonal? */ 728 if (ilu->nonzerosalongdiagonal) { 729 ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 730 } 731 ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 732 PetscLogObjectParent(pc,ilu->fact); 733 } else if (pc->flag != SAME_NONZERO_PATTERN) { 734 if (!ilu->reuseordering) { 735 /* compute a new ordering for the ILU */ 736 ierr = ISDestroy(ilu->row);CHKERRQ(ierr); 737 ierr = ISDestroy(ilu->col);CHKERRQ(ierr); 738 ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 739 if (ilu->row) PetscLogObjectParent(pc,ilu->row); 740 if (ilu->col) PetscLogObjectParent(pc,ilu->col); 741 /* Remove zeros along diagonal? */ 742 if (ilu->nonzerosalongdiagonal) { 743 ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 744 } 745 } 746 ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 747 ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 748 PetscLogObjectParent(pc,ilu->fact); 749 } 750 ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr); 751 } 752 PetscFunctionReturn(0); 753 } 754 755 #undef __FUNCT__ 756 #define __FUNCT__ "PCDestroy_ILU" 757 static PetscErrorCode PCDestroy_ILU(PC pc) 758 { 759 PC_ILU *ilu = (PC_ILU*)pc->data; 760 PetscErrorCode ierr; 761 762 PetscFunctionBegin; 763 ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 764 ierr = PetscStrfree(ilu->ordering);CHKERRQ(ierr); 765 ierr = PetscFree(ilu);CHKERRQ(ierr); 766 PetscFunctionReturn(0); 767 } 768 769 #undef __FUNCT__ 770 #define __FUNCT__ "PCApply_ILU" 771 static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 772 { 773 PC_ILU *ilu = (PC_ILU*)pc->data; 774 PetscErrorCode ierr; 775 776 PetscFunctionBegin; 777 ierr = MatSolve(ilu->fact,x,y);CHKERRQ(ierr); 778 PetscFunctionReturn(0); 779 } 780 781 #undef __FUNCT__ 782 #define __FUNCT__ "PCApplyTranspose_ILU" 783 static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 784 { 785 PC_ILU *ilu = (PC_ILU*)pc->data; 786 PetscErrorCode ierr; 787 788 PetscFunctionBegin; 789 ierr = MatSolveTranspose(ilu->fact,x,y);CHKERRQ(ierr); 790 PetscFunctionReturn(0); 791 } 792 793 #undef __FUNCT__ 794 #define __FUNCT__ "PCGetFactoredMatrix_ILU" 795 static PetscErrorCode PCGetFactoredMatrix_ILU(PC pc,Mat *mat) 796 { 797 PC_ILU *ilu = (PC_ILU*)pc->data; 798 799 PetscFunctionBegin; 800 if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 801 *mat = ilu->fact; 802 PetscFunctionReturn(0); 803 } 804 805 /*MC 806 PCILU - Incomplete factorization preconditioners. 807 808 Options Database Keys: 809 + -pc_ilu_levels <k> - number of levels of fill for ILU(k) 810 . -pc_ilu_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 811 its factorization (overwrites original matrix) 812 . -pc_ilu_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 813 . -pc_ilu_reuse_ordering - reuse ordering of factorized matrix from previous factorization 814 . -pc_ilu_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt 815 . -pc_ilu_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 816 . -pc_ilu_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 817 this decreases the chance of getting a zero pivot 818 . -pc_ilu_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 819 - -pc_ilu_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 820 than 1 the diagonal blocks are factored with partial pivoting (this increases the 821 stability of the ILU factorization 822 823 Level: beginner 824 825 Concepts: incomplete factorization 826 827 Notes: Only implemented for some matrix formats. Not implemented in parallel 828 829 For BAIJ matrices this implements a point block ILU 830 831 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 832 PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), PCILUSetUseDropTolerance(), 833 PCILUSetFill(), PCILUSetMatOrdering(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill(), 834 PCILUSetLevels(), PCILUSetUseInPlace(), PCILUSetAllowDiagonalFill(), PCILUSetPivotInBlocks(), 835 836 M*/ 837 838 EXTERN_C_BEGIN 839 #undef __FUNCT__ 840 #define __FUNCT__ "PCCreate_ILU" 841 PetscErrorCode PCCreate_ILU(PC pc) 842 { 843 PetscErrorCode ierr; 844 PC_ILU *ilu; 845 846 PetscFunctionBegin; 847 ierr = PetscNew(PC_ILU,&ilu);CHKERRQ(ierr); 848 PetscLogObjectMemory(pc,sizeof(PC_ILU)); 849 850 ilu->fact = 0; 851 ierr = MatFactorInfoInitialize(&ilu->info);CHKERRQ(ierr); 852 ilu->info.levels = 0; 853 ilu->info.fill = 1.0; 854 ilu->col = 0; 855 ilu->row = 0; 856 ilu->inplace = PETSC_FALSE; 857 ierr = PetscStrallocpy(MATORDERING_NATURAL,&ilu->ordering);CHKERRQ(ierr); 858 ilu->reuseordering = PETSC_FALSE; 859 ilu->usedt = PETSC_FALSE; 860 ilu->info.dt = PETSC_DEFAULT; 861 ilu->info.dtcount = PETSC_DEFAULT; 862 ilu->info.dtcol = PETSC_DEFAULT; 863 ilu->info.shiftnz = 0.0; 864 ilu->info.shiftpd = PETSC_FALSE; 865 ilu->info.shift_fraction = 0.0; 866 ilu->info.zeropivot = 1.e-12; 867 ilu->info.pivotinblocks = 1.0; 868 ilu->reusefill = PETSC_FALSE; 869 ilu->info.diagonal_fill = 0; 870 pc->data = (void*)ilu; 871 872 pc->ops->destroy = PCDestroy_ILU; 873 pc->ops->apply = PCApply_ILU; 874 pc->ops->applytranspose = PCApplyTranspose_ILU; 875 pc->ops->setup = PCSetUp_ILU; 876 pc->ops->setfromoptions = PCSetFromOptions_ILU; 877 pc->ops->getfactoredmatrix = PCGetFactoredMatrix_ILU; 878 pc->ops->view = PCView_ILU; 879 pc->ops->applyrichardson = 0; 880 881 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseDropTolerance_C","PCILUSetUseDropTolerance_ILU", 882 PCILUSetUseDropTolerance_ILU);CHKERRQ(ierr); 883 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetFill_C","PCILUSetFill_ILU", 884 PCILUSetFill_ILU);CHKERRQ(ierr); 885 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetMatOrdering_C","PCILUSetMatOrdering_ILU", 886 PCILUSetMatOrdering_ILU);CHKERRQ(ierr); 887 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetReuseOrdering_C","PCILUSetReuseOrdering_ILU", 888 PCILUSetReuseOrdering_ILU);CHKERRQ(ierr); 889 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUDTSetReuseFill_C","PCILUDTSetReuseFill_ILUDT", 890 PCILUDTSetReuseFill_ILUDT);CHKERRQ(ierr); 891 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetLevels_C","PCILUSetLevels_ILU", 892 PCILUSetLevels_ILU);CHKERRQ(ierr); 893 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseInPlace_C","PCILUSetUseInPlace_ILU", 894 PCILUSetUseInPlace_ILU);CHKERRQ(ierr); 895 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetAllowDiagonalFill_C","PCILUSetAllowDiagonalFill_ILU", 896 PCILUSetAllowDiagonalFill_ILU);CHKERRQ(ierr); 897 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetPivotInBlocks_C","PCILUSetPivotInBlocks_ILU", 898 PCILUSetPivotInBlocks_ILU);CHKERRQ(ierr); 899 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUReorderForNonzeroDiagonal_C","PCILUReorderForNonzeroDiagonal_ILU", 900 PCILUReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr); 901 PetscFunctionReturn(0); 902 } 903 EXTERN_C_END 904