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