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(((PetscObject)pc)->comm,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(pc,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) { 104 pc->mat = eis->A; 105 } 106 PetscFunctionReturn(0); 107 } 108 109 #undef __FUNCT__ 110 #define __FUNCT__ "PCReset_Eisenstat" 111 static PetscErrorCode PCReset_Eisenstat(PC pc) 112 { 113 PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 114 PetscErrorCode ierr; 115 116 PetscFunctionBegin; 117 ierr = VecDestroy(&eis->b[0]);CHKERRQ(ierr); 118 ierr = VecDestroy(&eis->b[1]);CHKERRQ(ierr); 119 ierr = MatDestroy(&eis->shell);CHKERRQ(ierr); 120 ierr = VecDestroy(&eis->diag);CHKERRQ(ierr); 121 PetscFunctionReturn(0); 122 } 123 124 #undef __FUNCT__ 125 #define __FUNCT__ "PCDestroy_Eisenstat" 126 static PetscErrorCode PCDestroy_Eisenstat(PC pc) 127 { 128 PetscErrorCode ierr; 129 130 PetscFunctionBegin; 131 ierr = PCReset_Eisenstat(pc);CHKERRQ(ierr); 132 ierr = PetscFree(pc->data);CHKERRQ(ierr); 133 PetscFunctionReturn(0); 134 } 135 136 #undef __FUNCT__ 137 #define __FUNCT__ "PCSetFromOptions_Eisenstat" 138 static PetscErrorCode PCSetFromOptions_Eisenstat(PC pc) 139 { 140 PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 141 PetscErrorCode ierr; 142 PetscBool flg = PETSC_FALSE; 143 144 PetscFunctionBegin; 145 ierr = PetscOptionsHead("Eisenstat SSOR options");CHKERRQ(ierr); 146 ierr = PetscOptionsReal("-pc_eisenstat_omega","Relaxation factor 0 < omega < 2","PCEisenstatSetOmega",eis->omega,&eis->omega,PETSC_NULL);CHKERRQ(ierr); 147 ierr = PetscOptionsBool("-pc_eisenstat_no_diagonal_scaling","Do not use standard diagonal scaling","PCEisenstatNoDiagonalScaling",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 148 if (flg) { 149 ierr = PCEisenstatNoDiagonalScaling(pc);CHKERRQ(ierr); 150 } 151 ierr = PetscOptionsTail();CHKERRQ(ierr); 152 PetscFunctionReturn(0); 153 } 154 155 #undef __FUNCT__ 156 #define __FUNCT__ "PCView_Eisenstat" 157 static PetscErrorCode PCView_Eisenstat(PC pc,PetscViewer viewer) 158 { 159 PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 160 PetscErrorCode ierr; 161 PetscBool iascii; 162 163 PetscFunctionBegin; 164 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 165 if (iascii) { 166 ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: omega = %G\n",eis->omega);CHKERRQ(ierr); 167 if (eis->usediag) { 168 ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Using diagonal scaling (default)\n");CHKERRQ(ierr); 169 } else { 170 ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Not using diagonal scaling\n");CHKERRQ(ierr); 171 } 172 } else { 173 SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for Eisenstat PC",((PetscObject)viewer)->type_name); 174 } 175 PetscFunctionReturn(0); 176 } 177 178 #undef __FUNCT__ 179 #define __FUNCT__ "PCSetUp_Eisenstat" 180 static PetscErrorCode PCSetUp_Eisenstat(PC pc) 181 { 182 PetscErrorCode ierr; 183 PetscInt M,N,m,n; 184 PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 185 186 PetscFunctionBegin; 187 if (!pc->setupcalled) { 188 ierr = MatGetSize(pc->mat,&M,&N);CHKERRQ(ierr); 189 ierr = MatGetLocalSize(pc->mat,&m,&n);CHKERRQ(ierr); 190 ierr = MatCreate(((PetscObject)pc)->comm,&eis->shell);CHKERRQ(ierr); 191 ierr = MatSetSizes(eis->shell,m,n,M,N);CHKERRQ(ierr); 192 ierr = MatSetType(eis->shell,MATSHELL);CHKERRQ(ierr); 193 ierr = MatSetUp(eis->shell);CHKERRQ(ierr); 194 ierr = MatShellSetContext(eis->shell,(void*)pc);CHKERRQ(ierr); 195 ierr = PetscLogObjectParent(pc,eis->shell);CHKERRQ(ierr); 196 ierr = MatShellSetOperation(eis->shell,MATOP_MULT,(void(*)(void))PCMult_Eisenstat);CHKERRQ(ierr); 197 } 198 if (!eis->usediag) PetscFunctionReturn(0); 199 if (!pc->setupcalled) { 200 ierr = MatGetVecs(pc->pmat,&eis->diag,0);CHKERRQ(ierr); 201 ierr = PetscLogObjectParent(pc,eis->diag);CHKERRQ(ierr); 202 } 203 ierr = MatGetDiagonal(pc->pmat,eis->diag);CHKERRQ(ierr); 204 PetscFunctionReturn(0); 205 } 206 207 /* --------------------------------------------------------------------*/ 208 209 EXTERN_C_BEGIN 210 #undef __FUNCT__ 211 #define __FUNCT__ "PCEisenstatSetOmega_Eisenstat" 212 PetscErrorCode PCEisenstatSetOmega_Eisenstat(PC pc,PetscReal omega) 213 { 214 PC_Eisenstat *eis; 215 216 PetscFunctionBegin; 217 if (omega >= 2.0 || omega <= 0.0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range"); 218 eis = (PC_Eisenstat*)pc->data; 219 eis->omega = omega; 220 PetscFunctionReturn(0); 221 } 222 EXTERN_C_END 223 224 EXTERN_C_BEGIN 225 #undef __FUNCT__ 226 #define __FUNCT__ "PCEisenstatNoDiagonalScaling_Eisenstat" 227 PetscErrorCode PCEisenstatNoDiagonalScaling_Eisenstat(PC pc) 228 { 229 PC_Eisenstat *eis; 230 231 PetscFunctionBegin; 232 eis = (PC_Eisenstat*)pc->data; 233 eis->usediag = PETSC_FALSE; 234 PetscFunctionReturn(0); 235 } 236 EXTERN_C_END 237 238 #undef __FUNCT__ 239 #define __FUNCT__ "PCEisenstatSetOmega" 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 .keywords: PC, Eisenstat, set, SOR, SSOR, relaxation, omega 266 267 .seealso: PCSORSetOmega() 268 @*/ 269 PetscErrorCode PCEisenstatSetOmega(PC pc,PetscReal omega) 270 { 271 PetscErrorCode ierr; 272 273 PetscFunctionBegin; 274 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 275 PetscValidLogicalCollectiveReal(pc,omega,2); 276 ierr = PetscTryMethod(pc,"PCEisenstatSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr); 277 PetscFunctionReturn(0); 278 } 279 280 #undef __FUNCT__ 281 #define __FUNCT__ "PCEisenstatNoDiagonalScaling" 282 /*@ 283 PCEisenstatNoDiagonalScaling - Causes the Eisenstat preconditioner 284 not to do additional diagonal preconditioning. For matrices with a constant 285 along the diagonal, this may save a small amount of work. 286 287 Logically Collective on PC 288 289 Input Parameter: 290 . pc - the preconditioner context 291 292 Options Database Key: 293 . -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling() 294 295 Level: intermediate 296 297 Note: 298 If you use the KPSSetDiagonalScaling() or -ksp_diagonal_scale option then you will 299 likley want to use this routine since it will save you some unneeded flops. 300 301 .keywords: PC, Eisenstat, use, diagonal, scaling, SSOR 302 303 .seealso: PCEisenstatSetOmega() 304 @*/ 305 PetscErrorCode PCEisenstatNoDiagonalScaling(PC pc) 306 { 307 PetscErrorCode ierr; 308 309 PetscFunctionBegin; 310 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 311 ierr = PetscTryMethod(pc,"PCEisenstatNoDiagonalScaling_C",(PC),(pc));CHKERRQ(ierr); 312 PetscFunctionReturn(0); 313 } 314 315 /* --------------------------------------------------------------------*/ 316 317 /*MC 318 PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel) 319 preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed. 320 321 Options Database Keys: 322 + -pc_eisenstat_omega <omega> - Sets omega 323 - -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling() 324 325 Level: beginner 326 327 Concepts: SOR, preconditioners, Gauss-Seidel, Eisenstat's trick 328 329 Notes: Only implemented for the SeqAIJ matrix format. 330 Not a true parallel SOR, in parallel this implementation corresponds to block 331 Jacobi with SOR on each block. 332 333 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 334 PCEisenstatNoDiagonalScaling(), PCEisenstatSetOmega(), PCSOR 335 M*/ 336 337 EXTERN_C_BEGIN 338 #undef __FUNCT__ 339 #define __FUNCT__ "PCCreate_Eisenstat" 340 PetscErrorCode PCCreate_Eisenstat(PC pc) 341 { 342 PetscErrorCode ierr; 343 PC_Eisenstat *eis; 344 345 PetscFunctionBegin; 346 ierr = PetscNewLog(pc,PC_Eisenstat,&eis);CHKERRQ(ierr); 347 348 pc->ops->apply = PCApply_Eisenstat; 349 pc->ops->presolve = PCPreSolve_Eisenstat; 350 pc->ops->postsolve = PCPostSolve_Eisenstat; 351 pc->ops->applyrichardson = 0; 352 pc->ops->setfromoptions = PCSetFromOptions_Eisenstat; 353 pc->ops->destroy = PCDestroy_Eisenstat; 354 pc->ops->reset = PCReset_Eisenstat; 355 pc->ops->view = PCView_Eisenstat; 356 pc->ops->setup = PCSetUp_Eisenstat; 357 358 pc->data = (void*)eis; 359 eis->omega = 1.0; 360 eis->b[0] = 0; 361 eis->b[1] = 0; 362 eis->diag = 0; 363 eis->usediag = PETSC_TRUE; 364 365 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatSetOmega_C","PCEisenstatSetOmega_Eisenstat", 366 PCEisenstatSetOmega_Eisenstat);CHKERRQ(ierr); 367 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C", 368 "PCEisenstatNoDiagonalScaling_Eisenstat", 369 PCEisenstatNoDiagonalScaling_Eisenstat);CHKERRQ(ierr); 370 PetscFunctionReturn(0); 371 } 372 EXTERN_C_END 373