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","PCFactorSetMatOrdering",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr); 151 if (flg) { 152 ierr = PCFactorSetMatOrdering(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 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__ "PCFactorSetUseInPlace_LU" 341 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_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__ "PCFactorSetMatOrdering_LU" 355 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrdering_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__ "PCFactorSetPivoting_LU" 370 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting_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__ "PCFactorSetPivotInBlocks_LU" 384 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_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__ "PCFactorReorderForNonzeroDiagonal" 398 /*@ 399 PCFactorReorderForNonzeroDiagonal - 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_factor_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 PCFactorReorderForNonzeroDiagonal(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,"PCFactorReorderForNonzeroDiagonal_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__ "PCFactorSetFill" 431 /*@ 432 PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 433 fill = number nonzeros in factor/number nonzeros in original matrix. 434 435 Collective on PC 436 437 Input Parameters: 438 + pc - the preconditioner context 439 - fill - amount of expected fill 440 441 Options Database Key: 442 . -pc_factor_fill <fill> - Sets fill amount 443 444 Level: intermediate 445 446 Note: 447 For sparse matrix factorizations it is difficult to predict how much 448 fill to expect. By running with the option -info PETSc will print the 449 actual amount of fill used; allowing you to set the value accurately for 450 future runs. Default PETSc uses a value of 5.0 451 452 .keywords: PC, set, factorization, direct, fill 453 454 @*/ 455 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill) 456 { 457 PetscErrorCode ierr,(*f)(PC,PetscReal); 458 459 PetscFunctionBegin; 460 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 461 if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0"); 462 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 463 if (f) { 464 ierr = (*f)(pc,fill);CHKERRQ(ierr); 465 } 466 PetscFunctionReturn(0); 467 } 468 469 #undef __FUNCT__ 470 #define __FUNCT__ "PCFactorSetUseInPlace" 471 /*@ 472 PCFactorSetUseInPlace - Tells the system to do an in-place factorization. 473 For dense matrices, this enables the solution of much larger problems. 474 For sparse matrices the factorization cannot be done truly in-place 475 so this does not save memory during the factorization, but after the matrix 476 is factored, the original unfactored matrix is freed, thus recovering that 477 space. 478 479 Collective on PC 480 481 Input Parameters: 482 . pc - the preconditioner context 483 484 Options Database Key: 485 . -pc_factor_in_place - Activates in-place factorization 486 487 Notes: 488 PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when 489 a different matrix is provided for the multiply and the preconditioner in 490 a call to KSPSetOperators(). 491 This is because the Krylov space methods require an application of the 492 matrix multiplication, which is not possible here because the matrix has 493 been factored in-place, replacing the original matrix. 494 495 Level: intermediate 496 497 .keywords: PC, set, factorization, direct, inplace, in-place, LU 498 499 .seealso: PCILUSetUseInPlace() 500 @*/ 501 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc) 502 { 503 PetscErrorCode ierr,(*f)(PC); 504 505 PetscFunctionBegin; 506 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 507 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr); 508 if (f) { 509 ierr = (*f)(pc);CHKERRQ(ierr); 510 } 511 PetscFunctionReturn(0); 512 } 513 514 #undef __FUNCT__ 515 #define __FUNCT__ "PCFactorSetMatOrdering" 516 /*@C 517 PCFactorSetMatOrdering - Sets the ordering routine (to reduce fill) to 518 be used in the LU factorization. 519 520 Collective on PC 521 522 Input Parameters: 523 + pc - the preconditioner context 524 - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 525 526 Options Database Key: 527 . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 528 529 Level: intermediate 530 531 Notes: nested dissection is used by default 532 533 .seealso: PCILUSetMatOrdering() 534 @*/ 535 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrdering(PC pc,MatOrderingType ordering) 536 { 537 PetscErrorCode ierr,(*f)(PC,MatOrderingType); 538 539 PetscFunctionBegin; 540 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 541 if (f) { 542 ierr = (*f)(pc,ordering);CHKERRQ(ierr); 543 } 544 PetscFunctionReturn(0); 545 } 546 547 #undef __FUNCT__ 548 #define __FUNCT__ "PCFactorSetPivoting" 549 /*@ 550 PCFactorSetPivoting - Determines when pivoting is done during LU. 551 For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 552 it is never done. For the Matlab and SuperLU factorization this is used. 553 554 Collective on PC 555 556 Input Parameters: 557 + pc - the preconditioner context 558 - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 559 560 Options Database Key: 561 . -pc_factor_pivoting <dtcol> 562 563 Level: intermediate 564 565 .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks() 566 @*/ 567 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC pc,PetscReal dtcol) 568 { 569 PetscErrorCode ierr,(*f)(PC,PetscReal); 570 571 PetscFunctionBegin; 572 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr); 573 if (f) { 574 ierr = (*f)(pc,dtcol);CHKERRQ(ierr); 575 } 576 PetscFunctionReturn(0); 577 } 578 579 #undef __FUNCT__ 580 #define __FUNCT__ "PCFactorSetPivotInBlocks" 581 /*@ 582 PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 583 with BAIJ or SBAIJ matrices 584 585 Collective on PC 586 587 Input Parameters: 588 + pc - the preconditioner context 589 - pivot - PETSC_TRUE or PETSC_FALSE 590 591 Options Database Key: 592 . -pc_factor_pivot_in_blocks <true,false> 593 594 Level: intermediate 595 596 .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting() 597 @*/ 598 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot) 599 { 600 PetscErrorCode ierr,(*f)(PC,PetscTruth); 601 602 PetscFunctionBegin; 603 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 604 if (f) { 605 ierr = (*f)(pc,pivot);CHKERRQ(ierr); 606 } 607 PetscFunctionReturn(0); 608 } 609 610 /* ------------------------------------------------------------------------ */ 611 612 /*MC 613 PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 614 615 Options Database Keys: 616 + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 617 . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 618 . -pc_factor_fill <fill> - Sets fill amount 619 . -pc_factor_in_place - Activates in-place factorization 620 . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 621 . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 622 stability of factorization. 623 . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 624 - -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 625 is optional with PETSC_TRUE being the default 626 627 Notes: Not all options work for all matrix formats 628 Run with -help to see additional options for particular matrix formats or factorization 629 algorithms 630 631 Level: beginner 632 633 Concepts: LU factorization, direct solver 634 635 Notes: Usually this will compute an "exact" solution in one iteration and does 636 not need a Krylov method (i.e. you can use -ksp_type preonly, or 637 KSPSetType(ksp,KSPPREONLY) for the Krylov method 638 639 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 640 PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCGetFactoredMatrix(), 641 PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrdering(), PCFactorSetPivoting(), 642 PCFactorSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd() 643 M*/ 644 645 EXTERN_C_BEGIN 646 #undef __FUNCT__ 647 #define __FUNCT__ "PCCreate_LU" 648 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc) 649 { 650 PetscErrorCode ierr; 651 PetscMPIInt size; 652 PC_LU *dir; 653 654 PetscFunctionBegin; 655 ierr = PetscNew(PC_LU,&dir);CHKERRQ(ierr); 656 ierr = PetscLogObjectMemory(pc,sizeof(PC_LU));CHKERRQ(ierr); 657 658 ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr); 659 dir->fact = 0; 660 dir->inplace = PETSC_FALSE; 661 dir->nonzerosalongdiagonal = PETSC_FALSE; 662 663 dir->info.fill = 5.0; 664 dir->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 665 dir->info.shiftnz = 0.0; 666 dir->info.zeropivot = 1.e-12; 667 dir->info.pivotinblocks = 1.0; 668 dir->info.shiftpd = 0.0; /* false */ 669 dir->info.shift_fraction = 0.0; 670 dir->col = 0; 671 dir->row = 0; 672 ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr); 673 if (size == 1) { 674 ierr = PetscStrallocpy(MATORDERING_ND,&dir->ordering);CHKERRQ(ierr); 675 } else { 676 ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr); 677 } 678 dir->reusefill = PETSC_FALSE; 679 dir->reuseordering = PETSC_FALSE; 680 pc->data = (void*)dir; 681 682 pc->ops->destroy = PCDestroy_LU; 683 pc->ops->apply = PCApply_LU; 684 pc->ops->applytranspose = PCApplyTranspose_LU; 685 pc->ops->setup = PCSetUp_LU; 686 pc->ops->setfromoptions = PCSetFromOptions_LU; 687 pc->ops->view = PCView_LU; 688 pc->ops->applyrichardson = 0; 689 pc->ops->getfactoredmatrix = PCGetFactoredMatrix_LU; 690 691 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_LU", 692 PCFactorSetZeroPivot_LU);CHKERRQ(ierr); 693 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_LU", 694 PCFactorSetShiftNonzero_LU);CHKERRQ(ierr); 695 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_LU", 696 PCFactorSetShiftPd_LU);CHKERRQ(ierr); 697 698 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_LU", 699 PCFactorSetFill_LU);CHKERRQ(ierr); 700 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU", 701 PCFactorSetUseInPlace_LU);CHKERRQ(ierr); 702 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrdering_C","PCFactorSetMatOrdering_LU", 703 PCFactorSetMatOrdering_LU);CHKERRQ(ierr); 704 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU", 705 PCFactorSetReuseOrdering_LU);CHKERRQ(ierr); 706 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU", 707 PCFactorSetReuseFill_LU);CHKERRQ(ierr); 708 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivoting_C","PCFactorSetPivoting_LU", 709 PCFactorSetPivoting_LU);CHKERRQ(ierr); 710 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_LU", 711 PCFactorSetPivotInBlocks_LU);CHKERRQ(ierr); 712 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU", 713 PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 714 PetscFunctionReturn(0); 715 } 716 EXTERN_C_END 717