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