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