DMPlexBuildFromCellSectionParallel#
Build distributed DMPLEX
topology from a list of vertices for each cell (common mesh generator output) allowing multiple celltypes
Synopsis#
#include "petscdmplex.h"
#include "petscdmplextransform.h"
PetscErrorCode DMPlexBuildFromCellSectionParallel(DM dm, PetscInt numCells, PetscInt numVertices, PetscInt NVertices, PetscSection cellSection, const PetscInt cells[], PetscSF *vertexSF, PetscInt **verticesAdjSaved)
Collective; No Fortran Support
Input Parameters#
dm - The
DM
numCells - The number of cells owned by this process
numVertices - The number of vertices to be owned by this process, or
PETSC_DECIDE
NVertices - The global number of vertices, or
PETSC_DETERMINE
cellSection - The
PetscSection
giving the number of vertices for each cell (layout of cells)cells - An array of the global vertex numbers for each cell
Output Parameters#
vertexSF - (Optional)
PetscSF
describing complete vertex ownershipverticesAdjSaved - (Optional) vertex adjacency array
Notes#
A triangle and quadrilateral sharing a face
2----------3
/ | |
/ | |
/ | |
0 0 | 1 |
\ | |
\ | |
\ | |
1----------4
would have input
numCells = 2, numVertices = 5
cells = [0 1 2 1 4 3 2]
which would result in the DMPLEX
4----------5
/ | |
/ | |
/ | |
2 0 | 1 |
\ | |
\ | |
\ | |
3----------6
Vertices are implicitly numbered consecutively 0,…,NVertices.
Each rank owns a chunk of numVertices consecutive vertices.
If numVertices is PETSC_DECIDE
, PETSc will distribute them as evenly as possible using PetscLayout.
If NVertices is PETSC_DETERMINE
and numVertices is PETSC_DECIDE, NVertices is computed by PETSc as the maximum vertex index in cells + 1.
If only NVertices is PETSC_DETERMINE
, it is computed as the sum of numVertices over all ranks.
The cell distribution is arbitrary non-overlapping, independent of the vertex distribution.
See Also#
DMPlex: Unstructured Grids, DM
, DMPLEX
, DMPlexBuildFromCellListParallel()
, DMPlexCreateFromCellSectionParallel()
, DMPlexBuildCoordinatesFromCellListParallel()
,
PetscSF
Level#
advanced
Location#
src/dm/impls/plex/plexcreate.c
Index of all DMPlex routines
Table of Contents for all manual pages
Index of all manual pages