1 #define PETSCKSP_DLL 2 3 /* 4 Defines a Eisenstat trick SSOR preconditioner. This uses about 5 %50 of the usual amount of floating point ops used for SSOR + Krylov 6 method. But it requires actually solving the preconditioned problem 7 with both left and right preconditioning. 8 */ 9 #include "private/pcimpl.h" /*I "petscpc.h" I*/ 10 11 typedef struct { 12 Mat shell,A; 13 Vec b,diag; /* temporary storage for true right hand side */ 14 PetscReal omega; 15 PetscBool usediag; /* indicates preconditioner should include diagonal scaling*/ 16 } PC_Eisenstat; 17 18 19 #undef __FUNCT__ 20 #define __FUNCT__ "PCMult_Eisenstat" 21 static PetscErrorCode PCMult_Eisenstat(Mat mat,Vec b,Vec x) 22 { 23 PetscErrorCode ierr; 24 PC pc; 25 PC_Eisenstat *eis; 26 27 PetscFunctionBegin; 28 ierr = MatShellGetContext(mat,(void **)&pc);CHKERRQ(ierr); 29 eis = (PC_Eisenstat*)pc->data; 30 ierr = MatSOR(eis->A,b,eis->omega,SOR_EISENSTAT,0.0,1,1,x);CHKERRQ(ierr); 31 PetscFunctionReturn(0); 32 } 33 34 #undef __FUNCT__ 35 #define __FUNCT__ "PCApply_Eisenstat" 36 static PetscErrorCode PCApply_Eisenstat(PC pc,Vec x,Vec y) 37 { 38 PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 39 PetscErrorCode ierr; 40 PetscBool hasop; 41 42 PetscFunctionBegin; 43 if (eis->usediag) { 44 ierr = MatHasOperation(pc->pmat,MATOP_MULT_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr); 45 if (hasop) { 46 ierr = MatMultDiagonalBlock(pc->pmat,x,y);CHKERRQ(ierr); 47 } else { 48 ierr = VecPointwiseMult(y,x,eis->diag);CHKERRQ(ierr); 49 } 50 } else {ierr = VecCopy(x,y);CHKERRQ(ierr);} 51 PetscFunctionReturn(0); 52 } 53 54 #undef __FUNCT__ 55 #define __FUNCT__ "PCPreSolve_Eisenstat" 56 static PetscErrorCode PCPreSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x) 57 { 58 PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 59 PetscBool nonzero; 60 PetscErrorCode ierr; 61 62 PetscFunctionBegin; 63 if (pc->mat != pc->pmat) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cannot have different mat and pmat"); 64 65 /* swap shell matrix and true matrix */ 66 eis->A = pc->mat; 67 pc->mat = eis->shell; 68 69 if (!eis->b) { 70 ierr = VecDuplicate(b,&eis->b);CHKERRQ(ierr); 71 ierr = PetscLogObjectParent(pc,eis->b);CHKERRQ(ierr); 72 } 73 74 75 /* if nonzero initial guess, modify x */ 76 ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 77 if (nonzero) { 78 ierr = VecCopy(x,eis->b);CHKERRQ(ierr); 79 ierr = MatSOR(eis->A,eis->b,eis->omega,SOR_APPLY_UPPER,0.0,1,1,x);CHKERRQ(ierr); 80 } 81 82 /* save true b, other option is to swap pointers */ 83 ierr = VecCopy(b,eis->b);CHKERRQ(ierr); 84 85 /* modify b by (L + D/omega)^{-1} */ 86 ierr = MatSOR(eis->A,eis->b,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),0.0,1,1,b);CHKERRQ(ierr); 87 PetscFunctionReturn(0); 88 } 89 90 #undef __FUNCT__ 91 #define __FUNCT__ "PCPostSolve_Eisenstat" 92 static PetscErrorCode PCPostSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x) 93 { 94 PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 95 PetscErrorCode ierr; 96 97 PetscFunctionBegin; 98 /* get back true b */ 99 ierr = VecCopy(eis->b,b);CHKERRQ(ierr); 100 101 /* modify x by (U + D/omega)^{-1} */ 102 ierr = VecCopy(x,eis->b);CHKERRQ(ierr); 103 ierr = MatSOR(eis->A,eis->b,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),0.0,1,1,x);CHKERRQ(ierr); 104 pc->mat = eis->A; 105 PetscFunctionReturn(0); 106 } 107 108 #undef __FUNCT__ 109 #define __FUNCT__ "PCReset_Eisenstat" 110 static PetscErrorCode PCReset_Eisenstat(PC pc) 111 { 112 PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 113 PetscErrorCode ierr; 114 115 PetscFunctionBegin; 116 if (eis->b) {ierr = VecDestroy(eis->b);CHKERRQ(ierr);} 117 if (eis->shell) {ierr = MatDestroy(eis->shell);CHKERRQ(ierr);} 118 if (eis->diag) {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 PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 127 PetscErrorCode ierr; 128 129 PetscFunctionBegin; 130 ierr = PCReset_Eisentat(pc);CHKERRQ(ierr); 131 ierr = PetscFree(pc->data);CHKERRQ(ierr); 132 PetscFunctionReturn(0); 133 } 134 135 #undef __FUNCT__ 136 #define __FUNCT__ "PCSetFromOptions_Eisenstat" 137 static PetscErrorCode PCSetFromOptions_Eisenstat(PC pc) 138 { 139 PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 140 PetscErrorCode ierr; 141 PetscBool flg = PETSC_FALSE; 142 143 PetscFunctionBegin; 144 ierr = PetscOptionsHead("Eisenstat SSOR options");CHKERRQ(ierr); 145 ierr = PetscOptionsReal("-pc_eisenstat_omega","Relaxation factor 0 < omega < 2","PCEisenstatSetOmega",eis->omega,&eis->omega,PETSC_NULL);CHKERRQ(ierr); 146 ierr = PetscOptionsBool("-pc_eisenstat_no_diagonal_scaling","Do not use standard diagonal scaling","PCEisenstatNoDiagonalScaling",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 147 if (flg) { 148 ierr = PCEisenstatNoDiagonalScaling(pc);CHKERRQ(ierr); 149 } 150 ierr = PetscOptionsTail();CHKERRQ(ierr); 151 PetscFunctionReturn(0); 152 } 153 154 #undef __FUNCT__ 155 #define __FUNCT__ "PCView_Eisenstat" 156 static PetscErrorCode PCView_Eisenstat(PC pc,PetscViewer viewer) 157 { 158 PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 159 PetscErrorCode ierr; 160 PetscBool iascii; 161 162 PetscFunctionBegin; 163 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 164 if (iascii) { 165 ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: omega = %G\n",eis->omega);CHKERRQ(ierr); 166 if (eis->usediag) { 167 ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Using diagonal scaling (default)\n");CHKERRQ(ierr); 168 } else { 169 ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Not using diagonal scaling\n");CHKERRQ(ierr); 170 } 171 } else { 172 SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for Eisenstat PC",((PetscObject)viewer)->type_name); 173 } 174 PetscFunctionReturn(0); 175 } 176 177 #undef __FUNCT__ 178 #define __FUNCT__ "PCSetUp_Eisenstat" 179 static PetscErrorCode PCSetUp_Eisenstat(PC pc) 180 { 181 PetscErrorCode ierr; 182 PetscInt M,N,m,n; 183 PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 184 185 PetscFunctionBegin; 186 if (!pc->setupcalled) { 187 ierr = MatGetSize(pc->mat,&M,&N);CHKERRQ(ierr); 188 ierr = MatGetLocalSize(pc->mat,&m,&n);CHKERRQ(ierr); 189 ierr = MatCreate(((PetscObject)pc)->comm,&eis->shell);CHKERRQ(ierr); 190 ierr = MatSetSizes(eis->shell,m,n,M,N);CHKERRQ(ierr); 191 ierr = MatSetType(eis->shell,MATSHELL);CHKERRQ(ierr); 192 ierr = MatShellSetContext(eis->shell,(void*)pc);CHKERRQ(ierr); 193 ierr = PetscLogObjectParent(pc,eis->shell);CHKERRQ(ierr); 194 ierr = MatShellSetOperation(eis->shell,MATOP_MULT,(void(*)(void))PCMult_Eisenstat);CHKERRQ(ierr); 195 } 196 if (!eis->usediag) PetscFunctionReturn(0); 197 if (!pc->setupcalled) { 198 ierr = MatGetVecs(pc->pmat,&eis->diag,0);CHKERRQ(ierr); 199 ierr = PetscLogObjectParent(pc,eis->diag);CHKERRQ(ierr); 200 } 201 ierr = MatGetDiagonal(pc->pmat,eis->diag);CHKERRQ(ierr); 202 PetscFunctionReturn(0); 203 } 204 205 /* --------------------------------------------------------------------*/ 206 207 EXTERN_C_BEGIN 208 #undef __FUNCT__ 209 #define __FUNCT__ "PCEisenstatSetOmega_Eisenstat" 210 PetscErrorCode PCEisenstatSetOmega_Eisenstat(PC pc,PetscReal omega) 211 { 212 PC_Eisenstat *eis; 213 214 PetscFunctionBegin; 215 if (omega >= 2.0 || omega <= 0.0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range"); 216 eis = (PC_Eisenstat*)pc->data; 217 eis->omega = omega; 218 PetscFunctionReturn(0); 219 } 220 EXTERN_C_END 221 222 EXTERN_C_BEGIN 223 #undef __FUNCT__ 224 #define __FUNCT__ "PCEisenstatNoDiagonalScaling_Eisenstat" 225 PetscErrorCode PCEisenstatNoDiagonalScaling_Eisenstat(PC pc) 226 { 227 PC_Eisenstat *eis; 228 229 PetscFunctionBegin; 230 eis = (PC_Eisenstat*)pc->data; 231 eis->usediag = PETSC_FALSE; 232 PetscFunctionReturn(0); 233 } 234 EXTERN_C_END 235 236 #undef __FUNCT__ 237 #define __FUNCT__ "PCEisenstatSetOmega" 238 /*@ 239 PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega, 240 to use with Eisenstat's trick (where omega = 1.0 by default). 241 242 Logically Collective on PC 243 244 Input Parameters: 245 + pc - the preconditioner context 246 - omega - relaxation coefficient (0 < omega < 2) 247 248 Options Database Key: 249 . -pc_eisenstat_omega <omega> - Sets omega 250 251 Notes: 252 The Eisenstat trick implementation of SSOR requires about 50% of the 253 usual amount of floating point operations used for SSOR + Krylov method; 254 however, the preconditioned problem must be solved with both left 255 and right preconditioning. 256 257 To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner, 258 which can be chosen with the database options 259 $ -pc_type sor -pc_sor_symmetric 260 261 Level: intermediate 262 263 .keywords: PC, Eisenstat, set, SOR, SSOR, relaxation, omega 264 265 .seealso: PCSORSetOmega() 266 @*/ 267 PetscErrorCode PCEisenstatSetOmega(PC pc,PetscReal omega) 268 { 269 PetscErrorCode ierr; 270 271 PetscFunctionBegin; 272 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 273 PetscValidLogicalCollectiveReal(pc,omega,2); 274 ierr = PetscTryMethod(pc,"PCEisenstatSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr); 275 PetscFunctionReturn(0); 276 } 277 278 #undef __FUNCT__ 279 #define __FUNCT__ "PCEisenstatNoDiagonalScaling" 280 /*@ 281 PCEisenstatNoDiagonalScaling - Causes the Eisenstat preconditioner 282 not to do additional diagonal preconditioning. For matrices with a constant 283 along the diagonal, this may save a small amount of work. 284 285 Logically Collective on PC 286 287 Input Parameter: 288 . pc - the preconditioner context 289 290 Options Database Key: 291 . -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling() 292 293 Level: intermediate 294 295 Note: 296 If you use the KPSSetDiagonalScaling() or -ksp_diagonal_scale option then you will 297 likley want to use this routine since it will save you some unneeded flops. 298 299 .keywords: PC, Eisenstat, use, diagonal, scaling, SSOR 300 301 .seealso: PCEisenstatSetOmega() 302 @*/ 303 PetscErrorCode PCEisenstatNoDiagonalScaling(PC pc) 304 { 305 PetscErrorCode ierr; 306 307 PetscFunctionBegin; 308 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 309 ierr = PetscTryMethod(pc,"PCEisenstatNoDiagonalScaling_C",(PC),(pc));CHKERRQ(ierr); 310 PetscFunctionReturn(0); 311 } 312 313 /* --------------------------------------------------------------------*/ 314 315 /*MC 316 PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel) 317 preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed. 318 319 Options Database Keys: 320 + -pc_eisenstat_omega <omega> - Sets omega 321 - -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling() 322 323 Level: beginner 324 325 Concepts: SOR, preconditioners, Gauss-Seidel, Eisenstat's trick 326 327 Notes: Only implemented for the SeqAIJ matrix format. 328 Not a true parallel SOR, in parallel this implementation corresponds to block 329 Jacobi with SOR on each block. 330 331 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 332 PCEisenstatNoDiagonalScaling(), PCEisenstatSetOmega(), PCSOR 333 M*/ 334 335 EXTERN_C_BEGIN 336 #undef __FUNCT__ 337 #define __FUNCT__ "PCCreate_Eisenstat" 338 PetscErrorCode PCCreate_Eisenstat(PC pc) 339 { 340 PetscErrorCode ierr; 341 PC_Eisenstat *eis; 342 343 PetscFunctionBegin; 344 ierr = PetscNewLog(pc,PC_Eisenstat,&eis);CHKERRQ(ierr); 345 346 pc->ops->apply = PCApply_Eisenstat; 347 pc->ops->presolve = PCPreSolve_Eisenstat; 348 pc->ops->postsolve = PCPostSolve_Eisenstat; 349 pc->ops->applyrichardson = 0; 350 pc->ops->setfromoptions = PCSetFromOptions_Eisenstat; 351 pc->ops->destroy = PCDestroy_Eisenstat; 352 pc->ops->reset = PCReset_Eisenstat; 353 pc->ops->view = PCView_Eisenstat; 354 pc->ops->setup = PCSetUp_Eisenstat; 355 356 pc->data = (void*)eis; 357 eis->omega = 1.0; 358 eis->b = 0; 359 eis->diag = 0; 360 eis->usediag = PETSC_TRUE; 361 362 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatSetOmega_C","PCEisenstatSetOmega_Eisenstat", 363 PCEisenstatSetOmega_Eisenstat);CHKERRQ(ierr); 364 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C", 365 "PCEisenstatNoDiagonalScaling_Eisenstat", 366 PCEisenstatNoDiagonalScaling_Eisenstat);CHKERRQ(ierr); 367 PetscFunctionReturn(0); 368 } 369 EXTERN_C_END 370