xref: /petsc/src/dm/impls/plex/plexgenerate.c (revision 589a23caa660d2a5f330cc8d1ed213e9cfaf51a7)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2d9deefdfSMatthew G. Knepley 
3d9deefdfSMatthew G. Knepley /*@C
496ca5757SLisandro Dalcin   DMPlexInvertCell - Flips cell orientations since Plex stores some of them internally with outward normals.
5d9deefdfSMatthew G. Knepley 
6d9deefdfSMatthew G. Knepley   Input Parameters:
796ca5757SLisandro Dalcin + cellType - The cell type
8d9deefdfSMatthew G. Knepley - cone - The incoming cone
9d9deefdfSMatthew G. Knepley 
10d9deefdfSMatthew G. Knepley   Output Parameter:
11d9deefdfSMatthew G. Knepley . cone - The inverted cone (in-place)
12d9deefdfSMatthew G. Knepley 
13d9deefdfSMatthew G. Knepley   Level: developer
14d9deefdfSMatthew G. Knepley 
15d9deefdfSMatthew G. Knepley .seealso: DMPlexGenerate()
16d9deefdfSMatthew G. Knepley @*/
1796ca5757SLisandro Dalcin PetscErrorCode DMPlexInvertCell(DMPolytopeType cellType, PetscInt cone[])
18d9deefdfSMatthew G. Knepley {
1996ca5757SLisandro Dalcin #define SWAPCONE(cone,i,j)  \
2096ca5757SLisandro Dalcin   do {                      \
2196ca5757SLisandro Dalcin     PetscInt _cone_tmp;     \
2296ca5757SLisandro Dalcin     _cone_tmp = (cone)[i];  \
2396ca5757SLisandro Dalcin     (cone)[i] = (cone)[j];  \
2496ca5757SLisandro Dalcin     (cone)[j] = _cone_tmp;  \
2596ca5757SLisandro Dalcin   } while (0)
26d9deefdfSMatthew G. Knepley 
27d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
2896ca5757SLisandro Dalcin   switch (cellType) {
2996ca5757SLisandro Dalcin   case DM_POLYTOPE_POINT:              break;
3096ca5757SLisandro Dalcin   case DM_POLYTOPE_SEGMENT:            break;
3196ca5757SLisandro Dalcin   case DM_POLYTOPE_POINT_PRISM_TENSOR: break;
3296ca5757SLisandro Dalcin   case DM_POLYTOPE_TRIANGLE:           break;
3396ca5757SLisandro Dalcin   case DM_POLYTOPE_QUADRILATERAL:      break;
3496ca5757SLisandro Dalcin   case DM_POLYTOPE_SEG_PRISM_TENSOR:   SWAPCONE(cone,2,3); break;
3596ca5757SLisandro Dalcin   case DM_POLYTOPE_TETRAHEDRON:        SWAPCONE(cone,0,1); break;
3696ca5757SLisandro Dalcin   case DM_POLYTOPE_HEXAHEDRON:         SWAPCONE(cone,1,3); break;
3796ca5757SLisandro Dalcin   case DM_POLYTOPE_TRI_PRISM:          SWAPCONE(cone,1,2); break;
3896ca5757SLisandro Dalcin   case DM_POLYTOPE_TRI_PRISM_TENSOR:   break;
3996ca5757SLisandro Dalcin   case DM_POLYTOPE_QUAD_PRISM_TENSOR:  break;
40d9deefdfSMatthew G. Knepley   default: break;
41d9deefdfSMatthew G. Knepley   }
42d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
4396ca5757SLisandro Dalcin #undef SWAPCONE
4496ca5757SLisandro Dalcin }
4596ca5757SLisandro Dalcin 
4696ca5757SLisandro Dalcin /*@C
4796ca5757SLisandro Dalcin   DMPlexReorderCell - Flips cell orientations since Plex stores some of them internally with outward normals.
4896ca5757SLisandro Dalcin 
4996ca5757SLisandro Dalcin   Input Parameters:
5096ca5757SLisandro Dalcin + dm - The DMPlex object
5196ca5757SLisandro Dalcin . cell - The cell
5296ca5757SLisandro Dalcin - cone - The incoming cone
5396ca5757SLisandro Dalcin 
5496ca5757SLisandro Dalcin   Output Parameter:
5596ca5757SLisandro Dalcin . cone - The reordered cone (in-place)
5696ca5757SLisandro Dalcin 
5796ca5757SLisandro Dalcin   Level: developer
5896ca5757SLisandro Dalcin 
5996ca5757SLisandro Dalcin .seealso: DMPlexGenerate()
6096ca5757SLisandro Dalcin @*/
6196ca5757SLisandro Dalcin PetscErrorCode DMPlexReorderCell(DM dm, PetscInt cell, PetscInt cone[])
6296ca5757SLisandro Dalcin {
6396ca5757SLisandro Dalcin   DMPolytopeType cellType;
6496ca5757SLisandro Dalcin   PetscErrorCode ierr;
6596ca5757SLisandro Dalcin 
6696ca5757SLisandro Dalcin   PetscFunctionBegin;
6796ca5757SLisandro Dalcin   ierr = DMPlexGetCellType(dm, cell, &cellType);CHKERRQ(ierr);
6896ca5757SLisandro Dalcin   ierr = DMPlexInvertCell(cellType, cone);CHKERRQ(ierr);
6996ca5757SLisandro Dalcin   PetscFunctionReturn(0);
70d9deefdfSMatthew G. Knepley }
71d9deefdfSMatthew G. Knepley 
72d9deefdfSMatthew G. Knepley 
7394ef8ddeSSatish Balay /*@C
74d9deefdfSMatthew G. Knepley   DMPlexTriangleSetOptions - Set the options used for the Triangle mesh generator
75d9deefdfSMatthew G. Knepley 
76d9deefdfSMatthew G. Knepley   Not Collective
77d9deefdfSMatthew G. Knepley 
78d9deefdfSMatthew G. Knepley   Inputs Parameters:
79d9deefdfSMatthew G. Knepley + dm - The DMPlex object
80d9deefdfSMatthew G. Knepley - opts - The command line options
81d9deefdfSMatthew G. Knepley 
82d9deefdfSMatthew G. Knepley   Level: developer
83d9deefdfSMatthew G. Knepley 
84d9deefdfSMatthew G. Knepley .seealso: DMPlexTetgenSetOptions(), DMPlexGenerate()
85d9deefdfSMatthew G. Knepley @*/
86d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexTriangleSetOptions(DM dm, const char *opts)
87d9deefdfSMatthew G. Knepley {
88d9deefdfSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex*) dm->data;
89d9deefdfSMatthew G. Knepley   PetscErrorCode ierr;
90d9deefdfSMatthew G. Knepley 
91d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
92d9deefdfSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
93d9deefdfSMatthew G. Knepley   PetscValidPointer(opts, 2);
94606ac3a5SMatthew G. Knepley   ierr = PetscFree(mesh->triangleOpts);CHKERRQ(ierr);
95606ac3a5SMatthew G. Knepley   ierr = PetscStrallocpy(opts, &mesh->triangleOpts);CHKERRQ(ierr);
96d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
97d9deefdfSMatthew G. Knepley }
98d9deefdfSMatthew G. Knepley 
9994ef8ddeSSatish Balay /*@C
100d9deefdfSMatthew G. Knepley   DMPlexTetgenSetOptions - Set the options used for the Tetgen mesh generator
101d9deefdfSMatthew G. Knepley 
102d9deefdfSMatthew G. Knepley   Not Collective
103d9deefdfSMatthew G. Knepley 
104d9deefdfSMatthew G. Knepley   Inputs Parameters:
105d9deefdfSMatthew G. Knepley + dm - The DMPlex object
106d9deefdfSMatthew G. Knepley - opts - The command line options
107d9deefdfSMatthew G. Knepley 
108d9deefdfSMatthew G. Knepley   Level: developer
109d9deefdfSMatthew G. Knepley 
110d9deefdfSMatthew G. Knepley .seealso: DMPlexTriangleSetOptions(), DMPlexGenerate()
111d9deefdfSMatthew G. Knepley @*/
112d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexTetgenSetOptions(DM dm, const char *opts)
113d9deefdfSMatthew G. Knepley {
114d9deefdfSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex*) dm->data;
115d9deefdfSMatthew G. Knepley   PetscErrorCode ierr;
116d9deefdfSMatthew G. Knepley 
117d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
118d9deefdfSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
119d9deefdfSMatthew G. Knepley   PetscValidPointer(opts, 2);
120606ac3a5SMatthew G. Knepley   ierr = PetscFree(mesh->tetgenOpts);CHKERRQ(ierr);
121606ac3a5SMatthew G. Knepley   ierr = PetscStrallocpy(opts, &mesh->tetgenOpts);CHKERRQ(ierr);
122d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
123d9deefdfSMatthew G. Knepley }
124d9deefdfSMatthew G. Knepley 
1253a074057SBarry Smith /*
1263a074057SBarry Smith    Contains the list of registered DMPlexGenerators routines
1273a074057SBarry Smith */
1283a074057SBarry Smith PetscFunctionList DMPlexGenerateList = NULL;
129d9deefdfSMatthew G. Knepley 
1303a074057SBarry Smith struct _n_PetscFunctionList {
1313a074057SBarry Smith   PetscErrorCode    (*generate)(DM, PetscBool, DM*);
1323a074057SBarry Smith   PetscErrorCode    (*refine)(DM,double*, DM*);
1333a074057SBarry Smith   char              *name;               /* string to identify routine */
1343a074057SBarry Smith   PetscInt          dim;
1353a074057SBarry Smith   PetscFunctionList next;                /* next pointer */
1363a074057SBarry Smith };
137d9deefdfSMatthew G. Knepley 
138d9deefdfSMatthew G. Knepley /*@C
139d9deefdfSMatthew G. Knepley   DMPlexGenerate - Generates a mesh.
140d9deefdfSMatthew G. Knepley 
141d9deefdfSMatthew G. Knepley   Not Collective
142d9deefdfSMatthew G. Knepley 
143d9deefdfSMatthew G. Knepley   Input Parameters:
144d9deefdfSMatthew G. Knepley + boundary - The DMPlex boundary object
145d9deefdfSMatthew G. Knepley . name - The mesh generation package name
146d9deefdfSMatthew G. Knepley - interpolate - Flag to create intermediate mesh elements
147d9deefdfSMatthew G. Knepley 
148d9deefdfSMatthew G. Knepley   Output Parameter:
149d9deefdfSMatthew G. Knepley . mesh - The DMPlex object
150d9deefdfSMatthew G. Knepley 
1512c4dacf2SBarry Smith   Options Database:
1522c4dacf2SBarry Smith +  -dm_plex_generate <name> - package to generate mesh, for example, triangle, ctetgen or tetgen
1532c4dacf2SBarry Smith -  -dm_plex_generator <name> - package to generate mesh, for example, triangle, ctetgen or tetgen (deprecated)
1542c4dacf2SBarry Smith 
155d9deefdfSMatthew G. Knepley   Level: intermediate
156d9deefdfSMatthew G. Knepley 
157d9deefdfSMatthew G. Knepley .seealso: DMPlexCreate(), DMRefine()
158d9deefdfSMatthew G. Knepley @*/
159d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexGenerate(DM boundary, const char name[], PetscBool interpolate, DM *mesh)
160d9deefdfSMatthew G. Knepley {
161d9deefdfSMatthew G. Knepley   PetscInt          dim;
162d9deefdfSMatthew G. Knepley   char              genname[1024];
1633a074057SBarry Smith   PetscBool         flg;
164d9deefdfSMatthew G. Knepley   PetscErrorCode    ierr;
1653a074057SBarry Smith   PetscFunctionList fl;
1662c4dacf2SBarry Smith   const char*       suggestions;
167d9deefdfSMatthew G. Knepley 
168d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
169d9deefdfSMatthew G. Knepley   PetscValidHeaderSpecific(boundary, DM_CLASSID, 1);
170d9deefdfSMatthew G. Knepley   PetscValidLogicalCollectiveBool(boundary, interpolate, 2);
171d9deefdfSMatthew G. Knepley   ierr = DMGetDimension(boundary, &dim);CHKERRQ(ierr);
172*589a23caSBarry Smith   ierr = PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_plex_generator", genname, sizeof(genname), &flg);CHKERRQ(ierr);
173d9deefdfSMatthew G. Knepley   if (flg) name = genname;
1742c4dacf2SBarry Smith   else {
175*589a23caSBarry Smith     ierr = PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_plex_generate", genname, sizeof(genname), &flg);CHKERRQ(ierr);
1762c4dacf2SBarry Smith     if (flg) name = genname;
1772c4dacf2SBarry Smith   }
1783a074057SBarry Smith 
1793a074057SBarry Smith   fl = DMPlexGenerateList;
180d9deefdfSMatthew G. Knepley   if (name) {
1813a074057SBarry Smith     while (fl) {
1823a074057SBarry Smith       ierr = PetscStrcmp(fl->name,name,&flg);CHKERRQ(ierr);
1833a074057SBarry Smith       if (flg) {
184367003a6SStefano Zampini         ierr = (*fl->generate)(boundary,interpolate,mesh);CHKERRQ(ierr);
1853a074057SBarry Smith         PetscFunctionReturn(0);
186d9deefdfSMatthew G. Knepley       }
1873a074057SBarry Smith       fl = fl->next;
188d9deefdfSMatthew G. Knepley     }
1892c4dacf2SBarry Smith     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);
1903a074057SBarry Smith   } else {
1913a074057SBarry Smith     while (fl) {
1923a074057SBarry Smith       if (boundary->dim == fl->dim) {
193367003a6SStefano Zampini         ierr = (*fl->generate)(boundary,interpolate,mesh);CHKERRQ(ierr);
1943a074057SBarry Smith         PetscFunctionReturn(0);
1953a074057SBarry Smith       }
1963a074057SBarry Smith       fl = fl->next;
1973a074057SBarry Smith     }
1982c4dacf2SBarry Smith     suggestions = "";
1992c4dacf2SBarry Smith     if (boundary->dim+1 == 2) suggestions = " You may need to add --download-triangle to your ./configure options";
2002c4dacf2SBarry Smith     else if (boundary->dim+1 == 3) suggestions = " You may need to add --download-ctetgen or --download-tetgen in your ./configure options";
2012c4dacf2SBarry Smith     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid generator of dimension %D registered%s",boundary->dim+1,suggestions);
2023a074057SBarry Smith   }
2033a074057SBarry Smith   PetscFunctionReturn(0);
2043a074057SBarry Smith }
2053a074057SBarry Smith 
2063a074057SBarry Smith /*@C
2073a074057SBarry Smith   DMPlexGenerateRegister -  Adds a grid generator to DMPlex
2083a074057SBarry Smith 
2093a074057SBarry Smith    Not Collective
2103a074057SBarry Smith 
2113a074057SBarry Smith    Input Parameters:
212cd921712SStefano Zampini +  name_solver - name of a new user-defined grid generator
2133a074057SBarry Smith .  fnc - generator function
2143a074057SBarry Smith .  rfnc - refinement function
2153a074057SBarry Smith -  dim - dimension of boundary of domain
2163a074057SBarry Smith 
2173a074057SBarry Smith    Notes:
2183a074057SBarry Smith    DMPlexGenerateRegister() may be called multiple times to add several user-defined solvers.
2193a074057SBarry Smith 
2203a074057SBarry Smith    Sample usage:
2213a074057SBarry Smith .vb
222cd921712SStefano Zampini    DMPlexGenerateRegister("my_generator",MyGeneratorCreate,MyGeneratorRefiner,dim);
2233a074057SBarry Smith .ve
2243a074057SBarry Smith 
2253a074057SBarry Smith    Then, your generator can be chosen with the procedural interface via
226cd921712SStefano Zampini $     DMPlexGenerate(dm,"my_generator",...)
2273a074057SBarry Smith    or at runtime via the option
228cd921712SStefano Zampini $     -dm_plex_generator my_generator
2293a074057SBarry Smith 
2303a074057SBarry Smith    Level: advanced
2313a074057SBarry Smith 
232cd921712SStefano Zampini .seealso: DMPlexGenerateRegisterAll(), DMPlexGenerate(), DMPlexGenerateRegisterDestroy()
2333a074057SBarry Smith 
2343a074057SBarry Smith @*/
2353a074057SBarry Smith PetscErrorCode  DMPlexGenerateRegister(const char sname[],PetscErrorCode (*fnc)(DM, PetscBool,DM*), PetscErrorCode (*rfnc)(DM, double*,DM*),PetscInt dim)
2363a074057SBarry Smith {
2373a074057SBarry Smith   PetscErrorCode    ierr;
2383a074057SBarry Smith   PetscFunctionList entry;
2393a074057SBarry Smith 
2403a074057SBarry Smith   PetscFunctionBegin;
2413a074057SBarry Smith   ierr            = PetscNew(&entry);CHKERRQ(ierr);
2423a074057SBarry Smith   ierr            = PetscStrallocpy(sname,&entry->name);CHKERRQ(ierr);
2433a074057SBarry Smith   entry->generate = fnc;
2443a074057SBarry Smith   entry->refine   = rfnc;
2453a074057SBarry Smith   entry->dim      = dim;
2463a074057SBarry Smith   entry->next     = NULL;
2473a074057SBarry Smith   if (!DMPlexGenerateList) DMPlexGenerateList = entry;
2483a074057SBarry Smith   else {
2493a074057SBarry Smith     PetscFunctionList fl = DMPlexGenerateList;
2503a074057SBarry Smith     while (fl->next) fl = fl->next;
2513a074057SBarry Smith     fl->next = entry;
2523a074057SBarry Smith   }
2533a074057SBarry Smith   PetscFunctionReturn(0);
2543a074057SBarry Smith }
2553a074057SBarry Smith 
2563a074057SBarry Smith extern PetscBool DMPlexGenerateRegisterAllCalled;
2573a074057SBarry Smith 
2583a074057SBarry Smith PetscErrorCode  DMPlexGenerateRegisterDestroy(void)
2593a074057SBarry Smith {
2603a074057SBarry Smith   PetscFunctionList next,fl;
2613a074057SBarry Smith   PetscErrorCode    ierr;
2623a074057SBarry Smith 
2633a074057SBarry Smith   PetscFunctionBegin;
2643a074057SBarry Smith   next = fl =  DMPlexGenerateList;
2653a074057SBarry Smith     while (next) {
2663a074057SBarry Smith     next = fl ? fl->next : NULL;
2673a074057SBarry Smith     if (fl) {ierr = PetscFree(fl->name);CHKERRQ(ierr);}
2683a074057SBarry Smith     ierr = PetscFree(fl);CHKERRQ(ierr);
2693a074057SBarry Smith     fl = next;
2703a074057SBarry Smith   }
2713a074057SBarry Smith   DMPlexGenerateList              = NULL;
2723a074057SBarry Smith   DMPlexGenerateRegisterAllCalled = PETSC_FALSE;
273d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
274d9deefdfSMatthew G. Knepley }
275