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