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