# TaoGetObjectiveAndGradient
Gets the combined objective function and gradient evaluation routine for the function to be optimized 
## Synopsis
```
#include "petsctao.h" 
PetscErrorCode TaoGetObjectiveAndGradient(Tao tao, Vec *g, PetscErrorCode (**func)(Tao, Vec, PetscReal *, Vec, void *), void **ctx)
```
Not collective


## Input Parameter

- ***tao -*** the Tao context



## Output Parameters

- ***g -*** the vector to internally hold the gradient computation
- ***func -*** the gradient function
- ***ctx -*** user-defined context for private data for the gradient evaluation routine



## Calling sequence of func
```none
func (Tao tao, Vec x, PetscReal *f, Vec g, void *ctx);
```


- ***x -*** input vector
- ***f -*** objective value (output)
- ***g -*** gradient value (output)
- ***ctx -*** [optional] user-defined function context





## See Also
 `Tao`, `TaoSolve()`, `TaoSetObjective()`, `TaoSetGradient()`, `TaoSetHessian()`, `TaoSetObjectiveAndGradient()`

## Level
beginner

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/tao/interface/taosolver_fg.c.html#TaoGetObjectiveAndGradient">src/tao/interface/taosolver_fg.c</A>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/tao/interface/taosolver_fg.c)


[Index of all Tao routines](index.md)  
[Table of Contents for all manual pages](/docs/manualpages/index.md)  
[Index of all manual pages](/docs/manualpages/singleindex.md)  
