PETSc version 3.15.5
DMDATSSetIJacobianLocal
set a local residual evaluation function 
Synopsis
#include "petscdmda.h" 
#include "petscts.h" 
PetscErrorCode DMDATSSetIJacobianLocal(DM dm,DMDATSIJacobianLocal func,void *ctx)
Logically Collective
Input Arguments
|  | dm | - DM to associate callback with | 
|  | func | - local residual evaluation | 
|  | ctx | - optional context for local residual evaluation | 
Calling sequence for func
func(DMDALocalInfo* info,PetscReal t,void* x,void *xdot,PetscScalar shift,Mat J,Mat B,void *ctx);
|  | info | - DMDALocalInfo defining the subdomain to evaluate the residual on | 
|  | t | - time at which to evaluate the jacobian | 
|  | x | - array of local state information | 
|  | xdot | - time derivative at this state | 
|  | shift | - see TSSetIJacobian() for the meaning of this parameter | 
|  | J | - Jacobian matrix | 
|  | B | - preconditioner matrix; often same as J | 
|  | ctx | - optional context passed above | 
See Also
 DMTSSetJacobian(), DMDATSSetIFunctionLocal(), DMDASNESSetJacobianLocal()
Level
beginner
Location
src/ts/utils/dmdats.c
Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages