1(ch_snes)= 2 3# SNES: Nonlinear Solvers 4 5The solution of large-scale nonlinear problems pervades many facets of 6computational science and demands robust and flexible solution 7strategies. The `SNES` library of PETSc provides a powerful suite of 8data-structure-neutral numerical routines for such problems. Built on 9top of the linear solvers and data structures discussed in preceding 10chapters, `SNES` enables the user to easily customize the nonlinear 11solvers according to the application at hand. Also, the `SNES` 12interface is *identical* for the uniprocess and parallel cases; the only 13difference in the parallel version is that each process typically forms 14only its local contribution to various matrices and vectors. 15 16The `SNES` class includes methods for solving systems of nonlinear 17equations of the form 18 19$$ 20\mathbf{F}(\mathbf{x}) = 0, 21$$ (fx0) 22 23where $\mathbf{F}: \, \Re^n \to \Re^n$. Newton-like methods provide the 24core of the package, including both line search and trust region 25techniques. A suite of nonlinear Krylov methods and methods based upon 26problem decomposition are also included. The solvers are discussed 27further in {any}`sec_nlsolvers`. Following the PETSc design 28philosophy, the interfaces to the various solvers are all virtually 29identical. In addition, the `SNES` software is completely flexible, so 30that the user can at runtime change any facet of the solution process. 31 32PETSc’s default method for solving the nonlinear equation is Newton’s 33method with line search, `SNESNEWTONLS`. The general form of the $n$-dimensional Newton’s method 34for solving {math:numref}`fx0` is 35 36$$ 37\mathbf{x}_{k+1} = \mathbf{x}_k - \mathbf{J}(\mathbf{x}_k)^{-1} \mathbf{F}(\mathbf{x}_k), \;\; k=0,1, \ldots, 38$$ (newton) 39 40where $\mathbf{x}_0$ is an initial approximation to the solution and 41$\mathbf{J}(\mathbf{x}_k) = \mathbf{F}'(\mathbf{x}_k)$, the Jacobian, is nonsingular at each 42iteration. In practice, the Newton iteration {math:numref}`newton` is 43implemented by the following two steps: 44 45$$ 46\begin{aligned} 471. & \text{(Approximately) solve} & \mathbf{J}(\mathbf{x}_k) \Delta \mathbf{x}_k &= -\mathbf{F}(\mathbf{x}_k). \\ 482. & \text{Update} & \mathbf{x}_{k+1} &\gets \mathbf{x}_k + \Delta \mathbf{x}_k. 49\end{aligned} 50$$ 51 52Other defect-correction algorithms can be implemented by using different 53choices for $J(\mathbf{x}_k)$. 54 55(sec_snesusage)= 56 57## Basic SNES Usage 58 59In the simplest usage of the nonlinear solvers, the user must merely 60provide a C, C++, Fortran, or Python routine to evaluate the nonlinear function 61{math:numref}`fx0`. The corresponding Jacobian matrix 62can be approximated with finite differences. For codes that are 63typically more efficient and accurate, the user can provide a routine to 64compute the Jacobian; details regarding these application-provided 65routines are discussed below. To provide an overview of the use of the 66nonlinear solvers, browse the concrete example in {ref}`ex1.c <snes-ex1>` or skip ahead to the discussion. 67 68(snes_ex1)= 69 70:::{admonition} Listing: `src/snes/tutorials/ex1.c` 71```{literalinclude} /../src/snes/tutorials/ex1.c 72:end-before: /*TEST 73``` 74::: 75 76To create a `SNES` solver, one must first call `SNESCreate()` as 77follows: 78 79``` 80SNESCreate(MPI_Comm comm,SNES *snes); 81``` 82 83The user must then set routines for evaluating the residual function {math:numref}`fx0` 84and, *possibly*, its associated Jacobian matrix, as 85discussed in the following sections. 86 87To choose a nonlinear solution method, the user can either call 88 89``` 90SNESSetType(SNES snes,SNESType method); 91``` 92 93or use the option `-snes_type <method>`, where details regarding the 94available methods are presented in {any}`sec_nlsolvers`. The 95application code can take complete control of the linear and nonlinear 96techniques used in the Newton-like method by calling 97 98``` 99SNESSetFromOptions(snes); 100``` 101 102This routine provides an interface to the PETSc options database, so 103that at runtime the user can select a particular nonlinear solver, set 104various parameters and customized routines (e.g., specialized line 105search variants), prescribe the convergence tolerance, and set 106monitoring routines. With this routine the user can also control all 107linear solver options in the `KSP`, and `PC` modules, as discussed 108in {any}`ch_ksp`. 109 110After having set these routines and options, the user solves the problem 111by calling 112 113``` 114SNESSolve(SNES snes,Vec b,Vec x); 115``` 116 117where `x` should be initialized to the initial guess before calling and contains the solution on return. 118In particular, to employ an initial guess of 119zero, the user should explicitly set this vector to zero by calling 120`VecZeroEntries(x)`. Finally, after solving the nonlinear system (or several 121systems), the user should destroy the `SNES` context with 122 123``` 124SNESDestroy(SNES *snes); 125``` 126 127(sec_snesfunction)= 128 129### Nonlinear Function Evaluation 130 131When solving a system of nonlinear equations, the user must provide a 132a residual function {math:numref}`fx0`, which is set using 133 134``` 135SNESSetFunction(SNES snes,Vec f,PetscErrorCode (*FormFunction)(SNES snes,Vec x,Vec f,void *ctx),void *ctx); 136``` 137 138The argument `f` is an optional vector for storing the solution; pass `NULL` to have the `SNES` allocate it for you. 139The argument `ctx` is an optional user-defined context, which can 140store any private, application-specific data required by the function 141evaluation routine; `NULL` should be used if such information is not 142needed. In C and C++, a user-defined context is merely a structure in 143which various objects can be stashed; in Fortran a user context can be 144an integer array that contains both parameters and pointers to PETSc 145objects. 146<a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex5.c.html">SNES Tutorial ex5</a> 147and 148<a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex5f90.F90.html">SNES Tutorial ex5f90</a> 149give examples of user-defined application contexts in C and Fortran, 150respectively. 151 152(sec_snesjacobian)= 153 154### Jacobian Evaluation 155 156The user may also specify a routine to form some approximation of the 157Jacobian matrix, `A`, at the current iterate, `x`, as is typically 158done with 159 160``` 161SNESSetJacobian(SNES snes,Mat Amat,Mat Pmat,PetscErrorCode (*FormJacobian)(SNES snes,Vec x,Mat A,Mat B,void *ctx),void *ctx); 162``` 163 164The arguments of the routine `FormJacobian()` are the current iterate, 165`x`; the (approximate) Jacobian matrix, `Amat`; the matrix from 166which the preconditioner is constructed, `Pmat` (which is usually the 167same as `Amat`); and an optional user-defined Jacobian context, 168`ctx`, for application-specific data. The `FormJacobian()` 169callback is only invoked if the solver requires it, always 170*after* `FormFunction()` has been called at the current iterate. 171 172Note that the `SNES` solvers 173are all data-structure neutral, so the full range of PETSc matrix 174formats (including “matrix-free” methods) can be used. 175{any}`ch_matrices` discusses information regarding 176available matrix formats and options, while {any}`sec_nlmatrixfree` focuses on matrix-free methods in 177`SNES`. We briefly touch on a few details of matrix usage that are 178particularly important for efficient use of the nonlinear solvers. 179 180A common usage paradigm is to assemble the problem Jacobian in the 181preconditioner storage `B`, rather than `A`. In the case where they 182are identical, as in many simulations, this makes no difference. 183However, it allows us to check the analytic Jacobian we construct in 184`FormJacobian()` by passing the `-snes_mf_operator` flag. This 185causes PETSc to approximate the Jacobian using finite differencing of 186the function evaluation (discussed in {any}`sec_fdmatrix`), 187and the analytic Jacobian becomes merely the preconditioner. Even if the 188analytic Jacobian is incorrect, it is likely that the finite difference 189approximation will converge, and thus this is an excellent method to 190verify the analytic Jacobian. Moreover, if the analytic Jacobian is 191incomplete (some terms are missing or approximate), 192`-snes_mf_operator` may be used to obtain the exact solution, where 193the Jacobian approximation has been transferred to the preconditioner. 194 195One such approximate Jacobian comes from “Picard linearization”, use `SNESSetPicard()`, which 196writes the nonlinear system as 197 198$$ 199\mathbf{F}(\mathbf{x}) := \mathbf{A}(\mathbf{x}) \mathbf{x} - \mathbf{b} = 0 200$$ 201 202where $\mathbf{A}(\mathbf{x})$ usually contains the lower-derivative parts of the 203equation. For example, the nonlinear diffusion problem 204 205$$ 206- \nabla\cdot(\kappa(u) \nabla u) = 0 207$$ 208 209would be linearized as 210 211$$ 212A(u) v \simeq -\nabla\cdot(\kappa(u) \nabla v). 213$$ 214 215Usually this linearization is simpler to implement than Newton and the 216linear problems are somewhat easier to solve. In addition to using 217`-snes_mf_operator` with this approximation to the Jacobian, the 218Picard iterative procedure can be performed by defining $\mathbf{J}(\mathbf{x})$ 219to be $\mathbf{A}(\mathbf{x})$. Sometimes this iteration exhibits better global 220convergence than Newton linearization. 221 222During successive calls to `FormJacobian()`, the user can either 223insert new matrix contexts or reuse old ones, depending on the 224application requirements. For many sparse matrix formats, reusing the 225old space (and merely changing the matrix elements) is more efficient; 226however, if the matrix nonzero structure completely changes, creating an 227entirely new matrix context may be preferable. Upon subsequent calls to 228the `FormJacobian()` routine, the user may wish to reinitialize the 229matrix entries to zero by calling `MatZeroEntries()`. See 230{any}`sec_othermat` for details on the reuse of the matrix 231context. 232 233The directory `$PETSC_DIR/src/snes/tutorials` provides a variety of 234examples. 235 236Sometimes a nonlinear solver may produce a step that is not within the domain 237of a given function, for example one with a negative pressure. When this occurs 238one can call `SNESSetFunctionDomainError()` or `SNESSetJacobianDomainError()` 239to indicate to `SNES` the step is not valid. One must also use `SNESGetConvergedReason()` 240and check the reason to confirm if the solver succeeded. See {any}`sec_vi` for how to 241provide `SNES` with bounds on the variables to solve (differential) variational inequalities 242and how to control properties of the line step computed. 243 244(sec_nlsolvers)= 245 246## The Nonlinear Solvers 247 248As summarized in Table {any}`tab-snesdefaults`, `SNES` includes 249several Newton-like nonlinear solvers based on line search techniques 250and trust region methods. Also provided are several nonlinear Krylov 251methods, as well as nonlinear methods involving decompositions of the 252problem. 253 254Each solver may have associated with it a set of options, which can be 255set with routines and options database commands provided for this 256purpose. A complete list can be found by consulting the manual pages or 257by running a program with the `-help` option; we discuss just a few in 258the sections below. 259 260```{eval-rst} 261.. list-table:: PETSc Nonlinear Solvers 262 :name: tab-snesdefaults 263 :header-rows: 1 264 265 * - Method 266 - SNESType 267 - Options Name 268 - Default Line Search 269 * - Line Search Newton 270 - ``SNESNEWTONLS`` 271 - ``newtonls`` 272 - ``SNESLINESEARCHBT`` 273 * - Trust region Newton 274 - ``SNESNEWTONTR`` 275 - ``newtontr`` 276 - — 277 * - Newton with Arc Length Continuation 278 - ``SNESNEWTONAL`` 279 - ``newtonal`` 280 - — 281 * - Nonlinear Richardson 282 - ``SNESNRICHARDSON`` 283 - ``nrichardson`` 284 - ``SNESLINESEARCHL2`` 285 * - Nonlinear CG 286 - ``SNESNCG`` 287 - ``ncg`` 288 - ``SNESLINESEARCHCP`` 289 * - Nonlinear GMRES 290 - ``SNESNGMRES`` 291 - ``ngmres`` 292 - ``SNESLINESEARCHL2`` 293 * - Quasi-Newton 294 - ``SNESQN`` 295 - ``qn`` 296 - see :any:`tab-qndefaults` 297 * - Full Approximation Scheme 298 - ``SNESFAS`` 299 - ``fas`` 300 - — 301 * - Nonlinear ASM 302 - ``SNESNASM`` 303 - ``nasm`` 304 - – 305 * - ASPIN 306 - ``SNESASPIN`` 307 - ``aspin`` 308 - ``SNESLINESEARCHBT`` 309 * - Nonlinear Gauss-Seidel 310 - ``SNESNGS`` 311 - ``ngs`` 312 - – 313 * - Anderson Mixing 314 - ``SNESANDERSON`` 315 - ``anderson`` 316 - – 317 * - Newton with constraints (1) 318 - ``SNESVINEWTONRSLS`` 319 - ``vinewtonrsls`` 320 - ``SNESLINESEARCHBT`` 321 * - Newton with constraints (2) 322 - ``SNESVINEWTONSSLS`` 323 - ``vinewtonssls`` 324 - ``SNESLINESEARCHBT`` 325 * - Multi-stage Smoothers 326 - ``SNESMS`` 327 - ``ms`` 328 - – 329 * - Composite 330 - ``SNESCOMPOSITE`` 331 - ``composite`` 332 - – 333 * - Linear solve only 334 - ``SNESKSPONLY`` 335 - ``ksponly`` 336 - – 337 * - Python Shell 338 - ``SNESPYTHON`` 339 - ``python`` 340 - – 341 * - Shell (user-defined) 342 - ``SNESSHELL`` 343 - ``shell`` 344 - – 345 346``` 347 348### Line Search Newton 349 350The method `SNESNEWTONLS` (`-snes_type newtonls`) provides a 351line search Newton method for solving systems of nonlinear equations. By 352default, this technique employs cubic backtracking 353{cite}`dennis:83`. Alternative line search techniques are 354listed in Table {any}`tab-linesearches`. 355 356```{eval-rst} 357.. table:: PETSc Line Search Methods 358 :name: tab-linesearches 359 360 ==================== =========================== ================ 361 **Line Search** **SNESLineSearchType** **Options Name** 362 ==================== =========================== ================ 363 Backtracking ``SNESLINESEARCHBT`` ``bt`` 364 (damped) step ``SNESLINESEARCHBASIC`` ``basic`` 365 identical to above ``SNESLINESEARCHNONE`` ``none`` 366 L2-norm Minimization ``SNESLINESEARCHL2`` ``l2`` 367 Critical point ``SNESLINESEARCHCP`` ``cp`` 368 Bisection ``SNESLINESEARCHBISECTION`` ``bisection`` 369 Shell ``SNESLINESEARCHSHELL`` ``shell`` 370 ==================== =========================== ================ 371``` 372 373Every `SNES` has a line search context of type `SNESLineSearch` that 374may be retrieved using 375 376``` 377SNESGetLineSearch(SNES snes,SNESLineSearch *ls);. 378``` 379 380There are several default options for the line searches. The order of 381polynomial approximation may be set with `-snes_linesearch_order` or 382 383``` 384SNESLineSearchSetOrder(SNESLineSearch ls, PetscInt order); 385``` 386 387for instance, 2 for quadratic or 3 for cubic. Sometimes, it may not be 388necessary to monitor the progress of the nonlinear iteration. In this 389case, `-snes_linesearch_norms` or 390 391``` 392SNESLineSearchSetComputeNorms(SNESLineSearch ls,PetscBool norms); 393``` 394 395may be used to turn off function, step, and solution norm computation at 396the end of the linesearch. 397 398The default line search for the line search Newton method, 399`SNESLINESEARCHBT` involves several parameters, which are set to 400defaults that are reasonable for many applications. The user can 401override the defaults by using the following options: 402 403- `-snes_linesearch_alpha <alpha>` 404- `-snes_linesearch_maxstep <max>` 405- `-snes_linesearch_minlambda <tol>` 406 407Besides the backtracking linesearch, there are `SNESLINESEARCHL2`, 408which uses a polynomial secant minimization of $||F(x)||_2$, and 409`SNESLINESEARCHCP`, which minimizes $F(x) \cdot Y$ where 410$Y$ is the search direction. These are both potentially iterative 411line searches, which may be used to find a better-fitted steplength in 412the case where a single secant search is not sufficient. The number of 413iterations may be set with `-snes_linesearch_max_it`. In addition, the 414convergence criteria of the iterative line searches may be set using 415function tolerances `-snes_linesearch_rtol` and 416`-snes_linesearch_atol`, and steplength tolerance 417`snes_linesearch_ltol`. 418 419For highly non-linear problems, the bisection line search `SNESLINESEARCHBISECTION` 420may prove useful due to its robustness. Similar to the critical point line search 421`SNESLINESEARCHCP`, it seeks to find the root of $F(x) \cdot Y$. 422While the latter does so through a secant method, the bisection line search 423does so by iteratively bisecting the step length interval. 424It works as follows (with $f(\lambda)=F(x-\lambda Y) \cdot Y / ||Y||$ for brevity): 425 4261. initialize: $j=1$, $\lambda_0 = \lambda_{\text{left}} = 0.0$, $\lambda_j = \lambda_{\text{right}} = \alpha$, compute $f(\lambda_0)$ and $f(\lambda_j)$ 427 4282. check whether there is a change of sign in the interval: $f(\lambda_{\text{left}}) f(\lambda_j) \leq 0$; if not accept the full step length $\lambda_1$ 429 4303. if there is a change of sign, enter iterative bisection procedure 431 432 1. check convergence/ exit criteria: 433 434 - absolute tolerance $f(\lambda_j) < \mathtt{atol}$ 435 - relative tolerance $f(\lambda_j) < \mathtt{rtol} \cdot f(\lambda_0)$ 436 - change of step length $\lambda_j - \lambda_{j-1} < \mathtt{ltol}$ 437 - number of iterations $j < \mathtt{max\_it}$ 438 439 2. if $j > 1$, determine direction of bisection 440 441 $$ 442 \begin{aligned}\lambda_{\text{left}} &= \begin{cases}\lambda_{\text{left}} &f(\lambda_{\text{left}}) f(\lambda_j) \leq 0\\\lambda_{j} &\text{else}\\ \end{cases}\\ \lambda_{\text{right}} &= \begin{cases} \lambda_j &f(\lambda_{\text{left}}) f(\lambda_j) \leq 0\\\lambda_{\text{right}} &\text{else}\\ \end{cases}\\\end{aligned} 443 $$ 444 445 3. bisect the interval: $\lambda_{j+1} = (\lambda_{\text{left}} + \lambda_{\text{right}})/2$, compute $f(\lambda_{j+1})$ 446 4. update variables for the next iteration: $\lambda_j \gets \lambda_{j+1}$, $f(\lambda_j) \gets f(\lambda_{j+1})$, $j \gets j+1$ 447 448Custom line search types may either be defined using 449`SNESLineSearchShell`, or by creating a custom user line search type 450in the model of the preexisting ones and register it using 451 452``` 453SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch));. 454``` 455 456### Trust Region Methods 457 458The trust region method in `SNES` for solving systems of nonlinear 459equations, `SNESNEWTONTR` (`-snes_type newtontr`), is similar to the one developed in the 460MINPACK project {cite}`more84`. Several parameters can be 461set to control the variation of the trust region size during the 462solution process. In particular, the user can control the initial trust 463region radius, computed by 464 465$$ 466\Delta = \Delta_0 \| F_0 \|_2, 467$$ 468 469by setting $\Delta_0$ via the option `-snes_tr_delta0 <delta0>`. 470 471### Newton with Arc Length Continuation 472 473The Newton method with arc length continuation reformulates the linearized system 474$K\delta \mathbf x = -\mathbf F(\mathbf x)$ by introducing the load parameter 475$\lambda$ and splitting the residual into two components, commonly 476corresponding to internal and external forces: 477 478$$ 479\mathbf F(x, \lambda) = \mathbf F^{\mathrm{int}}(\mathbf x) - \mathbf F^{\mathrm{ext}}(\mathbf x, \lambda) 480$$ 481 482Often, $\mathbf F^{\mathrm{ext}}(\mathbf x, \lambda)$ is linear in $\lambda$, 483which can be thought of as applying the external force in proportional load 484increments. By default, this is how the right-hand side vector is handled in the 485implemented method. Generally, however, $\mathbf F^{\mathrm{ext}}(\mathbf x, \lambda)$ 486may depend non-linearly on $\lambda$ or $\mathbf x$, or both. 487To accommodate this possibility, we provide the `SNESNewtonALGetLoadParameter()` 488function, which allows for the current value of $\lambda$ to be queried in the 489functions provided to `SNESSetFunction()` and `SNESSetJacobian()`. 490 491Additionally, we split the solution update into two components: 492 493$$ 494\delta \mathbf x = \delta s\delta\mathbf x^F + \delta\lambda\delta\mathbf x^Q, 495$$ 496 497where $\delta s = 1$ unless partial corrections are used (discussed more 498below). Each of $\delta \mathbf x^F$ and $\delta \mathbf x^Q$ are found via 499solving a linear system with the Jacobian $K$: 500 501- $\delta \mathbf x^F$ is the full Newton step for a given value of $\lambda$: $K \delta \mathbf x^F = -\mathbf F(\mathbf x, \lambda)$ 502- $\delta \mathbf x^Q$ is the variation in $\mathbf x$ with respect to $\lambda$, computed by $K \delta\mathbf x^Q = \mathbf Q(\mathbf x, \lambda)$, where $\mathbf Q(\mathbf x, \lambda) = -\partial \mathbf F (\mathbf x, \lambda) / \partial \lambda$ is the tangent load vector. 503 504Often, the tangent load vector $\mathbf Q$ is constant within a load increment, 505which corresponds to the case of proportional loading discussed above. By default, 506$\mathbf Q$ is the full right-hand-side vector, if one was provided. 507The user can also provide a function which computes $\mathbf Q$ to 508`SNESNewtonALSetFunction()`. This function should have the same signature as for 509`SNESSetFunction`, and the user should use `SNESNewtonALGetLoadParameter()` to get 510$\lambda$ if it is needed. 511 512**The Constraint Surface.** Considering the $n+1$ dimensional space of 513$\mathbf x$ and $\lambda$, we define the linearized equilibrium line to be 514the set of points for which the linearized equilibrium equations are satisfied. 515Given the previous iterative solution 516$\mathbf t^{(j-1)} = [\mathbf x^{(j-1)}, \lambda^{(j-1)}]$, 517this line is defined by the point $\mathbf t^{(j-1)} + [\delta\mathbf x^F, 0]$ and 518the vector $\mathbf t^Q [\delta\mathbf x^Q, 1]$. 519The arc length method seeks the intersection of this linearized equilibrium line 520with a quadratic constraint surface, defined by 521 522% math::L^2 = \|\Delta x\|^2 + \psi^2 (\Delta\lambda)^2, 523 524where $L$ is a user-provided step size corresponding to the radius of the 525constraint surface, $\Delta\mathbf x$ and $\Delta\lambda$ are the 526accumulated updates over the current load step, and $\psi^2$ is a 527user-provided consistency parameter determining the shape of the constraint surface. 528Generally, $\psi^2 > 0$ leads to a hyper-sphere constraint surface, while 529$\psi^2 = 0$ leads to a hyper-cylinder constraint surface. 530 531Since the solution will always fall on the constraint surface, the method will often 532require multiple incremental steps to fully solve the non-linear problem. 533This is necessary to accurately trace the equilibrium path. 534Importantly, this is fundamentally different from time stepping. 535While a similar process could be implemented as a `TS`, this method is 536particularly designed to be used as a SNES, either standalone or within a `TS`. 537 538To this end, by default, the load parameter is used such that the full external 539forces are applied at $\lambda = 1$, although we allow for the user to specify 540a different value via `-snes_newtonal_lambda_max`. 541To ensure that the solution corresponds exactly to the external force prescribed by 542the user, i.e. that the load parameter is exactly $\lambda_{max}$ at the end 543of the SNES solve, we clamp the value before computing the solution update. 544As such, the final increment will likely be a hybrid of arc length continuation and 545normal Newton iterations. 546 547**Choosing the Continuation Step.** For the first iteration from an equilibrium 548point, there is a single correct way to choose $\delta\lambda$, which follows 549from the constraint equations. Specifically the constraint equations yield the 550quadratic equation $a\delta\lambda^2 + b\delta\lambda + c = 0$, where 551 552$$ 553\begin{aligned} 554a &= \|\delta\mathbf x^Q\|^2 + \psi^2,\\ 555b &= 2\delta\mathbf x^Q\cdot (\Delta\mathbf x + \delta s\delta\mathbf x^F) + 2\psi^2 \Delta\lambda,\\ 556c &= \|\Delta\mathbf x + \delta s\delta\mathbf x^F\|^2 + \psi^2 \Delta\lambda^2 - L^2. 557\end{aligned} 558$$ 559 560Since in the first iteration, $\Delta\mathbf x = \delta\mathbf x^F = \mathbf 0$ and 561$\Delta\lambda = 0$, $b = 0$ and the equation simplifies to a pair of 562real roots: 563 564$$ 565\delta\lambda = \pm\frac{L}{\sqrt{\|\delta\mathbf x^Q\|^2 + \psi^2}}, 566$$ 567 568where the sign is positive for the first increment and is determined by the previous 569increment otherwise as 570 571$$ 572\text{sign}(\delta\lambda) = \text{sign}\big(\delta\mathbf x^Q \cdot (\Delta\mathbf x)_{i-1} + \psi^2(\Delta\lambda)_{i-1}\big), 573$$ 574 575where $(\Delta\mathbf x)_{i-1}$ and $(\Delta\lambda)_{i-1}$ are the 576accumulated updates over the previous load step. 577 578In subsequent iterations, there are different approaches to selecting 579$\delta\lambda$, all of which have trade-offs. 580The main difference is whether the iterative solution falls on the constraint 581surface at every iteration, or only when fully converged. 582This MR implements one of each of these approaches, set via 583`SNESNewtonALSetCorrectionType()` or 584`-snes_newtonal_correction_type <normal|exact>` on the command line. 585 586**Corrections in the Normal Hyperplane.** The `SNES_NEWTONAL_CORRECTION_NORMAL` 587option is simpler and computationally less expensive, but may fail to converge, as 588the constraint equation is not satisfied at every iteration. 589The update $\delta \lambda$ is chosen such that the update is within the 590normal hyper-surface to the quadratic constraint surface. 591Mathematically, that is 592 593$$ 594\delta \lambda = -\frac{\Delta \mathbf x \cdot \delta \mathbf x^F}{\Delta\mathbf x \cdot \delta\mathbf x^Q + \psi^2 \Delta\lambda}. 595$$ 596 597This implementation is based on {cite}`LeonPaulinoPereiraMenezesLages_2011`. 598 599**Exact Corrections.** The `SNES_NEWTONAL_CORRECTION_EXACT` option is far more 600complex, but ensures that the constraint is exactly satisfied at every Newton 601iteration. As such, it is generally more robust. 602By evaluating the intersection of constraint surface and equilibrium line at each 603iteration, $\delta\lambda$ is chosen as one of the roots of the above 604quadratic equation $a\delta\lambda^2 + b\delta\lambda + c = 0$. 605This method encounters issues, however, if the linearized equilibrium line and 606constraint surface do not intersect due to particularly large linearized error. 607In this case, the roots are complex. 608To continue progressing toward a solution, this method uses a partial correction by 609choosing $\delta s$ such that the quadratic equation has a single real root. 610Geometrically, this is selecting the point on the constraint surface closest to the 611linearized equilibrium line. See the code or {cite}`Ritto-CorreaCamotim2008` for a 612mathematical description of these partial corrections. 613 614### Nonlinear Krylov Methods 615 616A number of nonlinear Krylov methods are provided, including Nonlinear 617Richardson (`SNESNRICHARDSON`), nonlinear conjugate gradient (`SNESNCG`), nonlinear GMRES (`SNESNGMRES`), and Anderson Mixing (`SNESANDERSON`). These 618methods are described individually below. They are all instrumental to 619PETSc’s nonlinear preconditioning. 620 621**Nonlinear Richardson.** The nonlinear Richardson iteration, `SNESNRICHARDSON`, merely 622takes the form of a line search-damped fixed-point iteration of the form 623 624$$ 625\mathbf{x}_{k+1} = \mathbf{x}_k - \lambda \mathbf{F}(\mathbf{x}_k), \;\; k=0,1, \ldots, 626$$ 627 628where the default linesearch is `SNESLINESEARCHL2`. This simple solver 629is mostly useful as a nonlinear smoother, or to provide line search 630stabilization to an inner method. 631 632**Nonlinear Conjugate Gradients.** Nonlinear CG, `SNESNCG`, is equivalent to linear 633CG, but with the steplength determined by line search 634(`SNESLINESEARCHCP` by default). Five variants (Fletcher-Reed, 635Hestenes-Steifel, Polak-Ribiere-Polyak, Dai-Yuan, and Conjugate Descent) 636are implemented in PETSc and may be chosen using 637 638``` 639SNESNCGSetType(SNES snes, SNESNCGType btype); 640``` 641 642**Anderson Mixing and Nonlinear GMRES Methods.** Nonlinear GMRES (`SNESNGMRES`), and 643Anderson Mixing (`SNESANDERSON`) methods combine the last $m$ iterates, plus a new 644fixed-point iteration iterate, into an approximate residual-minimizing new iterate. 645 646All of the above methods have support for using a nonlinear preconditioner to compute the preliminary update step, rather than the default 647which is the nonlinear function's residual, \$ mathbf\{F}(mathbf\{x}\_k)\$. The different update is obtained by solving a nonlinear preconditioner nonlinear problem, which has its own 648`SNES` object that may be obtained with `SNESGetNPC()`. 649Quasi-Newton Methods 650^^^^^^^^^^^^^^^^^^^^ 651 652Quasi-Newton methods store iterative rank-one updates to the Jacobian 653instead of computing the Jacobian directly. Three limited-memory quasi-Newton 654methods are provided, L-BFGS, which are described in 655Table {any}`tab-qndefaults`. These all are encapsulated under 656`-snes_type qn` and may be changed with `snes_qn_type`. The default 657is L-BFGS, which provides symmetric updates to an approximate Jacobian. 658This iteration is similar to the line search Newton methods. 659 660The quasi-Newton methods support the use of a nonlinear preconditioner that can be obtained with `SNESGetNPC()` and then configured; or that can be configured with 661`SNES`, `KSP`, and `PC` options using the options database prefix `-npc_`. 662 663```{eval-rst} 664.. list-table:: PETSc quasi-Newton solvers 665 :name: tab-qndefaults 666 :header-rows: 1 667 668 * - QN Method 669 - ``SNESQNType`` 670 - Options Name 671 - Default Line Search 672 * - L-BFGS 673 - ``SNES_QN_LBFGS`` 674 - ``lbfgs`` 675 - ``SNESLINESEARCHCP`` 676 * - “Good” Broyden 677 - ``SNES_QN_BROYDEN`` 678 - ``broyden`` 679 - ``SNESLINESEARCHBASIC`` (or equivalently ``SNESLINESEARCHNONE`` 680 * - “Bad” Broyden 681 - ``SNES_QN_BADBROYDEN`` 682 - ``badbroyden`` 683 - ``SNESLINESEARCHL2`` 684``` 685 686One may also control the form of the initial Jacobian approximation with 687 688``` 689SNESQNSetScaleType(SNES snes, SNESQNScaleType stype); 690``` 691 692and the restart type with 693 694``` 695SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype); 696``` 697 698### The Full Approximation Scheme 699 700The Nonlinear Full Approximation Scheme (FAS) `SNESFAS`, is a nonlinear multigrid method. At 701each level, there is a recursive cycle control `SNES` instance, and 702either one or two nonlinear solvers that act as smoothers (up and down). Problems 703set up using the `SNES` `DMDA` interface are automatically 704coarsened. FAS, `SNESFAS`, differs slightly from linear multigrid `PCMG`, in that the hierarchy is 705constructed recursively. However, much of the interface is a one-to-one 706map. We describe the “get” operations here, and it can be assumed that 707each has a corresponding “set” operation. For instance, the number of 708levels in the hierarchy may be retrieved using 709 710``` 711SNESFASGetLevels(SNES snes, PetscInt *levels); 712``` 713 714There are four `SNESFAS` cycle types, `SNES_FAS_MULTIPLICATIVE`, 715`SNES_FAS_ADDITIVE`, `SNES_FAS_FULL`, and `SNES_FAS_KASKADE`. The 716type may be set with 717 718``` 719SNESFASSetType(SNES snes,SNESFASType fastype);. 720``` 721 722and the cycle type, 1 for V, 2 for W, may be set with 723 724``` 725SNESFASSetCycles(SNES snes, PetscInt cycles);. 726``` 727 728Much like the interface to `PCMG` described in {any}`sec_mg`, there are interfaces to recover the 729various levels’ cycles and smoothers. The level smoothers may be 730accessed with 731 732``` 733SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth); 734SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth); 735SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth); 736``` 737 738and the level cycles with 739 740``` 741SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes);. 742``` 743 744Also akin to `PCMG`, the restriction and prolongation at a level may 745be acquired with 746 747``` 748SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat); 749SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat); 750``` 751 752In addition, FAS requires special restriction for solution-like 753variables, called injection. This may be set with 754 755``` 756SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat);. 757``` 758 759The coarse solve context may be acquired with 760 761``` 762SNESFASGetCoarseSolve(SNES snes, SNES *smooth); 763``` 764 765### Nonlinear Additive Schwarz 766 767Nonlinear Additive Schwarz methods (NASM) take a number of local 768nonlinear subproblems, solves them independently in parallel, and 769combines those solutions into a new approximate solution. 770 771``` 772SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]); 773``` 774 775allows for the user to create these local subdomains. Problems set up 776using the `SNES` `DMDA` interface are automatically decomposed. To 777begin, the type of subdomain updates to the whole solution are limited 778to two types borrowed from `PCASM`: `PC_ASM_BASIC`, in which the 779overlapping updates added. `PC_ASM_RESTRICT` updates in a 780nonoverlapping fashion. This may be set with 781 782``` 783SNESNASMSetType(SNES snes,PCASMType type);. 784``` 785 786`SNESASPIN` is a helper `SNES` type that sets up a nonlinearly 787preconditioned Newton’s method using NASM as the preconditioner. 788 789## General Options 790 791This section discusses options and routines that apply to all `SNES` 792solvers and problem classes. In particular, we focus on convergence 793tests, monitoring routines, and tools for checking derivative 794computations. 795 796(sec_snesconvergence)= 797 798### Convergence Tests 799 800Convergence of the nonlinear solvers can be detected in a variety of 801ways; the user can even specify a customized test, as discussed below. 802Most of the nonlinear solvers use `SNESConvergenceTestDefault()`, 803however, `SNESNEWTONTR` uses a method-specific additional convergence 804test as well. The convergence tests involves several parameters, which 805are set by default to values that should be reasonable for a wide range 806of problems. The user can customize the parameters to the problem at 807hand by using some of the following routines and options. 808 809One method of convergence testing is to declare convergence when the 810norm of the change in the solution between successive iterations is less 811than some tolerance, `stol`. Convergence can also be determined based 812on the norm of the function. Such a test can use either the absolute 813size of the norm, `atol`, or its relative decrease, `rtol`, from an 814initial guess. The following routine sets these parameters, which are 815used in many of the default `SNES` convergence tests: 816 817``` 818SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol, PetscInt its,PetscInt fcts); 819``` 820 821This routine also sets the maximum numbers of allowable nonlinear 822iterations, `its`, and function evaluations, `fcts`. The 823corresponding options database commands for setting these parameters are: 824 825- `-snes_atol <atol>` 826- `-snes_rtol <rtol>` 827- `-snes_stol <stol>` 828- `-snes_max_it <its>` 829- `-snes_max_funcs <fcts>` (use `unlimited` for no maximum) 830 831A related routine is `SNESGetTolerances()`. `PETSC_CURRENT` may be used 832for any parameter to indicate the current value should be retained; use `PETSC_DETERMINE` to restore to the default value from when the object was created. 833 834Users can set their own customized convergence tests in `SNES` by 835using the command 836 837``` 838SNESSetConvergenceTest(SNES snes,PetscErrorCode (*test)(SNES snes,PetscInt it,PetscReal xnorm, PetscReal gnorm,PetscReal f,SNESConvergedReason reason, void *cctx),void *cctx,PetscErrorCode (*destroy)(void *cctx)); 839``` 840 841The final argument of the convergence test routine, `cctx`, denotes an 842optional user-defined context for private data. When solving systems of 843nonlinear equations, the arguments `xnorm`, `gnorm`, and `f` are 844the current iterate norm, current step norm, and function norm, 845respectively. `SNESConvergedReason` should be set positive for 846convergence and negative for divergence. See `include/petscsnes.h` for 847a list of values for `SNESConvergedReason`. 848 849(sec_snesmonitor)= 850 851### Convergence Monitoring 852 853By default the `SNES` solvers run silently without displaying 854information about the iterations. The user can initiate monitoring with 855the command 856 857``` 858SNESMonitorSet(SNES snes, PetscErrorCode (*mon)(SNES snes, PetscInt its, PetscReal norm, void* mctx), void *mctx, (PetscCtxDestroyFn *)*monitordestroy); 859``` 860 861The routine, `mon`, indicates a user-defined monitoring routine, where 862`its` and `mctx` respectively denote the iteration number and an 863optional user-defined context for private data for the monitor routine. 864The argument `norm` is the function norm. 865 866The routine set by `SNESMonitorSet()` is called once after every 867successful step computation within the nonlinear solver. Hence, the user 868can employ this routine for any application-specific computations that 869should be done after the solution update. The option `-snes_monitor` 870activates the default `SNES` monitor routine, 871`SNESMonitorDefault()`, while `-snes_monitor_lg_residualnorm` draws 872a simple line graph of the residual norm’s convergence. 873 874One can cancel hardwired monitoring routines for `SNES` at runtime 875with `-snes_monitor_cancel`. 876 877As the Newton method converges so that the residual norm is small, say 878$10^{-10}$, many of the final digits printed with the 879`-snes_monitor` option are meaningless. Worse, they are different on 880different machines; due to different round-off rules used by, say, the 881IBM RS6000 and the Sun SPARC. This makes testing between different 882machines difficult. The option `-snes_monitor_short` causes PETSc to 883print fewer of the digits of the residual norm as it gets smaller; thus 884on most of the machines it will always print the same numbers making 885cross-process testing easier. 886 887The routines 888 889``` 890SNESGetSolution(SNES snes,Vec *x); 891SNESGetFunction(SNES snes,Vec *r,void *ctx,int(**func)(SNES,Vec,Vec,void*)); 892``` 893 894return the solution vector and function vector from a `SNES` context. 895These routines are useful, for instance, if the convergence test 896requires some property of the solution or function other than those 897passed with routine arguments. 898 899(sec_snesderivs)= 900 901### Checking Accuracy of Derivatives 902 903Since hand-coding routines for Jacobian matrix evaluation can be error 904prone, `SNES` provides easy-to-use support for checking these matrices 905against finite difference versions. In the simplest form of comparison, 906users can employ the option `-snes_test_jacobian` to compare the 907matrices at several points. Although not exhaustive, this test will 908generally catch obvious problems. One can compare the elements of the 909two matrices by using the option `-snes_test_jacobian_view` , which 910causes the two matrices to be printed to the screen. 911 912Another means for verifying the correctness of a code for Jacobian 913computation is running the problem with either the finite difference or 914matrix-free variant, `-snes_fd` or `-snes_mf`; see {any}`sec_fdmatrix` or {any}`sec_nlmatrixfree`. 915If a 916problem converges well with these matrix approximations but not with a 917user-provided routine, the problem probably lies with the hand-coded 918matrix. See the note in {any}`sec_snesjacobian` about 919assembling your Jabobian in the "preconditioner" slot `Pmat`. 920 921The correctness of user provided `MATSHELL` Jacobians in general can be 922checked with `MatShellTestMultTranspose()` and `MatShellTestMult()`. 923 924The correctness of user provided `MATSHELL` Jacobians via `TSSetRHSJacobian()` 925can be checked with `TSRHSJacobianTestTranspose()` and `TSRHSJacobianTest()` 926that check the correction of the matrix-transpose vector product and the 927matrix-product. From the command line, these can be checked with 928 929- `-ts_rhs_jacobian_test_mult_transpose` 930- `-mat_shell_test_mult_transpose_view` 931- `-ts_rhs_jacobian_test_mult` 932- `-mat_shell_test_mult_view` 933 934## Inexact Newton-like Methods 935 936Since exact solution of the linear Newton systems within {math:numref}`newton` 937at each iteration can be costly, modifications 938are often introduced that significantly reduce these expenses and yet 939retain the rapid convergence of Newton’s method. Inexact or truncated 940Newton techniques approximately solve the linear systems using an 941iterative scheme. In comparison with using direct methods for solving 942the Newton systems, iterative methods have the virtue of requiring 943little space for matrix storage and potentially saving significant 944computational work. Within the class of inexact Newton methods, of 945particular interest are Newton-Krylov methods, where the subsidiary 946iterative technique for solving the Newton system is chosen from the 947class of Krylov subspace projection methods. Note that at runtime the 948user can set any of the linear solver options discussed in {any}`ch_ksp`, 949such as `-ksp_type <ksp_method>` and 950`-pc_type <pc_method>`, to set the Krylov subspace and preconditioner 951methods. 952 953Two levels of iterations occur for the inexact techniques, where during 954each global or outer Newton iteration a sequence of subsidiary inner 955iterations of a linear solver is performed. Appropriate control of the 956accuracy to which the subsidiary iterative method solves the Newton 957system at each global iteration is critical, since these inner 958iterations determine the asymptotic convergence rate for inexact Newton 959techniques. While the Newton systems must be solved well enough to 960retain fast local convergence of the Newton’s iterates, use of excessive 961inner iterations, particularly when $\| \mathbf{x}_k - \mathbf{x}_* \|$ is large, 962is neither necessary nor economical. Thus, the number of required inner 963iterations typically increases as the Newton process progresses, so that 964the truncated iterates approach the true Newton iterates. 965 966A sequence of nonnegative numbers $\{\eta_k\}$ can be used to 967indicate the variable convergence criterion. In this case, when solving 968a system of nonlinear equations, the update step of the Newton process 969remains unchanged, and direct solution of the linear system is replaced 970by iteration on the system until the residuals 971 972$$ 973\mathbf{r}_k^{(i)} = \mathbf{F}'(\mathbf{x}_k) \Delta \mathbf{x}_k + \mathbf{F}(\mathbf{x}_k) 974$$ 975 976satisfy 977 978$$ 979\frac{ \| \mathbf{r}_k^{(i)} \| }{ \| \mathbf{F}(\mathbf{x}_k) \| } \leq \eta_k \leq \eta < 1. 980$$ 981 982Here $\mathbf{x}_0$ is an initial approximation of the solution, and 983$\| \cdot \|$ denotes an arbitrary norm in $\Re^n$ . 984 985By default a constant relative convergence tolerance is used for solving 986the subsidiary linear systems within the Newton-like methods of 987`SNES`. When solving a system of nonlinear equations, one can instead 988employ the techniques of Eisenstat and Walker {cite}`ew96` 989to compute $\eta_k$ at each step of the nonlinear solver by using 990the option `-snes_ksp_ew` . In addition, by adding one’s own 991`KSP` convergence test (see {any}`sec_convergencetests`), one can easily create one’s own, 992problem-dependent, inner convergence tests. 993 994(sec_nlmatrixfree)= 995 996## Matrix-Free Methods 997 998The `SNES` class fully supports matrix-free methods. The matrices 999specified in the Jacobian evaluation routine need not be conventional 1000matrices; instead, they can point to the data required to implement a 1001particular matrix-free method. The matrix-free variant is allowed *only* 1002when the linear systems are solved by an iterative method in combination 1003with no preconditioning (`PCNONE` or `-pc_type` `none`), a 1004user-provided preconditioner matrix, or a user-provided preconditioner 1005shell (`PCSHELL`, discussed in {any}`sec_pc`); that 1006is, obviously matrix-free methods cannot be used with a direct solver, 1007approximate factorization, or other preconditioner which requires access 1008to explicit matrix entries. 1009 1010The user can create a matrix-free context for use within `SNES` with 1011the routine 1012 1013``` 1014MatCreateSNESMF(SNES snes,Mat *mat); 1015``` 1016 1017This routine creates the data structures needed for the matrix-vector 1018products that arise within Krylov space iterative 1019methods {cite}`brownsaad:90`. 1020The default `SNES` 1021matrix-free approximations can also be invoked with the command 1022`-snes_mf`. Or, one can retain the user-provided Jacobian 1023preconditioner, but replace the user-provided Jacobian matrix with the 1024default matrix-free variant with the option `-snes_mf_operator`. 1025 1026`MatCreateSNESMF()` uses 1027 1028``` 1029MatCreateMFFD(Vec x, Mat *mat); 1030``` 1031 1032which can also be used directly for users who need a matrix-free matrix but are not using `SNES`. 1033 1034The user can set one parameter to control the Jacobian-vector product 1035approximation with the command 1036 1037``` 1038MatMFFDSetFunctionError(Mat mat,PetscReal rerror); 1039``` 1040 1041The parameter `rerror` should be set to the square root of the 1042relative error in the function evaluations, $e_{rel}$; the default 1043is the square root of machine epsilon (about $10^{-8}$ in double 1044precision), which assumes that the functions are evaluated to full 1045floating-point precision accuracy. This parameter can also be set from 1046the options database with `-mat_mffd_err <err>` 1047 1048In addition, PETSc provides ways to register new routines to compute 1049the differencing parameter ($h$); see the manual page for 1050`MatMFFDSetType()` and `MatMFFDRegister()`. We currently provide two 1051default routines accessible via `-mat_mffd_type <ds or wp>`. For 1052the default approach there is one “tuning” parameter, set with 1053 1054``` 1055MatMFFDDSSetUmin(Mat mat,PetscReal umin); 1056``` 1057 1058This parameter, `umin` (or $u_{min}$), is a bit involved; its 1059default is $10^{-6}$ . Its command line form is `-mat_mffd_umin <umin>`. 1060 1061The Jacobian-vector product is approximated 1062via the formula 1063 1064$$ 1065F'(u) a \approx \frac{F(u + h*a) - F(u)}{h} 1066$$ 1067 1068where $h$ is computed via 1069 1070$$ 1071h = e_{\text{rel}} \cdot \begin{cases} 1072u^{T}a/\lVert a \rVert^2_2 & \text{if $|u^T a| > u_{\min} \lVert a \rVert_{1}$} \\ 1073u_{\min} \operatorname{sign}(u^{T}a) \lVert a \rVert_{1}/\lVert a\rVert^2_2 & \text{otherwise}. 1074\end{cases} 1075$$ 1076 1077This approach is taken from Brown and Saad 1078{cite}`brownsaad:90`. The second approach, taken from Walker and Pernice, 1079{cite}`pw98`, computes $h$ via 1080 1081$$ 1082\begin{aligned} 1083 h = \frac{\sqrt{1 + ||u||}e_{rel}}{||a||}\end{aligned} 1084$$ 1085 1086This has no tunable parameters, but note that inside the nonlinear solve 1087for the entire *linear* iterative process $u$ does not change 1088hence $\sqrt{1 + ||u||}$ need be computed only once. This 1089information may be set with the options 1090 1091``` 1092MatMFFDWPSetComputeNormU(Mat mat,PetscBool ); 1093``` 1094 1095or `-mat_mffd_compute_normu <true or false>`. This information is used 1096to eliminate the redundant computation of these parameters, therefore 1097reducing the number of collective operations and improving the 1098efficiency of the application code. This takes place automatically for the PETSc GMRES solver with left preconditioning. 1099 1100It is also possible to monitor the differencing parameters h that are 1101computed via the routines 1102 1103``` 1104MatMFFDSetHHistory(Mat,PetscScalar *,int); 1105MatMFFDResetHHistory(Mat,PetscScalar *,int); 1106MatMFFDGetH(Mat,PetscScalar *); 1107``` 1108 1109We include an explicit example of using matrix-free methods in {any}`ex3.c <snes_ex3>`. 1110Note that by using the option `-snes_mf` one can 1111easily convert any `SNES` code to use a matrix-free Newton-Krylov 1112method without a preconditioner. As shown in this example, 1113`SNESSetFromOptions()` must be called *after* `SNESSetJacobian()` to 1114enable runtime switching between the user-specified Jacobian and the 1115default `SNES` matrix-free form. 1116 1117(snes_ex3)= 1118 1119:::{admonition} Listing: `src/snes/tutorials/ex3.c` 1120```{literalinclude} /../src/snes/tutorials/ex3.c 1121:end-before: /*TEST 1122``` 1123::: 1124 1125Table {any}`tab-jacobians` summarizes the various matrix situations 1126that `SNES` supports. In particular, different linear system matrices 1127and preconditioning matrices are allowed, as well as both matrix-free 1128and application-provided preconditioners. If {any}`ex3.c <snes_ex3>` is run with 1129the options `-snes_mf` and `-user_precond` then it uses a 1130matrix-free application of the matrix-vector multiple and a user 1131provided matrix-free Jacobian. 1132 1133```{eval-rst} 1134.. list-table:: Jacobian Options 1135 :name: tab-jacobians 1136 1137 * - Matrix Use 1138 - Conventional Matrix Formats 1139 - Matrix-free versions 1140 * - Jacobian Matrix 1141 - Create matrix with ``MatCreate()``:math:`^*`. Assemble matrix with user-defined routine :math:`^\dagger` 1142 - Create matrix with ``MatCreateShell()``. Use ``MatShellSetOperation()`` to set various matrix actions, or use ``MatCreateMFFD()`` or ``MatCreateSNESMF()``. 1143 * - Preconditioning Matrix 1144 - Create matrix with ``MatCreate()``:math:`^*`. Assemble matrix with user-defined routine :math:`^\dagger` 1145 - Use ``SNESGetKSP()`` and ``KSPGetPC()`` to access the ``PC``, then use ``PCSetType(pc, PCSHELL)`` followed by ``PCShellSetApply()``. 1146``` 1147 1148$^*$ Use either the generic `MatCreate()` or a format-specific variant such as `MatCreateAIJ()`. 1149 1150$^\dagger$ Set user-defined matrix formation routine with `SNESSetJacobian()` or with a `DM` variant such as `DMDASNESSetJacobianLocal()` 1151 1152SNES also provides some less well-integrated code to apply matrix-free finite differencing using an automatically computed measurement of the 1153noise of the functions. This can be selected with `-snes_mf_version 2`; it does not use `MatCreateMFFD()` but has similar options that start with 1154`-snes_mf_` instead of `-mat_mffd_`. Note that this alternative prefix **only** works for version 2 differencing. 1155 1156(sec_fdmatrix)= 1157 1158## Finite Difference Jacobian Approximations 1159 1160PETSc provides some tools to help approximate the Jacobian matrices 1161efficiently via finite differences. These tools are intended for use in 1162certain situations where one is unable to compute Jacobian matrices 1163analytically, and matrix-free methods do not work well without a 1164preconditioner, due to very poor conditioning. The approximation 1165requires several steps: 1166 1167- First, one colors the columns of the (not yet built) Jacobian matrix, 1168 so that columns of the same color do not share any common rows. 1169- Next, one creates a `MatFDColoring` data structure that will be 1170 used later in actually computing the Jacobian. 1171- Finally, one tells the nonlinear solvers of `SNES` to use the 1172 `SNESComputeJacobianDefaultColor()` routine to compute the 1173 Jacobians. 1174 1175A code fragment that demonstrates this process is given below. 1176 1177``` 1178ISColoring iscoloring; 1179MatFDColoring fdcoloring; 1180MatColoring coloring; 1181 1182/* 1183 This initializes the nonzero structure of the Jacobian. This is artificial 1184 because clearly if we had a routine to compute the Jacobian we wouldn't 1185 need to use finite differences. 1186*/ 1187FormJacobian(snes,x, &J, &J, &user); 1188 1189/* 1190 Color the matrix, i.e. determine groups of columns that share no common 1191 rows. These columns in the Jacobian can all be computed simultaneously. 1192*/ 1193MatColoringCreate(J, &coloring); 1194MatColoringSetType(coloring,MATCOLORINGSL); 1195MatColoringSetFromOptions(coloring); 1196MatColoringApply(coloring, &iscoloring); 1197MatColoringDestroy(&coloring); 1198/* 1199 Create the data structure that SNESComputeJacobianDefaultColor() uses 1200 to compute the actual Jacobians via finite differences. 1201*/ 1202MatFDColoringCreate(J,iscoloring, &fdcoloring); 1203ISColoringDestroy(&iscoloring); 1204MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))FormFunction, &user); 1205MatFDColoringSetFromOptions(fdcoloring); 1206 1207/* 1208 Tell SNES to use the routine SNESComputeJacobianDefaultColor() 1209 to compute Jacobians. 1210*/ 1211SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,fdcoloring); 1212``` 1213 1214Of course, we are cheating a bit. If we do not have an analytic formula 1215for computing the Jacobian, then how do we know what its nonzero 1216structure is so that it may be colored? Determining the structure is 1217problem dependent, but fortunately, for most structured grid problems 1218(the class of problems for which PETSc was originally designed) if one 1219knows the stencil used for the nonlinear function one can usually fairly 1220easily obtain an estimate of the location of nonzeros in the matrix. 1221This is harder in the unstructured case, but one typically knows where the nonzero entries are from the mesh topology and distribution of degrees of freedom. 1222If using `DMPlex` ({any}`ch_unstructured`) for unstructured meshes, the nonzero locations will be identified in `DMCreateMatrix()` and the procedure above can be used. 1223Most external packages for unstructured meshes have similar functionality. 1224 1225One need not necessarily use a `MatColoring` object to determine a 1226coloring. For example, if a grid can be colored directly (without using 1227the associated matrix), then that coloring can be provided to 1228`MatFDColoringCreate()`. Note that the user must always preset the 1229nonzero structure in the matrix regardless of which coloring routine is 1230used. 1231 1232PETSc provides the following coloring algorithms, which can be selected using `MatColoringSetType()` or via the command line argument `-mat_coloring_type`. 1233 1234```{eval-rst} 1235.. list-table:: 1236 :header-rows: 1 1237 1238 * - Algorithm 1239 - ``MatColoringType`` 1240 - ``-mat_coloring_type`` 1241 - Parallel 1242 * - smallest-last :cite:`more84` 1243 - ``MATCOLORINGSL`` 1244 - ``sl`` 1245 - No 1246 * - largest-first :cite:`more84` 1247 - ``MATCOLORINGLF`` 1248 - ``lf`` 1249 - No 1250 * - incidence-degree :cite:`more84` 1251 - ``MATCOLORINGID`` 1252 - ``id`` 1253 - No 1254 * - Jones-Plassmann :cite:`jp:pcolor` 1255 - ``MATCOLORINGJP`` 1256 - ``jp`` 1257 - Yes 1258 * - Greedy 1259 - ``MATCOLORINGGREEDY`` 1260 - ``greedy`` 1261 - Yes 1262 * - Natural (1 color per column) 1263 - ``MATCOLORINGNATURAL`` 1264 - ``natural`` 1265 - Yes 1266 * - Power (:math:`A^k` followed by 1-coloring) 1267 - ``MATCOLORINGPOWER`` 1268 - ``power`` 1269 - Yes 1270``` 1271 1272As for the matrix-free computation of Jacobians ({any}`sec_nlmatrixfree`), two parameters affect the accuracy of the 1273finite difference Jacobian approximation. These are set with the command 1274 1275``` 1276MatFDColoringSetParameters(MatFDColoring fdcoloring,PetscReal rerror,PetscReal umin); 1277``` 1278 1279The parameter `rerror` is the square root of the relative error in the 1280function evaluations, $e_{rel}$; the default is the square root of 1281machine epsilon (about $10^{-8}$ in double precision), which 1282assumes that the functions are evaluated approximately to floating-point 1283precision accuracy. The second parameter, `umin`, is a bit more 1284involved; its default is $10e^{-6}$ . Column $i$ of the 1285Jacobian matrix (denoted by $F_{:i}$) is approximated by the 1286formula 1287 1288$$ 1289F'_{:i} \approx \frac{F(u + h*dx_{i}) - F(u)}{h} 1290$$ 1291 1292where $h$ is computed via: 1293 1294$$ 1295h = e_{\text{rel}} \cdot \begin{cases} 1296u_{i} & \text{if $|u_{i}| > u_{\min}$} \\ 1297u_{\min} \cdot \operatorname{sign}(u_{i}) & \text{otherwise}. 1298\end{cases} 1299$$ 1300 1301for `MATMFFD_DS` or: 1302 1303$$ 1304h = e_{\text{rel}} \sqrt(\|u\|) 1305$$ 1306 1307for `MATMFFD_WP` (default). These parameters may be set from the options 1308database with 1309 1310``` 1311-mat_fd_coloring_err <err> 1312-mat_fd_coloring_umin <umin> 1313-mat_fd_type <htype> 1314``` 1315 1316Note that `MatColoring` type `MATCOLORINGSL`, `MATCOLORINGLF`, and 1317`MATCOLORINGID` are sequential algorithms. `MATCOLORINGJP` and 1318`MATCOLORINGGREEDY` are parallel algorithms, although in practice they 1319may create more colors than the sequential algorithms. If one computes 1320the coloring `iscoloring` reasonably with a parallel algorithm or by 1321knowledge of the discretization, the routine `MatFDColoringCreate()` 1322is scalable. An example of this for 2D distributed arrays is given below 1323that uses the utility routine `DMCreateColoring()`. 1324 1325``` 1326DMCreateColoring(da,IS_COLORING_GHOSTED, &iscoloring); 1327MatFDColoringCreate(J,iscoloring, &fdcoloring); 1328MatFDColoringSetFromOptions(fdcoloring); 1329ISColoringDestroy( &iscoloring); 1330``` 1331 1332Note that the routine `MatFDColoringCreate()` currently is only 1333supported for the AIJ and BAIJ matrix formats. 1334 1335(sec_vi)= 1336 1337## Variational Inequalities 1338 1339`SNES` can also solve (differential) variational inequalities with box (bound) constraints. 1340These are nonlinear algebraic systems with additional inequality 1341constraints on some or all of the variables: 1342$L_i \le u_i \le H_i$. For example, the pressure variable cannot be negative. 1343Some, or all, of the lower bounds may be 1344negative infinity (indicated to PETSc with `SNES_VI_NINF`) and some, or 1345all, of the upper bounds may be infinity (indicated by `SNES_VI_INF`). 1346The commands 1347 1348``` 1349SNESVISetVariableBounds(SNES,Vec L,Vec H); 1350SNESVISetComputeVariableBounds(SNES snes, PetscErrorCode (*compute)(SNES,Vec,Vec)) 1351``` 1352 1353are used to indicate that one is solving a variational inequality. Problems with box constraints can be solved with 1354the reduced space, `SNESVINEWTONRSLS`, and semi-smooth `SNESVINEWTONSSLS` solvers. 1355 1356The 1357option `-snes_vi_monitor` turns on extra monitoring of the active set 1358associated with the bounds and `-snes_vi_type` allows selecting from 1359several VI solvers, the default is preferred. 1360 1361`SNESLineSearchSetPreCheck()` and `SNESLineSearchSetPostCheck()` can also be used to control properties 1362of the steps selected by `SNES`. 1363 1364(sec_snespc)= 1365 1366## Nonlinear Preconditioning 1367 1368The mathematical framework of nonlinear preconditioning is explained in detail in {cite}`bruneknepleysmithtu15`. 1369Nonlinear preconditioning in PETSc involves the use of an inner `SNES` 1370instance to define the step for an outer `SNES` instance. The inner 1371instance may be extracted using 1372 1373``` 1374SNESGetNPC(SNES snes,SNES *npc); 1375``` 1376 1377and passed run-time options using the `-npc_` prefix. Nonlinear 1378preconditioning comes in two flavors: left and right. The side may be 1379changed using `-snes_npc_side` or `SNESSetNPCSide()`. Left nonlinear 1380preconditioning redefines the nonlinear function as the action of the 1381nonlinear preconditioner $\mathbf{M}$; 1382 1383$$ 1384\mathbf{F}_{M}(x) = \mathbf{M}(\mathbf{x},\mathbf{b}) - \mathbf{x}. 1385$$ 1386 1387Right nonlinear preconditioning redefines the nonlinear function as the 1388function on the action of the nonlinear preconditioner; 1389 1390$$ 1391\mathbf{F}(\mathbf{M}(\mathbf{x},\mathbf{b})) = \mathbf{b}, 1392$$ 1393 1394which can be interpreted as putting the preconditioner into “striking 1395distance” of the solution by outer acceleration. 1396 1397In addition, basic patterns of solver composition are available with the 1398`SNESType` `SNESCOMPOSITE`. This allows for two or more `SNES` 1399instances to be combined additively or multiplicatively. By command 1400line, a set of `SNES` types may be given by comma separated list 1401argument to `-snes_composite_sneses`. There are additive 1402(`SNES_COMPOSITE_ADDITIVE`), additive with optimal damping 1403(`SNES_COMPOSITE_ADDITIVEOPTIMAL`), and multiplicative 1404(`SNES_COMPOSITE_MULTIPLICATIVE`) variants which may be set with 1405 1406``` 1407SNESCompositeSetType(SNES,SNESCompositeType); 1408``` 1409 1410New subsolvers may be added to the composite solver with 1411 1412``` 1413SNESCompositeAddSNES(SNES,SNESType); 1414``` 1415 1416and accessed with 1417 1418``` 1419SNESCompositeGetSNES(SNES,PetscInt,SNES *); 1420``` 1421 1422```{eval-rst} 1423.. bibliography:: /petsc.bib 1424 :filter: docname in docnames 1425``` 1426