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