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