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