:orphan:
# DMProjectFunction
This projects the given function into the function space provided by a `DM`, putting the coefficients in a global vector.
## Synopsis
```
#include "petscdm.h"
#include "petscdmlabel.h"
#include "petscds.h"
PetscErrorCode DMProjectFunction(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X)
```
Collective
## Input Parameters
- ***dm -*** The `DM`
- ***time -*** The time
- ***funcs -*** The coordinate functions to evaluate, one per field
- ***ctxs -*** Optional array of contexts to pass to each coordinate function. ctxs itself may be null.
- ***mode -*** The insertion mode for values
## Output Parameter
- ***X -*** vector
## Calling sequence of `funcs`
```none
PetscErrorCode funcs(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx);
```
- ***dim -*** The spatial dimension
- ***time -*** The time at which to sample
- ***x -*** The coordinates
- ***Nc -*** The number of components
- ***u -*** The output field values
- ***ctx -*** optional user-defined function context
## Developer Notes
This API is specific to only particular usage of `DM`
The notes need to provide some information about what has to be provided to the `DM` to be able to perform the computation.
## See Also
[](ch_dmbase), `DM`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabel()`, `DMComputeL2Diff()`
## Level
developer
## Location
src/dm/interface/dm.c
## Examples
src/dm/dt/dualspace/impls/lagrange/tutorials/ex2.c
src/snes/tutorials/ex12.c
src/snes/tutorials/ex62.c
src/snes/tutorials/ex63.c
src/snes/tutorials/ex69.c
src/snes/tutorials/ex71.c
src/snes/tutorials/ex75.c
src/snes/tutorials/ex76.c
src/snes/tutorials/ex77.c
src/tao/tutorials/ex1.c
src/tao/tutorials/ex2.c
---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/dm/interface/dm.c)
[Index of all DM routines](index.md)
[Table of Contents for all manual pages](/manualpages/index.md)
[Index of all manual pages](/manualpages/singleindex.md)