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],solvertype[64]; 117 MatSolverType stype; 118 PetscFList ordlist; 119 PetscReal tol; 120 121 PetscFunctionBegin; 122 ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 123 ierr = PetscOptionsHead("LU options");CHKERRQ(ierr); 124 ierr = PetscOptionsName("-pc_factor_in_place","Form LU in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr); 125 if (flg) { 126 ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr); 127 } 128 ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in LU/non-zeros in matrix","PCFactorSetFill",lu->info.fill,&lu->info.fill,0);CHKERRQ(ierr); 129 130 ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 131 if (flg) { 132 ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr); 133 } 134 ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",lu->info.shiftnz,&lu->info.shiftnz,0);CHKERRQ(ierr); 135 ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr); 136 if (flg) { 137 ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr); 138 } 139 ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr); 140 141 ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr); 142 if (flg) { 143 ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr); 144 } 145 ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr); 146 if (flg) { 147 ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr); 148 } 149 150 ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 151 ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in LU","PCFactorSetMatOrderingType",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr); 152 if (flg) { 153 ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 154 } 155 156 ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 157 if (flg) { 158 tol = PETSC_DECIDE; 159 ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 160 ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 161 } 162 163 ierr = PetscOptionsReal("-pc_factor_pivoting","Pivoting tolerance (used only for some factorization)","PCFactorSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);CHKERRQ(ierr); 164 165 flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 166 ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 167 if (set) { 168 ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 169 } 170 ierr = PetscOptionsTail();CHKERRQ(ierr); 171 PetscFunctionReturn(0); 172 } 173 174 #undef __FUNCT__ 175 #define __FUNCT__ "PCView_LU" 176 static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer) 177 { 178 PC_LU *lu = (PC_LU*)pc->data; 179 PetscErrorCode ierr; 180 PetscTruth iascii,isstring; 181 182 PetscFunctionBegin; 183 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 184 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 185 if (iascii) { 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 = PetscViewerASCIIPrintf(viewer," LU: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr); 194 ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 195 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 196 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 197 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 198 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 199 ierr = MatView(lu->fact,viewer);CHKERRQ(ierr); 200 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 201 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 202 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 203 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 204 } 205 if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 206 if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 207 } else if (isstring) { 208 ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 209 } else { 210 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name); 211 } 212 PetscFunctionReturn(0); 213 } 214 215 #undef __FUNCT__ 216 #define __FUNCT__ "PCFactorGetMatrix_LU" 217 static PetscErrorCode PCFactorGetMatrix_LU(PC pc,Mat *mat) 218 { 219 PC_LU *dir = (PC_LU*)pc->data; 220 221 PetscFunctionBegin; 222 if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 223 *mat = dir->fact; 224 PetscFunctionReturn(0); 225 } 226 227 #undef __FUNCT__ 228 #define __FUNCT__ "PCSetUp_LU" 229 static PetscErrorCode PCSetUp_LU(PC pc) 230 { 231 PetscErrorCode ierr; 232 PC_LU *dir = (PC_LU*)pc->data; 233 234 PetscFunctionBegin; 235 if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill; 236 237 if (dir->inplace) { 238 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 239 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 240 ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 241 if (dir->row) { 242 ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 243 ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 244 } 245 ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr); 246 dir->fact = pc->pmat; 247 } else { 248 MatInfo info; 249 if (!pc->setupcalled) { 250 ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 251 if (dir->nonzerosalongdiagonal) { 252 ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 253 } 254 if (dir->row) { 255 ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 256 ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 257 } 258 ierr = MatGetFactor(pc->pmat,"petsc",MAT_FACTOR_LU,&dir->fact);CHKERRQ(ierr); 259 ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 260 ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 261 dir->actualfill = info.fill_ratio_needed; 262 ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr); 263 } else if (pc->flag != SAME_NONZERO_PATTERN) { 264 if (!dir->reuseordering) { 265 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 266 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 267 ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 268 if (dir->nonzerosalongdiagonal) { 269 ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 270 } 271 if (dir->row) { 272 ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 273 ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 274 } 275 } 276 ierr = MatDestroy(dir->fact);CHKERRQ(ierr); 277 ierr = MatGetFactor(pc->pmat,"petsc",MAT_FACTOR_LU,&dir->fact);CHKERRQ(ierr); 278 ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 279 ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 280 dir->actualfill = info.fill_ratio_needed; 281 ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr); 282 } 283 ierr = MatLUFactorNumeric(pc->pmat,&dir->info,&dir->fact);CHKERRQ(ierr); 284 } 285 PetscFunctionReturn(0); 286 } 287 288 #undef __FUNCT__ 289 #define __FUNCT__ "PCDestroy_LU" 290 static PetscErrorCode PCDestroy_LU(PC pc) 291 { 292 PC_LU *dir = (PC_LU*)pc->data; 293 PetscErrorCode ierr; 294 295 PetscFunctionBegin; 296 if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);} 297 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 298 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 299 ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 300 ierr = PetscFree(dir);CHKERRQ(ierr); 301 PetscFunctionReturn(0); 302 } 303 304 #undef __FUNCT__ 305 #define __FUNCT__ "PCApply_LU" 306 static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 307 { 308 PC_LU *dir = (PC_LU*)pc->data; 309 PetscErrorCode ierr; 310 311 PetscFunctionBegin; 312 if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 313 else {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);} 314 PetscFunctionReturn(0); 315 } 316 317 #undef __FUNCT__ 318 #define __FUNCT__ "PCApplyTranspose_LU" 319 static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 320 { 321 PC_LU *dir = (PC_LU*)pc->data; 322 PetscErrorCode ierr; 323 324 PetscFunctionBegin; 325 if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 326 else {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);} 327 PetscFunctionReturn(0); 328 } 329 330 /* -----------------------------------------------------------------------------------*/ 331 332 EXTERN_C_BEGIN 333 #undef __FUNCT__ 334 #define __FUNCT__ "PCFactorSetFill_LU" 335 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_LU(PC pc,PetscReal fill) 336 { 337 PC_LU *dir; 338 339 PetscFunctionBegin; 340 dir = (PC_LU*)pc->data; 341 dir->info.fill = fill; 342 PetscFunctionReturn(0); 343 } 344 EXTERN_C_END 345 346 EXTERN_C_BEGIN 347 #undef __FUNCT__ 348 #define __FUNCT__ "PCFactorSetUseInPlace_LU" 349 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_LU(PC pc) 350 { 351 PC_LU *dir; 352 353 PetscFunctionBegin; 354 dir = (PC_LU*)pc->data; 355 dir->inplace = PETSC_TRUE; 356 PetscFunctionReturn(0); 357 } 358 EXTERN_C_END 359 360 EXTERN_C_BEGIN 361 #undef __FUNCT__ 362 #define __FUNCT__ "PCFactorSetMatOrderingType_LU" 363 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_LU(PC pc,MatOrderingType ordering) 364 { 365 PC_LU *dir = (PC_LU*)pc->data; 366 PetscErrorCode ierr; 367 368 PetscFunctionBegin; 369 ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 370 ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 371 PetscFunctionReturn(0); 372 } 373 EXTERN_C_END 374 375 EXTERN_C_BEGIN 376 #undef __FUNCT__ 377 #define __FUNCT__ "PCFactorSetPivoting_LU" 378 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting_LU(PC pc,PetscReal dtcol) 379 { 380 PC_LU *dir = (PC_LU*)pc->data; 381 382 PetscFunctionBegin; 383 if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %G must be between 0 and 1",dtcol); 384 dir->info.dtcol = dtcol; 385 PetscFunctionReturn(0); 386 } 387 EXTERN_C_END 388 389 EXTERN_C_BEGIN 390 #undef __FUNCT__ 391 #define __FUNCT__ "PCFactorSetPivotInBlocks_LU" 392 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_LU(PC pc,PetscTruth pivot) 393 { 394 PC_LU *dir = (PC_LU*)pc->data; 395 396 PetscFunctionBegin; 397 dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 398 PetscFunctionReturn(0); 399 } 400 EXTERN_C_END 401 402 /* -----------------------------------------------------------------------------------*/ 403 404 #undef __FUNCT__ 405 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal" 406 /*@ 407 PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 408 409 Collective on PC 410 411 Input Parameters: 412 + pc - the preconditioner context 413 - tol - diagonal entries smaller than this in absolute value are considered zero 414 415 Options Database Key: 416 . -pc_factor_nonzeros_along_diagonal 417 418 Level: intermediate 419 420 .keywords: PC, set, factorization, direct, fill 421 422 .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal() 423 @*/ 424 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 425 { 426 PetscErrorCode ierr,(*f)(PC,PetscReal); 427 428 PetscFunctionBegin; 429 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 430 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr); 431 if (f) { 432 ierr = (*f)(pc,rtol);CHKERRQ(ierr); 433 } 434 PetscFunctionReturn(0); 435 } 436 437 #undef __FUNCT__ 438 #define __FUNCT__ "PCFactorSetFill" 439 /*@ 440 PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 441 fill = number nonzeros in factor/number nonzeros in original matrix. 442 443 Collective on PC 444 445 Input Parameters: 446 + pc - the preconditioner context 447 - fill - amount of expected fill 448 449 Options Database Key: 450 . -pc_factor_fill <fill> - Sets fill amount 451 452 Level: intermediate 453 454 Note: 455 For sparse matrix factorizations it is difficult to predict how much 456 fill to expect. By running with the option -info PETSc will print the 457 actual amount of fill used; allowing you to set the value accurately for 458 future runs. Default PETSc uses a value of 5.0 459 460 .keywords: PC, set, factorization, direct, fill 461 462 @*/ 463 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(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,"PCFactorSetFill_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__ "PCFactorSetUseInPlace" 479 /*@ 480 PCFactorSetUseInPlace - 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_factor_in_place - Activates in-place factorization 494 495 Notes: 496 PCFactorSetUseInplace() 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 PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc) 510 { 511 PetscErrorCode ierr,(*f)(PC); 512 513 PetscFunctionBegin; 514 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 515 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_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__ "PCFactorSetMatOrderingType" 524 /*@C 525 PCFactorSetMatOrderingType - 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_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 536 537 Level: intermediate 538 539 Notes: nested dissection is used by default 540 541 For Cholesky and ICC and the SBAIJ format reorderings are not available, 542 since only the upper triangular part of the matrix is stored. You can use the 543 SeqAIJ format in this case to get reorderings. 544 545 @*/ 546 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,MatOrderingType ordering) 547 { 548 PetscErrorCode ierr,(*f)(PC,MatOrderingType); 549 550 PetscFunctionBegin; 551 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr); 552 if (f) { 553 ierr = (*f)(pc,ordering);CHKERRQ(ierr); 554 } 555 PetscFunctionReturn(0); 556 } 557 558 #undef __FUNCT__ 559 #define __FUNCT__ "PCFactorSetPivoting" 560 /*@ 561 PCFactorSetPivoting - Determines when pivoting is done during LU. 562 For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 563 it is never done. For the Matlab and SuperLU factorization this is used. 564 565 Collective on PC 566 567 Input Parameters: 568 + pc - the preconditioner context 569 - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 570 571 Options Database Key: 572 . -pc_factor_pivoting <dtcol> 573 574 Level: intermediate 575 576 .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks() 577 @*/ 578 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC pc,PetscReal dtcol) 579 { 580 PetscErrorCode ierr,(*f)(PC,PetscReal); 581 582 PetscFunctionBegin; 583 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr); 584 if (f) { 585 ierr = (*f)(pc,dtcol);CHKERRQ(ierr); 586 } 587 PetscFunctionReturn(0); 588 } 589 590 #undef __FUNCT__ 591 #define __FUNCT__ "PCFactorSetPivotInBlocks" 592 /*@ 593 PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 594 with BAIJ or SBAIJ matrices 595 596 Collective on PC 597 598 Input Parameters: 599 + pc - the preconditioner context 600 - pivot - PETSC_TRUE or PETSC_FALSE 601 602 Options Database Key: 603 . -pc_factor_pivot_in_blocks <true,false> 604 605 Level: intermediate 606 607 .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting() 608 @*/ 609 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot) 610 { 611 PetscErrorCode ierr,(*f)(PC,PetscTruth); 612 613 PetscFunctionBegin; 614 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 615 if (f) { 616 ierr = (*f)(pc,pivot);CHKERRQ(ierr); 617 } 618 PetscFunctionReturn(0); 619 } 620 621 /* ------------------------------------------------------------------------ */ 622 623 /*MC 624 PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 625 626 Options Database Keys: 627 + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 628 . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 629 . -pc_factor_fill <fill> - Sets fill amount 630 . -pc_factor_in_place - Activates in-place factorization 631 . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 632 . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 633 stability of factorization. 634 . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 635 . -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 636 is optional with PETSC_TRUE being the default 637 - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 638 diagonal. 639 640 Notes: Not all options work for all matrix formats 641 Run with -help to see additional options for particular matrix formats or factorization 642 algorithms 643 644 Level: beginner 645 646 Concepts: LU factorization, direct solver 647 648 Notes: Usually this will compute an "exact" solution in one iteration and does 649 not need a Krylov method (i.e. you can use -ksp_type preonly, or 650 KSPSetType(ksp,KSPPREONLY) for the Krylov method 651 652 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 653 PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 654 PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetPivoting(), 655 PCFactorSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(), PCFactorReorderForNonzeroDiagonal() 656 M*/ 657 658 EXTERN_C_BEGIN 659 #undef __FUNCT__ 660 #define __FUNCT__ "PCCreate_LU" 661 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc) 662 { 663 PetscErrorCode ierr; 664 PetscMPIInt size; 665 PC_LU *dir; 666 667 PetscFunctionBegin; 668 ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr); 669 670 ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr); 671 dir->fact = 0; 672 dir->inplace = PETSC_FALSE; 673 dir->nonzerosalongdiagonal = PETSC_FALSE; 674 675 dir->info.fill = 5.0; 676 dir->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 677 dir->info.shiftnz = 0.0; 678 dir->info.zeropivot = 1.e-12; 679 dir->info.pivotinblocks = 1.0; 680 dir->info.shiftpd = 0.0; /* false */ 681 dir->info.shift_fraction = 0.0; 682 dir->col = 0; 683 dir->row = 0; 684 ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 685 if (size == 1) { 686 ierr = PetscStrallocpy(MATORDERING_ND,&dir->ordering);CHKERRQ(ierr); 687 } else { 688 ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr); 689 } 690 dir->reusefill = PETSC_FALSE; 691 dir->reuseordering = PETSC_FALSE; 692 pc->data = (void*)dir; 693 694 pc->ops->destroy = PCDestroy_LU; 695 pc->ops->apply = PCApply_LU; 696 pc->ops->applytranspose = PCApplyTranspose_LU; 697 pc->ops->setup = PCSetUp_LU; 698 pc->ops->setfromoptions = PCSetFromOptions_LU; 699 pc->ops->view = PCView_LU; 700 pc->ops->applyrichardson = 0; 701 pc->ops->getfactoredmatrix = PCFactorGetMatrix_LU; 702 703 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_LU", 704 PCFactorSetZeroPivot_LU);CHKERRQ(ierr); 705 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_LU", 706 PCFactorSetShiftNonzero_LU);CHKERRQ(ierr); 707 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_LU", 708 PCFactorSetShiftPd_LU);CHKERRQ(ierr); 709 710 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_LU", 711 PCFactorSetFill_LU);CHKERRQ(ierr); 712 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU", 713 PCFactorSetUseInPlace_LU);CHKERRQ(ierr); 714 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_LU", 715 PCFactorSetMatOrderingType_LU);CHKERRQ(ierr); 716 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU", 717 PCFactorSetReuseOrdering_LU);CHKERRQ(ierr); 718 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU", 719 PCFactorSetReuseFill_LU);CHKERRQ(ierr); 720 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivoting_C","PCFactorSetPivoting_LU", 721 PCFactorSetPivoting_LU);CHKERRQ(ierr); 722 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_LU", 723 PCFactorSetPivotInBlocks_LU);CHKERRQ(ierr); 724 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU", 725 PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 726 PetscFunctionReturn(0); 727 } 728 EXTERN_C_END 729