xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision db397164376a084390f1f7c2c738f890136c9c53)
1331830f3SMatthew G. Knepley #define PETSCDM_DLL
2331830f3SMatthew G. Knepley #include <petsc-private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3331830f3SMatthew G. Knepley 
4331830f3SMatthew G. Knepley #undef __FUNCT__
5331830f3SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateGmsh"
6331830f3SMatthew G. Knepley /*@
7331830f3SMatthew G. Knepley   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file.
8331830f3SMatthew G. Knepley 
9331830f3SMatthew G. Knepley   Collective on comm
10331830f3SMatthew G. Knepley 
11331830f3SMatthew G. Knepley   Input Parameters:
12331830f3SMatthew G. Knepley + comm  - The MPI communicator
13331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file
14331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
15331830f3SMatthew G. Knepley 
16331830f3SMatthew G. Knepley   Output Parameter:
17331830f3SMatthew G. Knepley . dm  - The DM object representing the mesh
18331830f3SMatthew G. Knepley 
19331830f3SMatthew G. Knepley   Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format
20331830f3SMatthew G. Knepley 
21331830f3SMatthew G. Knepley   Level: beginner
22331830f3SMatthew G. Knepley 
23331830f3SMatthew G. Knepley .keywords: mesh,Gmsh
24331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate()
25331830f3SMatthew G. Knepley @*/
26331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
27331830f3SMatthew G. Knepley {
28331830f3SMatthew G. Knepley   FILE          *fd;
29331830f3SMatthew G. Knepley   PetscSection   coordSection;
30331830f3SMatthew G. Knepley   Vec            coordinates;
31331830f3SMatthew G. Knepley   PetscScalar   *coords, *coordsIn = NULL;
32*db397164SMichael Lange   PetscInt       dim = 0, tdim = 0, coordSize, c, v, d;
33*db397164SMichael Lange   int            numVertices = 0, numCells = 0, trueNumCells = 0, numFacets = 0, snum;
34261b4668SMatthew G. Knepley   long           fpos = 0;
35331830f3SMatthew G. Knepley   PetscMPIInt    num_proc, rank;
36331830f3SMatthew G. Knepley   char           line[PETSC_MAX_PATH_LEN];
37331830f3SMatthew G. Knepley   PetscBool      match;
38331830f3SMatthew G. Knepley   PetscErrorCode ierr;
39331830f3SMatthew G. Knepley 
40331830f3SMatthew G. Knepley   PetscFunctionBegin;
41331830f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
42331830f3SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
43331830f3SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
44331830f3SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
45331830f3SMatthew G. Knepley   if (!rank) {
46331830f3SMatthew G. Knepley     PetscBool match;
47331830f3SMatthew G. Knepley     int       fileType, dataSize;
48331830f3SMatthew G. Knepley 
49331830f3SMatthew G. Knepley     ierr = PetscViewerASCIIGetPointer(viewer, &fd);CHKERRQ(ierr);
50331830f3SMatthew G. Knepley     /* Read header */
51331830f3SMatthew G. Knepley     fgets(line, PETSC_MAX_PATH_LEN, fd);
52331830f3SMatthew G. Knepley     ierr = PetscStrncmp(line, "$MeshFormat\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
53331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
54331830f3SMatthew G. Knepley     snum = fscanf(fd, "2.2 %d %d\n", &fileType, &dataSize);CHKERRQ(snum != 2);
55331830f3SMatthew G. Knepley     if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
56331830f3SMatthew G. Knepley     if (dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data size %d is not valid for a Gmsh file", dataSize);
57331830f3SMatthew G. Knepley     fgets(line, PETSC_MAX_PATH_LEN, fd);
58331830f3SMatthew G. Knepley     ierr = PetscStrncmp(line, "$EndMeshFormat\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
59331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
60331830f3SMatthew G. Knepley     /* Read vertices */
61331830f3SMatthew G. Knepley     fgets(line, PETSC_MAX_PATH_LEN, fd);
62331830f3SMatthew G. Knepley     ierr = PetscStrncmp(line, "$Nodes\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
63331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
64331830f3SMatthew G. Knepley     snum = fscanf(fd, "%d\n", &numVertices);CHKERRQ(snum != 1);
65331830f3SMatthew G. Knepley     ierr = PetscMalloc(numVertices*3 * sizeof(PetscScalar), &coordsIn);CHKERRQ(ierr);
66331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
67331830f3SMatthew G. Knepley       double x, y, z;
68331830f3SMatthew G. Knepley       int    i;
69331830f3SMatthew G. Knepley 
70331830f3SMatthew G. Knepley       snum = fscanf(fd, "%d %lg %lg %lg\n", &i, &x, &y, &z);CHKERRQ(snum != 4);
71331830f3SMatthew G. Knepley       coordsIn[v*3+0] = x; coordsIn[v*3+1] = y; coordsIn[v*3+2] = z;
72331830f3SMatthew G. Knepley       if (i != v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1);
73331830f3SMatthew G. Knepley     }
74331830f3SMatthew G. Knepley     fgets(line, PETSC_MAX_PATH_LEN, fd);
75331830f3SMatthew G. Knepley     ierr = PetscStrncmp(line, "$EndNodes\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
76331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
77331830f3SMatthew G. Knepley     /* Read cells */
78331830f3SMatthew G. Knepley     fgets(line, PETSC_MAX_PATH_LEN, fd);
79331830f3SMatthew G. Knepley     ierr = PetscStrncmp(line, "$Elements\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
80331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
81331830f3SMatthew G. Knepley     snum = fscanf(fd, "%d\n", &numCells);CHKERRQ(snum != 1);
82331830f3SMatthew G. Knepley   }
83*db397164SMichael Lange 
84331830f3SMatthew G. Knepley   if (!rank) {
85331830f3SMatthew G. Knepley     fpos = ftell(fd);
86*db397164SMichael Lange     /* The Gmsh format disguises facets as elements, so we have to run through all "element" entries
87*db397164SMichael Lange        to get the correct numCells and decide the topological dimension of the mesh */
88*db397164SMichael Lange     trueNumCells = 0;
89*db397164SMichael Lange     for (c = 0; c < numCells; ++c) {
90*db397164SMichael Lange       PetscInt numCorners, t;
91*db397164SMichael Lange       int      cone[8], i, cellType, numTags, tag;
92*db397164SMichael Lange       snum = fscanf(fd, "%d %d %d", &i, &cellType, &numTags);CHKERRQ(snum != 3);
93*db397164SMichael Lange       if (numTags) for (t = 0; t < numTags; ++t) {snum = fscanf(fd, "%d", &tag);CHKERRQ(snum != 1);}
94*db397164SMichael Lange       switch (cellType) {
95*db397164SMichael Lange       case 1: /* 2-node line */
96*db397164SMichael Lange         dim = 1;
97*db397164SMichael Lange         numCorners = 2;
98*db397164SMichael Lange         snum = fscanf(fd, "%d %d\n", &cone[0], &cone[1]);CHKERRQ(snum != numCorners);
99*db397164SMichael Lange         break;
100*db397164SMichael Lange       case 2: /* 3-node triangle */
101*db397164SMichael Lange         dim = 2;
102*db397164SMichael Lange         numCorners = 3;
103*db397164SMichael Lange         snum = fscanf(fd, "%d %d %d\n", &cone[0], &cone[1], &cone[2]);CHKERRQ(snum != numCorners);
104*db397164SMichael Lange         break;
105*db397164SMichael Lange       case 3: /* 4-node quadrangle */
106*db397164SMichael Lange         dim = 2;
107*db397164SMichael Lange         numCorners = 4;
108*db397164SMichael Lange         snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners);
109*db397164SMichael Lange         break;
110*db397164SMichael Lange       case 4: /* 4-node tetrahedron */
111*db397164SMichael Lange         dim  = 3;
112*db397164SMichael Lange         numCorners = 4;
113*db397164SMichael Lange         snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners);
114*db397164SMichael Lange         break;
115*db397164SMichael Lange       case 5: /* 8-node hexahedron */
116*db397164SMichael Lange         dim = 3;
117*db397164SMichael Lange         numCorners = 8;
118*db397164SMichael Lange         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);
119*db397164SMichael Lange         break;
120*db397164SMichael Lange       default:
121*db397164SMichael Lange         SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
122*db397164SMichael Lange       }
123*db397164SMichael Lange       if (dim > tdim) {
124*db397164SMichael Lange         tdim = dim;
125*db397164SMichael Lange         trueNumCells = 0;
126*db397164SMichael Lange       }
127*db397164SMichael Lange       trueNumCells++;
128*db397164SMichael Lange     }
129*db397164SMichael Lange   }
130*db397164SMichael Lange   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
131*db397164SMichael Lange   numFacets = numCells - trueNumCells;
132*db397164SMichael Lange   if (!rank) {
133*db397164SMichael Lange     ierr = fseek(fd, fpos, SEEK_SET);CHKERRQ(ierr);
134331830f3SMatthew G. Knepley     for (c = 0; c < numCells; ++c) {
13521328a2fSMatthew G. Knepley       PetscInt numCorners, t;
13621328a2fSMatthew G. Knepley       int      cone[8], i, cellType, numTags, tag;
137331830f3SMatthew G. Knepley 
138331830f3SMatthew G. Knepley       snum = fscanf(fd, "%d %d %d", &i, &cellType, &numTags);CHKERRQ(snum != 3);
139331830f3SMatthew G. Knepley       if (numTags) for (t = 0; t < numTags; ++t) {snum = fscanf(fd, "%d", &tag);CHKERRQ(snum != 1);}
140331830f3SMatthew G. Knepley       switch (cellType) {
141331830f3SMatthew G. Knepley       case 1: /* 2-node line */
142331830f3SMatthew G. Knepley         dim = 1;
143331830f3SMatthew G. Knepley         numCorners = 2;
144331830f3SMatthew G. Knepley         snum = fscanf(fd, "%d %d\n", &cone[0], &cone[1]);CHKERRQ(snum != numCorners);
145331830f3SMatthew G. Knepley         break;
146331830f3SMatthew G. Knepley       case 2: /* 3-node triangle */
147331830f3SMatthew G. Knepley         dim = 2;
148331830f3SMatthew G. Knepley         numCorners = 3;
149331830f3SMatthew G. Knepley         snum = fscanf(fd, "%d %d %d\n", &cone[0], &cone[1], &cone[2]);CHKERRQ(snum != numCorners);
150331830f3SMatthew G. Knepley         break;
151331830f3SMatthew G. Knepley       case 3: /* 4-node quadrangle */
152331830f3SMatthew G. Knepley         dim = 2;
153331830f3SMatthew G. Knepley         numCorners = 4;
154331830f3SMatthew G. Knepley         snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners);
155331830f3SMatthew G. Knepley         break;
156331830f3SMatthew G. Knepley       case 4: /* 4-node tetrahedron */
157331830f3SMatthew G. Knepley         dim  = 3;
158331830f3SMatthew G. Knepley         numCorners = 4;
159331830f3SMatthew G. Knepley         snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners);
160331830f3SMatthew G. Knepley         break;
161331830f3SMatthew G. Knepley       case 5: /* 8-node hexahedron */
162331830f3SMatthew G. Knepley         dim = 3;
163331830f3SMatthew G. Knepley         numCorners = 8;
164331830f3SMatthew G. Knepley         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);
165331830f3SMatthew G. Knepley         break;
166331830f3SMatthew G. Knepley       default:
167331830f3SMatthew G. Knepley         SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
168331830f3SMatthew G. Knepley       }
169*db397164SMichael Lange       if (dim == tdim) {
170*db397164SMichael Lange         ierr = DMPlexSetConeSize(*dm, c-numFacets, numCorners);CHKERRQ(ierr);
171*db397164SMichael Lange       }
172331830f3SMatthew G. Knepley       if (i != c+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell number %d should be %d", i, c+1);
173331830f3SMatthew G. Knepley     }
174331830f3SMatthew G. Knepley   }
175331830f3SMatthew G. Knepley   ierr = DMSetUp(*dm);CHKERRQ(ierr);
176331830f3SMatthew G. Knepley   if (!rank) {
177331830f3SMatthew G. Knepley     ierr = fseek(fd, fpos, SEEK_SET);CHKERRQ(ierr);
178331830f3SMatthew G. Knepley     for (c = 0; c < numCells; ++c) {
17921328a2fSMatthew G. Knepley       PetscInt pcone[8], numCorners, corner, t;
18021328a2fSMatthew G. Knepley       int      cone[8], i, cellType, numTags, tag;
181331830f3SMatthew G. Knepley 
182331830f3SMatthew G. Knepley       snum = fscanf(fd, "%d %d %d", &i, &cellType, &numTags);CHKERRQ(snum != 3);
183331830f3SMatthew G. Knepley       if (numTags) for (t = 0; t < numTags; ++t) {snum = fscanf(fd, "%d", &tag);CHKERRQ(snum != 1);}
184331830f3SMatthew G. Knepley       switch (cellType) {
185331830f3SMatthew G. Knepley       case 1: /* 2-node line */
186331830f3SMatthew G. Knepley         dim = 1;
187331830f3SMatthew G. Knepley         numCorners = 2;
188331830f3SMatthew G. Knepley         snum = fscanf(fd, "%d %d\n", &cone[0], &cone[1]);CHKERRQ(snum != numCorners);
189331830f3SMatthew G. Knepley         break;
190331830f3SMatthew G. Knepley       case 2: /* 3-node triangle */
191331830f3SMatthew G. Knepley         dim = 2;
192331830f3SMatthew G. Knepley         numCorners = 3;
193331830f3SMatthew G. Knepley         snum = fscanf(fd, "%d %d %d\n", &cone[0], &cone[1], &cone[2]);CHKERRQ(snum != numCorners);
194331830f3SMatthew G. Knepley         break;
195331830f3SMatthew G. Knepley       case 3: /* 4-node quadrangle */
196331830f3SMatthew G. Knepley         dim = 2;
197331830f3SMatthew G. Knepley         numCorners = 4;
198331830f3SMatthew G. Knepley         snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners);
199331830f3SMatthew G. Knepley         break;
200331830f3SMatthew G. Knepley       case 4: /* 4-node tetrahedron */
201331830f3SMatthew G. Knepley         dim  = 3;
202331830f3SMatthew G. Knepley         numCorners = 4;
203331830f3SMatthew G. Knepley         snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners);
2047083e7f6SMatthew G. Knepley         ierr = DMPlexInvertCell(dim, numCorners, cone);CHKERRQ(ierr);
205331830f3SMatthew G. Knepley         break;
206331830f3SMatthew G. Knepley       case 5: /* 8-node hexahedron */
207331830f3SMatthew G. Knepley         dim = 3;
208331830f3SMatthew G. Knepley         numCorners = 8;
209331830f3SMatthew G. Knepley         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);
2107083e7f6SMatthew G. Knepley         ierr = DMPlexInvertCell(dim, numCorners, cone);CHKERRQ(ierr);
211331830f3SMatthew G. Knepley         break;
212331830f3SMatthew G. Knepley       default:
213331830f3SMatthew G. Knepley         SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
214331830f3SMatthew G. Knepley       }
215*db397164SMichael Lange       if (dim == tdim) {
216*db397164SMichael Lange         for (corner = 0; corner < numCorners; ++corner) pcone[corner] = cone[corner] + trueNumCells-1;
217*db397164SMichael Lange         ierr = DMPlexSetCone(*dm, c-numFacets, pcone);CHKERRQ(ierr);
218*db397164SMichael Lange       }
219331830f3SMatthew G. Knepley       if (i != c+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell number %d should be %d", i, c+1);
220331830f3SMatthew G. Knepley     }
221331830f3SMatthew G. Knepley     fgets(line, PETSC_MAX_PATH_LEN, fd);
222331830f3SMatthew G. Knepley     ierr = PetscStrncmp(line, "$EndElements\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
223331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
224331830f3SMatthew G. Knepley   }
2256228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
226331830f3SMatthew G. Knepley   ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr);
227331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
228331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
229331830f3SMatthew G. Knepley   if (interpolate) {
230331830f3SMatthew G. Knepley     DM idm;
231331830f3SMatthew G. Knepley 
232331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
233331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
234331830f3SMatthew G. Knepley     *dm  = idm;
235331830f3SMatthew G. Knepley   }
236331830f3SMatthew G. Knepley   /* Read coordinates */
237979c4b60SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
238331830f3SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
239331830f3SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
240331830f3SMatthew G. Knepley   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
241331830f3SMatthew G. Knepley   for (v = numCells; v < numCells+numVertices; ++v) {
242331830f3SMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
243331830f3SMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
244331830f3SMatthew G. Knepley   }
245331830f3SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
246331830f3SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
247331830f3SMatthew G. Knepley   ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr);
248331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
249331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
250331830f3SMatthew G. Knepley   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
251331830f3SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
252331830f3SMatthew G. Knepley   if (!rank) {
253331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
254331830f3SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
255331830f3SMatthew G. Knepley         coords[v*dim+d] = coordsIn[v*3+d];
256331830f3SMatthew G. Knepley       }
257331830f3SMatthew G. Knepley     }
258331830f3SMatthew G. Knepley   }
259331830f3SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
260331830f3SMatthew G. Knepley   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
261331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
262331830f3SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
263331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
264331830f3SMatthew G. Knepley }
265*db397164SMichael Lange 
266*db397164SMichael Lange #undef __FUNCT__
267*db397164SMichael Lange #define __FUNCT__ "DMPlexView_Ascii"
268*db397164SMichael Lange PetscErrorCode DMPlexView_Ascii(DM dm, PetscViewer viewer)
269*db397164SMichael Lange {
270*db397164SMichael Lange 
271*db397164SMichael Lange }
272