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