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 int 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 int PCSetUp_Jacobi(PC pc) 96 { 97 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 98 Vec diag,diagsqrt; 99 int ierr,n,i,zeroflag=0; 100 PetscScalar *x; 101 102 PetscFunctionBegin; 103 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 166 PetscFunctionReturn(0); 167 } 168 /* -------------------------------------------------------------------------- */ 169 /* 170 PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the 171 inverse of the square root of the diagonal entries of the matrix. This 172 is used for symmetric application of the Jacobi preconditioner. 173 174 Input Parameter: 175 . pc - the preconditioner context 176 */ 177 #undef __FUNCT__ 178 #define __FUNCT__ "PCSetUp_Jacobi_Symmetric" 179 static int PCSetUp_Jacobi_Symmetric(PC pc) 180 { 181 int ierr; 182 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 183 184 PetscFunctionBegin; 185 186 ierr = MatGetVecs(pc->pmat,&jac->diagsqrt,0);CHKERRQ(ierr); 187 PetscLogObjectParent(pc,jac->diagsqrt); 188 ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 189 PetscFunctionReturn(0); 190 } 191 /* -------------------------------------------------------------------------- */ 192 /* 193 PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the 194 inverse of the diagonal entries of the matrix. This is used for left of 195 right application of the Jacobi preconditioner. 196 197 Input Parameter: 198 . pc - the preconditioner context 199 */ 200 #undef __FUNCT__ 201 #define __FUNCT__ "PCSetUp_Jacobi_NonSymmetric" 202 static int PCSetUp_Jacobi_NonSymmetric(PC pc) 203 { 204 int ierr; 205 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 206 207 PetscFunctionBegin; 208 209 ierr = MatGetVecs(pc->pmat,&jac->diag,0);CHKERRQ(ierr); 210 PetscLogObjectParent(pc,jac->diag); 211 ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 212 PetscFunctionReturn(0); 213 } 214 /* -------------------------------------------------------------------------- */ 215 /* 216 PCApply_Jacobi - Applies the Jacobi preconditioner to a vector. 217 218 Input Parameters: 219 . pc - the preconditioner context 220 . x - input vector 221 222 Output Parameter: 223 . y - output vector 224 225 Application Interface Routine: PCApply() 226 */ 227 #undef __FUNCT__ 228 #define __FUNCT__ "PCApply_Jacobi" 229 static int PCApply_Jacobi(PC pc,Vec x,Vec y) 230 { 231 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 232 int ierr; 233 234 PetscFunctionBegin; 235 if (!jac->diag) { 236 ierr = PCSetUp_Jacobi_NonSymmetric(pc);CHKERRQ(ierr); 237 } 238 ierr = VecPointwiseMult(x,jac->diag,y);CHKERRQ(ierr); 239 PetscFunctionReturn(0); 240 } 241 /* -------------------------------------------------------------------------- */ 242 /* 243 PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a 244 symmetric preconditioner to a vector. 245 246 Input Parameters: 247 . pc - the preconditioner context 248 . x - input vector 249 250 Output Parameter: 251 . y - output vector 252 253 Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight() 254 */ 255 #undef __FUNCT__ 256 #define __FUNCT__ "PCApplySymmetricLeftOrRight_Jacobi" 257 static int PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y) 258 { 259 int ierr; 260 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 261 262 PetscFunctionBegin; 263 if (!jac->diagsqrt) { 264 ierr = PCSetUp_Jacobi_Symmetric(pc);CHKERRQ(ierr); 265 } 266 VecPointwiseMult(x,jac->diagsqrt,y); 267 PetscFunctionReturn(0); 268 } 269 /* -------------------------------------------------------------------------- */ 270 /* 271 PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner 272 that was created with PCCreate_Jacobi(). 273 274 Input Parameter: 275 . pc - the preconditioner context 276 277 Application Interface Routine: PCDestroy() 278 */ 279 #undef __FUNCT__ 280 #define __FUNCT__ "PCDestroy_Jacobi" 281 static int PCDestroy_Jacobi(PC pc) 282 { 283 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 284 int ierr; 285 286 PetscFunctionBegin; 287 if (jac->diag) {ierr = VecDestroy(jac->diag);CHKERRQ(ierr);} 288 if (jac->diagsqrt) {ierr = VecDestroy(jac->diagsqrt);CHKERRQ(ierr);} 289 290 /* 291 Free the private data structure that was hanging off the PC 292 */ 293 ierr = PetscFree(jac);CHKERRQ(ierr); 294 PetscFunctionReturn(0); 295 } 296 297 #undef __FUNCT__ 298 #define __FUNCT__ "PCSetFromOptions_Jacobi" 299 static int PCSetFromOptions_Jacobi(PC pc) 300 { 301 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 302 int ierr; 303 304 PetscFunctionBegin; 305 ierr = PetscOptionsHead("Jacobi options");CHKERRQ(ierr); 306 ierr = PetscOptionsLogical("-pc_jacobi_rowmax","Use row maximums for diagonal","PCJacobiSetUseRowMax",jac->userowmax, 307 &jac->userowmax,PETSC_NULL);CHKERRQ(ierr); 308 ierr = PetscOptionsTail();CHKERRQ(ierr); 309 PetscFunctionReturn(0); 310 } 311 312 /* -------------------------------------------------------------------------- */ 313 /* 314 PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi, 315 and sets this as the private data within the generic preconditioning 316 context, PC, that was created within PCCreate(). 317 318 Input Parameter: 319 . pc - the preconditioner context 320 321 Application Interface Routine: PCCreate() 322 */ 323 324 /*MC 325 PCJacobi - Jacobi (i.e. diagonal scaling preconditioning) 326 327 Options Database Key: 328 . -pc_jacobi_rowmax - use the maximum absolute value in each row as the scaling factor, 329 rather than the diagonal 330 331 Level: beginner 332 333 Concepts: Jacobi, diagonal scaling, preconditioners 334 335 Notes: By using KSPSetPreconditionerSide(ksp,PC_SYMMETRIC) or -ksp_symmetric_pc you 336 can scale each side of the matrix by the squareroot of the diagonal entries. 337 338 Zero entries along the diagonal are replaced with the value 1.0 339 340 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 341 PCJacobiSetUseRowMax(), 342 M*/ 343 344 EXTERN_C_BEGIN 345 #undef __FUNCT__ 346 #define __FUNCT__ "PCCreate_Jacobi" 347 int PCCreate_Jacobi(PC pc) 348 { 349 PC_Jacobi *jac; 350 int ierr; 351 352 PetscFunctionBegin; 353 354 /* 355 Creates the private data structure for this preconditioner and 356 attach it to the PC object. 357 */ 358 ierr = PetscNew(PC_Jacobi,&jac);CHKERRQ(ierr); 359 pc->data = (void*)jac; 360 361 /* 362 Logs the memory usage; this is not needed but allows PETSc to 363 monitor how much memory is being used for various purposes. 364 */ 365 PetscLogObjectMemory(pc,sizeof(PC_Jacobi)); 366 367 /* 368 Initialize the pointers to vectors to ZERO; these will be used to store 369 diagonal entries of the matrix for fast preconditioner application. 370 */ 371 jac->diag = 0; 372 jac->diagsqrt = 0; 373 jac->userowmax = PETSC_FALSE; 374 375 /* 376 Set the pointers for the functions that are provided above. 377 Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 378 are called, they will automatically call these functions. Note we 379 choose not to provide a couple of these functions since they are 380 not needed. 381 */ 382 pc->ops->apply = PCApply_Jacobi; 383 pc->ops->applytranspose = PCApply_Jacobi; 384 pc->ops->setup = PCSetUp_Jacobi; 385 pc->ops->destroy = PCDestroy_Jacobi; 386 pc->ops->setfromoptions = PCSetFromOptions_Jacobi; 387 pc->ops->view = 0; 388 pc->ops->applyrichardson = 0; 389 pc->ops->applysymmetricleft = PCApplySymmetricLeftOrRight_Jacobi; 390 pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi; 391 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseRowMax_C","PCJacobiSetUseRowMax_Jacobi", 392 PCJacobiSetUseRowMax_Jacobi);CHKERRQ(ierr); 393 PetscFunctionReturn(0); 394 } 395 EXTERN_C_END 396 397 #undef __FUNCT__ 398 #define __FUNCT__ "PCJacobiSetUseRowMax" 399 /*@ 400 PCJacobiSetUseRowMax - Causes the Jacobi preconditioner to use the 401 maximum entry in each row as the diagonal preconditioner, instead of 402 the diagonal entry 403 404 Collective on PC 405 406 Input Parameters: 407 . pc - the preconditioner context 408 409 410 Options Database Key: 411 . -pc_jacobi_rowmax 412 413 Level: intermediate 414 415 Concepts: Jacobi preconditioner 416 417 @*/ 418 int PCJacobiSetUseRowMax(PC pc) 419 { 420 int ierr,(*f)(PC); 421 422 PetscFunctionBegin; 423 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 424 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCJacobiSetRowMax_C",(void (**)(void))&f);CHKERRQ(ierr); 425 if (f) { 426 ierr = (*f)(pc);CHKERRQ(ierr); 427 } 428 PetscFunctionReturn(0); 429 } 430 431