petsc-3.13.6 2020-09-29
Report Typos and Errors

TaoSetJacobianResidualRoutine

Sets the function to compute the least-squares residual Jacobian as well as the location to store the matrix.

Synopsis

#include "petsctao.h" 
PetscErrorCode TaoSetJacobianResidualRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, 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 preconditioner, can be same as J
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,void *ctx);

tao - the Tao context
x - input vector
J - Jacobian matrix
Jpre - preconditioning matrix, usually the same as J
ctx - [optional] user-defined Jacobian context

Level

intermediate

Location

src/tao/interface/taosolver_hj.c

Examples

src/tao/leastsquares/tutorials/chwirut1.c.html
src/tao/leastsquares/tutorials/cs1.c.html
src/tao/leastsquares/tutorials/tomography.c.html

Index of all Tao routines
Table of Contents for all manual pages
Index of all manual pages