petsc-3.13.6 2020-09-29
Report Typos and Errors

MATMFFD_DS

the code for compute the "h" used in the finite difference matrix-free matrix vector product. This code implements the strategy in Dennis and Schnabel, "Numerical Methods for Unconstrained Optimization and Nonlinear Equations".

Options Database Keys

-mat_mffd_umin <umin> see MatMFFDDSSetUmin() -

Notes

Requires 2 norms and 1 inner product, but they are computed together so only one parallel collective operation is needed. See MATMFFD_WP for a method (with GMRES) that requires NO collective operations.

Formula used

F'(u)*a = [F(u+h*a) - F(u)]/h where h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise where error_rel = square root of relative error in function evaluation umin = minimum iterate parameter

See Also

MATMFFD, MatCreateMFFD(), MatCreateSNESMF(), MATMFFD_WP, MatMFFDDSSetUmin()

Level

intermediate

Location

src/mat/impls/mffd/mffddef.c

Examples

src/snes/tutorials/ex14.c.html

Index of all Mat routines
Table of Contents for all manual pages
Index of all manual pages