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