1 #define PETSCKSP_DLL 2 3 /* 4 Defines a ILU factorization preconditioner for any Mat implementation 5 */ 6 #include "../src/ksp/pc/impls/factor/ilu/ilu.h" /*I "petscpc.h" I*/ 7 8 /* ------------------------------------------------------------------------------------------*/ 9 EXTERN_C_BEGIN 10 #undef __FUNCT__ 11 #define __FUNCT__ "PCFactorSetReuseFill_ILU" 12 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_ILU(PC pc,PetscTruth flag) 13 { 14 PC_ILU *lu = (PC_ILU*)pc->data; 15 16 PetscFunctionBegin; 17 lu->reusefill = flag; 18 PetscFunctionReturn(0); 19 } 20 EXTERN_C_END 21 22 EXTERN_C_BEGIN 23 #undef __FUNCT__ 24 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_ILU" 25 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z) 26 { 27 PC_ILU *ilu = (PC_ILU*)pc->data; 28 29 PetscFunctionBegin; 30 ilu->nonzerosalongdiagonal = PETSC_TRUE; 31 if (z == PETSC_DECIDE) { 32 ilu->nonzerosalongdiagonaltol = 1.e-10; 33 } else { 34 ilu->nonzerosalongdiagonaltol = z; 35 } 36 PetscFunctionReturn(0); 37 } 38 EXTERN_C_END 39 40 #undef __FUNCT__ 41 #define __FUNCT__ "PCDestroy_ILU_Internal" 42 PetscErrorCode PCDestroy_ILU_Internal(PC pc) 43 { 44 PC_ILU *ilu = (PC_ILU*)pc->data; 45 PetscErrorCode ierr; 46 47 PetscFunctionBegin; 48 if (!ilu->inplace && ((PC_Factor*)ilu)->fact) {ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr);} 49 if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 50 if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 51 PetscFunctionReturn(0); 52 } 53 54 EXTERN_C_BEGIN 55 #undef __FUNCT__ 56 #define __FUNCT__ "PCFactorSetDropTolerance_ILU" 57 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 58 { 59 PC_ILU *ilu = (PC_ILU*)pc->data; 60 61 PetscFunctionBegin; 62 if (pc->setupcalled && (((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) { 63 SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cannot change drop tolerance after using PC"); 64 } 65 ((PC_Factor*)ilu)->info.dt = dt; 66 ((PC_Factor*)ilu)->info.dtcol = dtcol; 67 ((PC_Factor*)ilu)->info.dtcount = dtcount; 68 ((PC_Factor*)ilu)->info.usedt = 1.0; 69 PetscFunctionReturn(0); 70 } 71 EXTERN_C_END 72 73 EXTERN_C_BEGIN 74 #undef __FUNCT__ 75 #define __FUNCT__ "PCFactorSetReuseOrdering_ILU" 76 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_ILU(PC pc,PetscTruth flag) 77 { 78 PC_ILU *ilu = (PC_ILU*)pc->data; 79 80 PetscFunctionBegin; 81 ilu->reuseordering = flag; 82 PetscFunctionReturn(0); 83 } 84 EXTERN_C_END 85 86 EXTERN_C_BEGIN 87 #undef __FUNCT__ 88 #define __FUNCT__ "PCFactorSetUseInPlace_ILU" 89 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_ILU(PC pc) 90 { 91 PC_ILU *dir = (PC_ILU*)pc->data; 92 93 PetscFunctionBegin; 94 dir->inplace = PETSC_TRUE; 95 PetscFunctionReturn(0); 96 } 97 EXTERN_C_END 98 99 #undef __FUNCT__ 100 #define __FUNCT__ "PCSetFromOptions_ILU" 101 static PetscErrorCode PCSetFromOptions_ILU(PC pc) 102 { 103 PetscErrorCode ierr; 104 PetscInt itmp; 105 PetscTruth flg; 106 /* PetscReal dt[3]; */ 107 PC_ILU *ilu = (PC_ILU*)pc->data; 108 PetscReal tol; 109 110 PetscFunctionBegin; 111 ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr); 112 ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr); 113 114 ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);CHKERRQ(ierr); 115 if (flg) ((PC_Factor*)ilu)->info.levels = itmp; 116 flg = PETSC_FALSE; 117 ierr = PetscOptionsTruth("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 118 ((PC_Factor*)ilu)->info.diagonal_fill = (double) flg; 119 /* 120 dt[0] = ((PC_Factor*)ilu)->info.dt; 121 dt[1] = ((PC_Factor*)ilu)->info.dtcol; 122 dt[2] = ((PC_Factor*)ilu)->info.dtcount; 123 124 PetscInt dtmax = 3; 125 ierr = PetscOptionsRealArray("-pc_factor_drop_tolerance,","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr); 126 if (flg) { 127 ierr = PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr); 128 } 129 */ 130 ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 131 if (flg) { 132 tol = PETSC_DECIDE; 133 ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 134 ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 135 } 136 137 ierr = PetscOptionsTail();CHKERRQ(ierr); 138 PetscFunctionReturn(0); 139 } 140 141 #undef __FUNCT__ 142 #define __FUNCT__ "PCView_ILU" 143 static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer) 144 { 145 PC_ILU *ilu = (PC_ILU*)pc->data; 146 PetscErrorCode ierr; 147 PetscTruth iascii; 148 149 PetscFunctionBegin; 150 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 151 if (iascii) { 152 if (ilu->inplace) { 153 ierr = PetscViewerASCIIPrintf(viewer," ILU: in-place factorization\n");CHKERRQ(ierr); 154 } else { 155 ierr = PetscViewerASCIIPrintf(viewer," ILU: out-of-place factorization\n");CHKERRQ(ierr); 156 } 157 158 if (ilu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," ILU: Reusing fill from past factorization\n");CHKERRQ(ierr);} 159 if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," ILU: Reusing reordering from past factorization\n");CHKERRQ(ierr);} 160 } 161 ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 162 PetscFunctionReturn(0); 163 } 164 165 #undef __FUNCT__ 166 #define __FUNCT__ "PCSetUp_ILU" 167 static PetscErrorCode PCSetUp_ILU(PC pc) 168 { 169 PetscErrorCode ierr; 170 PC_ILU *ilu = (PC_ILU*)pc->data; 171 MatInfo info; 172 PetscTruth flg; 173 174 PetscFunctionBegin; 175 /* ugly hack to change default, since it is not support by some matrix types */ 176 if (((PC_Factor*)ilu)->info.shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 177 ierr = PetscTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg); 178 if (!flg) { 179 ierr = PetscTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&flg); 180 if (!flg) { 181 ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_INBLOCKS; 182 PetscInfo(pc,"Changing shift type from NONZERO to INBLOCKS because block matrices do not support NONZERO"); 183 } 184 } 185 } 186 187 if (ilu->inplace) { 188 CHKMEMQ; 189 if (!pc->setupcalled) { 190 191 /* In-place factorization only makes sense with the natural ordering, 192 so we only need to get the ordering once, even if nonzero structure changes */ 193 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 194 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 195 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 196 } 197 198 /* In place ILU only makes sense with fill factor of 1.0 because 199 cannot have levels of fill */ 200 ((PC_Factor*)ilu)->info.fill = 1.0; 201 ((PC_Factor*)ilu)->info.diagonal_fill = 0.0; 202 ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);CHKMEMQ; 203 ((PC_Factor*)ilu)->fact = pc->pmat; 204 } else { 205 if (!pc->setupcalled) { 206 /* first time in so compute reordering and symbolic factorization */ 207 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 208 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 209 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 210 /* Remove zeros along diagonal? */ 211 if (ilu->nonzerosalongdiagonal) { 212 ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 213 } 214 if (!((PC_Factor*)ilu)->fact){ 215 ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 216 } 217 ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 218 ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 219 ilu->actualfill = info.fill_ratio_needed; 220 ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 221 } else if (pc->flag != SAME_NONZERO_PATTERN) { 222 if (!ilu->reuseordering) { 223 /* compute a new ordering for the ILU */ 224 ierr = ISDestroy(ilu->row);CHKERRQ(ierr); 225 ierr = ISDestroy(ilu->col);CHKERRQ(ierr); 226 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 227 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 228 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 229 /* Remove zeros along diagonal? */ 230 if (ilu->nonzerosalongdiagonal) { 231 ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 232 } 233 } 234 ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 235 ierr = MatGetFactor(pc->pmat,MATSOLVERPETSC,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 236 ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 237 ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 238 ilu->actualfill = info.fill_ratio_needed; 239 ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 240 } 241 CHKMEMQ; 242 ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 243 CHKMEMQ; 244 } 245 PetscFunctionReturn(0); 246 } 247 248 #undef __FUNCT__ 249 #define __FUNCT__ "PCDestroy_ILU" 250 static PetscErrorCode PCDestroy_ILU(PC pc) 251 { 252 PC_ILU *ilu = (PC_ILU*)pc->data; 253 PetscErrorCode ierr; 254 255 PetscFunctionBegin; 256 ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 257 ierr = PetscFree(((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr); 258 ierr = PetscFree(((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 259 ierr = PetscFree(ilu);CHKERRQ(ierr); 260 PetscFunctionReturn(0); 261 } 262 263 #undef __FUNCT__ 264 #define __FUNCT__ "PCApply_ILU" 265 static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 266 { 267 PC_ILU *ilu = (PC_ILU*)pc->data; 268 PetscErrorCode ierr; 269 270 PetscFunctionBegin; 271 ierr = MatSolve(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr); 272 PetscFunctionReturn(0); 273 } 274 275 #undef __FUNCT__ 276 #define __FUNCT__ "PCApplyTranspose_ILU" 277 static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 278 { 279 PC_ILU *ilu = (PC_ILU*)pc->data; 280 PetscErrorCode ierr; 281 282 PetscFunctionBegin; 283 ierr = MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr); 284 PetscFunctionReturn(0); 285 } 286 287 /*MC 288 PCILU - Incomplete factorization preconditioners. 289 290 Options Database Keys: 291 + -pc_factor_levels <k> - number of levels of fill for ILU(k) 292 . -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 293 its factorization (overwrites original matrix) 294 . -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 295 . -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization 296 . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 297 . -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 298 this decreases the chance of getting a zero pivot 299 . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 300 . -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 301 than 1 the diagonal blocks are factored with partial pivoting (this increases the 302 stability of the ILU factorization 303 304 Level: beginner 305 306 Concepts: incomplete factorization 307 308 Notes: Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU) 309 310 For BAIJ matrices this implements a point block ILU 311 312 References: 313 T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving 314 self-adjoint elliptic difference equations. SIAM J. Numer. Anal., 5:559--573, 1968. 315 316 T.A. Oliphant. An implicit numerical method for solving two-dimensional time-dependent dif- 317 fusion problems. Quart. Appl. Math., 19:221--229, 1961. 318 319 Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST 320 http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical 321 Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 322 Science and Engineering, Kluwer, pp. 167--202. 323 324 325 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 326 PCFactorSetZeroPivot(), PCFactorSetShiftSetType(), PCFactorSetAmount(), 327 PCFactorSetDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 328 PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks() 329 330 M*/ 331 332 EXTERN_C_BEGIN 333 #undef __FUNCT__ 334 #define __FUNCT__ "PCCreate_ILU" 335 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ILU(PC pc) 336 { 337 PetscErrorCode ierr; 338 PC_ILU *ilu; 339 340 PetscFunctionBegin; 341 ierr = PetscNewLog(pc,PC_ILU,&ilu);CHKERRQ(ierr); 342 343 ((PC_Factor*)ilu)->fact = 0; 344 ierr = MatFactorInfoInitialize(&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 345 ((PC_Factor*)ilu)->factortype = MAT_FACTOR_ILU; 346 ((PC_Factor*)ilu)->info.levels = 0.; 347 ((PC_Factor*)ilu)->info.fill = 1.0; 348 ilu->col = 0; 349 ilu->row = 0; 350 ilu->inplace = PETSC_FALSE; 351 ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr); 352 ierr = PetscStrallocpy(MATORDERINGNATURAL,&((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 353 ilu->reuseordering = PETSC_FALSE; 354 ((PC_Factor*)ilu)->info.dt = PETSC_DEFAULT; 355 ((PC_Factor*)ilu)->info.dtcount = PETSC_DEFAULT; 356 ((PC_Factor*)ilu)->info.dtcol = PETSC_DEFAULT; 357 ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_NONZERO; 358 ((PC_Factor*)ilu)->info.shiftamount = 1.e-12; 359 ((PC_Factor*)ilu)->info.zeropivot = 1.e-12; 360 ((PC_Factor*)ilu)->info.pivotinblocks = 1.0; 361 ilu->reusefill = PETSC_FALSE; 362 ((PC_Factor*)ilu)->info.diagonal_fill = 0.; 363 pc->data = (void*)ilu; 364 365 pc->ops->destroy = PCDestroy_ILU; 366 pc->ops->apply = PCApply_ILU; 367 pc->ops->applytranspose = PCApplyTranspose_ILU; 368 pc->ops->setup = PCSetUp_ILU; 369 pc->ops->setfromoptions = PCSetFromOptions_ILU; 370 pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 371 pc->ops->view = PCView_ILU; 372 pc->ops->applyrichardson = 0; 373 374 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 375 PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 376 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor", 377 PCFactorSetShiftType_Factor);CHKERRQ(ierr); 378 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor", 379 PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 380 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 381 PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 382 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor", 383 PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 384 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetDropTolerance_C","PCFactorSetDropTolerance_ILU", 385 PCFactorSetDropTolerance_ILU);CHKERRQ(ierr); 386 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 387 PCFactorSetFill_Factor);CHKERRQ(ierr); 388 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 389 PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 390 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU", 391 PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr); 392 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU", 393 PCFactorSetReuseFill_ILU);CHKERRQ(ierr); 394 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor", 395 PCFactorSetLevels_Factor);CHKERRQ(ierr); 396 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU", 397 PCFactorSetUseInPlace_ILU);CHKERRQ(ierr); 398 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_Factor", 399 PCFactorSetAllowDiagonalFill_Factor);CHKERRQ(ierr); 400 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor", 401 PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr); 402 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU", 403 PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr); 404 PetscFunctionReturn(0); 405 } 406 EXTERN_C_END 407