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