xref: /petsc/src/dm/impls/plex/plexgenerate.c (revision 1ed6e3ff8437baa411029a28c2b08f047df9ad9a)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
2d9deefdfSMatthew G. Knepley 
3*cc4c1da9SBarry Smith /*@
4a1cb98faSBarry Smith   DMPlexInvertCell - Flips cell orientations since `DMPLEX` stores some of them internally with outward normals.
5d9deefdfSMatthew G. Knepley 
6d9deefdfSMatthew G. Knepley   Input Parameters:
796ca5757SLisandro Dalcin + cellType - The cell type
8d9deefdfSMatthew G. Knepley - cone     - The incoming cone
9d9deefdfSMatthew G. Knepley 
10d9deefdfSMatthew G. Knepley   Output Parameter:
11d9deefdfSMatthew G. Knepley . cone - The inverted cone (in-place)
12d9deefdfSMatthew G. Knepley 
13d9deefdfSMatthew G. Knepley   Level: developer
14d9deefdfSMatthew G. Knepley 
151cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPolytopeType`, `DMPlexGenerate()`
16d9deefdfSMatthew G. Knepley @*/
DMPlexInvertCell(DMPolytopeType cellType,PetscInt cone[])17d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexInvertCell(DMPolytopeType cellType, PetscInt cone[])
18d71ae5a4SJacob Faibussowitsch {
1996ca5757SLisandro Dalcin #define SWAPCONE(cone, i, j) \
2096ca5757SLisandro Dalcin   do { \
2196ca5757SLisandro Dalcin     PetscInt _cone_tmp; \
2296ca5757SLisandro Dalcin     _cone_tmp = (cone)[i]; \
2396ca5757SLisandro Dalcin     (cone)[i] = (cone)[j]; \
2496ca5757SLisandro Dalcin     (cone)[j] = _cone_tmp; \
2596ca5757SLisandro Dalcin   } while (0)
26d9deefdfSMatthew G. Knepley 
27d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
2896ca5757SLisandro Dalcin   switch (cellType) {
29d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
30d71ae5a4SJacob Faibussowitsch     break;
31d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEGMENT:
32d71ae5a4SJacob Faibussowitsch     break;
33d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
34d71ae5a4SJacob Faibussowitsch     break;
35d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
36d71ae5a4SJacob Faibussowitsch     break;
37d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
38d71ae5a4SJacob Faibussowitsch     break;
39d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
40d71ae5a4SJacob Faibussowitsch     SWAPCONE(cone, 2, 3);
41d71ae5a4SJacob Faibussowitsch     break;
42d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
43d71ae5a4SJacob Faibussowitsch     SWAPCONE(cone, 0, 1);
44d71ae5a4SJacob Faibussowitsch     break;
45d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
46d71ae5a4SJacob Faibussowitsch     SWAPCONE(cone, 1, 3);
47d71ae5a4SJacob Faibussowitsch     break;
48d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
49d71ae5a4SJacob Faibussowitsch     SWAPCONE(cone, 1, 2);
50d71ae5a4SJacob Faibussowitsch     break;
51d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
52d71ae5a4SJacob Faibussowitsch     break;
53d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
54d71ae5a4SJacob Faibussowitsch     break;
55d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
56d71ae5a4SJacob Faibussowitsch     SWAPCONE(cone, 1, 3);
57d71ae5a4SJacob Faibussowitsch     break;
58d71ae5a4SJacob Faibussowitsch   default:
59d71ae5a4SJacob Faibussowitsch     break;
60d9deefdfSMatthew G. Knepley   }
613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6296ca5757SLisandro Dalcin #undef SWAPCONE
6396ca5757SLisandro Dalcin }
6496ca5757SLisandro Dalcin 
65*cc4c1da9SBarry Smith /*@
66a1cb98faSBarry Smith   DMPlexReorderCell - Flips cell orientations since `DMPLEX` stores some of them internally with outward normals.
6796ca5757SLisandro Dalcin 
6896ca5757SLisandro Dalcin   Input Parameters:
69a1cb98faSBarry Smith + dm   - The `DMPLEX` object
7096ca5757SLisandro Dalcin . cell - The cell
7196ca5757SLisandro Dalcin - cone - The incoming cone
7296ca5757SLisandro Dalcin 
7396ca5757SLisandro Dalcin   Output Parameter:
7496ca5757SLisandro Dalcin . cone - The reordered cone (in-place)
7596ca5757SLisandro Dalcin 
7696ca5757SLisandro Dalcin   Level: developer
7796ca5757SLisandro Dalcin 
781cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPolytopeType`, `DMPlexGenerate()`
7996ca5757SLisandro Dalcin @*/
DMPlexReorderCell(DM dm,PetscInt cell,PetscInt cone[])80d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexReorderCell(DM dm, PetscInt cell, PetscInt cone[])
81d71ae5a4SJacob Faibussowitsch {
8296ca5757SLisandro Dalcin   DMPolytopeType cellType;
8396ca5757SLisandro Dalcin 
8496ca5757SLisandro Dalcin   PetscFunctionBegin;
859566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cell, &cellType));
869566063dSJacob Faibussowitsch   PetscCall(DMPlexInvertCell(cellType, cone));
873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
88d9deefdfSMatthew G. Knepley }
89d9deefdfSMatthew G. Knepley 
9094ef8ddeSSatish Balay /*@C
91d9deefdfSMatthew G. Knepley   DMPlexTriangleSetOptions - Set the options used for the Triangle mesh generator
92d9deefdfSMatthew G. Knepley 
93d9deefdfSMatthew G. Knepley   Not Collective
94d9deefdfSMatthew G. Knepley 
9560225df5SJacob Faibussowitsch   Input Parameters:
96a1cb98faSBarry Smith + dm   - The `DMPLEX` object
97d9deefdfSMatthew G. Knepley - opts - The command line options
98d9deefdfSMatthew G. Knepley 
99d9deefdfSMatthew G. Knepley   Level: developer
100d9deefdfSMatthew G. Knepley 
1011cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTetgenSetOptions()`, `DMPlexGenerate()`
102d9deefdfSMatthew G. Knepley @*/
DMPlexTriangleSetOptions(DM dm,const char * opts)103d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexTriangleSetOptions(DM dm, const char *opts)
104d71ae5a4SJacob Faibussowitsch {
105d9deefdfSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
106d9deefdfSMatthew G. Knepley 
107d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
108d9deefdfSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1094f572ea9SToby Isaac   PetscAssertPointer(opts, 2);
1109566063dSJacob Faibussowitsch   PetscCall(PetscFree(mesh->triangleOpts));
1119566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(opts, &mesh->triangleOpts));
1123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
113d9deefdfSMatthew G. Knepley }
114d9deefdfSMatthew G. Knepley 
11594ef8ddeSSatish Balay /*@C
116d9deefdfSMatthew G. Knepley   DMPlexTetgenSetOptions - Set the options used for the Tetgen mesh generator
117d9deefdfSMatthew G. Knepley 
118d9deefdfSMatthew G. Knepley   Not Collective
119d9deefdfSMatthew G. Knepley 
12060225df5SJacob Faibussowitsch   Input Parameters:
121a1cb98faSBarry Smith + dm   - The `DMPLEX` object
122d9deefdfSMatthew G. Knepley - opts - The command line options
123d9deefdfSMatthew G. Knepley 
124d9deefdfSMatthew G. Knepley   Level: developer
125d9deefdfSMatthew G. Knepley 
1261cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTriangleSetOptions()`, `DMPlexGenerate()`
127d9deefdfSMatthew G. Knepley @*/
DMPlexTetgenSetOptions(DM dm,const char * opts)128d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexTetgenSetOptions(DM dm, const char *opts)
129d71ae5a4SJacob Faibussowitsch {
130d9deefdfSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
131d9deefdfSMatthew G. Knepley 
132d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
133d9deefdfSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1344f572ea9SToby Isaac   PetscAssertPointer(opts, 2);
1359566063dSJacob Faibussowitsch   PetscCall(PetscFree(mesh->tetgenOpts));
1369566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(opts, &mesh->tetgenOpts));
1373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
138d9deefdfSMatthew G. Knepley }
139d9deefdfSMatthew G. Knepley 
140*cc4c1da9SBarry Smith /*@
141d9deefdfSMatthew G. Knepley   DMPlexGenerate - Generates a mesh.
142d9deefdfSMatthew G. Knepley 
143d9deefdfSMatthew G. Knepley   Not Collective
144d9deefdfSMatthew G. Knepley 
145d9deefdfSMatthew G. Knepley   Input Parameters:
146a1cb98faSBarry Smith + boundary    - The `DMPLEX` boundary object
147d9deefdfSMatthew G. Knepley . name        - The mesh generation package name
148d9deefdfSMatthew G. Knepley - interpolate - Flag to create intermediate mesh elements
149d9deefdfSMatthew G. Knepley 
150d9deefdfSMatthew G. Knepley   Output Parameter:
151a1cb98faSBarry Smith . mesh - The `DMPLEX` object
152d9deefdfSMatthew G. Knepley 
153a1cb98faSBarry Smith   Options Database Keys:
1542c4dacf2SBarry Smith + -dm_plex_generate <name> - package to generate mesh, for example, triangle, ctetgen or tetgen
155c0517cd5SMatthew G. Knepley - -dm_generator <name>     - package to generate mesh, for example, triangle, ctetgen or tetgen
1562c4dacf2SBarry Smith 
157d9deefdfSMatthew G. Knepley   Level: intermediate
158d9deefdfSMatthew G. Knepley 
1591cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMRefine()`
160d9deefdfSMatthew G. Knepley @*/
DMPlexGenerate(DM boundary,const char name[],PetscBool interpolate,DM * mesh)161d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGenerate(DM boundary, const char name[], PetscBool interpolate, DM *mesh)
162d71ae5a4SJacob Faibussowitsch {
163c0517cd5SMatthew G. Knepley   DMGeneratorFunctionList fl;
1643f2a96e3SMatthew G. Knepley   char                    genname[PETSC_MAX_PATH_LEN];
1653f2a96e3SMatthew G. Knepley   const char             *suggestions;
166d9deefdfSMatthew G. Knepley   PetscInt                dim;
1673a074057SBarry Smith   PetscBool               flg;
168d9deefdfSMatthew G. Knepley 
169d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
170d9deefdfSMatthew G. Knepley   PetscValidHeaderSpecific(boundary, DM_CLASSID, 1);
171064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(boundary, interpolate, 3);
1729566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(boundary, &dim));
1739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(((PetscObject)boundary)->options, ((PetscObject)boundary)->prefix, "-dm_generator", genname, sizeof(genname), &flg));
174d9deefdfSMatthew G. Knepley   if (flg) name = genname;
1752c4dacf2SBarry Smith   else {
1769566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetString(((PetscObject)boundary)->options, ((PetscObject)boundary)->prefix, "-dm_plex_generate", genname, sizeof(genname), &flg));
1772c4dacf2SBarry Smith     if (flg) name = genname;
1782c4dacf2SBarry Smith   }
1793a074057SBarry Smith 
180c0517cd5SMatthew G. Knepley   fl = DMGenerateList;
181d9deefdfSMatthew G. Knepley   if (name) {
1823a074057SBarry Smith     while (fl) {
1839566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(fl->name, name, &flg));
1843a074057SBarry Smith       if (flg) {
1859566063dSJacob Faibussowitsch         PetscCall((*fl->generate)(boundary, interpolate, mesh));
1863ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
187d9deefdfSMatthew G. Knepley       }
1883a074057SBarry Smith       fl = fl->next;
189d9deefdfSMatthew G. Knepley     }
19098921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Grid generator %s not registered; you may need to add --download-%s to your ./configure options", name, name);
1913a074057SBarry Smith   } else {
1923a074057SBarry Smith     while (fl) {
1933a074057SBarry Smith       if (boundary->dim == fl->dim) {
1949566063dSJacob Faibussowitsch         PetscCall((*fl->generate)(boundary, interpolate, mesh));
1953ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
1963a074057SBarry Smith       }
1973a074057SBarry Smith       fl = fl->next;
1983a074057SBarry Smith     }
1992c4dacf2SBarry Smith     suggestions = "";
2002c4dacf2SBarry Smith     if (boundary->dim + 1 == 2) suggestions = " You may need to add --download-triangle to your ./configure options";
2012c4dacf2SBarry Smith     else if (boundary->dim + 1 == 3) suggestions = " You may need to add --download-ctetgen or --download-tetgen in your ./configure options";
20263a3b9bcSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No grid generator of dimension %" PetscInt_FMT " registered%s", boundary->dim + 1, suggestions);
2033a074057SBarry Smith   }
2043a074057SBarry Smith }
205