petsc-3.6.1 2015-08-06
Report Typos and Errors

DMPlexCreateFromCellList

This takes as input common mesh generator output, a list of the vertices for each cell, and produces a DM

Synopsis

#include "petscdmplex.h"   
PetscErrorCode DMPlexCreateFromCellList(MPI_Comm comm, PetscInt dim, PetscInt numCells, PetscInt numVertices, PetscInt numCorners, PetscBool interpolate, const int cells[], PetscInt spaceDim, const double vertexCoords[], DM *dm)

Input Parameters

comm - The communicator
dim - The topological dimension of the mesh
numCells - The number of cells
numVertices - The number of vertices
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 vertices 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

Note: Two triangles sharing a face


       2
     / | \
    /  |  \
   /   |   \
  0  0 | 1  3
   \   |   /
    \  |  /
     \ | /
       1
would have input
 numCells = 2, numVertices = 4
 cells = [0 1 2  1 3 2]

which would result in the DMPlex

       4
     / | \
    /  |  \
   /   |   \
  2  0 | 1  5
   \   |   /
    \  |  /
     \ | /
       3

See Also

DMPlexCreateFromDAG(), DMPlexCreate()

Level:beginner
Location:
src/dm/impls/plex/plexcreate.c
Index of all DM routines
Table of Contents for all manual pages
Index of all manual pages