:orphan: # DMFieldEvaluate Evaluate the field and its derivatives on a set of points ## Synopsis ``` #include "petscdmfield.h" #include "petscdmfield.h" PetscErrorCode DMFieldEvaluate(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H) ``` Collective ## Input Parameters - ***field -*** The `DMField` object - ***points -*** The points at which to evaluate the field. Should have size d x n, where d is the coordinate dimension of the manifold and n is the number of points - ***datatype -*** The PetscDataType of the output arrays: either `PETSC_REAL` or `PETSC_SCALAR`. If the field is complex and datatype is `PETSC_REAL`, the real part of the field is returned. ## Output Parameters - ***B -*** pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field. If B is not NULL, the values of the field are written in this array, varying first by component, then by point. - ***D -*** pointer to data of size d * c * n * sizeof(datatype). If D is not NULL, the values of the field's spatial derivatives are written in this array, varying first by the partial derivative component, then by field component, then by point. - ***H -*** pointer to data of size d * d * c * n * sizeof(datatype). If H is not NULL, the values of the field's second spatial derivatives are written in this array, varying first by the second partial derivative component, then by field component, then by point. ## See Also `DMField`, `DMFieldGetDM()`, `DMFieldGetNumComponents()`, `DMFieldEvaluateFE()`, `DMFieldEvaluateFV()`, `PetscDataType` ## Level intermediate ## Location src/dm/field/interface/dmfield.c ## Examples src/dm/field/tutorials/ex1.c
## Implementations DMFieldEvaluate_DA in src/dm/field/impls/da/dmfieldda.c
DMFieldEvaluate_DS in src/dm/field/impls/ds/dmfieldds.c
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/dm/field/interface/dmfield.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)