xref: /petsc/src/dm/impls/plex/plexvtk.c (revision e6ccafaeb117de91eefc647fc1341c9a4ed99e5b)
1552f7358SJed Brown #define PETSCDM_DLL
2552f7358SJed Brown #include <petsc-private/pleximpl.h>    /*I   "petscdmplex.h"   I*/
3552f7358SJed Brown #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h>
4552f7358SJed Brown 
5552f7358SJed Brown #undef __FUNCT__
6552f7358SJed Brown #define __FUNCT__ "DMPlexVTKGetCellType"
70adebc6cSBarry Smith PetscErrorCode DMPlexVTKGetCellType(DM dm, PetscInt dim, PetscInt corners, PetscInt *cellType)
80adebc6cSBarry Smith {
9552f7358SJed Brown   PetscFunctionBegin;
10552f7358SJed Brown   *cellType = -1;
11552f7358SJed Brown   switch (dim) {
12552f7358SJed Brown   case 0:
13552f7358SJed Brown     switch (corners) {
14552f7358SJed Brown     case 1:
15552f7358SJed Brown       *cellType = 1; /* VTK_VERTEX */
16552f7358SJed Brown       break;
17552f7358SJed Brown     default:
18552f7358SJed Brown       break;
19552f7358SJed Brown     }
20552f7358SJed Brown     break;
21552f7358SJed Brown   case 1:
22552f7358SJed Brown     switch (corners) {
23552f7358SJed Brown     case 2:
24552f7358SJed Brown       *cellType = 3; /* VTK_LINE */
25552f7358SJed Brown       break;
26552f7358SJed Brown     case 3:
27552f7358SJed Brown       *cellType = 21; /* VTK_QUADRATIC_EDGE */
28552f7358SJed Brown       break;
29552f7358SJed Brown     default:
30552f7358SJed Brown       break;
31552f7358SJed Brown     }
32552f7358SJed Brown     break;
33552f7358SJed Brown   case 2:
34552f7358SJed Brown     switch (corners) {
35552f7358SJed Brown     case 3:
36552f7358SJed Brown       *cellType = 5; /* VTK_TRIANGLE */
37552f7358SJed Brown       break;
38552f7358SJed Brown     case 4:
39552f7358SJed Brown       *cellType = 9; /* VTK_QUAD */
40552f7358SJed Brown       break;
41552f7358SJed Brown     case 6:
42552f7358SJed Brown       *cellType = 22; /* VTK_QUADRATIC_TRIANGLE */
43552f7358SJed Brown       break;
44552f7358SJed Brown     case 9:
45552f7358SJed Brown       *cellType = 23; /* VTK_QUADRATIC_QUAD */
46552f7358SJed Brown       break;
47552f7358SJed Brown     default:
48552f7358SJed Brown       break;
49552f7358SJed Brown     }
50552f7358SJed Brown     break;
51552f7358SJed Brown   case 3:
52552f7358SJed Brown     switch (corners) {
53552f7358SJed Brown     case 4:
54552f7358SJed Brown       *cellType = 10; /* VTK_TETRA */
55552f7358SJed Brown       break;
56552f7358SJed Brown     case 8:
57552f7358SJed Brown       *cellType = 12; /* VTK_HEXAHEDRON */
58552f7358SJed Brown       break;
59552f7358SJed Brown     case 10:
60552f7358SJed Brown       *cellType = 24; /* VTK_QUADRATIC_TETRA */
61552f7358SJed Brown       break;
62552f7358SJed Brown     case 27:
63552f7358SJed Brown       *cellType = 29; /* VTK_QUADRATIC_HEXAHEDRON */
64552f7358SJed Brown       break;
65552f7358SJed Brown     default:
66552f7358SJed Brown       break;
67552f7358SJed Brown     }
68552f7358SJed Brown   }
69552f7358SJed Brown   PetscFunctionReturn(0);
70552f7358SJed Brown }
71552f7358SJed Brown 
72552f7358SJed Brown #undef __FUNCT__
73552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteCells_ASCII"
74552f7358SJed Brown PetscErrorCode DMPlexVTKWriteCells_ASCII(DM dm, FILE *fp, PetscInt *totalCells)
75552f7358SJed Brown {
76552f7358SJed Brown   MPI_Comm       comm = ((PetscObject) dm)->comm;
77552f7358SJed Brown   IS             globalVertexNumbers;
78552f7358SJed Brown   const PetscInt *gvertex;
79552f7358SJed Brown   PetscInt       dim;
80552f7358SJed Brown   PetscInt       numCorners = 0, totCorners = 0, maxCorners, *corners;
81552f7358SJed Brown   PetscInt       numCells   = 0, totCells   = 0, maxCells, cellHeight;
82552f7358SJed Brown   PetscInt       numLabelCells, cMax, cStart, cEnd, c, vStart, vEnd, v;
83552f7358SJed Brown   PetscMPIInt    numProcs, rank, proc, tag;
84552f7358SJed Brown   PetscBool      hasLabel;
85552f7358SJed Brown   PetscErrorCode ierr;
86552f7358SJed Brown 
87552f7358SJed Brown   PetscFunctionBegin;
88552f7358SJed Brown   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
89552f7358SJed Brown   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
90552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
91552f7358SJed Brown   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
92552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
93552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
94552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
95770b213bSMatthew G Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
968865f1eaSKarl Rupp   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
97552f7358SJed Brown   ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
988865f1eaSKarl Rupp 
99552f7358SJed Brown   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
100552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
101552f7358SJed Brown     PetscInt *closure = PETSC_NULL;
102552f7358SJed Brown     PetscInt closureSize;
103552f7358SJed Brown 
104552f7358SJed Brown     if (hasLabel) {
105552f7358SJed Brown       PetscInt value;
106552f7358SJed Brown 
107552f7358SJed Brown       ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
108552f7358SJed Brown       if (value != 1) continue;
109552f7358SJed Brown     }
110552f7358SJed Brown     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
111552f7358SJed Brown     for (v = 0; v < closureSize*2; v += 2) {
112552f7358SJed Brown       if ((closure[v] >= vStart) && (closure[v] < vEnd)) ++numCorners;
113552f7358SJed Brown     }
114552f7358SJed Brown     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
115552f7358SJed Brown     ++numCells;
116552f7358SJed Brown   }
117552f7358SJed Brown   maxCells = numCells;
118552f7358SJed Brown   ierr     = MPI_Reduce(&numCells, &totCells, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRQ(ierr);
119552f7358SJed Brown   ierr     = MPI_Reduce(&numCells, &maxCells, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr);
120552f7358SJed Brown   ierr     = MPI_Reduce(&numCorners, &totCorners, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRQ(ierr);
121552f7358SJed Brown   ierr     = MPI_Reduce(&numCorners, &maxCorners, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr);
122552f7358SJed Brown   ierr     = DMPlexGetVertexNumbering(dm, &globalVertexNumbers);CHKERRQ(ierr);
123552f7358SJed Brown   ierr     = ISGetIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr);
124552f7358SJed Brown   ierr     = PetscMalloc(maxCells * sizeof(PetscInt), &corners);CHKERRQ(ierr);
125552f7358SJed Brown   ierr     = PetscFPrintf(comm, fp, "CELLS %d %d\n", totCells, totCorners+totCells);CHKERRQ(ierr);
126552f7358SJed Brown   if (!rank) {
127552f7358SJed Brown     PetscInt *remoteVertices;
128552f7358SJed Brown 
129552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
130552f7358SJed Brown       PetscInt *closure = PETSC_NULL;
131552f7358SJed Brown       PetscInt closureSize, nC = 0;
132552f7358SJed Brown 
133552f7358SJed Brown       if (hasLabel) {
134552f7358SJed Brown         PetscInt value;
135552f7358SJed Brown 
136552f7358SJed Brown         ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
137552f7358SJed Brown         if (value != 1) continue;
138552f7358SJed Brown       }
139552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
140552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
141552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
142552f7358SJed Brown           const PetscInt gv = gvertex[closure[v] - vStart];
143552f7358SJed Brown           closure[nC++] = gv < 0 ? -(gv+1) : gv;
144552f7358SJed Brown         }
145552f7358SJed Brown       }
146552f7358SJed Brown       corners[numCells++] = nC;
147552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "%d ", nC);CHKERRQ(ierr);
148552f7358SJed Brown       for (v = 0; v < nC; ++v) {
149552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, " %d", closure[v]);CHKERRQ(ierr);
150552f7358SJed Brown       }
151552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
152552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
153552f7358SJed Brown     }
154552f7358SJed Brown     if (numProcs > 1) {ierr = PetscMalloc((maxCorners+maxCells) * sizeof(PetscInt), &remoteVertices);CHKERRQ(ierr);}
155552f7358SJed Brown     for (proc = 1; proc < numProcs; ++proc) {
156552f7358SJed Brown       MPI_Status status;
157552f7358SJed Brown 
158552f7358SJed Brown       ierr = MPI_Recv(&numCorners, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
159552f7358SJed Brown       ierr = MPI_Recv(remoteVertices, numCorners, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
160552f7358SJed Brown       for (c = 0; c < numCorners;) {
161552f7358SJed Brown         PetscInt nC = remoteVertices[c++];
162552f7358SJed Brown 
163552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "%d ", nC);CHKERRQ(ierr);
164552f7358SJed Brown         for (v = 0; v < nC; ++v, ++c) {
165552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " %d", remoteVertices[c]);CHKERRQ(ierr);
166552f7358SJed Brown         }
167552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
168552f7358SJed Brown       }
169552f7358SJed Brown     }
170552f7358SJed Brown     if (numProcs > 1) {ierr = PetscFree(remoteVertices);CHKERRQ(ierr);}
171552f7358SJed Brown   } else {
172552f7358SJed Brown     PetscInt *localVertices, numSend = numCells+numCorners, k = 0;
173552f7358SJed Brown 
174552f7358SJed Brown     ierr = PetscMalloc(numSend * sizeof(PetscInt), &localVertices);CHKERRQ(ierr);
175552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
176552f7358SJed Brown       PetscInt *closure = PETSC_NULL;
177552f7358SJed Brown       PetscInt closureSize, nC = 0;
178552f7358SJed Brown 
179552f7358SJed Brown       if (hasLabel) {
180552f7358SJed Brown         PetscInt value;
181552f7358SJed Brown 
182552f7358SJed Brown         ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
183552f7358SJed Brown         if (value != 1) continue;
184552f7358SJed Brown       }
185552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
186552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
187552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
188552f7358SJed Brown           const PetscInt gv = gvertex[closure[v] - vStart];
189552f7358SJed Brown           closure[nC++] = gv < 0 ? -(gv+1) : gv;
190552f7358SJed Brown         }
191552f7358SJed Brown       }
192552f7358SJed Brown       corners[numCells++] = nC;
193552f7358SJed Brown       localVertices[k++]  = nC;
194552f7358SJed Brown       for (v = 0; v < nC; ++v, ++k) {
195552f7358SJed Brown         localVertices[k] = closure[v];
196552f7358SJed Brown       }
197552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
198552f7358SJed Brown     }
199552f7358SJed Brown     if (k != numSend) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid number of vertices to send %d should be %d", k, numSend);
200552f7358SJed Brown     ierr = MPI_Send(&numSend, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
201552f7358SJed Brown     ierr = MPI_Send(localVertices, numSend, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
202552f7358SJed Brown     ierr = PetscFree(localVertices);CHKERRQ(ierr);
203552f7358SJed Brown   }
204552f7358SJed Brown   ierr = ISRestoreIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr);
205552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "CELL_TYPES %d\n", totCells);CHKERRQ(ierr);
206552f7358SJed Brown   if (!rank) {
207552f7358SJed Brown     PetscInt cellType;
208552f7358SJed Brown 
209552f7358SJed Brown     for (c = 0; c < numCells; ++c) {
210552f7358SJed Brown       ierr = DMPlexVTKGetCellType(dm, dim, corners[c], &cellType);CHKERRQ(ierr);
211552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "%d\n", cellType);CHKERRQ(ierr);
212552f7358SJed Brown     }
213552f7358SJed Brown     for (proc = 1; proc < numProcs; ++proc) {
214552f7358SJed Brown       MPI_Status status;
215552f7358SJed Brown 
216552f7358SJed Brown       ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
217552f7358SJed Brown       ierr = MPI_Recv(corners, numCells, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
218552f7358SJed Brown       for (c = 0; c < numCells; ++c) {
219552f7358SJed Brown         ierr = DMPlexVTKGetCellType(dm, dim, corners[c], &cellType);CHKERRQ(ierr);
220552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "%d\n", cellType);CHKERRQ(ierr);
221552f7358SJed Brown       }
222552f7358SJed Brown     }
223552f7358SJed Brown   } else {
224552f7358SJed Brown     ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
225552f7358SJed Brown     ierr = MPI_Send(corners, numCells, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
226552f7358SJed Brown   }
227552f7358SJed Brown   ierr        = PetscFree(corners);CHKERRQ(ierr);
228552f7358SJed Brown   *totalCells = totCells;
229552f7358SJed Brown   PetscFunctionReturn(0);
230552f7358SJed Brown }
231552f7358SJed Brown 
232552f7358SJed Brown #undef __FUNCT__
233552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteSection_ASCII"
2340adebc6cSBarry Smith PetscErrorCode DMPlexVTKWriteSection_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec v, FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale)
2350adebc6cSBarry Smith {
236552f7358SJed Brown   MPI_Comm           comm    = ((PetscObject) dm)->comm;
237552f7358SJed Brown   const MPI_Datatype mpiType = MPIU_SCALAR;
238552f7358SJed Brown   PetscScalar        *array;
239552f7358SJed Brown   PetscInt           numDof = 0, maxDof;
240552f7358SJed Brown   PetscInt           numLabelCells, cellHeight, cMax, cStart, cEnd, numLabelVertices, vMax, vStart, vEnd, pStart, pEnd, p;
241552f7358SJed Brown   PetscMPIInt        numProcs, rank, proc, tag;
242552f7358SJed Brown   PetscBool          hasLabel;
243552f7358SJed Brown   PetscErrorCode     ierr;
244552f7358SJed Brown 
245552f7358SJed Brown   PetscFunctionBegin;
246552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
247552f7358SJed Brown   PetscValidHeaderSpecific(v,VEC_CLASSID,4);
248552f7358SJed Brown   if (precision < 0) precision = 6;
249552f7358SJed Brown   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
250552f7358SJed Brown   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
251552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
252552f7358SJed Brown   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
253552f7358SJed Brown   /* VTK only wants the values at cells or vertices */
254552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
255552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
256552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
257770b213bSMatthew G Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, PETSC_NULL, PETSC_NULL, &vMax);CHKERRQ(ierr);
2588865f1eaSKarl Rupp   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
2598865f1eaSKarl Rupp   if (vMax >= 0) vEnd = PetscMin(vEnd, vMax);
260552f7358SJed Brown   pStart   = PetscMax(PetscMin(cStart, vStart), pStart);
261552f7358SJed Brown   pEnd     = PetscMin(PetscMax(cEnd,   vEnd),   pEnd);
262552f7358SJed Brown   ierr     = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
263552f7358SJed Brown   ierr     = DMPlexGetStratumSize(dm, "vtk", 2, &numLabelVertices);CHKERRQ(ierr);
264552f7358SJed Brown   hasLabel = numLabelCells > 0 || numLabelVertices > 0 ? PETSC_TRUE : PETSC_FALSE;
265552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
266552f7358SJed Brown     /* Reject points not either cells or vertices */
267552f7358SJed Brown     if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
268552f7358SJed Brown     if (hasLabel) {
269552f7358SJed Brown       PetscInt value;
270552f7358SJed Brown 
271552f7358SJed Brown       if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
272552f7358SJed Brown           ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
273552f7358SJed Brown         ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
274552f7358SJed Brown         if (value != 1) continue;
275552f7358SJed Brown       }
276552f7358SJed Brown     }
277552f7358SJed Brown     ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr);
2788865f1eaSKarl Rupp     if (numDof) break;
279552f7358SJed Brown   }
280552f7358SJed Brown   ierr = MPI_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
281552f7358SJed Brown   enforceDof = PetscMax(enforceDof, maxDof);
282552f7358SJed Brown   ierr = VecGetArray(v, &array);CHKERRQ(ierr);
283552f7358SJed Brown   if (!rank) {
284552f7358SJed Brown     char formatString[8];
285552f7358SJed Brown 
286552f7358SJed Brown     ierr = PetscSNPrintf(formatString, 8, "%%.%de", precision);CHKERRQ(ierr);
287552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
288552f7358SJed Brown       /* Here we lose a way to filter points by keeping them out of the Numbering */
289552f7358SJed Brown       PetscInt dof, off, goff, d;
290552f7358SJed Brown 
291552f7358SJed Brown       /* Reject points not either cells or vertices */
292552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
293552f7358SJed Brown       if (hasLabel) {
294552f7358SJed Brown         PetscInt value;
295552f7358SJed Brown 
296552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
297552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
298552f7358SJed Brown           ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
299552f7358SJed Brown           if (value != 1) continue;
300552f7358SJed Brown         }
301552f7358SJed Brown       }
302552f7358SJed Brown       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
303552f7358SJed Brown       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
304552f7358SJed Brown       ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
305552f7358SJed Brown       if (dof && goff >= 0) {
306552f7358SJed Brown         for (d = 0; d < dof; d++) {
307552f7358SJed Brown           if (d > 0) {
308552f7358SJed Brown             ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr);
309552f7358SJed Brown           }
310552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(array[off+d])*scale);CHKERRQ(ierr);
311552f7358SJed Brown         }
312552f7358SJed Brown         for (d = dof; d < enforceDof; d++) {
313552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr);
314552f7358SJed Brown         }
315552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
316552f7358SJed Brown       }
317552f7358SJed Brown     }
318552f7358SJed Brown     for (proc = 1; proc < numProcs; ++proc) {
319552f7358SJed Brown       PetscScalar *remoteValues;
320552f7358SJed Brown       PetscInt    size, d;
321552f7358SJed Brown       MPI_Status  status;
322552f7358SJed Brown 
323552f7358SJed Brown       ierr = MPI_Recv(&size, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
324552f7358SJed Brown       ierr = PetscMalloc(size * sizeof(PetscScalar), &remoteValues);CHKERRQ(ierr);
325552f7358SJed Brown       ierr = MPI_Recv(remoteValues, size, mpiType, proc, tag, comm, &status);CHKERRQ(ierr);
326552f7358SJed Brown       for (p = 0; p < size/maxDof; ++p) {
327552f7358SJed Brown         for (d = 0; d < maxDof; ++d) {
328552f7358SJed Brown           if (d > 0) {
329552f7358SJed Brown             ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr);
330552f7358SJed Brown           }
331552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(remoteValues[p*maxDof+d])*scale);CHKERRQ(ierr);
332552f7358SJed Brown         }
333552f7358SJed Brown         for (d = maxDof; d < enforceDof; ++d) {
334552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr);
335552f7358SJed Brown         }
336552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
337552f7358SJed Brown       }
338552f7358SJed Brown       ierr = PetscFree(remoteValues);CHKERRQ(ierr);
339552f7358SJed Brown     }
340552f7358SJed Brown   } else {
341552f7358SJed Brown     PetscScalar *localValues;
342552f7358SJed Brown     PetscInt    size, k = 0;
343552f7358SJed Brown 
344552f7358SJed Brown     ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
345552f7358SJed Brown     ierr = PetscMalloc(size * sizeof(PetscScalar), &localValues);CHKERRQ(ierr);
346552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
347552f7358SJed Brown       PetscInt dof, off, goff, d;
348552f7358SJed Brown 
349552f7358SJed Brown       /* Reject points not either cells or vertices */
350552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
351552f7358SJed Brown       if (hasLabel) {
352552f7358SJed Brown         PetscInt value;
353552f7358SJed Brown 
354552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
355552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
356552f7358SJed Brown           ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
357552f7358SJed Brown           if (value != 1) continue;
358552f7358SJed Brown         }
359552f7358SJed Brown       }
360552f7358SJed Brown       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
361552f7358SJed Brown       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
362552f7358SJed Brown       ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
363552f7358SJed Brown       if (goff >= 0) {
364552f7358SJed Brown         for (d = 0; d < dof; ++d) {
365552f7358SJed Brown           localValues[k++] = array[off+d];
366552f7358SJed Brown         }
367552f7358SJed Brown       }
368552f7358SJed Brown     }
369552f7358SJed Brown     ierr = MPI_Send(&k, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
370552f7358SJed Brown     ierr = MPI_Send(localValues, k, mpiType, 0, tag, comm);CHKERRQ(ierr);
371552f7358SJed Brown     ierr = PetscFree(localValues);CHKERRQ(ierr);
372552f7358SJed Brown   }
373552f7358SJed Brown   ierr = VecRestoreArray(v, &array);CHKERRQ(ierr);
374552f7358SJed Brown   PetscFunctionReturn(0);
375552f7358SJed Brown }
376552f7358SJed Brown 
377552f7358SJed Brown #undef __FUNCT__
378552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteField_ASCII"
379552f7358SJed Brown PetscErrorCode DMPlexVTKWriteField_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec field, const char name[], FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale)
380552f7358SJed Brown {
381552f7358SJed Brown   MPI_Comm       comm   = ((PetscObject) dm)->comm;
382552f7358SJed Brown   PetscInt       numDof = 0, maxDof;
383552f7358SJed Brown   PetscInt       pStart, pEnd, p;
384552f7358SJed Brown   PetscErrorCode ierr;
385552f7358SJed Brown 
386552f7358SJed Brown   PetscFunctionBegin;
387552f7358SJed Brown   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
388552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
389552f7358SJed Brown     ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr);
3908865f1eaSKarl Rupp     if (numDof) break;
391552f7358SJed Brown   }
392552f7358SJed Brown   numDof = PetscMax(numDof, enforceDof);
393552f7358SJed Brown   ierr = MPI_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, ((PetscObject) dm)->comm);CHKERRQ(ierr);
394552f7358SJed Brown   if (!name) name = "Unknown";
395552f7358SJed Brown   if (maxDof == 3) {
396552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "VECTORS %s double\n", name);CHKERRQ(ierr);
397552f7358SJed Brown   } else {
398552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "SCALARS %s double %d\n", name, maxDof);CHKERRQ(ierr);
399552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr);
400552f7358SJed Brown   }
401552f7358SJed Brown   ierr = DMPlexVTKWriteSection_ASCII(dm, section, globalSection, field, fp, enforceDof, precision, scale);CHKERRQ(ierr);
402552f7358SJed Brown   PetscFunctionReturn(0);
403552f7358SJed Brown }
404552f7358SJed Brown 
405552f7358SJed Brown #undef __FUNCT__
406552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll_ASCII"
407552f7358SJed Brown static PetscErrorCode DMPlexVTKWriteAll_ASCII(DM dm, PetscViewer viewer)
408552f7358SJed Brown {
409552f7358SJed Brown   MPI_Comm                 comm = ((PetscObject) dm)->comm;
410552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*) viewer->data;
411552f7358SJed Brown   FILE                     *fp;
412552f7358SJed Brown   PetscViewerVTKObjectLink link;
413552f7358SJed Brown   PetscSection             coordSection, globalCoordSection;
414552f7358SJed Brown   PetscLayout              vLayout;
415552f7358SJed Brown   Vec                      coordinates;
416552f7358SJed Brown   PetscReal                lengthScale;
417552f7358SJed Brown   PetscInt                 vMax, totVertices, totCells;
418552f7358SJed Brown   PetscBool                hasPoint = PETSC_FALSE, hasCell = PETSC_FALSE;
419552f7358SJed Brown   PetscErrorCode           ierr;
420552f7358SJed Brown 
421552f7358SJed Brown   PetscFunctionBegin;
422552f7358SJed Brown   ierr = PetscFOpen(comm, vtk->filename, "wb", &fp);CHKERRQ(ierr);
423552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "# vtk DataFile Version 2.0\n");CHKERRQ(ierr);
424552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "Simplicial Mesh Example\n");CHKERRQ(ierr);
425552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "ASCII\n");CHKERRQ(ierr);
426552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "DATASET UNSTRUCTURED_GRID\n");CHKERRQ(ierr);
427552f7358SJed Brown   /* Vertices */
428552f7358SJed Brown   ierr = DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);CHKERRQ(ierr);
429552f7358SJed Brown   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
430552f7358SJed Brown   ierr = PetscSectionCreateGlobalSection(coordSection, dm->sf, PETSC_FALSE, &globalCoordSection);CHKERRQ(ierr);
431552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
432770b213bSMatthew G Knepley   ierr = DMPlexGetHybridBounds(dm, PETSC_NULL, PETSC_NULL, PETSC_NULL, &vMax);CHKERRQ(ierr);
433552f7358SJed Brown   if (vMax >= 0) {
434552f7358SJed Brown     PetscInt pStart, pEnd, p, localSize = 0;
435552f7358SJed Brown 
436552f7358SJed Brown     ierr = PetscSectionGetChart(globalCoordSection, &pStart, &pEnd);CHKERRQ(ierr);
437552f7358SJed Brown     pEnd = PetscMin(pEnd, vMax);
438552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
439552f7358SJed Brown       PetscInt dof;
440552f7358SJed Brown 
441552f7358SJed Brown       ierr = PetscSectionGetDof(globalCoordSection, p, &dof);CHKERRQ(ierr);
4428865f1eaSKarl Rupp       if (dof > 0) ++localSize;
443552f7358SJed Brown     }
444552f7358SJed Brown     ierr = PetscLayoutCreate(((PetscObject) dm)->comm, &vLayout);CHKERRQ(ierr);
445552f7358SJed Brown     ierr = PetscLayoutSetLocalSize(vLayout, localSize);CHKERRQ(ierr);
446552f7358SJed Brown     ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr);
447552f7358SJed Brown     ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr);
448552f7358SJed Brown   } else {
449552f7358SJed Brown     ierr = PetscSectionGetPointLayout(((PetscObject) dm)->comm, globalCoordSection, &vLayout);CHKERRQ(ierr);
450552f7358SJed Brown   }
451552f7358SJed Brown   ierr = PetscLayoutGetSize(vLayout, &totVertices);CHKERRQ(ierr);
452552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "POINTS %d double\n", totVertices);CHKERRQ(ierr);
453552f7358SJed Brown   ierr = DMPlexVTKWriteSection_ASCII(dm, coordSection, globalCoordSection, coordinates, fp, 3, PETSC_DETERMINE, lengthScale);CHKERRQ(ierr);
454552f7358SJed Brown   /* Cells */
455552f7358SJed Brown   ierr = DMPlexVTKWriteCells_ASCII(dm, fp, &totCells);CHKERRQ(ierr);
456552f7358SJed Brown   /* Vertex fields */
457552f7358SJed Brown   for (link = vtk->link; link; link = link->next) {
458552f7358SJed Brown     if ((link->ft == PETSC_VTK_POINT_FIELD) || (link->ft == PETSC_VTK_POINT_VECTOR_FIELD)) hasPoint = PETSC_TRUE;
459552f7358SJed Brown     if ((link->ft == PETSC_VTK_CELL_FIELD)  || (link->ft == PETSC_VTK_CELL_VECTOR_FIELD))  hasCell  = PETSC_TRUE;
460552f7358SJed Brown   }
461552f7358SJed Brown   if (hasPoint) {
46222d28d08SBarry Smith     ierr = PetscFPrintf(comm, fp, "POINT_DATA %d\n", totVertices);CHKERRQ(ierr);
463552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
464552f7358SJed Brown       Vec          X = (Vec) link->vec;
465552f7358SJed Brown       DM           dmX;
466839dd189SMatthew G Knepley       PetscSection section, globalSection, newSection = PETSC_NULL;
467552f7358SJed Brown       const char   *name;
468552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
469552f7358SJed Brown 
470552f7358SJed Brown       if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
471552f7358SJed Brown       if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) enforceDof = 3;
472552f7358SJed Brown       ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr);
473552f7358SJed Brown       ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
474552f7358SJed Brown       if (dmX) {
475088580cfSMatthew G Knepley         DMLabel  subpointMap, subpointMapX;
476839dd189SMatthew G Knepley         PetscInt dim, dimX, pStart, pEnd, qStart, qEnd;
477839dd189SMatthew G Knepley 
47822d28d08SBarry Smith         ierr = DMGetDefaultSection(dmX, &section);CHKERRQ(ierr);
479839dd189SMatthew G Knepley         /* Here is where we check whether dmX is a submesh of dm */
480839dd189SMatthew G Knepley         ierr = DMPlexGetDimension(dm,  &dim);CHKERRQ(ierr);
481839dd189SMatthew G Knepley         ierr = DMPlexGetDimension(dmX, &dimX);CHKERRQ(ierr);
482839dd189SMatthew G Knepley         ierr = DMPlexGetChart(dm,  &pStart, &pEnd);CHKERRQ(ierr);
483839dd189SMatthew G Knepley         ierr = DMPlexGetChart(dmX, &qStart, &qEnd);CHKERRQ(ierr);
484839dd189SMatthew G Knepley         ierr = DMPlexGetSubpointMap(dm,  &subpointMap);CHKERRQ(ierr);
485839dd189SMatthew G Knepley         ierr = DMPlexGetSubpointMap(dmX, &subpointMapX);CHKERRQ(ierr);
486839dd189SMatthew G Knepley         if (((dim != dimX) || ((pEnd-pStart) < (qEnd-qStart))) && subpointMap && !subpointMapX) {
487839dd189SMatthew G Knepley           const PetscInt *ind;
488*e6ccafaeSMatthew G Knepley           IS              subpointIS;
489839dd189SMatthew G Knepley           PetscInt        n, q;
490839dd189SMatthew G Knepley 
491839dd189SMatthew G Knepley           ierr = PetscPrintf(PETSC_COMM_SELF, "Making translation PetscSection\n");CHKERRQ(ierr);
492839dd189SMatthew G Knepley           ierr = PetscSectionGetChart(section, &qStart, &qEnd);CHKERRQ(ierr);
493*e6ccafaeSMatthew G Knepley           ierr = DMPlexCreateSubpointIS(dm, &subpointIS);CHKERRQ(ierr);
494*e6ccafaeSMatthew G Knepley           ierr = ISGetLocalSize(subpointIS, &n);CHKERRQ(ierr);
495*e6ccafaeSMatthew G Knepley           ierr = ISGetIndices(subpointIS, &ind);CHKERRQ(ierr);
496839dd189SMatthew G Knepley           ierr = PetscSectionCreate(comm, &newSection);CHKERRQ(ierr);
497839dd189SMatthew G Knepley           ierr = PetscSectionSetChart(newSection, pStart, pEnd);CHKERRQ(ierr);
498839dd189SMatthew G Knepley           for (q = qStart; q < qEnd; ++q) {
499839dd189SMatthew G Knepley             PetscInt dof, off, p;
500839dd189SMatthew G Knepley 
501839dd189SMatthew G Knepley             ierr = PetscSectionGetDof(section, q, &dof);CHKERRQ(ierr);
502839dd189SMatthew G Knepley             if (dof) {
503839dd189SMatthew G Knepley               ierr = PetscFindInt(q, n, ind, &p);CHKERRQ(ierr);
504839dd189SMatthew G Knepley               if (p >= pStart) {
505839dd189SMatthew G Knepley                 ierr = PetscSectionSetDof(newSection, p, dof);CHKERRQ(ierr);
506839dd189SMatthew G Knepley                 ierr = PetscSectionGetOffset(section, q, &off);CHKERRQ(ierr);
507839dd189SMatthew G Knepley                 ierr = PetscSectionSetOffset(newSection, p, off);CHKERRQ(ierr);
508839dd189SMatthew G Knepley               }
509839dd189SMatthew G Knepley             }
510839dd189SMatthew G Knepley           }
511*e6ccafaeSMatthew G Knepley           ierr = ISRestoreIndices(subpointIS, &ind);CHKERRQ(ierr);
512*e6ccafaeSMatthew G Knepley           ierr = ISDestroy(&subpointIS);CHKERRQ(ierr);
513839dd189SMatthew G Knepley           /* No need to setup section */
514839dd189SMatthew G Knepley           section = newSection;
515839dd189SMatthew G Knepley         }
516552f7358SJed Brown       } else {
517552f7358SJed Brown         PetscContainer c;
518552f7358SJed Brown 
519552f7358SJed Brown         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &c);CHKERRQ(ierr);
520552f7358SJed Brown         if (!c) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
521552f7358SJed Brown         ierr = PetscContainerGetPointer(c, (void**) &section);CHKERRQ(ierr);
522552f7358SJed Brown       }
523552f7358SJed Brown       if (!section) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
524552f7358SJed Brown       ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, &globalSection);CHKERRQ(ierr);
525552f7358SJed Brown       ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr);
526552f7358SJed Brown       ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr);
527839dd189SMatthew G Knepley       if (newSection) {ierr = PetscSectionDestroy(&newSection);CHKERRQ(ierr);}
528552f7358SJed Brown     }
529552f7358SJed Brown   }
530552f7358SJed Brown   /* Cell Fields */
531552f7358SJed Brown   if (hasCell) {
53222d28d08SBarry Smith     ierr = PetscFPrintf(comm, fp, "CELL_DATA %d\n", totCells);CHKERRQ(ierr);
533552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
534552f7358SJed Brown       Vec          X = (Vec) link->vec;
535552f7358SJed Brown       DM           dmX;
536552f7358SJed Brown       PetscSection section, globalSection;
537552f7358SJed Brown       const char   *name;
538552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
539552f7358SJed Brown 
540552f7358SJed Brown       if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
541552f7358SJed Brown       if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) enforceDof = 3;
542552f7358SJed Brown       ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr);
543552f7358SJed Brown       ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
544552f7358SJed Brown       if (dmX) {
54522d28d08SBarry Smith         ierr = DMGetDefaultSection(dmX, &section);CHKERRQ(ierr);
546552f7358SJed Brown       } else {
547552f7358SJed Brown         PetscContainer c;
548552f7358SJed Brown 
549552f7358SJed Brown         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &c);CHKERRQ(ierr);
550552f7358SJed Brown         if (!c) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
551552f7358SJed Brown         ierr = PetscContainerGetPointer(c, (void**) &section);CHKERRQ(ierr);
552552f7358SJed Brown       }
553552f7358SJed Brown       if (!section) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
554552f7358SJed Brown       ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, &globalSection);CHKERRQ(ierr);
555552f7358SJed Brown       ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr);
556552f7358SJed Brown       ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr);
557552f7358SJed Brown     }
558552f7358SJed Brown   }
559552f7358SJed Brown   /* Cleanup */
560552f7358SJed Brown   ierr = PetscSectionDestroy(&globalCoordSection);CHKERRQ(ierr);
561552f7358SJed Brown   ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr);
562552f7358SJed Brown   ierr = PetscFClose(comm, fp);CHKERRQ(ierr);
563552f7358SJed Brown   PetscFunctionReturn(0);
564552f7358SJed Brown }
565552f7358SJed Brown 
566552f7358SJed Brown #undef __FUNCT__
567552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll"
568552f7358SJed Brown /*@C
569552f7358SJed Brown   DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
570552f7358SJed Brown 
571552f7358SJed Brown   Collective
572552f7358SJed Brown 
573552f7358SJed Brown   Input Arguments:
574552f7358SJed Brown + odm - The DMPlex specifying the mesh, passed as a PetscObject
575552f7358SJed Brown - viewer - viewer of type VTK
576552f7358SJed Brown 
577552f7358SJed Brown   Level: developer
578552f7358SJed Brown 
579552f7358SJed Brown   Note:
580552f7358SJed Brown   This function is a callback used by the VTK viewer to actually write the file.
581552f7358SJed Brown   The reason for this odd model is that the VTK file format does not provide any way to write one field at a time.
582552f7358SJed Brown   Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
583552f7358SJed Brown 
584552f7358SJed Brown .seealso: PETSCVIEWERVTK
585552f7358SJed Brown @*/
586552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer)
587552f7358SJed Brown {
588552f7358SJed Brown   DM             dm = (DM) odm;
589552f7358SJed Brown   PetscBool      isvtk;
590552f7358SJed Brown   PetscErrorCode ierr;
591552f7358SJed Brown 
592552f7358SJed Brown   PetscFunctionBegin;
593552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
594552f7358SJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
595552f7358SJed Brown   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr);
596552f7358SJed Brown   if (!isvtk) SETERRQ1(((PetscObject) viewer)->comm, PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name);
597552f7358SJed Brown   switch (viewer->format) {
598552f7358SJed Brown   case PETSC_VIEWER_ASCII_VTK:
599552f7358SJed Brown     ierr = DMPlexVTKWriteAll_ASCII(dm, viewer);CHKERRQ(ierr);
600552f7358SJed Brown     break;
601552f7358SJed Brown   case PETSC_VIEWER_VTK_VTU:
602552f7358SJed Brown     ierr = DMPlexVTKWriteAll_VTU(dm, viewer);CHKERRQ(ierr);
603552f7358SJed Brown     break;
604552f7358SJed Brown   default: SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]);
605552f7358SJed Brown   }
606552f7358SJed Brown   PetscFunctionReturn(0);
607552f7358SJed Brown }
608