:orphan:
# DMPlexCreateFromCellListPetsc
Create `DMPLEX` from a list of vertices for each cell (common mesh generator output), but only process 0 takes in the input 
## Synopsis
```
#include "petscdmplex.h"   
#include "petscdmplex.h"   
PetscErrorCode DMPlexCreateFromCellListPetsc(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const PetscInt cells[], PetscInt spaceDim, const PetscReal vertexCoords[], DM *dm)
```
Collective


## Input Parameters

- ***comm -*** The communicator
- ***dim -*** The topological dimension of the mesh
- ***numCells -*** The number of cells, only on process 0
- ***numVertices -*** The number of vertices owned by this process, or `PETSC_DECIDE`, only on process 0
- ***numCorners -*** The number of vertices for each cell, only on process 0
- ***interpolate -*** Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
- ***cells -*** An array of numCells*numCorners numbers, the vertices for each cell, only on process 0
- ***spaceDim -*** The spatial dimension used for coordinates
- ***vertexCoords -*** An array of numVertices*spaceDim numbers, the coordinates of each vertex, only on process 0



## Output Parameter

- ***dm -*** The `DM`, which only has points on process 0





## Notes
This function is just a convenient sequence of `DMCreate()`, `DMSetType()`, `DMSetDimension()`, `DMPlexBuildFromCellList()`,
`DMPlexInterpolate()`, `DMPlexBuildCoordinatesFromCellList()`

See `DMPlexBuildFromCellList()` for an example and details about the topology-related parameters.
See `DMPlexBuildCoordinatesFromCellList()` for details about the geometry-related parameters.
See `DMPlexCreateFromCellListParallelPetsc()` for parallel input


## See Also
 [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateFromCellListParallelPetsc()`, `DMPlexBuildFromCellList()`, `DMPlexBuildCoordinatesFromCellList()`, `DMPlexCreateFromDAG()`, `DMPlexCreate()`

## Level
intermediate

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/dm/impls/plex/plexcreate.c.html#DMPlexCreateFromCellListPetsc">src/dm/impls/plex/plexcreate.c</A>

## Examples
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/dm/tutorials/ex20.c.html">src/dm/tutorials/ex20.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/tao/tutorials/ex3.c.html">src/tao/tutorials/ex3.c</A><BR>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/dm/impls/plex/plexcreate.c)


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