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 CHKMEMQ; 215 Mat fact=((PC_Factor*)ilu)->fact; 216 if (!fact){ 217 ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 218 } 219 CHKMEMQ; 220 ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 221 ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 222 ilu->actualfill = info.fill_ratio_needed; 223 ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 224 } else if (pc->flag != SAME_NONZERO_PATTERN) { 225 if (!ilu->reuseordering) { 226 /* compute a new ordering for the ILU */ 227 ierr = ISDestroy(ilu->row);CHKERRQ(ierr); 228 ierr = ISDestroy(ilu->col);CHKERRQ(ierr); 229 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 230 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 231 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 232 /* Remove zeros along diagonal? */ 233 if (ilu->nonzerosalongdiagonal) { 234 ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 235 } 236 } 237 ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 238 ierr = MatGetFactor(pc->pmat,MATSOLVERPETSC,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 239 ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 240 ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 241 ilu->actualfill = info.fill_ratio_needed; 242 ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 243 } 244 CHKMEMQ; 245 ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 246 CHKMEMQ; 247 } 248 PetscFunctionReturn(0); 249 } 250 251 #undef __FUNCT__ 252 #define __FUNCT__ "PCDestroy_ILU" 253 static PetscErrorCode PCDestroy_ILU(PC pc) 254 { 255 PC_ILU *ilu = (PC_ILU*)pc->data; 256 PetscErrorCode ierr; 257 258 PetscFunctionBegin; 259 ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 260 ierr = PetscFree(((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr); 261 ierr = PetscFree(((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 262 ierr = PetscFree(ilu);CHKERRQ(ierr); 263 PetscFunctionReturn(0); 264 } 265 266 #undef __FUNCT__ 267 #define __FUNCT__ "PCApply_ILU" 268 static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 269 { 270 PC_ILU *ilu = (PC_ILU*)pc->data; 271 PetscErrorCode ierr; 272 273 PetscFunctionBegin; 274 ierr = MatSolve(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr); 275 PetscFunctionReturn(0); 276 } 277 278 #undef __FUNCT__ 279 #define __FUNCT__ "PCApplyTranspose_ILU" 280 static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 281 { 282 PC_ILU *ilu = (PC_ILU*)pc->data; 283 PetscErrorCode ierr; 284 285 PetscFunctionBegin; 286 ierr = MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr); 287 PetscFunctionReturn(0); 288 } 289 290 /*MC 291 PCILU - Incomplete factorization preconditioners. 292 293 Options Database Keys: 294 + -pc_factor_levels <k> - number of levels of fill for ILU(k) 295 . -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 296 its factorization (overwrites original matrix) 297 . -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 298 . -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization 299 . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 300 . -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 301 this decreases the chance of getting a zero pivot 302 . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 303 . -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 304 than 1 the diagonal blocks are factored with partial pivoting (this increases the 305 stability of the ILU factorization 306 307 Level: beginner 308 309 Concepts: incomplete factorization 310 311 Notes: Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU) 312 313 For BAIJ matrices this implements a point block ILU 314 315 References: 316 T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving 317 self-adjoint elliptic difference equations. SIAM J. Numer. Anal., 5:559--573, 1968. 318 319 T.A. Oliphant. An implicit numerical method for solving two-dimensional time-dependent dif- 320 fusion problems. Quart. Appl. Math., 19:221--229, 1961. 321 322 Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST 323 http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical 324 Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 325 Science and Engineering, Kluwer, pp. 167--202. 326 327 328 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 329 PCFactorSetZeroPivot(), PCFactorSetShiftSetType(), PCFactorSetAmount(), 330 PCFactorSetDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 331 PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks() 332 333 M*/ 334 335 EXTERN_C_BEGIN 336 #undef __FUNCT__ 337 #define __FUNCT__ "PCCreate_ILU" 338 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ILU(PC pc) 339 { 340 PetscErrorCode ierr; 341 PC_ILU *ilu; 342 343 PetscFunctionBegin; 344 ierr = PetscNewLog(pc,PC_ILU,&ilu);CHKERRQ(ierr); 345 346 ((PC_Factor*)ilu)->fact = 0; 347 ierr = MatFactorInfoInitialize(&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 348 ((PC_Factor*)ilu)->factortype = MAT_FACTOR_ILU; 349 ((PC_Factor*)ilu)->info.levels = 0.; 350 ((PC_Factor*)ilu)->info.fill = 1.0; 351 ilu->col = 0; 352 ilu->row = 0; 353 ilu->inplace = PETSC_FALSE; 354 ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr); 355 ierr = PetscStrallocpy(MATORDERINGNATURAL,&((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 356 ilu->reuseordering = PETSC_FALSE; 357 ((PC_Factor*)ilu)->info.dt = PETSC_DEFAULT; 358 ((PC_Factor*)ilu)->info.dtcount = PETSC_DEFAULT; 359 ((PC_Factor*)ilu)->info.dtcol = PETSC_DEFAULT; 360 ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_NONZERO; 361 ((PC_Factor*)ilu)->info.shiftamount = 1.e-12; 362 ((PC_Factor*)ilu)->info.zeropivot = 1.e-12; 363 ((PC_Factor*)ilu)->info.pivotinblocks = 1.0; 364 ilu->reusefill = PETSC_FALSE; 365 ((PC_Factor*)ilu)->info.diagonal_fill = 0.; 366 pc->data = (void*)ilu; 367 368 pc->ops->destroy = PCDestroy_ILU; 369 pc->ops->apply = PCApply_ILU; 370 pc->ops->applytranspose = PCApplyTranspose_ILU; 371 pc->ops->setup = PCSetUp_ILU; 372 pc->ops->setfromoptions = PCSetFromOptions_ILU; 373 pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 374 pc->ops->view = PCView_ILU; 375 pc->ops->applyrichardson = 0; 376 377 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 378 PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 379 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor", 380 PCFactorSetShiftType_Factor);CHKERRQ(ierr); 381 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor", 382 PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 383 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 384 PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 385 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor", 386 PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 387 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetDropTolerance_C","PCFactorSetDropTolerance_ILU", 388 PCFactorSetDropTolerance_ILU);CHKERRQ(ierr); 389 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 390 PCFactorSetFill_Factor);CHKERRQ(ierr); 391 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 392 PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 393 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU", 394 PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr); 395 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU", 396 PCFactorSetReuseFill_ILU);CHKERRQ(ierr); 397 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor", 398 PCFactorSetLevels_Factor);CHKERRQ(ierr); 399 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU", 400 PCFactorSetUseInPlace_ILU);CHKERRQ(ierr); 401 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_Factor", 402 PCFactorSetAllowDiagonalFill_Factor);CHKERRQ(ierr); 403 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor", 404 PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr); 405 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU", 406 PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr); 407 PetscFunctionReturn(0); 408 } 409 EXTERN_C_END 410