TSGetIJacobianP#
Gets the function that computes the Jacobian of \( F\) w.r.t. the parameters \(p\) where \(F(Udot,U,p,t) = G(U,p,t) \), as well as the location to store the matrix.
Synopsis#
#include <petscts.h>
PetscErrorCode TSGetIJacobianP(TS ts, Mat *Amat, PetscErrorCode (**func)(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal shift, Mat A, void *ctx), void **ctx)
Logically Collective
Input Parameter#
Output Parameters#
Amat - JacobianP matrix
func - the function that computes the JacobianP
ctx - [optional] user-defined function context
Calling sequence of func
#
ts - the
TS
contextt - current timestep
U - input vector (current ODE solution)
Udot - time derivative of state vector
shift - shift to apply, see the note in
TSSetIJacobian()
A - output matrix
ctx - [optional] user-defined function context
Note#
Amat
has the same number of rows and the same row parallel layout as u
, Amat
has the same number of columns and parallel layout as p
See Also#
TS: Scalable ODE and DAE Solvers, TSSetRHSJacobianP()
, TS
, TSSetIJacobianP()
, TSGetRHSJacobianP()
Level#
intermediate
Location#
src/ts/interface/sensitivity/tssen.c
Index of all Sensitivity routines
Table of Contents for all manual pages
Index of all manual pages