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