petsc-3.3-p7 2013-05-11

SNESSetPicard

Use SNES to solve the semilinear-system A(x) x = b(x) via a Picard type iteration (Picard linearization)

Synopsis

#include "petscsnes.h"  
#include "petscdmshell.h" 
#include "petscsys.h" 
PetscErrorCode  SNESSetPicard(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),Mat jmat, Mat mat, PetscErrorCode (*mfunc)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
Logically Collective on SNES

Input Parameters

snes - the SNES context
r - vector to store function value
func - function evaluation routine
jmat - normally the same as mat but you can pass another matrix for which you compute the Jacobian of A(x) x - b(x) (see jmat below)
mat - matrix to store A
mfunc - function to compute matrix value
ctx - [optional] user-defined context for private data for the function evaluation routine (may be PETSC_NULL)

Calling sequence of func

   func (SNES snes,Vec x,Vec f,void *ctx);

f - function vector
ctx - optional user-defined function context

Calling sequence of mfunc

    mfunc (SNES snes,Vec x,Mat *jmat,Mat *mat,int *flag,void *ctx);

x - input vector
jmat - Form Jacobian matrix of A(x) x - b(x) if available, not there is really no reason to use it in this way since then you can just use SNESSetJacobian(), normally just pass mat in this location
mat - form A(x) matrix
flag - flag indicating information about the preconditioner matrix structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER
ctx - [optional] user-defined Jacobian context

Notes

One can call SNESSetPicard() or SNESSetFunction() (and possibly SNESSetJacobian()) but cannot call both

    Solves the equation A(x) x = b(x) via the defect correction algorithm A(x^{n}) (x^{n+1} - x^{n}) = b(x^{n}) - A(x^{n})x^{n}
    Note that when an exact solver is used this corresponds to the "classic" Picard A(x^{n}) x^{n+1} = b(x^{n}) iteration.

Run with -snes_mf_operator to solve the system with Newton's method using A(x^{n}) to construct the preconditioner.

We implement the defect correction form of the Picard iteration because it converges much more generally when inexact linear solvers are used then the direct Picard iteration A(x^n) x^{n+1} = b(x^n)

There is some controversity over the definition of a Picard iteration for nonlinear systems but almost everyone agrees that it involves a linear solve and some believe it is the iteration A(x^{n}) x^{n+1} = b(x^{n}) hence we use the name Picard. If anyone has an authoritative reference that defines the Picard iteration different please contact us at [email protected] and we'll have an entirely new argument :-).

Keywords

SNES, nonlinear, set, function

See Also

SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESGetPicard(), SNESLineSearchPreCheckPicard()

Level:beginner
Location:
src/snes/interface/snes.c
Index of all SNES routines
Table of Contents for all manual pages
Index of all manual pages

Examples

src/snes/examples/tutorials/ex15.c.html