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