DMPlexCreateDoublet#
Creates a mesh of two cells of the specified type, optionally with later refinement.
Synopsis#
#include "petscdmplex.h"
#include "petscdmplex.h"
PetscErrorCode DMPlexCreateDoublet(MPI_Comm comm, PetscInt dim, PetscBool simplex, PetscBool interpolate, PetscReal refinementLimit, DM *newdm)
Collective
Input Parameters#
comm - The communicator for the
DMobjectdim - The spatial dimension
simplex - Flag for simplicial cells, otherwise they are tensor product cells
interpolate - Flag to create intermediate mesh pieces (edges, faces)
refinementLimit - A nonzero number indicates the largest admissible volume for a refined cell
Output Parameter#
newdm - The
DMobject
See Also#
DMPlex: Unstructured Grids, DM, DMPLEX, DMSetType(), DMCreate()
Level#
beginner
Location#
src/dm/impls/plex/plexcreate.c
Index of all DMPlex routines
Table of Contents for all manual pages
Index of all manual pages