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