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