1 2 #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/ 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "PCFactorSetUpMatSolverPackage" 6 /*@ 7 PCFactorSetUpMatSolverPackage - Can be called after KSPSetOperators() or PCSetOperators(), causes MatGetFactor() to be called so then one may 8 set the options for that particular factorization object. 9 10 Input Parameter: 11 . pc - the preconditioner context 12 13 Notes: After you have called this function (which has to be after the KSPSetOperators() or PCSetOperators()) you can call PCFactorGetMatrix() and then set factor options on that matrix. 14 15 .seealso: PCFactorSetMatSolverPackage(), PCFactorGetMatrix() 16 17 Level: intermediate 18 19 @*/ 20 PetscErrorCode PCFactorSetUpMatSolverPackage(PC pc) 21 { 22 PetscErrorCode ierr; 23 24 PetscFunctionBegin; 25 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 26 ierr = PetscTryMethod(pc,"PCFactorSetUpMatSolverPackage_C",(PC),(pc));CHKERRQ(ierr); 27 PetscFunctionReturn(0); 28 } 29 30 #undef __FUNCT__ 31 #define __FUNCT__ "PCFactorSetZeroPivot" 32 /*@ 33 PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 34 35 Logically Collective on PC 36 37 Input Parameters: 38 + pc - the preconditioner context 39 - zero - all pivots smaller than this will be considered zero 40 41 Options Database Key: 42 . -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot 43 44 Level: intermediate 45 46 .keywords: PC, set, factorization, direct, fill 47 48 .seealso: PCFactorSetShiftType(), PCFactorSetShiftAmount() 49 @*/ 50 PetscErrorCode PCFactorSetZeroPivot(PC pc,PetscReal zero) 51 { 52 PetscErrorCode ierr; 53 54 PetscFunctionBegin; 55 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 56 PetscValidLogicalCollectiveReal(pc,zero,2); 57 ierr = PetscTryMethod(pc,"PCFactorSetZeroPivot_C",(PC,PetscReal),(pc,zero));CHKERRQ(ierr); 58 PetscFunctionReturn(0); 59 } 60 61 #undef __FUNCT__ 62 #define __FUNCT__ "PCFactorSetShiftType" 63 /*@ 64 PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during 65 numerical factorization, thus the matrix has nonzero pivots 66 67 Logically Collective on PC 68 69 Input Parameters: 70 + pc - the preconditioner context 71 - shifttype - type of shift; one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO, MAT_SHIFT_POSITIVE_DEFINITE, MAT_SHIFT_INBLOCKS 72 73 Options Database Key: 74 . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types 75 76 Level: intermediate 77 78 .keywords: PC, set, factorization, 79 80 .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftAmount() 81 @*/ 82 PetscErrorCode PCFactorSetShiftType(PC pc,MatFactorShiftType shifttype) 83 { 84 PetscErrorCode ierr; 85 86 PetscFunctionBegin; 87 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 88 PetscValidLogicalCollectiveEnum(pc,shifttype,2); 89 ierr = PetscTryMethod(pc,"PCFactorSetShiftType_C",(PC,MatFactorShiftType),(pc,shifttype));CHKERRQ(ierr); 90 PetscFunctionReturn(0); 91 } 92 93 #undef __FUNCT__ 94 #define __FUNCT__ "PCFactorSetShiftAmount" 95 /*@ 96 PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during 97 numerical factorization, thus the matrix has nonzero pivots 98 99 Logically Collective on PC 100 101 Input Parameters: 102 + pc - the preconditioner context 103 - shiftamount - amount of shift 104 105 Options Database Key: 106 . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 107 108 Level: intermediate 109 110 .keywords: PC, set, factorization, 111 112 .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftType() 113 @*/ 114 PetscErrorCode PCFactorSetShiftAmount(PC pc,PetscReal shiftamount) 115 { 116 PetscErrorCode ierr; 117 118 PetscFunctionBegin; 119 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 120 PetscValidLogicalCollectiveReal(pc,shiftamount,2); 121 ierr = PetscTryMethod(pc,"PCFactorSetShiftAmount_C",(PC,PetscReal),(pc,shiftamount));CHKERRQ(ierr); 122 PetscFunctionReturn(0); 123 } 124 125 #undef __FUNCT__ 126 #define __FUNCT__ "PCFactorSetDropTolerance" 127 /* 128 PCFactorSetDropTolerance - The preconditioner will use an ILU 129 based on a drop tolerance. (Under development) 130 131 Logically Collective on PC 132 133 Input Parameters: 134 + pc - the preconditioner context 135 . dt - the drop tolerance, try from 1.e-10 to .1 136 . dtcol - tolerance for column pivot, good values [0.1 to 0.01] 137 - maxrowcount - the max number of nonzeros allowed in a row, best value 138 depends on the number of nonzeros in row of original matrix 139 140 Options Database Key: 141 . -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance 142 143 Level: intermediate 144 145 There are NO default values for the 3 parameters, you must set them with reasonable values for your 146 matrix. We don't know how to compute reasonable values. 147 148 .keywords: PC, levels, reordering, factorization, incomplete, ILU 149 */ 150 PetscErrorCode PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount) 151 { 152 PetscErrorCode ierr; 153 154 PetscFunctionBegin; 155 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 156 PetscValidLogicalCollectiveReal(pc,dtcol,2); 157 PetscValidLogicalCollectiveInt(pc,maxrowcount,3); 158 ierr = PetscTryMethod(pc,"PCFactorSetDropTolerance_C",(PC,PetscReal,PetscReal,PetscInt),(pc,dt,dtcol,maxrowcount));CHKERRQ(ierr); 159 PetscFunctionReturn(0); 160 } 161 162 #undef __FUNCT__ 163 #define __FUNCT__ "PCFactorGetLevels" 164 /*@ 165 PCFactorGetLevels - Gets the number of levels of fill to use. 166 167 Logically Collective on PC 168 169 Input Parameters: 170 . pc - the preconditioner context 171 172 Output Parameter: 173 . levels - number of levels of fill 174 175 Level: intermediate 176 177 .keywords: PC, levels, fill, factorization, incomplete, ILU 178 @*/ 179 PetscErrorCode PCFactorGetLevels(PC pc,PetscInt *levels) 180 { 181 PetscErrorCode ierr; 182 183 PetscFunctionBegin; 184 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 185 ierr = PetscUseMethod(pc,"PCFactorGetLevels_C",(PC,PetscInt*),(pc,levels));CHKERRQ(ierr); 186 PetscFunctionReturn(0); 187 } 188 189 #undef __FUNCT__ 190 #define __FUNCT__ "PCFactorSetLevels" 191 /*@ 192 PCFactorSetLevels - Sets the number of levels of fill to use. 193 194 Logically Collective on PC 195 196 Input Parameters: 197 + pc - the preconditioner context 198 - levels - number of levels of fill 199 200 Options Database Key: 201 . -pc_factor_levels <levels> - Sets fill level 202 203 Level: intermediate 204 205 .keywords: PC, levels, fill, factorization, incomplete, ILU 206 @*/ 207 PetscErrorCode PCFactorSetLevels(PC pc,PetscInt levels) 208 { 209 PetscErrorCode ierr; 210 211 PetscFunctionBegin; 212 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 213 if (levels < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 214 PetscValidLogicalCollectiveInt(pc,levels,2); 215 ierr = PetscTryMethod(pc,"PCFactorSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr); 216 PetscFunctionReturn(0); 217 } 218 219 #undef __FUNCT__ 220 #define __FUNCT__ "PCFactorSetAllowDiagonalFill" 221 /*@ 222 PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 223 treated as level 0 fill even if there is no non-zero location. 224 225 Logically Collective on PC 226 227 Input Parameters: 228 + pc - the preconditioner context 229 - flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off 230 231 Options Database Key: 232 . -pc_factor_diagonal_fill 233 234 Notes: 235 Does not apply with 0 fill. 236 237 Level: intermediate 238 239 .keywords: PC, levels, fill, factorization, incomplete, ILU 240 241 .seealso: PCFactorGetAllowDiagonalFill() 242 243 @*/ 244 PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc,PetscBool flg) 245 { 246 PetscErrorCode ierr; 247 248 PetscFunctionBegin; 249 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 250 ierr = PetscTryMethod(pc,"PCFactorSetAllowDiagonalFill_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 251 PetscFunctionReturn(0); 252 } 253 254 #undef __FUNCT__ 255 #define __FUNCT__ "PCFactorGetAllowDiagonalFill" 256 /*@ 257 PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are 258 treated as level 0 fill even if there is no non-zero location. 259 260 Logically Collective on PC 261 262 Input Parameter: 263 . pc - the preconditioner context 264 265 Output Parameter: 266 . flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off 267 268 Options Database Key: 269 . -pc_factor_diagonal_fill 270 271 Notes: 272 Does not apply with 0 fill. 273 274 Level: intermediate 275 276 .keywords: PC, levels, fill, factorization, incomplete, ILU 277 278 .seealso: PCFactorSetAllowDiagonalFill() 279 280 @*/ 281 PetscErrorCode PCFactorGetAllowDiagonalFill(PC pc,PetscBool *flg) 282 { 283 PetscErrorCode ierr; 284 285 PetscFunctionBegin; 286 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 287 ierr = PetscUseMethod(pc,"PCFactorGetAllowDiagonalFill_C",(PC,PetscBool*),(pc,flg));CHKERRQ(ierr); 288 PetscFunctionReturn(0); 289 } 290 291 #undef __FUNCT__ 292 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal" 293 /*@ 294 PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 295 296 Logically Collective on PC 297 298 Input Parameters: 299 + pc - the preconditioner context 300 - tol - diagonal entries smaller than this in absolute value are considered zero 301 302 Options Database Key: 303 . -pc_factor_nonzeros_along_diagonal <tol> 304 305 Level: intermediate 306 307 .keywords: PC, set, factorization, direct, fill 308 309 .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal() 310 @*/ 311 PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 312 { 313 PetscErrorCode ierr; 314 315 PetscFunctionBegin; 316 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 317 PetscValidLogicalCollectiveReal(pc,rtol,2); 318 ierr = PetscTryMethod(pc,"PCFactorReorderForNonzeroDiagonal_C",(PC,PetscReal),(pc,rtol));CHKERRQ(ierr); 319 PetscFunctionReturn(0); 320 } 321 322 #undef __FUNCT__ 323 #define __FUNCT__ "PCFactorSetMatSolverPackage" 324 /*@C 325 PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization 326 327 Logically Collective on PC 328 329 Input Parameters: 330 + pc - the preconditioner context 331 - stype - for example, superlu, superlu_dist 332 333 Options Database Key: 334 . -pc_factor_mat_solver_package <stype> - petsc, superlu, superlu_dist, mumps, cusparse 335 336 Level: intermediate 337 338 Note: 339 By default this will use the PETSc factorization if it exists 340 341 342 .keywords: PC, set, factorization, direct, fill 343 344 .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage() 345 346 @*/ 347 PetscErrorCode PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype) 348 { 349 PetscErrorCode ierr; 350 351 PetscFunctionBegin; 352 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 353 ierr = PetscTryMethod(pc,"PCFactorSetMatSolverPackage_C",(PC,const MatSolverPackage),(pc,stype));CHKERRQ(ierr); 354 PetscFunctionReturn(0); 355 } 356 357 #undef __FUNCT__ 358 #define __FUNCT__ "PCFactorGetMatSolverPackage" 359 /*@C 360 PCFactorGetMatSolverPackage - gets the software that is used to perform the factorization 361 362 Not Collective 363 364 Input Parameter: 365 . pc - the preconditioner context 366 367 Output Parameter: 368 . stype - for example, superlu, superlu_dist (NULL if the PC does not have a solver package) 369 370 Level: intermediate 371 372 373 .keywords: PC, set, factorization, direct, fill 374 375 .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage() 376 377 @*/ 378 PetscErrorCode PCFactorGetMatSolverPackage(PC pc,const MatSolverPackage *stype) 379 { 380 PetscErrorCode ierr,(*f)(PC,const MatSolverPackage*); 381 382 PetscFunctionBegin; 383 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 384 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",&f);CHKERRQ(ierr); 385 if (f) { 386 ierr = (*f)(pc,stype);CHKERRQ(ierr); 387 } else { 388 *stype = NULL; 389 } 390 PetscFunctionReturn(0); 391 } 392 393 #undef __FUNCT__ 394 #define __FUNCT__ "PCFactorSetFill" 395 /*@ 396 PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 397 fill = number nonzeros in factor/number nonzeros in original matrix. 398 399 Not Collective, each process can expect a different amount of fill 400 401 Input Parameters: 402 + pc - the preconditioner context 403 - fill - amount of expected fill 404 405 Options Database Key: 406 . -pc_factor_fill <fill> - Sets fill amount 407 408 Level: intermediate 409 410 Note: 411 For sparse matrix factorizations it is difficult to predict how much 412 fill to expect. By running with the option -info PETSc will print the 413 actual amount of fill used; allowing you to set the value accurately for 414 future runs. Default PETSc uses a value of 5.0 415 416 This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with PCFactorSetLevels() or -pc_factor_levels. 417 418 419 .keywords: PC, set, factorization, direct, fill 420 421 @*/ 422 PetscErrorCode PCFactorSetFill(PC pc,PetscReal fill) 423 { 424 PetscErrorCode ierr; 425 426 PetscFunctionBegin; 427 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 428 if (fill < 1.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0"); 429 ierr = PetscTryMethod(pc,"PCFactorSetFill_C",(PC,PetscReal),(pc,fill));CHKERRQ(ierr); 430 PetscFunctionReturn(0); 431 } 432 433 #undef __FUNCT__ 434 #define __FUNCT__ "PCFactorSetUseInPlace" 435 /*@ 436 PCFactorSetUseInPlace - Tells the system to do an in-place factorization. 437 For dense matrices, this enables the solution of much larger problems. 438 For sparse matrices the factorization cannot be done truly in-place 439 so this does not save memory during the factorization, but after the matrix 440 is factored, the original unfactored matrix is freed, thus recovering that 441 space. 442 443 Logically Collective on PC 444 445 Input Parameters: 446 + pc - the preconditioner context 447 - flg - PETSC_TRUE to enable, PETSC_FALSE to disable 448 449 Options Database Key: 450 . -pc_factor_in_place <true,false>- Activate/deactivate in-place factorization 451 452 Notes: 453 PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when 454 a different matrix is provided for the multiply and the preconditioner in 455 a call to KSPSetOperators(). 456 This is because the Krylov space methods require an application of the 457 matrix multiplication, which is not possible here because the matrix has 458 been factored in-place, replacing the original matrix. 459 460 Level: intermediate 461 462 .keywords: PC, set, factorization, direct, inplace, in-place, LU 463 464 .seealso: PCFactorGetUseInPlace() 465 @*/ 466 PetscErrorCode PCFactorSetUseInPlace(PC pc,PetscBool flg) 467 { 468 PetscErrorCode ierr; 469 470 PetscFunctionBegin; 471 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 472 ierr = PetscTryMethod(pc,"PCFactorSetUseInPlace_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 473 PetscFunctionReturn(0); 474 } 475 476 #undef __FUNCT__ 477 #define __FUNCT__ "PCFactorGetUseInPlace" 478 /*@ 479 PCFactorGetUseInPlace - Determines if an in-place factorization is being used. 480 481 Logically Collective on PC 482 483 Input Parameter: 484 . pc - the preconditioner context 485 486 Output Parameter: 487 . flg - PETSC_TRUE to enable, PETSC_FALSE to disable 488 489 Level: intermediate 490 491 .keywords: PC, set, factorization, direct, inplace, in-place, LU 492 493 .seealso: PCFactorSetUseInPlace() 494 @*/ 495 PetscErrorCode PCFactorGetUseInPlace(PC pc,PetscBool *flg) 496 { 497 PetscErrorCode ierr; 498 499 PetscFunctionBegin; 500 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 501 ierr = PetscUseMethod(pc,"PCFactorGetUseInPlace_C",(PC,PetscBool*),(pc,flg));CHKERRQ(ierr); 502 PetscFunctionReturn(0); 503 } 504 505 #undef __FUNCT__ 506 #define __FUNCT__ "PCFactorSetMatOrderingType" 507 /*@C 508 PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 509 be used in the LU factorization. 510 511 Logically Collective on PC 512 513 Input Parameters: 514 + pc - the preconditioner context 515 - ordering - the matrix ordering name, for example, MATORDERINGND or MATORDERINGRCM 516 517 Options Database Key: 518 . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 519 520 Level: intermediate 521 522 Notes: nested dissection is used by default 523 524 For Cholesky and ICC and the SBAIJ format reorderings are not available, 525 since only the upper triangular part of the matrix is stored. You can use the 526 SeqAIJ format in this case to get reorderings. 527 528 @*/ 529 PetscErrorCode PCFactorSetMatOrderingType(PC pc,MatOrderingType ordering) 530 { 531 PetscErrorCode ierr; 532 533 PetscFunctionBegin; 534 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 535 ierr = PetscTryMethod(pc,"PCFactorSetMatOrderingType_C",(PC,MatOrderingType),(pc,ordering));CHKERRQ(ierr); 536 PetscFunctionReturn(0); 537 } 538 539 #undef __FUNCT__ 540 #define __FUNCT__ "PCFactorSetColumnPivot" 541 /*@ 542 PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization. 543 For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 544 it is never done. For the MATLAB and SuperLU factorization this is used. 545 546 Logically Collective on PC 547 548 Input Parameters: 549 + pc - the preconditioner context 550 - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 551 552 Options Database Key: 553 . -pc_factor_pivoting <dtcol> 554 555 Level: intermediate 556 557 .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks() 558 @*/ 559 PetscErrorCode PCFactorSetColumnPivot(PC pc,PetscReal dtcol) 560 { 561 PetscErrorCode ierr; 562 563 PetscFunctionBegin; 564 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 565 PetscValidLogicalCollectiveReal(pc,dtcol,2); 566 ierr = PetscTryMethod(pc,"PCFactorSetColumnPivot_C",(PC,PetscReal),(pc,dtcol));CHKERRQ(ierr); 567 PetscFunctionReturn(0); 568 } 569 570 #undef __FUNCT__ 571 #define __FUNCT__ "PCFactorSetPivotInBlocks" 572 /*@ 573 PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 574 with BAIJ or SBAIJ matrices 575 576 Logically Collective on PC 577 578 Input Parameters: 579 + pc - the preconditioner context 580 - pivot - PETSC_TRUE or PETSC_FALSE 581 582 Options Database Key: 583 . -pc_factor_pivot_in_blocks <true,false> 584 585 Level: intermediate 586 587 .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot() 588 @*/ 589 PetscErrorCode PCFactorSetPivotInBlocks(PC pc,PetscBool pivot) 590 { 591 PetscErrorCode ierr; 592 593 PetscFunctionBegin; 594 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 595 PetscValidLogicalCollectiveBool(pc,pivot,2); 596 ierr = PetscTryMethod(pc,"PCFactorSetPivotInBlocks_C",(PC,PetscBool),(pc,pivot));CHKERRQ(ierr); 597 PetscFunctionReturn(0); 598 } 599 600 #undef __FUNCT__ 601 #define __FUNCT__ "PCFactorSetReuseFill" 602 /*@ 603 PCFactorSetReuseFill - When matrices with same different nonzero structure are factored, 604 this causes later ones to use the fill ratio computed in the initial factorization. 605 606 Logically Collective on PC 607 608 Input Parameters: 609 + pc - the preconditioner context 610 - flag - PETSC_TRUE to reuse else PETSC_FALSE 611 612 Options Database Key: 613 . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 614 615 Level: intermediate 616 617 .keywords: PC, levels, reordering, factorization, incomplete, Cholesky 618 619 .seealso: PCFactorSetReuseOrdering() 620 @*/ 621 PetscErrorCode PCFactorSetReuseFill(PC pc,PetscBool flag) 622 { 623 PetscErrorCode ierr; 624 625 PetscFunctionBegin; 626 PetscValidHeaderSpecific(pc,PC_CLASSID,2); 627 PetscValidLogicalCollectiveBool(pc,flag,2); 628 ierr = PetscTryMethod(pc,"PCFactorSetReuseFill_C",(PC,PetscBool),(pc,flag));CHKERRQ(ierr); 629 PetscFunctionReturn(0); 630 } 631