xref: /petsc/src/dm/impls/plex/plexgenerate.c (revision 98921bda46e76d7aaed9e0138c5ff9d0ce93f355)
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 
7294ef8ddeSSatish Balay /*@C
73d9deefdfSMatthew G. Knepley   DMPlexTriangleSetOptions - Set the options used for the Triangle mesh generator
74d9deefdfSMatthew G. Knepley 
75d9deefdfSMatthew G. Knepley   Not Collective
76d9deefdfSMatthew G. Knepley 
77d9deefdfSMatthew G. Knepley   Inputs Parameters:
78d9deefdfSMatthew G. Knepley + dm - The DMPlex object
79d9deefdfSMatthew G. Knepley - opts - The command line options
80d9deefdfSMatthew G. Knepley 
81d9deefdfSMatthew G. Knepley   Level: developer
82d9deefdfSMatthew G. Knepley 
83d9deefdfSMatthew G. Knepley .seealso: DMPlexTetgenSetOptions(), DMPlexGenerate()
84d9deefdfSMatthew G. Knepley @*/
85d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexTriangleSetOptions(DM dm, const char *opts)
86d9deefdfSMatthew G. Knepley {
87d9deefdfSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex*) dm->data;
88d9deefdfSMatthew G. Knepley   PetscErrorCode ierr;
89d9deefdfSMatthew G. Knepley 
90d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
91d9deefdfSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
92d9deefdfSMatthew G. Knepley   PetscValidPointer(opts, 2);
93606ac3a5SMatthew G. Knepley   ierr = PetscFree(mesh->triangleOpts);CHKERRQ(ierr);
94606ac3a5SMatthew G. Knepley   ierr = PetscStrallocpy(opts, &mesh->triangleOpts);CHKERRQ(ierr);
95d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
96d9deefdfSMatthew G. Knepley }
97d9deefdfSMatthew G. Knepley 
9894ef8ddeSSatish Balay /*@C
99d9deefdfSMatthew G. Knepley   DMPlexTetgenSetOptions - Set the options used for the Tetgen mesh generator
100d9deefdfSMatthew G. Knepley 
101d9deefdfSMatthew G. Knepley   Not Collective
102d9deefdfSMatthew G. Knepley 
103d9deefdfSMatthew G. Knepley   Inputs Parameters:
104d9deefdfSMatthew G. Knepley + dm - The DMPlex object
105d9deefdfSMatthew G. Knepley - opts - The command line options
106d9deefdfSMatthew G. Knepley 
107d9deefdfSMatthew G. Knepley   Level: developer
108d9deefdfSMatthew G. Knepley 
109d9deefdfSMatthew G. Knepley .seealso: DMPlexTriangleSetOptions(), DMPlexGenerate()
110d9deefdfSMatthew G. Knepley @*/
111d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexTetgenSetOptions(DM dm, const char *opts)
112d9deefdfSMatthew G. Knepley {
113d9deefdfSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex*) dm->data;
114d9deefdfSMatthew G. Knepley   PetscErrorCode ierr;
115d9deefdfSMatthew G. Knepley 
116d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
117d9deefdfSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
118d9deefdfSMatthew G. Knepley   PetscValidPointer(opts, 2);
119606ac3a5SMatthew G. Knepley   ierr = PetscFree(mesh->tetgenOpts);CHKERRQ(ierr);
120606ac3a5SMatthew G. Knepley   ierr = PetscStrallocpy(opts, &mesh->tetgenOpts);CHKERRQ(ierr);
121d9deefdfSMatthew G. Knepley   PetscFunctionReturn(0);
122d9deefdfSMatthew G. Knepley }
123d9deefdfSMatthew G. Knepley 
124d9deefdfSMatthew G. Knepley /*@C
125d9deefdfSMatthew G. Knepley   DMPlexGenerate - Generates a mesh.
126d9deefdfSMatthew G. Knepley 
127d9deefdfSMatthew G. Knepley   Not Collective
128d9deefdfSMatthew G. Knepley 
129d9deefdfSMatthew G. Knepley   Input Parameters:
130d9deefdfSMatthew G. Knepley + boundary - The DMPlex boundary object
131d9deefdfSMatthew G. Knepley . name - The mesh generation package name
132d9deefdfSMatthew G. Knepley - interpolate - Flag to create intermediate mesh elements
133d9deefdfSMatthew G. Knepley 
134d9deefdfSMatthew G. Knepley   Output Parameter:
135d9deefdfSMatthew G. Knepley . mesh - The DMPlex object
136d9deefdfSMatthew G. Knepley 
1372c4dacf2SBarry Smith   Options Database:
1382c4dacf2SBarry Smith +  -dm_plex_generate <name> - package to generate mesh, for example, triangle, ctetgen or tetgen
139c0517cd5SMatthew G. Knepley -  -dm_generator <name> - package to generate mesh, for example, triangle, ctetgen or tetgen
1402c4dacf2SBarry Smith 
141d9deefdfSMatthew G. Knepley   Level: intermediate
142d9deefdfSMatthew G. Knepley 
143d9deefdfSMatthew G. Knepley .seealso: DMPlexCreate(), DMRefine()
144d9deefdfSMatthew G. Knepley @*/
145d9deefdfSMatthew G. Knepley PetscErrorCode DMPlexGenerate(DM boundary, const char name[], PetscBool interpolate, DM *mesh)
146d9deefdfSMatthew G. Knepley {
147c0517cd5SMatthew G. Knepley   DMGeneratorFunctionList fl;
1483f2a96e3SMatthew G. Knepley   char                    genname[PETSC_MAX_PATH_LEN];
1493f2a96e3SMatthew G. Knepley   const char             *suggestions;
150d9deefdfSMatthew G. Knepley   PetscInt                dim;
1513a074057SBarry Smith   PetscBool               flg;
152d9deefdfSMatthew G. Knepley   PetscErrorCode          ierr;
153d9deefdfSMatthew G. Knepley 
154d9deefdfSMatthew G. Knepley   PetscFunctionBegin;
155d9deefdfSMatthew G. Knepley   PetscValidHeaderSpecific(boundary, DM_CLASSID, 1);
156064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveBool(boundary, interpolate, 3);
157d9deefdfSMatthew G. Knepley   ierr = DMGetDimension(boundary, &dim);CHKERRQ(ierr);
158c0517cd5SMatthew G. Knepley   ierr = PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_generator", genname, sizeof(genname), &flg);CHKERRQ(ierr);
159d9deefdfSMatthew G. Knepley   if (flg) name = genname;
1602c4dacf2SBarry Smith   else {
161589a23caSBarry Smith     ierr = PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_plex_generate", genname, sizeof(genname), &flg);CHKERRQ(ierr);
1622c4dacf2SBarry Smith     if (flg) name = genname;
1632c4dacf2SBarry Smith   }
1643a074057SBarry Smith 
165c0517cd5SMatthew G. Knepley   fl = DMGenerateList;
166d9deefdfSMatthew G. Knepley   if (name) {
1673a074057SBarry Smith     while (fl) {
1683a074057SBarry Smith       ierr = PetscStrcmp(fl->name,name,&flg);CHKERRQ(ierr);
1693a074057SBarry Smith       if (flg) {
170367003a6SStefano Zampini         ierr = (*fl->generate)(boundary,interpolate,mesh);CHKERRQ(ierr);
1713a074057SBarry Smith         PetscFunctionReturn(0);
172d9deefdfSMatthew G. Knepley       }
1733a074057SBarry Smith       fl = fl->next;
174d9deefdfSMatthew G. Knepley     }
175*98921bdaSJacob Faibussowitsch     SETERRQ(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);
1763a074057SBarry Smith   } else {
1773a074057SBarry Smith     while (fl) {
1783a074057SBarry Smith       if (boundary->dim == fl->dim) {
179367003a6SStefano Zampini         ierr = (*fl->generate)(boundary,interpolate,mesh);CHKERRQ(ierr);
1803a074057SBarry Smith         PetscFunctionReturn(0);
1813a074057SBarry Smith       }
1823a074057SBarry Smith       fl = fl->next;
1833a074057SBarry Smith     }
1842c4dacf2SBarry Smith     suggestions = "";
1852c4dacf2SBarry Smith     if (boundary->dim+1 == 2) suggestions = " You may need to add --download-triangle to your ./configure options";
1862c4dacf2SBarry Smith     else if (boundary->dim+1 == 3) suggestions = " You may need to add --download-ctetgen or --download-tetgen in your ./configure options";
187*98921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid generator of dimension %D registered%s",boundary->dim+1,suggestions);
1883a074057SBarry Smith   }
1893a074057SBarry Smith }
190