xref: /petsc/src/dm/impls/plex/plexgenerate.c (revision 4e278199b78715991f5c71ebbd945c1489263e6c)
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   PetscErrorCode ierr;
65 
66   PetscFunctionBegin;
67   ierr = DMPlexGetCellType(dm, cell, &cellType);CHKERRQ(ierr);
68   ierr = DMPlexInvertCell(cellType, cone);CHKERRQ(ierr);
69   PetscFunctionReturn(0);
70 }
71 
72 /*@C
73   DMPlexTriangleSetOptions - Set the options used for the Triangle mesh generator
74 
75   Not Collective
76 
77   Inputs Parameters:
78 + dm - The DMPlex object
79 - opts - The command line options
80 
81   Level: developer
82 
83 .seealso: DMPlexTetgenSetOptions(), DMPlexGenerate()
84 @*/
85 PetscErrorCode DMPlexTriangleSetOptions(DM dm, const char *opts)
86 {
87   DM_Plex       *mesh = (DM_Plex*) dm->data;
88   PetscErrorCode ierr;
89 
90   PetscFunctionBegin;
91   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
92   PetscValidPointer(opts, 2);
93   ierr = PetscFree(mesh->triangleOpts);CHKERRQ(ierr);
94   ierr = PetscStrallocpy(opts, &mesh->triangleOpts);CHKERRQ(ierr);
95   PetscFunctionReturn(0);
96 }
97 
98 /*@C
99   DMPlexTetgenSetOptions - Set the options used for the Tetgen mesh generator
100 
101   Not Collective
102 
103   Inputs Parameters:
104 + dm - The DMPlex object
105 - opts - The command line options
106 
107   Level: developer
108 
109 .seealso: DMPlexTriangleSetOptions(), DMPlexGenerate()
110 @*/
111 PetscErrorCode DMPlexTetgenSetOptions(DM dm, const char *opts)
112 {
113   DM_Plex       *mesh = (DM_Plex*) dm->data;
114   PetscErrorCode ierr;
115 
116   PetscFunctionBegin;
117   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
118   PetscValidPointer(opts, 2);
119   ierr = PetscFree(mesh->tetgenOpts);CHKERRQ(ierr);
120   ierr = PetscStrallocpy(opts, &mesh->tetgenOpts);CHKERRQ(ierr);
121   PetscFunctionReturn(0);
122 }
123 
124 /*
125    Contains the list of registered DMPlexGenerators routines
126 */
127 PlexGeneratorFunctionList DMPlexGenerateList = NULL;
128 
129 /*@C
130   DMPlexGenerate - Generates a mesh.
131 
132   Not Collective
133 
134   Input Parameters:
135 + boundary - The DMPlex boundary object
136 . name - The mesh generation package name
137 - interpolate - Flag to create intermediate mesh elements
138 
139   Output Parameter:
140 . mesh - The DMPlex object
141 
142   Options Database:
143 +  -dm_plex_generate <name> - package to generate mesh, for example, triangle, ctetgen or tetgen
144 -  -dm_plex_generator <name> - package to generate mesh, for example, triangle, ctetgen or tetgen (deprecated)
145 
146   Level: intermediate
147 
148 .seealso: DMPlexCreate(), DMRefine()
149 @*/
150 PetscErrorCode DMPlexGenerate(DM boundary, const char name[], PetscBool interpolate, DM *mesh)
151 {
152   PlexGeneratorFunctionList fl;
153   char                      genname[PETSC_MAX_PATH_LEN];
154   const char               *suggestions;
155   PetscInt                  dim;
156   PetscBool                 flg;
157   PetscErrorCode            ierr;
158 
159   PetscFunctionBegin;
160   PetscValidHeaderSpecific(boundary, DM_CLASSID, 1);
161   PetscValidLogicalCollectiveBool(boundary, interpolate, 3);
162   ierr = DMGetDimension(boundary, &dim);CHKERRQ(ierr);
163   ierr = PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_plex_generator", genname, sizeof(genname), &flg);CHKERRQ(ierr);
164   if (flg) name = genname;
165   else {
166     ierr = PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_plex_generate", genname, sizeof(genname), &flg);CHKERRQ(ierr);
167     if (flg) name = genname;
168   }
169 
170   fl = DMPlexGenerateList;
171   if (name) {
172     while (fl) {
173       ierr = PetscStrcmp(fl->name,name,&flg);CHKERRQ(ierr);
174       if (flg) {
175         ierr = (*fl->generate)(boundary,interpolate,mesh);CHKERRQ(ierr);
176         PetscFunctionReturn(0);
177       }
178       fl = fl->next;
179     }
180     SETERRQ2(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);
181   } else {
182     while (fl) {
183       if (boundary->dim == fl->dim) {
184         ierr = (*fl->generate)(boundary,interpolate,mesh);CHKERRQ(ierr);
185         PetscFunctionReturn(0);
186       }
187       fl = fl->next;
188     }
189     suggestions = "";
190     if (boundary->dim+1 == 2) suggestions = " You may need to add --download-triangle to your ./configure options";
191     else if (boundary->dim+1 == 3) suggestions = " You may need to add --download-ctetgen or --download-tetgen in your ./configure options";
192     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid generator of dimension %D registered%s",boundary->dim+1,suggestions);
193   }
194 }
195 
196 /*@C
197   DMPlexGenerateRegister -  Adds a grid generator to DMPlex
198 
199    Not Collective
200 
201    Input Parameters:
202 +  name_solver - name of a new user-defined grid generator
203 .  fnc - generator function
204 .  rfnc - refinement function
205 .  alfnc - adapt by label function
206 -  dim - dimension of boundary of domain
207 
208    Notes:
209    DMPlexGenerateRegister() may be called multiple times to add several user-defined solvers.
210 
211    Sample usage:
212 .vb
213    DMPlexGenerateRegister("my_generator",MyGeneratorCreate,MyGeneratorRefiner,MyGeneratorAdaptor,dim);
214 .ve
215 
216    Then, your generator can be chosen with the procedural interface via
217 $     DMPlexGenerate(dm,"my_generator",...)
218    or at runtime via the option
219 $     -dm_plex_generator my_generator
220 
221    Level: advanced
222 
223 .seealso: DMPlexGenerateRegisterAll(), DMPlexGenerate(), DMPlexGenerateRegisterDestroy()
224 
225 @*/
226 PetscErrorCode  DMPlexGenerateRegister(const char sname[], PetscErrorCode (*fnc)(DM, PetscBool, DM*), PetscErrorCode (*rfnc)(DM, PetscReal*, DM*), PetscErrorCode (*alfnc)(DM, DMLabel, DM*), PetscInt dim)
227 {
228   PlexGeneratorFunctionList entry;
229   PetscErrorCode            ierr;
230 
231   PetscFunctionBegin;
232   ierr = PetscNew(&entry);CHKERRQ(ierr);
233   ierr = PetscStrallocpy(sname,&entry->name);CHKERRQ(ierr);
234   entry->generate   = fnc;
235   entry->refine     = rfnc;
236   entry->adaptlabel = alfnc;
237   entry->dim        = dim;
238   entry->next       = NULL;
239   if (!DMPlexGenerateList) DMPlexGenerateList = entry;
240   else {
241     PlexGeneratorFunctionList fl = DMPlexGenerateList;
242     while (fl->next) fl = fl->next;
243     fl->next = entry;
244   }
245   PetscFunctionReturn(0);
246 }
247 
248 extern PetscBool DMPlexGenerateRegisterAllCalled;
249 
250 PetscErrorCode  DMPlexGenerateRegisterDestroy(void)
251 {
252   PlexGeneratorFunctionList next,fl;
253   PetscErrorCode            ierr;
254 
255   PetscFunctionBegin;
256   next = fl =  DMPlexGenerateList;
257     while (next) {
258     next = fl ? fl->next : NULL;
259     if (fl) {ierr = PetscFree(fl->name);CHKERRQ(ierr);}
260     ierr = PetscFree(fl);CHKERRQ(ierr);
261     fl = next;
262   }
263   DMPlexGenerateList              = NULL;
264   DMPlexGenerateRegisterAllCalled = PETSC_FALSE;
265   PetscFunctionReturn(0);
266 }
267