DMDASNESSetJacobianLocal#
set a local Jacobian evaluation function for use with DMDA
Synopsis#
#include "petscdmda.h"
#include "petscsnes.h"
PetscErrorCode DMDASNESSetJacobianLocal(DM dm, PetscErrorCode (*func)(DMDALocalInfo *info, void *x, Mat J, Mat M, void *ctx), void *ctx)
Logically Collective
Input Parameters#
dm -
DM
to associate callback withfunc - local Jacobian evaluation function
ctx - optional context for local Jacobian evaluation
Calling sequence of func
#
info -
DMDALocalInfo
defining the subdomain to evaluate the Jacobian atx - dimensional pointer to state at which to evaluate Jacobian (e.g. PetscScalar *x or **x or ***x)
J -
Mat
object for the JacobianM -
Mat
object for the Jacobian preconditioner matrix, oftenJ
ctx - optional context passed above
Note#
The J
and M
matrices are created internally by DMCreateMatrix()
See Also#
SNES: Nonlinear Solvers, DMDA
, DMDASNESSetFunctionLocal()
, DMSNESSetJacobian()
, DMDACreate1d()
, DMDACreate2d()
, DMDACreate3d()
Level#
beginner
Location#
Examples#
src/snes/tutorials/ex5.c
src/snes/tutorials/ex4.c
src/snes/tutorials/ex55.c
src/snes/tutorials/ex9.c
src/snes/tutorials/ex48.c
src/snes/tutorials/ex15.c
src/snes/tutorials/ex16.c
Index of all SNES routines
Table of Contents for all manual pages
Index of all manual pages