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