xref: /petsc/src/dm/impls/plex/plexgenerate.c (revision 3f02e49b19195914bf17f317a25cb39636853415)
1 #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
2 
3 /*@
4   DMPlexInvertCell - Flips cell orientations since `DMPLEX` stores some of them internally with outward normals.
5 
6   Input Parameters:
7 + cellType - The cell type
8 - cone     - The incoming cone
9 
10   Output Parameter:
11 . cone - The inverted cone (in-place)
12 
13   Level: developer
14 
15 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPolytopeType`, `DMPlexGenerate()`
16 @*/
17 PetscErrorCode DMPlexInvertCell(DMPolytopeType cellType, PetscInt cone[])
18 {
19 #define SWAPCONE(cone, i, j) \
20   do { \
21     PetscInt _cone_tmp; \
22     _cone_tmp = (cone)[i]; \
23     (cone)[i] = (cone)[j]; \
24     (cone)[j] = _cone_tmp; \
25   } while (0)
26 
27   PetscFunctionBegin;
28   switch (cellType) {
29   case DM_POLYTOPE_POINT:
30     break;
31   case DM_POLYTOPE_SEGMENT:
32     break;
33   case DM_POLYTOPE_POINT_PRISM_TENSOR:
34     break;
35   case DM_POLYTOPE_TRIANGLE:
36     break;
37   case DM_POLYTOPE_QUADRILATERAL:
38     break;
39   case DM_POLYTOPE_SEG_PRISM_TENSOR:
40     SWAPCONE(cone, 2, 3);
41     break;
42   case DM_POLYTOPE_TETRAHEDRON:
43     SWAPCONE(cone, 0, 1);
44     break;
45   case DM_POLYTOPE_HEXAHEDRON:
46     SWAPCONE(cone, 1, 3);
47     break;
48   case DM_POLYTOPE_TRI_PRISM:
49     SWAPCONE(cone, 1, 2);
50     break;
51   case DM_POLYTOPE_TRI_PRISM_TENSOR:
52     break;
53   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
54     break;
55   case DM_POLYTOPE_PYRAMID:
56     SWAPCONE(cone, 1, 3);
57     break;
58   default:
59     break;
60   }
61   PetscFunctionReturn(PETSC_SUCCESS);
62 #undef SWAPCONE
63 }
64 
65 /*@
66   DMPlexReorderCell - Flips cell orientations since `DMPLEX` stores some of them internally with outward normals.
67 
68   Input Parameters:
69 + dm   - The `DMPLEX` object
70 . cell - The cell
71 - cone - The incoming cone
72 
73   Output Parameter:
74 . cone - The reordered cone (in-place)
75 
76   Level: developer
77 
78 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPolytopeType`, `DMPlexGenerate()`
79 @*/
80 PetscErrorCode DMPlexReorderCell(DM dm, PetscInt cell, PetscInt cone[])
81 {
82   DMPolytopeType cellType;
83 
84   PetscFunctionBegin;
85   PetscCall(DMPlexGetCellType(dm, cell, &cellType));
86   PetscCall(DMPlexInvertCell(cellType, cone));
87   PetscFunctionReturn(PETSC_SUCCESS);
88 }
89 
90 /*@C
91   DMPlexTriangleSetOptions - Set the options used for the Triangle mesh generator
92 
93   Not Collective
94 
95   Input Parameters:
96 + dm   - The `DMPLEX` object
97 - opts - The command line options
98 
99   Level: developer
100 
101 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTetgenSetOptions()`, `DMPlexGenerate()`
102 @*/
103 PetscErrorCode DMPlexTriangleSetOptions(DM dm, const char *opts)
104 {
105   DM_Plex *mesh = (DM_Plex *)dm->data;
106 
107   PetscFunctionBegin;
108   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
109   PetscAssertPointer(opts, 2);
110   PetscCall(PetscFree(mesh->triangleOpts));
111   PetscCall(PetscStrallocpy(opts, &mesh->triangleOpts));
112   PetscFunctionReturn(PETSC_SUCCESS);
113 }
114 
115 /*@C
116   DMPlexTetgenSetOptions - Set the options used for the Tetgen mesh generator
117 
118   Not Collective
119 
120   Input Parameters:
121 + dm   - The `DMPLEX` object
122 - opts - The command line options
123 
124   Level: developer
125 
126 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTriangleSetOptions()`, `DMPlexGenerate()`
127 @*/
128 PetscErrorCode DMPlexTetgenSetOptions(DM dm, const char *opts)
129 {
130   DM_Plex *mesh = (DM_Plex *)dm->data;
131 
132   PetscFunctionBegin;
133   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
134   PetscAssertPointer(opts, 2);
135   PetscCall(PetscFree(mesh->tetgenOpts));
136   PetscCall(PetscStrallocpy(opts, &mesh->tetgenOpts));
137   PetscFunctionReturn(PETSC_SUCCESS);
138 }
139 
140 /*@
141   DMPlexGenerate - Generates a mesh.
142 
143   Not Collective
144 
145   Input Parameters:
146 + boundary    - The `DMPLEX` boundary object
147 . name        - The mesh generation package name
148 - interpolate - Flag to create intermediate mesh elements
149 
150   Output Parameter:
151 . mesh - The `DMPLEX` object
152 
153   Options Database Keys:
154 + -dm_plex_generate <name> - package to generate mesh, for example, triangle, ctetgen or tetgen
155 - -dm_generator <name>     - package to generate mesh, for example, triangle, ctetgen or tetgen
156 
157   Level: intermediate
158 
159 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMRefine()`
160 @*/
161 PetscErrorCode DMPlexGenerate(DM boundary, const char name[], PetscBool interpolate, DM *mesh)
162 {
163   DMGeneratorFunctionList fl;
164   char                    genname[PETSC_MAX_PATH_LEN];
165   const char             *suggestions;
166   PetscInt                dim;
167   PetscBool               flg;
168 
169   PetscFunctionBegin;
170   PetscValidHeaderSpecific(boundary, DM_CLASSID, 1);
171   PetscValidLogicalCollectiveBool(boundary, interpolate, 3);
172   PetscCall(DMGetDimension(boundary, &dim));
173   PetscCall(PetscOptionsGetString(((PetscObject)boundary)->options, ((PetscObject)boundary)->prefix, "-dm_generator", genname, sizeof(genname), &flg));
174   if (flg) name = genname;
175   else {
176     PetscCall(PetscOptionsGetString(((PetscObject)boundary)->options, ((PetscObject)boundary)->prefix, "-dm_plex_generate", genname, sizeof(genname), &flg));
177     if (flg) name = genname;
178   }
179 
180   fl = DMGenerateList;
181   if (name) {
182     while (fl) {
183       PetscCall(PetscStrcmp(fl->name, name, &flg));
184       if (flg) {
185         PetscCall((*fl->generate)(boundary, interpolate, mesh));
186         PetscFunctionReturn(PETSC_SUCCESS);
187       }
188       fl = fl->next;
189     }
190     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);
191   } else {
192     while (fl) {
193       if (boundary->dim == fl->dim) {
194         PetscCall((*fl->generate)(boundary, interpolate, mesh));
195         PetscFunctionReturn(PETSC_SUCCESS);
196       }
197       fl = fl->next;
198     }
199     suggestions = "";
200     if (boundary->dim + 1 == 2) suggestions = " You may need to add --download-triangle to your ./configure options";
201     else if (boundary->dim + 1 == 3) suggestions = " You may need to add --download-ctetgen or --download-tetgen in your ./configure options";
202     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No grid generator of dimension %" PetscInt_FMT " registered%s", boundary->dim + 1, suggestions);
203   }
204 }
205