:orphan: # DMCreateInjection Gets injection matrix between two `DM` objects. This is an operator that applied to a vector obtained with `DMCreateGlobalVector()` on the fine grid maps the values to a vector on the vector on the coarse `DM` by simply selecting the values on the coarse grid points. This compares to the operator obtained by `DMCreateRestriction()` or the transpose of the operator obtained by `DMCreateInterpolation()` that uses a "local weighted average" of the values around the coarse grid point as the coarse grid value. ## Synopsis ``` #include "petscdm.h" #include "petscdmlabel.h" #include "petscds.h" PetscErrorCode DMCreateInjection(DM dac, DM daf, Mat *mat) ``` Collective ## Input Parameters - ***dac -*** the `DM` object - ***daf -*** the second, finer `DM` object ## Output Parameter - ***mat -*** the injection ## Note For `DMDA` objects this only works for "uniform refinement", that is the refined mesh was obtained `DMRefine()` or the coarse mesh was obtained by `DMCoarsen()`. The coordinates set into the `DMDA` are completely ignored in computing the injection. ## See Also [](ch_dmbase), `DM`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateColoring()`, `DMCreateMatrix()`, `DMCreateMassMatrix()`, `DMCreateInterpolation()`, `DMCreateRestriction()`, `MatRestrict()`, `MatInterpolate()` ## Level developer ## Location src/dm/interface/dm.c ## Implementations DMCreateInjection_DA in src/dm/impls/da/dainterp.c
DMCreateInjection_pforest in src/dm/impls/forest/p4est/pforest.h
DMCreateInjection_Moab in src/dm/impls/moab/dmmbmg.cxx
DMCreateInjection_Plex in src/dm/impls/plex/plex.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)