xref: /petsc/doc/manual/ksp.md (revision ac84dfd5778759083efa0c46d3820bac8a11500e)
1(ch_ksp)=
2
3# KSP: Linear System Solvers
4
5The `KSP` object is the heart of PETSc, because it provides uniform
6and efficient access to all of the package’s linear system solvers,
7including parallel and sequential, direct and iterative. `KSP` is
8intended for solving systems of the form
9
10$$
11A x = b,
12$$ (eq_axeqb)
13
14where $A$ denotes the matrix representation of a linear operator,
15$b$ is the right-hand-side vector, and $x$ is the solution
16vector. `KSP` uses the same calling sequence for both direct and
17iterative solution of a linear system. In addition, particular solution
18techniques and their associated options can be selected at runtime.
19
20The combination of a Krylov subspace method and a preconditioner is at
21the center of most modern numerical codes for the iterative solution of
22linear systems. Many textbooks (e.g. {cite}`fgn` {cite}`vandervorst2003`, or {cite}`saad2003`) provide an
23overview of the theory of such methods.
24The `KSP` package, discussed in
25{any}`sec_ksp`, provides many popular Krylov subspace
26iterative methods; the `PC` module, described in
27{any}`sec_pc`, includes a variety of preconditioners.
28
29(sec_usingksp)=
30
31## Using KSP
32
33To solve a linear system with `KSP`, one must first create a solver
34context with the command
35
36```
37KSPCreate(MPI_Comm comm,KSP *ksp);
38```
39
40Here `comm` is the MPI communicator and `ksp` is the newly formed
41solver context. Before actually solving a linear system with `KSP`,
42the user must call the following routine to set the matrices associated
43with the linear system:
44
45```
46KSPSetOperators(KSP ksp,Mat Amat,Mat Pmat);
47```
48
49The argument `Amat`, representing the matrix that defines the linear
50system, is a symbolic placeholder for any kind of matrix or operator. In
51particular, `KSP` *does* support matrix-free methods. The routine
52`MatCreateShell()` in {any}`sec_matrixfree`
53provides further information regarding matrix-free methods. Typically,
54the matrix from which the preconditioner is to be constructed, `Pmat`,
55is the same as the matrix that defines the linear system, `Amat`;
56however, occasionally these matrices differ (for instance, when a
57matrix used to compute the preconditioner is obtained from a lower order method than that
58employed to form the linear system matrix).
59
60Much of the power of `KSP` can be accessed through the single routine
61
62```
63KSPSetFromOptions(KSP ksp);
64```
65
66This routine accepts the option `-help` as well as any of
67the `KSP` and `PC` options discussed below. To solve a linear
68system, one sets the right hand size and solution vectors using the
69command
70
71```
72KSPSolve(KSP ksp,Vec b,Vec x);
73```
74
75where `b` and `x` respectively denote the right-hand side and
76solution vectors. On return, the iteration number at which the iterative
77process stopped can be obtained using
78
79```
80KSPGetIterationNumber(KSP ksp, PetscInt *its);
81```
82
83Note that this does not state that the method converged at this
84iteration: it can also have reached the maximum number of iterations, or
85have diverged.
86
87{any}`sec_convergencetests` gives more details
88regarding convergence testing. Note that multiple linear solves can be
89performed by the same `KSP` context. Once the `KSP` context is no
90longer needed, it should be destroyed with the command
91
92```
93KSPDestroy(KSP *ksp);
94```
95
96The above procedure is sufficient for general use of the `KSP`
97package. One additional step is required for users who wish to customize
98certain preconditioners (e.g., see {any}`sec_bjacobi`) or
99to log certain performance data using the PETSc profiling facilities (as
100discussed in {any}`ch_profiling`). In this case, the user can
101optionally explicitly call
102
103```
104KSPSetUp(KSP ksp);
105```
106
107before calling `KSPSolve()` to perform any setup required for the
108linear solvers. The explicit call of this routine enables the separate
109profiling of any computations performed during the set up phase, such
110as incomplete factorization for the ILU preconditioner.
111
112The default solver within `KSP` is restarted GMRES, `KSPGMRES`, preconditioned for
113the uniprocess case with ILU(0), and for the multiprocess case with the
114block Jacobi method (with one block per process, each of which is solved
115with ILU(0)). A variety of other solvers and options are also available.
116To allow application programmers to set any of the preconditioner or
117Krylov subspace options directly within the code, we provide routines
118that extract the `PC` and `KSP` contexts,
119
120```
121KSPGetPC(KSP ksp,PC *pc);
122```
123
124The application programmer can then directly call any of the `PC` or
125`KSP` routines to modify the corresponding default options.
126
127To solve a linear system with a direct solver (supported by
128PETSc for sequential matrices, and by several external solvers through
129PETSc interfaces, see {any}`sec_externalsol`) one may use
130the options `-ksp_type` `preonly` (or the equivalent `-ksp_type` `none`)
131`-pc_type` `lu` or `-pc_type` `cholesky` (see below).
132
133By default, if a direct solver is used, the factorization is *not* done
134in-place. This approach prevents the user from the unexpected surprise
135of having a corrupted matrix after a linear solve. The routine
136`PCFactorSetUseInPlace()`, discussed below, causes factorization to be
137done in-place.
138
139## Solving Successive Linear Systems
140
141When solving multiple linear systems of the same size with the same
142method, several options are available. To solve successive linear
143systems having the *same* matrix from which to construct the preconditioner (i.e., the same data
144structure with exactly the same matrix elements) but different
145right-hand-side vectors, the user should simply call `KSPSolve()`
146multiple times. The preconditioner setup operations (e.g., factorization
147for ILU) will be done during the first call to `KSPSolve()` only; such
148operations will *not* be repeated for successive solves.
149
150To solve successive linear systems that have *different* matrix values, because you
151have changed the matrix values in the `Mat` objects you passed to `KSPSetOperators()`,
152still simply call `KPSSolve()`. In this case the preconditioner will be recomputed
153automatically. Use the option `-ksp_reuse_preconditioner true`, or call
154`KSPSetReusePreconditioner()`, to reuse the previously computed preconditioner.
155For many problems, if the matrix changes values only slightly, reusing the
156old preconditioner can be more efficient.
157
158If you wish to reuse the `KSP` with a different sized matrix and vectors, you must
159call `KSPReset()` before calling `KSPSetOperators()` with the new matrix.
160
161(sec_ksp)=
162
163## Krylov Methods
164
165The Krylov subspace methods accept a number of options, many of which
166are discussed below. First, to set the Krylov subspace method that is to
167be used, one calls the command
168
169```
170KSPSetType(KSP ksp,KSPType method);
171```
172
173The type can be one of `KSPRICHARDSON`, `KSPCHEBYSHEV`, `KSPCG`,
174`KSPGMRES`, `KSPTCQMR`, `KSPBCGS`, `KSPCGS`, `KSPTFQMR`,
175`KSPCR`, `KSPLSQR`, `KSPBICG`, `KSPPREONLY` (or the equivalent `KSPNONE`), or others; see
176{any}`tab-kspdefaults` or the `KSPType` man page for more.
177The `KSP` method can also be set with the options database command
178`-ksp_type`, followed by one of the options `richardson`,
179`chebyshev`, `cg`, `gmres`, `tcqmr`, `bcgs`, `cgs`,
180`tfqmr`, `cr`, `lsqr`, `bicg`, `preonly` (or the equivalent `none`), or others (see
181{any}`tab-kspdefaults` or the `KSPType` man page). There are
182method-specific options. For instance, for the Richardson, Chebyshev, and
183GMRES methods:
184
185```
186KSPRichardsonSetScale(KSP ksp,PetscReal scale);
187KSPChebyshevSetEigenvalues(KSP ksp,PetscReal emax,PetscReal emin);
188KSPGMRESSetRestart(KSP ksp,PetscInt max_steps);
189```
190
191The default parameter values are
192`scale=1.0, emax=0.01, emin=100.0`, and `max_steps=30`. The
193GMRES restart and Richardson damping factor can also be set with the
194options `-ksp_gmres_restart <n>` and
195`-ksp_richardson_scale <factor>`.
196
197The default technique for orthogonalization of the Krylov vectors in
198GMRES is the unmodified (classical) Gram-Schmidt method, which can be
199set with
200
201```
202KSPGMRESSetOrthogonalization(KSP ksp,KSPGMRESClassicalGramSchmidtOrthogonalization);
203```
204
205or the options database command `-ksp_gmres_classicalgramschmidt`. By
206default this will *not* use iterative refinement to improve the
207stability of the orthogonalization. This can be changed with the option
208
209```
210KSPGMRESSetCGSRefinementType(KSP ksp,KSPGMRESCGSRefinementType type)
211```
212
213or via the options database with
214
215```
216-ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always>
217```
218
219The values for `KSPGMRESCGSRefinementType()` are
220`KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_IFNEEDED`
221and `KSP_GMRES_CGS_REFINE_ALWAYS`.
222
223One can also use modified Gram-Schmidt, by using the orthogonalization
224routine `KSPGMRESModifiedGramSchmidtOrthogonalization()` or by using
225the command line option `-ksp_gmres_modifiedgramschmidt`.
226
227For the conjugate gradient method with complex numbers, there are two
228slightly different algorithms depending on whether the matrix is
229Hermitian symmetric or truly symmetric (the default is to assume that it
230is Hermitian symmetric). To indicate that it is symmetric, one uses the
231command
232
233```
234KSPCGSetType(ksp,KSP_CG_SYMMETRIC);
235```
236
237Note that this option is not valid for all matrices.
238
239Some `KSP` types do not support preconditioning. For instance,
240the CGLS algorithm does not involve a preconditioner; any preconditioner
241set to work with the `KSP` object is ignored if `KSPCGLS` was
242selected.
243
244By default, `KSP` assumes an initial guess of zero by zeroing the
245initial value for the solution vector that is given; this zeroing is
246done at the call to `KSPSolve()`. To use a nonzero initial guess, the
247user *must* call
248
249```
250KSPSetInitialGuessNonzero(KSP ksp,PetscBool flg);
251```
252
253(sec_ksppc)=
254
255### Preconditioning within KSP
256
257Since the rate of convergence of Krylov projection methods for a
258particular linear system is strongly dependent on its spectrum,
259preconditioning is typically used to alter the spectrum and hence
260accelerate the convergence rate of iterative techniques. Preconditioning
261can be applied to the system {eq}`eq_axeqb` by
262
263$$
264(M_L^{-1} A M_R^{-1}) \, (M_R x) = M_L^{-1} b,
265$$ (eq_prec)
266
267where $M_L$ and $M_R$ indicate preconditioning matrices (or,
268matrices from which the preconditioner is to be constructed). If
269$M_L = I$ in {eq}`eq_prec`, right preconditioning
270results, and the residual of {eq}`eq_axeqb`,
271
272$$
273r \equiv b - Ax = b - A M_R^{-1} \, M_R x,
274$$
275
276is preserved. In contrast, the residual is altered for left
277($M_R = I$) and symmetric preconditioning, as given by
278
279$$
280r_L \equiv M_L^{-1} b - M_L^{-1} A x = M_L^{-1} r.
281$$
282
283By default, most KSP implementations use left preconditioning. Some more
284naturally use other options, though. For instance, `KSPQCG` defaults
285to use symmetric preconditioning and `KSPFGMRES` uses right
286preconditioning by default. Right preconditioning can be activated for
287some methods by using the options database command
288`-ksp_pc_side right` or calling the routine
289
290```
291KSPSetPCSide(ksp,PC_RIGHT);
292```
293
294Attempting to use right preconditioning for a method that does not
295currently support it results in an error message of the form
296
297```none
298KSPSetUp_Richardson:No right preconditioning for KSPRICHARDSON
299```
300
301```{eval-rst}
302.. list-table:: KSP Objects
303  :name: tab-kspdefaults
304  :header-rows: 1
305
306  * - Method
307    - KSPType
308    - Options Database
309  * - Richardson
310    - ``KSPRICHARDSON``
311    - ``richardson``
312  * - Chebyshev
313    - ``KSPCHEBYSHEV``
314    - ``chebyshev``
315  * - Conjugate Gradient :cite:`hs:52`
316    - ``KSPCG``
317    - ``cg``
318  * - Pipelined Conjugate Gradients :cite:`ghyselsvanroose2014`
319    - ``KSPPIPECG``
320    - ``pipecg``
321  * - Pipelined Conjugate Gradients (Gropp)
322    - ``KSPGROPPCG``
323    - ``groppcg``
324  * - Pipelined Conjugate Gradients with Residual Replacement
325    - ``KSPPIPECGRR``
326    - ``pipecgrr``
327  * - Conjugate Gradients for the Normal Equations
328    - ``KSPCGNE``
329    - ``cgne``
330  * - Flexible Conjugate Gradients :cite:`flexiblecg`
331    - ``KSPFCG``
332    - ``fcg``
333  * -  Pipelined, Flexible Conjugate Gradients :cite:`sananschneppmay2016`
334    - ``KSPPIPEFCG``
335    - ``pipefcg``
336  * - Conjugate Gradients for Least Squares
337    - ``KSPCGLS``
338    - ``cgls``
339  * - Conjugate Gradients with Constraint (1)
340    - ``KSPNASH``
341    - ``nash``
342  * - Conjugate Gradients with Constraint (2)
343    - ``KSPSTCG``
344    - ``stcg``
345  * - Conjugate Gradients with Constraint (3)
346    - ``KSPGLTR``
347    - ``gltr``
348  * - Conjugate Gradients with Constraint (4)
349    - ``KSPQCG``
350    - ``qcg``
351  * - BiConjugate Gradient
352    - ``KSPBICG``
353    - ``bicg``
354  * - BiCGSTAB :cite:`v:92`
355    - ``KSPBCGS``
356    - ``bcgs``
357  * - Improved BiCGSTAB
358    - ``KSPIBCGS``
359    - ``ibcgs``
360  * - QMRCGSTAB :cite:`chan1994qmrcgs`
361    - ``KSPQMRCGS``
362    - ``qmrcgs``
363  * - Flexible BiCGSTAB
364    - ``KSPFBCGS``
365    - ``fbcgs``
366  * - Flexible BiCGSTAB (variant)
367    - ``KSPFBCGSR``
368    - ``fbcgsr``
369  * - Enhanced BiCGSTAB(L)
370    - ``KSPBCGSL``
371    - ``bcgsl``
372  * - Minimal Residual Method :cite:`paige.saunders:solution`
373    - ``KSPMINRES``
374    - ``minres``
375  * - Generalized Minimal Residual :cite:`saad.schultz:gmres`
376    - ``KSPGMRES``
377    - ``gmres``
378  * - Flexible Generalized Minimal Residual :cite:`saad1993`
379    - ``KSPFGMRES``
380    - ``fgmres``
381  * - Deflated Generalized Minimal Residual
382    - ``KSPDGMRES``
383    - ``dgmres``
384  * - Pipelined Generalized Minimal Residual :cite:`ghyselsashbymeerbergenvanroose2013`
385    - ``KSPPGMRES``
386    - ``pgmres``
387  * - Pipelined, Flexible Generalized Minimal Residual :cite:`sananschneppmay2016`
388    - ``KSPPIPEFGMRES``
389    - ``pipefgmres``
390  * - Generalized Minimal Residual with Accelerated Restart
391    - ``KSPLGMRES``
392    - ``lgmres``
393  * - Conjugate Residual :cite:`eisenstat1983variational`
394    - ``KSPCR``
395    - ``cr``
396  * - Generalized Conjugate Residual
397    - ``KSPGCR``
398    - ``gcr``
399  * - Pipelined Conjugate Residual
400    - ``KSPPIPECR``
401    - ``pipecr``
402  * - Pipelined, Flexible Conjugate Residual :cite:`sananschneppmay2016`
403    - ``KSPPIPEGCR``
404    - ``pipegcr``
405  * - FETI-DP
406    - ``KSPFETIDP``
407    - ``fetidp``
408  * - Conjugate Gradient Squared :cite:`so:89`
409    - ``KSPCGS``
410    - ``cgs``
411  * - Transpose-Free Quasi-Minimal Residual (1) :cite:`f:93`
412    - ``KSPTFQMR``
413    - ``tfqmr``
414  * - Transpose-Free Quasi-Minimal Residual (2)
415    - ``KSPTCQMR``
416    - ``tcqmr``
417  * - Least Squares Method
418    - ``KSPLSQR``
419    - ``lsqr``
420  * - Symmetric LQ Method :cite:`paige.saunders:solution`
421    - ``KSPSYMMLQ``
422    - ``symmlq``
423  * - TSIRM
424    - ``KSPTSIRM``
425    - ``tsirm``
426  * - Python Shell
427    - ``KSPPYTHON``
428    - ``python``
429  * - Shell for no ``KSP`` method
430    - ``KSPNONE``
431    - ``none``
432
433```
434
435Note: the bi-conjugate gradient method requires application of both the
436matrix and its transpose plus the preconditioner and its transpose.
437Currently not all matrices and preconditioners provide this support and
438thus the `KSPBICG` cannot always be used.
439
440Note: PETSc implements the FETI-DP (Finite Element Tearing and
441Interconnecting Dual-Primal) method as an implementation of `KSP` since it recasts the
442original problem into a constrained minimization one with Lagrange
443multipliers. The only matrix type supported is `MATIS`. Support for
444saddle point problems is provided. See the man page for `KSPFETIDP` for
445further details.
446
447(sec_convergencetests)=
448
449### Convergence Tests
450
451The default convergence test, `KSPConvergedDefault()`, uses the \$ l_2 \$ norm of the preconditioned \$ B(b - A x) \$ or unconditioned residual \$ b - Ax\$, depending on the `KSPType` and the value of `KSPNormType` set with `KSPSetNormType`. For `KSPCG` and `KSPGMRES` the default is the norm of the preconditioned residual.
452The preconditioned residual is used by default for
453convergence testing of all left-preconditioned `KSP` methods. For the
454conjugate gradient, Richardson, and Chebyshev methods the true residual
455can be used by the options database command
456`-ksp_norm_type unpreconditioned` or by calling the routine
457
458```
459KSPSetNormType(ksp, KSP_NORM_UNPRECONDITIONED);
460```
461
462`KSPCG` also supports using the natural norm induced by the symmetric positive-definite
463matrix that defines the linear system with the options database command `-ksp_norm_type natural` or by calling the routine
464
465```
466KSPSetNormType(ksp, KSP_NORM_NATURAL);
467```
468
469Convergence (or divergence) is decided
470by three quantities: the decrease of the residual norm relative to the
471norm of the right-hand side, `rtol`, the absolute size of the residual
472norm, `atol`, and the relative increase in the residual, `dtol`.
473Convergence is detected at iteration $k$ if
474
475$$
476\| r_k \|_2 < {\rm max} ( \text{rtol} * \| b \|_2, \text{atol}),
477$$
478
479where $r_k = b - A x_k$. Divergence is detected if
480
481$$
482\| r_k \|_2 > \text{dtol} * \| b \|_2.
483$$
484
485These parameters, as well as the maximum number of allowable iterations,
486can be set with the routine
487
488```
489KSPSetTolerances(KSP ksp,PetscReal rtol,PetscReal atol,PetscReal dtol,PetscInt maxits);
490```
491
492The user can retain the current value of any of these parameters by
493specifying `PETSC_CURRENT` as the corresponding tolerance; the
494defaults are `rtol=1e-5`, `atol=1e-50`, `dtol=1e5`, and
495`maxits=1e4`. Using `PETSC_DETERMINE` will set the parameters back to their
496initial values when the object's type was set. These parameters can also be set from the options
497database with the commands `-ksp_rtol` `<rtol>`, `-ksp_atol`
498`<atol>`, `-ksp_divtol` `<dtol>`, and `-ksp_max_it` `<its>`.
499
500In addition to providing an interface to a simple convergence test,
501`KSP` allows the application programmer the flexibility to provide
502customized convergence-testing routines. The user can specify a
503customized routine with the command
504
505```
506KSPSetConvergenceTest(KSP ksp,PetscErrorCode (*test)(KSP ksp,PetscInt it,PetscReal rnorm, KSPConvergedReason *reason,void *ctx),void *ctx,PetscErrorCode (*destroy)(void *ctx));
507```
508
509The final routine argument, `ctx`, is an optional context for private
510data for the user-defined convergence routine, `test`. Other `test`
511routine arguments are the iteration number, `it`, and the residual’s
512norm, `rnorm`. The routine for detecting convergence,
513`test`, should set `reason` to positive for convergence, 0 for no
514convergence, and negative for failure to converge. A full list of
515possible values is given in the `KSPConvergedReason` manual page.
516You can use `KSPGetConvergedReason()` after
517`KSPSolve()` to see why convergence/divergence was detected.
518
519(sec_kspmonitor)=
520
521### Convergence Monitoring
522
523By default, the Krylov solvers, `KSPSolve()`, run silently without displaying
524information about the iterations. The user can indicate that the norms
525of the residuals should be displayed at each iteration by using `-ksp_monitor` with
526the options database. To display the residual norms in a graphical
527window (running under X Windows), one should use
528`-ksp_monitor draw::draw_lg`. Application programmers can also
529provide their own routines to perform the monitoring by using the
530command
531
532```
533KSPMonitorSet(KSP ksp, PetscErrorCode (*mon)(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx), void *ctx, (PetscCtxDestroyFn *)mondestroy);
534```
535
536The final routine argument, `ctx`, is an optional context for private
537data for the user-defined monitoring routine, `mon`. Other `mon`
538routine arguments are the iteration number (`it`) and the residual’s
539norm (`rnorm`), as discussed above in {any}`sec_convergencetests`.
540A helpful routine within user-defined
541monitors is `PetscObjectGetComm((PetscObject)ksp,MPI_Comm *comm)`,
542which returns in `comm` the MPI communicator for the `KSP` context.
543See {any}`sec_writing` for more discussion of the use of
544MPI communicators within PETSc.
545
546Many monitoring routines are supplied with PETSc, including
547
548```
549KSPMonitorResidual(KSP,PetscInt,PetscReal, void *);
550KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
551KSPMonitorTrueResidual(KSP,PetscInt,PetscReal, void *);
552```
553
554The default monitor simply prints an estimate of a norm of
555the residual at each iteration. The routine
556`KSPMonitorSingularValue()` is appropriate only for use with the
557conjugate gradient method or GMRES, since it prints estimates of the
558extreme singular values of the preconditioned operator at each
559iteration computed via the Lanczos or Arnoldi algorithms.
560
561Since `KSPMonitorTrueResidual()` prints the true
562residual at each iteration by actually computing the residual using the
563formula $r = b - Ax$, the routine is slow and should be used only
564for testing or convergence studies, not for timing. These `KSPSolve()` monitors may
565be accessed with the command line options `-ksp_monitor`,
566`-ksp_monitor_singular_value`, and `-ksp_monitor_true_residual`.
567
568To employ the default graphical monitor, one should use the command
569`-ksp_monitor draw::draw_lg`.
570
571One can cancel hardwired monitoring routines for KSP at runtime with
572`-ksp_monitor_cancel`.
573
574### Understanding the Operator’s Spectrum
575
576Since the convergence of Krylov subspace methods depends strongly on the
577spectrum (eigenvalues) of the preconditioned operator, PETSc has
578specific routines for eigenvalue approximation via the Arnoldi or
579Lanczos iteration. First, before the linear solve one must call
580
581```
582KSPSetComputeEigenvalues(ksp,PETSC_TRUE);
583```
584
585Then after the `KSP` solve one calls
586
587```
588KSPComputeEigenvalues(KSP ksp,PetscInt n,PetscReal *realpart,PetscReal *complexpart,PetscInt *neig);
589```
590
591Here, `n` is the size of the two arrays and the eigenvalues are
592inserted into those two arrays. `neig` is the number of eigenvalues
593computed; this number depends on the size of the Krylov space generated
594during the linear system solution, for GMRES it is never larger than the
595`restart` parameter. There is an additional routine
596
597```
598KSPComputeEigenvaluesExplicitly(KSP ksp, PetscInt n,PetscReal *realpart,PetscReal *complexpart);
599```
600
601that is useful only for very small problems. It explicitly computes the
602full representation of the preconditioned operator and calls LAPACK to
603compute its eigenvalues. It should be only used for matrices of size up
604to a couple hundred. The `PetscDrawSP*()` routines are very useful for
605drawing scatter plots of the eigenvalues.
606
607The eigenvalues may also be computed and displayed graphically with the
608options data base commands `-ksp_view_eigenvalues draw` and
609`-ksp_view_eigenvalues_explicit draw`. Or they can be dumped to the
610screen in ASCII text via `-ksp_view_eigenvalues` and
611`-ksp_view_eigenvalues_explicit`.
612
613(sec_flexibleksp)=
614
615### Flexible Krylov Methods
616
617Standard Krylov methods require that the preconditioner be a linear operator, thus, for example, a standard `KSP` method
618cannot use a `KSP` in its preconditioner, as is common in the Block-Jacobi method `PCBJACOBI`, for example.
619Flexible Krylov methods are a subset of methods that allow (with modest additional requirements
620on memory) the preconditioner to be nonlinear. For example, they can be used with the `PCKSP` preconditioner.
621The flexible `KSP` methods have the label "Flexible" in {any}`tab-kspdefaults`.
622
623One can use `KSPMonitorDynamicTolerance()` to control the tolerances used by inner `KSP` solvers in `PCKSP`, `PCBJACOBI`, and `PCDEFLATION`.
624
625In addition to supporting `PCKSP`, the flexible methods support `KSP*SetModifyPC()`, for example, `KSPFGMRESSetModifyPC()`, these functions
626allow the user to provide a callback function that changes the preconditioner at each Krylov iteration. Its calling sequence is as follows.
627
628```
629PetscErrorCode f(KSP ksp,PetscInt total_its,PetscInt its_since_restart,PetscReal res_norm,void *ctx);
630```
631
632(sec_pipelineksp)=
633
634### Pipelined Krylov Methods
635
636Standard Krylov methods have one or more global reductions resulting from the computations of inner products or norms in each iteration.
637These reductions need to block until all MPI processes have received the results. For a large number of MPI processes (this number is machine dependent
638but can be above 10,000 processes) this synchronization is very time consuming and can significantly slow the computation. Pipelined Krylov
639methods overlap the reduction operations with local computations (generally the application of the matrix-vector products and precondtiioners)
640thus effectively "hiding" the time of the reductions. In addition, they may reduce the number of global synchronizations by rearranging the
641computations in a way that some of them can be collapsed, e.g., two or more calls to `MPI_Allreduce()` may be combined into one call.
642The pipeline `KSP` methods have the label "Pipeline" in {any}`tab-kspdefaults`.
643
644Special configuration of MPI may be necessary for reductions to make asynchronous progress, which is important for
645performance of pipelined methods. See {any}`doc_faq_pipelined` for details.
646
647### Other KSP Options
648
649To obtain the solution vector and right-hand side from a `KSP`
650context, one uses
651
652```
653KSPGetSolution(KSP ksp,Vec *x);
654KSPGetRhs(KSP ksp,Vec *rhs);
655```
656
657During the iterative process the solution may not yet have been
658calculated or it may be stored in a different location. To access the
659approximate solution during the iterative process, one uses the command
660
661```
662KSPBuildSolution(KSP ksp,Vec w,Vec *v);
663```
664
665where the solution is returned in `v`. The user can optionally provide
666a vector in `w` as the location to store the vector; however, if `w`
667is `NULL`, space allocated by PETSc in the `KSP` context is used.
668One should not destroy this vector. For certain `KSP` methods (e.g.,
669GMRES), the construction of the solution is expensive, while for many
670others it doesn’t even require a vector copy.
671
672Access to the residual is done in a similar way with the command
673
674```
675KSPBuildResidual(KSP ksp,Vec t,Vec w,Vec *v);
676```
677
678Again, for GMRES and certain other methods this is an expensive
679operation.
680
681(sec_pc)=
682
683## Preconditioners
684
685As discussed in {any}`sec_ksppc`, Krylov subspace methods
686are typically used in conjunction with a preconditioner. To employ a
687particular preconditioning method, the user can either select it from
688the options database using input of the form `-pc_type <methodname>`
689or set the method with the command
690
691```
692PCSetType(PC pc,PCType method);
693```
694
695In {any}`tab-pcdefaults` we summarize the basic
696preconditioning methods supported in PETSc. See the `PCType` manual
697page for a complete list.
698
699The `PCSHELL` preconditioner allows users to provide their own
700specific, application-provided custom preconditioner.
701
702The direct
703preconditioner, `PCLU` , is, in fact, a direct solver for the linear
704system that uses LU factorization. `PCLU` is included as a
705preconditioner so that PETSc has a consistent interface among direct and
706iterative linear solvers.
707
708PETSc provides several domain decomposition methods/preconditioners including
709`PCASM`, `PCGASM`, `PCBDDC`, and `PCHPDDM`. In addition PETSc provides
710multiple multigrid solvers/preconditioners including `PCMG`, `PCGAMG`, `PCHYPRE`,
711and `PCML`. See further discussion below.
712
713```{eval-rst}
714.. list-table:: PETSc Preconditioners (partial list)
715   :name: tab-pcdefaults
716   :header-rows: 1
717
718   * - Method
719     - PCType
720     - Options Database
721   * - Jacobi
722     - ``PCJACOBI``
723     - ``jacobi``
724   * - Block Jacobi
725     - ``PCBJACOBI``
726     - ``bjacobi``
727   * - SOR (and SSOR)
728     - ``PCSOR``
729     - ``sor``
730   * - SOR with Eisenstat trick
731     - ``PCEISENSTAT``
732     - ``eisenstat``
733   * - Incomplete Cholesky
734     - ``PCICC``
735     - ``icc``
736   * - Incomplete LU
737     - ``PCILU``
738     - ``ilu``
739   * - Additive Schwarz
740     - ``PCASM``
741     - ``asm``
742   * - Generalized Additive Schwarz
743     - ``PCGASM``
744     - ``gasm``
745   * - Algebraic Multigrid
746     - ``PCGAMG``
747     - ``gamg``
748   * - Balancing Domain Decomposition by Constraints
749     - ``PCBDDC``
750     - ``bddc``
751   * - Linear solver
752     - ``PCKSP``
753     - ``ksp``
754   * - Combination of preconditioners
755     - ``PCCOMPOSITE``
756     - ``composite``
757   * - LU
758     - ``PCLU``
759     - ``lu``
760   * - Cholesky
761     - ``PCCHOLESKY``
762     - ``cholesky``
763   * - No preconditioning
764     - ``PCNONE``
765     - ``none``
766   * - Shell for user-defined ``PC``
767     - ``PCSHELL``
768     - ``shell``
769```
770
771Each preconditioner may have associated with it a set of options, which
772can be set with routines and options database commands provided for this
773purpose. Such routine names and commands are all of the form
774`PC<TYPE><Option>` and `-pc_<type>_<option> [value]`. A complete
775list can be found by consulting the `PCType` manual page; we discuss
776just a few in the sections below.
777
778(sec_ilu_icc)=
779
780### ILU and ICC Preconditioners
781
782Some of the options for ILU preconditioner are
783
784```
785PCFactorSetLevels(PC pc,PetscInt levels);
786PCFactorSetReuseOrdering(PC pc,PetscBool flag);
787PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount);
788PCFactorSetReuseFill(PC pc,PetscBool flag);
789PCFactorSetUseInPlace(PC pc,PetscBool flg);
790PCFactorSetAllowDiagonalFill(PC pc,PetscBool flg);
791```
792
793When repeatedly solving linear systems with the same `KSP` context,
794one can reuse some information computed during the first linear solve.
795In particular, `PCFactorSetReuseOrdering()` causes the ordering (for
796example, set with `-pc_factor_mat_ordering_type` `order`) computed
797in the first factorization to be reused for later factorizations.
798`PCFactorSetUseInPlace()` is often used with `PCASM` or
799`PCBJACOBI` when zero fill is used, since it reuses the matrix space
800to store the incomplete factorization it saves memory and copying time.
801Note that in-place factorization is not appropriate with any ordering
802besides natural and cannot be used with the drop tolerance
803factorization. These options may be set in the database with
804
805- `-pc_factor_levels <levels>`
806- `-pc_factor_reuse_ordering`
807- `-pc_factor_reuse_fill`
808- `-pc_factor_in_place`
809- `-pc_factor_nonzeros_along_diagonal`
810- `-pc_factor_diagonal_fill`
811
812See {any}`sec_symbolfactor` for information on
813preallocation of memory for anticipated fill during factorization. By
814alleviating the considerable overhead for dynamic memory allocation,
815such tuning can significantly enhance performance.
816
817PETSc supports incomplete factorization preconditioners
818for several matrix types for sequential matrices (for example
819`MATSEQAIJ`, `MATSEQBAIJ`, and `MATSEQSBAIJ`).
820
821### SOR and SSOR Preconditioners
822
823PETSc provides only a sequential SOR preconditioner; it can only be
824used with sequential matrices or as the subblock preconditioner when
825using block Jacobi or ASM preconditioning (see below).
826
827The options for SOR preconditioning with `PCSOR` are
828
829```
830PCSORSetOmega(PC pc,PetscReal omega);
831PCSORSetIterations(PC pc,PetscInt its,PetscInt lits);
832PCSORSetSymmetric(PC pc,MatSORType type);
833```
834
835The first of these commands sets the relaxation factor for successive
836over (under) relaxation. The second command sets the number of inner
837iterations `its` and local iterations `lits` (the number of
838smoothing sweeps on a process before doing a ghost point update from the
839other processes) to use between steps of the Krylov space method. The
840total number of SOR sweeps is given by `its*lits`. The third command
841sets the kind of SOR sweep, where the argument `type` can be one of
842`SOR_FORWARD_SWEEP`, `SOR_BACKWARD_SWEEP` or
843`SOR_SYMMETRIC_SWEEP`, the default being `SOR_FORWARD_SWEEP`.
844Setting the type to be `SOR_SYMMETRIC_SWEEP` produces the SSOR method.
845In addition, each process can locally and independently perform the
846specified variant of SOR with the types `SOR_LOCAL_FORWARD_SWEEP`,
847`SOR_LOCAL_BACKWARD_SWEEP`, and `SOR_LOCAL_SYMMETRIC_SWEEP`. These
848variants can also be set with the options `-pc_sor_omega <omega>`,
849`-pc_sor_its <its>`, `-pc_sor_lits <lits>`, `-pc_sor_backward`,
850`-pc_sor_symmetric`, `-pc_sor_local_forward`,
851`-pc_sor_local_backward`, and `-pc_sor_local_symmetric`.
852
853The Eisenstat trick {cite}`eisenstat81` for SSOR
854preconditioning can be employed with the method `PCEISENSTAT`
855(`-pc_type` `eisenstat`). By using both left and right
856preconditioning of the linear system, this variant of SSOR requires
857about half of the floating-point operations for conventional SSOR. The
858option `-pc_eisenstat_no_diagonal_scaling` (or the routine
859`PCEisenstatSetNoDiagonalScaling()`) turns off diagonal scaling in
860conjunction with Eisenstat SSOR method, while the option
861`-pc_eisenstat_omega <omega>` (or the routine
862`PCEisenstatSetOmega(PC pc,PetscReal omega)`) sets the SSOR relaxation
863coefficient, `omega`, as discussed above.
864
865(sec_factorization)=
866
867### LU Factorization
868
869The LU preconditioner provides several options. The first, given by the
870command
871
872```
873PCFactorSetUseInPlace(PC pc,PetscBool flg);
874```
875
876causes the factorization to be performed in-place and hence destroys the
877original matrix. The options database variant of this command is
878`-pc_factor_in_place`. Another direct preconditioner option is
879selecting the ordering of equations with the command
880`-pc_factor_mat_ordering_type <ordering>`. The possible orderings are
881
882- `MATORDERINGNATURAL` - Natural
883- `MATORDERINGND` - Nested Dissection
884- `MATORDERING1WD` - One-way Dissection
885- `MATORDERINGRCM` - Reverse Cuthill-McKee
886- `MATORDERINGQMD` - Quotient Minimum Degree
887
888These orderings can also be set through the options database by
889specifying one of the following: `-pc_factor_mat_ordering_type`
890`natural`, or `nd`, or `1wd`, or `rcm`, or `qmd`. In addition,
891see `MatGetOrdering()`, discussed in {any}`sec_matfactor`.
892
893The sparse LU factorization provided in PETSc does not perform pivoting
894for numerical stability (since they are designed to preserve nonzero
895structure), and thus occasionally an LU factorization will fail with a
896zero pivot when, in fact, the matrix is non-singular. The option
897`-pc_factor_nonzeros_along_diagonal <tol>` will often help eliminate
898the zero pivot, by preprocessing the column ordering to remove small
899values from the diagonal. Here, `tol` is an optional tolerance to
900decide if a value is nonzero; by default it is `1.e-10`.
901
902In addition, {any}`sec_symbolfactor` provides information
903on preallocation of memory for anticipated fill during factorization.
904Such tuning can significantly enhance performance, since it eliminates
905the considerable overhead for dynamic memory allocation.
906
907(sec_bjacobi)=
908
909### Block Jacobi and Overlapping Additive Schwarz Preconditioners
910
911The block Jacobi and overlapping additive Schwarz (domain decomposition) methods in PETSc are
912supported in parallel; however, only the uniprocess version of the block
913Gauss-Seidel method is available. By default, the PETSc
914implementations of these methods employ ILU(0) factorization on each
915individual block (that is, the default solver on each subblock is
916`PCType=PCILU`, `KSPType=KSPPREONLY` (or equivalently `KSPType=KSPNONE`); the user can set alternative
917linear solvers via the options `-sub_ksp_type` and `-sub_pc_type`.
918In fact, all of the `KSP` and `PC` options can be applied to the
919subproblems by inserting the prefix `-sub_` at the beginning of the
920option name. These options database commands set the particular options
921for *all* of the blocks within the global problem. In addition, the
922routines
923
924```
925PCBJacobiGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **subksp);
926PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **subksp);
927```
928
929extract the `KSP` context for each local block. The argument
930`n_local` is the number of blocks on the calling process, and
931`first_local` indicates the global number of the first block on the
932process. The blocks are numbered successively by processes from zero
933through $b_g-1$, where $b_g$ is the number of global blocks.
934The array of `KSP` contexts for the local blocks is given by
935`subksp`. This mechanism enables the user to set different solvers for
936the various blocks. To set the appropriate data structures, the user
937*must* explicitly call `KSPSetUp()` before calling
938`PCBJacobiGetSubKSP()` or `PCASMGetSubKSP(`). For further details,
939see
940<a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex7.c.html">KSP Tutorial ex7</a>
941or
942<a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex8.c.html">KSP Tutorial ex8</a>.
943
944The block Jacobi, block Gauss-Seidel, and additive Schwarz
945preconditioners allow the user to set the number of blocks into which
946the problem is divided. The options database commands to set this value
947are `-pc_bjacobi_blocks` `n` and `-pc_bgs_blocks` `n`, and,
948within a program, the corresponding routines are
949
950```
951PCBJacobiSetTotalBlocks(PC pc,PetscInt blocks,PetscInt *size);
952PCASMSetTotalSubdomains(PC pc,PetscInt n,IS *is,IS *islocal);
953PCASMSetType(PC pc,PCASMType type);
954```
955
956The optional argument `size` is an array indicating the size of each
957block. Currently, for certain parallel matrix formats, only a single
958block per process is supported. However, the `MATMPIAIJ` and
959`MATMPIBAIJ` formats support the use of general blocks as long as no
960blocks are shared among processes. The `is` argument contains the
961index sets that define the subdomains.
962
963The object `PCASMType` is one of `PC_ASM_BASIC`,
964`PC_ASM_INTERPOLATE`, `PC_ASM_RESTRICT`, or `PC_ASM_NONE` and may
965also be set with the options database `-pc_asm_type` `[basic`,
966`interpolate`, `restrict`, `none]`. The type `PC_ASM_BASIC` (or
967`-pc_asm_type` `basic`) corresponds to the standard additive Schwarz
968method that uses the full restriction and interpolation operators. The
969type `PC_ASM_RESTRICT` (or `-pc_asm_type` `restrict`) uses a full
970restriction operator, but during the interpolation process ignores the
971off-process values. Similarly, `PC_ASM_INTERPOLATE` (or
972`-pc_asm_type` `interpolate`) uses a limited restriction process in
973conjunction with a full interpolation, while `PC_ASM_NONE` (or
974`-pc_asm_type` `none`) ignores off-process values for both
975restriction and interpolation. The ASM types with limited restriction or
976interpolation were suggested by Xiao-Chuan Cai and Marcus Sarkis
977{cite}`cs99`. `PC_ASM_RESTRICT` is the PETSc default, as
978it saves substantial communication and for many problems has the added
979benefit of requiring fewer iterations for convergence than the standard
980additive Schwarz method.
981
982The user can also set the number of blocks and sizes on a per-process
983basis with the commands
984
985```
986PCBJacobiSetLocalBlocks(PC pc,PetscInt blocks,PetscInt *size);
987PCASMSetLocalSubdomains(PC pc,PetscInt N,IS *is,IS *islocal);
988```
989
990For the ASM preconditioner one can use the following command to set the
991overlap to compute in constructing the subdomains.
992
993```
994PCASMSetOverlap(PC pc,PetscInt overlap);
995```
996
997The overlap defaults to 1, so if one desires that no additional overlap
998be computed beyond what may have been set with a call to
999`PCASMSetTotalSubdomains()` or `PCASMSetLocalSubdomains()`, then
1000`overlap` must be set to be 0. In particular, if one does *not*
1001explicitly set the subdomains in an application code, then all overlap
1002would be computed internally by PETSc, and using an overlap of 0 would
1003result in an ASM variant that is equivalent to the block Jacobi
1004preconditioner. Note that one can define initial index sets `is` with
1005*any* overlap via `PCASMSetTotalSubdomains()` or
1006`PCASMSetLocalSubdomains()`; the routine `PCASMSetOverlap()` merely
1007allows PETSc to extend that overlap further if desired.
1008
1009`PCGASM` is a generalization of `PCASM` that allows
1010the user to specify subdomains that span multiple MPI processes. This can be
1011useful for problems where small subdomains result in poor convergence.
1012To be effective, the multi-processor subproblems must be solved using a
1013sufficiently strong subsolver, such as `PCLU`, for which `SuperLU_DIST` or a
1014similar parallel direct solver could be used; other choices may include
1015a multigrid solver on the subdomains.
1016
1017The interface for `PCGASM` is similar to that of `PCASM`. In
1018particular, `PCGASMType` is one of `PC_GASM_BASIC`,
1019`PC_GASM_INTERPOLATE`, `PC_GASM_RESTRICT`, `PC_GASM_NONE`. These
1020options have the same meaning as with `PCASM` and may also be set with
1021the options database `-pc_gasm_type` `[basic`, `interpolate`,
1022`restrict`, `none]`.
1023
1024Unlike `PCASM`, however, `PCGASM` allows the user to define
1025subdomains that span multiple MPI processes. The simplest way to do this is
1026using a call to `PCGASMSetTotalSubdomains(PC pc,PetscInt N)` with
1027the total number of subdomains `N` that is smaller than the MPI
1028communicator `size`. In this case `PCGASM` will coalesce `size/N`
1029consecutive single-rank subdomains into a single multi-rank subdomain.
1030The single-rank subdomains contain the degrees of freedom corresponding
1031to the locally-owned rows of the `PCGASM` matrix used to compute the preconditioner –
1032these are the subdomains `PCASM` and `PCGASM` use by default.
1033
1034Each of the multirank subdomain subproblems is defined on the
1035subcommunicator that contains the coalesced `PCGASM` processes. In general
1036this might not result in a very good subproblem if the single-rank
1037problems corresponding to the coalesced processes are not very strongly
1038connected. In the future this will be addressed with a hierarchical
1039partitioner that generates well-connected coarse subdomains first before
1040subpartitioning them into the single-rank subdomains.
1041
1042In the meantime the user can provide his or her own multi-rank
1043subdomains by calling `PCGASMSetSubdomains(PC,IS[],IS[])` where each
1044of the `IS` objects on the list defines the inner (without the
1045overlap) or the outer (including the overlap) subdomain on the
1046subcommunicator of the `IS` object. A helper subroutine
1047`PCGASMCreateSubdomains2D()` is similar to PCASM’s but is capable of
1048constructing multi-rank subdomains that can be then used with
1049`PCGASMSetSubdomains()`. An alternative way of creating multi-rank
1050subdomains is by using the underlying `DM` object, if it is capable of
1051generating such decompositions via `DMCreateDomainDecomposition()`.
1052Ordinarily the decomposition specified by the user via
1053`PCGASMSetSubdomains()` takes precedence, unless
1054`PCGASMSetUseDMSubdomains()` instructs `PCGASM` to prefer
1055`DM`-created decompositions.
1056
1057Currently there is no support for increasing the overlap of multi-rank
1058subdomains via `PCGASMSetOverlap()` – this functionality works only
1059for subdomains that fit within a single MPI process, exactly as in
1060`PCASM`.
1061
1062Examples of the described `PCGASM` usage can be found in
1063<a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex62.c.html">KSP Tutorial ex62</a>.
1064In particular, `runex62_superlu_dist` illustrates the use of
1065`SuperLU_DIST` as the subdomain solver on coalesced multi-rank
1066subdomains. The `runex62_2D_*` examples illustrate the use of
1067`PCGASMCreateSubdomains2D()`.
1068
1069(sec_amg)=
1070
1071### Algebraic Multigrid (AMG) Preconditioners
1072
1073PETSc has a native algebraic multigrid preconditioner `PCGAMG` –
1074*gamg* – and interfaces to three external AMG packages: *hypre*, *ML*
1075and *AMGx* (CUDA platforms only) that can be downloaded in the
1076configuration phase (e.g., `--download-hypre` ) and used by
1077specifying that command line parameter (e.g., `-pc_type hypre`).
1078*Hypre* is relatively monolithic in that a PETSc matrix is converted into a hypre
1079matrix, and then *hypre* is called to solve the entire problem. *ML* is more
1080modular because PETSc only has *ML* generate the coarse grid spaces
1081(columns of the prolongation operator), which is the core of an AMG method,
1082and then constructs a `PCMG` with Galerkin coarse grid operator
1083construction. `PCGAMG` is designed from the beginning to be modular, to
1084allow for new components to be added easily and also populates a
1085multigrid preconditioner `PCMG` so generic multigrid parameters are
1086used (see {any}`sec_mg`). PETSc provides a fully supported (smoothed) aggregation AMG, but supports the addition of new methods
1087(`-pc_type gamg -pc_gamg_type agg` or `PCSetType(pc,PCGAMG)` and
1088`PCGAMGSetType(pc, PCGAMGAGG)`. Examples of extension are reference implementations of
1089a classical AMG method (`-pc_gamg_type classical`), a (2D) hybrid geometric
1090AMG method (`-pc_gamg_type geo`) that are not supported. A 2.5D AMG method DofColumns
1091{cite}`isaacstadlerghattas2015` supports 2D coarsenings extruded in the third dimension. `PCGAMG` does require the use
1092of `MATAIJ` matrices. For instance, `MATBAIJ` matrices are not supported. One
1093can use `MATAIJ` instead of `MATBAIJ` without changing any code other than the
1094constructor (or the `-mat_type` from the command line). For instance,
1095`MatSetValuesBlocked` works with `MATAIJ` matrices.
1096
1097**Important parameters for PCGAMGAGG**
1098
1099- Control the generation of the coarse grid
1100
1101  > - `-pc_gamg_aggressive_coarsening` \<n:int:1> Use aggressive coarsening on the finest n levels to construct the coarser mesh.
1102  >   See `PCGAMGAGGSetNSmooths()`. The larger value produces a faster preconditioner to create and solve, but the convergence may be slower.
1103  > - `-pc_gamg_low_memory_threshold_filter` \<bool:false> Filter small matrix entries before coarsening the mesh.
1104  >   See `PCGAMGSetLowMemoryFilter()`.
1105  > - `-pc_gamg_threshold` \<tol:real:0.0> The threshold of small values to drop when `-pc_gamg_low_memory_threshold_filter` is used. A
1106  >   negative value means keeping even the locations with 0.0. See `PCGAMGSetThreshold()`
1107  > - `-pc_gamg_threshold_scale` \<v>:real:1.0> Set a scale factor applied to each coarser level when `-pc_gamg_low_memory_threshold_filter` is used.
1108  >   See `PCGAMGSetThresholdScale()`.
1109  > - `-pc_gamg_mat_coarsen_type` \<mis|hem|misk:misk> Algorithm used to coarsen the matrix graph. See `MatCoarsenSetType()`.
1110  > - `-pc_gamg_mat_coarsen_max_it` \<it:int:4> Maximum HEM iterations to use. See `MatCoarsenSetMaximumIterations()`.
1111  > - `-pc_gamg_aggressive_mis_k` \<k:int:2> k distance in MIS coarsening (>2 is 'aggressive') to use in coarsening.
1112  >   See `PCGAMGMISkSetAggressive()`. The larger value produces a preconditioner that is faster to create and solve with but the convergence may be slower.
1113  >   This option and the previous option work to determine how aggressively the grids are coarsened.
1114  > - `-pc_gamg_mis_k_minimum_degree_ordering` \<bool:false> Use a minimum degree ordering in the greedy MIS algorithm used to coarsen.
1115  >   See `PCGAMGMISkSetMinDegreeOrdering()`
1116
1117- Control the generation of the prolongation for `PCGAMGAGG`
1118
1119  > - `-pc_gamg_agg_nsmooths` \<n:int:1> Number of smoothing steps to be used in constructing the prolongation. For symmetric problems,
1120  >   generally, one or more is best. For some strongly nonsymmetric problems, 0 may be best. See `PCGAMGSetNSmooths()`.
1121
1122- Control the amount of parallelism on the levels
1123
1124  > - `-pc_gamg_process_eq_limit` \<n:int:50> Sets the minimum number of equations allowed per process when coarsening (otherwise, fewer MPI processes
1125  >   are used for the coarser mesh). A larger value will cause the coarser problems to be run on fewer MPI processes, resulting
1126  >   in less communication and possibly a faster time to solution. See `PCGAMGSetProcEqLim()`.
1127  >
1128  > - `-pc_gamg_rank_reduction_factors` \<rn,rn-1,...,r1:int> Set a schedule for MPI rank reduction on coarse grids. `See PCGAMGSetRankReductionFactors()`
1129  >   This overrides the lessening of processes that would arise from `-pc_gamg_process_eq_limit`.
1130  >
1131  > - `-pc_gamg_repartition` \<bool:false> Run a partitioner on each coarser mesh generated rather than using the default partition arising from the
1132  >   finer mesh. See `PCGAMGSetRepartition()`. This increases the preconditioner setup time but will result in less time per
1133  >   iteration of the solver.
1134  >
1135  > - `-pc_gamg_parallel_coarse_grid_solver` \<bool:false> Allow the coarse grid solve to run in parallel, depending on the value of `-pc_gamg_coarse_eq_limit`.
1136  >   See `PCGAMGSetParallelCoarseGridSolve()`. If the coarse grid problem is large then this can
1137  >   improve the time to solution.
1138  >
1139  >   - `-pc_gamg_coarse_eq_limit` \<n:int:50> Sets the minimum number of equations allowed per process on the coarsest level when coarsening
1140  >     (otherwise fewer MPI processes will be used). A larger value will cause the coarse problems to be run on fewer MPI processes.
1141  >     This only applies if `-pc_gamg_parallel_coarse_grid_solver` is set to true. See `PCGAMGSetCoarseEqLim()`.
1142
1143- Control the smoothers
1144
1145  > - `-pc_mg_levels` \<n:int> Set the maximum number of levels to use.
1146  > - `-mg_levels_ksp_type` \<KSPType:chebyshev> If `KSPCHEBYSHEV` or `KSPRICHARDSON` is not used, then the Krylov
1147  >   method for the entire multigrid solve has to be a flexible method such as `KSPFGMRES`. Generally, the
1148  >   stronger the Krylov method the faster the convergence, but with more cost per iteration. See `KSPSetType()`.
1149  > - `-mg_levels_ksp_max_it` \<its:int:2> Sets the number of iterations to run the smoother on each level. Generally, the more iterations
1150  >   , the faster the convergence, but with more cost per multigrid iteration. See `PCMGSetNumberSmooth()`.
1151  > - `-mg_levels_ksp_xxx` Sets options for the `KSP` in the smoother on the levels.
1152  > - `-mg_levels_pc_type` \<PCType:jacobi> Sets the smoother to use on each level. See `PCSetType()`. Generally, the
1153  >   stronger the preconditioner the faster the convergence, but with more cost per iteration.
1154  > - `-mg_levels_pc_xxx` Sets options for the `PC` in the smoother on the levels.
1155  > - `-mg_coarse_ksp_type` \<KSPType:none> Sets the solver `KSPType` to use on the coarsest level.
1156  > - `-mg_coarse_pc_type` \<PCType:lu> Sets the solver `PCType` to use on the coarsest level.
1157  > - `-pc_gamg_asm_use_agg` \<bool:false> Use `PCASM` as the smoother on each level with the aggregates defined by the coarsening process are
1158  >   the subdomains. This option automatically switches the smoother on the levels to be `PCASM`.
1159  > - `-mg_levels_pc_asm_overlap` \<n:int:0> Use non-zero overlap with `-pc_gamg_asm_use_agg`. See `PCASMSetOverlap()`.
1160
1161- Control the multigrid algorithm
1162
1163  > - `-pc_mg_type` \<additive|multiplicative|full|kaskade:multiplicative> The type of multigrid to use. Usually, multiplicative is the fastest.
1164  > - `-pc_mg_cycle_type` \<v|w:v> Use V- or W-cycle with `-pc_mg_type` `multiplicative`
1165
1166`PCGAMG` provides unsmoothed aggregation (`-pc_gamg_agg_nsmooths 0`) and
1167smoothed aggregation (`-pc_gamg_agg_nsmooths 1` or
1168`PCGAMGSetNSmooths(pc,1)`). Smoothed aggregation (SA), {cite}`vanek1996algebraic`, {cite}`vanek2001convergence`, is recommended
1169for symmetric positive definite systems. Unsmoothed aggregation can be
1170useful for asymmetric problems and problems where the highest eigenestimates are problematic. If poor convergence rates are observed using
1171the smoothed version, one can test unsmoothed aggregation.
1172
1173**Eigenvalue estimates:** The parameters for the KSP eigen estimator,
1174used for SA, can be set with `-pc_gamg_esteig_ksp_max_it` and
1175`-pc_gamg_esteig_ksp_type`. For example, CG generally converges to the
1176highest eigenvalue faster than GMRES (the default for KSP) if your problem
1177is symmetric positive definite. One can specify CG with
1178`-pc_gamg_esteig_ksp_type cg`. The default for
1179`-pc_gamg_esteig_ksp_max_it` is 10, which we have found is pretty safe
1180with a (default) safety factor of 1.1. One can specify the range of real
1181eigenvalues in the same way as with Chebyshev KSP solvers
1182(smoothers), with `-pc_gamg_eigenvalues <emin,emax>`. GAMG sets the MG
1183smoother type to chebyshev by default. By default, GAMG uses its eigen
1184estimate, if it has one, for Chebyshev smoothers if the smoother uses
1185Jacobi preconditioning. This can be overridden with
1186`-pc_gamg_use_sa_esteig  <true,false>`.
1187
1188AMG methods require knowledge of the number of degrees of freedom per
1189vertex; the default is one (a scalar problem). Vector problems like
1190elasticity should set the block size of the matrix appropriately with
1191`-mat_block_size bs` or `MatSetBlockSize(mat,bs)`. Equations must be
1192ordered in “vertex-major” ordering (e.g.,
1193$x_1,y_1,z_1,x_2,y_2,...$).
1194
1195**Near null space:** Smoothed aggregation requires an explicit
1196representation of the (near) null space of the operator for optimal
1197performance. One can provide an orthonormal set of null space vectors
1198with `MatSetNearNullSpace()`. The vector of all ones is the default
1199for each variable given by the block size (e.g., the translational rigid
1200body modes). For elasticity, where rotational rigid body modes are
1201required to complete the near null-space you can use
1202`MatNullSpaceCreateRigidBody()` to create the null space vectors and
1203then `MatSetNearNullSpace()`.
1204
1205**Coarse grid data model:** The GAMG framework provides for reducing the
1206number of active processes on coarse grids to reduce communication costs
1207when there is not enough parallelism to keep relative communication
1208costs down. Most AMG solvers reduce to just one active process on the
1209coarsest grid (the PETSc MG framework also supports redundantly solving
1210the coarse grid on all processes to reduce communication
1211costs potentially). However, this forcing to one process can be overridden if one
1212wishes to use a parallel coarse grid solver. GAMG generalizes this by
1213reducing the active number of processes on other coarse grids.
1214GAMG will select the number of active processors by fitting the desired
1215number of equations per process (set with
1216`-pc_gamg_process_eq_limit <50>,`) at each level given that size of
1217each level. If $P_i < P$ processors are desired on a level
1218$i$, then the first $P_i$ processes are populated with the grid
1219and the remaining are empty on that grid. One can, and probably should,
1220repartition the coarse grids with `-pc_gamg_repartition <true>`,
1221otherwise an integer process reduction factor ($q$) is selected
1222and the equations on the first $q$ processes are moved to process
12230, and so on. As mentioned, multigrid generally coarsens the problem
1224until it is small enough to be solved with an exact solver (e.g., LU or
1225SVD) in a relatively short time. GAMG will stop coarsening when the
1226number of the equation on a grid falls below the threshold given by
1227`-pc_gamg_coarse_eq_limit <50>,`.
1228
1229**Coarse grid parameters:** There are several options to provide
1230parameters to the coarsening algorithm and parallel data layout. Run a
1231code using `PCGAMG` with `-help` to get a full listing of GAMG
1232parameters with short descriptions. The rate of coarsening is
1233critical in AMG performance – too slow coarsening will result in an
1234overly expensive solver per iteration and too fast coarsening will
1235result in decrease in the convergence rate. `-pc_gamg_threshold <-1>`
1236and `-pc_gamg_aggressive_coarsening <N>` are the primary parameters that
1237control coarsening rates, which is very important for AMG performance. A
1238greedy maximal independent set (MIS) algorithm is used in coarsening.
1239Squaring the graph implements MIS-2; the root vertex in an
1240aggregate is more than two edges away from another root vertex instead
1241of more than one in MIS. The threshold parameter sets a normalized
1242threshold for which edges are removed from the MIS graph, thereby
1243coarsening slower. Zero will keep all non-zero edges, a negative number
1244will keep zero edges, and a positive number will drop small edges. Typical
1245finite threshold values are in the range of $0.01 - 0.05$. There
1246are additional parameters for changing the weights on coarse grids.
1247
1248The parallel MIS algorithms require symmetric weights/matrices. Thus `PCGAMG`
1249will automatically make the graph symmetric if it is not symmetric. Since this
1250has additional cost, users should indicate the symmetry of the matrices they
1251provide by calling
1252
1253```
1254MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE (or PETSC_FALSE))
1255```
1256
1257or
1258
1259```
1260MatSetOption(mat,MAT_STRUCTURALLY_SYMMETRIC,PETSC_TRUE (or PETSC_FALSE)).
1261```
1262
1263If they know that the matrix will always have symmetry despite future changes
1264to the matrix (with, for example, `MatSetValues()`) then they should also call
1265
1266```
1267MatSetOption(mat,MAT_SYMMETRY_ETERNAL,PETSC_TRUE (or PETSC_FALSE))
1268```
1269
1270or
1271
1272```
1273MatSetOption(mat,MAT_STRUCTURAL_SYMMETRY_ETERNAL,PETSC_TRUE (or PETSC_FALSE)).
1274```
1275
1276Using this information allows the algorithm to skip unnecessary computations.
1277
1278**Troubleshooting algebraic multigrid methods:** If `PCGAMG`, *ML*, *AMGx* or
1279*hypre* does not perform well; the first thing to try is one of the other
1280methods. Often, the default parameters or just the strengths of different
1281algorithms can fix performance problems or provide useful information to
1282guide further debugging. There are several sources of poor performance
1283of AMG solvers and often special purpose methods must be developed to
1284achieve the full potential of multigrid. To name just a few sources of
1285performance degradation that may not be fixed with parameters in PETSc
1286currently: non-elliptic operators, curl/curl operators, highly stretched
1287grids or highly anisotropic problems, large jumps in material
1288coefficients with complex geometry (AMG is particularly well suited to
1289jumps in coefficients, but it is not a perfect solution), highly
1290incompressible elasticity, not to mention ill-posed problems and many
1291others. For Grad-Div and Curl-Curl operators, you may want to try the
1292Auxiliary-space Maxwell Solver (AMS,
1293`-pc_type hypre -pc_hypre_type ams`) or the Auxiliary-space Divergence
1294Solver (ADS, `-pc_type hypre -pc_hypre_type ads`) solvers. These
1295solvers need some additional information on the underlying mesh;
1296specifically, AMS needs the discrete gradient operator, which can be
1297specified via `PCHYPRESetDiscreteGradient()`. In addition to the
1298discrete gradient, ADS also needs the specification of the discrete curl
1299operator, which can be set using `PCHYPRESetDiscreteCurl()`.
1300
1301**I am converging slowly, what do I do?** AMG methods are sensitive to
1302coarsening rates and methods; for GAMG use `-pc_gamg_threshold <x>`
1303or `PCGAMGSetThreshold()` to regulate coarsening rates; higher values decrease
1304the coarsening rate. A high threshold (e.g., $x=0.08$) will result in an
1305expensive but potentially powerful preconditioner, and a low threshold
1306(e.g., $x=0.0$) will result in faster coarsening, fewer levels,
1307cheaper solves, and generally worse convergence rates.
1308
1309Aggressive_coarsening is the second mechanism for
1310increasing the coarsening rate and thereby decreasing the cost of the
1311coarse grids and generally decreasing the solver convergence rate.
1312Use `-pc_gamg_aggressive_coarsening <N>`, or
1313`PCGAMGSetAggressiveLevels(pc,N)`, to aggressively coarsen the graph on the finest N
1314levels. The default is $N=1$. There are two options for aggressive coarsening: 1) the
1315default, square graph: use $A^T A$ in the MIS coarsening algorithm and 2) coarsen with MIS-2 (instead
1316of the default of MIS-1). Use `-pc_gamg_aggressive_square_graph false`
1317to use MIS-k coarsening and `-pc_gamg_aggressive_mis_k k` to select
1318the level of MIS other than the default $k=2$.
1319The square graph approach seems to coarsen slower, which results in
1320larger coarse grids and is more expensive, but generally improves
1321the convergence rate.
1322If the coarse grids are expensive to compute, and use a lot of memory,
1323using MIS-2 is a good alternative (setting MIS-1 effectively turns
1324aggressive coarsening off). Note that MIS-3 is also supported.
1325
1326One can run with `-info :pc` and grep for `PCGAMG` to get statistics on
1327each level, which can be used to see if you are coarsening at an
1328appropriate rate. With smoothed aggregation, you generally want to coarse
1329at about a rate of 3:1 in each dimension. Coarsening too slowly will
1330result in large numbers of non-zeros per row on coarse grids (this is
1331reported). The number of non-zeros can go up very high, say about 300
1332(times the degrees of freedom per vertex) on a 3D hex mesh. One can also
1333look at the grid complexity, which is also reported (the ratio of the
1334total number of matrix entries for all levels to the number of matrix
1335entries on the fine level). Grid complexity should be well under 2.0 and
1336preferably around $1.3$ or lower. If convergence is poor and the
1337Galerkin coarse grid construction is much smaller than the time for each
1338solve, one can safely decrease the coarsening rate.
1339`-pc_gamg_threshold` $-1.0$ is the simplest and most robust
1340option and is recommended if poor convergence rates are observed, at
1341least until the source of the problem is discovered. In conclusion, decreasing the coarsening rate (increasing the
1342threshold) should be tried if convergence is slow.
1343
1344**A note on Chebyshev smoothers.** Chebyshev solvers are attractive as
1345multigrid smoothers because they can target a specific interval of the
1346spectrum, which is the purpose of a smoother. The spectral bounds for
1347Chebyshev solvers are simple to compute because they rely on the highest
1348eigenvalue of your (diagonally preconditioned) operator, which is
1349conceptually simple to compute. However, if this highest eigenvalue
1350estimate is not accurate (too low), the solvers can fail with an
1351indefinite preconditioner message. One can run with `-info` and grep
1352for `PCGAMG` to get these estimates or use `-ksp_view`. These highest
1353eigenvalues are generally between 1.5-3.0. For symmetric positive
1354definite systems, CG is a better eigenvalue estimator
1355`-mg_levels_esteig_ksp_type cg`. Bad Eigen estimates often cause indefinite matrix messages. Explicitly damped Jacobi or Krylov
1356smoothers can provide an alternative to Chebyshev, and *hypre* has
1357alternative smoothers.
1358
1359**Now, am I solving alright? Can I expect better?** If you find that you
1360are getting nearly one digit in reduction of the residual per iteration
1361and are using a modest number of point smoothing steps (e.g., 1-4
1362iterations of SOR), then you may be fairly close to textbook multigrid
1363efficiency. However, you also need to check the setup costs. This can be
1364determined by running with `-log_view` and check that the time for the
1365Galerkin coarse grid construction (`MatPtAP()`) is not (much) more than
1366the time spent in each solve (`KSPSolve()`). If the `MatPtAP()` time is
1367too large, then one can increase the coarsening rate by decreasing the
1368threshold and using aggressive coarsening
1369(`-pc_gamg_aggressive_coarsening <N>`, squares the graph on the finest N
1370levels). Likewise, if your `MatPtAP()` time is short and your convergence
1371If the rate is not ideal, you could decrease the coarsening rate.
1372
1373PETSc’s AMG solver is a framework for developers to
1374easily add AMG capabilities, like new AMG methods or an AMG component
1375like a matrix triple product. Contact us directly if you are interested
1376in contributing.
1377
1378Using algebraic multigrid as a "standalone" solver is possible but not recommended, as it does not accelerate it with a Krylov method.
1379Use a `KSPType` of `KSPRICHARDSON`
1380(or equivalently `-ksp_type richardson`) to achieve this. Using `KSPPREONLY` will not work since it only applies a single multigrid cycle.
1381
1382#### Adaptive Interpolation
1383
1384**Interpolation** transfers a function from the coarse space to the fine space. We would like this process to be accurate for the functions resolved by the coarse grid, in particular the approximate solution computed there. By default, we create these matrices using local interpolation of the fine grid dual basis functions in the coarse basis. However, an adaptive procedure can optimize the coefficients of the interpolator to reproduce pairs of coarse/fine functions which should approximate the lowest modes of the generalized eigenproblem
1385
1386$$
1387A x = \lambda M x
1388$$
1389
1390where $A$ is the system matrix and $M$ is the smoother. Note that for defect-correction MG, the interpolated solution from the coarse space need not be as accurate as the fine solution, for the same reason that updates in iterative refinement can be less accurate. However, in FAS or in the final interpolation step for each level of Full Multigrid, we must have interpolation as accurate as the fine solution since we are moving the entire solution itself.
1391
1392**Injection** should accurately transfer the fine solution to the coarse grid. Accuracy here means that the action of a coarse dual function on either should produce approximately the same result. In the structured grid case, this means that we just use the same values on coarse points. This can result in aliasing.
1393
1394**Restriction** is intended to transfer the fine residual to the coarse space. Here we use averaging (often the transpose of the interpolation operation) to damp out the fine space contributions. Thus, it is less accurate than injection, but avoids aliasing of the high modes.
1395
1396For a multigrid cycle, the interpolator $P$ is intended to accurately reproduce "smooth" functions from the coarse space in the fine space, keeping the energy of the interpolant about the same. For the Laplacian on a structured mesh, it is easy to determine what these low-frequency functions are. They are the Fourier modes. However an arbitrary operator $A$ will have different coarse modes that we want to resolve accurately on the fine grid, so that our coarse solve produces a good guess for the fine problem. How do we make sure that our interpolator $P$ can do this?
1397
1398We first must decide what we mean by accurate interpolation of some functions. Suppose we know the continuum function $f$ that we care about, and we are only interested in a finite element description of discrete functions. Then the coarse function representing $f$ is given by
1399
1400$$
1401f^C = \sum_i f^C_i \phi^C_i,
1402$$
1403
1404and similarly the fine grid form is
1405
1406$$
1407f^F = \sum_i f^F_i \phi^F_i.
1408$$
1409
1410Now we would like the interpolant of the coarse representer to the fine grid to be as close as possible to the fine representer in a least squares sense, meaning we want to solve the minimization problem
1411
1412$$
1413\min_{P} \| f^F - P f^C \|_2
1414$$
1415
1416Now we can express $P$ as a matrix by looking at the matrix elements $P_{ij} = \phi^F_i P \phi^C_j$. Then we have
1417
1418$$
1419\begin{aligned}
1420  &\phi^F_i f^F - \phi^F_i P f^C \\
1421= &f^F_i - \sum_j P_{ij} f^C_j
1422\end{aligned}
1423$$
1424
1425so that our discrete optimization problem is
1426
1427$$
1428\min_{P_{ij}} \| f^F_i - \sum_j P_{ij} f^C_j \|_2
1429$$
1430
1431and we will treat each row of the interpolator as a separate optimization problem. We could allow an arbitrary sparsity pattern, or try to determine adaptively, as is done in sparse approximate inverse preconditioning. However, we know the supports of the basis functions in finite elements, and thus the naive sparsity pattern from local interpolation can be used.
1432
1433We note here that the BAMG framework of Brannick et al. {cite}`brandtbrannickkahllivshits2011` does not use fine and coarse functions spaces, but rather a fine point/coarse point division which we will not employ here. Our general PETSc routine should work for both since the input would be the checking set (fine basis coefficients or fine space points) and the approximation set (coarse basis coefficients in the support or coarse points in the sparsity pattern).
1434
1435We can easily solve the above problem using QR factorization. However, there are many smooth functions from the coarse space that we want interpolated accurately, and a single $f$ would not constrain the values $P_{ij}`$ well. Therefore, we will use several functions $\{f_k\}$ in our minimization,
1436
1437$$
1438\begin{aligned}
1439  &\min_{P_{ij}} \sum_k w_k \| f^{F,k}_i - \sum_j P_{ij} f^{C,k}_j \|_2 \\
1440= &\min_{P_{ij}} \sum_k \| \sqrt{w_k} f^{F,k}_i - \sqrt{w_k} \sum_j P_{ij} f^{C,k}_j \|_2 \\
1441= &\min_{P_{ij}} \| W^{1/2} \mathbf{f}^{F}_i - W^{1/2} \mathbf{f}^{C} p_i \|_2
1442\end{aligned}
1443$$
1444
1445where
1446
1447$$
1448\begin{aligned}
1449W         &= \begin{pmatrix} w_0 & & \\ & \ddots & \\ & & w_K \end{pmatrix} \\
1450\mathbf{f}^{F}_i &= \begin{pmatrix} f^{F,0}_i \\ \vdots \\ f^{F,K}_i \end{pmatrix} \\
1451\mathbf{f}^{C}   &= \begin{pmatrix} f^{C,0}_0 & \cdots & f^{C,0}_n \\ \vdots & \ddots &  \vdots \\ f^{C,K}_0 & \cdots & f^{C,K}_n \end{pmatrix} \\
1452p_i       &= \begin{pmatrix} P_{i0} \\ \vdots \\ P_{in} \end{pmatrix}
1453\end{aligned}
1454$$
1455
1456or alternatively
1457
1458$$
1459\begin{aligned}
1460[W]_{kk}     &= w_k \\
1461[f^{F}_i]_k  &= f^{F,k}_i \\
1462[f^{C}]_{kj} &= f^{C,k}_j \\
1463[p_i]_j      &= P_{ij}
1464\end{aligned}
1465$$
1466
1467We thus have a standard least-squares problem
1468
1469$$
1470\min_{P_{ij}} \| b - A x \|_2
1471$$
1472
1473where
1474
1475$$
1476\begin{aligned}
1477A &= W^{1/2} f^{C} \\
1478b &= W^{1/2} f^{F}_i \\
1479x &= p_i
1480\end{aligned}
1481$$
1482
1483which can be solved using LAPACK.
1484
1485We will typically perform this optimization on a multigrid level $l$ when the change in eigenvalue from level $l+1$ is relatively large, meaning
1486
1487$$
1488\frac{|\lambda_l - \lambda_{l+1}|}{|\lambda_l|}.
1489$$
1490
1491This indicates that the generalized eigenvector associated with that eigenvalue was not adequately represented by $P^l_{l+1}`$, and the interpolator should be recomputed.
1492
1493```{raw} html
1494<hr>
1495```
1496
1497### Balancing Domain Decomposition by Constraints
1498
1499PETSc provides the Balancing Domain Decomposition by Constraints (`PCBDDC`)
1500method for preconditioning parallel finite element problems stored in
1501unassembled format (see `MATIS`). `PCBDDC` is a 2-level non-overlapping
1502domain decomposition method which can be easily adapted to different
1503problems and discretizations by means of few user customizations. The
1504application of the preconditioner to a vector consists in the static
1505condensation of the residual at the interior of the subdomains by means
1506of local Dirichlet solves, followed by an additive combination of Neumann
1507local corrections and the solution of a global coupled coarse problem.
1508Command line options for the underlying `KSP` objects are prefixed by
1509`-pc_bddc_dirichlet`, `-pc_bddc_neumann`, and `-pc_bddc_coarse`
1510respectively.
1511
1512The implementation supports any kind of linear system, and
1513assumes a one-to-one mapping between subdomains and MPI processes.
1514Complex numbers are supported as well. For non-symmetric problems, use
1515the runtime option `-pc_bddc_symmetric 0`.
1516
1517Unlike conventional non-overlapping methods that iterates just on the
1518degrees of freedom at the interface between subdomain, `PCBDDC`
1519iterates on the whole set of degrees of freedom, allowing the use of
1520approximate subdomain solvers. When using approximate solvers, the
1521command line switches `-pc_bddc_dirichlet_approximate` and/or
1522`-pc_bddc_neumann_approximate` should be used to inform `PCBDDC`. If
1523any of the local problems is singular, the nullspace of the local
1524operator should be attached to the local matrix via
1525`MatSetNullSpace()`.
1526
1527At the basis of the method there’s the analysis of the connected
1528components of the interface for the detection of vertices, edges and
1529faces equivalence classes. Additional information on the degrees of
1530freedom can be supplied to `PCBDDC` by using the following functions:
1531
1532- `PCBDDCSetDofsSplitting()`
1533- `PCBDDCSetLocalAdjacencyGraph()`
1534- `PCBDDCSetPrimalVerticesLocalIS()`
1535- `PCBDDCSetNeumannBoundaries()`
1536- `PCBDDCSetDirichletBoundaries()`
1537- `PCBDDCSetNeumannBoundariesLocal()`
1538- `PCBDDCSetDirichletBoundariesLocal()`
1539
1540Crucial for the convergence of the iterative process is the
1541specification of the primal constraints to be imposed at the interface
1542between subdomains. `PCBDDC` uses by default vertex continuities and
1543edge arithmetic averages, which are enough for the three-dimensional
1544Poisson problem with constant coefficients. The user can switch on and
1545off the usage of vertices, edges or face constraints by using the
1546command line switches `-pc_bddc_use_vertices`, `-pc_bddc_use_edges`,
1547`-pc_bddc_use_faces`. A customization of the constraints is available
1548by attaching a `MatNullSpace` object to the  matrix used to compute the preconditioner via
1549`MatSetNearNullSpace()`. The vectors of the `MatNullSpace` object
1550should represent the constraints in the form of quadrature rules;
1551quadrature rules for different classes of the interface can be listed in
1552the same vector. The number of vectors of the `MatNullSpace` object
1553corresponds to the maximum number of constraints that can be imposed for
1554each class. Once all the quadrature rules for a given interface class
1555have been extracted, an SVD operation is performed to retain the
1556non-singular modes. As an example, the rigid body modes represent an
1557effective choice for elasticity, even in the almost incompressible case.
1558For particular problems, e.g. edge-based discretization with Nedelec
1559elements, a user defined change of basis of the degrees of freedom can
1560be beneficial for `PCBDDC`; use `PCBDDCSetChangeOfBasisMat()` to
1561customize the change of basis.
1562
1563The `PCBDDC` method is usually robust with respect to jumps in the material
1564parameters aligned with the interface; for PDEs with more than one
1565material parameter you may also consider to use the so-called deluxe
1566scaling, available via the command line switch
1567`-pc_bddc_use_deluxe_scaling`. Other scalings are available, see
1568`PCISSetSubdomainScalingFactor()`,
1569`PCISSetSubdomainDiagonalScaling()` or
1570`PCISSetUseStiffnessScaling()`. However, the convergence properties of
1571the `PCBDDC` method degrades in presence of large jumps in the material
1572coefficients not aligned with the interface; for such cases, PETSc has
1573the capability of adaptively computing the primal constraints. Adaptive
1574selection of constraints could be requested by specifying a threshold
1575value at command line by using `-pc_bddc_adaptive_threshold x`. Valid
1576values for the threshold `x` ranges from 1 to infinity, with smaller
1577values corresponding to more robust preconditioners. For SPD problems in
15782D, or in 3D with only face degrees of freedom (like in the case of
1579Raviart-Thomas or Brezzi-Douglas-Marini elements), such a threshold is a
1580very accurate estimator of the condition number of the resulting
1581preconditioned operator. Since the adaptive selection of constraints for
1582`PCBDDC` methods is still an active topic of research, its implementation is
1583currently limited to SPD problems; moreover, because the technique
1584requires the explicit knowledge of the local Schur complements, it needs
1585the external package MUMPS.
1586
1587When solving problems decomposed in thousands of subdomains or more, the
1588solution of the `PCBDDC` coarse problem could become a bottleneck; in order
1589to overcome this issue, the user could either consider to solve the
1590parallel coarse problem on a subset of the communicator associated with
1591`PCBDDC` by using the command line switch
1592`-pc_bddc_coarse_redistribute`, or instead use a multilevel approach.
1593The latter can be requested by specifying the number of requested level
1594at command line (`-pc_bddc_levels`) or by using `PCBDDCSetLevels()`.
1595An additional parameter (see `PCBDDCSetCoarseningRatio()`) controls
1596the number of subdomains that will be generated at the next level; the
1597larger the coarsening ratio, the lower the number of coarser subdomains.
1598
1599For further details, see the example
1600<a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex59.c">KSP Tutorial ex59</a>
1601and the online documentation for `PCBDDC`.
1602
1603### Shell Preconditioners
1604
1605The shell preconditioner simply uses an application-provided routine to
1606implement the preconditioner. That is, it allows users to write or wrap their
1607own custom preconditioners as a `PC` and use it with `KSP`, etc.
1608
1609To provide a custom preconditioner application, use
1610
1611```
1612PCShellSetApply(PC pc,PetscErrorCode (*apply)(PC,Vec,Vec));
1613```
1614
1615Often a preconditioner needs access to an application-provided data
1616structured. For this, one should use
1617
1618```
1619PCShellSetContext(PC pc,void *ctx);
1620```
1621
1622to set this data structure and
1623
1624```
1625PCShellGetContext(PC pc,void *ctx);
1626```
1627
1628to retrieve it in `apply`. The three routine arguments of `apply()`
1629are the `PC`, the input vector, and the output vector, respectively.
1630
1631For a preconditioner that requires some sort of “setup” before being
1632used, that requires a new setup every time the operator is changed, one
1633can provide a routine that is called every time the operator is changed
1634(usually via `KSPSetOperators()`).
1635
1636```
1637PCShellSetSetUp(PC pc,PetscErrorCode (*setup)(PC));
1638```
1639
1640The argument to the `setup` routine is the same `PC` object which
1641can be used to obtain the operators with `PCGetOperators()` and the
1642application-provided data structure that was set with
1643`PCShellSetContext()`.
1644
1645(sec_combining_pcs)=
1646
1647### Combining Preconditioners
1648
1649The `PC` type `PCCOMPOSITE` allows one to form new preconditioners
1650by combining already-defined preconditioners and solvers. Combining
1651preconditioners usually requires some experimentation to find a
1652combination of preconditioners that works better than any single method.
1653It is a tricky business and is not recommended until your application
1654code is complete and running and you are trying to improve performance.
1655In many cases using a single preconditioner is better than a
1656combination; an exception is the multigrid/multilevel preconditioners
1657(solvers) that are always combinations of some sort, see {any}`sec_mg`.
1658
1659Let $B_1$ and $B_2$ represent the application of two
1660preconditioners of type `type1` and `type2`. The preconditioner
1661$B = B_1 + B_2$ can be obtained with
1662
1663```
1664PCSetType(pc,PCCOMPOSITE);
1665PCCompositeAddPCType(pc,type1);
1666PCCompositeAddPCType(pc,type2);
1667```
1668
1669Any number of preconditioners may added in this way.
1670
1671This way of combining preconditioners is called additive, since the
1672actions of the preconditioners are added together. This is the default
1673behavior. An alternative can be set with the option
1674
1675```
1676PCCompositeSetType(pc,PC_COMPOSITE_MULTIPLICATIVE);
1677```
1678
1679In this form the new residual is updated after the application of each
1680preconditioner and the next preconditioner applied to the next residual.
1681For example, with two composed preconditioners: $B_1$ and
1682$B_2$; $y = B x$ is obtained from
1683
1684$$
1685\begin{aligned}
1686y    = B_1 x \\
1687w_1  = x - A y \\
1688y    = y + B_2 w_1\end{aligned}
1689$$
1690
1691Loosely, this corresponds to a Gauss-Seidel iteration, while additive
1692corresponds to a Jacobi iteration.
1693
1694Under most circumstances, the multiplicative form requires one-half the
1695number of iterations as the additive form; however, the multiplicative
1696form does require the application of $A$ inside the
1697preconditioner.
1698
1699In the multiplicative version, the calculation of the residual inside
1700the preconditioner can be done in two ways: using the original linear
1701system matrix or using the matrix used to build the preconditioners
1702$B_1$, $B_2$, etc. By default it uses the “preconditioner
1703matrix”, to use the `Amat` matrix use the option
1704
1705```
1706PCSetUseAmat(PC pc);
1707```
1708
1709The individual preconditioners can be accessed (in order to set options)
1710via
1711
1712```
1713PCCompositeGetPC(PC pc,PetscInt count,PC *subpc);
1714```
1715
1716For example, to set the first sub preconditioners to use ILU(1)
1717
1718```
1719PC subpc;
1720PCCompositeGetPC(pc,0,&subpc);
1721PCFactorSetFill(subpc,1);
1722```
1723
1724One can also change the operator that is used to construct a particular
1725`PC` in the composite `PC` calling `PCSetOperators()` on the obtained `PC`.
1726`PCFIELDSPLIT`, {any}`sec_block_matrices`, provides an alternative approach to defining composite preconditioners
1727with a variety of pre-defined compositions.
1728
1729These various options can also be set via the options database. For
1730example, `-pc_type` `composite` `-pc_composite_pcs` `jacobi,ilu`
1731causes the composite preconditioner to be used with two preconditioners:
1732Jacobi and ILU. The option `-pc_composite_type` `multiplicative`
1733initiates the multiplicative version of the algorithm, while
1734`-pc_composite_type` `additive` the additive version. Using the
1735`Amat` matrix is obtained with the option `-pc_use_amat`. One sets
1736options for the sub-preconditioners with the extra prefix `-sub_N_`
1737where `N` is the number of the sub-preconditioner. For example,
1738`-sub_0_pc_ifactor_fill` `0`.
1739
1740PETSc also allows a preconditioner to be a complete `KSPSolve()` linear solver. This
1741is achieved with the `PCKSP` type.
1742
1743```
1744PCSetType(PC pc,PCKSP);
1745PCKSPGetKSP(pc,&ksp);
1746 /* set any KSP/PC options */
1747```
1748
1749From the command line one can use 5 iterations of biCG-stab with ILU(0)
1750preconditioning as the preconditioner with
1751`-pc_type ksp -ksp_pc_type ilu -ksp_ksp_max_it 5 -ksp_ksp_type bcgs`.
1752
1753By default the inner `KSP` solver uses the outer preconditioner
1754matrix, `Pmat`, as the matrix to be solved in the linear system; to
1755use the matrix that defines the linear system, `Amat` use the option
1756
1757```
1758PCSetUseAmat(PC pc);
1759```
1760
1761or at the command line with `-pc_use_amat`.
1762
1763Naturally, one can use a `PCKSP` preconditioner inside a composite
1764preconditioner. For example,
1765`-pc_type composite -pc_composite_pcs ilu,ksp -sub_1_pc_type jacobi -sub_1_ksp_max_it 10`
1766uses two preconditioners: ILU(0) and 10 iterations of GMRES with Jacobi
1767preconditioning. However, it is not clear whether one would ever wish to
1768do such a thing.
1769
1770(sec_mg)=
1771
1772### Multigrid Preconditioners
1773
1774A large suite of routines is available for using geometric multigrid as
1775a preconditioner [^id3]. In the `PC` framework, the user is required to
1776provide the coarse grid solver, smoothers, restriction and interpolation
1777operators, and code to calculate residuals. The `PC` package allows
1778these components to be encapsulated within a PETSc-compliant
1779preconditioner. We fully support both matrix-free and matrix-based
1780multigrid solvers.
1781
1782A multigrid preconditioner is created with the four commands
1783
1784```
1785KSPCreate(MPI_Comm comm,KSP *ksp);
1786KSPGetPC(KSP ksp,PC *pc);
1787PCSetType(PC pc,PCMG);
1788PCMGSetLevels(pc,PetscInt levels,MPI_Comm *comms);
1789```
1790
1791A large number of parameters affect the multigrid behavior. The command
1792
1793```
1794PCMGSetType(PC pc,PCMGType mode);
1795```
1796
1797indicates which form of multigrid to apply {cite}`1sbg`.
1798
1799For standard V or W-cycle multigrids, one sets the `mode` to be
1800`PC_MG_MULTIPLICATIVE`; for the additive form (which in certain cases
1801reduces to the BPX method, or additive multilevel Schwarz, or multilevel
1802diagonal scaling), one uses `PC_MG_ADDITIVE` as the `mode`. For a
1803variant of full multigrid, one can use `PC_MG_FULL`, and for the
1804Kaskade algorithm `PC_MG_KASKADE`. For the multiplicative and full
1805multigrid options, one can use a W-cycle by calling
1806
1807```
1808PCMGSetCycleType(PC pc,PCMGCycleType ctype);
1809```
1810
1811with a value of `PC_MG_CYCLE_W` for `ctype`. The commands above can
1812also be set from the options database. The option names are
1813`-pc_mg_type [multiplicative, additive, full, kaskade]`, and
1814`-pc_mg_cycle_type` `<ctype>`.
1815
1816The user can control the amount of smoothing by configuring the solvers
1817on the levels. By default, the up and down smoothers are identical. If
1818separate configuration of up and down smooths is required, it can be
1819requested with the option `-pc_mg_distinct_smoothup` or the routine
1820
1821```
1822PCMGSetDistinctSmoothUp(PC pc);
1823```
1824
1825The multigrid routines, which determine the solvers and
1826interpolation/restriction operators that are used, are mandatory. To set
1827the coarse grid solver, one must call
1828
1829```
1830PCMGGetCoarseSolve(PC pc,KSP *ksp);
1831```
1832
1833and set the appropriate options in `ksp`. Similarly, the smoothers are
1834controlled by first calling
1835
1836```
1837PCMGGetSmoother(PC pc,PetscInt level,KSP *ksp);
1838```
1839
1840and then setting the various options in the `ksp.` For example,
1841
1842```
1843PCMGGetSmoother(pc,1,&ksp);
1844KSPSetOperators(ksp,A1,A1);
1845```
1846
1847sets the matrix that defines the smoother on level 1 of the multigrid.
1848While
1849
1850```
1851PCMGGetSmoother(pc,1,&ksp);
1852KSPGetPC(ksp,&pc);
1853PCSetType(pc,PCSOR);
1854```
1855
1856sets SOR as the smoother to use on level 1.
1857
1858To use a different pre- or postsmoother, one should call the following
1859routines instead.
1860
1861```
1862PCMGGetSmootherUp(PC pc,PetscInt level,KSP *upksp);
1863PCMGGetSmootherDown(PC pc,PetscInt level,KSP *downksp);
1864```
1865
1866Use
1867
1868```
1869PCMGSetInterpolation(PC pc,PetscInt level,Mat P);
1870```
1871
1872and
1873
1874```
1875PCMGSetRestriction(PC pc,PetscInt level,Mat R);
1876```
1877
1878to define the intergrid transfer operations. If only one of these is
1879set, its transpose will be used for the other.
1880
1881It is possible for these interpolation operations to be matrix-free (see
1882{any}`sec_matrixfree`); One should then make
1883sure that these operations are defined for the (matrix-free) matrices
1884passed in. Note that this system is arranged so that if the
1885interpolation is the transpose of the restriction, you can pass the same
1886`mat` argument to both `PCMGSetRestriction()` and
1887`PCMGSetInterpolation()`.
1888
1889On each level except the coarsest, one must also set the routine to
1890compute the residual. The following command suffices:
1891
1892```
1893PCMGSetResidual(PC pc,PetscInt level,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat);
1894```
1895
1896The `residual()` function normally does not need to be set if one’s
1897operator is stored in `Mat` format. In certain circumstances, where it
1898is much cheaper to calculate the residual directly, rather than through
1899the usual formula $b - Ax$, the user may wish to provide an
1900alternative.
1901
1902Finally, the user may provide three work vectors for each level (except
1903on the finest, where only the residual work vector is required). The
1904work vectors are set with the commands
1905
1906```
1907PCMGSetRhs(PC pc,PetscInt level,Vec b);
1908PCMGSetX(PC pc,PetscInt level,Vec x);
1909PCMGSetR(PC pc,PetscInt level,Vec r);
1910```
1911
1912The `PC` references these vectors, so you should call `VecDestroy()`
1913when you are finished with them. If any of these vectors are not
1914provided, the preconditioner will allocate them.
1915
1916One can control the `KSP` and `PC` options used on the various
1917levels (as well as the coarse grid) using the prefix `mg_levels_`
1918(`mg_coarse_` for the coarse grid). For example,
1919`-mg_levels_ksp_type cg` will cause the CG method to be used as the
1920Krylov method for each level. Or
1921`-mg_levels_pc_type ilu -mg_levels_pc_factor_levels 2` will cause the
1922ILU preconditioner to be used on each level with two levels of fill in
1923the incomplete factorization.
1924
1925(sec_block_matrices)=
1926
1927## Solving Block Matrices with PCFIELDSPLIT
1928
1929Block matrices represent an important class of problems in numerical
1930linear algebra and offer the possibility of far more efficient iterative
1931solvers than just treating the entire matrix as a black box. In this
1932section, we use the common linear algebra definition of block matrices, where matrices are divided into a small, problem-size independent (two,
1933three, or so) number of very large blocks. These blocks arise naturally
1934from the underlying physics or discretization of the problem, such as the velocity and pressure. Under a certain numbering of
1935unknowns, the matrix can be written as
1936
1937$$
1938\left( \begin{array}{cccc}
1939A_{00}   & A_{01} & A_{02} & A_{03} \\
1940A_{10}   & A_{11} & A_{12} & A_{13} \\
1941A_{20}   & A_{21} & A_{22} & A_{23} \\
1942A_{30}   & A_{31} & A_{32} & A_{33} \\
1943\end{array} \right),
1944$$
1945
1946where each $A_{ij}$ is an entire block. The matrices on a parallel computer are not explicitly stored this way. Instead, each process will
1947own some rows of $A_{0*}$, $A_{1*}$ etc. On a
1948process, the blocks may be stored in one block followed by another
1949
1950$$
1951\left( \begin{array}{ccccccc}
1952A_{{00}_{00}}   & A_{{00}_{01}} & A_{{00}_{02}} & ... & A_{{01}_{00}} & A_{{01}_{01}} & ...  \\
1953A_{{00}_{10}}   & A_{{00}_{11}} & A_{{00}_{12}} & ... & A_{{01}_{10}} & A_{{01}_{11}} & ... \\
1954A_{{00}_{20}}   & A_{{00}_{21}} & A_{{00}_{22}} & ... & A_{{01}_{20}} & A_{{01}_{21}}  & ...\\
1955... \\
1956A_{{10}_{00}}   & A_{{10}_{01}} & A_{{10}_{02}} & ... & A_{{11}_{00}} & A_{{11}_{01}}  & ... \\
1957A_{{10}_{10}}   & A_{{10}_{11}} & A_{{10}_{12}} & ... & A_{{11}_{10}} & A_{{11}_{11}}  & ... \\
1958... \\
1959\end{array} \right)
1960$$
1961
1962or interlaced, for example, with four blocks
1963
1964$$
1965\left( \begin{array}{ccccc}
1966A_{{00}_{00}}   & A_{{01}_{00}} &  A_{{00}_{01}} & A_{{01}_{01}} &  ... \\
1967A_{{10}_{00}}   & A_{{11}_{00}} &  A_{{10}_{01}} & A_{{11}_{01}} &  ... \\
1968A_{{00}_{10}}   & A_{{01}_{10}} & A_{{00}_{11}} & A_{{01}_{11}} & ...\\
1969A_{{10}_{10}}   & A_{{11}_{10}} & A_{{10}_{11}} & A_{{11}_{11}} & ...\\
1970...
1971\end{array} \right).
1972$$
1973
1974Note that for interlaced storage, the number of rows/columns of each
1975block must be the same size. Matrices obtained with `DMCreateMatrix()`
1976where the `DM` is a `DMDA` are always stored interlaced. Block
1977matrices can also be stored using the `MATNEST` format, which holds
1978separate assembled blocks. Each of these nested matrices is itself
1979distributed in parallel. It is more efficient to use `MATNEST` with
1980the methods described in this section because there are fewer copies and
1981better formats (e.g., `MATBAIJ` or `MATSBAIJ`) can be used for the
1982components, but it is not possible to use many other methods with
1983`MATNEST`. See {any}`sec_matnest` for more on assembling
1984block matrices without depending on a specific matrix format.
1985
1986The PETSc `PCFIELDSPLIT` preconditioner implements the
1987“block” solvers in PETSc, {cite}`elman2008tcp`. There are three ways to provide the
1988information that defines the blocks. If the matrices are stored as
1989interlaced then `PCFieldSplitSetFields()` can be called repeatedly to
1990indicate which fields belong to each block. More generally
1991`PCFieldSplitSetIS()` can be used to indicate exactly which
1992rows/columns of the matrix belong to a particular block (field). You can provide
1993names for each block with these routines; if you do not, they are numbered from 0. With these two approaches, the blocks may
1994overlap (though they generally will not overlap). If only one block is defined,
1995then the complement of the matrices is used to define the other block.
1996Finally, the option `-pc_fieldsplit_detect_saddle_point` causes two
1997diagonal blocks to be found, one associated with all rows/columns that
1998have zeros on the diagonals and the rest.
1999
2000**Important parameters for PCFIELDSPLIT**
2001
2002- Control the fields used
2003
2004  - `-pc_fieldsplit_detect_saddle_point` \<bool:false> Generate two fields, the first consists of all rows with a nonzero on the diagonal, and the second will be all rows
2005    with zero on the diagonal. See `PCFieldSplitSetDetectSaddlePoint()`.
2006
2007  - `-pc_fieldsplit_dm_splits` \<bool:true> Use the `DM` attached to the preconditioner to determine the fields. See `PCFieldSplitSetDMSplits()` and
2008    `DMCreateFieldDecomposition()`.
2009
2010  - `-pc_fieldsplit_%d_fields` \<f1,f2,...:int> Use f1, f2, .. to define field `d`. The `fn` are in the range of 0, ..., bs-1 where bs is the block size
2011    of the matrix or set with `PCFieldSplitSetBlockSize()`. See `PCFieldSplitSetFields()`.
2012
2013    - `-pc_fieldsplit_default` \<bool:true> Automatically add any fields needed that have not been supplied explicitly by `-pc_fieldsplit_%d_fields`.
2014
2015  - `DMFieldsplitSetIS()` Provide the `IS` that defines a particular field.
2016
2017- Control the type of the block preconditioner
2018
2019  - `-pc_fieldsplit_type` \<additive|multiplicative|symmetric_multiplicative|schur|gkb:multiplicative> The order in which the field solves are applied.
2020    For symmetric problems where `KSPCG` is used `symmetric_multiplicative` must be used instead of `multiplicative`. `additive` is the least expensive
2021    to apply but provides the worst convergence. `schur` requires either a good preconditioner for the Schur complement or a naturally well-conditioned
2022    Schur complement, but when it works well can be extremely effective. See `PCFieldSplitSetType()`. `gkb` is for symmetric saddle-point problems (the lower-right
2023    the block is zero).
2024
2025  - `-pc_fieldsplit_diag_use_amat` \<bool:false> Use the first matrix that is passed to `KSPSetJacobian()` to construct the block-diagonal sub-matrices used in the algorithms,
2026    by default, the second matrix is used.
2027
2028  - Options for Schur preconditioner: `-pc_fieldsplit_type`
2029    `schur`
2030
2031    - `-pc_fieldsplit_schur_fact_type` \<diag|lower|upper|full:diag> See `PCFieldSplitSetSchurFactType()`. `full` reduces the iterations but each iteration requires additional
2032      field solves.
2033
2034    - `-pc_fieldsplit_schur_precondition` \<self|selfp|user|a11|full:user> How the Schur complement is preconditioned. See `PCFieldSplitSetSchurPre()`.
2035
2036      - `-fieldsplit_1_mat_schur_complement_ainv_type` \<diag|lump:diag> Use the lumped diagonal of $A_{00}$ when `-pc_fieldsplit_schur_precondition`
2037        `selfp` is used.
2038
2039    - `-pc_fieldsplit_schur_scale` \<real:-1.0> Controls the sign flip of S for `-pc_fieldsplit_schur_fact_type` `diag`.
2040      See `PCFieldSplitSetSchurScale()`
2041
2042    - `fieldsplit_1_xxx` controls the solver for the Schur complement system.
2043      If a `DM` provided the fields, use the second field name set in the `DM` instead of 1.
2044
2045      - `-fieldsplit_1_pc_type` `lsc` `-fieldsplit_1_lsc_pc_xxx` use
2046        the least squares commutators {cite}`elmanhowleshadidshuttleworthtuminaro2006` {cite}`silvester2001efficient`
2047        preconditioner for the Schur complement with any preconditioner for the least-squares matrix, see `PCLSC`.
2048        If a `DM` provided the fields, use the second field name set in the `DM` instead of 1.
2049
2050    - `-fieldsplit_upper_xxx` Set options for the solver in the upper solver when `-pc_fieldsplit_schur_fact_type`
2051      `upper` or `full` is used. Defaults to
2052      using the solver as provided with `-fieldsplit_0_xxx`.
2053
2054    - `-fieldsplit_1_inner_xxx` Set the options for the solver inside the application of the Schur complement;
2055      defaults to using the solver as provided with `-fieldsplit_0_xxx`. If a `DM` provides the fields use the name of the second field name set in the `DM` instead of 1.
2056
2057  - Options for GKB preconditioner: `-pc_fieldsplit_type` gkb
2058
2059    - `-pc_fieldsplit_gkb_tol` \<real:1e-5> See `PCFieldSplitSetGKBTol()`.
2060    - `-pc_fieldsplit_gkb_delay` \<int:5> See `PCFieldSplitSetGKBDelay()`.
2061    - `-pc_fieldsplit_gkb_nu` \<real:1.0> See `PCFieldSplitSetGKBNu()`.
2062    - `-pc_fieldsplit_gkb_maxit` \<int:100> See `PCFieldSplitSetGKBMaxit()`.
2063    - `-pc_fieldsplit_gkb_monitor` \<bool:false> Monitor the convergence of the inner solver.
2064
2065- Options for additive and multiplication field solvers:
2066
2067   - `-fieldsplit_%d_xxx` Set options for the solver for field number `d`. For example, `-fieldsplit_0_pc_type`
2068    `jacobi`. When the fields are obtained from a `DM` use the
2069    field name instead of `d`.
2070
2071For simplicity, we restrict our matrices to two-by-two blocks in the rest of the section. So the matrix is
2072
2073$$
2074\left( \begin{array}{cc}
2075A_{00}   & A_{01} \\
2076A_{10}   & A_{11} \\
2077\end{array} \right).
2078$$
2079
2080On occasion, the user may provide another matrix that is used to
2081construct parts of the preconditioner
2082
2083$$
2084\left( \begin{array}{cc}
2085Ap_{00}   & Ap_{01} \\
2086Ap_{10}   & Ap_{11} \\
2087\end{array} \right).
2088$$
2089
2090For notational simplicity define $\text{ksp}(A,Ap)$ to mean
2091approximately solving a linear system using `KSP` with the operator
2092$A$ and preconditioner built from matrix $Ap$.
2093
2094For matrices defined with any number of blocks, there are three “block”
2095algorithms available: block Jacobi,
2096
2097$$
2098\left( \begin{array}{cc}
2099  \text{ksp}(A_{00},Ap_{00})   & 0 \\
2100  0   & \text{ksp}(A_{11},Ap_{11}) \\
2101\end{array} \right)
2102$$
2103
2104block Gauss-Seidel,
2105
2106$$
2107\left( \begin{array}{cc}
2108I   & 0 \\
21090 & A^{-1}_{11} \\
2110\end{array} \right)
2111\left( \begin{array}{cc}
2112I   & 0 \\
2113-A_{10} & I \\
2114\end{array} \right)
2115\left( \begin{array}{cc}
2116A^{-1}_{00}   & 0 \\
21170 & I \\
2118\end{array} \right)
2119$$
2120
2121which is implemented [^id4] as
2122
2123$$
2124\left( \begin{array}{cc}
2125I   & 0 \\
2126  0 & \text{ksp}(A_{11},Ap_{11}) \\
2127\end{array} \right)
2128$$
2129
2130$$
2131\left[
2132\left( \begin{array}{cc}
21330   & 0 \\
21340 & I \\
2135\end{array} \right) +
2136\left( \begin{array}{cc}
2137I   & 0 \\
2138-A_{10} & -A_{11} \\
2139\end{array} \right)
2140\left( \begin{array}{cc}
2141I   & 0 \\
21420 & 0 \\
2143\end{array} \right)
2144\right]
2145$$
2146
2147$$
2148\left( \begin{array}{cc}
2149  \text{ksp}(A_{00},Ap_{00})   & 0 \\
21500 & I \\
2151\end{array} \right)
2152$$
2153
2154and symmetric block Gauss-Seidel
2155
2156$$
2157\left( \begin{array}{cc}
2158A_{00}^{-1}   & 0 \\
21590 & I \\
2160\end{array} \right)
2161\left( \begin{array}{cc}
2162I   & -A_{01} \\
21630 & I \\
2164\end{array} \right)
2165\left( \begin{array}{cc}
2166A_{00}   & 0 \\
21670 & A_{11}^{-1} \\
2168\end{array} \right)
2169\left( \begin{array}{cc}
2170I   & 0 \\
2171-A_{10} & I \\
2172\end{array} \right)
2173\left( \begin{array}{cc}
2174A_{00}^{-1}   & 0 \\
21750 & I \\
2176\end{array} \right).
2177$$
2178
2179These can be accessed with
2180`-pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative>`
2181or the function `PCFieldSplitSetType()`. The option prefixes for the
2182internal KSPs are given by `-fieldsplit_name_`.
2183
2184By default blocks $A_{00}, A_{01}$ and so on are extracted out of
2185`Pmat`, the matrix that the `KSP` uses to build the preconditioner,
2186and not out of `Amat` (i.e., $A$ itself). As discussed above, in
2187{any}`sec_combining_pcs`, however, it is
2188possible to use `Amat` instead of `Pmat` by calling
2189`PCSetUseAmat(pc)` or using `-pc_use_amat` on the command line.
2190Alternatively, you can have `PCFIELDSPLIT` extract the diagonal blocks
2191$A_{00}, A_{11}$ etc. out of `Amat` by calling
2192`PCFieldSplitSetDiagUseAmat(pc,PETSC_TRUE)` or supplying command-line
2193argument `-pc_fieldsplit_diag_use_amat`. Similarly,
2194`PCFieldSplitSetOffDiagUseAmat(pc,{PETSC_TRUE`) or
2195`-pc_fieldsplit_off_diag_use_amat` will cause the off-diagonal blocks
2196$A_{01},A_{10}$ etc. to be extracted out of `Amat`.
2197
2198For two-by-two blocks only, there is another family of solvers based on
2199Schur complements. The inverse of the Schur complement factorization is
2200
2201$$
2202\left[
2203\left( \begin{array}{cc}
2204I   & 0 \\
2205A_{10}A_{00}^{-1} & I \\
2206\end{array} \right)
2207\left( \begin{array}{cc}
2208A_{00}  & 0 \\
22090 & S \\
2210\end{array} \right)
2211\left( \begin{array}{cc}
2212I   & A_{00}^{-1} A_{01} \\
22130 & I \\
2214\end{array} \right)
2215\right]^{-1} =
2216$$
2217
2218$$
2219\left( \begin{array}{cc}
2220I   & A_{00}^{-1} A_{01} \\
22210 & I \\
2222\end{array} \right)^{-1}
2223\left( \begin{array}{cc}
2224A_{00}^{-1}  & 0 \\
22250 & S^{-1} \\
2226\end{array} \right)
2227\left( \begin{array}{cc}
2228I   & 0 \\
2229A_{10}A_{00}^{-1} & I \\
2230\end{array} \right)^{-1} =
2231$$
2232
2233$$
2234\left( \begin{array}{cc}
2235I   & -A_{00}^{-1} A_{01} \\
22360 & I \\
2237\end{array} \right)
2238\left( \begin{array}{cc}
2239A_{00}^{-1}  & 0 \\
22400 & S^{-1} \\
2241\end{array} \right)
2242\left( \begin{array}{cc}
2243I   & 0 \\
2244-A_{10}A_{00}^{-1} & I \\
2245\end{array} \right) =
2246$$
2247
2248$$
2249\left( \begin{array}{cc}
2250A_{00}^{-1}   & 0 \\
22510 & I \\
2252\end{array} \right)
2253\left( \begin{array}{cc}
2254I   & -A_{01} \\
22550 & I \\
2256\end{array} \right)
2257\left( \begin{array}{cc}
2258A_{00} & 0 \\
22590 & S^{-1} \\
2260\end{array} \right)
2261\left( \begin{array}{cc}
2262I   & 0 \\
2263-A_{10} & I \\
2264\end{array} \right)
2265\left( \begin{array}{cc}
2266A_{00}^{-1}   & 0 \\
22670 & I \\
2268\end{array} \right).
2269$$
2270
2271The preconditioner is accessed with `-pc_fieldsplit_type` `schur` and is
2272implemented as
2273
2274$$
2275\left( \begin{array}{cc}
2276  \text{ksp}(A_{00},Ap_{00})   & 0 \\
22770 & I \\
2278\end{array} \right)
2279\left( \begin{array}{cc}
2280I   & -A_{01} \\
22810 & I \\
2282\end{array} \right)
2283$$
2284
2285$$
2286\left( \begin{array}{cc}
2287I  & 0 \\
2288  0 & \text{ksp}(\hat{S},\hat{S}p) \\
2289\end{array} \right)
2290\left( \begin{array}{cc}
2291I   & 0 \\
2292  -A_{10} \text{ksp}(A_{00},Ap_{00}) & I \\
2293\end{array} \right).
2294$$
2295
2296Where
2297$\hat{S} = A_{11} - A_{10} \text{ksp}(A_{00},Ap_{00}) A_{01}$ is
2298the approximate Schur complement.
2299
2300There are several variants of the Schur complement preconditioner
2301obtained by dropping some of the terms; these can be obtained with
2302`-pc_fieldsplit_schur_fact_type <diag,lower,upper,full>` or the
2303function `PCFieldSplitSetSchurFactType()`. Note that the `diag` form
2304uses the preconditioner
2305
2306$$
2307\left( \begin{array}{cc}
2308  \text{ksp}(A_{00},Ap_{00})   & 0 \\
2309  0 & -\text{ksp}(\hat{S},\hat{S}p) \\
2310\end{array} \right).
2311$$
2312
2313This is done to ensure the preconditioner is positive definite for a
2314a common class of problems, saddle points with a positive definite
2315$A_{00}$: for these, the Schur complement is negative definite.
2316
2317The effectiveness of the Schur complement preconditioner depends on the
2318availability of a good preconditioner $\hat Sp$ for the Schur
2319complement matrix. In general, you are responsible for supplying
2320$\hat Sp$ via
2321`PCFieldSplitSetSchurPre(pc,PC_FIELDSPLIT_SCHUR_PRE_USER,Sp)`.
2322Without a good problem-specific $\hat Sp$, you can use
2323some built-in options.
2324
2325Using `-pc_fieldsplit_schur_precondition user` on the command line
2326activates the matrix supplied programmatically, as explained above.
2327
2328With `-pc_fieldsplit_schur_precondition a11` (default)
2329$\hat Sp = A_{11}$ is used to build a preconditioner for
2330$\hat S$.
2331
2332Otherwise, `-pc_fieldsplit_schur_precondition self` will set
2333$\hat Sp = \hat S$ and use the Schur complement matrix itself to
2334build the preconditioner.
2335
2336The problem with the last approach is that $\hat S$ is used in
2337the unassembled, matrix-free form, and many preconditioners (e.g., ILU)
2338cannot be built out of such matrices. Instead, you can *assemble* an
2339approximation to $\hat S$ by inverting $A_{00}$, but only
2340approximately, to ensure the sparsity of $\hat Sp$ as much
2341as possible. Specifically, using
2342`-pc_fieldsplit_schur_precondition selfp` will assemble
2343$\hat Sp = A_{11} - A_{10} \text{inv}(A_{00}) A_{01}$.
2344
2345By default $\text{inv}(A_{00})$ is the inverse of the diagonal of
2346$A_{00}$, but using
2347`-fieldsplit_1_mat_schur_complement_ainv_type lump` will lump
2348$A_{00}$ first. Using
2349`-fieldsplit_1_mat_schur_complement_ainv_type blockdiag` will use the
2350inverse of the block diagonal of $A_{00}$. Option
2351`-mat_schur_complement_ainv_type` applies to any matrix of
2352`MatSchurComplement` type and here it is used with the prefix
2353`-fieldsplit_1` of the linear system in the second split.
2354
2355Finally, you can use the `PCLSC` preconditioner for the Schur
2356complement with `-pc_fieldsplit_type schur -fieldsplit_1_pc_type lsc`.
2357This uses for the preconditioner to $\hat{S}$ the operator
2358
2359$$
2360\text{ksp}(A_{10} A_{01},A_{10} A_{01}) A_{10} A_{00} A_{01} \text{ksp}(A_{10} A_{01},A_{10} A_{01})
2361$$
2362
2363Which, of course, introduces two additional inner solves for each
2364application of the Schur complement. The options prefix for this inner
2365`KSP` is `-fieldsplit_1_lsc_`. Instead of constructing the matrix
2366$A_{10} A_{01}$, users can provide their own matrix. This is
2367done by attaching the matrix/matrices to the $Sp$ matrix they
2368provide with
2369
2370```
2371PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L);
2372PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp);
2373```
2374
2375(sec_singular)=
2376
2377## Solving Singular Systems
2378
2379Sometimes one is required to solver singular linear systems. In this
2380case, the system matrix has a nontrivial null space. For example, the
2381discretization of the Laplacian operator with Neumann boundary
2382conditions has a null space of the constant functions. PETSc has tools
2383to help solve these systems. This approach is only guaranteed to work for left preconditioning (see `KSPSetPCSide()`); for example it
2384may not work in some situations with `KSPFGMRES`.
2385
2386First, one must know what the null space is and store it using an
2387orthonormal basis in an array of PETSc Vecs. The constant functions can
2388be handled separately, since they are such a common case. Create a
2389`MatNullSpace` object with the command
2390
2391```
2392MatNullSpaceCreate(MPI_Comm,PetscBool hasconstants,PetscInt dim,Vec *basis,MatNullSpace *nsp);
2393```
2394
2395Here, `dim` is the number of vectors in `basis` and `hasconstants`
2396indicates if the null space contains the constant functions. If the null
2397space contains the constant functions you do not need to include it in
2398the `basis` vectors you provide, nor in the count `dim`.
2399
2400One then tells the `KSP` object you are using what the null space is
2401with the call
2402
2403```
2404MatSetNullSpace(Mat Amat,MatNullSpace nsp);
2405```
2406
2407The `Amat` should be the *first* matrix argument used with
2408`KSPSetOperators()`, `SNESSetJacobian()`, or `TSSetIJacobian()`.
2409The PETSc solvers will now
2410handle the null space during the solution process.
2411
2412If the right-hand side of linear system is not in the range of `Amat`, that is it is not
2413orthogonal to the null space of `Amat` transpose, then the residual
2414norm of the Krylov iteration will not converge to zero; it will converge to a non-zero value while the
2415solution is converging to the least squares solution of the linear system. One can, if one desires,
2416apply `MatNullSpaceRemove()` with the null space of `Amat` transpose to the right-hand side before calling
2417`KSPSolve()`. Then the residual norm will converge to zero.
2418
2419If one chooses a direct solver (or an incomplete factorization) it may
2420still detect a zero pivot. You can run with the additional options or
2421`-pc_factor_shift_type NONZERO`
2422`-pc_factor_shift_amount  <dampingfactor>` to prevent the zero pivot.
2423A good choice for the `dampingfactor` is 1.e-10.
2424
2425If the matrix is non-symmetric and you wish to solve the transposed linear system
2426you must provide the null space of the transposed matrix with `MatSetTransposeNullSpace()`.
2427
2428(sec_externalsol)=
2429
2430## Using External Linear Solvers
2431
2432PETSc interfaces to several external linear solvers (also see {any}`acknowledgements`).
2433To use these solvers, one may:
2434
24351. Run `configure` with the additional options
2436   `--download-packagename` e.g. `--download-superlu_dist`
2437   `--download-parmetis` (SuperLU_DIST needs ParMetis) or
2438   `--download-mumps` `--download-scalapack` (MUMPS requires
2439   ScaLAPACK).
24402. Build the PETSc libraries.
24413. Use the runtime option: `-ksp_type preonly` (or equivalently `-ksp_type none`) `-pc_type <pctype>`
2442   `-pc_factor_mat_solver_type <packagename>`. For eg:
2443   `-ksp_type preonly` `-pc_type lu`
2444   `-pc_factor_mat_solver_type superlu_dist`.
2445
2446```{eval-rst}
2447.. list-table:: Options for External Solvers
2448   :name: tab-externaloptions
2449   :header-rows: 1
2450
2451   * - MatType
2452     - PCType
2453     - MatSolverType
2454     - Package
2455   * - ``seqaij``
2456     - ``lu``
2457     - ``MATSOLVERESSL``
2458     - ``essl``
2459   * - ``seqaij``
2460     - ``lu``
2461     - ``MATSOLVERLUSOL``
2462     -  ``lusol``
2463   * - ``seqaij``
2464     - ``lu``
2465     - ``MATSOLVERMATLAB``
2466     - ``matlab``
2467   * - ``aij``
2468     - ``lu``
2469     - ``MATSOLVERMUMPS``
2470     - ``mumps``
2471   * - ``aij``
2472     - ``cholesky``
2473     - -
2474     - -
2475   * - ``sbaij``
2476     - ``cholesky``
2477     - -
2478     - -
2479   * - ``seqaij``
2480     - ``lu``
2481     - ``MATSOLVERSUPERLU``
2482     - ``superlu``
2483   * - ``aij``
2484     - ``lu``
2485     - ``MATSOLVERSUPERLU_DIST``
2486     - ``superlu_dist``
2487   * - ``seqaij``
2488     - ``lu``
2489     - ``MATSOLVERUMFPACK``
2490     - ``umfpack``
2491   * - ``seqaij``
2492     - ``cholesky``
2493     - ``MATSOLVERCHOLMOD``
2494     - ``cholmod``
2495   * - ``seqaij``
2496     - ``lu``
2497     - ``MATSOLVERKLU``
2498     -  ``klu``
2499   * - ``dense``
2500     - ``lu``
2501     - ``MATSOLVERELEMENTAL``
2502     -  ``elemental``
2503   * - ``dense``
2504     - ``cholesky``
2505     - -
2506     - -
2507   * - ``seqaij``
2508     - ``lu``
2509     - ``MATSOLVERMKL_PARDISO``
2510     - ``mkl_pardiso``
2511   * - ``aij``
2512     - ``lu``
2513     - ``MATSOLVERMKL_CPARDISO``
2514     - ``mkl_cpardiso``
2515   * - ``aij``
2516     - ``lu``
2517     - ``MATSOLVERPASTIX``
2518     -  ``pastix``
2519   * - ``aij``
2520     - ``cholesky``
2521     - ``MATSOLVERBAS``
2522     -  ``bas``
2523   * - ``aijcusparse``
2524     - ``lu``
2525     - ``MATSOLVERCUSPARSE``
2526     - ``cusparse``
2527   * - ``aijcusparse``
2528     - ``cholesky``
2529     -  -
2530     -  -
2531   * - ``aij``
2532     - ``lu``, ``cholesky``
2533     - ``MATSOLVERPETSC``
2534     - ``petsc``
2535   * - ``baij``
2536     - -
2537     - -
2538     - -
2539   * - ``aijcrl``
2540     - -
2541     - -
2542     - -
2543   * - ``aijperm``
2544     - -
2545     - -
2546     - -
2547   * - ``seqdense``
2548     - -
2549     - -
2550     - -
2551   * - ``aij``
2552     - -
2553     - -
2554     - -
2555   * - ``baij``
2556     - -
2557     - -
2558     - -
2559   * - ``aijcrl``
2560     - -
2561     - -
2562     - -
2563   * - ``aijperm``
2564     - -
2565     - -
2566     - -
2567   * - ``seqdense``
2568     - -
2569     - -
2570     - -
2571```
2572
2573The default and available input options for each external software can
2574be found by specifying `-help` at runtime.
2575
2576As an alternative to using runtime flags to employ these external
2577packages, procedural calls are provided for some packages. For example,
2578the following procedural calls are equivalent to runtime options
2579`-ksp_type preonly` (or equivalently `-ksp_type none`) `-pc_type lu`
2580`-pc_factor_mat_solver_type mumps` `-mat_mumps_icntl_7 3`:
2581
2582```
2583KSPSetType(ksp, KSPPREONLY); // (or equivalently KSPSetType(ksp, KSPNONE))
2584KSPGetPC(ksp, &pc);
2585PCSetType(pc, PCLU);
2586PCFactorSetMatSolverType(pc, MATSOLVERMUMPS);
2587PCFactorSetUpMatSolverType(pc);
2588PCFactorGetMatrix(pc, &F);
2589icntl=7; ival = 3;
2590MatMumpsSetIcntl(F, icntl, ival);
2591```
2592
2593One can also create matrices with the appropriate capabilities by
2594calling `MatCreate()` followed by `MatSetType()` specifying the
2595desired matrix type from {any}`tab-externaloptions`. These
2596matrix types inherit capabilities from their PETSc matrix parents:
2597`MATSEQAIJ`, `MATMPIAIJ`, etc. As a result, the preallocation routines
2598`MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, etc.
2599and any other type specific routines of the base class are supported.
2600One can also call `MatConvert()` inplace to convert the matrix to and
2601from its base class without performing an expensive data copy.
2602`MatConvert()` cannot be called on matrices that have already been
2603factored.
2604
2605In {any}`tab-externaloptions`, the base class `aij` refers
2606to the fact that inheritance is based on `MATSEQAIJ` when constructed
2607with a single process communicator, and from `MATMPIAIJ` otherwise.
2608The same holds for `baij` and `sbaij`. For codes that are intended
2609to be run as both a single process or with multiple processes, depending
2610on the `mpiexec` command, it is recommended that both sets of
2611preallocation routines are called for these communicator morphing types.
2612The call for the incorrect type will simply be ignored without any harm
2613or message.
2614
2615(sec_pcmpi)=
2616
2617## Using PETSc's MPI parallel linear solvers from a non-MPI program
2618
2619Using PETSc's MPI linear solver server it is possible to use multiple MPI processes to solve a
2620a linear system when the application code, including the matrix generation, is run on a single
2621MPI process (with or without OpenMP). The application code must be built with MPI and must call
2622`PetscInitialize()` at the very beginning of the program and end with `PetscFinalize()`. The
2623application code may utilize OpenMP.
2624The code may create multiple matrices and `KSP` objects and call `KSPSolve()`, similarly the
2625code may utilize the `SNES` nonlinear solvers, the `TS` ODE integrators, and the `Tao` optimization algorithms
2626which use `KSP`.
2627
2628The program must then be launched using the standard approaches for launching MPI programs with the additional
2629PETSc option `-mpi_linear_solver_server`. The linear solves are controlled via the options database in the usual manner (using any options prefix
2630you may have provided via `KSPSetOptionsPrefix()`, for example `-ksp_type cg -ksp_monitor -pc_type bjacobi -ksp_view`. The solver options cannot be set via
2631the functional interface, for example `KSPSetType()` etc.
2632
2633The option `-mpi_linear_solver_server_view` will print
2634a summary of all the systems solved by the MPI linear solver server when the program completes. By default the linear solver server
2635will only parallelize the linear solve to the extent that it believes is appropriate to obtain speedup for the parallel solve, for example, if the
2636matrix has 1,000 rows and columns the solution will not be parallelized by default. One can use the option `-mpi_linear_solver_server_minimum_count_per_rank 5000`
2637to cause the linear solver server to allow as few as 5,000 unknowns per MPI process in the parallel solve.
2638
2639See `PCMPI`, `PCMPIServerBegin()`, and `PCMPIServerEnd()` for more details on the solvers.
2640
2641For help when anything goes wrong with the MPI linear solver server see `PCMPIServerBegin()`.
2642
2643Amdahl's law makes clear that parallelizing only a portion of a numerical code can only provide a limited improvement
2644in the computation time; thus it is crucial to understand what phases of a computation must be parallelized (via MPI, OpenMP, or some other model)
2645to ensure a useful increase in performance. One of the crucial phases is likely the generation of the matrix entries; the
2646use of `MatSetPreallocationCOO()` and `MatSetValuesCOO()` in an OpenMP code allows parallelizing the generation of the matrix.
2647
2648See {any}`sec_pcmpi_study` for a study of the use of `PCMPI` on a specific PETSc application.
2649
2650```{rubric} Footnotes
2651```
2652
2653[^id3]: See {any}`sec_amg` for information on using algebraic multigrid.
2654
2655[^id4]: This may seem an odd way to implement since it involves the "extra"
2656    multiply by $-A_{11}$. The reason is this is implemented this
2657    way is that this approach works for any number of blocks that may
2658    overlap.
2659
2660```{rubric} References
2661```
2662
2663```{eval-rst}
2664.. bibliography:: /petsc.bib
2665   :filter: docname in docnames
2666```
2667