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__ "PCFactorReorderForNonzeroDiagonal_LU" 66 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_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__ "PCFactorSetReuseOrdering_LU" 84 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_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__ "PCFactorSetReuseFill_LU" 98 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_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_factor_in_place","Form LU in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr); 124 if (flg) { 125 ierr = PCFactorSetUseInPlace(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_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr); 141 if (flg) { 142 ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr); 143 } 144 ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr); 145 if (flg) { 146 ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr); 147 } 148 149 ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 150 ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in LU","PCFactorSetMatOrderingType",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr); 151 if (flg) { 152 ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 153 } 154 155 ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 156 if (flg) { 157 tol = PETSC_DECIDE; 158 ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 159 ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 160 } 161 162 ierr = PetscOptionsReal("-pc_factor_pivoting","Pivoting tolerance (used only for some factorization)","PCFactorSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);CHKERRQ(ierr); 163 164 flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 165 ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 166 if (set) { 167 ierr = PCFactorSetPivotInBlocks(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 186 if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");CHKERRQ(ierr);} 187 else {ierr = PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");CHKERRQ(ierr);} 188 ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr); 189 ierr = PetscViewerASCIIPrintf(viewer," LU: tolerance for zero pivot %G\n",lu->info.zeropivot);CHKERRQ(ierr); 190 if (lu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," LU: using Manteuffel shift\n");CHKERRQ(ierr);} 191 if (lu->fact) { 192 ierr = PetscViewerASCIIPrintf(viewer," LU: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr); 193 ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 194 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 195 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 196 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 197 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 198 ierr = MatView(lu->fact,viewer);CHKERRQ(ierr); 199 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 200 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 201 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 202 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 203 } 204 if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 205 if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 206 } else if (isstring) { 207 ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 208 } else { 209 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name); 210 } 211 PetscFunctionReturn(0); 212 } 213 214 #undef __FUNCT__ 215 #define __FUNCT__ "PCFactorGetMatrix_LU" 216 static PetscErrorCode PCFactorGetMatrix_LU(PC pc,Mat *mat) 217 { 218 PC_LU *dir = (PC_LU*)pc->data; 219 220 PetscFunctionBegin; 221 if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 222 *mat = dir->fact; 223 PetscFunctionReturn(0); 224 } 225 226 #undef __FUNCT__ 227 #define __FUNCT__ "PCSetUp_LU" 228 static PetscErrorCode PCSetUp_LU(PC pc) 229 { 230 PetscErrorCode ierr; 231 PC_LU *dir = (PC_LU*)pc->data; 232 233 PetscFunctionBegin; 234 if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill; 235 236 if (dir->inplace) { 237 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 238 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 239 ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 240 if (dir->row) { 241 ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 242 ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 243 } 244 ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr); 245 dir->fact = pc->pmat; 246 } else { 247 MatInfo info; 248 if (!pc->setupcalled) { 249 ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 250 if (dir->nonzerosalongdiagonal) { 251 ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 252 } 253 if (dir->row) { 254 ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 255 ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 256 } 257 ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 258 ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 259 dir->actualfill = info.fill_ratio_needed; 260 ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr); 261 } else if (pc->flag != SAME_NONZERO_PATTERN) { 262 if (!dir->reuseordering) { 263 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 264 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 265 ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 266 if (dir->nonzerosalongdiagonal) { 267 ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 268 } 269 if (dir->row) { 270 ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 271 ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 272 } 273 } 274 ierr = MatDestroy(dir->fact);CHKERRQ(ierr); 275 ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 276 ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 277 dir->actualfill = info.fill_ratio_needed; 278 ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr); 279 } 280 ierr = MatLUFactorNumeric(pc->pmat,&dir->info,&dir->fact);CHKERRQ(ierr); 281 } 282 PetscFunctionReturn(0); 283 } 284 285 #undef __FUNCT__ 286 #define __FUNCT__ "PCDestroy_LU" 287 static PetscErrorCode PCDestroy_LU(PC pc) 288 { 289 PC_LU *dir = (PC_LU*)pc->data; 290 PetscErrorCode ierr; 291 292 PetscFunctionBegin; 293 if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);} 294 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 295 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 296 ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 297 ierr = PetscFree(dir);CHKERRQ(ierr); 298 PetscFunctionReturn(0); 299 } 300 301 #undef __FUNCT__ 302 #define __FUNCT__ "PCApply_LU" 303 static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 304 { 305 PC_LU *dir = (PC_LU*)pc->data; 306 PetscErrorCode ierr; 307 308 PetscFunctionBegin; 309 if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 310 else {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);} 311 PetscFunctionReturn(0); 312 } 313 314 #undef __FUNCT__ 315 #define __FUNCT__ "PCApplyTranspose_LU" 316 static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 317 { 318 PC_LU *dir = (PC_LU*)pc->data; 319 PetscErrorCode ierr; 320 321 PetscFunctionBegin; 322 if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 323 else {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);} 324 PetscFunctionReturn(0); 325 } 326 327 /* -----------------------------------------------------------------------------------*/ 328 329 EXTERN_C_BEGIN 330 #undef __FUNCT__ 331 #define __FUNCT__ "PCFactorSetFill_LU" 332 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_LU(PC pc,PetscReal fill) 333 { 334 PC_LU *dir; 335 336 PetscFunctionBegin; 337 dir = (PC_LU*)pc->data; 338 dir->info.fill = fill; 339 PetscFunctionReturn(0); 340 } 341 EXTERN_C_END 342 343 EXTERN_C_BEGIN 344 #undef __FUNCT__ 345 #define __FUNCT__ "PCFactorSetUseInPlace_LU" 346 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_LU(PC pc) 347 { 348 PC_LU *dir; 349 350 PetscFunctionBegin; 351 dir = (PC_LU*)pc->data; 352 dir->inplace = PETSC_TRUE; 353 PetscFunctionReturn(0); 354 } 355 EXTERN_C_END 356 357 EXTERN_C_BEGIN 358 #undef __FUNCT__ 359 #define __FUNCT__ "PCFactorSetMatOrderingType_LU" 360 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_LU(PC pc,MatOrderingType ordering) 361 { 362 PC_LU *dir = (PC_LU*)pc->data; 363 PetscErrorCode ierr; 364 365 PetscFunctionBegin; 366 ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 367 ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 368 PetscFunctionReturn(0); 369 } 370 EXTERN_C_END 371 372 EXTERN_C_BEGIN 373 #undef __FUNCT__ 374 #define __FUNCT__ "PCFactorSetPivoting_LU" 375 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting_LU(PC pc,PetscReal dtcol) 376 { 377 PC_LU *dir = (PC_LU*)pc->data; 378 379 PetscFunctionBegin; 380 if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %G must be between 0 and 1",dtcol); 381 dir->info.dtcol = dtcol; 382 PetscFunctionReturn(0); 383 } 384 EXTERN_C_END 385 386 EXTERN_C_BEGIN 387 #undef __FUNCT__ 388 #define __FUNCT__ "PCFactorSetPivotInBlocks_LU" 389 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_LU(PC pc,PetscTruth pivot) 390 { 391 PC_LU *dir = (PC_LU*)pc->data; 392 393 PetscFunctionBegin; 394 dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 395 PetscFunctionReturn(0); 396 } 397 EXTERN_C_END 398 399 /* -----------------------------------------------------------------------------------*/ 400 401 #undef __FUNCT__ 402 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal" 403 /*@ 404 PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 405 406 Collective on PC 407 408 Input Parameters: 409 + pc - the preconditioner context 410 - tol - diagonal entries smaller than this in absolute value are considered zero 411 412 Options Database Key: 413 . -pc_factor_nonzeros_along_diagonal 414 415 Level: intermediate 416 417 .keywords: PC, set, factorization, direct, fill 418 419 .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal() 420 @*/ 421 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 422 { 423 PetscErrorCode ierr,(*f)(PC,PetscReal); 424 425 PetscFunctionBegin; 426 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 427 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr); 428 if (f) { 429 ierr = (*f)(pc,rtol);CHKERRQ(ierr); 430 } 431 PetscFunctionReturn(0); 432 } 433 434 #undef __FUNCT__ 435 #define __FUNCT__ "PCFactorSetFill" 436 /*@ 437 PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 438 fill = number nonzeros in factor/number nonzeros in original matrix. 439 440 Collective on PC 441 442 Input Parameters: 443 + pc - the preconditioner context 444 - fill - amount of expected fill 445 446 Options Database Key: 447 . -pc_factor_fill <fill> - Sets fill amount 448 449 Level: intermediate 450 451 Note: 452 For sparse matrix factorizations it is difficult to predict how much 453 fill to expect. By running with the option -info PETSc will print the 454 actual amount of fill used; allowing you to set the value accurately for 455 future runs. Default PETSc uses a value of 5.0 456 457 .keywords: PC, set, factorization, direct, fill 458 459 @*/ 460 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill) 461 { 462 PetscErrorCode ierr,(*f)(PC,PetscReal); 463 464 PetscFunctionBegin; 465 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 466 if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0"); 467 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 468 if (f) { 469 ierr = (*f)(pc,fill);CHKERRQ(ierr); 470 } 471 PetscFunctionReturn(0); 472 } 473 474 #undef __FUNCT__ 475 #define __FUNCT__ "PCFactorSetUseInPlace" 476 /*@ 477 PCFactorSetUseInPlace - Tells the system to do an in-place factorization. 478 For dense matrices, this enables the solution of much larger problems. 479 For sparse matrices the factorization cannot be done truly in-place 480 so this does not save memory during the factorization, but after the matrix 481 is factored, the original unfactored matrix is freed, thus recovering that 482 space. 483 484 Collective on PC 485 486 Input Parameters: 487 . pc - the preconditioner context 488 489 Options Database Key: 490 . -pc_factor_in_place - Activates in-place factorization 491 492 Notes: 493 PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when 494 a different matrix is provided for the multiply and the preconditioner in 495 a call to KSPSetOperators(). 496 This is because the Krylov space methods require an application of the 497 matrix multiplication, which is not possible here because the matrix has 498 been factored in-place, replacing the original matrix. 499 500 Level: intermediate 501 502 .keywords: PC, set, factorization, direct, inplace, in-place, LU 503 504 .seealso: PCILUSetUseInPlace() 505 @*/ 506 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc) 507 { 508 PetscErrorCode ierr,(*f)(PC); 509 510 PetscFunctionBegin; 511 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 512 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr); 513 if (f) { 514 ierr = (*f)(pc);CHKERRQ(ierr); 515 } 516 PetscFunctionReturn(0); 517 } 518 519 #undef __FUNCT__ 520 #define __FUNCT__ "PCFactorSetMatOrderingType" 521 /*@C 522 PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 523 be used in the LU factorization. 524 525 Collective on PC 526 527 Input Parameters: 528 + pc - the preconditioner context 529 - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 530 531 Options Database Key: 532 . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 533 534 Level: intermediate 535 536 Notes: nested dissection is used by default 537 538 @*/ 539 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,MatOrderingType ordering) 540 { 541 PetscErrorCode ierr,(*f)(PC,MatOrderingType); 542 543 PetscFunctionBegin; 544 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr); 545 if (f) { 546 ierr = (*f)(pc,ordering);CHKERRQ(ierr); 547 } 548 PetscFunctionReturn(0); 549 } 550 551 #undef __FUNCT__ 552 #define __FUNCT__ "PCFactorSetPivoting" 553 /*@ 554 PCFactorSetPivoting - Determines when pivoting is done during LU. 555 For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 556 it is never done. For the Matlab and SuperLU factorization this is used. 557 558 Collective on PC 559 560 Input Parameters: 561 + pc - the preconditioner context 562 - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 563 564 Options Database Key: 565 . -pc_factor_pivoting <dtcol> 566 567 Level: intermediate 568 569 .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks() 570 @*/ 571 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC pc,PetscReal dtcol) 572 { 573 PetscErrorCode ierr,(*f)(PC,PetscReal); 574 575 PetscFunctionBegin; 576 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr); 577 if (f) { 578 ierr = (*f)(pc,dtcol);CHKERRQ(ierr); 579 } 580 PetscFunctionReturn(0); 581 } 582 583 #undef __FUNCT__ 584 #define __FUNCT__ "PCFactorSetPivotInBlocks" 585 /*@ 586 PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 587 with BAIJ or SBAIJ matrices 588 589 Collective on PC 590 591 Input Parameters: 592 + pc - the preconditioner context 593 - pivot - PETSC_TRUE or PETSC_FALSE 594 595 Options Database Key: 596 . -pc_factor_pivot_in_blocks <true,false> 597 598 Level: intermediate 599 600 .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting() 601 @*/ 602 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot) 603 { 604 PetscErrorCode ierr,(*f)(PC,PetscTruth); 605 606 PetscFunctionBegin; 607 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 608 if (f) { 609 ierr = (*f)(pc,pivot);CHKERRQ(ierr); 610 } 611 PetscFunctionReturn(0); 612 } 613 614 /* ------------------------------------------------------------------------ */ 615 616 /*MC 617 PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 618 619 Options Database Keys: 620 + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 621 . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 622 . -pc_factor_fill <fill> - Sets fill amount 623 . -pc_factor_in_place - Activates in-place factorization 624 . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 625 . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 626 stability of factorization. 627 . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 628 . -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 629 is optional with PETSC_TRUE being the default 630 - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 631 diagonal. 632 633 Notes: Not all options work for all matrix formats 634 Run with -help to see additional options for particular matrix formats or factorization 635 algorithms 636 637 Level: beginner 638 639 Concepts: LU factorization, direct solver 640 641 Notes: Usually this will compute an "exact" solution in one iteration and does 642 not need a Krylov method (i.e. you can use -ksp_type preonly, or 643 KSPSetType(ksp,KSPPREONLY) for the Krylov method 644 645 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 646 PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 647 PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetPivoting(), 648 PCFactorSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(), PCFactorReorderForNonzeroDiagonal() 649 M*/ 650 651 EXTERN_C_BEGIN 652 #undef __FUNCT__ 653 #define __FUNCT__ "PCCreate_LU" 654 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc) 655 { 656 PetscErrorCode ierr; 657 PetscMPIInt size; 658 PC_LU *dir; 659 660 PetscFunctionBegin; 661 ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr); 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 = 0.0; /* 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 = PCFactorGetMatrix_LU; 695 696 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_LU", 697 PCFactorSetZeroPivot_LU);CHKERRQ(ierr); 698 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_LU", 699 PCFactorSetShiftNonzero_LU);CHKERRQ(ierr); 700 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_LU", 701 PCFactorSetShiftPd_LU);CHKERRQ(ierr); 702 703 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_LU", 704 PCFactorSetFill_LU);CHKERRQ(ierr); 705 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU", 706 PCFactorSetUseInPlace_LU);CHKERRQ(ierr); 707 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_LU", 708 PCFactorSetMatOrderingType_LU);CHKERRQ(ierr); 709 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU", 710 PCFactorSetReuseOrdering_LU);CHKERRQ(ierr); 711 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU", 712 PCFactorSetReuseFill_LU);CHKERRQ(ierr); 713 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivoting_C","PCFactorSetPivoting_LU", 714 PCFactorSetPivoting_LU);CHKERRQ(ierr); 715 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_LU", 716 PCFactorSetPivotInBlocks_LU);CHKERRQ(ierr); 717 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU", 718 PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 719 PetscFunctionReturn(0); 720 } 721 EXTERN_C_END 722