xref: /petsc/src/dm/impls/plex/plexgenerate.c (revision 2da392cc7c10228af19ad9843ce5155178acb644)
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 PetscFunctionList DMPlexGenerateList = NULL;
129 
130 struct _n_PetscFunctionList {
131   PetscErrorCode    (*generate)(DM, PetscBool, DM*);
132   PetscErrorCode    (*refine)(DM,double*, DM*);
133   char              *name;               /* string to identify routine */
134   PetscInt          dim;
135   PetscFunctionList next;                /* next pointer */
136 };
137 
138 /*@C
139   DMPlexGenerate - Generates a mesh.
140 
141   Not Collective
142 
143   Input Parameters:
144 + boundary - The DMPlex boundary object
145 . name - The mesh generation package name
146 - interpolate - Flag to create intermediate mesh elements
147 
148   Output Parameter:
149 . mesh - The DMPlex object
150 
151   Options Database:
152 +  -dm_plex_generate <name> - package to generate mesh, for example, triangle, ctetgen or tetgen
153 -  -dm_plex_generator <name> - package to generate mesh, for example, triangle, ctetgen or tetgen (deprecated)
154 
155   Level: intermediate
156 
157 .seealso: DMPlexCreate(), DMRefine()
158 @*/
159 PetscErrorCode DMPlexGenerate(DM boundary, const char name[], PetscBool interpolate, DM *mesh)
160 {
161   PetscInt          dim;
162   char              genname[1024];
163   PetscBool         flg;
164   PetscErrorCode    ierr;
165   PetscFunctionList fl;
166   const char*       suggestions;
167 
168   PetscFunctionBegin;
169   PetscValidHeaderSpecific(boundary, DM_CLASSID, 1);
170   PetscValidLogicalCollectiveBool(boundary, interpolate, 2);
171   ierr = DMGetDimension(boundary, &dim);CHKERRQ(ierr);
172   ierr = PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_plex_generator", genname, sizeof(genname), &flg);CHKERRQ(ierr);
173   if (flg) name = genname;
174   else {
175     ierr = PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_plex_generate", genname, sizeof(genname), &flg);CHKERRQ(ierr);
176     if (flg) name = genname;
177   }
178 
179   fl = DMPlexGenerateList;
180   if (name) {
181     while (fl) {
182       ierr = PetscStrcmp(fl->name,name,&flg);CHKERRQ(ierr);
183       if (flg) {
184         ierr = (*fl->generate)(boundary,interpolate,mesh);CHKERRQ(ierr);
185         PetscFunctionReturn(0);
186       }
187       fl = fl->next;
188     }
189     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);
190   } else {
191     while (fl) {
192       if (boundary->dim == fl->dim) {
193         ierr = (*fl->generate)(boundary,interpolate,mesh);CHKERRQ(ierr);
194         PetscFunctionReturn(0);
195       }
196       fl = fl->next;
197     }
198     suggestions = "";
199     if (boundary->dim+1 == 2) suggestions = " You may need to add --download-triangle to your ./configure options";
200     else if (boundary->dim+1 == 3) suggestions = " You may need to add --download-ctetgen or --download-tetgen in your ./configure options";
201     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid generator of dimension %D registered%s",boundary->dim+1,suggestions);
202   }
203   PetscFunctionReturn(0);
204 }
205 
206 /*@C
207   DMPlexGenerateRegister -  Adds a grid generator to DMPlex
208 
209    Not Collective
210 
211    Input Parameters:
212 +  name_solver - name of a new user-defined grid generator
213 .  fnc - generator function
214 .  rfnc - refinement function
215 -  dim - dimension of boundary of domain
216 
217    Notes:
218    DMPlexGenerateRegister() may be called multiple times to add several user-defined solvers.
219 
220    Sample usage:
221 .vb
222    DMPlexGenerateRegister("my_generator",MyGeneratorCreate,MyGeneratorRefiner,dim);
223 .ve
224 
225    Then, your generator can be chosen with the procedural interface via
226 $     DMPlexGenerate(dm,"my_generator",...)
227    or at runtime via the option
228 $     -dm_plex_generator my_generator
229 
230    Level: advanced
231 
232 .seealso: DMPlexGenerateRegisterAll(), DMPlexGenerate(), DMPlexGenerateRegisterDestroy()
233 
234 @*/
235 PetscErrorCode  DMPlexGenerateRegister(const char sname[],PetscErrorCode (*fnc)(DM, PetscBool,DM*), PetscErrorCode (*rfnc)(DM, double*,DM*),PetscInt dim)
236 {
237   PetscErrorCode    ierr;
238   PetscFunctionList entry;
239 
240   PetscFunctionBegin;
241   ierr            = PetscNew(&entry);CHKERRQ(ierr);
242   ierr            = PetscStrallocpy(sname,&entry->name);CHKERRQ(ierr);
243   entry->generate = fnc;
244   entry->refine   = rfnc;
245   entry->dim      = dim;
246   entry->next     = NULL;
247   if (!DMPlexGenerateList) DMPlexGenerateList = entry;
248   else {
249     PetscFunctionList fl = DMPlexGenerateList;
250     while (fl->next) fl = fl->next;
251     fl->next = entry;
252   }
253   PetscFunctionReturn(0);
254 }
255 
256 extern PetscBool DMPlexGenerateRegisterAllCalled;
257 
258 PetscErrorCode  DMPlexGenerateRegisterDestroy(void)
259 {
260   PetscFunctionList next,fl;
261   PetscErrorCode    ierr;
262 
263   PetscFunctionBegin;
264   next = fl =  DMPlexGenerateList;
265     while (next) {
266     next = fl ? fl->next : NULL;
267     if (fl) {ierr = PetscFree(fl->name);CHKERRQ(ierr);}
268     ierr = PetscFree(fl);CHKERRQ(ierr);
269     fl = next;
270   }
271   DMPlexGenerateList              = NULL;
272   DMPlexGenerateRegisterAllCalled = PETSC_FALSE;
273   PetscFunctionReturn(0);
274 }
275