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