1 #define PETSCKSP_DLL 2 3 /* 4 Defines a ILU factorization preconditioner for any Mat implementation 5 */ 6 #include "private/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__ "PCFactorSetReuseFill_ILU" 14 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_ILU(PC pc,PetscTruth flag) 15 { 16 PC_ILU *lu; 17 18 PetscFunctionBegin; 19 lu = (PC_ILU*)pc->data; 20 lu->reusefill = flag; 21 PetscFunctionReturn(0); 22 } 23 EXTERN_C_END 24 25 EXTERN_C_BEGIN 26 #undef __FUNCT__ 27 #define __FUNCT__ "PCFactorSetZeroPivot_ILU" 28 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_ILU(PC pc,PetscReal z) 29 { 30 PC_ILU *ilu; 31 32 PetscFunctionBegin; 33 ilu = (PC_ILU*)pc->data; 34 ilu->info.zeropivot = z; 35 PetscFunctionReturn(0); 36 } 37 EXTERN_C_END 38 39 EXTERN_C_BEGIN 40 #undef __FUNCT__ 41 #define __FUNCT__ "PCFactorSetShiftNonzero_ILU" 42 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_ILU(PC pc,PetscReal shift) 43 { 44 PC_ILU *dir; 45 46 PetscFunctionBegin; 47 dir = (PC_ILU*)pc->data; 48 if (shift == (PetscReal) PETSC_DECIDE) { 49 dir->info.shiftnz = 1.e-12; 50 } else { 51 dir->info.shiftnz = shift; 52 } 53 PetscFunctionReturn(0); 54 } 55 EXTERN_C_END 56 57 EXTERN_C_BEGIN 58 #undef __FUNCT__ 59 #define __FUNCT__ "PCFactorSetShiftPd_ILU" 60 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_ILU(PC pc,PetscTruth shift) 61 { 62 PC_ILU *dir; 63 64 PetscFunctionBegin; 65 dir = (PC_ILU*)pc->data; 66 if (shift) { 67 dir->info.shift_fraction = 0.0; 68 dir->info.shiftpd = 1.0; 69 } else { 70 dir->info.shiftpd = 0.0; 71 } 72 PetscFunctionReturn(0); 73 } 74 EXTERN_C_END 75 76 EXTERN_C_BEGIN 77 #undef __FUNCT__ 78 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_ILU" 79 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z) 80 { 81 PC_ILU *ilu = (PC_ILU*)pc->data; 82 83 PetscFunctionBegin; 84 ilu->nonzerosalongdiagonal = PETSC_TRUE; 85 if (z == PETSC_DECIDE) { 86 ilu->nonzerosalongdiagonaltol = 1.e-10; 87 } else { 88 ilu->nonzerosalongdiagonaltol = z; 89 } 90 PetscFunctionReturn(0); 91 } 92 EXTERN_C_END 93 94 #undef __FUNCT__ 95 #define __FUNCT__ "PCDestroy_ILU_Internal" 96 PetscErrorCode PCDestroy_ILU_Internal(PC pc) 97 { 98 PC_ILU *ilu = (PC_ILU*)pc->data; 99 PetscErrorCode ierr; 100 101 PetscFunctionBegin; 102 if (!ilu->inplace && ilu->fact) {ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);} 103 if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 104 if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 105 PetscFunctionReturn(0); 106 } 107 108 EXTERN_C_BEGIN 109 #undef __FUNCT__ 110 #define __FUNCT__ "PCFactorSetUseDropTolerance_ILU" 111 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 112 { 113 PC_ILU *ilu; 114 PetscErrorCode ierr; 115 116 PetscFunctionBegin; 117 ilu = (PC_ILU*)pc->data; 118 if (pc->setupcalled && (!ilu->usedt || ilu->info.dt != dt || ilu->info.dtcol != dtcol || ilu->info.dtcount != dtcount)) { 119 pc->setupcalled = 0; 120 ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 121 } 122 ilu->usedt = PETSC_TRUE; 123 ilu->info.dt = dt; 124 ilu->info.dtcol = dtcol; 125 ilu->info.dtcount = dtcount; 126 ilu->info.fill = PETSC_DEFAULT; 127 PetscFunctionReturn(0); 128 } 129 EXTERN_C_END 130 131 EXTERN_C_BEGIN 132 #undef __FUNCT__ 133 #define __FUNCT__ "PCFactorSetFill_ILU" 134 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_ILU(PC pc,PetscReal fill) 135 { 136 PC_ILU *dir; 137 138 PetscFunctionBegin; 139 dir = (PC_ILU*)pc->data; 140 dir->info.fill = fill; 141 PetscFunctionReturn(0); 142 } 143 EXTERN_C_END 144 145 EXTERN_C_BEGIN 146 #undef __FUNCT__ 147 #define __FUNCT__ "PCFactorSetMatOrdering_ILU" 148 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrdering_ILU(PC pc,MatOrderingType ordering) 149 { 150 PC_ILU *dir = (PC_ILU*)pc->data; 151 PetscErrorCode ierr; 152 PetscTruth flg; 153 154 PetscFunctionBegin; 155 if (!pc->setupcalled) { 156 ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 157 ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 158 } else { 159 ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr); 160 if (!flg) { 161 pc->setupcalled = 0; 162 ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 163 ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 164 /* free the data structures, then create them again */ 165 ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 166 } 167 } 168 PetscFunctionReturn(0); 169 } 170 EXTERN_C_END 171 172 EXTERN_C_BEGIN 173 #undef __FUNCT__ 174 #define __FUNCT__ "PCFactorSetReuseOrdering_ILU" 175 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_ILU(PC pc,PetscTruth flag) 176 { 177 PC_ILU *ilu; 178 179 PetscFunctionBegin; 180 ilu = (PC_ILU*)pc->data; 181 ilu->reuseordering = flag; 182 PetscFunctionReturn(0); 183 } 184 EXTERN_C_END 185 186 EXTERN_C_BEGIN 187 #undef __FUNCT__ 188 #define __FUNCT__ "PCFactorSetLevels_ILU" 189 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels_ILU(PC pc,PetscInt levels) 190 { 191 PC_ILU *ilu; 192 PetscErrorCode ierr; 193 194 PetscFunctionBegin; 195 ilu = (PC_ILU*)pc->data; 196 197 if (!pc->setupcalled) { 198 ilu->info.levels = levels; 199 } else if (ilu->usedt || ilu->info.levels != levels) { 200 ilu->info.levels = levels; 201 pc->setupcalled = 0; 202 ilu->usedt = PETSC_FALSE; 203 ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 204 } 205 PetscFunctionReturn(0); 206 } 207 EXTERN_C_END 208 209 EXTERN_C_BEGIN 210 #undef __FUNCT__ 211 #define __FUNCT__ "PCFactorSetUseInPlace_ILU" 212 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_ILU(PC pc) 213 { 214 PC_ILU *dir; 215 216 PetscFunctionBegin; 217 dir = (PC_ILU*)pc->data; 218 dir->inplace = PETSC_TRUE; 219 PetscFunctionReturn(0); 220 } 221 EXTERN_C_END 222 223 EXTERN_C_BEGIN 224 #undef __FUNCT__ 225 #define __FUNCT__ "PCFactorSetAllowDiagonalFill_ILU" 226 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill_ILU(PC pc) 227 { 228 PC_ILU *dir; 229 230 PetscFunctionBegin; 231 dir = (PC_ILU*)pc->data; 232 dir->info.diagonal_fill = 1; 233 PetscFunctionReturn(0); 234 } 235 EXTERN_C_END 236 237 #undef __FUNCT__ 238 #define __FUNCT__ "PCFactorSetUseDropTolerance" 239 /*@ 240 PCFactorSetUseDropTolerance - The preconditioner will use an ILU 241 based on a drop tolerance. 242 243 Collective on PC 244 245 Input Parameters: 246 + pc - the preconditioner context 247 . dt - the drop tolerance, try from 1.e-10 to .1 248 . dtcol - tolerance for column pivot, good values [0.1 to 0.01] 249 - maxrowcount - the max number of nonzeros allowed in a row, best value 250 depends on the number of nonzeros in row of original matrix 251 252 Options Database Key: 253 . -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance 254 255 Level: intermediate 256 257 Notes: 258 This uses the iludt() code of Saad's SPARSKIT package 259 260 .keywords: PC, levels, reordering, factorization, incomplete, ILU 261 @*/ 262 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount) 263 { 264 PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt); 265 266 PetscFunctionBegin; 267 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 268 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 269 if (f) { 270 ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr); 271 } 272 PetscFunctionReturn(0); 273 } 274 275 #undef __FUNCT__ 276 #define __FUNCT__ "PCFactorSetLevels" 277 /*@ 278 PCFactorSetLevels - Sets the number of levels of fill to use. 279 280 Collective on PC 281 282 Input Parameters: 283 + pc - the preconditioner context 284 - levels - number of levels of fill 285 286 Options Database Key: 287 . -pc_factor_levels <levels> - Sets fill level 288 289 Level: intermediate 290 291 .keywords: PC, levels, fill, factorization, incomplete, ILU 292 @*/ 293 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC pc,PetscInt levels) 294 { 295 PetscErrorCode ierr,(*f)(PC,PetscInt); 296 297 PetscFunctionBegin; 298 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 299 if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 300 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr); 301 if (f) { 302 ierr = (*f)(pc,levels);CHKERRQ(ierr); 303 } 304 PetscFunctionReturn(0); 305 } 306 307 #undef __FUNCT__ 308 #define __FUNCT__ "PCFactorSetAllowDiagonalFill" 309 /*@ 310 PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 311 treated as level 0 fill even if there is no non-zero location. 312 313 Collective on PC 314 315 Input Parameters: 316 + pc - the preconditioner context 317 318 Options Database Key: 319 . -pc_factor_diagonal_fill 320 321 Notes: 322 Does not apply with 0 fill. 323 324 Level: intermediate 325 326 .keywords: PC, levels, fill, factorization, incomplete, ILU 327 @*/ 328 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC pc) 329 { 330 PetscErrorCode ierr,(*f)(PC); 331 332 PetscFunctionBegin; 333 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 334 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr); 335 if (f) { 336 ierr = (*f)(pc);CHKERRQ(ierr); 337 } 338 PetscFunctionReturn(0); 339 } 340 341 /* ------------------------------------------------------------------------------------------*/ 342 343 EXTERN_C_BEGIN 344 #undef __FUNCT__ 345 #define __FUNCT__ "PCFactorSetPivotInBlocks_ILU" 346 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_ILU(PC pc,PetscTruth pivot) 347 { 348 PC_ILU *dir = (PC_ILU*)pc->data; 349 350 PetscFunctionBegin; 351 dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 352 PetscFunctionReturn(0); 353 } 354 EXTERN_C_END 355 356 #undef __FUNCT__ 357 #define __FUNCT__ "PCSetFromOptions_ILU" 358 static PetscErrorCode PCSetFromOptions_ILU(PC pc) 359 { 360 PetscErrorCode ierr; 361 PetscInt dtmax = 3,itmp; 362 PetscTruth flg,set; 363 PetscReal dt[3]; 364 char tname[256]; 365 PC_ILU *ilu = (PC_ILU*)pc->data; 366 PetscFList ordlist; 367 PetscReal tol; 368 369 PetscFunctionBegin; 370 ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 371 ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr); 372 ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)ilu->info.levels,&itmp,&flg);CHKERRQ(ierr); 373 if (flg) ilu->info.levels = itmp; 374 ierr = PetscOptionsName("-pc_factor_in_place","do factorization in place","PCFactorSetUseInPlace",&ilu->inplace);CHKERRQ(ierr); 375 ierr = PetscOptionsName("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",&flg);CHKERRQ(ierr); 376 ilu->info.diagonal_fill = (double) flg; 377 ierr = PetscOptionsName("-pc_factor_reuse_fill","Reuse fill ratio from previous factorization","PCFactorSetReuseFill",&ilu->reusefill);CHKERRQ(ierr); 378 ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse previous reordering","PCFactorSetReuseOrdering",&ilu->reuseordering);CHKERRQ(ierr); 379 380 ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 381 if (flg) { 382 ierr = PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr); 383 } 384 ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",ilu->info.shiftnz,&ilu->info.shiftnz,0);CHKERRQ(ierr); 385 386 ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr); 387 if (flg) { 388 ierr = PetscOptionsInt("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",(PetscInt)ilu->info.shiftpd,&itmp,&flg); CHKERRQ(ierr); 389 if (flg && !itmp) { 390 ierr = PCFactorSetShiftPd(pc,PETSC_FALSE);CHKERRQ(ierr); 391 } else { 392 ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr); 393 } 394 } 395 ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",ilu->info.zeropivot,&ilu->info.zeropivot,0);CHKERRQ(ierr); 396 397 dt[0] = ilu->info.dt; 398 dt[1] = ilu->info.dtcol; 399 dt[2] = ilu->info.dtcount; 400 ierr = PetscOptionsRealArray("-pc_factor_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetUseDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr); 401 if (flg) { 402 ierr = PCFactorSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr); 403 } 404 ierr = PetscOptionsReal("-pc_factor_fill","Expected fill in factorization","PCFactorSetFill",ilu->info.fill,&ilu->info.fill,&flg);CHKERRQ(ierr); 405 ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 406 if (flg) { 407 tol = PETSC_DECIDE; 408 ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 409 ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 410 } 411 412 ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 413 ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCFactorSetMatOrdering",ordlist,ilu->ordering,tname,256,&flg);CHKERRQ(ierr); 414 if (flg) { 415 ierr = PCFactorSetMatOrdering(pc,tname);CHKERRQ(ierr); 416 } 417 flg = ilu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 418 ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 419 if (set) { 420 ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 421 } 422 ierr = PetscOptionsTail();CHKERRQ(ierr); 423 PetscFunctionReturn(0); 424 } 425 426 #undef __FUNCT__ 427 #define __FUNCT__ "PCView_ILU" 428 static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer) 429 { 430 PC_ILU *ilu = (PC_ILU*)pc->data; 431 PetscErrorCode ierr; 432 PetscTruth isstring,iascii; 433 434 PetscFunctionBegin; 435 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 436 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 437 if (iascii) { 438 if (ilu->usedt) { 439 ierr = PetscViewerASCIIPrintf(viewer," ILU: drop tolerance %G\n",ilu->info.dt);CHKERRQ(ierr); 440 ierr = PetscViewerASCIIPrintf(viewer," ILU: max nonzeros per row %D\n",(PetscInt)ilu->info.dtcount);CHKERRQ(ierr); 441 ierr = PetscViewerASCIIPrintf(viewer," ILU: column permutation tolerance %G\n",ilu->info.dtcol);CHKERRQ(ierr); 442 } else if (ilu->info.levels == 1) { 443 ierr = PetscViewerASCIIPrintf(viewer," ILU: %D level of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr); 444 } else { 445 ierr = PetscViewerASCIIPrintf(viewer," ILU: %D levels of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr); 446 } 447 ierr = PetscViewerASCIIPrintf(viewer," ILU: max fill ratio allocated %G\n",ilu->info.fill);CHKERRQ(ierr); 448 ierr = PetscViewerASCIIPrintf(viewer," ILU: tolerance for zero pivot %G\n",ilu->info.zeropivot);CHKERRQ(ierr); 449 if (ilu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," ILU: using Manteuffel shift\n");CHKERRQ(ierr);} 450 if (ilu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," in-place factorization\n");CHKERRQ(ierr);} 451 else {ierr = PetscViewerASCIIPrintf(viewer," out-of-place factorization\n");CHKERRQ(ierr);} 452 ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",ilu->ordering);CHKERRQ(ierr); 453 if (ilu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 454 if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 455 if (ilu->fact) { 456 ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 457 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 458 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 459 ierr = MatView(ilu->fact,viewer);CHKERRQ(ierr); 460 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 461 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 462 } 463 } else if (isstring) { 464 ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)ilu->info.levels,ilu->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 465 } else { 466 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name); 467 } 468 PetscFunctionReturn(0); 469 } 470 471 #undef __FUNCT__ 472 #define __FUNCT__ "PCSetUp_ILU" 473 static PetscErrorCode PCSetUp_ILU(PC pc) 474 { 475 PetscErrorCode ierr; 476 PC_ILU *ilu = (PC_ILU*)pc->data; 477 478 PetscFunctionBegin; 479 if (ilu->inplace) { 480 if (!pc->setupcalled) { 481 482 /* In-place factorization only makes sense with the natural ordering, 483 so we only need to get the ordering once, even if nonzero structure changes */ 484 ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 485 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 486 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 487 } 488 489 /* In place ILU only makes sense with fill factor of 1.0 because 490 cannot have levels of fill */ 491 ilu->info.fill = 1.0; 492 ilu->info.diagonal_fill = 0; 493 ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&ilu->info);CHKERRQ(ierr); 494 ilu->fact = pc->pmat; 495 } else if (ilu->usedt) { 496 if (!pc->setupcalled) { 497 ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 498 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 499 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 500 ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 501 ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 502 } else if (pc->flag != SAME_NONZERO_PATTERN) { 503 ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 504 if (!ilu->reuseordering) { 505 if (ilu->row) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 506 if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 507 ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 508 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 509 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 510 } 511 ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 512 ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 513 } else if (!ilu->reusefill) { 514 ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 515 ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 516 ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 517 } else { 518 ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr); 519 } 520 } else { 521 if (!pc->setupcalled) { 522 /* first time in so compute reordering and symbolic factorization */ 523 ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 524 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 525 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 526 /* Remove zeros along diagonal? */ 527 if (ilu->nonzerosalongdiagonal) { 528 ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 529 } 530 ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 531 ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 532 } else if (pc->flag != SAME_NONZERO_PATTERN) { 533 if (!ilu->reuseordering) { 534 /* compute a new ordering for the ILU */ 535 ierr = ISDestroy(ilu->row);CHKERRQ(ierr); 536 ierr = ISDestroy(ilu->col);CHKERRQ(ierr); 537 ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 538 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 539 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 540 /* Remove zeros along diagonal? */ 541 if (ilu->nonzerosalongdiagonal) { 542 ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 543 } 544 } 545 ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 546 ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 547 ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 548 } 549 ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr); 550 } 551 PetscFunctionReturn(0); 552 } 553 554 #undef __FUNCT__ 555 #define __FUNCT__ "PCDestroy_ILU" 556 static PetscErrorCode PCDestroy_ILU(PC pc) 557 { 558 PC_ILU *ilu = (PC_ILU*)pc->data; 559 PetscErrorCode ierr; 560 561 PetscFunctionBegin; 562 ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 563 ierr = PetscStrfree(ilu->ordering);CHKERRQ(ierr); 564 ierr = PetscFree(ilu);CHKERRQ(ierr); 565 PetscFunctionReturn(0); 566 } 567 568 #undef __FUNCT__ 569 #define __FUNCT__ "PCApply_ILU" 570 static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 571 { 572 PC_ILU *ilu = (PC_ILU*)pc->data; 573 PetscErrorCode ierr; 574 575 PetscFunctionBegin; 576 ierr = MatSolve(ilu->fact,x,y);CHKERRQ(ierr); 577 PetscFunctionReturn(0); 578 } 579 580 #undef __FUNCT__ 581 #define __FUNCT__ "PCApplyTranspose_ILU" 582 static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 583 { 584 PC_ILU *ilu = (PC_ILU*)pc->data; 585 PetscErrorCode ierr; 586 587 PetscFunctionBegin; 588 ierr = MatSolveTranspose(ilu->fact,x,y);CHKERRQ(ierr); 589 PetscFunctionReturn(0); 590 } 591 592 #undef __FUNCT__ 593 #define __FUNCT__ "PCGetFactoredMatrix_ILU" 594 static PetscErrorCode PCGetFactoredMatrix_ILU(PC pc,Mat *mat) 595 { 596 PC_ILU *ilu = (PC_ILU*)pc->data; 597 598 PetscFunctionBegin; 599 if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 600 *mat = ilu->fact; 601 PetscFunctionReturn(0); 602 } 603 604 /*MC 605 PCILU - Incomplete factorization preconditioners. 606 607 Options Database Keys: 608 + -pc_factor_levels <k> - number of levels of fill for ILU(k) 609 . -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 610 its factorization (overwrites original matrix) 611 . -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 612 . -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization 613 . -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt 614 . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 615 . -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 616 this decreases the chance of getting a zero pivot 617 . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 618 . -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 619 than 1 the diagonal blocks are factored with partial pivoting (this increases the 620 stability of the ILU factorization 621 . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 622 - -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 623 is optional with PETSC_TRUE being the default 624 625 Level: beginner 626 627 Concepts: incomplete factorization 628 629 Notes: Only implemented for some matrix formats. Not implemented in parallel 630 631 For BAIJ matrices this implements a point block ILU 632 633 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 634 PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), PCFactorSetUseDropTolerance(), 635 PCFactorSetFill(), PCFactorSetMatOrdering(), PCFactorSetReuseOrdering(), 636 PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(), 637 PCFactorSetShiftNonzero(),PCFactorSetShiftPd() 638 639 M*/ 640 641 EXTERN_C_BEGIN 642 #undef __FUNCT__ 643 #define __FUNCT__ "PCCreate_ILU" 644 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ILU(PC pc) 645 { 646 PetscErrorCode ierr; 647 PC_ILU *ilu; 648 649 PetscFunctionBegin; 650 ierr = PetscNew(PC_ILU,&ilu);CHKERRQ(ierr); 651 ierr = PetscLogObjectMemory(pc,sizeof(PC_ILU));CHKERRQ(ierr); 652 653 ilu->fact = 0; 654 ierr = MatFactorInfoInitialize(&ilu->info);CHKERRQ(ierr); 655 ilu->info.levels = 0; 656 ilu->info.fill = 1.0; 657 ilu->col = 0; 658 ilu->row = 0; 659 ilu->inplace = PETSC_FALSE; 660 ierr = PetscStrallocpy(MATORDERING_NATURAL,&ilu->ordering);CHKERRQ(ierr); 661 ilu->reuseordering = PETSC_FALSE; 662 ilu->usedt = PETSC_FALSE; 663 ilu->info.dt = PETSC_DEFAULT; 664 ilu->info.dtcount = PETSC_DEFAULT; 665 ilu->info.dtcol = PETSC_DEFAULT; 666 ilu->info.shiftnz = 0.0; 667 ilu->info.shiftpd = 0.0; /* false */ 668 ilu->info.shift_fraction = 0.0; 669 ilu->info.zeropivot = 1.e-12; 670 ilu->info.pivotinblocks = 1.0; 671 ilu->reusefill = PETSC_FALSE; 672 ilu->info.diagonal_fill = 0; 673 pc->data = (void*)ilu; 674 675 pc->ops->destroy = PCDestroy_ILU; 676 pc->ops->apply = PCApply_ILU; 677 pc->ops->applytranspose = PCApplyTranspose_ILU; 678 pc->ops->setup = PCSetUp_ILU; 679 pc->ops->setfromoptions = PCSetFromOptions_ILU; 680 pc->ops->getfactoredmatrix = PCGetFactoredMatrix_ILU; 681 pc->ops->view = PCView_ILU; 682 pc->ops->applyrichardson = 0; 683 684 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_ILU", 685 PCFactorSetZeroPivot_ILU);CHKERRQ(ierr); 686 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_ILU", 687 PCFactorSetShiftNonzero_ILU);CHKERRQ(ierr); 688 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_ILU", 689 PCFactorSetShiftPd_ILU);CHKERRQ(ierr); 690 691 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseDropTolerance_C","PCFactorSetUseDropTolerance_ILU", 692 PCFactorSetUseDropTolerance_ILU);CHKERRQ(ierr); 693 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_ILU", 694 PCFactorSetFill_ILU);CHKERRQ(ierr); 695 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrdering_C","PCFactorSetMatOrdering_ILU", 696 PCFactorSetMatOrdering_ILU);CHKERRQ(ierr); 697 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU", 698 PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr); 699 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU", 700 PCFactorSetReuseFill_ILU);CHKERRQ(ierr); 701 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_ILU", 702 PCFactorSetLevels_ILU);CHKERRQ(ierr); 703 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU", 704 PCFactorSetUseInPlace_ILU);CHKERRQ(ierr); 705 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_ILU", 706 PCFactorSetAllowDiagonalFill_ILU);CHKERRQ(ierr); 707 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_ILU", 708 PCFactorSetPivotInBlocks_ILU);CHKERRQ(ierr); 709 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU", 710 PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr); 711 PetscFunctionReturn(0); 712 } 713 EXTERN_C_END 714