1 2 /* -------------------------------------------------------------------- 3 4 This file implements a Jacobi preconditioner for matrices that use 5 the Mat interface (various matrix formats). Actually, the only 6 matrix operation used here is MatGetDiagonal(), which extracts 7 diagonal elements of the preconditioning matrix. 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 "src/ksp/pc/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 58 of the preconditioner matrix */ 59 Vec diagsqrt; /* vector containing the reciprocals of the square roots of 60 the diagonal elements of the preconditioner matrix (used 61 only for symmetric preconditioner application) */ 62 PetscTruth userowmax; 63 } PC_Jacobi; 64 65 EXTERN_C_BEGIN 66 #undef __FUNCT__ 67 #define __FUNCT__ "PCJacobiSetUseRowMax_Jacobi" 68 PetscErrorCode PCJacobiSetUseRowMax_Jacobi(PC pc) 69 { 70 PC_Jacobi *j; 71 72 PetscFunctionBegin; 73 j = (PC_Jacobi*)pc->data; 74 j->userowmax = PETSC_TRUE; 75 PetscFunctionReturn(0); 76 } 77 EXTERN_C_END 78 79 /* -------------------------------------------------------------------------- */ 80 /* 81 PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner 82 by setting data structures and options. 83 84 Input Parameter: 85 . pc - the preconditioner context 86 87 Application Interface Routine: PCSetUp() 88 89 Notes: 90 The interface routine PCSetUp() is not usually called directly by 91 the user, but instead is called by PCApply() if necessary. 92 */ 93 #undef __FUNCT__ 94 #define __FUNCT__ "PCSetUp_Jacobi" 95 static PetscErrorCode PCSetUp_Jacobi(PC pc) 96 { 97 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 98 Vec diag,diagsqrt; 99 PetscErrorCode ierr; 100 PetscInt n,i,zeroflag=0; 101 PetscScalar *x; 102 103 PetscFunctionBegin; 104 /* 105 For most preconditioners the code would begin here something like 106 107 if (pc->setupcalled == 0) { allocate space the first time this is ever called 108 ierr = MatGetVecs(pc->mat,&jac->diag);CHKERRQ(ierr); 109 PetscLogObjectParent(pc,jac->diag); 110 } 111 112 But for this preconditioner we want to support use of both the matrix' diagonal 113 elements (for left or right preconditioning) and square root of diagonal elements 114 (for symmetric preconditioning). Hence we do not allocate space here, since we 115 don't know at this point which will be needed (diag and/or diagsqrt) until the user 116 applies the preconditioner, and we don't want to allocate BOTH unless we need 117 them both. Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric() 118 and PCSetUp_Jacobi_Symmetric(), respectively. 119 */ 120 121 /* 122 Here we set up the preconditioner; that is, we copy the diagonal values from 123 the matrix and put them into a format to make them quick to apply as a preconditioner. 124 */ 125 diag = jac->diag; 126 diagsqrt = jac->diagsqrt; 127 128 if (diag) { 129 if (jac->userowmax) { 130 ierr = MatGetRowMax(pc->pmat,diag);CHKERRQ(ierr); 131 } else { 132 ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr); 133 } 134 ierr = VecReciprocal(diag);CHKERRQ(ierr); 135 ierr = VecGetLocalSize(diag,&n);CHKERRQ(ierr); 136 ierr = VecGetArray(diag,&x);CHKERRQ(ierr); 137 for (i=0; i<n; i++) { 138 if (x[i] == 0.0) { 139 x[i] = 1.0; 140 zeroflag = 1; 141 } 142 } 143 ierr = VecRestoreArray(diag,&x);CHKERRQ(ierr); 144 } 145 if (diagsqrt) { 146 if (jac->userowmax) { 147 ierr = MatGetRowMax(pc->pmat,diagsqrt);CHKERRQ(ierr); 148 } else { 149 ierr = MatGetDiagonal(pc->pmat,diagsqrt);CHKERRQ(ierr); 150 } 151 ierr = VecGetLocalSize(diagsqrt,&n);CHKERRQ(ierr); 152 ierr = VecGetArray(diagsqrt,&x);CHKERRQ(ierr); 153 for (i=0; i<n; i++) { 154 if (x[i] != 0.0) x[i] = 1.0/sqrt(PetscAbsScalar(x[i])); 155 else { 156 x[i] = 1.0; 157 zeroflag = 1; 158 } 159 } 160 ierr = VecRestoreArray(diagsqrt,&x);CHKERRQ(ierr); 161 } 162 if (zeroflag) { 163 PetscLogInfo(pc,"PCSetUp_Jacobi:Zero detected in diagonal of matrix, using 1 at those locations\n"); 164 } 165 PetscFunctionReturn(0); 166 } 167 /* -------------------------------------------------------------------------- */ 168 /* 169 PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the 170 inverse of the square root of the diagonal entries of the matrix. This 171 is used for symmetric application of the Jacobi preconditioner. 172 173 Input Parameter: 174 . pc - the preconditioner context 175 */ 176 #undef __FUNCT__ 177 #define __FUNCT__ "PCSetUp_Jacobi_Symmetric" 178 static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc) 179 { 180 PetscErrorCode ierr; 181 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 182 183 PetscFunctionBegin; 184 ierr = MatGetVecs(pc->pmat,&jac->diagsqrt,0);CHKERRQ(ierr); 185 PetscLogObjectParent(pc,jac->diagsqrt); 186 ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 187 PetscFunctionReturn(0); 188 } 189 /* -------------------------------------------------------------------------- */ 190 /* 191 PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the 192 inverse of the diagonal entries of the matrix. This is used for left of 193 right application of the Jacobi preconditioner. 194 195 Input Parameter: 196 . pc - the preconditioner context 197 */ 198 #undef __FUNCT__ 199 #define __FUNCT__ "PCSetUp_Jacobi_NonSymmetric" 200 static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc) 201 { 202 PetscErrorCode ierr; 203 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 204 205 PetscFunctionBegin; 206 ierr = MatGetVecs(pc->pmat,&jac->diag,0);CHKERRQ(ierr); 207 PetscLogObjectParent(pc,jac->diag); 208 ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 209 PetscFunctionReturn(0); 210 } 211 /* -------------------------------------------------------------------------- */ 212 /* 213 PCApply_Jacobi - Applies the Jacobi preconditioner to a vector. 214 215 Input Parameters: 216 . pc - the preconditioner context 217 . x - input vector 218 219 Output Parameter: 220 . y - output vector 221 222 Application Interface Routine: PCApply() 223 */ 224 #undef __FUNCT__ 225 #define __FUNCT__ "PCApply_Jacobi" 226 static PetscErrorCode PCApply_Jacobi(PC pc,Vec x,Vec y) 227 { 228 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 229 PetscErrorCode ierr; 230 231 PetscFunctionBegin; 232 if (!jac->diag) { 233 ierr = PCSetUp_Jacobi_NonSymmetric(pc);CHKERRQ(ierr); 234 } 235 ierr = VecPointwiseMult(x,jac->diag,y);CHKERRQ(ierr); 236 PetscFunctionReturn(0); 237 } 238 /* -------------------------------------------------------------------------- */ 239 /* 240 PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a 241 symmetric preconditioner to a vector. 242 243 Input Parameters: 244 . pc - the preconditioner context 245 . x - input vector 246 247 Output Parameter: 248 . y - output vector 249 250 Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight() 251 */ 252 #undef __FUNCT__ 253 #define __FUNCT__ "PCApplySymmetricLeftOrRight_Jacobi" 254 static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y) 255 { 256 PetscErrorCode ierr; 257 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 258 259 PetscFunctionBegin; 260 if (!jac->diagsqrt) { 261 ierr = PCSetUp_Jacobi_Symmetric(pc);CHKERRQ(ierr); 262 } 263 VecPointwiseMult(x,jac->diagsqrt,y); 264 PetscFunctionReturn(0); 265 } 266 /* -------------------------------------------------------------------------- */ 267 /* 268 PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner 269 that was created with PCCreate_Jacobi(). 270 271 Input Parameter: 272 . pc - the preconditioner context 273 274 Application Interface Routine: PCDestroy() 275 */ 276 #undef __FUNCT__ 277 #define __FUNCT__ "PCDestroy_Jacobi" 278 static PetscErrorCode PCDestroy_Jacobi(PC pc) 279 { 280 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 281 PetscErrorCode ierr; 282 283 PetscFunctionBegin; 284 if (jac->diag) {ierr = VecDestroy(jac->diag);CHKERRQ(ierr);} 285 if (jac->diagsqrt) {ierr = VecDestroy(jac->diagsqrt);CHKERRQ(ierr);} 286 287 /* 288 Free the private data structure that was hanging off the PC 289 */ 290 ierr = PetscFree(jac);CHKERRQ(ierr); 291 PetscFunctionReturn(0); 292 } 293 294 #undef __FUNCT__ 295 #define __FUNCT__ "PCSetFromOptions_Jacobi" 296 static PetscErrorCode PCSetFromOptions_Jacobi(PC pc) 297 { 298 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 299 PetscErrorCode ierr; 300 301 PetscFunctionBegin; 302 ierr = PetscOptionsHead("Jacobi options");CHKERRQ(ierr); 303 ierr = PetscOptionsLogical("-pc_jacobi_rowmax","Use row maximums for diagonal","PCJacobiSetUseRowMax",jac->userowmax, 304 &jac->userowmax,PETSC_NULL);CHKERRQ(ierr); 305 ierr = PetscOptionsTail();CHKERRQ(ierr); 306 PetscFunctionReturn(0); 307 } 308 309 /* -------------------------------------------------------------------------- */ 310 /* 311 PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi, 312 and sets this as the private data within the generic preconditioning 313 context, PC, that was created within PCCreate(). 314 315 Input Parameter: 316 . pc - the preconditioner context 317 318 Application Interface Routine: PCCreate() 319 */ 320 321 /*MC 322 PCJacobi - Jacobi (i.e. diagonal scaling preconditioning) 323 324 Options Database Key: 325 . -pc_jacobi_rowmax - use the maximum absolute value in each row as the scaling factor, 326 rather than the diagonal 327 328 Level: beginner 329 330 Concepts: Jacobi, diagonal scaling, preconditioners 331 332 Notes: By using KSPSetPreconditionerSide(ksp,PC_SYMMETRIC) or -ksp_symmetric_pc you 333 can scale each side of the matrix by the squareroot of the diagonal entries. 334 335 Zero entries along the diagonal are replaced with the value 1.0 336 337 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 338 PCJacobiSetUseRowMax(), 339 M*/ 340 341 EXTERN_C_BEGIN 342 #undef __FUNCT__ 343 #define __FUNCT__ "PCCreate_Jacobi" 344 PetscErrorCode PCCreate_Jacobi(PC pc) 345 { 346 PC_Jacobi *jac; 347 PetscErrorCode ierr; 348 349 PetscFunctionBegin; 350 /* 351 Creates the private data structure for this preconditioner and 352 attach it to the PC object. 353 */ 354 ierr = PetscNew(PC_Jacobi,&jac);CHKERRQ(ierr); 355 pc->data = (void*)jac; 356 357 /* 358 Logs the memory usage; this is not needed but allows PETSc to 359 monitor how much memory is being used for various purposes. 360 */ 361 PetscLogObjectMemory(pc,sizeof(PC_Jacobi)); 362 363 /* 364 Initialize the pointers to vectors to ZERO; these will be used to store 365 diagonal entries of the matrix for fast preconditioner application. 366 */ 367 jac->diag = 0; 368 jac->diagsqrt = 0; 369 jac->userowmax = PETSC_FALSE; 370 371 /* 372 Set the pointers for the functions that are provided above. 373 Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 374 are called, they will automatically call these functions. Note we 375 choose not to provide a couple of these functions since they are 376 not needed. 377 */ 378 pc->ops->apply = PCApply_Jacobi; 379 pc->ops->applytranspose = PCApply_Jacobi; 380 pc->ops->setup = PCSetUp_Jacobi; 381 pc->ops->destroy = PCDestroy_Jacobi; 382 pc->ops->setfromoptions = PCSetFromOptions_Jacobi; 383 pc->ops->view = 0; 384 pc->ops->applyrichardson = 0; 385 pc->ops->applysymmetricleft = PCApplySymmetricLeftOrRight_Jacobi; 386 pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi; 387 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseRowMax_C","PCJacobiSetUseRowMax_Jacobi", 388 PCJacobiSetUseRowMax_Jacobi);CHKERRQ(ierr); 389 PetscFunctionReturn(0); 390 } 391 EXTERN_C_END 392 393 #undef __FUNCT__ 394 #define __FUNCT__ "PCJacobiSetUseRowMax" 395 /*@ 396 PCJacobiSetUseRowMax - Causes the Jacobi preconditioner to use the 397 maximum entry in each row as the diagonal preconditioner, instead of 398 the diagonal entry 399 400 Collective on PC 401 402 Input Parameters: 403 . pc - the preconditioner context 404 405 406 Options Database Key: 407 . -pc_jacobi_rowmax 408 409 Level: intermediate 410 411 Concepts: Jacobi preconditioner 412 413 @*/ 414 PetscErrorCode PCJacobiSetUseRowMax(PC pc) 415 { 416 PetscErrorCode ierr,(*f)(PC); 417 418 PetscFunctionBegin; 419 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 420 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCJacobiSetRowMax_C",(void (**)(void))&f);CHKERRQ(ierr); 421 if (f) { 422 ierr = (*f)(pc);CHKERRQ(ierr); 423 } 424 PetscFunctionReturn(0); 425 } 426 427