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 "include/private/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__ "PCFactorSetMatOrderingType_ILU" 148 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_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 There are NO default values for the 3 parameters, you must set them with reasonable values for your 261 matrix. We don't know how to compute reasonable values. 262 263 .keywords: PC, levels, reordering, factorization, incomplete, ILU 264 @*/ 265 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount) 266 { 267 PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt); 268 269 PetscFunctionBegin; 270 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 271 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 272 if (f) { 273 ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr); 274 } 275 PetscFunctionReturn(0); 276 } 277 278 #undef __FUNCT__ 279 #define __FUNCT__ "PCFactorSetLevels" 280 /*@ 281 PCFactorSetLevels - Sets the number of levels of fill to use. 282 283 Collective on PC 284 285 Input Parameters: 286 + pc - the preconditioner context 287 - levels - number of levels of fill 288 289 Options Database Key: 290 . -pc_factor_levels <levels> - Sets fill level 291 292 Level: intermediate 293 294 .keywords: PC, levels, fill, factorization, incomplete, ILU 295 @*/ 296 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC pc,PetscInt levels) 297 { 298 PetscErrorCode ierr,(*f)(PC,PetscInt); 299 300 PetscFunctionBegin; 301 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 302 if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 303 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr); 304 if (f) { 305 ierr = (*f)(pc,levels);CHKERRQ(ierr); 306 } 307 PetscFunctionReturn(0); 308 } 309 310 #undef __FUNCT__ 311 #define __FUNCT__ "PCFactorSetAllowDiagonalFill" 312 /*@ 313 PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 314 treated as level 0 fill even if there is no non-zero location. 315 316 Collective on PC 317 318 Input Parameters: 319 + pc - the preconditioner context 320 321 Options Database Key: 322 . -pc_factor_diagonal_fill 323 324 Notes: 325 Does not apply with 0 fill. 326 327 Level: intermediate 328 329 .keywords: PC, levels, fill, factorization, incomplete, ILU 330 @*/ 331 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC pc) 332 { 333 PetscErrorCode ierr,(*f)(PC); 334 335 PetscFunctionBegin; 336 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 337 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr); 338 if (f) { 339 ierr = (*f)(pc);CHKERRQ(ierr); 340 } 341 PetscFunctionReturn(0); 342 } 343 344 /* ------------------------------------------------------------------------------------------*/ 345 346 EXTERN_C_BEGIN 347 #undef __FUNCT__ 348 #define __FUNCT__ "PCFactorSetPivotInBlocks_ILU" 349 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_ILU(PC pc,PetscTruth pivot) 350 { 351 PC_ILU *dir = (PC_ILU*)pc->data; 352 353 PetscFunctionBegin; 354 dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 355 PetscFunctionReturn(0); 356 } 357 EXTERN_C_END 358 359 #undef __FUNCT__ 360 #define __FUNCT__ "PCSetFromOptions_ILU" 361 static PetscErrorCode PCSetFromOptions_ILU(PC pc) 362 { 363 PetscErrorCode ierr; 364 PetscInt dtmax = 3,itmp; 365 PetscTruth flg,set; 366 PetscReal dt[3]; 367 char tname[256]; 368 PC_ILU *ilu = (PC_ILU*)pc->data; 369 PetscFList ordlist; 370 PetscReal tol; 371 372 PetscFunctionBegin; 373 ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 374 ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr); 375 ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)ilu->info.levels,&itmp,&flg);CHKERRQ(ierr); 376 if (flg) ilu->info.levels = itmp; 377 ierr = PetscOptionsName("-pc_factor_in_place","do factorization in place","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr); 378 if (flg) ilu->inplace = PETSC_TRUE; 379 ierr = PetscOptionsName("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",&flg);CHKERRQ(ierr); 380 ilu->info.diagonal_fill = (double) flg; 381 ierr = PetscOptionsName("-pc_factor_reuse_fill","Reuse fill ratio from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr); 382 if (flg) ilu->reusefill = PETSC_TRUE; 383 ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse previous reordering","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr); 384 if (flg) ilu->reuseordering = PETSC_TRUE; 385 ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 386 if (flg) { 387 ierr = PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr); 388 } 389 ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",ilu->info.shiftnz,&ilu->info.shiftnz,0);CHKERRQ(ierr); 390 391 ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr); 392 if (flg) { 393 ierr = PetscOptionsInt("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",(PetscInt)ilu->info.shiftpd,&itmp,&flg); CHKERRQ(ierr); 394 if (flg && !itmp) { 395 ierr = PCFactorSetShiftPd(pc,PETSC_FALSE);CHKERRQ(ierr); 396 } else { 397 ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr); 398 } 399 } 400 ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",ilu->info.zeropivot,&ilu->info.zeropivot,0);CHKERRQ(ierr); 401 402 dt[0] = ilu->info.dt; 403 dt[1] = ilu->info.dtcol; 404 dt[2] = ilu->info.dtcount; 405 ierr = PetscOptionsRealArray("-pc_factor_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetUseDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr); 406 if (flg) { 407 ierr = PCFactorSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr); 408 } 409 ierr = PetscOptionsReal("-pc_factor_fill","Expected fill in factorization","PCFactorSetFill",ilu->info.fill,&ilu->info.fill,&flg);CHKERRQ(ierr); 410 ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 411 if (flg) { 412 tol = PETSC_DECIDE; 413 ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 414 ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 415 } 416 417 ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 418 ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCFactorSetMatOrderingType",ordlist,ilu->ordering,tname,256,&flg);CHKERRQ(ierr); 419 if (flg) { 420 ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 421 } 422 flg = ilu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 423 ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 424 if (set) { 425 ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 426 } 427 ierr = PetscOptionsTail();CHKERRQ(ierr); 428 PetscFunctionReturn(0); 429 } 430 431 #undef __FUNCT__ 432 #define __FUNCT__ "PCView_ILU" 433 static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer) 434 { 435 PC_ILU *ilu = (PC_ILU*)pc->data; 436 PetscErrorCode ierr; 437 PetscTruth isstring,iascii; 438 439 PetscFunctionBegin; 440 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 441 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 442 if (iascii) { 443 if (ilu->usedt) { 444 ierr = PetscViewerASCIIPrintf(viewer," ILU: drop tolerance %G\n",ilu->info.dt);CHKERRQ(ierr); 445 ierr = PetscViewerASCIIPrintf(viewer," ILU: max nonzeros per row %D\n",(PetscInt)ilu->info.dtcount);CHKERRQ(ierr); 446 ierr = PetscViewerASCIIPrintf(viewer," ILU: column permutation tolerance %G\n",ilu->info.dtcol);CHKERRQ(ierr); 447 } else if (ilu->info.levels == 1) { 448 ierr = PetscViewerASCIIPrintf(viewer," ILU: %D level of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr); 449 } else { 450 ierr = PetscViewerASCIIPrintf(viewer," ILU: %D levels of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr); 451 } 452 ierr = PetscViewerASCIIPrintf(viewer," ILU: factor fill ratio allocated %G\n",ilu->info.fill);CHKERRQ(ierr); 453 ierr = PetscViewerASCIIPrintf(viewer," ILU: tolerance for zero pivot %G\n",ilu->info.zeropivot);CHKERRQ(ierr); 454 if (ilu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," ILU: using Manteuffel shift\n");CHKERRQ(ierr);} 455 if (ilu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," in-place factorization\n");CHKERRQ(ierr);} 456 else {ierr = PetscViewerASCIIPrintf(viewer," out-of-place factorization\n");CHKERRQ(ierr);} 457 ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",ilu->ordering);CHKERRQ(ierr); 458 if (ilu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 459 if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 460 if (ilu->fact) { 461 ierr = PetscViewerASCIIPrintf(viewer," ILU: factor fill ratio needed %G\n",ilu->actualfill);CHKERRQ(ierr); 462 ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 463 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 464 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 465 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 466 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 467 ierr = MatView(ilu->fact,viewer);CHKERRQ(ierr); 468 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 469 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 470 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 471 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 472 } 473 } else if (isstring) { 474 ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)ilu->info.levels,ilu->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 475 } else { 476 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name); 477 } 478 PetscFunctionReturn(0); 479 } 480 481 #undef __FUNCT__ 482 #define __FUNCT__ "PCSetUp_ILU" 483 static PetscErrorCode PCSetUp_ILU(PC pc) 484 { 485 PetscErrorCode ierr; 486 PC_ILU *ilu = (PC_ILU*)pc->data; 487 MatInfo info; 488 489 PetscFunctionBegin; 490 if (ilu->inplace) { 491 if (!pc->setupcalled) { 492 493 /* In-place factorization only makes sense with the natural ordering, 494 so we only need to get the ordering once, even if nonzero structure changes */ 495 ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 496 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 497 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 498 } 499 500 /* In place ILU only makes sense with fill factor of 1.0 because 501 cannot have levels of fill */ 502 ilu->info.fill = 1.0; 503 ilu->info.diagonal_fill = 0; 504 ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&ilu->info);CHKERRQ(ierr); 505 ilu->fact = pc->pmat; 506 } else if (ilu->usedt) { 507 if (!pc->setupcalled) { 508 ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 509 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 510 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 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 (pc->flag != SAME_NONZERO_PATTERN) { 514 ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 515 if (!ilu->reuseordering) { 516 if (ilu->row) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 517 if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 518 ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 519 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 520 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 521 } 522 ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 523 ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 524 } else if (!ilu->reusefill) { 525 ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 526 ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 527 ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 528 } else { 529 ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr); 530 } 531 } else { 532 if (!pc->setupcalled) { 533 /* first time in so compute reordering and symbolic factorization */ 534 ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 535 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 536 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 537 /* Remove zeros along diagonal? */ 538 if (ilu->nonzerosalongdiagonal) { 539 ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 540 } 541 ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 542 ierr = MatGetInfo(ilu->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 543 ilu->actualfill = info.fill_ratio_needed; 544 ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 545 } else if (pc->flag != SAME_NONZERO_PATTERN) { 546 if (!ilu->reuseordering) { 547 /* compute a new ordering for the ILU */ 548 ierr = ISDestroy(ilu->row);CHKERRQ(ierr); 549 ierr = ISDestroy(ilu->col);CHKERRQ(ierr); 550 ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 551 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 552 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 553 /* Remove zeros along diagonal? */ 554 if (ilu->nonzerosalongdiagonal) { 555 ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 556 } 557 } 558 ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 559 ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 560 ierr = MatGetInfo(ilu->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 561 ilu->actualfill = info.fill_ratio_needed; 562 ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 563 } 564 ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr); 565 } 566 PetscFunctionReturn(0); 567 } 568 569 #undef __FUNCT__ 570 #define __FUNCT__ "PCDestroy_ILU" 571 static PetscErrorCode PCDestroy_ILU(PC pc) 572 { 573 PC_ILU *ilu = (PC_ILU*)pc->data; 574 PetscErrorCode ierr; 575 576 PetscFunctionBegin; 577 ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 578 ierr = PetscStrfree(ilu->ordering);CHKERRQ(ierr); 579 ierr = PetscFree(ilu);CHKERRQ(ierr); 580 PetscFunctionReturn(0); 581 } 582 583 #undef __FUNCT__ 584 #define __FUNCT__ "PCApply_ILU" 585 static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 586 { 587 PC_ILU *ilu = (PC_ILU*)pc->data; 588 PetscErrorCode ierr; 589 590 PetscFunctionBegin; 591 ierr = MatSolve(ilu->fact,x,y);CHKERRQ(ierr); 592 PetscFunctionReturn(0); 593 } 594 595 #undef __FUNCT__ 596 #define __FUNCT__ "PCApplyTranspose_ILU" 597 static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 598 { 599 PC_ILU *ilu = (PC_ILU*)pc->data; 600 PetscErrorCode ierr; 601 602 PetscFunctionBegin; 603 ierr = MatSolveTranspose(ilu->fact,x,y);CHKERRQ(ierr); 604 PetscFunctionReturn(0); 605 } 606 607 #undef __FUNCT__ 608 #define __FUNCT__ "PCGetFactoredMatrix_ILU" 609 static PetscErrorCode PCGetFactoredMatrix_ILU(PC pc,Mat *mat) 610 { 611 PC_ILU *ilu = (PC_ILU*)pc->data; 612 613 PetscFunctionBegin; 614 if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 615 *mat = ilu->fact; 616 PetscFunctionReturn(0); 617 } 618 619 /*MC 620 PCILU - Incomplete factorization preconditioners. 621 622 Options Database Keys: 623 + -pc_factor_levels <k> - number of levels of fill for ILU(k) 624 . -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 625 its factorization (overwrites original matrix) 626 . -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 627 . -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization 628 . -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt 629 . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 630 . -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 631 this decreases the chance of getting a zero pivot 632 . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 633 . -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 634 than 1 the diagonal blocks are factored with partial pivoting (this increases the 635 stability of the ILU factorization 636 . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 637 - -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 638 is optional with PETSC_TRUE being the default 639 640 Level: beginner 641 642 Concepts: incomplete factorization 643 644 Notes: Only implemented for some matrix formats. Not implemented in parallel 645 646 For BAIJ matrices this implements a point block ILU 647 648 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 649 PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), PCFactorSetUseDropTolerance(), 650 PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 651 PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(), 652 PCFactorSetShiftNonzero(),PCFactorSetShiftPd() 653 654 M*/ 655 656 EXTERN_C_BEGIN 657 #undef __FUNCT__ 658 #define __FUNCT__ "PCCreate_ILU" 659 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ILU(PC pc) 660 { 661 PetscErrorCode ierr; 662 PC_ILU *ilu; 663 664 PetscFunctionBegin; 665 ierr = PetscNewLog(pc,PC_ILU,&ilu);CHKERRQ(ierr); 666 667 ilu->fact = 0; 668 ierr = MatFactorInfoInitialize(&ilu->info);CHKERRQ(ierr); 669 ilu->info.levels = 0; 670 ilu->info.fill = 1.0; 671 ilu->col = 0; 672 ilu->row = 0; 673 ilu->inplace = PETSC_FALSE; 674 ierr = PetscStrallocpy(MATORDERING_NATURAL,&ilu->ordering);CHKERRQ(ierr); 675 ilu->reuseordering = PETSC_FALSE; 676 ilu->usedt = PETSC_FALSE; 677 ilu->info.dt = PETSC_DEFAULT; 678 ilu->info.dtcount = PETSC_DEFAULT; 679 ilu->info.dtcol = PETSC_DEFAULT; 680 ilu->info.shiftnz = 0.0; 681 ilu->info.shiftpd = 0.0; /* false */ 682 ilu->info.shift_fraction = 0.0; 683 ilu->info.zeropivot = 1.e-12; 684 ilu->info.pivotinblocks = 1.0; 685 ilu->reusefill = PETSC_FALSE; 686 ilu->info.diagonal_fill = 0; 687 pc->data = (void*)ilu; 688 689 pc->ops->destroy = PCDestroy_ILU; 690 pc->ops->apply = PCApply_ILU; 691 pc->ops->applytranspose = PCApplyTranspose_ILU; 692 pc->ops->setup = PCSetUp_ILU; 693 pc->ops->setfromoptions = PCSetFromOptions_ILU; 694 pc->ops->getfactoredmatrix = PCGetFactoredMatrix_ILU; 695 pc->ops->view = PCView_ILU; 696 pc->ops->applyrichardson = 0; 697 698 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_ILU", 699 PCFactorSetZeroPivot_ILU);CHKERRQ(ierr); 700 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_ILU", 701 PCFactorSetShiftNonzero_ILU);CHKERRQ(ierr); 702 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_ILU", 703 PCFactorSetShiftPd_ILU);CHKERRQ(ierr); 704 705 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseDropTolerance_C","PCFactorSetUseDropTolerance_ILU", 706 PCFactorSetUseDropTolerance_ILU);CHKERRQ(ierr); 707 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_ILU", 708 PCFactorSetFill_ILU);CHKERRQ(ierr); 709 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_ILU", 710 PCFactorSetMatOrderingType_ILU);CHKERRQ(ierr); 711 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU", 712 PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr); 713 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU", 714 PCFactorSetReuseFill_ILU);CHKERRQ(ierr); 715 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_ILU", 716 PCFactorSetLevels_ILU);CHKERRQ(ierr); 717 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU", 718 PCFactorSetUseInPlace_ILU);CHKERRQ(ierr); 719 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_ILU", 720 PCFactorSetAllowDiagonalFill_ILU);CHKERRQ(ierr); 721 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_ILU", 722 PCFactorSetPivotInBlocks_ILU);CHKERRQ(ierr); 723 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU", 724 PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr); 725 PetscFunctionReturn(0); 726 } 727 EXTERN_C_END 728