xref: /petsc/src/dm/impls/plex/plexgenerate.c (revision 3a074057bf25314dc21e177b1a2e45e66b0fa3b1)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2d9deefdfSMatthew G. Knepley 
3d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexInvertCell_Internal(PetscInt dim, PetscInt numCorners, PetscInt cone[])
4d9deefdfSMatthew G. Knepley {
5d9deefdfSMatthew G. Knepley   int tmpc;
6d9deefdfSMatthew G. Knepley 
7d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
8d9deefdfSMatthew G. Knepley   if (dim != 3) PetscFunctionReturn(0);
9d9deefdfSMatthew G. Knepley   switch (numCorners) {
10d9deefdfSMatthew G. Knepley   case 4:
11d9deefdfSMatthew G. Knepley     tmpc    = cone[0];
12d9deefdfSMatthew G. Knepley     cone[0] = cone[1];
13d9deefdfSMatthew G. Knepley     cone[1] = tmpc;
14d9deefdfSMatthew G. Knepley     break;
15d9deefdfSMatthew G. Knepley   case 8:
16d9deefdfSMatthew G. Knepley     tmpc    = cone[1];
17d9deefdfSMatthew G. Knepley     cone[1] = cone[3];
18d9deefdfSMatthew G. Knepley     cone[3] = tmpc;
19d9deefdfSMatthew G. Knepley     break;
20d9deefdfSMatthew G. Knepley   default: break;
21d9deefdfSMatthew G. Knepley   }
22d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
23d9deefdfSMatthew G. Knepley }
24d9deefdfSMatthew G. Knepley 
25d9deefdfSMatthew G. Knepley /*@C
26d9deefdfSMatthew G. Knepley   DMPlexInvertCell - This flips tetrahedron and hexahedron orientation since Plex stores them internally with outward normals. Other cells are left untouched.
27d9deefdfSMatthew G. Knepley 
28d9deefdfSMatthew G. Knepley   Input Parameters:
29d9deefdfSMatthew G. Knepley + numCorners - The number of vertices in a cell
30d9deefdfSMatthew G. Knepley - cone - The incoming cone
31d9deefdfSMatthew G. Knepley 
32d9deefdfSMatthew G. Knepley   Output Parameter:
33d9deefdfSMatthew G. Knepley . cone - The inverted cone (in-place)
34d9deefdfSMatthew G. Knepley 
35d9deefdfSMatthew G. Knepley   Level: developer
36d9deefdfSMatthew G. Knepley 
37d9deefdfSMatthew G. Knepley .seealso: DMPlexGenerate()
38d9deefdfSMatthew G. Knepley @*/
39d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexInvertCell(PetscInt dim, PetscInt numCorners, int cone[])
40d9deefdfSMatthew G. Knepley {
41d9deefdfSMatthew G. Knepley   int tmpc;
42d9deefdfSMatthew G. Knepley 
43d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
44d9deefdfSMatthew G. Knepley   if (dim != 3) PetscFunctionReturn(0);
45d9deefdfSMatthew G. Knepley   switch (numCorners) {
46d9deefdfSMatthew G. Knepley   case 4:
47d9deefdfSMatthew G. Knepley     tmpc    = cone[0];
48d9deefdfSMatthew G. Knepley     cone[0] = cone[1];
49d9deefdfSMatthew G. Knepley     cone[1] = tmpc;
50d9deefdfSMatthew G. Knepley     break;
51d9deefdfSMatthew G. Knepley   case 8:
52d9deefdfSMatthew G. Knepley     tmpc    = cone[1];
53d9deefdfSMatthew G. Knepley     cone[1] = cone[3];
54d9deefdfSMatthew G. Knepley     cone[3] = tmpc;
55d9deefdfSMatthew G. Knepley     break;
56d9deefdfSMatthew G. Knepley   default: break;
57d9deefdfSMatthew G. Knepley   }
58d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
59d9deefdfSMatthew G. Knepley }
60d9deefdfSMatthew G. Knepley 
61d9deefdfSMatthew G. Knepley 
6294ef8ddeSSatish Balay /*@C
63d9deefdfSMatthew G. Knepley   DMPlexTriangleSetOptions - Set the options used for the Triangle mesh generator
64d9deefdfSMatthew G. Knepley 
65d9deefdfSMatthew G. Knepley   Not Collective
66d9deefdfSMatthew G. Knepley 
67d9deefdfSMatthew G. Knepley   Inputs Parameters:
68d9deefdfSMatthew G. Knepley + dm - The DMPlex object
69d9deefdfSMatthew G. Knepley - opts - The command line options
70d9deefdfSMatthew G. Knepley 
71d9deefdfSMatthew G. Knepley   Level: developer
72d9deefdfSMatthew G. Knepley 
73d9deefdfSMatthew G. Knepley .keywords: mesh, points
74d9deefdfSMatthew G. Knepley .seealso: DMPlexTetgenSetOptions(), DMPlexGenerate()
75d9deefdfSMatthew G. Knepley @*/
76d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexTriangleSetOptions(DM dm, const char *opts)
77d9deefdfSMatthew G. Knepley {
78d9deefdfSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex*) dm->data;
79d9deefdfSMatthew G. Knepley   PetscErrorCode ierr;
80d9deefdfSMatthew G. Knepley 
81d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
82d9deefdfSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
83d9deefdfSMatthew G. Knepley   PetscValidPointer(opts, 2);
84606ac3a5SMatthew G. Knepley   ierr = PetscFree(mesh->triangleOpts);CHKERRQ(ierr);
85606ac3a5SMatthew G. Knepley   ierr = PetscStrallocpy(opts, &mesh->triangleOpts);CHKERRQ(ierr);
86d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
87d9deefdfSMatthew G. Knepley }
88d9deefdfSMatthew G. Knepley 
8994ef8ddeSSatish Balay /*@C
90d9deefdfSMatthew G. Knepley   DMPlexTetgenSetOptions - Set the options used for the Tetgen mesh generator
91d9deefdfSMatthew G. Knepley 
92d9deefdfSMatthew G. Knepley   Not Collective
93d9deefdfSMatthew G. Knepley 
94d9deefdfSMatthew G. Knepley   Inputs Parameters:
95d9deefdfSMatthew G. Knepley + dm - The DMPlex object
96d9deefdfSMatthew G. Knepley - opts - The command line options
97d9deefdfSMatthew G. Knepley 
98d9deefdfSMatthew G. Knepley   Level: developer
99d9deefdfSMatthew G. Knepley 
100d9deefdfSMatthew G. Knepley .keywords: mesh, points
101d9deefdfSMatthew G. Knepley .seealso: DMPlexTriangleSetOptions(), DMPlexGenerate()
102d9deefdfSMatthew G. Knepley @*/
103d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexTetgenSetOptions(DM dm, const char *opts)
104d9deefdfSMatthew G. Knepley {
105d9deefdfSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex*) dm->data;
106d9deefdfSMatthew G. Knepley   PetscErrorCode ierr;
107d9deefdfSMatthew G. Knepley 
108d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
109d9deefdfSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
110d9deefdfSMatthew G. Knepley   PetscValidPointer(opts, 2);
111606ac3a5SMatthew G. Knepley   ierr = PetscFree(mesh->tetgenOpts);CHKERRQ(ierr);
112606ac3a5SMatthew G. Knepley   ierr = PetscStrallocpy(opts, &mesh->tetgenOpts);CHKERRQ(ierr);
113d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
114d9deefdfSMatthew G. Knepley }
115d9deefdfSMatthew G. Knepley 
116*3a074057SBarry Smith /*
117*3a074057SBarry Smith    Contains the list of registered DMPlexGenerators routines
118*3a074057SBarry Smith */
119*3a074057SBarry Smith PetscFunctionList DMPlexGenerateList = NULL;
120d9deefdfSMatthew G. Knepley 
121*3a074057SBarry Smith struct _n_PetscFunctionList {
122*3a074057SBarry Smith   PetscErrorCode    (*generate)(DM, PetscBool, DM*);
123*3a074057SBarry Smith   PetscErrorCode    (*refine)(DM,double*, DM*);
124*3a074057SBarry Smith   char              *name;               /* string to identify routine */
125*3a074057SBarry Smith   PetscInt          dim;
126*3a074057SBarry Smith   PetscFunctionList next;                /* next pointer */
127*3a074057SBarry Smith };
128d9deefdfSMatthew G. Knepley 
129d9deefdfSMatthew G. Knepley /*@C
130d9deefdfSMatthew G. Knepley   DMPlexGenerate - Generates a mesh.
131d9deefdfSMatthew G. Knepley 
132d9deefdfSMatthew G. Knepley   Not Collective
133d9deefdfSMatthew G. Knepley 
134d9deefdfSMatthew G. Knepley   Input Parameters:
135d9deefdfSMatthew G. Knepley + boundary - The DMPlex boundary object
136d9deefdfSMatthew G. Knepley . name - The mesh generation package name
137d9deefdfSMatthew G. Knepley - interpolate - Flag to create intermediate mesh elements
138d9deefdfSMatthew G. Knepley 
139d9deefdfSMatthew G. Knepley   Output Parameter:
140d9deefdfSMatthew G. Knepley . mesh - The DMPlex object
141d9deefdfSMatthew G. Knepley 
142d9deefdfSMatthew G. Knepley   Level: intermediate
143d9deefdfSMatthew G. Knepley 
144d9deefdfSMatthew G. Knepley .keywords: mesh, elements
145d9deefdfSMatthew G. Knepley .seealso: DMPlexCreate(), DMRefine()
146d9deefdfSMatthew G. Knepley @*/
147d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexGenerate(DM boundary, const char name[], PetscBool interpolate, DM *mesh)
148d9deefdfSMatthew G. Knepley {
149d9deefdfSMatthew G. Knepley   PetscInt          dim;
150d9deefdfSMatthew G. Knepley   char              genname[1024];
151*3a074057SBarry Smith   PetscBool         flg;
152d9deefdfSMatthew G. Knepley   PetscErrorCode    ierr;
153*3a074057SBarry Smith   PetscFunctionList fl;
154d9deefdfSMatthew G. Knepley 
155d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
156d9deefdfSMatthew G. Knepley   PetscValidHeaderSpecific(boundary, DM_CLASSID, 1);
157d9deefdfSMatthew G. Knepley   PetscValidLogicalCollectiveBool(boundary, interpolate, 2);
158d9deefdfSMatthew G. Knepley   ierr = DMGetDimension(boundary, &dim);CHKERRQ(ierr);
159c5929fdfSBarry Smith   ierr = PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_plex_generator", genname, 1024, &flg);CHKERRQ(ierr);
160d9deefdfSMatthew G. Knepley   if (flg) name = genname;
161*3a074057SBarry Smith 
162*3a074057SBarry Smith   fl = DMPlexGenerateList;
163d9deefdfSMatthew G. Knepley   if (name) {
164*3a074057SBarry Smith     while (fl) {
165*3a074057SBarry Smith       ierr = PetscStrcmp(fl->name,name,&flg);CHKERRQ(ierr);
166*3a074057SBarry Smith       if (flg) {
167*3a074057SBarry Smith         (*fl->generate)(boundary,interpolate,mesh);CHKERRQ(ierr);
168*3a074057SBarry Smith         PetscFunctionReturn(0);
169d9deefdfSMatthew G. Knepley       }
170*3a074057SBarry Smith       fl = fl->next;
171d9deefdfSMatthew G. Knepley     }
172*3a074057SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Grid generator %g not registered",name);CHKERRQ(ierr);
173*3a074057SBarry Smith   } else {
174*3a074057SBarry Smith     while (fl) {
175*3a074057SBarry Smith       if (boundary->dim == fl->dim) {
176*3a074057SBarry Smith         (*fl->generate)(boundary,interpolate,mesh);CHKERRQ(ierr);
177*3a074057SBarry Smith         PetscFunctionReturn(0);
178*3a074057SBarry Smith       }
179*3a074057SBarry Smith       fl = fl->next;
180*3a074057SBarry Smith     }
181*3a074057SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid generator of dimension %D registered",boundary->dim);CHKERRQ(ierr);
182*3a074057SBarry Smith   }
183*3a074057SBarry Smith   PetscFunctionReturn(0);
184*3a074057SBarry Smith }
185*3a074057SBarry Smith 
186*3a074057SBarry Smith /*@C
187*3a074057SBarry Smith   DMPlexGenerateRegister -  Adds a grid generator to DMPlex
188*3a074057SBarry Smith 
189*3a074057SBarry Smith    Not Collective
190*3a074057SBarry Smith 
191*3a074057SBarry Smith    Input Parameters:
192*3a074057SBarry Smith +  name_solver - name of a new user-defined solver
193*3a074057SBarry Smith .  fnc - generator function
194*3a074057SBarry Smith .  rfnc - refinement function
195*3a074057SBarry Smith -  dim - dimension of boundary of domain
196*3a074057SBarry Smith 
197*3a074057SBarry Smith    Notes:
198*3a074057SBarry Smith    DMPlexGenerateRegister() may be called multiple times to add several user-defined solvers.
199*3a074057SBarry Smith 
200*3a074057SBarry Smith    Sample usage:
201*3a074057SBarry Smith .vb
202*3a074057SBarry Smith    DMPlexGenerateRegister("my_generator",MyGeneratorCreate);
203*3a074057SBarry Smith .ve
204*3a074057SBarry Smith 
205*3a074057SBarry Smith    Then, your generator can be chosen with the procedural interface via
206*3a074057SBarry Smith $     DMPlexGenerateSetType(ksp,"my_generator")
207*3a074057SBarry Smith    or at runtime via the option
208*3a074057SBarry Smith $     -dmplex_generate_type my_generator
209*3a074057SBarry Smith 
210*3a074057SBarry Smith    Level: advanced
211*3a074057SBarry Smith 
212*3a074057SBarry Smith .keywords: DMPlexGenerate, register
213*3a074057SBarry Smith 
214*3a074057SBarry Smith .seealso: DMPlexGenerateRegisterAll(), DMPlexGenerateRegisterDestroy()
215*3a074057SBarry Smith 
216*3a074057SBarry Smith @*/
217*3a074057SBarry Smith PetscErrorCode  DMPlexGenerateRegister(const char sname[],PetscErrorCode (*fnc)(DM, PetscBool,DM*), PetscErrorCode (*rfnc)(DM, double*,DM*),PetscInt dim)
218*3a074057SBarry Smith {
219*3a074057SBarry Smith   PetscErrorCode    ierr;
220*3a074057SBarry Smith   PetscFunctionList entry;
221*3a074057SBarry Smith 
222*3a074057SBarry Smith   PetscFunctionBegin;
223*3a074057SBarry Smith   ierr            = PetscNew(&entry);CHKERRQ(ierr);
224*3a074057SBarry Smith   ierr            = PetscStrallocpy(sname,&entry->name);CHKERRQ(ierr);
225*3a074057SBarry Smith   entry->generate = fnc;
226*3a074057SBarry Smith   entry->refine   = rfnc;
227*3a074057SBarry Smith   entry->dim      = dim;
228*3a074057SBarry Smith   entry->next     = NULL;
229*3a074057SBarry Smith   if (!DMPlexGenerateList) DMPlexGenerateList = entry;
230*3a074057SBarry Smith   else {
231*3a074057SBarry Smith     PetscFunctionList fl = DMPlexGenerateList;
232*3a074057SBarry Smith     while (fl->next) fl = fl->next;
233*3a074057SBarry Smith     fl->next = entry;
234*3a074057SBarry Smith   }
235*3a074057SBarry Smith   PetscFunctionReturn(0);
236*3a074057SBarry Smith }
237*3a074057SBarry Smith 
238*3a074057SBarry Smith extern PetscBool DMPlexGenerateRegisterAllCalled;
239*3a074057SBarry Smith 
240*3a074057SBarry Smith PetscErrorCode  DMPlexGenerateRegisterDestroy(void)
241*3a074057SBarry Smith {
242*3a074057SBarry Smith   PetscFunctionList next,fl;
243*3a074057SBarry Smith   PetscErrorCode    ierr;
244*3a074057SBarry Smith 
245*3a074057SBarry Smith   PetscFunctionBegin;
246*3a074057SBarry Smith   next = fl =  DMPlexGenerateList;
247*3a074057SBarry Smith     while (next) {
248*3a074057SBarry Smith     next = fl ? fl->next : NULL;
249*3a074057SBarry Smith     if (fl) {ierr = PetscFree(fl->name);CHKERRQ(ierr);}
250*3a074057SBarry Smith     ierr = PetscFree(fl);CHKERRQ(ierr);
251*3a074057SBarry Smith     fl = next;
252*3a074057SBarry Smith   }
253*3a074057SBarry Smith   DMPlexGenerateList              = NULL;
254*3a074057SBarry Smith   DMPlexGenerateRegisterAllCalled = PETSC_FALSE;
255d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
256d9deefdfSMatthew G. Knepley }
257