DMPlexCreateHypercubicMesh#
Creates a periodic mesh on the tensor product of unit intervals using only vertices and edges.
Synopsis#
#include "petscdmplex.h"
#include "petscdmplextransform.h"
PetscErrorCode DMPlexCreateHypercubicMesh(MPI_Comm comm, PetscInt dim, const PetscInt edges[], const PetscReal lower[], const PetscReal upper[], PetscInt overlap, DM *dm)
Collective
Input Parameters#
comm - The communicator for the
DM
objectdim - The spatial dimension
edges - Number of edges per dimension, or
NULL
for (1,) in 1D and (2, 2) in 2D and (1, 1, 1) in 3Dlower - The lower left corner, or
NULL
for (0, 0, 0)upper - The upper right corner, or
NULL
for (1, 1, 1)overlap - The number of vertices in each direction to include in the overlap (default is 1)
Output Parameter#
dm - The DM object
Note#
If you want to customize this mesh using options, you just need to
DMCreate(comm, &dm);
DMSetType(dm, DMPLEX);
DMSetFromOptions(dm);
and use the options on the DMSetFromOptions()
page.
The vertices are numbered is lexicographic order, and the dim edges exiting a vertex in the positive orthant are number consecutively,
18--0-19--2-20--4-18
| | | |
13 15 17 13
| | | |
24-12-25-14-26-16-24
| | | |
7 9 11 7
| | | |
21--6-22--8-23-10-21
| | | |
1 3 5 1
| | | |
18--0-19--2-20--4-18
See Also#
DMSetFromOptions()
, DMPlexCreateFromFile()
, DMPlexCreateHexCylinderMesh()
, DMSetType()
, DMCreate()
Level#
beginner
Location#
Implementations#
DMPlexCreateHypercubicMesh_Internal() in src/dm/impls/plex/plexcreate.c
Index of all DMPlex routines
Table of Contents for all manual pages
Index of all manual pages