DMDASNESSetJacobianLocalVec#
set a local Jacobian evaluation function that operates on a local vector with DMDA
Synopsis#
#include "petscdmda.h" 
#include "petscsnes.h" 
PetscErrorCode DMDASNESSetJacobianLocalVec(DM dm, PetscErrorCode (*func)(DMDALocalInfo *info, Vec x, Mat J, Mat M, void *), void *ctx)
Logically Collective
Input Parameters#
- dm - - DMto associate callback with
- func - local Jacobian evaluation 
- ctx - optional context for local Jacobian evaluation 
Calling sequence of func#
- info - - DMDALocalInfodefining the subdomain to evaluate the Jacobian at
- x - state vector at which to evaluate Jacobian 
- J - the Jacobian 
- M - approximate Jacobian from which the preconditioner will be computed, often - J
- ctx - optional context passed above 
See Also#
SNES: Nonlinear Solvers, DMDA, DMDASNESSetJacobianLocal(), DMDASNESSetFunctionLocalVec(), DMSNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
Level#
beginner
Location#
Examples#
Index of all SNES routines
Table of Contents for all manual pages
Index of all manual pages