petsc-3.14.6 2021-03-30
Report Typos and Errors

DMPlexCreateFromCellListParallelPetsc

Create distributed DMPLEX from a list of vertices for each cell (common mesh generator output)

Synopsis

#include "petscdmplex.h"   
#include "petscdmplex.h"   
PetscErrorCode DMPlexCreateFromCellListParallelPetsc(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt NVertices, PetscInt numCorners, PetscBool interpolate, const PetscInt cells[], PetscInt spaceDim, const PetscReal vertexCoords[], PetscSF *vertexSF, DM *dm)

Input Parameters

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

Output Parameter

dm - The DM
vertexSF - (Optional) SF describing complete vertex ownership

Notes

This function is just a convenient sequence of DMCreate(), DMSetType(), DMSetDimension(), DMPlexBuildFromCellListParallel(), DMPlexInterpolate(), DMPlexBuildCoordinatesFromCellListParallel()

See DMPlexBuildFromCellListParallel() for an example and details about the topology-related parameters. See DMPlexBuildCoordinatesFromCellListParallel() for details about the geometry-related parameters.

See Also

DMPlexCreateFromCellListPetsc(), DMPlexBuildFromCellListParallel(), DMPlexBuildCoordinatesFromCellListParallel(), DMPlexCreateFromDAG(), DMPlexCreate()

Level

intermediate

Location

src/dm/impls/plex/plexcreate.c
Index of all DMPLEX routines
Table of Contents for all manual pages
Index of all manual pages