Actual source code: plexgenerate.c
1: #include <petsc/private/dmpleximpl.h>
3: /*@C
4: DMPlexInvertCell - Flips cell orientations since Plex stores some of them internally with outward normals.
6: Input Parameters:
7: + cellType - The cell type
8: - cone - The incoming cone
10: Output Parameter:
11: . cone - The inverted cone (in-place)
13: Level: developer
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)
27: switch (cellType) {
28: case DM_POLYTOPE_POINT: break;
29: case DM_POLYTOPE_SEGMENT: break;
30: case DM_POLYTOPE_POINT_PRISM_TENSOR: break;
31: case DM_POLYTOPE_TRIANGLE: break;
32: case DM_POLYTOPE_QUADRILATERAL: break;
33: case DM_POLYTOPE_SEG_PRISM_TENSOR: SWAPCONE(cone,2,3); break;
34: case DM_POLYTOPE_TETRAHEDRON: SWAPCONE(cone,0,1); break;
35: case DM_POLYTOPE_HEXAHEDRON: SWAPCONE(cone,1,3); break;
36: case DM_POLYTOPE_TRI_PRISM: SWAPCONE(cone,1,2); break;
37: case DM_POLYTOPE_TRI_PRISM_TENSOR: break;
38: case DM_POLYTOPE_QUAD_PRISM_TENSOR: break;
39: default: break;
40: }
41: return 0;
42: #undef SWAPCONE
43: }
45: /*@C
46: DMPlexReorderCell - Flips cell orientations since Plex stores some of them internally with outward normals.
48: Input Parameters:
49: + dm - The DMPlex object
50: . cell - The cell
51: - cone - The incoming cone
53: Output Parameter:
54: . cone - The reordered cone (in-place)
56: Level: developer
58: .seealso: DMPlexGenerate()
59: @*/
60: PetscErrorCode DMPlexReorderCell(DM dm, PetscInt cell, PetscInt cone[])
61: {
62: DMPolytopeType cellType;
64: DMPlexGetCellType(dm, cell, &cellType);
65: DMPlexInvertCell(cellType, cone);
66: return 0;
67: }
69: /*@C
70: DMPlexTriangleSetOptions - Set the options used for the Triangle mesh generator
72: Not Collective
74: Inputs Parameters:
75: + dm - The DMPlex object
76: - opts - The command line options
78: Level: developer
80: .seealso: DMPlexTetgenSetOptions(), DMPlexGenerate()
81: @*/
82: PetscErrorCode DMPlexTriangleSetOptions(DM dm, const char *opts)
83: {
84: DM_Plex *mesh = (DM_Plex*) dm->data;
88: PetscFree(mesh->triangleOpts);
89: PetscStrallocpy(opts, &mesh->triangleOpts);
90: return 0;
91: }
93: /*@C
94: DMPlexTetgenSetOptions - Set the options used for the Tetgen mesh generator
96: Not Collective
98: Inputs Parameters:
99: + dm - The DMPlex object
100: - opts - The command line options
102: Level: developer
104: .seealso: DMPlexTriangleSetOptions(), DMPlexGenerate()
105: @*/
106: PetscErrorCode DMPlexTetgenSetOptions(DM dm, const char *opts)
107: {
108: DM_Plex *mesh = (DM_Plex*) dm->data;
112: PetscFree(mesh->tetgenOpts);
113: PetscStrallocpy(opts, &mesh->tetgenOpts);
114: return 0;
115: }
117: /*@C
118: DMPlexGenerate - Generates a mesh.
120: Not Collective
122: Input Parameters:
123: + boundary - The DMPlex boundary object
124: . name - The mesh generation package name
125: - interpolate - Flag to create intermediate mesh elements
127: Output Parameter:
128: . mesh - The DMPlex object
130: Options Database:
131: + -dm_plex_generate <name> - package to generate mesh, for example, triangle, ctetgen or tetgen
132: - -dm_generator <name> - package to generate mesh, for example, triangle, ctetgen or tetgen
134: Level: intermediate
136: .seealso: DMPlexCreate(), DMRefine()
137: @*/
138: PetscErrorCode DMPlexGenerate(DM boundary, const char name[], PetscBool interpolate, DM *mesh)
139: {
140: DMGeneratorFunctionList fl;
141: char genname[PETSC_MAX_PATH_LEN];
142: const char *suggestions;
143: PetscInt dim;
144: PetscBool flg;
148: DMGetDimension(boundary, &dim);
149: PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_generator", genname, sizeof(genname), &flg);
150: if (flg) name = genname;
151: else {
152: PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_plex_generate", genname, sizeof(genname), &flg);
153: if (flg) name = genname;
154: }
156: fl = DMGenerateList;
157: if (name) {
158: while (fl) {
159: PetscStrcmp(fl->name,name,&flg);
160: if (flg) {
161: (*fl->generate)(boundary,interpolate,mesh);
162: return 0;
163: }
164: fl = fl->next;
165: }
166: 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);
167: } else {
168: while (fl) {
169: if (boundary->dim == fl->dim) {
170: (*fl->generate)(boundary,interpolate,mesh);
171: return 0;
172: }
173: fl = fl->next;
174: }
175: suggestions = "";
176: if (boundary->dim+1 == 2) suggestions = " You may need to add --download-triangle to your ./configure options";
177: else if (boundary->dim+1 == 3) suggestions = " You may need to add --download-ctetgen or --download-tetgen in your ./configure options";
178: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid generator of dimension %D registered%s",boundary->dim+1,suggestions);
179: }
180: }