Actual source code: plexgmsh.c
petsc-dev 2014-02-02
1: #define PETSCDM_DLL
2: #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/
6: /*@
7: DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file.
9: Collective on comm
11: Input Parameters:
12: + comm - The MPI communicator
13: . viewer - The Viewer associated with a Gmsh file
14: - interpolate - Create faces and edges in the mesh
16: Output Parameter:
17: . dm - The DM object representing the mesh
19: Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format
21: Level: beginner
23: .keywords: mesh,Gmsh
24: .seealso: DMPLEX, DMCreate()
25: @*/
26: PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
27: {
28: FILE *fd;
29: PetscSection coordSection;
30: Vec coordinates;
31: PetscScalar *coords, *coordsIn = NULL;
32: PetscInt dim = 0, coordSize, c, v, d;
33: int numVertices = 0, numCells = 0, snum;
34: long fpos = 0;
35: PetscMPIInt num_proc, rank;
36: char line[PETSC_MAX_PATH_LEN];
37: PetscBool match;
41: MPI_Comm_rank(comm, &rank);
42: MPI_Comm_size(comm, &num_proc);
43: DMCreate(comm, dm);
44: DMSetType(*dm, DMPLEX);
45: if (!rank) {
46: PetscBool match;
47: int fileType, dataSize;
49: PetscViewerASCIIGetPointer(viewer, &fd);
50: /* Read header */
51: fgets(line, PETSC_MAX_PATH_LEN, fd);
52: PetscStrncmp(line, "$MeshFormat\n", PETSC_MAX_PATH_LEN, &match);
53: if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
54: snum = fscanf(fd, "2.2 %d %d\n", &fileType, &dataSize);CHKERRQ(snum != 2);
55: if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
56: if (dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data size %d is not valid for a Gmsh file", dataSize);
57: fgets(line, PETSC_MAX_PATH_LEN, fd);
58: PetscStrncmp(line, "$EndMeshFormat\n", PETSC_MAX_PATH_LEN, &match);
59: if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
60: /* Read vertices */
61: fgets(line, PETSC_MAX_PATH_LEN, fd);
62: PetscStrncmp(line, "$Nodes\n", PETSC_MAX_PATH_LEN, &match);
63: if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
64: snum = fscanf(fd, "%d\n", &numVertices);CHKERRQ(snum != 1);
65: PetscMalloc(numVertices*3 * sizeof(PetscScalar), &coordsIn);
66: for (v = 0; v < numVertices; ++v) {
67: double x, y, z;
68: int i;
70: snum = fscanf(fd, "%d %lg %lg %lg\n", &i, &x, &y, &z);CHKERRQ(snum != 4);
71: coordsIn[v*3+0] = x; coordsIn[v*3+1] = y; coordsIn[v*3+2] = z;
72: if (i != v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1);
73: }
74: fgets(line, PETSC_MAX_PATH_LEN, fd);
75: PetscStrncmp(line, "$EndNodes\n", PETSC_MAX_PATH_LEN, &match);
76: if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
77: /* Read cells */
78: fgets(line, PETSC_MAX_PATH_LEN, fd);
79: PetscStrncmp(line, "$Elements\n", PETSC_MAX_PATH_LEN, &match);
80: if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
81: snum = fscanf(fd, "%d\n", &numCells);CHKERRQ(snum != 1);
82: }
83: DMPlexSetChart(*dm, 0, numCells+numVertices);
84: if (!rank) {
85: fpos = ftell(fd);
86: for (c = 0; c < numCells; ++c) {
87: PetscInt numCorners, t;
88: int cone[8], i, cellType, numTags, tag;
90: snum = fscanf(fd, "%d %d %d", &i, &cellType, &numTags);CHKERRQ(snum != 3);
91: if (numTags) for (t = 0; t < numTags; ++t) {snum = fscanf(fd, "%d", &tag);CHKERRQ(snum != 1);}
92: switch (cellType) {
93: case 1: /* 2-node line */
94: dim = 1;
95: numCorners = 2;
96: snum = fscanf(fd, "%d %d\n", &cone[0], &cone[1]);CHKERRQ(snum != numCorners);
97: break;
98: case 2: /* 3-node triangle */
99: dim = 2;
100: numCorners = 3;
101: snum = fscanf(fd, "%d %d %d\n", &cone[0], &cone[1], &cone[2]);CHKERRQ(snum != numCorners);
102: break;
103: case 3: /* 4-node quadrangle */
104: dim = 2;
105: numCorners = 4;
106: snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners);
107: break;
108: case 4: /* 4-node tetrahedron */
109: dim = 3;
110: numCorners = 4;
111: snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners);
112: break;
113: case 5: /* 8-node hexahedron */
114: dim = 3;
115: numCorners = 8;
116: snum = fscanf(fd, "%d %d %d %d %d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3], &cone[4], &cone[5], &cone[6], &cone[7]);CHKERRQ(snum != numCorners);
117: break;
118: default:
119: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
120: }
121: DMPlexSetConeSize(*dm, c, numCorners);
122: if (i != c+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell number %d should be %d", i, c+1);
123: }
124: }
125: DMSetUp(*dm);
126: if (!rank) {
127: fseek(fd, fpos, SEEK_SET);
128: for (c = 0; c < numCells; ++c) {
129: PetscInt pcone[8], numCorners, corner, t;
130: int cone[8], i, cellType, numTags, tag;
132: snum = fscanf(fd, "%d %d %d", &i, &cellType, &numTags);CHKERRQ(snum != 3);
133: if (numTags) for (t = 0; t < numTags; ++t) {snum = fscanf(fd, "%d", &tag);CHKERRQ(snum != 1);}
134: switch (cellType) {
135: case 1: /* 2-node line */
136: dim = 1;
137: numCorners = 2;
138: snum = fscanf(fd, "%d %d\n", &cone[0], &cone[1]);CHKERRQ(snum != numCorners);
139: break;
140: case 2: /* 3-node triangle */
141: dim = 2;
142: numCorners = 3;
143: snum = fscanf(fd, "%d %d %d\n", &cone[0], &cone[1], &cone[2]);CHKERRQ(snum != numCorners);
144: break;
145: case 3: /* 4-node quadrangle */
146: dim = 2;
147: numCorners = 4;
148: snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners);
149: break;
150: case 4: /* 4-node tetrahedron */
151: dim = 3;
152: numCorners = 4;
153: snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners);
154: DMPlexInvertCell(dim, numCorners, cone);
155: break;
156: case 5: /* 8-node hexahedron */
157: dim = 3;
158: numCorners = 8;
159: snum = fscanf(fd, "%d %d %d %d %d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3], &cone[4], &cone[5], &cone[6], &cone[7]);CHKERRQ(snum != numCorners);
160: DMPlexInvertCell(dim, numCorners, cone);
161: break;
162: default:
163: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
164: }
165: for (corner = 0; corner < numCorners; ++corner) pcone[corner] = cone[corner] + numCells-1;
166: DMPlexSetCone(*dm, c, pcone);
167: if (i != c+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell number %d should be %d", i, c+1);
168: }
169: fgets(line, PETSC_MAX_PATH_LEN, fd);
170: PetscStrncmp(line, "$EndElements\n", PETSC_MAX_PATH_LEN, &match);
171: if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
172: }
173: MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);
174: DMPlexSetDimension(*dm, dim);
175: DMPlexSymmetrize(*dm);
176: DMPlexStratify(*dm);
177: if (interpolate) {
178: DM idm;
180: DMPlexInterpolate(*dm, &idm);
181: DMDestroy(dm);
182: *dm = idm;
183: }
184: /* Read coordinates */
185: DMGetCoordinateSection(*dm, &coordSection);
186: PetscSectionSetNumFields(coordSection, 1);
187: PetscSectionSetFieldComponents(coordSection, 0, dim);
188: PetscSectionSetChart(coordSection, numCells, numCells + numVertices);
189: for (v = numCells; v < numCells+numVertices; ++v) {
190: PetscSectionSetDof(coordSection, v, dim);
191: PetscSectionSetFieldDof(coordSection, v, 0, dim);
192: }
193: PetscSectionSetUp(coordSection);
194: PetscSectionGetStorageSize(coordSection, &coordSize);
195: VecCreate(comm, &coordinates);
196: PetscObjectSetName((PetscObject) coordinates, "coordinates");
197: VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);
198: VecSetType(coordinates, VECSTANDARD);
199: VecGetArray(coordinates, &coords);
200: if (!rank) {
201: for (v = 0; v < numVertices; ++v) {
202: for (d = 0; d < dim; ++d) {
203: coords[v*dim+d] = coordsIn[v*3+d];
204: }
205: }
206: }
207: VecRestoreArray(coordinates, &coords);
208: PetscFree(coordsIn);
209: DMSetCoordinatesLocal(*dm, coordinates);
210: VecDestroy(&coordinates);
211: return(0);
212: }