MatCreateMFFD#

Creates a matrix-free matrix of type MATMFFD that uses finite differences on a provided function to approximately multiply a vector by the matrix (Jacobian) . See also MatCreateSNESMF()

Synopsis#

#include "petscmat.h"   
PetscErrorCode MatCreateMFFD(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, Mat *J)

Collective

Input Parameters#

  • comm - MPI communicator

  • m - number of local rows (or PETSC_DECIDE to have calculated if M is given) This value should be the same as the local size used in creating the y vector for the matrix-vector product y = Ax.

  • n - This value should be the same as the local size used in creating the x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have calculated if N is given) For square matrices n is almost always m.

  • M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)

  • N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)

Output Parameter#

  • J - the matrix-free matrix

Options Database Keys#

  • -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS)

  • -mat_mffd_err - square root of estimated relative error in function evaluation

  • -mat_mffd_period - how often h is recomputed, defaults to 1, every time

  • -mat_mffd_check_positivity - possibly decrease h until U + h*a has only positive values

  • -mat_mffd_umin - Sets umin (for default PETSc routine that computes h only)

  • -mat_mffd_complex - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing (requires real valued functions but that PETSc be configured for complex numbers)

  • -snes_mf - use the finite difference based matrix-free matrix with SNESSolve() and no preconditioner

  • -snes_mf_operator - use the finite difference based matrix-free matrix with SNESSolve() but construct a preconditioner using the matrix passed as pmat to SNESSetJacobian().

Notes#

Use MatMFFDSetFunction() to provide the function that will be differenced to compute the matrix-vector product.

The matrix-free matrix context contains the function pointers and work space for performing finite difference approximations of Jacobian-vector products, F’(u)*a,

The default code uses the following approach to compute h

     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

To have SNES use the matrix-free finite difference matrix-vector product and not provide a separate matrix from which to compute the preconditioner (the pmat argument SNESSetJacobian()), then simply call SNESSetJacobian() with NULL for the matrices and MatMFFDComputeJacobian(). Or use the options database option -snes_mf

The user can set error_rel via MatMFFDSetFunctionError() and umin via MatMFFDDSSetUmin().

Use MATSHELL or MatCreateShell() to provide your own custom matrix-vector operation.

See Also#

Matrices, Mat, MATMFFD, MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), MatCreateShell(), MATSHELL, MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian()

Level#

advanced

Location#

src/mat/impls/mffd/mffd.c

Examples#

src/snes/tutorials/ex22.c

Implementations#

MatCreateMFFD_DS() in src/mat/impls/mffd/mffddef.c
MatCreateMFFD_WP() in src/mat/impls/mffd/wp.c


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