#include "petscdmda.h" #include "petscts.h" PetscErrorCode DMDATSSetRHSJacobianLocal(DM dm,DMDATSRHSJacobianLocal func,void *ctx)Logically Collective
dm | - DM to associate callback with | |
func | - local RHS Jacobian evaluation routine | |
ctx | - optional context for local jacobian evaluation |
func(DMDALocalInfo* info,PetscReal t,void* x,Mat J,Mat B,MatStructure *flg,void *ctx);
info | - DMDALocalInfo defining the subdomain to evaluate the residual on | |
t | - time at which to evaluate residual | |
x | - array of local state information | |
J | - Jacobian matrix | |
B | - preconditioner matrix; often same as J | |
flg | - flag indicating information about the preconditioner matrix structure (same as flag in KSPSetOperators()) | |
ctx | - optional context passed above |
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