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 #include "src/ksp/pc/impls/factor/factor.h" /*I "petscpc.h" I*/ 9 10 typedef struct { 11 PC_Factor hdr; 12 PetscReal actualfill; /* actual fill in factor */ 13 PetscTruth inplace; /* flag indicating in-place factorization */ 14 IS row,col; /* index sets used for reordering */ 15 PetscTruth reuseordering; /* reuses previous reordering computed */ 16 PetscTruth reusefill; /* reuse fill from previous Cholesky */ 17 } PC_Cholesky; 18 19 EXTERN_C_BEGIN 20 #undef __FUNCT__ 21 #define __FUNCT__ "PCFactorSetMatSolverPackage_Cholesky" 22 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage_Cholesky(PC pc,const MatSolverPackage stype) 23 { 24 PetscErrorCode ierr; 25 PC_Cholesky *cholesky = (PC_Cholesky*)pc->data; 26 27 PetscFunctionBegin; 28 if (((PC_Factor*)cholesky)->fact) { 29 const MatSolverPackage ltype; 30 PetscTruth flg; 31 ierr = MatFactorGetSolverPackage(((PC_Factor*)cholesky)->fact,<ype);CHKERRQ(ierr); 32 ierr = PetscStrcmp(stype,ltype,&flg);CHKERRQ(ierr); 33 if (!flg) { 34 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used"); 35 } 36 } 37 ierr = PetscStrfree(((PC_Factor*)cholesky)->solvertype);CHKERRQ(ierr); 38 ierr = PetscStrallocpy(stype,&((PC_Factor*)cholesky)->solvertype);CHKERRQ(ierr); 39 PetscFunctionReturn(0); 40 } 41 EXTERN_C_END 42 43 EXTERN_C_BEGIN 44 #undef __FUNCT__ 45 #define __FUNCT__ "PCFactorSetZeroPivot_Cholesky" 46 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_Cholesky(PC pc,PetscReal z) 47 { 48 PC_Cholesky *ch = (PC_Cholesky*)pc->data; 49 50 PetscFunctionBegin; 51 ((PC_Factor*)ch)->info.zeropivot = z; 52 PetscFunctionReturn(0); 53 } 54 EXTERN_C_END 55 56 EXTERN_C_BEGIN 57 #undef __FUNCT__ 58 #define __FUNCT__ "PCFactorSetShiftNonzero_Cholesky" 59 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_Cholesky(PC pc,PetscReal shift) 60 { 61 PC_Cholesky *dir; 62 63 PetscFunctionBegin; 64 dir = (PC_Cholesky*)pc->data; 65 if (shift == (PetscReal) PETSC_DECIDE) { 66 ((PC_Factor*)dir)->info.shiftnz = 1.e-12; 67 } else { 68 ((PC_Factor*)dir)->info.shiftnz = shift; 69 } 70 PetscFunctionReturn(0); 71 } 72 EXTERN_C_END 73 74 EXTERN_C_BEGIN 75 #undef __FUNCT__ 76 #define __FUNCT__ "PCFactorSetShiftPd_Cholesky" 77 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_Cholesky(PC pc,PetscTruth shift) 78 { 79 PC_Cholesky *dir; 80 81 PetscFunctionBegin; 82 dir = (PC_Cholesky*)pc->data; 83 if (shift) { 84 ((PC_Factor*)dir)->info.shiftpd = 1.0; 85 } else { 86 ((PC_Factor*)dir)->info.shiftpd = 0.0; 87 } 88 PetscFunctionReturn(0); 89 } 90 EXTERN_C_END 91 92 EXTERN_C_BEGIN 93 #undef __FUNCT__ 94 #define __FUNCT__ "PCFactorSetReuseOrdering_Cholesky" 95 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_Cholesky(PC pc,PetscTruth flag) 96 { 97 PC_Cholesky *lu; 98 99 PetscFunctionBegin; 100 lu = (PC_Cholesky*)pc->data; 101 lu->reuseordering = flag; 102 PetscFunctionReturn(0); 103 } 104 EXTERN_C_END 105 106 EXTERN_C_BEGIN 107 #undef __FUNCT__ 108 #define __FUNCT__ "PCFactorSetReuseFill_Cholesky" 109 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_Cholesky(PC pc,PetscTruth flag) 110 { 111 PC_Cholesky *lu; 112 113 PetscFunctionBegin; 114 lu = (PC_Cholesky*)pc->data; 115 lu->reusefill = flag; 116 PetscFunctionReturn(0); 117 } 118 EXTERN_C_END 119 120 #undef __FUNCT__ 121 #define __FUNCT__ "PCSetFromOptions_Cholesky" 122 static PetscErrorCode PCSetFromOptions_Cholesky(PC pc) 123 { 124 PC_Cholesky *lu = (PC_Cholesky*)pc->data; 125 PetscErrorCode ierr; 126 PetscTruth flg; 127 char tname[256], solvertype[64]; 128 PetscFList ordlist; 129 130 PetscFunctionBegin; 131 ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 132 ierr = PetscOptionsHead("Cholesky options");CHKERRQ(ierr); 133 ierr = PetscOptionsName("-pc_factor_in_place","Form Cholesky in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr); 134 if (flg) { 135 ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr); 136 } 137 ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in Cholesky/non-zeros in matrix","PCFactorSetFill",((PC_Factor*)lu)->info.fill,&((PC_Factor*)lu)->info.fill,0);CHKERRQ(ierr); 138 139 ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr); 140 if (flg) { 141 ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr); 142 } 143 ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr); 144 if (flg) { 145 ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr); 146 } 147 148 ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 149 ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in Cholesky","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)lu)->ordering,tname,256,&flg);CHKERRQ(ierr); 150 if (flg) { 151 ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 152 } 153 154 /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */ 155 ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific Cholesky solver to use","MatGetFactor",((PC_Factor*)lu)->solvertype,solvertype,64,&flg);CHKERRQ(ierr); 156 if (flg) { 157 ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr); 158 } 159 160 ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 161 if (flg) { 162 ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr); 163 } 164 ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",((PC_Factor*)lu)->info.shiftnz,&((PC_Factor*)lu)->info.shiftnz,0);CHKERRQ(ierr); 165 ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr); 166 if (flg) { 167 ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr); 168 } 169 ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)lu)->info.zeropivot,&((PC_Factor*)lu)->info.zeropivot,0);CHKERRQ(ierr); 170 171 ierr = PetscOptionsTail();CHKERRQ(ierr); 172 PetscFunctionReturn(0); 173 } 174 175 #undef __FUNCT__ 176 #define __FUNCT__ "PCView_Cholesky" 177 static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer) 178 { 179 PC_Cholesky *lu = (PC_Cholesky*)pc->data; 180 PetscErrorCode ierr; 181 PetscTruth iascii,isstring; 182 183 PetscFunctionBegin; 184 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 185 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 186 if (iascii) { 187 188 if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," Cholesky: in-place factorization\n");CHKERRQ(ierr);} 189 else {ierr = PetscViewerASCIIPrintf(viewer," Cholesky: out-of-place factorization\n");CHKERRQ(ierr);} 190 ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",((PC_Factor*)lu)->ordering);CHKERRQ(ierr); 191 if (((PC_Factor*)lu)->fact) { 192 ierr = PetscViewerASCIIPrintf(viewer," Cholesky: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr); 193 ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 194 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 195 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 196 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 197 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 198 ierr = MatView(((PC_Factor*)lu)->fact,viewer);CHKERRQ(ierr); 199 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 200 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 201 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 202 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 203 } 204 if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 205 if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 206 } else if (isstring) { 207 ierr = PetscViewerStringSPrintf(viewer," order=%s",((PC_Factor*)lu)->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 208 } else { 209 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCCHOLESKY",((PetscObject)viewer)->type_name); 210 } 211 PetscFunctionReturn(0); 212 } 213 214 #undef __FUNCT__ 215 #define __FUNCT__ "PCFactorGetMatrix_Cholesky" 216 static PetscErrorCode PCFactorGetMatrix_Cholesky(PC pc,Mat *mat) 217 { 218 PC_Cholesky *dir = (PC_Cholesky*)pc->data; 219 220 PetscFunctionBegin; 221 if (!((PC_Factor*)dir)->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 222 *mat = ((PC_Factor*)dir)->fact; 223 PetscFunctionReturn(0); 224 } 225 226 #undef __FUNCT__ 227 #define __FUNCT__ "PCSetUp_Cholesky" 228 static PetscErrorCode PCSetUp_Cholesky(PC pc) 229 { 230 PetscErrorCode ierr; 231 PetscTruth flg; 232 PC_Cholesky *dir = (PC_Cholesky*)pc->data; 233 234 PetscFunctionBegin; 235 if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 236 237 if (dir->inplace) { 238 if (dir->row && dir->col && (dir->row != dir->col)) { 239 ierr = ISDestroy(dir->row);CHKERRQ(ierr); 240 dir->row = 0; 241 } 242 if (dir->col) { 243 ierr = ISDestroy(dir->col);CHKERRQ(ierr); 244 dir->col = 0; 245 } 246 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 247 if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 248 ierr = ISDestroy(dir->col);CHKERRQ(ierr); 249 dir->col=0; 250 } 251 if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);} 252 ierr = MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 253 ((PC_Factor*)dir)->fact = pc->pmat; 254 } else { 255 MatInfo info; 256 if (!pc->setupcalled) { 257 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 258 /* check if dir->row == dir->col */ 259 ierr = ISEqual(dir->row,dir->col,&flg);CHKERRQ(ierr); 260 if (!flg) SETERRQ(PETSC_ERR_ARG_INCOMP,"row and column permutations must equal"); 261 ierr = ISDestroy(dir->col);CHKERRQ(ierr); /* only pass one ordering into CholeskyFactor */ 262 dir->col=0; 263 264 ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg);CHKERRQ(ierr); 265 if (flg) { 266 PetscReal tol = 1.e-10; 267 ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr); 268 ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 269 } 270 if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);} 271 ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 272 ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 273 ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 274 dir->actualfill = info.fill_ratio_needed; 275 ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 276 } else if (pc->flag != SAME_NONZERO_PATTERN) { 277 if (!dir->reuseordering) { 278 if (dir->row && dir->col && (dir->row != dir->col)) { 279 ierr = ISDestroy(dir->row);CHKERRQ(ierr); 280 dir->row = 0; 281 } 282 if (dir->col) { 283 ierr = ISDestroy(dir->col);CHKERRQ(ierr); 284 dir->col =0; 285 } 286 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 287 if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 288 ierr = ISDestroy(dir->col);CHKERRQ(ierr); 289 dir->col=0; 290 } 291 ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg);CHKERRQ(ierr); 292 if (flg) { 293 PetscReal tol = 1.e-10; 294 ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr); 295 ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 296 } 297 if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);} 298 } 299 ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 300 ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 301 ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 302 ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 303 dir->actualfill = info.fill_ratio_needed; 304 ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 305 } 306 ierr = MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 307 } 308 PetscFunctionReturn(0); 309 } 310 311 #undef __FUNCT__ 312 #define __FUNCT__ "PCDestroy_Cholesky" 313 static PetscErrorCode PCDestroy_Cholesky(PC pc) 314 { 315 PC_Cholesky *dir = (PC_Cholesky*)pc->data; 316 PetscErrorCode ierr; 317 318 PetscFunctionBegin; 319 if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 320 if (dir->row) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 321 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 322 ierr = PetscStrfree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 323 ierr = PetscStrfree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 324 ierr = PetscFree(dir);CHKERRQ(ierr); 325 PetscFunctionReturn(0); 326 } 327 328 #undef __FUNCT__ 329 #define __FUNCT__ "PCApply_Cholesky" 330 static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y) 331 { 332 PC_Cholesky *dir = (PC_Cholesky*)pc->data; 333 PetscErrorCode ierr; 334 335 PetscFunctionBegin; 336 if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 337 else {ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 338 PetscFunctionReturn(0); 339 } 340 341 #undef __FUNCT__ 342 #define __FUNCT__ "PCApplyTranspose_Cholesky" 343 static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y) 344 { 345 PC_Cholesky *dir = (PC_Cholesky*)pc->data; 346 PetscErrorCode ierr; 347 348 PetscFunctionBegin; 349 if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 350 else {ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 351 PetscFunctionReturn(0); 352 } 353 354 /* -----------------------------------------------------------------------------------*/ 355 356 EXTERN_C_BEGIN 357 #undef __FUNCT__ 358 #define __FUNCT__ "PCFactorSetFill_Cholesky" 359 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_Cholesky(PC pc,PetscReal fill) 360 { 361 PC_Cholesky *dir; 362 363 PetscFunctionBegin; 364 dir = (PC_Cholesky*)pc->data; 365 ((PC_Factor*)dir)->info.fill = fill; 366 PetscFunctionReturn(0); 367 } 368 EXTERN_C_END 369 370 EXTERN_C_BEGIN 371 #undef __FUNCT__ 372 #define __FUNCT__ "PCFactorSetUseInPlace_Cholesky" 373 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_Cholesky(PC pc) 374 { 375 PC_Cholesky *dir; 376 377 PetscFunctionBegin; 378 dir = (PC_Cholesky*)pc->data; 379 dir->inplace = PETSC_TRUE; 380 PetscFunctionReturn(0); 381 } 382 EXTERN_C_END 383 384 EXTERN_C_BEGIN 385 #undef __FUNCT__ 386 #define __FUNCT__ "PCFactorSetMatOrderingType_Cholesky" 387 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_Cholesky(PC pc,const MatOrderingType ordering) 388 { 389 PC_Cholesky *dir = (PC_Cholesky*)pc->data; 390 PetscErrorCode ierr; 391 392 PetscFunctionBegin; 393 ierr = PetscStrfree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 394 ierr = PetscStrallocpy(ordering,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 395 PetscFunctionReturn(0); 396 } 397 EXTERN_C_END 398 399 /* -----------------------------------------------------------------------------------*/ 400 401 #undef __FUNCT__ 402 #define __FUNCT__ "PCFactorSetReuseOrdering" 403 /*@ 404 PCFactorSetReuseOrdering - When similar matrices are factored, this 405 causes the ordering computed in the first factor to be used for all 406 following factors. 407 408 Collective on PC 409 410 Input Parameters: 411 + pc - the preconditioner context 412 - flag - PETSC_TRUE to reuse else PETSC_FALSE 413 414 Options Database Key: 415 . -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 416 417 Level: intermediate 418 419 .keywords: PC, levels, reordering, factorization, incomplete, LU 420 421 .seealso: PCFactorSetReuseFill() 422 @*/ 423 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering(PC pc,PetscTruth flag) 424 { 425 PetscErrorCode ierr,(*f)(PC,PetscTruth); 426 427 PetscFunctionBegin; 428 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 429 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 430 if (f) { 431 ierr = (*f)(pc,flag);CHKERRQ(ierr); 432 } 433 PetscFunctionReturn(0); 434 } 435 436 #undef __FUNCT__ 437 #define __FUNCT__ "PCFactorSetReuseFill" 438 /*@ 439 PCFactorSetReuseFill - When matrices with same different nonzero structure are factored, 440 this causes later ones to use the fill ratio computed in the initial factorization. 441 442 Collective on PC 443 444 Input Parameters: 445 + pc - the preconditioner context 446 - flag - PETSC_TRUE to reuse else PETSC_FALSE 447 448 Options Database Key: 449 . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 450 451 Level: intermediate 452 453 .keywords: PC, levels, reordering, factorization, incomplete, Cholesky 454 455 .seealso: PCFactorSetReuseOrdering() 456 @*/ 457 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC pc,PetscTruth flag) 458 { 459 PetscErrorCode ierr,(*f)(PC,PetscTruth); 460 461 PetscFunctionBegin; 462 PetscValidHeaderSpecific(pc,PC_COOKIE,2); 463 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr); 464 if (f) { 465 ierr = (*f)(pc,flag);CHKERRQ(ierr); 466 } 467 PetscFunctionReturn(0); 468 } 469 470 /*MC 471 PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 472 473 Options Database Keys: 474 + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 475 . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles 476 . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 477 . -pc_factor_fill <fill> - Sets fill amount 478 . -pc_factor_in_place - Activates in-place factorization 479 . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 480 . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 481 - -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 482 is optional with PETSC_TRUE being the default 483 484 Notes: Not all options work for all matrix formats 485 486 Level: beginner 487 488 Concepts: Cholesky factorization, direct solver 489 490 Notes: Usually this will compute an "exact" solution in one iteration and does 491 not need a Krylov method (i.e. you can use -ksp_type preonly, or 492 KSPSetType(ksp,KSPPREONLY) for the Krylov method 493 494 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 495 PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 496 PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), 497 PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd() 498 499 M*/ 500 501 EXTERN_C_BEGIN 502 #undef __FUNCT__ 503 #define __FUNCT__ "PCCreate_Cholesky" 504 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Cholesky(PC pc) 505 { 506 PetscErrorCode ierr; 507 PC_Cholesky *dir; 508 509 PetscFunctionBegin; 510 ierr = PetscNewLog(pc,PC_Cholesky,&dir);CHKERRQ(ierr); 511 512 ((PC_Factor*)dir)->fact = 0; 513 dir->inplace = PETSC_FALSE; 514 ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 515 ((PC_Factor*)dir)->info.fill = 5.0; 516 ((PC_Factor*)dir)->info.shiftnz = 0.0; 517 ((PC_Factor*)dir)->info.shiftpd = 0.0; /* false */ 518 ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 519 dir->col = 0; 520 dir->row = 0; 521 ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 522 ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 523 dir->reusefill = PETSC_FALSE; 524 dir->reuseordering = PETSC_FALSE; 525 pc->data = (void*)dir; 526 527 pc->ops->destroy = PCDestroy_Cholesky; 528 pc->ops->apply = PCApply_Cholesky; 529 pc->ops->applytranspose = PCApplyTranspose_Cholesky; 530 pc->ops->setup = PCSetUp_Cholesky; 531 pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 532 pc->ops->view = PCView_Cholesky; 533 pc->ops->applyrichardson = 0; 534 pc->ops->getfactoredmatrix = PCFactorGetMatrix_Cholesky; 535 536 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Cholesky", 537 PCFactorSetMatSolverPackage_Cholesky);CHKERRQ(ierr); 538 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Cholesky", 539 PCFactorSetZeroPivot_Cholesky);CHKERRQ(ierr); 540 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Cholesky", 541 PCFactorSetShiftNonzero_Cholesky);CHKERRQ(ierr); 542 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Cholesky", 543 PCFactorSetShiftPd_Cholesky);CHKERRQ(ierr); 544 545 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Cholesky", 546 PCFactorSetFill_Cholesky);CHKERRQ(ierr); 547 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_Cholesky", 548 PCFactorSetUseInPlace_Cholesky);CHKERRQ(ierr); 549 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Cholesky", 550 PCFactorSetMatOrderingType_Cholesky);CHKERRQ(ierr); 551 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_Cholesky", 552 PCFactorSetReuseOrdering_Cholesky);CHKERRQ(ierr); 553 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_Cholesky", 554 PCFactorSetReuseFill_Cholesky);CHKERRQ(ierr); 555 PetscFunctionReturn(0); 556 } 557 EXTERN_C_END 558