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