DMPlexCreateFromCellSectionParallel#
Create distributed DMPLEX
from a list of vertices for each cell (common mesh generator output) and supports multiple celltypes
Synopsis#
#include "petscdmplex.h"
#include "petscdmplextransform.h"
PetscErrorCode DMPlexCreateFromCellSectionParallel(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt NVertices, PetscSection cellSection, PetscBool interpolate, const PetscInt cells[], PetscInt spaceDim, const PetscReal vertexCoords[], PetscSF *vertexSF, PetscInt **verticesAdj, DM *dm)
Collective
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
cellSection - The
PetscSection
giving the number of vertices for each cell (layout of cells)interpolate - Flag indicating that intermediate mesh entities (faces, edges) should be created automatically
cells - An array of 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 Parameters#
Notes#
This function is just a convenient sequence of DMCreate()
, DMSetType()
, DMSetDimension()
,
DMPlexBuildFromCellSectionParallel()
, DMPlexInterpolate()
, DMPlexBuildCoordinatesFromCellListParallel()
See DMPlexBuildFromCellSectionParallel()
for an example and details about the topology-related parameters.
See DMPlexBuildCoordinatesFromCellListParallel()
for details about the geometry-related parameters.
See Also#
DMPlex: Unstructured Grids, DM
, DMPLEX
, 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