PETSc version 3.15.5
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 on Tao
Input Parameters
|  | tao | - the Tao context | 
|  | J | - Matrix used for the jacobian | 
|  | Jpre | - Matrix that will be used operated on by PETSc 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
   func(Tao tao,Vec x,Mat J,Mat Jpre,Mat Jinv,void *ctx);
|  | tao | - the Tao  context | 
|  | x | - input vector | 
|  | J | - Jacobian matrix | 
|  | Jpre | - preconditioner matrix, usually the same as J | 
|  | Jinv | - inverse of J | 
|  | ctx | - [optional] user-defined Jacobian context | 
See Also
 TaoComputeJacobianState(), TaoSetJacobianDesignRoutine(), TaoSetStateDesignIS()
Level
intermediate
Location
src/tao/interface/taosolver_hj.c
Examples
src/tao/pde_constrained/tutorials/elliptic.c.html
src/tao/pde_constrained/tutorials/parabolic.c.html
src/tao/pde_constrained/tutorials/hyperbolic.c.html
Index of all Tao routines
Table of Contents for all manual pages
Index of all manual pages