xref: /petsc/src/dm/impls/plex/plexvtk.c (revision b2566f29c2b6470df769aa9f7deb9e2726b0959e)
1552f7358SJed Brown #define PETSCDM_DLL
2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.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 {
7682f516ccSBarry Smith   MPI_Comm       comm;
7702281ff3SMatthew G. Knepley   DMLabel        label;
7802281ff3SMatthew G. Knepley   IS             globalVertexNumbers = NULL;
79552f7358SJed Brown   const PetscInt *gvertex;
80552f7358SJed Brown   PetscInt       dim;
81552f7358SJed Brown   PetscInt       numCorners = 0, totCorners = 0, maxCorners, *corners;
82552f7358SJed Brown   PetscInt       numCells   = 0, totCells   = 0, maxCells, cellHeight;
83452b7979SMatthew G. Knepley   PetscInt       numLabelCells, maxLabelCells, cMax, cStart, cEnd, c, vStart, vEnd, v;
84552f7358SJed Brown   PetscMPIInt    numProcs, rank, proc, tag;
85552f7358SJed Brown   PetscErrorCode ierr;
86552f7358SJed Brown 
87552f7358SJed Brown   PetscFunctionBegin;
8882f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
89552f7358SJed Brown   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
90552f7358SJed Brown   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
91552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
92c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
93552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
94552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
95552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
960298fd71SBarry Smith   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
978865f1eaSKarl Rupp   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
9802281ff3SMatthew G. Knepley   ierr = DMPlexGetLabel(dm, "vtk", &label);CHKERRQ(ierr);
99552f7358SJed Brown   ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
100*b2566f29SBarry Smith   ierr = MPIU_Allreduce(&numLabelCells, &maxLabelCells, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
101452b7979SMatthew G. Knepley   if (!maxLabelCells) label = NULL;
102552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
1030298fd71SBarry Smith     PetscInt *closure = NULL;
10402281ff3SMatthew G. Knepley     PetscInt closureSize, value;
105552f7358SJed Brown 
10602281ff3SMatthew G. Knepley     if (label) {
10702281ff3SMatthew G. Knepley       ierr = DMLabelGetValue(label, 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);
124785e854fSJed Brown   ierr     = PetscMalloc1(maxCells, &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;
1280dcf9c70SMatthew G. Knepley     int      *vertices;
129552f7358SJed Brown 
130785e854fSJed Brown     ierr = PetscMalloc1(maxCorners, &vertices);CHKERRQ(ierr);
131552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
1320298fd71SBarry Smith       PetscInt *closure = NULL;
13302281ff3SMatthew G. Knepley       PetscInt closureSize, value, nC = 0;
134552f7358SJed Brown 
13502281ff3SMatthew G. Knepley       if (label) {
13602281ff3SMatthew G. Knepley         ierr = DMLabelGetValue(label, 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];
1430dcf9c70SMatthew G. Knepley           vertices[nC++] = gv < 0 ? -(gv+1) : gv;
144552f7358SJed Brown         }
145552f7358SJed Brown       }
1460dcf9c70SMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
147552f7358SJed Brown       corners[numCells++] = nC;
148552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "%d ", nC);CHKERRQ(ierr);
1490dcf9c70SMatthew G. Knepley       ierr = DMPlexInvertCell(dim, nC, vertices);CHKERRQ(ierr);
150552f7358SJed Brown       for (v = 0; v < nC; ++v) {
1510dcf9c70SMatthew G. Knepley         ierr = PetscFPrintf(comm, fp, " %d", vertices[v]);CHKERRQ(ierr);
152552f7358SJed Brown       }
153552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
154552f7358SJed Brown     }
155854ce69bSBarry Smith     if (numProcs > 1) {ierr = PetscMalloc1(maxCorners+maxCells, &remoteVertices);CHKERRQ(ierr);}
156552f7358SJed Brown     for (proc = 1; proc < numProcs; ++proc) {
157552f7358SJed Brown       MPI_Status status;
158552f7358SJed Brown 
159552f7358SJed Brown       ierr = MPI_Recv(&numCorners, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
160552f7358SJed Brown       ierr = MPI_Recv(remoteVertices, numCorners, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
161552f7358SJed Brown       for (c = 0; c < numCorners;) {
162552f7358SJed Brown         PetscInt nC = remoteVertices[c++];
163552f7358SJed Brown 
164552f7358SJed Brown         for (v = 0; v < nC; ++v, ++c) {
1650dcf9c70SMatthew G. Knepley           vertices[v] = remoteVertices[c];
1660dcf9c70SMatthew G. Knepley         }
1670dcf9c70SMatthew G. Knepley         ierr = DMPlexInvertCell(dim, nC, vertices);CHKERRQ(ierr);
1680dcf9c70SMatthew G. Knepley         ierr = PetscFPrintf(comm, fp, "%d ", nC);CHKERRQ(ierr);
1690dcf9c70SMatthew G. Knepley         for (v = 0; v < nC; ++v) {
1700dcf9c70SMatthew G. Knepley           ierr = PetscFPrintf(comm, fp, " %d", vertices[v]);CHKERRQ(ierr);
171552f7358SJed Brown         }
172552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
173552f7358SJed Brown       }
174552f7358SJed Brown     }
175552f7358SJed Brown     if (numProcs > 1) {ierr = PetscFree(remoteVertices);CHKERRQ(ierr);}
1760dcf9c70SMatthew G. Knepley     ierr = PetscFree(vertices);CHKERRQ(ierr);
177552f7358SJed Brown   } else {
178552f7358SJed Brown     PetscInt *localVertices, numSend = numCells+numCorners, k = 0;
179552f7358SJed Brown 
180785e854fSJed Brown     ierr = PetscMalloc1(numSend, &localVertices);CHKERRQ(ierr);
181552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
1820298fd71SBarry Smith       PetscInt *closure = NULL;
18302281ff3SMatthew G. Knepley       PetscInt closureSize, value, nC = 0;
184552f7358SJed Brown 
18502281ff3SMatthew G. Knepley       if (label) {
18602281ff3SMatthew G. Knepley         ierr = DMLabelGetValue(label, c, &value);CHKERRQ(ierr);
187552f7358SJed Brown         if (value != 1) continue;
188552f7358SJed Brown       }
189552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
190552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
191552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
192552f7358SJed Brown           const PetscInt gv = gvertex[closure[v] - vStart];
193552f7358SJed Brown           closure[nC++] = gv < 0 ? -(gv+1) : gv;
194552f7358SJed Brown         }
195552f7358SJed Brown       }
196552f7358SJed Brown       corners[numCells++] = nC;
197552f7358SJed Brown       localVertices[k++]  = nC;
198552f7358SJed Brown       for (v = 0; v < nC; ++v, ++k) {
199552f7358SJed Brown         localVertices[k] = closure[v];
200552f7358SJed Brown       }
201552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
202552f7358SJed Brown     }
203552f7358SJed Brown     if (k != numSend) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid number of vertices to send %d should be %d", k, numSend);
204552f7358SJed Brown     ierr = MPI_Send(&numSend, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
205552f7358SJed Brown     ierr = MPI_Send(localVertices, numSend, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
206552f7358SJed Brown     ierr = PetscFree(localVertices);CHKERRQ(ierr);
207552f7358SJed Brown   }
208552f7358SJed Brown   ierr = ISRestoreIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr);
209552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "CELL_TYPES %d\n", totCells);CHKERRQ(ierr);
210552f7358SJed Brown   if (!rank) {
211552f7358SJed Brown     PetscInt cellType;
212552f7358SJed Brown 
213552f7358SJed Brown     for (c = 0; c < numCells; ++c) {
214552f7358SJed Brown       ierr = DMPlexVTKGetCellType(dm, dim, corners[c], &cellType);CHKERRQ(ierr);
215552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "%d\n", cellType);CHKERRQ(ierr);
216552f7358SJed Brown     }
217552f7358SJed Brown     for (proc = 1; proc < numProcs; ++proc) {
218552f7358SJed Brown       MPI_Status status;
219552f7358SJed Brown 
220552f7358SJed Brown       ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
221552f7358SJed Brown       ierr = MPI_Recv(corners, numCells, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
222552f7358SJed Brown       for (c = 0; c < numCells; ++c) {
223552f7358SJed Brown         ierr = DMPlexVTKGetCellType(dm, dim, corners[c], &cellType);CHKERRQ(ierr);
224552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "%d\n", cellType);CHKERRQ(ierr);
225552f7358SJed Brown       }
226552f7358SJed Brown     }
227552f7358SJed Brown   } else {
228552f7358SJed Brown     ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
229552f7358SJed Brown     ierr = MPI_Send(corners, numCells, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
230552f7358SJed Brown   }
231552f7358SJed Brown   ierr        = PetscFree(corners);CHKERRQ(ierr);
232552f7358SJed Brown   *totalCells = totCells;
233552f7358SJed Brown   PetscFunctionReturn(0);
234552f7358SJed Brown }
235552f7358SJed Brown 
236552f7358SJed Brown #undef __FUNCT__
237e8f6f0f6SMatthew G. Knepley #define __FUNCT__ "DMPlexVTKWritePartition_ASCII"
238e8f6f0f6SMatthew G. Knepley PetscErrorCode DMPlexVTKWritePartition_ASCII(DM dm, FILE *fp)
239e8f6f0f6SMatthew G. Knepley {
240e8f6f0f6SMatthew G. Knepley   MPI_Comm       comm;
241e8f6f0f6SMatthew G. Knepley   PetscInt       numCells = 0, cellHeight;
242e8f6f0f6SMatthew G. Knepley   PetscInt       numLabelCells, cMax, cStart, cEnd, c;
243e8f6f0f6SMatthew G. Knepley   PetscMPIInt    numProcs, rank, proc, tag;
244e8f6f0f6SMatthew G. Knepley   PetscBool      hasLabel;
245e8f6f0f6SMatthew G. Knepley   PetscErrorCode ierr;
246e8f6f0f6SMatthew G. Knepley 
247e8f6f0f6SMatthew G. Knepley   PetscFunctionBegin;
248e8f6f0f6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
249e8f6f0f6SMatthew G. Knepley   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
250e8f6f0f6SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
251e8f6f0f6SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
252e8f6f0f6SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
253e8f6f0f6SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
254e8f6f0f6SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
255e8f6f0f6SMatthew G. Knepley   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
256e8f6f0f6SMatthew G. Knepley   ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
257e8f6f0f6SMatthew G. Knepley   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
258e8f6f0f6SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
259e8f6f0f6SMatthew G. Knepley     if (hasLabel) {
260e8f6f0f6SMatthew G. Knepley       PetscInt value;
261e8f6f0f6SMatthew G. Knepley 
262e8f6f0f6SMatthew G. Knepley       ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
263e8f6f0f6SMatthew G. Knepley       if (value != 1) continue;
264e8f6f0f6SMatthew G. Knepley     }
265e8f6f0f6SMatthew G. Knepley     ++numCells;
266e8f6f0f6SMatthew G. Knepley   }
267e8f6f0f6SMatthew G. Knepley   if (!rank) {
268e8f6f0f6SMatthew G. Knepley     for (c = 0; c < numCells; ++c) {ierr = PetscFPrintf(comm, fp, "%d\n", rank);CHKERRQ(ierr);}
269e8f6f0f6SMatthew G. Knepley     for (proc = 1; proc < numProcs; ++proc) {
270e8f6f0f6SMatthew G. Knepley       MPI_Status status;
271e8f6f0f6SMatthew G. Knepley 
272e8f6f0f6SMatthew G. Knepley       ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
273e8f6f0f6SMatthew G. Knepley       for (c = 0; c < numCells; ++c) {ierr = PetscFPrintf(comm, fp, "%d\n", proc);CHKERRQ(ierr);}
274e8f6f0f6SMatthew G. Knepley     }
275e8f6f0f6SMatthew G. Knepley   } else {
276e8f6f0f6SMatthew G. Knepley     ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
277e8f6f0f6SMatthew G. Knepley   }
278e8f6f0f6SMatthew G. Knepley   PetscFunctionReturn(0);
279e8f6f0f6SMatthew G. Knepley }
280e8f6f0f6SMatthew G. Knepley 
281e8f6f0f6SMatthew G. Knepley #undef __FUNCT__
282552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteSection_ASCII"
2830adebc6cSBarry Smith PetscErrorCode DMPlexVTKWriteSection_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec v, FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale)
2840adebc6cSBarry Smith {
28582f516ccSBarry Smith   MPI_Comm           comm;
286552f7358SJed Brown   const MPI_Datatype mpiType = MPIU_SCALAR;
287552f7358SJed Brown   PetscScalar        *array;
288552f7358SJed Brown   PetscInt           numDof = 0, maxDof;
289552f7358SJed Brown   PetscInt           numLabelCells, cellHeight, cMax, cStart, cEnd, numLabelVertices, vMax, vStart, vEnd, pStart, pEnd, p;
290552f7358SJed Brown   PetscMPIInt        numProcs, rank, proc, tag;
291552f7358SJed Brown   PetscBool          hasLabel;
292552f7358SJed Brown   PetscErrorCode     ierr;
293552f7358SJed Brown 
294552f7358SJed Brown   PetscFunctionBegin;
29582f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
296552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
297552f7358SJed Brown   PetscValidHeaderSpecific(v,VEC_CLASSID,4);
298552f7358SJed Brown   if (precision < 0) precision = 6;
299552f7358SJed Brown   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
300552f7358SJed Brown   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
301552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
302552f7358SJed Brown   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
303552f7358SJed Brown   /* VTK only wants the values at cells or vertices */
304552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
305552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
306552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
3070298fd71SBarry Smith   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, &vMax);CHKERRQ(ierr);
3088865f1eaSKarl Rupp   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
3098865f1eaSKarl Rupp   if (vMax >= 0) vEnd = PetscMin(vEnd, vMax);
310552f7358SJed Brown   pStart   = PetscMax(PetscMin(cStart, vStart), pStart);
311552f7358SJed Brown   pEnd     = PetscMin(PetscMax(cEnd,   vEnd),   pEnd);
312552f7358SJed Brown   ierr     = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
313552f7358SJed Brown   ierr     = DMPlexGetStratumSize(dm, "vtk", 2, &numLabelVertices);CHKERRQ(ierr);
314552f7358SJed Brown   hasLabel = numLabelCells > 0 || numLabelVertices > 0 ? PETSC_TRUE : PETSC_FALSE;
315552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
316552f7358SJed Brown     /* Reject points not either cells or vertices */
317552f7358SJed Brown     if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
318552f7358SJed Brown     if (hasLabel) {
319552f7358SJed Brown       PetscInt value;
320552f7358SJed Brown 
321552f7358SJed Brown       if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
322552f7358SJed Brown           ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
323552f7358SJed Brown         ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
324552f7358SJed Brown         if (value != 1) continue;
325552f7358SJed Brown       }
326552f7358SJed Brown     }
327552f7358SJed Brown     ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr);
3288865f1eaSKarl Rupp     if (numDof) break;
329552f7358SJed Brown   }
330*b2566f29SBarry Smith   ierr = MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
331552f7358SJed Brown   enforceDof = PetscMax(enforceDof, maxDof);
332552f7358SJed Brown   ierr = VecGetArray(v, &array);CHKERRQ(ierr);
333552f7358SJed Brown   if (!rank) {
334552f7358SJed Brown     char formatString[8];
335552f7358SJed Brown 
336552f7358SJed Brown     ierr = PetscSNPrintf(formatString, 8, "%%.%de", precision);CHKERRQ(ierr);
337552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
338552f7358SJed Brown       /* Here we lose a way to filter points by keeping them out of the Numbering */
339552f7358SJed Brown       PetscInt dof, off, goff, d;
340552f7358SJed Brown 
341552f7358SJed Brown       /* Reject points not either cells or vertices */
342552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
343552f7358SJed Brown       if (hasLabel) {
344552f7358SJed Brown         PetscInt value;
345552f7358SJed Brown 
346552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
347552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
348552f7358SJed Brown           ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
349552f7358SJed Brown           if (value != 1) continue;
350552f7358SJed Brown         }
351552f7358SJed Brown       }
352552f7358SJed Brown       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
353552f7358SJed Brown       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
354552f7358SJed Brown       ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
355552f7358SJed Brown       if (dof && goff >= 0) {
356552f7358SJed Brown         for (d = 0; d < dof; d++) {
357552f7358SJed Brown           if (d > 0) {
358552f7358SJed Brown             ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr);
359552f7358SJed Brown           }
360552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(array[off+d])*scale);CHKERRQ(ierr);
361552f7358SJed Brown         }
362552f7358SJed Brown         for (d = dof; d < enforceDof; d++) {
363552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr);
364552f7358SJed Brown         }
365552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
366552f7358SJed Brown       }
367552f7358SJed Brown     }
368552f7358SJed Brown     for (proc = 1; proc < numProcs; ++proc) {
369552f7358SJed Brown       PetscScalar *remoteValues;
370d892089bSMatthew G. Knepley       PetscInt    size = 0, d;
371552f7358SJed Brown       MPI_Status  status;
372552f7358SJed Brown 
373552f7358SJed Brown       ierr = MPI_Recv(&size, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
374785e854fSJed Brown       ierr = PetscMalloc1(size, &remoteValues);CHKERRQ(ierr);
375552f7358SJed Brown       ierr = MPI_Recv(remoteValues, size, mpiType, proc, tag, comm, &status);CHKERRQ(ierr);
376552f7358SJed Brown       for (p = 0; p < size/maxDof; ++p) {
377552f7358SJed Brown         for (d = 0; d < maxDof; ++d) {
378552f7358SJed Brown           if (d > 0) {
379552f7358SJed Brown             ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr);
380552f7358SJed Brown           }
381552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(remoteValues[p*maxDof+d])*scale);CHKERRQ(ierr);
382552f7358SJed Brown         }
383552f7358SJed Brown         for (d = maxDof; d < enforceDof; ++d) {
384552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr);
385552f7358SJed Brown         }
386552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
387552f7358SJed Brown       }
388552f7358SJed Brown       ierr = PetscFree(remoteValues);CHKERRQ(ierr);
389552f7358SJed Brown     }
390552f7358SJed Brown   } else {
391552f7358SJed Brown     PetscScalar *localValues;
392552f7358SJed Brown     PetscInt    size, k = 0;
393552f7358SJed Brown 
394552f7358SJed Brown     ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
395785e854fSJed Brown     ierr = PetscMalloc1(size, &localValues);CHKERRQ(ierr);
396552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
397552f7358SJed Brown       PetscInt dof, off, goff, d;
398552f7358SJed Brown 
399552f7358SJed Brown       /* Reject points not either cells or vertices */
400552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
401552f7358SJed Brown       if (hasLabel) {
402552f7358SJed Brown         PetscInt value;
403552f7358SJed Brown 
404552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
405552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
406552f7358SJed Brown           ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
407552f7358SJed Brown           if (value != 1) continue;
408552f7358SJed Brown         }
409552f7358SJed Brown       }
410552f7358SJed Brown       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
411552f7358SJed Brown       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
412552f7358SJed Brown       ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
413552f7358SJed Brown       if (goff >= 0) {
414552f7358SJed Brown         for (d = 0; d < dof; ++d) {
415552f7358SJed Brown           localValues[k++] = array[off+d];
416552f7358SJed Brown         }
417552f7358SJed Brown       }
418552f7358SJed Brown     }
419552f7358SJed Brown     ierr = MPI_Send(&k, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
420552f7358SJed Brown     ierr = MPI_Send(localValues, k, mpiType, 0, tag, comm);CHKERRQ(ierr);
421552f7358SJed Brown     ierr = PetscFree(localValues);CHKERRQ(ierr);
422552f7358SJed Brown   }
423552f7358SJed Brown   ierr = VecRestoreArray(v, &array);CHKERRQ(ierr);
424552f7358SJed Brown   PetscFunctionReturn(0);
425552f7358SJed Brown }
426552f7358SJed Brown 
427552f7358SJed Brown #undef __FUNCT__
428552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteField_ASCII"
429552f7358SJed Brown PetscErrorCode DMPlexVTKWriteField_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec field, const char name[], FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale)
430552f7358SJed Brown {
43182f516ccSBarry Smith   MPI_Comm       comm;
432552f7358SJed Brown   PetscInt       numDof = 0, maxDof;
433552f7358SJed Brown   PetscInt       pStart, pEnd, p;
434552f7358SJed Brown   PetscErrorCode ierr;
435552f7358SJed Brown 
436552f7358SJed Brown   PetscFunctionBegin;
43782f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
438552f7358SJed Brown   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
439552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
440552f7358SJed Brown     ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr);
4418865f1eaSKarl Rupp     if (numDof) break;
442552f7358SJed Brown   }
443552f7358SJed Brown   numDof = PetscMax(numDof, enforceDof);
444*b2566f29SBarry Smith   ierr = MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
445552f7358SJed Brown   if (!name) name = "Unknown";
446552f7358SJed Brown   if (maxDof == 3) {
447552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "VECTORS %s double\n", name);CHKERRQ(ierr);
448552f7358SJed Brown   } else {
449552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "SCALARS %s double %d\n", name, maxDof);CHKERRQ(ierr);
450552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr);
451552f7358SJed Brown   }
452552f7358SJed Brown   ierr = DMPlexVTKWriteSection_ASCII(dm, section, globalSection, field, fp, enforceDof, precision, scale);CHKERRQ(ierr);
453552f7358SJed Brown   PetscFunctionReturn(0);
454552f7358SJed Brown }
455552f7358SJed Brown 
456552f7358SJed Brown #undef __FUNCT__
457552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll_ASCII"
458552f7358SJed Brown static PetscErrorCode DMPlexVTKWriteAll_ASCII(DM dm, PetscViewer viewer)
459552f7358SJed Brown {
46082f516ccSBarry Smith   MPI_Comm                 comm;
461552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*) viewer->data;
462552f7358SJed Brown   FILE                     *fp;
463552f7358SJed Brown   PetscViewerVTKObjectLink link;
464552f7358SJed Brown   PetscSection             coordSection, globalCoordSection;
465552f7358SJed Brown   PetscLayout              vLayout;
466552f7358SJed Brown   Vec                      coordinates;
467552f7358SJed Brown   PetscReal                lengthScale;
468552f7358SJed Brown   PetscInt                 vMax, totVertices, totCells;
469e8f6f0f6SMatthew G. Knepley   PetscBool                hasPoint = PETSC_FALSE, hasCell = PETSC_FALSE, writePartition = PETSC_FALSE;
470552f7358SJed Brown   PetscErrorCode           ierr;
471552f7358SJed Brown 
472552f7358SJed Brown   PetscFunctionBegin;
47382f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
474552f7358SJed Brown   ierr = PetscFOpen(comm, vtk->filename, "wb", &fp);CHKERRQ(ierr);
475552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "# vtk DataFile Version 2.0\n");CHKERRQ(ierr);
476552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "Simplicial Mesh Example\n");CHKERRQ(ierr);
477552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "ASCII\n");CHKERRQ(ierr);
478552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "DATASET UNSTRUCTURED_GRID\n");CHKERRQ(ierr);
479552f7358SJed Brown   /* Vertices */
480552f7358SJed Brown   ierr = DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);CHKERRQ(ierr);
48169d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
48215b58121SMatthew G. Knepley   ierr = PetscSectionCreateGlobalSection(coordSection, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalCoordSection);CHKERRQ(ierr);
483552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
4840298fd71SBarry Smith   ierr = DMPlexGetHybridBounds(dm, NULL, NULL, NULL, &vMax);CHKERRQ(ierr);
485552f7358SJed Brown   if (vMax >= 0) {
486552f7358SJed Brown     PetscInt pStart, pEnd, p, localSize = 0;
487552f7358SJed Brown 
488552f7358SJed Brown     ierr = PetscSectionGetChart(globalCoordSection, &pStart, &pEnd);CHKERRQ(ierr);
489552f7358SJed Brown     pEnd = PetscMin(pEnd, vMax);
490552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
491552f7358SJed Brown       PetscInt dof;
492552f7358SJed Brown 
493552f7358SJed Brown       ierr = PetscSectionGetDof(globalCoordSection, p, &dof);CHKERRQ(ierr);
4948865f1eaSKarl Rupp       if (dof > 0) ++localSize;
495552f7358SJed Brown     }
49682f516ccSBarry Smith     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &vLayout);CHKERRQ(ierr);
497552f7358SJed Brown     ierr = PetscLayoutSetLocalSize(vLayout, localSize);CHKERRQ(ierr);
498552f7358SJed Brown     ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr);
499552f7358SJed Brown     ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr);
500552f7358SJed Brown   } else {
50182f516ccSBarry Smith     ierr = PetscSectionGetPointLayout(PetscObjectComm((PetscObject)dm), globalCoordSection, &vLayout);CHKERRQ(ierr);
502552f7358SJed Brown   }
503552f7358SJed Brown   ierr = PetscLayoutGetSize(vLayout, &totVertices);CHKERRQ(ierr);
504552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "POINTS %d double\n", totVertices);CHKERRQ(ierr);
505552f7358SJed Brown   ierr = DMPlexVTKWriteSection_ASCII(dm, coordSection, globalCoordSection, coordinates, fp, 3, PETSC_DETERMINE, lengthScale);CHKERRQ(ierr);
506552f7358SJed Brown   /* Cells */
507552f7358SJed Brown   ierr = DMPlexVTKWriteCells_ASCII(dm, fp, &totCells);CHKERRQ(ierr);
508552f7358SJed Brown   /* Vertex fields */
509552f7358SJed Brown   for (link = vtk->link; link; link = link->next) {
510552f7358SJed Brown     if ((link->ft == PETSC_VTK_POINT_FIELD) || (link->ft == PETSC_VTK_POINT_VECTOR_FIELD)) hasPoint = PETSC_TRUE;
511552f7358SJed Brown     if ((link->ft == PETSC_VTK_CELL_FIELD)  || (link->ft == PETSC_VTK_CELL_VECTOR_FIELD))  hasCell  = PETSC_TRUE;
512552f7358SJed Brown   }
513552f7358SJed Brown   if (hasPoint) {
51422d28d08SBarry Smith     ierr = PetscFPrintf(comm, fp, "POINT_DATA %d\n", totVertices);CHKERRQ(ierr);
515552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
516552f7358SJed Brown       Vec          X = (Vec) link->vec;
517552f7358SJed Brown       DM           dmX;
5180298fd71SBarry Smith       PetscSection section, globalSection, newSection = NULL;
519552f7358SJed Brown       const char   *name;
520552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
521552f7358SJed Brown 
522552f7358SJed Brown       if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
523552f7358SJed Brown       if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) enforceDof = 3;
524552f7358SJed Brown       ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr);
525552f7358SJed Brown       ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
526552f7358SJed Brown       if (dmX) {
527088580cfSMatthew G Knepley         DMLabel  subpointMap, subpointMapX;
528839dd189SMatthew G Knepley         PetscInt dim, dimX, pStart, pEnd, qStart, qEnd;
529839dd189SMatthew G Knepley 
53022d28d08SBarry Smith         ierr = DMGetDefaultSection(dmX, &section);CHKERRQ(ierr);
531839dd189SMatthew G Knepley         /* Here is where we check whether dmX is a submesh of dm */
532c73cfb54SMatthew G. Knepley         ierr = DMGetDimension(dm,  &dim);CHKERRQ(ierr);
533c73cfb54SMatthew G. Knepley         ierr = DMGetDimension(dmX, &dimX);CHKERRQ(ierr);
534839dd189SMatthew G Knepley         ierr = DMPlexGetChart(dm,  &pStart, &pEnd);CHKERRQ(ierr);
535839dd189SMatthew G Knepley         ierr = DMPlexGetChart(dmX, &qStart, &qEnd);CHKERRQ(ierr);
536839dd189SMatthew G Knepley         ierr = DMPlexGetSubpointMap(dm,  &subpointMap);CHKERRQ(ierr);
537839dd189SMatthew G Knepley         ierr = DMPlexGetSubpointMap(dmX, &subpointMapX);CHKERRQ(ierr);
538839dd189SMatthew G Knepley         if (((dim != dimX) || ((pEnd-pStart) < (qEnd-qStart))) && subpointMap && !subpointMapX) {
539434e42b5SMatthew G. Knepley           const PetscInt *ind = NULL;
540e6ccafaeSMatthew G Knepley           IS              subpointIS;
541434e42b5SMatthew G. Knepley           PetscInt        n = 0, q;
542839dd189SMatthew G Knepley 
543839dd189SMatthew G Knepley           ierr = PetscSectionGetChart(section, &qStart, &qEnd);CHKERRQ(ierr);
544e6ccafaeSMatthew G Knepley           ierr = DMPlexCreateSubpointIS(dm, &subpointIS);CHKERRQ(ierr);
545434e42b5SMatthew G. Knepley           if (subpointIS) {
546e6ccafaeSMatthew G Knepley             ierr = ISGetLocalSize(subpointIS, &n);CHKERRQ(ierr);
547e6ccafaeSMatthew G Knepley             ierr = ISGetIndices(subpointIS, &ind);CHKERRQ(ierr);
548434e42b5SMatthew G. Knepley           }
549839dd189SMatthew G Knepley           ierr = PetscSectionCreate(comm, &newSection);CHKERRQ(ierr);
550839dd189SMatthew G Knepley           ierr = PetscSectionSetChart(newSection, pStart, pEnd);CHKERRQ(ierr);
551839dd189SMatthew G Knepley           for (q = qStart; q < qEnd; ++q) {
552839dd189SMatthew G Knepley             PetscInt dof, off, p;
553839dd189SMatthew G Knepley 
554839dd189SMatthew G Knepley             ierr = PetscSectionGetDof(section, q, &dof);CHKERRQ(ierr);
555839dd189SMatthew G Knepley             if (dof) {
556839dd189SMatthew G Knepley               ierr = PetscFindInt(q, n, ind, &p);CHKERRQ(ierr);
557839dd189SMatthew G Knepley               if (p >= pStart) {
558839dd189SMatthew G Knepley                 ierr = PetscSectionSetDof(newSection, p, dof);CHKERRQ(ierr);
559839dd189SMatthew G Knepley                 ierr = PetscSectionGetOffset(section, q, &off);CHKERRQ(ierr);
560839dd189SMatthew G Knepley                 ierr = PetscSectionSetOffset(newSection, p, off);CHKERRQ(ierr);
561839dd189SMatthew G Knepley               }
562839dd189SMatthew G Knepley             }
563839dd189SMatthew G Knepley           }
564434e42b5SMatthew G. Knepley           if (subpointIS) {
565e6ccafaeSMatthew G Knepley             ierr = ISRestoreIndices(subpointIS, &ind);CHKERRQ(ierr);
566e6ccafaeSMatthew G Knepley             ierr = ISDestroy(&subpointIS);CHKERRQ(ierr);
567434e42b5SMatthew G. Knepley           }
568839dd189SMatthew G Knepley           /* No need to setup section */
569839dd189SMatthew G Knepley           section = newSection;
570839dd189SMatthew G Knepley         }
571552f7358SJed Brown       } else {
572426ae2f1SMatthew G Knepley         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
573426ae2f1SMatthew G Knepley         if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
574552f7358SJed Brown       }
57582f516ccSBarry Smith       if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
57615b58121SMatthew G. Knepley       ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection);CHKERRQ(ierr);
577552f7358SJed Brown       ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr);
578552f7358SJed Brown       ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr);
579839dd189SMatthew G Knepley       if (newSection) {ierr = PetscSectionDestroy(&newSection);CHKERRQ(ierr);}
580552f7358SJed Brown     }
581552f7358SJed Brown   }
582552f7358SJed Brown   /* Cell Fields */
583e8f6f0f6SMatthew G. Knepley   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_view_partition", &writePartition, NULL);CHKERRQ(ierr);
584e8f6f0f6SMatthew G. Knepley   if (hasCell || writePartition) {
58522d28d08SBarry Smith     ierr = PetscFPrintf(comm, fp, "CELL_DATA %d\n", totCells);CHKERRQ(ierr);
586552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
587552f7358SJed Brown       Vec          X = (Vec) link->vec;
588552f7358SJed Brown       DM           dmX;
589552f7358SJed Brown       PetscSection section, globalSection;
590552f7358SJed Brown       const char   *name;
591552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
592552f7358SJed Brown 
593552f7358SJed Brown       if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
594552f7358SJed Brown       if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) enforceDof = 3;
595552f7358SJed Brown       ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr);
596552f7358SJed Brown       ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
597552f7358SJed Brown       if (dmX) {
59822d28d08SBarry Smith         ierr = DMGetDefaultSection(dmX, &section);CHKERRQ(ierr);
599552f7358SJed Brown       } else {
600552f7358SJed Brown         PetscContainer c;
601552f7358SJed Brown 
602552f7358SJed Brown         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &c);CHKERRQ(ierr);
60382f516ccSBarry Smith         if (!c) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
604552f7358SJed Brown         ierr = PetscContainerGetPointer(c, (void**) &section);CHKERRQ(ierr);
605552f7358SJed Brown       }
60682f516ccSBarry Smith       if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
60715b58121SMatthew G. Knepley       ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection);CHKERRQ(ierr);
608552f7358SJed Brown       ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr);
609552f7358SJed Brown       ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr);
610552f7358SJed Brown     }
611e8f6f0f6SMatthew G. Knepley     if (writePartition) {
612e8f6f0f6SMatthew G. Knepley       ierr = PetscFPrintf(comm, fp, "SCALARS partition int 1\n");CHKERRQ(ierr);
613e8f6f0f6SMatthew G. Knepley       ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr);
614e8f6f0f6SMatthew G. Knepley       ierr = DMPlexVTKWritePartition_ASCII(dm, fp);CHKERRQ(ierr);
615e8f6f0f6SMatthew G. Knepley     }
616552f7358SJed Brown   }
617552f7358SJed Brown   /* Cleanup */
618552f7358SJed Brown   ierr = PetscSectionDestroy(&globalCoordSection);CHKERRQ(ierr);
619552f7358SJed Brown   ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr);
620552f7358SJed Brown   ierr = PetscFClose(comm, fp);CHKERRQ(ierr);
621552f7358SJed Brown   PetscFunctionReturn(0);
622552f7358SJed Brown }
623552f7358SJed Brown 
624552f7358SJed Brown #undef __FUNCT__
625552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll"
626552f7358SJed Brown /*@C
627552f7358SJed Brown   DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
628552f7358SJed Brown 
629552f7358SJed Brown   Collective
630552f7358SJed Brown 
631552f7358SJed Brown   Input Arguments:
632552f7358SJed Brown + odm - The DMPlex specifying the mesh, passed as a PetscObject
633552f7358SJed Brown - viewer - viewer of type VTK
634552f7358SJed Brown 
635552f7358SJed Brown   Level: developer
636552f7358SJed Brown 
637552f7358SJed Brown   Note:
638552f7358SJed Brown   This function is a callback used by the VTK viewer to actually write the file.
639552f7358SJed 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.
640552f7358SJed Brown   Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
641552f7358SJed Brown 
642552f7358SJed Brown .seealso: PETSCVIEWERVTK
643552f7358SJed Brown @*/
644552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer)
645552f7358SJed Brown {
646552f7358SJed Brown   DM             dm = (DM) odm;
647552f7358SJed Brown   PetscBool      isvtk;
648552f7358SJed Brown   PetscErrorCode ierr;
649552f7358SJed Brown 
650552f7358SJed Brown   PetscFunctionBegin;
651552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
652552f7358SJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
653552f7358SJed Brown   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr);
65482f516ccSBarry Smith   if (!isvtk) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name);
655552f7358SJed Brown   switch (viewer->format) {
656552f7358SJed Brown   case PETSC_VIEWER_ASCII_VTK:
657552f7358SJed Brown     ierr = DMPlexVTKWriteAll_ASCII(dm, viewer);CHKERRQ(ierr);
658552f7358SJed Brown     break;
659552f7358SJed Brown   case PETSC_VIEWER_VTK_VTU:
660552f7358SJed Brown     ierr = DMPlexVTKWriteAll_VTU(dm, viewer);CHKERRQ(ierr);
661552f7358SJed Brown     break;
66282f516ccSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]);
663552f7358SJed Brown   }
664552f7358SJed Brown   PetscFunctionReturn(0);
665552f7358SJed Brown }
666