1 2 /* -------------------------------------------------------------------- 3 4 This file implements a Jacobi preconditioner in PETSc as part of PC. 5 You can use this as a starting point for implementing your own 6 preconditioner that is not provided with PETSc. (You might also consider 7 just using PCSHELL) 8 9 The following basic routines are required for each preconditioner. 10 PCCreate_XXX() - Creates a preconditioner context 11 PCSetFromOptions_XXX() - Sets runtime options 12 PCApply_XXX() - Applies the preconditioner 13 PCDestroy_XXX() - Destroys the preconditioner context 14 where the suffix "_XXX" denotes a particular implementation, in 15 this case we use _Jacobi (e.g., PCCreate_Jacobi, PCApply_Jacobi). 16 These routines are actually called via the common user interface 17 routines PCCreate(), PCSetFromOptions(), PCApply(), and PCDestroy(), 18 so the application code interface remains identical for all 19 preconditioners. 20 21 Another key routine is: 22 PCSetUp_XXX() - Prepares for the use of a preconditioner 23 by setting data structures and options. The interface routine PCSetUp() 24 is not usually called directly by the user, but instead is called by 25 PCApply() if necessary. 26 27 Additional basic routines are: 28 PCView_XXX() - Prints details of runtime options that 29 have actually been used. 30 These are called by application codes via the interface routines 31 PCView(). 32 33 The various types of solvers (preconditioners, Krylov subspace methods, 34 nonlinear solvers, timesteppers) are all organized similarly, so the 35 above description applies to these categories also. One exception is 36 that the analogues of PCApply() for these components are KSPSolve(), 37 SNESSolve(), and TSSolve(). 38 39 Additional optional functionality unique to preconditioners is left and 40 right symmetric preconditioner application via PCApplySymmetricLeft() 41 and PCApplySymmetricRight(). The Jacobi implementation is 42 PCApplySymmetricLeftOrRight_Jacobi(). 43 44 -------------------------------------------------------------------- */ 45 46 /* 47 Include files needed for the Jacobi preconditioner: 48 pcimpl.h - private include file intended for use by all preconditioners 49 */ 50 51 #include <petsc-private/pcimpl.h> /*I "petscpc.h" I*/ 52 53 /* 54 Private context (data structure) for the Jacobi preconditioner. 55 */ 56 typedef struct { 57 Vec diag; /* vector containing the reciprocals of the diagonal elements of the preconditioner matrix */ 58 Vec diagsqrt; /* vector containing the reciprocals of the square roots of 59 the diagonal elements of the preconditioner matrix (used 60 only for symmetric preconditioner application) */ 61 PetscBool userowmax; 62 PetscBool userowsum; 63 PetscBool useabs; /* use the absolute values of the diagonal entries */ 64 } PC_Jacobi; 65 66 EXTERN_C_BEGIN 67 #undef __FUNCT__ 68 #define __FUNCT__ "PCJacobiSetUseRowMax_Jacobi" 69 PetscErrorCode PCJacobiSetUseRowMax_Jacobi(PC pc) 70 { 71 PC_Jacobi *j; 72 73 PetscFunctionBegin; 74 j = (PC_Jacobi*)pc->data; 75 j->userowmax = PETSC_TRUE; 76 PetscFunctionReturn(0); 77 } 78 EXTERN_C_END 79 80 EXTERN_C_BEGIN 81 #undef __FUNCT__ 82 #define __FUNCT__ "PCJacobiSetUseRowSum_Jacobi" 83 PetscErrorCode PCJacobiSetUseRowSum_Jacobi(PC pc) 84 { 85 PC_Jacobi *j; 86 87 PetscFunctionBegin; 88 j = (PC_Jacobi*)pc->data; 89 j->userowsum = PETSC_TRUE; 90 PetscFunctionReturn(0); 91 } 92 EXTERN_C_END 93 94 EXTERN_C_BEGIN 95 #undef __FUNCT__ 96 #define __FUNCT__ "PCJacobiSetUseAbs_Jacobi" 97 PetscErrorCode PCJacobiSetUseAbs_Jacobi(PC pc) 98 { 99 PC_Jacobi *j; 100 101 PetscFunctionBegin; 102 j = (PC_Jacobi*)pc->data; 103 j->useabs = PETSC_TRUE; 104 PetscFunctionReturn(0); 105 } 106 EXTERN_C_END 107 108 /* -------------------------------------------------------------------------- */ 109 /* 110 PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner 111 by setting data structures and options. 112 113 Input Parameter: 114 . pc - the preconditioner context 115 116 Application Interface Routine: PCSetUp() 117 118 Notes: 119 The interface routine PCSetUp() is not usually called directly by 120 the user, but instead is called by PCApply() if necessary. 121 */ 122 #undef __FUNCT__ 123 #define __FUNCT__ "PCSetUp_Jacobi" 124 static PetscErrorCode PCSetUp_Jacobi(PC pc) 125 { 126 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 127 Vec diag,diagsqrt; 128 PetscErrorCode ierr; 129 PetscInt n,i; 130 PetscScalar *x; 131 PetscBool zeroflag = PETSC_FALSE; 132 133 PetscFunctionBegin; 134 /* 135 For most preconditioners the code would begin here something like 136 137 if (pc->setupcalled == 0) { allocate space the first time this is ever called 138 ierr = MatGetVecs(pc->mat,&jac->diag);CHKERRQ(ierr); 139 PetscLogObjectParent(pc,jac->diag); 140 } 141 142 But for this preconditioner we want to support use of both the matrix' diagonal 143 elements (for left or right preconditioning) and square root of diagonal elements 144 (for symmetric preconditioning). Hence we do not allocate space here, since we 145 don't know at this point which will be needed (diag and/or diagsqrt) until the user 146 applies the preconditioner, and we don't want to allocate BOTH unless we need 147 them both. Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric() 148 and PCSetUp_Jacobi_Symmetric(), respectively. 149 */ 150 151 /* 152 Here we set up the preconditioner; that is, we copy the diagonal values from 153 the matrix and put them into a format to make them quick to apply as a preconditioner. 154 */ 155 diag = jac->diag; 156 diagsqrt = jac->diagsqrt; 157 158 if (diag) { 159 if (jac->userowmax) { 160 ierr = MatGetRowMaxAbs(pc->pmat,diag,PETSC_NULL);CHKERRQ(ierr); 161 } else if (jac->userowsum) { 162 ierr = MatGetRowSum(pc->pmat,diag);CHKERRQ(ierr); 163 } else { 164 ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr); 165 } 166 ierr = VecReciprocal(diag);CHKERRQ(ierr); 167 ierr = VecGetLocalSize(diag,&n);CHKERRQ(ierr); 168 ierr = VecGetArray(diag,&x);CHKERRQ(ierr); 169 if (jac->useabs) { 170 for (i=0; i<n; i++) x[i] = PetscAbsScalar(x[i]); 171 } 172 for (i=0; i<n; i++) { 173 if (x[i] == 0.0) { 174 x[i] = 1.0; 175 zeroflag = PETSC_TRUE; 176 } 177 } 178 ierr = VecRestoreArray(diag,&x);CHKERRQ(ierr); 179 } 180 if (diagsqrt) { 181 if (jac->userowmax) { 182 ierr = MatGetRowMaxAbs(pc->pmat,diagsqrt,PETSC_NULL);CHKERRQ(ierr); 183 } else if (jac->userowsum) { 184 ierr = MatGetRowSum(pc->pmat,diagsqrt);CHKERRQ(ierr); 185 } else { 186 ierr = MatGetDiagonal(pc->pmat,diagsqrt);CHKERRQ(ierr); 187 } 188 ierr = VecGetLocalSize(diagsqrt,&n);CHKERRQ(ierr); 189 ierr = VecGetArray(diagsqrt,&x);CHKERRQ(ierr); 190 for (i=0; i<n; i++) { 191 if (x[i] != 0.0) x[i] = 1.0/PetscSqrtReal(PetscAbsScalar(x[i])); 192 else { 193 x[i] = 1.0; 194 zeroflag = PETSC_TRUE; 195 } 196 } 197 ierr = VecRestoreArray(diagsqrt,&x);CHKERRQ(ierr); 198 } 199 if (zeroflag) { 200 ierr = PetscInfo(pc,"Zero detected in diagonal of matrix, using 1 at those locations\n");CHKERRQ(ierr); 201 } 202 PetscFunctionReturn(0); 203 } 204 /* -------------------------------------------------------------------------- */ 205 /* 206 PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the 207 inverse of the square root of the diagonal entries of the matrix. This 208 is used for symmetric application of the Jacobi preconditioner. 209 210 Input Parameter: 211 . pc - the preconditioner context 212 */ 213 #undef __FUNCT__ 214 #define __FUNCT__ "PCSetUp_Jacobi_Symmetric" 215 static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc) 216 { 217 PetscErrorCode ierr; 218 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 219 220 PetscFunctionBegin; 221 ierr = MatGetVecs(pc->pmat,&jac->diagsqrt,0);CHKERRQ(ierr); 222 ierr = PetscLogObjectParent(pc,jac->diagsqrt);CHKERRQ(ierr); 223 ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 224 PetscFunctionReturn(0); 225 } 226 /* -------------------------------------------------------------------------- */ 227 /* 228 PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the 229 inverse of the diagonal entries of the matrix. This is used for left of 230 right application of the Jacobi preconditioner. 231 232 Input Parameter: 233 . pc - the preconditioner context 234 */ 235 #undef __FUNCT__ 236 #define __FUNCT__ "PCSetUp_Jacobi_NonSymmetric" 237 static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc) 238 { 239 PetscErrorCode ierr; 240 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 241 242 PetscFunctionBegin; 243 ierr = MatGetVecs(pc->pmat,&jac->diag,0);CHKERRQ(ierr); 244 ierr = PetscLogObjectParent(pc,jac->diag);CHKERRQ(ierr); 245 ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 246 PetscFunctionReturn(0); 247 } 248 /* -------------------------------------------------------------------------- */ 249 /* 250 PCApply_Jacobi - Applies the Jacobi preconditioner to a vector. 251 252 Input Parameters: 253 . pc - the preconditioner context 254 . x - input vector 255 256 Output Parameter: 257 . y - output vector 258 259 Application Interface Routine: PCApply() 260 */ 261 #undef __FUNCT__ 262 #define __FUNCT__ "PCApply_Jacobi" 263 static PetscErrorCode PCApply_Jacobi(PC pc,Vec x,Vec y) 264 { 265 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 266 PetscErrorCode ierr; 267 268 PetscFunctionBegin; 269 if (!jac->diag) { 270 ierr = PCSetUp_Jacobi_NonSymmetric(pc);CHKERRQ(ierr); 271 } 272 ierr = VecPointwiseMult(y,x,jac->diag);CHKERRQ(ierr); 273 PetscFunctionReturn(0); 274 } 275 /* -------------------------------------------------------------------------- */ 276 /* 277 PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a 278 symmetric preconditioner to a vector. 279 280 Input Parameters: 281 . pc - the preconditioner context 282 . x - input vector 283 284 Output Parameter: 285 . y - output vector 286 287 Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight() 288 */ 289 #undef __FUNCT__ 290 #define __FUNCT__ "PCApplySymmetricLeftOrRight_Jacobi" 291 static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y) 292 { 293 PetscErrorCode ierr; 294 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 295 296 PetscFunctionBegin; 297 if (!jac->diagsqrt) { 298 ierr = PCSetUp_Jacobi_Symmetric(pc);CHKERRQ(ierr); 299 } 300 VecPointwiseMult(y,x,jac->diagsqrt); 301 PetscFunctionReturn(0); 302 } 303 /* -------------------------------------------------------------------------- */ 304 #undef __FUNCT__ 305 #define __FUNCT__ "PCReset_Jacobi" 306 static PetscErrorCode PCReset_Jacobi(PC pc) 307 { 308 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 309 PetscErrorCode ierr; 310 311 PetscFunctionBegin; 312 ierr = VecDestroy(&jac->diag);CHKERRQ(ierr); 313 ierr = VecDestroy(&jac->diagsqrt);CHKERRQ(ierr); 314 PetscFunctionReturn(0); 315 } 316 317 /* 318 PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner 319 that was created with PCCreate_Jacobi(). 320 321 Input Parameter: 322 . pc - the preconditioner context 323 324 Application Interface Routine: PCDestroy() 325 */ 326 #undef __FUNCT__ 327 #define __FUNCT__ "PCDestroy_Jacobi" 328 static PetscErrorCode PCDestroy_Jacobi(PC pc) 329 { 330 PetscErrorCode ierr; 331 332 PetscFunctionBegin; 333 ierr = PCReset_Jacobi(pc);CHKERRQ(ierr); 334 335 /* 336 Free the private data structure that was hanging off the PC 337 */ 338 ierr = PetscFree(pc->data);CHKERRQ(ierr); 339 PetscFunctionReturn(0); 340 } 341 342 #undef __FUNCT__ 343 #define __FUNCT__ "PCSetFromOptions_Jacobi" 344 static PetscErrorCode PCSetFromOptions_Jacobi(PC pc) 345 { 346 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 347 PetscErrorCode ierr; 348 349 PetscFunctionBegin; 350 ierr = PetscOptionsHead("Jacobi options");CHKERRQ(ierr); 351 ierr = PetscOptionsBool("-pc_jacobi_rowmax","Use row maximums for diagonal","PCJacobiSetUseRowMax",jac->userowmax, 352 &jac->userowmax,PETSC_NULL);CHKERRQ(ierr); 353 ierr = PetscOptionsBool("-pc_jacobi_rowsum","Use row sums for diagonal","PCJacobiSetUseRowSum",jac->userowsum, 354 &jac->userowsum,PETSC_NULL);CHKERRQ(ierr); 355 ierr = PetscOptionsBool("-pc_jacobi_abs","Use absolute values of diagaonal entries","PCJacobiSetUseAbs",jac->useabs, 356 &jac->useabs,PETSC_NULL);CHKERRQ(ierr); 357 ierr = PetscOptionsTail();CHKERRQ(ierr); 358 PetscFunctionReturn(0); 359 } 360 361 /* -------------------------------------------------------------------------- */ 362 /* 363 PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi, 364 and sets this as the private data within the generic preconditioning 365 context, PC, that was created within PCCreate(). 366 367 Input Parameter: 368 . pc - the preconditioner context 369 370 Application Interface Routine: PCCreate() 371 */ 372 373 /*MC 374 PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning) 375 376 Options Database Key: 377 + -pc_jacobi_rowmax - use the maximum absolute value in each row as the scaling factor, 378 rather than the diagonal 379 . -pc_jacobi_rowsum - use the row sums as the scaling factor, 380 rather than the diagonal 381 - -pc_jacobi_abs - use the absolute value of the diagaonl entry 382 383 Level: beginner 384 385 Concepts: Jacobi, diagonal scaling, preconditioners 386 387 Notes: By using KSPSetPCSide(ksp,PC_SYMMETRIC) or -ksp_pc_side symmetric 388 can scale each side of the matrix by the squareroot of the diagonal entries. 389 390 Zero entries along the diagonal are replaced with the value 1.0 391 392 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 393 PCJacobiSetUseRowMax(), PCJacobiSetUseRowSum(), PCJacobiSetUseAbs() 394 M*/ 395 396 EXTERN_C_BEGIN 397 #undef __FUNCT__ 398 #define __FUNCT__ "PCCreate_Jacobi" 399 PetscErrorCode PCCreate_Jacobi(PC pc) 400 { 401 PC_Jacobi *jac; 402 PetscErrorCode ierr; 403 404 PetscFunctionBegin; 405 /* 406 Creates the private data structure for this preconditioner and 407 attach it to the PC object. 408 */ 409 ierr = PetscNewLog(pc,PC_Jacobi,&jac);CHKERRQ(ierr); 410 pc->data = (void*)jac; 411 412 /* 413 Initialize the pointers to vectors to ZERO; these will be used to store 414 diagonal entries of the matrix for fast preconditioner application. 415 */ 416 jac->diag = 0; 417 jac->diagsqrt = 0; 418 jac->userowmax = PETSC_FALSE; 419 jac->userowsum = PETSC_FALSE; 420 jac->useabs = PETSC_FALSE; 421 422 /* 423 Set the pointers for the functions that are provided above. 424 Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 425 are called, they will automatically call these functions. Note we 426 choose not to provide a couple of these functions since they are 427 not needed. 428 */ 429 pc->ops->apply = PCApply_Jacobi; 430 pc->ops->applytranspose = PCApply_Jacobi; 431 pc->ops->setup = PCSetUp_Jacobi; 432 pc->ops->reset = PCReset_Jacobi; 433 pc->ops->destroy = PCDestroy_Jacobi; 434 pc->ops->setfromoptions = PCSetFromOptions_Jacobi; 435 pc->ops->view = 0; 436 pc->ops->applyrichardson = 0; 437 pc->ops->applysymmetricleft = PCApplySymmetricLeftOrRight_Jacobi; 438 pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi; 439 440 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseRowMax_C","PCJacobiSetUseRowMax_Jacobi",PCJacobiSetUseRowMax_Jacobi);CHKERRQ(ierr); 441 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseRowSum_C","PCJacobiSetUseRowSum_Jacobi",PCJacobiSetUseRowSum_Jacobi);CHKERRQ(ierr); 442 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseAbs_C","PCJacobiSetUseAbs_Jacobi",PCJacobiSetUseAbs_Jacobi);CHKERRQ(ierr); 443 PetscFunctionReturn(0); 444 } 445 EXTERN_C_END 446 447 448 #undef __FUNCT__ 449 #define __FUNCT__ "PCJacobiSetUseAbs" 450 /*@ 451 PCJacobiSetUseAbs - Causes the Jacobi preconditioner to use the 452 absolute value of the diagonal to for the preconditioner 453 454 Logically Collective on PC 455 456 Input Parameters: 457 . pc - the preconditioner context 458 459 460 Options Database Key: 461 . -pc_jacobi_abs 462 463 Level: intermediate 464 465 Concepts: Jacobi preconditioner 466 467 .seealso: PCJacobiaUseRowMax(), PCJacobiaUseRowSum() 468 469 @*/ 470 PetscErrorCode PCJacobiSetUseAbs(PC pc) 471 { 472 PetscErrorCode ierr; 473 474 PetscFunctionBegin; 475 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 476 ierr = PetscTryMethod(pc,"PCJacobiSetUseAbs_C",(PC),(pc));CHKERRQ(ierr); 477 PetscFunctionReturn(0); 478 } 479 480 #undef __FUNCT__ 481 #define __FUNCT__ "PCJacobiSetUseRowMax" 482 /*@ 483 PCJacobiSetUseRowMax - Causes the Jacobi preconditioner to use the 484 maximum entry in each row as the diagonal preconditioner, instead of 485 the diagonal entry 486 487 Logically Collective on PC 488 489 Input Parameters: 490 . pc - the preconditioner context 491 492 493 Options Database Key: 494 . -pc_jacobi_rowmax 495 496 Level: intermediate 497 498 Concepts: Jacobi preconditioner 499 500 .seealso: PCJacobiaUseAbs() 501 @*/ 502 PetscErrorCode PCJacobiSetUseRowMax(PC pc) 503 { 504 PetscErrorCode ierr; 505 506 PetscFunctionBegin; 507 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 508 ierr = PetscTryMethod(pc,"PCJacobiSetUseRowMax_C",(PC),(pc));CHKERRQ(ierr); 509 PetscFunctionReturn(0); 510 } 511 512 #undef __FUNCT__ 513 #define __FUNCT__ "PCJacobiSetUseRowSum" 514 /*@ 515 PCJacobiSetUseRowSum - Causes the Jacobi preconditioner to use the 516 sum of each row as the diagonal preconditioner, instead of 517 the diagonal entry 518 519 Logical Collective on PC 520 521 Input Parameters: 522 . pc - the preconditioner context 523 524 525 Options Database Key: 526 . -pc_jacobi_rowsum 527 528 Level: intermediate 529 530 Concepts: Jacobi preconditioner 531 532 .seealso: PCJacobiaUseAbs(), PCJacobiaUseRowSum() 533 @*/ 534 PetscErrorCode PCJacobiSetUseRowSum(PC pc) 535 { 536 PetscErrorCode ierr; 537 538 PetscFunctionBegin; 539 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 540 ierr = PetscTryMethod(pc,"PCJacobiSetUseRowSum_C",(PC),(pc));CHKERRQ(ierr); 541 PetscFunctionReturn(0); 542 } 543 544