xref: /petsc/src/dm/impls/plex/plexgenerate.c (revision 367003a68d4e38a62ba2a0620cd4e6f42aa373fd)
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 
1163a074057SBarry Smith /*
1173a074057SBarry Smith    Contains the list of registered DMPlexGenerators routines
1183a074057SBarry Smith */
1193a074057SBarry Smith PetscFunctionList DMPlexGenerateList = NULL;
120d9deefdfSMatthew G. Knepley 
1213a074057SBarry Smith struct _n_PetscFunctionList {
1223a074057SBarry Smith   PetscErrorCode    (*generate)(DM, PetscBool, DM*);
1233a074057SBarry Smith   PetscErrorCode    (*refine)(DM,double*, DM*);
1243a074057SBarry Smith   char              *name;               /* string to identify routine */
1253a074057SBarry Smith   PetscInt          dim;
1263a074057SBarry Smith   PetscFunctionList next;                /* next pointer */
1273a074057SBarry 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];
1513a074057SBarry Smith   PetscBool         flg;
152d9deefdfSMatthew G. Knepley   PetscErrorCode    ierr;
1533a074057SBarry 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;
1613a074057SBarry Smith 
1623a074057SBarry Smith   fl = DMPlexGenerateList;
163d9deefdfSMatthew G. Knepley   if (name) {
1643a074057SBarry Smith     while (fl) {
1653a074057SBarry Smith       ierr = PetscStrcmp(fl->name,name,&flg);CHKERRQ(ierr);
1663a074057SBarry Smith       if (flg) {
167*367003a6SStefano Zampini         ierr = (*fl->generate)(boundary,interpolate,mesh);CHKERRQ(ierr);
1683a074057SBarry Smith         PetscFunctionReturn(0);
169d9deefdfSMatthew G. Knepley       }
1703a074057SBarry Smith       fl = fl->next;
171d9deefdfSMatthew G. Knepley     }
172*367003a6SStefano Zampini     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Grid generator %g not registered",name);
1733a074057SBarry Smith   } else {
1743a074057SBarry Smith     while (fl) {
1753a074057SBarry Smith       if (boundary->dim == fl->dim) {
176*367003a6SStefano Zampini         ierr = (*fl->generate)(boundary,interpolate,mesh);CHKERRQ(ierr);
1773a074057SBarry Smith         PetscFunctionReturn(0);
1783a074057SBarry Smith       }
1793a074057SBarry Smith       fl = fl->next;
1803a074057SBarry Smith     }
181*367003a6SStefano Zampini     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid generator of dimension %D registered",boundary->dim);
1823a074057SBarry Smith   }
1833a074057SBarry Smith   PetscFunctionReturn(0);
1843a074057SBarry Smith }
1853a074057SBarry Smith 
1863a074057SBarry Smith /*@C
1873a074057SBarry Smith   DMPlexGenerateRegister -  Adds a grid generator to DMPlex
1883a074057SBarry Smith 
1893a074057SBarry Smith    Not Collective
1903a074057SBarry Smith 
1913a074057SBarry Smith    Input Parameters:
192cd921712SStefano Zampini +  name_solver - name of a new user-defined grid generator
1933a074057SBarry Smith .  fnc - generator function
1943a074057SBarry Smith .  rfnc - refinement function
1953a074057SBarry Smith -  dim - dimension of boundary of domain
1963a074057SBarry Smith 
1973a074057SBarry Smith    Notes:
1983a074057SBarry Smith    DMPlexGenerateRegister() may be called multiple times to add several user-defined solvers.
1993a074057SBarry Smith 
2003a074057SBarry Smith    Sample usage:
2013a074057SBarry Smith .vb
202cd921712SStefano Zampini    DMPlexGenerateRegister("my_generator",MyGeneratorCreate,MyGeneratorRefiner,dim);
2033a074057SBarry Smith .ve
2043a074057SBarry Smith 
2053a074057SBarry Smith    Then, your generator can be chosen with the procedural interface via
206cd921712SStefano Zampini $     DMPlexGenerate(dm,"my_generator",...)
2073a074057SBarry Smith    or at runtime via the option
208cd921712SStefano Zampini $     -dm_plex_generator my_generator
2093a074057SBarry Smith 
2103a074057SBarry Smith    Level: advanced
2113a074057SBarry Smith 
2123a074057SBarry Smith .keywords: DMPlexGenerate, register
2133a074057SBarry Smith 
214cd921712SStefano Zampini .seealso: DMPlexGenerateRegisterAll(), DMPlexGenerate(), DMPlexGenerateRegisterDestroy()
2153a074057SBarry Smith 
2163a074057SBarry Smith @*/
2173a074057SBarry Smith PetscErrorCode  DMPlexGenerateRegister(const char sname[],PetscErrorCode (*fnc)(DM, PetscBool,DM*), PetscErrorCode (*rfnc)(DM, double*,DM*),PetscInt dim)
2183a074057SBarry Smith {
2193a074057SBarry Smith   PetscErrorCode    ierr;
2203a074057SBarry Smith   PetscFunctionList entry;
2213a074057SBarry Smith 
2223a074057SBarry Smith   PetscFunctionBegin;
2233a074057SBarry Smith   ierr            = PetscNew(&entry);CHKERRQ(ierr);
2243a074057SBarry Smith   ierr            = PetscStrallocpy(sname,&entry->name);CHKERRQ(ierr);
2253a074057SBarry Smith   entry->generate = fnc;
2263a074057SBarry Smith   entry->refine   = rfnc;
2273a074057SBarry Smith   entry->dim      = dim;
2283a074057SBarry Smith   entry->next     = NULL;
2293a074057SBarry Smith   if (!DMPlexGenerateList) DMPlexGenerateList = entry;
2303a074057SBarry Smith   else {
2313a074057SBarry Smith     PetscFunctionList fl = DMPlexGenerateList;
2323a074057SBarry Smith     while (fl->next) fl = fl->next;
2333a074057SBarry Smith     fl->next = entry;
2343a074057SBarry Smith   }
2353a074057SBarry Smith   PetscFunctionReturn(0);
2363a074057SBarry Smith }
2373a074057SBarry Smith 
2383a074057SBarry Smith extern PetscBool DMPlexGenerateRegisterAllCalled;
2393a074057SBarry Smith 
2403a074057SBarry Smith PetscErrorCode  DMPlexGenerateRegisterDestroy(void)
2413a074057SBarry Smith {
2423a074057SBarry Smith   PetscFunctionList next,fl;
2433a074057SBarry Smith   PetscErrorCode    ierr;
2443a074057SBarry Smith 
2453a074057SBarry Smith   PetscFunctionBegin;
2463a074057SBarry Smith   next = fl =  DMPlexGenerateList;
2473a074057SBarry Smith     while (next) {
2483a074057SBarry Smith     next = fl ? fl->next : NULL;
2493a074057SBarry Smith     if (fl) {ierr = PetscFree(fl->name);CHKERRQ(ierr);}
2503a074057SBarry Smith     ierr = PetscFree(fl);CHKERRQ(ierr);
2513a074057SBarry Smith     fl = next;
2523a074057SBarry Smith   }
2533a074057SBarry Smith   DMPlexGenerateList              = NULL;
2543a074057SBarry Smith   DMPlexGenerateRegisterAllCalled = PETSC_FALSE;
255d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
256d9deefdfSMatthew G. Knepley }
257