:orphan: # TaoSetJacobianStateRoutine Sets the function to compute the Jacobian (and its inverse) of the constraint function with respect to the state variables. Used only for PDE-constrained optimization. ## Synopsis ``` #include "petsctao.h" PetscErrorCode TaoSetJacobianStateRoutine(Tao tao, Mat J, Mat Jpre, Mat Jinv, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, Mat, void *), void *ctx) ``` Logically Collective ## Input Parameters - ***tao -*** the `Tao` context - ***J -*** Matrix used for the Jacobian - ***Jpre -*** Matrix that will be used to construct the preconditioner, can be same as `J`. Only used if `Jinv` is `NULL` - ***Jinv -*** [optional] Matrix used to apply the inverse of the state Jacobian. Use `NULL` to default to PETSc `KSP` solvers to apply the inverse. - ***func -*** Jacobian evaluation routine - ***ctx -*** [optional] user-defined context for private data for the Jacobian evaluation routine (may be `NULL`) ## Calling sequence of `func` ```none PetscErrorCode func(Tao tao, Vec x, Mat J, Mat Jpre, Mat Jinv, void *ctx); ``` - ***tao -*** the `Tao` context - ***x -*** input vector - ***J -*** Jacobian matrix - ***Jpre -*** matrix used to construct the preconditioner, usually the same as `J` - ***Jinv -*** inverse of `J` - ***ctx -*** [optional] user-defined Jacobian context ## See Also [](ch_tao), `Tao`, `TaoComputeJacobianState()`, `TaoSetJacobianDesignRoutine()`, `TaoSetStateDesignIS()` ## Level intermediate ## Location src/tao/interface/taosolver_hj.c ## Examples src/tao/pde_constrained/tutorials/elliptic.c
src/tao/pde_constrained/tutorials/hyperbolic.c
src/tao/pde_constrained/tutorials/parabolic.c
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/tao/interface/taosolver_hj.c) [Index of all Tao routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)