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