xref: /petsc/src/snes/mf/snesmfj.c (revision 0decc0a32d11c9048faef123c7738214ee57a3d4)
1 #define PETSCSNES_DLL
2 
3 #include "include/private/snesimpl.h"  /*I  "petscsnes.h" I*/
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "MatMFFDComputeJacobian"
7 /*@
8    MatMFFDComputeJacobian - Tells the matrix-free Jacobian object the new location at which
9        Jacobian matrix vector products will be computed at, i.e. J(x) * a.
10 
11    Collective on SNES
12 
13    Input Parameters:
14 +   snes - the nonlinear solver context
15 .   x - the point at which the Jacobian vector products will be performed
16 .   jac - the matrix-free Jacobian object
17 .   B - either the same as jac or another matrix type (ignored)
18 .   flag - not relevent for matrix-free form
19 -   dummy - the user context (ignored)
20 
21    Level: developer
22 
23    Notes:
24      This can be passed into SNESSetJacobian() when using a completely matrix-free solver,
25      that is the B matrix is also the same matrix operator. This is used when you select
26      -mat_mffd but rarely used directly by users.
27 
28 .seealso: MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD,
29           MatMFFDSetHHistory(),
30           MatMFFDKSPMonitor(), MatMFFDSetFunctionError(), MatMFFDCreate(), SNESSetJacobian()
31 
32 @*/
33 PetscErrorCode PETSCSNES_DLLEXPORT MatMFFDComputeJacobian(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure *flag,void *dummy)
34 {
35   PetscErrorCode ierr;
36   PetscFunctionBegin;
37   ierr = MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
38   ierr = MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
39   PetscFunctionReturn(0);
40 }
41 
42 #undef __FUNCT__
43 #define __FUNCT__ "MatCreateSNESMF"
44 /*@
45    MatCreateSNESMF - Creates a matrix-free matrix context for use with
46    a SNES solver.  This matrix can be used as the Jacobian argument for
47    the routine SNESSetJacobian().
48 
49    Collective on SNES and Vec
50 
51    Input Parameters:
52 +  snes - the SNES context
53 -  x - vector where SNES solution is to be stored.
54 
55    Output Parameter:
56 .  J - the matrix-free matrix
57 
58    Level: advanced
59 
60    Notes:
61    The matrix-free matrix context merely contains the function pointers
62    and work space for performing finite difference approximations of
63    Jacobian-vector products, F'(u)*a,
64 
65    The default code uses the following approach to compute h
66 
67 .vb
68      F'(u)*a = [F(u+h*a) - F(u)]/h where
69      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
70        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
71  where
72      error_rel = square root of relative error in function evaluation
73      umin = minimum iterate parameter
74 .ve
75    (see MATSNESMF_WP or MATSNESMF_DS)
76 
77    The user can set the error_rel via MatSNESMFSetFunctionError() and
78    umin via MatSNESMFDefaultSetUmin(); see the nonlinear solvers chapter
79    of the users manual for details.
80 
81    The user should call MatDestroy() when finished with the matrix-free
82    matrix context.
83 
84    Options Database Keys:
85 +  -mat_mffd_err <error_rel> - Sets error_rel
86 +  -mat_mffd_type - wp or ds (see MATSNESMF_WP or MATSNESMF_DS)
87 .  -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only)
88 -  -mat_mffd_ksp_monitor - KSP monitor routine that prints differencing h
89 
90 .keywords: SNES, default, matrix-free, create, matrix
91 
92 .seealso: MatDestroy(), MatSNESMFSetFunctionError(), MatSNESMFDefaultSetUmin()
93           MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), MatCreateMF(),
94           MatSNESMFGetH(),MatSNESMFKSPMonitor(), MatSNESMFRegisterDynamic), MatSNESMFComputeJacobian()
95 
96 @*/
97 PetscErrorCode PETSCSNES_DLLEXPORT MatCreateSNESMF(SNES snes,Vec x,Mat *J)
98 {
99   PetscErrorCode ierr;
100 
101   PetscFunctionBegin;
102   ierr = MatCreateMFFD(x,J);CHKERRQ(ierr);
103   ierr = MatMFFDSetFunction(*J,snes->vec_sol,(PetscErrorCode (*)(void*, _p_Vec*, _p_Vec*))SNESComputeFunction,snes);CHKERRQ(ierr);
104   PetscFunctionReturn(0);
105 }
106 
107 
108 
109 
110 
111 
112