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