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 PetscFList ordlist; 118 PetscReal tol; 119 120 PetscFunctionBegin; 121 ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 122 ierr = PetscOptionsHead("LU options");CHKERRQ(ierr); 123 ierr = PetscOptionsName("-pc_factor_in_place","Form LU in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr); 124 if (flg) { 125 ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr); 126 } 127 ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in LU/non-zeros in matrix","PCFactorSetFill",lu->info.fill,&lu->info.fill,0);CHKERRQ(ierr); 128 129 ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 130 if (flg) { 131 ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr); 132 } 133 ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",lu->info.shiftnz,&lu->info.shiftnz,0);CHKERRQ(ierr); 134 ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr); 135 if (flg) { 136 ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr); 137 } 138 ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr); 139 140 ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr); 141 if (flg) { 142 ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr); 143 } 144 ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr); 145 if (flg) { 146 ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr); 147 } 148 149 ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 150 ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in LU","PCFactorSetMatOrderingType",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr); 151 if (flg) { 152 ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 153 } 154 155 /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */ 156 ierr = PetscOptionsString("-pc_mat_solver_type","Specific LU solver to use","MatGetFactor",lu->solvertype,solvertype,64,&flg);CHKERRQ(ierr); 157 if (flg) { 158 ierr = PetscStrallocpy(solvertype,&lu->solvertype);CHKERRQ(ierr); 159 } 160 161 ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 162 if (flg) { 163 tol = PETSC_DECIDE; 164 ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 165 ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 166 } 167 168 ierr = PetscOptionsReal("-pc_factor_pivoting","Pivoting tolerance (used only for some factorization)","PCFactorSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);CHKERRQ(ierr); 169 170 flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 171 ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 172 if (set) { 173 ierr = PCFactorSetPivotInBlocks(pc,flg);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 = MatGetFactor(pc->pmat,dir->solvertype,MAT_FACTOR_LU,&dir->fact);CHKERRQ(ierr); 264 ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 265 ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 266 dir->actualfill = info.fill_ratio_needed; 267 ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr); 268 } else if (pc->flag != SAME_NONZERO_PATTERN) { 269 if (!dir->reuseordering) { 270 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 271 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 272 ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 273 if (dir->nonzerosalongdiagonal) { 274 ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 275 } 276 if (dir->row) { 277 ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 278 ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 279 } 280 } 281 ierr = MatDestroy(dir->fact);CHKERRQ(ierr); 282 ierr = MatGetFactor(pc->pmat,dir->solvertype,MAT_FACTOR_LU,&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 = PetscStrfree(dir->solvertype);CHKERRQ(ierr); 306 ierr = PetscFree(dir);CHKERRQ(ierr); 307 PetscFunctionReturn(0); 308 } 309 310 #undef __FUNCT__ 311 #define __FUNCT__ "PCApply_LU" 312 static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 313 { 314 PC_LU *dir = (PC_LU*)pc->data; 315 PetscErrorCode ierr; 316 317 PetscFunctionBegin; 318 if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 319 else {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);} 320 PetscFunctionReturn(0); 321 } 322 323 #undef __FUNCT__ 324 #define __FUNCT__ "PCApplyTranspose_LU" 325 static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 326 { 327 PC_LU *dir = (PC_LU*)pc->data; 328 PetscErrorCode ierr; 329 330 PetscFunctionBegin; 331 if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 332 else {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);} 333 PetscFunctionReturn(0); 334 } 335 336 /* -----------------------------------------------------------------------------------*/ 337 338 EXTERN_C_BEGIN 339 #undef __FUNCT__ 340 #define __FUNCT__ "PCFactorSetFill_LU" 341 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_LU(PC pc,PetscReal fill) 342 { 343 PC_LU *dir; 344 345 PetscFunctionBegin; 346 dir = (PC_LU*)pc->data; 347 dir->info.fill = fill; 348 PetscFunctionReturn(0); 349 } 350 EXTERN_C_END 351 352 EXTERN_C_BEGIN 353 #undef __FUNCT__ 354 #define __FUNCT__ "PCFactorSetUseInPlace_LU" 355 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_LU(PC pc) 356 { 357 PC_LU *dir; 358 359 PetscFunctionBegin; 360 dir = (PC_LU*)pc->data; 361 dir->inplace = PETSC_TRUE; 362 PetscFunctionReturn(0); 363 } 364 EXTERN_C_END 365 366 EXTERN_C_BEGIN 367 #undef __FUNCT__ 368 #define __FUNCT__ "PCFactorSetMatOrderingType_LU" 369 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_LU(PC pc,MatOrderingType ordering) 370 { 371 PC_LU *dir = (PC_LU*)pc->data; 372 PetscErrorCode ierr; 373 374 PetscFunctionBegin; 375 ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 376 ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 377 PetscFunctionReturn(0); 378 } 379 EXTERN_C_END 380 381 EXTERN_C_BEGIN 382 #undef __FUNCT__ 383 #define __FUNCT__ "PCFactorSetPivoting_LU" 384 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting_LU(PC pc,PetscReal dtcol) 385 { 386 PC_LU *dir = (PC_LU*)pc->data; 387 388 PetscFunctionBegin; 389 if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %G must be between 0 and 1",dtcol); 390 dir->info.dtcol = dtcol; 391 PetscFunctionReturn(0); 392 } 393 EXTERN_C_END 394 395 EXTERN_C_BEGIN 396 #undef __FUNCT__ 397 #define __FUNCT__ "PCFactorSetPivotInBlocks_LU" 398 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_LU(PC pc,PetscTruth pivot) 399 { 400 PC_LU *dir = (PC_LU*)pc->data; 401 402 PetscFunctionBegin; 403 dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 404 PetscFunctionReturn(0); 405 } 406 EXTERN_C_END 407 408 /* -----------------------------------------------------------------------------------*/ 409 410 #undef __FUNCT__ 411 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal" 412 /*@ 413 PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 414 415 Collective on PC 416 417 Input Parameters: 418 + pc - the preconditioner context 419 - tol - diagonal entries smaller than this in absolute value are considered zero 420 421 Options Database Key: 422 . -pc_factor_nonzeros_along_diagonal 423 424 Level: intermediate 425 426 .keywords: PC, set, factorization, direct, fill 427 428 .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal() 429 @*/ 430 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 431 { 432 PetscErrorCode ierr,(*f)(PC,PetscReal); 433 434 PetscFunctionBegin; 435 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 436 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr); 437 if (f) { 438 ierr = (*f)(pc,rtol);CHKERRQ(ierr); 439 } 440 PetscFunctionReturn(0); 441 } 442 443 #undef __FUNCT__ 444 #define __FUNCT__ "PCFactorSetFill" 445 /*@ 446 PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 447 fill = number nonzeros in factor/number nonzeros in original matrix. 448 449 Collective on PC 450 451 Input Parameters: 452 + pc - the preconditioner context 453 - fill - amount of expected fill 454 455 Options Database Key: 456 . -pc_factor_fill <fill> - Sets fill amount 457 458 Level: intermediate 459 460 Note: 461 For sparse matrix factorizations it is difficult to predict how much 462 fill to expect. By running with the option -info PETSc will print the 463 actual amount of fill used; allowing you to set the value accurately for 464 future runs. Default PETSc uses a value of 5.0 465 466 .keywords: PC, set, factorization, direct, fill 467 468 @*/ 469 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill) 470 { 471 PetscErrorCode ierr,(*f)(PC,PetscReal); 472 473 PetscFunctionBegin; 474 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 475 if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0"); 476 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 477 if (f) { 478 ierr = (*f)(pc,fill);CHKERRQ(ierr); 479 } 480 PetscFunctionReturn(0); 481 } 482 483 #undef __FUNCT__ 484 #define __FUNCT__ "PCFactorSetUseInPlace" 485 /*@ 486 PCFactorSetUseInPlace - Tells the system to do an in-place factorization. 487 For dense matrices, this enables the solution of much larger problems. 488 For sparse matrices the factorization cannot be done truly in-place 489 so this does not save memory during the factorization, but after the matrix 490 is factored, the original unfactored matrix is freed, thus recovering that 491 space. 492 493 Collective on PC 494 495 Input Parameters: 496 . pc - the preconditioner context 497 498 Options Database Key: 499 . -pc_factor_in_place - Activates in-place factorization 500 501 Notes: 502 PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when 503 a different matrix is provided for the multiply and the preconditioner in 504 a call to KSPSetOperators(). 505 This is because the Krylov space methods require an application of the 506 matrix multiplication, which is not possible here because the matrix has 507 been factored in-place, replacing the original matrix. 508 509 Level: intermediate 510 511 .keywords: PC, set, factorization, direct, inplace, in-place, LU 512 513 .seealso: PCILUSetUseInPlace() 514 @*/ 515 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc) 516 { 517 PetscErrorCode ierr,(*f)(PC); 518 519 PetscFunctionBegin; 520 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 521 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr); 522 if (f) { 523 ierr = (*f)(pc);CHKERRQ(ierr); 524 } 525 PetscFunctionReturn(0); 526 } 527 528 #undef __FUNCT__ 529 #define __FUNCT__ "PCFactorSetMatOrderingType" 530 /*@C 531 PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 532 be used in the LU factorization. 533 534 Collective on PC 535 536 Input Parameters: 537 + pc - the preconditioner context 538 - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 539 540 Options Database Key: 541 . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 542 543 Level: intermediate 544 545 Notes: nested dissection is used by default 546 547 For Cholesky and ICC and the SBAIJ format reorderings are not available, 548 since only the upper triangular part of the matrix is stored. You can use the 549 SeqAIJ format in this case to get reorderings. 550 551 @*/ 552 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,MatOrderingType ordering) 553 { 554 PetscErrorCode ierr,(*f)(PC,MatOrderingType); 555 556 PetscFunctionBegin; 557 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr); 558 if (f) { 559 ierr = (*f)(pc,ordering);CHKERRQ(ierr); 560 } 561 PetscFunctionReturn(0); 562 } 563 564 #undef __FUNCT__ 565 #define __FUNCT__ "PCFactorSetPivoting" 566 /*@ 567 PCFactorSetPivoting - Determines when pivoting is done during LU. 568 For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 569 it is never done. For the Matlab and SuperLU factorization this is used. 570 571 Collective on PC 572 573 Input Parameters: 574 + pc - the preconditioner context 575 - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 576 577 Options Database Key: 578 . -pc_factor_pivoting <dtcol> 579 580 Level: intermediate 581 582 .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks() 583 @*/ 584 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC pc,PetscReal dtcol) 585 { 586 PetscErrorCode ierr,(*f)(PC,PetscReal); 587 588 PetscFunctionBegin; 589 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr); 590 if (f) { 591 ierr = (*f)(pc,dtcol);CHKERRQ(ierr); 592 } 593 PetscFunctionReturn(0); 594 } 595 596 #undef __FUNCT__ 597 #define __FUNCT__ "PCFactorSetPivotInBlocks" 598 /*@ 599 PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 600 with BAIJ or SBAIJ matrices 601 602 Collective on PC 603 604 Input Parameters: 605 + pc - the preconditioner context 606 - pivot - PETSC_TRUE or PETSC_FALSE 607 608 Options Database Key: 609 . -pc_factor_pivot_in_blocks <true,false> 610 611 Level: intermediate 612 613 .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting() 614 @*/ 615 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot) 616 { 617 PetscErrorCode ierr,(*f)(PC,PetscTruth); 618 619 PetscFunctionBegin; 620 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 621 if (f) { 622 ierr = (*f)(pc,pivot);CHKERRQ(ierr); 623 } 624 PetscFunctionReturn(0); 625 } 626 627 /* ------------------------------------------------------------------------ */ 628 629 /*MC 630 PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 631 632 Options Database Keys: 633 + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 634 . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 635 . -pc_factor_fill <fill> - Sets fill amount 636 . -pc_factor_in_place - Activates in-place factorization 637 . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 638 . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 639 stability of factorization. 640 . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 641 . -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 642 is optional with PETSC_TRUE being the default 643 - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 644 diagonal. 645 646 Notes: Not all options work for all matrix formats 647 Run with -help to see additional options for particular matrix formats or factorization 648 algorithms 649 650 Level: beginner 651 652 Concepts: LU factorization, direct solver 653 654 Notes: Usually this will compute an "exact" solution in one iteration and does 655 not need a Krylov method (i.e. you can use -ksp_type preonly, or 656 KSPSetType(ksp,KSPPREONLY) for the Krylov method 657 658 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 659 PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 660 PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetPivoting(), 661 PCFactorSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(), PCFactorReorderForNonzeroDiagonal() 662 M*/ 663 664 EXTERN_C_BEGIN 665 #undef __FUNCT__ 666 #define __FUNCT__ "PCCreate_LU" 667 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc) 668 { 669 PetscErrorCode ierr; 670 PetscMPIInt size; 671 PC_LU *dir; 672 673 PetscFunctionBegin; 674 ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr); 675 676 ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr); 677 dir->fact = 0; 678 dir->inplace = PETSC_FALSE; 679 dir->nonzerosalongdiagonal = PETSC_FALSE; 680 681 dir->info.fill = 5.0; 682 dir->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 683 dir->info.shiftnz = 0.0; 684 dir->info.zeropivot = 1.e-12; 685 dir->info.pivotinblocks = 1.0; 686 dir->info.shiftpd = 0.0; /* false */ 687 dir->info.shift_fraction = 0.0; 688 dir->col = 0; 689 dir->row = 0; 690 691 ierr = PetscStrallocpy("petsc",&dir->solvertype);CHKERRQ(ierr); 692 ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 693 if (size == 1) { 694 ierr = PetscStrallocpy(MATORDERING_ND,&dir->ordering);CHKERRQ(ierr); 695 } else { 696 ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr); 697 } 698 dir->reusefill = PETSC_FALSE; 699 dir->reuseordering = PETSC_FALSE; 700 pc->data = (void*)dir; 701 702 pc->ops->destroy = PCDestroy_LU; 703 pc->ops->apply = PCApply_LU; 704 pc->ops->applytranspose = PCApplyTranspose_LU; 705 pc->ops->setup = PCSetUp_LU; 706 pc->ops->setfromoptions = PCSetFromOptions_LU; 707 pc->ops->view = PCView_LU; 708 pc->ops->applyrichardson = 0; 709 pc->ops->getfactoredmatrix = PCFactorGetMatrix_LU; 710 711 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_LU", 712 PCFactorSetZeroPivot_LU);CHKERRQ(ierr); 713 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_LU", 714 PCFactorSetShiftNonzero_LU);CHKERRQ(ierr); 715 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_LU", 716 PCFactorSetShiftPd_LU);CHKERRQ(ierr); 717 718 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_LU", 719 PCFactorSetFill_LU);CHKERRQ(ierr); 720 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU", 721 PCFactorSetUseInPlace_LU);CHKERRQ(ierr); 722 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_LU", 723 PCFactorSetMatOrderingType_LU);CHKERRQ(ierr); 724 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU", 725 PCFactorSetReuseOrdering_LU);CHKERRQ(ierr); 726 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU", 727 PCFactorSetReuseFill_LU);CHKERRQ(ierr); 728 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivoting_C","PCFactorSetPivoting_LU", 729 PCFactorSetPivoting_LU);CHKERRQ(ierr); 730 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_LU", 731 PCFactorSetPivotInBlocks_LU);CHKERRQ(ierr); 732 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU", 733 PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 734 PetscFunctionReturn(0); 735 } 736 EXTERN_C_END 737