1 /* 2 Defines a direct factorization preconditioner for any Mat implementation 3 Note: this need not be consided a preconditioner since it supplies 4 a direct solver. 5 */ 6 #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 7 8 typedef struct { 9 Mat fact; /* factored matrix */ 10 PetscReal actualfill; /* actual fill in factor */ 11 PetscTruth inplace; /* flag indicating in-place factorization */ 12 IS row,col; /* index sets used for reordering */ 13 MatOrderingType ordering; /* matrix ordering */ 14 PetscTruth reuseordering; /* reuses previous reordering computed */ 15 PetscTruth reusefill; /* reuse fill from previous LU */ 16 MatFactorInfo info; 17 PetscTruth nonzerosalongdiagonal; 18 PetscReal nonzerosalongdiagonaltol; 19 } PC_LU; 20 21 22 EXTERN_C_BEGIN 23 #undef __FUNCT__ 24 #define __FUNCT__ "PCLUReorderForNonzeroDiagonal_LU" 25 PetscErrorCode PCLUReorderForNonzeroDiagonal_LU(PC pc,PetscReal z) 26 { 27 PC_LU *lu = (PC_LU*)pc->data; 28 29 PetscFunctionBegin; 30 lu->nonzerosalongdiagonal = PETSC_TRUE; 31 if (z == PETSC_DECIDE) { 32 lu->nonzerosalongdiagonaltol = 1.e-10; 33 } else { 34 lu->nonzerosalongdiagonaltol = z; 35 } 36 PetscFunctionReturn(0); 37 } 38 EXTERN_C_END 39 40 EXTERN_C_BEGIN 41 #undef __FUNCT__ 42 #define __FUNCT__ "PCLUSetSetZeroPivot" 43 PetscErrorCode PCLUSetZeroPivot_LU(PC pc,PetscReal z) 44 { 45 PC_LU *lu; 46 47 PetscFunctionBegin; 48 lu = (PC_LU*)pc->data; 49 lu->info.zeropivot = z; 50 PetscFunctionReturn(0); 51 } 52 EXTERN_C_END 53 54 EXTERN_C_BEGIN 55 #undef __FUNCT__ 56 #define __FUNCT__ "PCLUSetReuseOrdering_LU" 57 PetscErrorCode PCLUSetReuseOrdering_LU(PC pc,PetscTruth flag) 58 { 59 PC_LU *lu; 60 61 PetscFunctionBegin; 62 lu = (PC_LU*)pc->data; 63 lu->reuseordering = flag; 64 PetscFunctionReturn(0); 65 } 66 EXTERN_C_END 67 68 EXTERN_C_BEGIN 69 #undef __FUNCT__ 70 #define __FUNCT__ "PCLUSetReuseFill_LU" 71 PetscErrorCode PCLUSetReuseFill_LU(PC pc,PetscTruth flag) 72 { 73 PC_LU *lu; 74 75 PetscFunctionBegin; 76 lu = (PC_LU*)pc->data; 77 lu->reusefill = flag; 78 PetscFunctionReturn(0); 79 } 80 EXTERN_C_END 81 82 EXTERN_C_BEGIN 83 #undef __FUNCT__ 84 #define __FUNCT__ "PCLUSetShift_LU" 85 PetscErrorCode PCLUSetShift_LU(PC pc,PetscTruth shift) 86 { 87 PC_LU *dir; 88 89 PetscFunctionBegin; 90 dir = (PC_LU*)pc->data; 91 dir->info.shift = shift; 92 if (shift) dir->info.shift_fraction = 0.0; 93 PetscFunctionReturn(0); 94 } 95 EXTERN_C_END 96 97 #undef __FUNCT__ 98 #define __FUNCT__ "PCSetFromOptions_LU" 99 static PetscErrorCode PCSetFromOptions_LU(PC pc) 100 { 101 PC_LU *lu = (PC_LU*)pc->data; 102 PetscErrorCode ierr; 103 PetscTruth flg,set; 104 char tname[256]; 105 PetscFList ordlist; 106 PetscReal tol; 107 108 PetscFunctionBegin; 109 ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 110 ierr = PetscOptionsHead("LU options");CHKERRQ(ierr); 111 ierr = PetscOptionsName("-pc_lu_in_place","Form LU in the same memory as the matrix","PCLUSetUseInPlace",&flg);CHKERRQ(ierr); 112 if (flg) { 113 ierr = PCLUSetUseInPlace(pc);CHKERRQ(ierr); 114 } 115 ierr = PetscOptionsReal("-pc_lu_fill","Expected non-zeros in LU/non-zeros in matrix","PCLUSetFill",lu->info.fill,&lu->info.fill,0);CHKERRQ(ierr); 116 117 ierr = PetscOptionsName("-pc_lu_damping","Damping added to diagonal","PCLUSetDamping",&flg);CHKERRQ(ierr); 118 if (flg) { 119 ierr = PCLUSetDamping(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr); 120 } 121 ierr = PetscOptionsReal("-pc_lu_damping","Damping added to diagonal","PCLUSetDamping",lu->info.damping,&lu->info.damping,0);CHKERRQ(ierr); 122 ierr = PetscOptionsName("-pc_lu_shift","Manteuffel shift applied to diagonal","PCLUSetShift",&flg);CHKERRQ(ierr); 123 if (flg) { 124 ierr = PCLUSetShift(pc,PETSC_TRUE);CHKERRQ(ierr); 125 } 126 ierr = PetscOptionsReal("-pc_lu_zeropivot","Pivot is considered zero if less than","PCLUSetSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr); 127 128 ierr = PetscOptionsName("-pc_lu_reuse_fill","Use fill from previous factorization","PCLUSetReuseFill",&flg);CHKERRQ(ierr); 129 if (flg) { 130 ierr = PCLUSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr); 131 } 132 ierr = PetscOptionsName("-pc_lu_reuse_ordering","Reuse ordering from previous factorization","PCLUSetReuseOrdering",&flg);CHKERRQ(ierr); 133 if (flg) { 134 ierr = PCLUSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr); 135 } 136 137 ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 138 ierr = PetscOptionsList("-pc_lu_mat_ordering_type","Reordering to reduce nonzeros in LU","PCLUSetMatOrdering",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr); 139 if (flg) { 140 ierr = PCLUSetMatOrdering(pc,tname);CHKERRQ(ierr); 141 } 142 143 ierr = PetscOptionsName("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCLUReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 144 if (flg) { 145 tol = PETSC_DECIDE; 146 ierr = PetscOptionsReal("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCLUReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 147 ierr = PCLUReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 148 } 149 150 ierr = PetscOptionsReal("-pc_lu_pivoting","Pivoting tolerance (used only for some factorization)","PCLUSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);CHKERRQ(ierr); 151 152 flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 153 ierr = PetscOptionsLogical("-pc_lu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCLUSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 154 if (set) { 155 ierr = PCLUSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 156 } 157 ierr = PetscOptionsTail();CHKERRQ(ierr); 158 PetscFunctionReturn(0); 159 } 160 161 #undef __FUNCT__ 162 #define __FUNCT__ "PCView_LU" 163 static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer) 164 { 165 PC_LU *lu = (PC_LU*)pc->data; 166 PetscErrorCode ierr; 167 PetscTruth iascii,isstring; 168 169 PetscFunctionBegin; 170 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 171 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 172 if (iascii) { 173 MatInfo info; 174 175 if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");CHKERRQ(ierr);} 176 else {ierr = PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");CHKERRQ(ierr);} 177 ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr); 178 ierr = PetscViewerASCIIPrintf(viewer," LU: tolerance for zero pivot %g\n",lu->info.zeropivot);CHKERRQ(ierr); 179 if (lu->info.shift) {ierr = PetscViewerASCIIPrintf(viewer," LU: using Manteuffel shift\n");CHKERRQ(ierr);} 180 if (lu->fact) { 181 ierr = MatGetInfo(lu->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 182 ierr = PetscViewerASCIIPrintf(viewer," LU nonzeros %g\n",info.nz_used);CHKERRQ(ierr); 183 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_FACTOR_INFO);CHKERRQ(ierr); 184 ierr = MatView(lu->fact,viewer);CHKERRQ(ierr); 185 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 186 } 187 if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 188 if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 189 } else if (isstring) { 190 ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 191 } else { 192 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name); 193 } 194 PetscFunctionReturn(0); 195 } 196 197 #undef __FUNCT__ 198 #define __FUNCT__ "PCGetFactoredMatrix_LU" 199 static PetscErrorCode PCGetFactoredMatrix_LU(PC pc,Mat *mat) 200 { 201 PC_LU *dir = (PC_LU*)pc->data; 202 203 PetscFunctionBegin; 204 if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 205 *mat = dir->fact; 206 PetscFunctionReturn(0); 207 } 208 209 #undef __FUNCT__ 210 #define __FUNCT__ "PCSetUp_LU" 211 static PetscErrorCode PCSetUp_LU(PC pc) 212 { 213 PetscErrorCode ierr; 214 PC_LU *dir = (PC_LU*)pc->data; 215 216 PetscFunctionBegin; 217 if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill; 218 219 if (dir->inplace) { 220 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 221 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 222 ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 223 if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);} 224 ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr); 225 dir->fact = pc->pmat; 226 } else { 227 MatInfo info; 228 if (!pc->setupcalled) { 229 ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 230 if (dir->nonzerosalongdiagonal) { 231 ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 232 } 233 if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);} 234 ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 235 ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 236 dir->actualfill = info.fill_ratio_needed; 237 PetscLogObjectParent(pc,dir->fact); 238 } else if (pc->flag != SAME_NONZERO_PATTERN) { 239 if (!dir->reuseordering) { 240 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 241 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 242 ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 243 if (dir->nonzerosalongdiagonal) { 244 ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 245 } 246 if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);} 247 } 248 ierr = MatDestroy(dir->fact);CHKERRQ(ierr); 249 ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 250 ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 251 dir->actualfill = info.fill_ratio_needed; 252 PetscLogObjectParent(pc,dir->fact); 253 } 254 ierr = MatLUFactorNumeric(pc->pmat,&dir->info,&dir->fact);CHKERRQ(ierr); 255 } 256 PetscFunctionReturn(0); 257 } 258 259 #undef __FUNCT__ 260 #define __FUNCT__ "PCDestroy_LU" 261 static PetscErrorCode PCDestroy_LU(PC pc) 262 { 263 PC_LU *dir = (PC_LU*)pc->data; 264 PetscErrorCode ierr; 265 266 PetscFunctionBegin; 267 if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);} 268 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 269 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 270 ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 271 ierr = PetscFree(dir);CHKERRQ(ierr); 272 PetscFunctionReturn(0); 273 } 274 275 #undef __FUNCT__ 276 #define __FUNCT__ "PCApply_LU" 277 static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 278 { 279 PC_LU *dir = (PC_LU*)pc->data; 280 PetscErrorCode ierr; 281 282 PetscFunctionBegin; 283 if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 284 else {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);} 285 PetscFunctionReturn(0); 286 } 287 288 #undef __FUNCT__ 289 #define __FUNCT__ "PCApplyTranspose_LU" 290 static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 291 { 292 PC_LU *dir = (PC_LU*)pc->data; 293 PetscErrorCode ierr; 294 295 PetscFunctionBegin; 296 if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 297 else {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);} 298 PetscFunctionReturn(0); 299 } 300 301 /* -----------------------------------------------------------------------------------*/ 302 303 EXTERN_C_BEGIN 304 #undef __FUNCT__ 305 #define __FUNCT__ "PCLUSetFill_LU" 306 PetscErrorCode PCLUSetFill_LU(PC pc,PetscReal fill) 307 { 308 PC_LU *dir; 309 310 PetscFunctionBegin; 311 dir = (PC_LU*)pc->data; 312 dir->info.fill = fill; 313 PetscFunctionReturn(0); 314 } 315 EXTERN_C_END 316 317 EXTERN_C_BEGIN 318 #undef __FUNCT__ 319 #define __FUNCT__ "PCLUSetDamping_LU" 320 PetscErrorCode PCLUSetDamping_LU(PC pc,PetscReal damping) 321 { 322 PC_LU *dir; 323 324 PetscFunctionBegin; 325 dir = (PC_LU*)pc->data; 326 if (damping == (PetscReal) PETSC_DECIDE) { 327 dir->info.damping = 1.e-12; 328 } else { 329 dir->info.damping = damping; 330 } 331 PetscFunctionReturn(0); 332 } 333 EXTERN_C_END 334 335 EXTERN_C_BEGIN 336 #undef __FUNCT__ 337 #define __FUNCT__ "PCLUSetUseInPlace_LU" 338 PetscErrorCode PCLUSetUseInPlace_LU(PC pc) 339 { 340 PC_LU *dir; 341 342 PetscFunctionBegin; 343 dir = (PC_LU*)pc->data; 344 dir->inplace = PETSC_TRUE; 345 PetscFunctionReturn(0); 346 } 347 EXTERN_C_END 348 349 EXTERN_C_BEGIN 350 #undef __FUNCT__ 351 #define __FUNCT__ "PCLUSetMatOrdering_LU" 352 PetscErrorCode PCLUSetMatOrdering_LU(PC pc,MatOrderingType ordering) 353 { 354 PC_LU *dir = (PC_LU*)pc->data; 355 PetscErrorCode ierr; 356 357 PetscFunctionBegin; 358 ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 359 ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 360 PetscFunctionReturn(0); 361 } 362 EXTERN_C_END 363 364 EXTERN_C_BEGIN 365 #undef __FUNCT__ 366 #define __FUNCT__ "PCLUSetPivoting_LU" 367 PetscErrorCode PCLUSetPivoting_LU(PC pc,PetscReal dtcol) 368 { 369 PC_LU *dir = (PC_LU*)pc->data; 370 371 PetscFunctionBegin; 372 if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %g must be between 0 and 1",dtcol); 373 dir->info.dtcol = dtcol; 374 PetscFunctionReturn(0); 375 } 376 EXTERN_C_END 377 378 EXTERN_C_BEGIN 379 #undef __FUNCT__ 380 #define __FUNCT__ "PCLUSetPivotInBlocks_LU" 381 PetscErrorCode PCLUSetPivotInBlocks_LU(PC pc,PetscTruth pivot) 382 { 383 PC_LU *dir = (PC_LU*)pc->data; 384 385 PetscFunctionBegin; 386 dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 387 PetscFunctionReturn(0); 388 } 389 EXTERN_C_END 390 391 /* -----------------------------------------------------------------------------------*/ 392 393 #undef __FUNCT__ 394 #define __FUNCT__ "PCLUReorderForNonzeroDiagonal" 395 /*@ 396 PCLUReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 397 398 Collective on PC 399 400 Input Parameters: 401 + pc - the preconditioner context 402 - tol - diagonal entries smaller than this in absolute value are considered zero 403 404 Options Database Key: 405 . -pc_lu_nonzeros_along_diagonal 406 407 Level: intermediate 408 409 .keywords: PC, set, factorization, direct, fill 410 411 .seealso: PCLUSetFill(), PCLUSetDamp(), PCILUSetZeroPivot(), MatReorderForNonzeroDiagonal() 412 @*/ 413 PetscErrorCode PCLUReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 414 { 415 PetscErrorCode ierr,(*f)(PC,PetscReal); 416 417 PetscFunctionBegin; 418 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 419 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr); 420 if (f) { 421 ierr = (*f)(pc,rtol);CHKERRQ(ierr); 422 } 423 PetscFunctionReturn(0); 424 } 425 426 #undef __FUNCT__ 427 #define __FUNCT__ "PCLUSetZeroPivot" 428 /*@ 429 PCLUSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 430 431 Collective on PC 432 433 Input Parameters: 434 + pc - the preconditioner context 435 - zero - all pivots smaller than this will be considered zero 436 437 Options Database Key: 438 . -pc_lu_zeropivot <zero> - Sets the zero pivot size 439 440 Level: intermediate 441 442 .keywords: PC, set, factorization, direct, fill 443 444 .seealso: PCLUSetFill(), PCLUSetDamp(), PCILUSetZeroPivot() 445 @*/ 446 PetscErrorCode PCLUSetZeroPivot(PC pc,PetscReal zero) 447 { 448 PetscErrorCode ierr,(*f)(PC,PetscReal); 449 450 PetscFunctionBegin; 451 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 452 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetZeroPivot_C",(void (**)(void))&f);CHKERRQ(ierr); 453 if (f) { 454 ierr = (*f)(pc,zero);CHKERRQ(ierr); 455 } 456 PetscFunctionReturn(0); 457 } 458 459 #undef __FUNCT__ 460 #define __FUNCT__ "PCLUSetShift" 461 /*@ 462 PCLUSetShift - specify whether to use Manteuffel shifting of LU. 463 If an LU factorisation breaks down because of nonpositive pivots, 464 adding sufficient identity to the diagonal will remedy this. 465 Setting this causes a bisection method to find the minimum shift that 466 will lead to a well-defined LU. 467 468 Input parameters: 469 + pc - the preconditioner context 470 - shifting - PETSC_TRUE to set shift else PETSC_FALSE 471 472 Options Database Key: 473 . -pc_lu_shift - Activate PCLUSetShift() 474 475 Level: intermediate 476 477 .keywords: PC, indefinite, factorization 478 479 .seealso: PCLUSetDamping(), PCILUSetShift() 480 @*/ 481 PetscErrorCode PCLUSetShift(PC pc,PetscTruth shifting) 482 { 483 PetscErrorCode ierr,(*f)(PC,PetscTruth); 484 485 PetscFunctionBegin; 486 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 487 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetShift_C",(void (**)(void))&f);CHKERRQ(ierr); 488 if (f) { 489 ierr = (*f)(pc,shifting);CHKERRQ(ierr); 490 } 491 PetscFunctionReturn(0); 492 } 493 494 495 #undef __FUNCT__ 496 #define __FUNCT__ "PCLUSetReuseOrdering" 497 /*@ 498 PCLUSetReuseOrdering - When similar matrices are factored, this 499 causes the ordering computed in the first factor to be used for all 500 following factors; applies to both fill and drop tolerance LUs. 501 502 Collective on PC 503 504 Input Parameters: 505 + pc - the preconditioner context 506 - flag - PETSC_TRUE to reuse else PETSC_FALSE 507 508 Options Database Key: 509 . -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering() 510 511 Level: intermediate 512 513 .keywords: PC, levels, reordering, factorization, incomplete, LU 514 515 .seealso: PCLUSetReuseFill(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill() 516 @*/ 517 PetscErrorCode PCLUSetReuseOrdering(PC pc,PetscTruth flag) 518 { 519 PetscErrorCode ierr,(*f)(PC,PetscTruth); 520 521 PetscFunctionBegin; 522 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 523 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 524 if (f) { 525 ierr = (*f)(pc,flag);CHKERRQ(ierr); 526 } 527 PetscFunctionReturn(0); 528 } 529 530 #undef __FUNCT__ 531 #define __FUNCT__ "PCLUSetReuseFill" 532 /*@ 533 PCLUSetReuseFill - When matrices with same nonzero structure are LU factored, 534 this causes later ones to use the fill computed in the initial factorization. 535 536 Collective on PC 537 538 Input Parameters: 539 + pc - the preconditioner context 540 - flag - PETSC_TRUE to reuse else PETSC_FALSE 541 542 Options Database Key: 543 . -pc_lu_reuse_fill - Activates PCLUSetReuseFill() 544 545 Level: intermediate 546 547 .keywords: PC, levels, reordering, factorization, incomplete, LU 548 549 .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCILUDTSetReuseFill() 550 @*/ 551 PetscErrorCode PCLUSetReuseFill(PC pc,PetscTruth flag) 552 { 553 PetscErrorCode ierr,(*f)(PC,PetscTruth); 554 555 PetscFunctionBegin; 556 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 557 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr); 558 if (f) { 559 ierr = (*f)(pc,flag);CHKERRQ(ierr); 560 } 561 PetscFunctionReturn(0); 562 } 563 564 #undef __FUNCT__ 565 #define __FUNCT__ "PCLUSetFill" 566 /*@ 567 PCLUSetFill - Indicate the amount of fill you expect in the factored matrix, 568 fill = number nonzeros in factor/number nonzeros in original matrix. 569 570 Collective on PC 571 572 Input Parameters: 573 + pc - the preconditioner context 574 - fill - amount of expected fill 575 576 Options Database Key: 577 . -pc_lu_fill <fill> - Sets fill amount 578 579 Level: intermediate 580 581 Note: 582 For sparse matrix factorizations it is difficult to predict how much 583 fill to expect. By running with the option -log_info PETSc will print the 584 actual amount of fill used; allowing you to set the value accurately for 585 future runs. Default PETSc uses a value of 5.0 586 587 .keywords: PC, set, factorization, direct, fill 588 589 .seealso: PCILUSetFill() 590 @*/ 591 PetscErrorCode PCLUSetFill(PC pc,PetscReal fill) 592 { 593 PetscErrorCode ierr,(*f)(PC,PetscReal); 594 595 PetscFunctionBegin; 596 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 597 if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0"); 598 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 599 if (f) { 600 ierr = (*f)(pc,fill);CHKERRQ(ierr); 601 } 602 PetscFunctionReturn(0); 603 } 604 605 #undef __FUNCT__ 606 #define __FUNCT__ "PCLUSetDamping" 607 /*@ 608 PCLUSetDamping - adds this quantity to the diagonal of the matrix during the 609 LU numerical factorization 610 611 Collective on PC 612 613 Input Parameters: 614 + pc - the preconditioner context 615 - damping - amount of damping (use PETSC_DECIDE for default of 1.e-12) 616 617 Options Database Key: 618 . -pc_lu_damping <damping> - Sets damping amount or PETSC_DECIDE for the default 619 620 Note: If 0.0 is given, then no damping is used. If a diagonal element is classified as a zero 621 pivot, then the damping is doubled until this is alleviated. 622 623 Level: intermediate 624 625 .keywords: PC, set, factorization, direct, fill 626 .seealso: PCILUSetFill(), PCILUSetDamp() 627 @*/ 628 PetscErrorCode PCLUSetDamping(PC pc,PetscReal damping) 629 { 630 PetscErrorCode ierr,(*f)(PC,PetscReal); 631 632 PetscFunctionBegin; 633 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 634 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 635 if (f) { 636 ierr = (*f)(pc,damping);CHKERRQ(ierr); 637 } 638 PetscFunctionReturn(0); 639 } 640 641 #undef __FUNCT__ 642 #define __FUNCT__ "PCLUSetUseInPlace" 643 /*@ 644 PCLUSetUseInPlace - Tells the system to do an in-place factorization. 645 For dense matrices, this enables the solution of much larger problems. 646 For sparse matrices the factorization cannot be done truly in-place 647 so this does not save memory during the factorization, but after the matrix 648 is factored, the original unfactored matrix is freed, thus recovering that 649 space. 650 651 Collective on PC 652 653 Input Parameters: 654 . pc - the preconditioner context 655 656 Options Database Key: 657 . -pc_lu_in_place - Activates in-place factorization 658 659 Notes: 660 PCLUSetUseInplace() can only be used with the KSP method KSPPREONLY or when 661 a different matrix is provided for the multiply and the preconditioner in 662 a call to KSPSetOperators(). 663 This is because the Krylov space methods require an application of the 664 matrix multiplication, which is not possible here because the matrix has 665 been factored in-place, replacing the original matrix. 666 667 Level: intermediate 668 669 .keywords: PC, set, factorization, direct, inplace, in-place, LU 670 671 .seealso: PCILUSetUseInPlace() 672 @*/ 673 PetscErrorCode PCLUSetUseInPlace(PC pc) 674 { 675 PetscErrorCode ierr,(*f)(PC); 676 677 PetscFunctionBegin; 678 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 679 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr); 680 if (f) { 681 ierr = (*f)(pc);CHKERRQ(ierr); 682 } 683 PetscFunctionReturn(0); 684 } 685 686 #undef __FUNCT__ 687 #define __FUNCT__ "PCLUSetMatOrdering" 688 /*@C 689 PCLUSetMatOrdering - Sets the ordering routine (to reduce fill) to 690 be used in the LU factorization. 691 692 Collective on PC 693 694 Input Parameters: 695 + pc - the preconditioner context 696 - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 697 698 Options Database Key: 699 . -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine 700 701 Level: intermediate 702 703 Notes: nested dissection is used by default 704 705 .seealso: PCILUSetMatOrdering() 706 @*/ 707 PetscErrorCode PCLUSetMatOrdering(PC pc,MatOrderingType ordering) 708 { 709 PetscErrorCode ierr,(*f)(PC,MatOrderingType); 710 711 PetscFunctionBegin; 712 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 713 if (f) { 714 ierr = (*f)(pc,ordering);CHKERRQ(ierr); 715 } 716 PetscFunctionReturn(0); 717 } 718 719 #undef __FUNCT__ 720 #define __FUNCT__ "PCLUSetPivoting" 721 /*@ 722 PCLUSetPivoting - Determines when pivoting is done during LU. 723 For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 724 it is never done. For the Matlab and SuperLU factorization this is used. 725 726 Collective on PC 727 728 Input Parameters: 729 + pc - the preconditioner context 730 - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 731 732 Options Database Key: 733 . -pc_lu_pivoting <dtcol> 734 735 Level: intermediate 736 737 .seealso: PCILUSetMatOrdering(), PCLUSetPivotInBlocks() 738 @*/ 739 PetscErrorCode PCLUSetPivoting(PC pc,PetscReal dtcol) 740 { 741 PetscErrorCode ierr,(*f)(PC,PetscReal); 742 743 PetscFunctionBegin; 744 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr); 745 if (f) { 746 ierr = (*f)(pc,dtcol);CHKERRQ(ierr); 747 } 748 PetscFunctionReturn(0); 749 } 750 751 #undef __FUNCT__ 752 #define __FUNCT__ "PCLUSetPivotInBlocks" 753 /*@ 754 PCLUSetPivotInBlocks - Determines if pivoting is done while factoring each block 755 with BAIJ or SBAIJ matrices 756 757 Collective on PC 758 759 Input Parameters: 760 + pc - the preconditioner context 761 - pivot - PETSC_TRUE or PETSC_FALSE 762 763 Options Database Key: 764 . -pc_lu_pivot_in_blocks <true,false> 765 766 Level: intermediate 767 768 .seealso: PCILUSetMatOrdering(), PCLUSetPivoting() 769 @*/ 770 PetscErrorCode PCLUSetPivotInBlocks(PC pc,PetscTruth pivot) 771 { 772 PetscErrorCode ierr,(*f)(PC,PetscTruth); 773 774 PetscFunctionBegin; 775 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 776 if (f) { 777 ierr = (*f)(pc,pivot);CHKERRQ(ierr); 778 } 779 PetscFunctionReturn(0); 780 } 781 782 /* ------------------------------------------------------------------------ */ 783 784 /*MC 785 PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 786 787 Options Database Keys: 788 + -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering() 789 . -pc_lu_reuse_fill - Activates PCLUSetReuseFill() 790 . -pc_lu_fill <fill> - Sets fill amount 791 . -pc_lu_damping <damping> - Sets damping amount 792 . -pc_lu_in_place - Activates in-place factorization 793 . -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine 794 . -pc_lu_pivoting <dtcol> 795 - -pc_lu_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 796 stability of factorization. 797 798 Notes: Not all options work for all matrix formats 799 Run with -help to see additional options for particular matrix formats or factorization 800 algorithms 801 802 Level: beginner 803 804 Concepts: LU factorization, direct solver 805 806 Notes: Usually this will compute an "exact" solution in one iteration and does 807 not need a Krylov method (i.e. you can use -ksp_type preonly, or 808 KSPSetType(ksp,KSPPREONLY) for the Krylov method 809 810 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 811 PCILU, PCCHOLESKY, PCICC, PCLUSetReuseOrdering(), PCLUSetReuseFill(), PCGetFactoredMatrix(), 812 PCLUSetFill(), PCLUSetDamping(), PCLUSetUseInPlace(), PCLUSetMatOrdering(), PCLUSetPivoting(), 813 PCLUSetPivotingInBlocks() 814 M*/ 815 816 EXTERN_C_BEGIN 817 #undef __FUNCT__ 818 #define __FUNCT__ "PCCreate_LU" 819 PetscErrorCode PCCreate_LU(PC pc) 820 { 821 PetscErrorCode ierr; 822 PetscMPIInt size; 823 PC_LU *dir; 824 825 PetscFunctionBegin; 826 ierr = PetscNew(PC_LU,&dir);CHKERRQ(ierr); 827 PetscLogObjectMemory(pc,sizeof(PC_LU)); 828 829 ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr); 830 dir->fact = 0; 831 dir->inplace = PETSC_FALSE; 832 dir->nonzerosalongdiagonal = PETSC_FALSE; 833 834 dir->info.fill = 5.0; 835 dir->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 836 dir->info.damping = 0.0; 837 dir->info.zeropivot = 1.e-12; 838 dir->info.pivotinblocks = 1.0; 839 dir->info.shift = PETSC_FALSE; 840 dir->info.shift_fraction = 0.0; 841 dir->col = 0; 842 dir->row = 0; 843 ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr); 844 if (size == 1) { 845 ierr = PetscStrallocpy(MATORDERING_ND,&dir->ordering);CHKERRQ(ierr); 846 } else { 847 ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr); 848 } 849 dir->reusefill = PETSC_FALSE; 850 dir->reuseordering = PETSC_FALSE; 851 pc->data = (void*)dir; 852 853 pc->ops->destroy = PCDestroy_LU; 854 pc->ops->apply = PCApply_LU; 855 pc->ops->applytranspose = PCApplyTranspose_LU; 856 pc->ops->setup = PCSetUp_LU; 857 pc->ops->setfromoptions = PCSetFromOptions_LU; 858 pc->ops->view = PCView_LU; 859 pc->ops->applyrichardson = 0; 860 pc->ops->getfactoredmatrix = PCGetFactoredMatrix_LU; 861 862 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetFill_C","PCLUSetFill_LU", 863 PCLUSetFill_LU);CHKERRQ(ierr); 864 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetDamping_C","PCLUSetDamping_LU", 865 PCLUSetDamping_LU);CHKERRQ(ierr); 866 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetShift_C","PCLUSetShift_LU", 867 PCLUSetShift_LU);CHKERRQ(ierr); 868 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetUseInPlace_C","PCLUSetUseInPlace_LU", 869 PCLUSetUseInPlace_LU);CHKERRQ(ierr); 870 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetMatOrdering_C","PCLUSetMatOrdering_LU", 871 PCLUSetMatOrdering_LU);CHKERRQ(ierr); 872 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseOrdering_C","PCLUSetReuseOrdering_LU", 873 PCLUSetReuseOrdering_LU);CHKERRQ(ierr); 874 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseFill_C","PCLUSetReuseFill_LU", 875 PCLUSetReuseFill_LU);CHKERRQ(ierr); 876 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivoting_C","PCLUSetPivoting_LU", 877 PCLUSetPivoting_LU);CHKERRQ(ierr); 878 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivotInBlocks_C","PCLUSetPivotInBlocks_LU", 879 PCLUSetPivotInBlocks_LU);CHKERRQ(ierr); 880 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetZeroPivot_C","PCLUSetZeroPivot_LU", 881 PCLUSetZeroPivot_LU);CHKERRQ(ierr); 882 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUReorderForNonzeroDiagonal_C","PCLUReorderForNonzeroDiagonal_LU", 883 PCLUReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 884 PetscFunctionReturn(0); 885 } 886 EXTERN_C_END 887