1 2 /* 3 Defines a Eisenstat trick SSOR preconditioner. This uses about 4 %50 of the usual amount of floating point ops used for SSOR + Krylov 5 method. But it requires actually solving the preconditioned problem 6 with both left and right preconditioning. 7 */ 8 #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 9 10 typedef struct { 11 Mat shell, A; 12 Vec b[2], diag; /* temporary storage for true right hand side */ 13 PetscReal omega; 14 PetscBool usediag; /* indicates preconditioner should include diagonal scaling*/ 15 } PC_Eisenstat; 16 17 static PetscErrorCode PCMult_Eisenstat(Mat mat, Vec b, Vec x) { 18 PC pc; 19 PC_Eisenstat *eis; 20 21 PetscFunctionBegin; 22 PetscCall(MatShellGetContext(mat, &pc)); 23 eis = (PC_Eisenstat *)pc->data; 24 PetscCall(MatSOR(eis->A, b, eis->omega, SOR_EISENSTAT, 0.0, 1, 1, x)); 25 PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason)); 26 PetscFunctionReturn(0); 27 } 28 29 static PetscErrorCode PCNorm_Eisenstat(Mat mat, NormType type, PetscReal *nrm) { 30 PC pc; 31 PC_Eisenstat *eis; 32 33 PetscFunctionBegin; 34 PetscCall(MatShellGetContext(mat, &pc)); 35 eis = (PC_Eisenstat *)pc->data; 36 PetscCall(MatNorm(eis->A, type, nrm)); 37 PetscFunctionReturn(0); 38 } 39 40 static PetscErrorCode PCApply_Eisenstat(PC pc, Vec x, Vec y) { 41 PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 42 PetscBool hasop; 43 44 PetscFunctionBegin; 45 if (eis->usediag) { 46 PetscCall(MatHasOperation(pc->pmat, MATOP_MULT_DIAGONAL_BLOCK, &hasop)); 47 if (hasop) { 48 PetscCall(MatMultDiagonalBlock(pc->pmat, x, y)); 49 } else { 50 PetscCall(VecPointwiseMult(y, x, eis->diag)); 51 } 52 } else PetscCall(VecCopy(x, y)); 53 PetscFunctionReturn(0); 54 } 55 56 static PetscErrorCode PCApplyTranspose_Eisenstat(PC pc, Vec x, Vec y) { 57 PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 58 PetscBool hasop, set, sym; 59 60 PetscFunctionBegin; 61 PetscCall(MatIsSymmetricKnown(eis->A, &set, &sym)); 62 PetscCheck(set && sym, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Can only apply transpose of Eisenstat if matrix is symmetric"); 63 if (eis->usediag) { 64 PetscCall(MatHasOperation(pc->pmat, MATOP_MULT_DIAGONAL_BLOCK, &hasop)); 65 if (hasop) { 66 PetscCall(MatMultDiagonalBlock(pc->pmat, x, y)); 67 } else { 68 PetscCall(VecPointwiseMult(y, x, eis->diag)); 69 } 70 } else PetscCall(VecCopy(x, y)); 71 PetscFunctionReturn(0); 72 } 73 74 static PetscErrorCode PCPreSolve_Eisenstat(PC pc, KSP ksp, Vec b, Vec x) { 75 PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 76 PetscBool nonzero; 77 78 PetscFunctionBegin; 79 if (pc->presolvedone < 2) { 80 PetscCheck(pc->mat == pc->pmat, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot have different mat and pmat"); 81 /* swap shell matrix and true matrix */ 82 eis->A = pc->mat; 83 pc->mat = eis->shell; 84 } 85 86 if (!eis->b[pc->presolvedone - 1]) { PetscCall(VecDuplicate(b, &eis->b[pc->presolvedone - 1])); } 87 88 /* if nonzero initial guess, modify x */ 89 PetscCall(KSPGetInitialGuessNonzero(ksp, &nonzero)); 90 if (nonzero) { 91 PetscCall(VecCopy(x, eis->b[pc->presolvedone - 1])); 92 PetscCall(MatSOR(eis->A, eis->b[pc->presolvedone - 1], eis->omega, SOR_APPLY_UPPER, 0.0, 1, 1, x)); 93 PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason)); 94 } 95 96 /* save true b, other option is to swap pointers */ 97 PetscCall(VecCopy(b, eis->b[pc->presolvedone - 1])); 98 99 /* modify b by (L + D/omega)^{-1} */ 100 PetscCall(MatSOR(eis->A, eis->b[pc->presolvedone - 1], eis->omega, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP), 0.0, 1, 1, b)); 101 PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason)); 102 PetscFunctionReturn(0); 103 } 104 105 static PetscErrorCode PCPostSolve_Eisenstat(PC pc, KSP ksp, Vec b, Vec x) { 106 PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 107 108 PetscFunctionBegin; 109 /* get back true b */ 110 PetscCall(VecCopy(eis->b[pc->presolvedone], b)); 111 112 /* modify x by (U + D/omega)^{-1} */ 113 PetscCall(VecCopy(x, eis->b[pc->presolvedone])); 114 PetscCall(MatSOR(eis->A, eis->b[pc->presolvedone], eis->omega, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP), 0.0, 1, 1, x)); 115 PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason)); 116 if (!pc->presolvedone) pc->mat = eis->A; 117 PetscFunctionReturn(0); 118 } 119 120 static PetscErrorCode PCReset_Eisenstat(PC pc) { 121 PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 122 123 PetscFunctionBegin; 124 PetscCall(VecDestroy(&eis->b[0])); 125 PetscCall(VecDestroy(&eis->b[1])); 126 PetscCall(MatDestroy(&eis->shell)); 127 PetscCall(VecDestroy(&eis->diag)); 128 PetscFunctionReturn(0); 129 } 130 131 static PetscErrorCode PCDestroy_Eisenstat(PC pc) { 132 PetscFunctionBegin; 133 PetscCall(PCReset_Eisenstat(pc)); 134 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetOmega_C", NULL)); 135 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetNoDiagonalScaling_C", NULL)); 136 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetOmega_C", NULL)); 137 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetNoDiagonalScaling_C", NULL)); 138 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", NULL)); 139 PetscCall(PetscFree(pc->data)); 140 PetscFunctionReturn(0); 141 } 142 143 static PetscErrorCode PCSetFromOptions_Eisenstat(PC pc, PetscOptionItems *PetscOptionsObject) { 144 PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 145 PetscBool set, flg; 146 147 PetscFunctionBegin; 148 PetscOptionsHeadBegin(PetscOptionsObject, "Eisenstat SSOR options"); 149 PetscCall(PetscOptionsReal("-pc_eisenstat_omega", "Relaxation factor 0 < omega < 2", "PCEisenstatSetOmega", eis->omega, &eis->omega, NULL)); 150 PetscCall(PetscOptionsBool("-pc_eisenstat_no_diagonal_scaling", "Do not use standard diagonal scaling", "PCEisenstatSetNoDiagonalScaling", eis->usediag ? PETSC_FALSE : PETSC_TRUE, &flg, &set)); 151 if (set) PetscCall(PCEisenstatSetNoDiagonalScaling(pc, flg)); 152 PetscOptionsHeadEnd(); 153 PetscFunctionReturn(0); 154 } 155 156 static PetscErrorCode PCView_Eisenstat(PC pc, PetscViewer viewer) { 157 PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 158 PetscBool iascii; 159 160 PetscFunctionBegin; 161 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 162 if (iascii) { 163 PetscCall(PetscViewerASCIIPrintf(viewer, " omega = %g\n", (double)eis->omega)); 164 if (eis->usediag) { 165 PetscCall(PetscViewerASCIIPrintf(viewer, " Using diagonal scaling (default)\n")); 166 } else { 167 PetscCall(PetscViewerASCIIPrintf(viewer, " Not using diagonal scaling\n")); 168 } 169 } 170 PetscFunctionReturn(0); 171 } 172 173 static PetscErrorCode PCSetUp_Eisenstat(PC pc) { 174 PetscInt M, N, m, n; 175 PetscBool set, sym; 176 PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 177 178 PetscFunctionBegin; 179 if (!pc->setupcalled) { 180 PetscCall(MatGetSize(pc->mat, &M, &N)); 181 PetscCall(MatGetLocalSize(pc->mat, &m, &n)); 182 PetscCall(MatIsSymmetricKnown(pc->mat, &set, &sym)); 183 PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &eis->shell)); 184 PetscCall(MatSetSizes(eis->shell, m, n, M, N)); 185 PetscCall(MatSetType(eis->shell, MATSHELL)); 186 PetscCall(MatSetUp(eis->shell)); 187 PetscCall(MatShellSetContext(eis->shell, pc)); 188 PetscCall(MatShellSetOperation(eis->shell, MATOP_MULT, (void (*)(void))PCMult_Eisenstat)); 189 if (set && sym) PetscCall(MatShellSetOperation(eis->shell, MATOP_MULT_TRANSPOSE, (void (*)(void))PCMult_Eisenstat)); 190 PetscCall(MatShellSetOperation(eis->shell, MATOP_NORM, (void (*)(void))PCNorm_Eisenstat)); 191 } 192 if (!eis->usediag) PetscFunctionReturn(0); 193 if (!pc->setupcalled) { PetscCall(MatCreateVecs(pc->pmat, &eis->diag, NULL)); } 194 PetscCall(MatGetDiagonal(pc->pmat, eis->diag)); 195 PetscFunctionReturn(0); 196 } 197 198 /* --------------------------------------------------------------------*/ 199 200 static PetscErrorCode PCEisenstatSetOmega_Eisenstat(PC pc, PetscReal omega) { 201 PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 202 203 PetscFunctionBegin; 204 PetscCheck(omega > 0.0 && omega < 2.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Relaxation out of range"); 205 eis->omega = omega; 206 PetscFunctionReturn(0); 207 } 208 209 static PetscErrorCode PCEisenstatSetNoDiagonalScaling_Eisenstat(PC pc, PetscBool flg) { 210 PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 211 212 PetscFunctionBegin; 213 eis->usediag = flg; 214 PetscFunctionReturn(0); 215 } 216 217 static PetscErrorCode PCEisenstatGetOmega_Eisenstat(PC pc, PetscReal *omega) { 218 PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 219 220 PetscFunctionBegin; 221 *omega = eis->omega; 222 PetscFunctionReturn(0); 223 } 224 225 static PetscErrorCode PCEisenstatGetNoDiagonalScaling_Eisenstat(PC pc, PetscBool *flg) { 226 PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 227 228 PetscFunctionBegin; 229 *flg = eis->usediag; 230 PetscFunctionReturn(0); 231 } 232 233 /*@ 234 PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega, 235 to use with Eisenstat's trick (where omega = 1.0 by default) 236 237 Logically Collective on pc 238 239 Input Parameters: 240 + pc - the preconditioner context 241 - omega - relaxation coefficient (0 < omega < 2) 242 243 Options Database Key: 244 . -pc_eisenstat_omega <omega> - Sets omega 245 246 Notes: 247 The Eisenstat trick implementation of SSOR requires about 50% of the 248 usual amount of floating point operations used for SSOR + Krylov method; 249 however, the preconditioned problem must be solved with both left 250 and right preconditioning. 251 252 To use SSOR without the Eisenstat trick, employ the `PCSOR` preconditioner, 253 which can be chosen with the database options 254 $ -pc_type sor -pc_sor_symmetric 255 256 Level: intermediate 257 258 .seealso: `PCSORSetOmega()`, `PCEISENSTAT` 259 @*/ 260 PetscErrorCode PCEisenstatSetOmega(PC pc, PetscReal omega) { 261 PetscFunctionBegin; 262 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 263 PetscValidLogicalCollectiveReal(pc, omega, 2); 264 PetscTryMethod(pc, "PCEisenstatSetOmega_C", (PC, PetscReal), (pc, omega)); 265 PetscFunctionReturn(0); 266 } 267 268 /*@ 269 PCEisenstatSetNoDiagonalScaling - Causes the Eisenstat preconditioner, `PCEISENSTAT` 270 not to do additional diagonal preconditioning. For matrices with a constant 271 along the diagonal, this may save a small amount of work. 272 273 Logically Collective on pc 274 275 Input Parameters: 276 + pc - the preconditioner context 277 - flg - `PETSC_TRUE` turns off diagonal scaling inside the algorithm 278 279 Options Database Key: 280 . -pc_eisenstat_no_diagonal_scaling - Activates `PCEisenstatSetNoDiagonalScaling()` 281 282 Level: intermediate 283 284 Note: 285 If you use the `KSPSetDiagonalScaling()` or -ksp_diagonal_scale option then you will 286 likely want to use this routine since it will save you some unneeded flops. 287 288 .seealso: `PCEisenstatSetOmega()`, `PCEISENSTAT` 289 @*/ 290 PetscErrorCode PCEisenstatSetNoDiagonalScaling(PC pc, PetscBool flg) { 291 PetscFunctionBegin; 292 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 293 PetscTryMethod(pc, "PCEisenstatSetNoDiagonalScaling_C", (PC, PetscBool), (pc, flg)); 294 PetscFunctionReturn(0); 295 } 296 297 /*@ 298 PCEisenstatGetOmega - Gets the SSOR relaxation coefficient, omega, 299 to use with Eisenstat's trick (where omega = 1.0 by default). 300 301 Logically Collective on pc 302 303 Input Parameter: 304 . pc - the preconditioner context 305 306 Output Parameter: 307 . omega - relaxation coefficient (0 < omega < 2) 308 309 Options Database Key: 310 . -pc_eisenstat_omega <omega> - Sets omega 311 312 Notes: 313 The Eisenstat trick implementation of SSOR requires about 50% of the 314 usual amount of floating point operations used for SSOR + Krylov method; 315 however, the preconditioned problem must be solved with both left 316 and right preconditioning. 317 318 To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner, 319 which can be chosen with the database options 320 $ -pc_type sor -pc_sor_symmetric 321 322 Level: intermediate 323 324 .seealso: `PCEISENSTAT`, `PCSORGetOmega()`, `PCEisenstatSetOmega()` 325 @*/ 326 PetscErrorCode PCEisenstatGetOmega(PC pc, PetscReal *omega) { 327 PetscFunctionBegin; 328 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 329 PetscUseMethod(pc, "PCEisenstatGetOmega_C", (PC, PetscReal *), (pc, omega)); 330 PetscFunctionReturn(0); 331 } 332 333 /*@ 334 PCEisenstatGetNoDiagonalScaling - Tells if the Eisenstat preconditioner 335 not to do additional diagonal preconditioning. For matrices with a constant 336 along the diagonal, this may save a small amount of work. 337 338 Logically Collective on pc 339 340 Input Parameter: 341 . pc - the preconditioner context 342 343 Output Parameter: 344 . flg - `PETSC_TRUE` means there is no diagonal scaling applied 345 346 Options Database Key: 347 . -pc_eisenstat_no_diagonal_scaling - Activates `PCEisenstatSetNoDiagonalScaling()` 348 349 Level: intermediate 350 351 Note: 352 If you use the KSPSetDiagonalScaling() or -ksp_diagonal_scale option then you will 353 likely want to use this routine since it will save you some unneeded flops. 354 355 .seealso: , `PCEISENSTAT`, `PCEisenstatGetOmega()` 356 @*/ 357 PetscErrorCode PCEisenstatGetNoDiagonalScaling(PC pc, PetscBool *flg) { 358 PetscFunctionBegin; 359 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 360 PetscUseMethod(pc, "PCEisenstatGetNoDiagonalScaling_C", (PC, PetscBool *), (pc, flg)); 361 PetscFunctionReturn(0); 362 } 363 364 static PetscErrorCode PCPreSolveChangeRHS_Eisenstat(PC pc, PetscBool *change) { 365 PetscFunctionBegin; 366 *change = PETSC_TRUE; 367 PetscFunctionReturn(0); 368 } 369 370 /*MC 371 PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel) 372 preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed. 373 374 Options Database Keys: 375 + -pc_eisenstat_omega <omega> - Sets omega 376 - -pc_eisenstat_no_diagonal_scaling - Activates `PCEisenstatSetNoDiagonalScaling()` 377 378 Level: beginner 379 380 Notes: 381 Only implemented for the `MATAIJ` matrix format. 382 383 Not a true parallel SOR, in parallel this implementation corresponds to block Jacobi with SOR on each block. 384 385 Developer Note: 386 Since this algorithm runs the Krylov method on a transformed linear system the implementation provides `PCPreSolve()` and `PCPostSolve()` 387 routines that `KSP` uses to set up the transformed linear system. 388 389 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCEisenstatGetOmega()`, 390 `PCEisenstatSetNoDiagonalScaling()`, `PCEisenstatSetOmega()`, `PCSOR` 391 M*/ 392 393 PETSC_EXTERN PetscErrorCode PCCreate_Eisenstat(PC pc) { 394 PC_Eisenstat *eis; 395 396 PetscFunctionBegin; 397 PetscCall(PetscNew(&eis)); 398 399 pc->ops->apply = PCApply_Eisenstat; 400 pc->ops->applytranspose = PCApplyTranspose_Eisenstat; 401 pc->ops->presolve = PCPreSolve_Eisenstat; 402 pc->ops->postsolve = PCPostSolve_Eisenstat; 403 pc->ops->applyrichardson = NULL; 404 pc->ops->setfromoptions = PCSetFromOptions_Eisenstat; 405 pc->ops->destroy = PCDestroy_Eisenstat; 406 pc->ops->reset = PCReset_Eisenstat; 407 pc->ops->view = PCView_Eisenstat; 408 pc->ops->setup = PCSetUp_Eisenstat; 409 410 pc->data = eis; 411 eis->omega = 1.0; 412 eis->b[0] = NULL; 413 eis->b[1] = NULL; 414 eis->diag = NULL; 415 eis->usediag = PETSC_TRUE; 416 417 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetOmega_C", PCEisenstatSetOmega_Eisenstat)); 418 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetNoDiagonalScaling_C", PCEisenstatSetNoDiagonalScaling_Eisenstat)); 419 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetOmega_C", PCEisenstatGetOmega_Eisenstat)); 420 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetNoDiagonalScaling_C", PCEisenstatGetNoDiagonalScaling_Eisenstat)); 421 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", PCPreSolveChangeRHS_Eisenstat)); 422 PetscFunctionReturn(0); 423 } 424