xref: /petsc/src/dm/impls/plex/plexvtk.c (revision e630c359805a7909577b2c5515bf56b45bba1d70)
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 
5de0a4605SMatthew G. Knepley PetscErrorCode DMPlexVTKGetCellType_Internal(DM dm, PetscInt dim, PetscInt corners, PetscInt *cellType)
60adebc6cSBarry Smith {
7552f7358SJed Brown   PetscFunctionBegin;
8552f7358SJed Brown   *cellType = -1;
9552f7358SJed Brown   switch (dim) {
10552f7358SJed Brown   case 0:
11552f7358SJed Brown     switch (corners) {
12552f7358SJed Brown     case 1:
13552f7358SJed Brown       *cellType = 1; /* VTK_VERTEX */
14552f7358SJed Brown       break;
15552f7358SJed Brown     default:
16552f7358SJed Brown       break;
17552f7358SJed Brown     }
18552f7358SJed Brown     break;
19552f7358SJed Brown   case 1:
20552f7358SJed Brown     switch (corners) {
21552f7358SJed Brown     case 2:
22552f7358SJed Brown       *cellType = 3; /* VTK_LINE */
23552f7358SJed Brown       break;
24552f7358SJed Brown     case 3:
25552f7358SJed Brown       *cellType = 21; /* VTK_QUADRATIC_EDGE */
26552f7358SJed Brown       break;
27552f7358SJed Brown     default:
28552f7358SJed Brown       break;
29552f7358SJed Brown     }
30552f7358SJed Brown     break;
31552f7358SJed Brown   case 2:
32552f7358SJed Brown     switch (corners) {
33552f7358SJed Brown     case 3:
34552f7358SJed Brown       *cellType = 5; /* VTK_TRIANGLE */
35552f7358SJed Brown       break;
36552f7358SJed Brown     case 4:
37552f7358SJed Brown       *cellType = 9; /* VTK_QUAD */
38552f7358SJed Brown       break;
39552f7358SJed Brown     case 6:
40552f7358SJed Brown       *cellType = 22; /* VTK_QUADRATIC_TRIANGLE */
41552f7358SJed Brown       break;
42552f7358SJed Brown     case 9:
43552f7358SJed Brown       *cellType = 23; /* VTK_QUADRATIC_QUAD */
44552f7358SJed Brown       break;
45552f7358SJed Brown     default:
46552f7358SJed Brown       break;
47552f7358SJed Brown     }
48552f7358SJed Brown     break;
49552f7358SJed Brown   case 3:
50552f7358SJed Brown     switch (corners) {
51552f7358SJed Brown     case 4:
52552f7358SJed Brown       *cellType = 10; /* VTK_TETRA */
53552f7358SJed Brown       break;
542f029394SStefano Zampini     case 6:
552f029394SStefano Zampini       *cellType = 13; /* VTK_WEDGE */
562f029394SStefano Zampini       break;
57552f7358SJed Brown     case 8:
58552f7358SJed Brown       *cellType = 12; /* VTK_HEXAHEDRON */
59552f7358SJed Brown       break;
60552f7358SJed Brown     case 10:
61552f7358SJed Brown       *cellType = 24; /* VTK_QUADRATIC_TETRA */
62552f7358SJed Brown       break;
63552f7358SJed Brown     case 27:
64552f7358SJed Brown       *cellType = 29; /* VTK_QUADRATIC_HEXAHEDRON */
65552f7358SJed Brown       break;
66552f7358SJed Brown     default:
67552f7358SJed Brown       break;
68552f7358SJed Brown     }
69552f7358SJed Brown   }
70552f7358SJed Brown   PetscFunctionReturn(0);
71552f7358SJed Brown }
72552f7358SJed Brown 
737508e736SMatthew G. Knepley static PetscErrorCode DMPlexVTKWriteCells_ASCII(DM dm, FILE *fp, PetscInt *totalCells)
74552f7358SJed Brown {
7582f516ccSBarry Smith   MPI_Comm       comm;
7602281ff3SMatthew G. Knepley   DMLabel        label;
7702281ff3SMatthew G. Knepley   IS             globalVertexNumbers = NULL;
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;
822f029394SStefano Zampini   PetscInt       numLabelCells, maxLabelCells, cStart, cEnd, c, vStart, vEnd, v;
839852e123SBarry Smith   PetscMPIInt    size, rank, proc, tag;
84552f7358SJed Brown   PetscErrorCode ierr;
85552f7358SJed Brown 
86552f7358SJed Brown   PetscFunctionBegin;
8782f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
88552f7358SJed Brown   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
899852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
90552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
91c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(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);
95c58f1c22SToby Isaac   ierr = DMGetLabel(dm, "vtk", &label);CHKERRQ(ierr);
96c58f1c22SToby Isaac   ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
97b2566f29SBarry Smith   ierr = MPIU_Allreduce(&numLabelCells, &maxLabelCells, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
98452b7979SMatthew G. Knepley   if (!maxLabelCells) label = NULL;
99552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
1000298fd71SBarry Smith     PetscInt *closure = NULL;
10102281ff3SMatthew G. Knepley     PetscInt closureSize, value;
102552f7358SJed Brown 
10302281ff3SMatthew G. Knepley     if (label) {
10402281ff3SMatthew G. Knepley       ierr = DMLabelGetValue(label, c, &value);CHKERRQ(ierr);
105552f7358SJed Brown       if (value != 1) continue;
106552f7358SJed Brown     }
107552f7358SJed Brown     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
108552f7358SJed Brown     for (v = 0; v < closureSize*2; v += 2) {
109552f7358SJed Brown       if ((closure[v] >= vStart) && (closure[v] < vEnd)) ++numCorners;
110552f7358SJed Brown     }
111552f7358SJed Brown     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
112552f7358SJed Brown     ++numCells;
113552f7358SJed Brown   }
114552f7358SJed Brown   maxCells = numCells;
115552f7358SJed Brown   ierr     = MPI_Reduce(&numCells, &totCells, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRQ(ierr);
116552f7358SJed Brown   ierr     = MPI_Reduce(&numCells, &maxCells, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr);
117552f7358SJed Brown   ierr     = MPI_Reduce(&numCorners, &totCorners, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRQ(ierr);
118552f7358SJed Brown   ierr     = MPI_Reduce(&numCorners, &maxCorners, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr);
119552f7358SJed Brown   ierr     = DMPlexGetVertexNumbering(dm, &globalVertexNumbers);CHKERRQ(ierr);
120552f7358SJed Brown   ierr     = ISGetIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr);
121785e854fSJed Brown   ierr     = PetscMalloc1(maxCells, &corners);CHKERRQ(ierr);
1226b8858a1SBarry Smith   ierr     = PetscFPrintf(comm, fp, "CELLS %D %D\n", totCells, totCorners+totCells);CHKERRQ(ierr);
123552f7358SJed Brown   if (!rank) {
124552f7358SJed Brown     PetscInt *remoteVertices;
1250dcf9c70SMatthew G. Knepley     int      *vertices;
126552f7358SJed Brown 
127785e854fSJed Brown     ierr = PetscMalloc1(maxCorners, &vertices);CHKERRQ(ierr);
128552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
1290298fd71SBarry Smith       PetscInt *closure = NULL;
13002281ff3SMatthew G. Knepley       PetscInt closureSize, value, nC = 0;
131552f7358SJed Brown 
13202281ff3SMatthew G. Knepley       if (label) {
13302281ff3SMatthew G. Knepley         ierr = DMLabelGetValue(label, c, &value);CHKERRQ(ierr);
134552f7358SJed Brown         if (value != 1) continue;
135552f7358SJed Brown       }
136552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
137552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
138552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
139552f7358SJed Brown           const PetscInt gv = gvertex[closure[v] - vStart];
1400dcf9c70SMatthew G. Knepley           vertices[nC++] = gv < 0 ? -(gv+1) : gv;
141552f7358SJed Brown         }
142552f7358SJed Brown       }
1430dcf9c70SMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
144552f7358SJed Brown       corners[numCells++] = nC;
1456b8858a1SBarry Smith       ierr = PetscFPrintf(comm, fp, "%D ", nC);CHKERRQ(ierr);
1460dcf9c70SMatthew G. Knepley       ierr = DMPlexInvertCell(dim, nC, vertices);CHKERRQ(ierr);
147552f7358SJed Brown       for (v = 0; v < nC; ++v) {
148532e6d4bSStefano Zampini         ierr = PetscFPrintf(comm, fp, " %d", vertices[v]);CHKERRQ(ierr);
149552f7358SJed Brown       }
150552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
151552f7358SJed Brown     }
1529852e123SBarry Smith     if (size > 1) {ierr = PetscMalloc1(maxCorners+maxCells, &remoteVertices);CHKERRQ(ierr);}
1539852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
154552f7358SJed Brown       MPI_Status status;
155552f7358SJed Brown 
156552f7358SJed Brown       ierr = MPI_Recv(&numCorners, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
157552f7358SJed Brown       ierr = MPI_Recv(remoteVertices, numCorners, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
158552f7358SJed Brown       for (c = 0; c < numCorners;) {
159552f7358SJed Brown         PetscInt nC = remoteVertices[c++];
160552f7358SJed Brown 
161552f7358SJed Brown         for (v = 0; v < nC; ++v, ++c) {
1620dcf9c70SMatthew G. Knepley           vertices[v] = remoteVertices[c];
1630dcf9c70SMatthew G. Knepley         }
1640dcf9c70SMatthew G. Knepley         ierr = DMPlexInvertCell(dim, nC, vertices);CHKERRQ(ierr);
1656b8858a1SBarry Smith         ierr = PetscFPrintf(comm, fp, "%D ", nC);CHKERRQ(ierr);
1660dcf9c70SMatthew G. Knepley         for (v = 0; v < nC; ++v) {
167532e6d4bSStefano Zampini           ierr = PetscFPrintf(comm, fp, " %d", vertices[v]);CHKERRQ(ierr);
168552f7358SJed Brown         }
169552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
170552f7358SJed Brown       }
171552f7358SJed Brown     }
1729852e123SBarry Smith     if (size > 1) {ierr = PetscFree(remoteVertices);CHKERRQ(ierr);}
1730dcf9c70SMatthew G. Knepley     ierr = PetscFree(vertices);CHKERRQ(ierr);
174552f7358SJed Brown   } else {
175552f7358SJed Brown     PetscInt *localVertices, numSend = numCells+numCorners, k = 0;
176552f7358SJed Brown 
177785e854fSJed Brown     ierr = PetscMalloc1(numSend, &localVertices);CHKERRQ(ierr);
178552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
1790298fd71SBarry Smith       PetscInt *closure = NULL;
18002281ff3SMatthew G. Knepley       PetscInt closureSize, value, nC = 0;
181552f7358SJed Brown 
18202281ff3SMatthew G. Knepley       if (label) {
18302281ff3SMatthew G. Knepley         ierr = DMLabelGetValue(label, c, &value);CHKERRQ(ierr);
184552f7358SJed Brown         if (value != 1) continue;
185552f7358SJed Brown       }
186552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
187552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
188552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
189552f7358SJed Brown           const PetscInt gv = gvertex[closure[v] - vStart];
190552f7358SJed Brown           closure[nC++] = gv < 0 ? -(gv+1) : gv;
191552f7358SJed Brown         }
192552f7358SJed Brown       }
193552f7358SJed Brown       corners[numCells++] = nC;
194552f7358SJed Brown       localVertices[k++]  = nC;
195552f7358SJed Brown       for (v = 0; v < nC; ++v, ++k) {
196552f7358SJed Brown         localVertices[k] = closure[v];
197552f7358SJed Brown       }
198552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
199552f7358SJed Brown     }
2006b8858a1SBarry Smith     if (k != numSend) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid number of vertices to send %D should be %D", k, numSend);
201552f7358SJed Brown     ierr = MPI_Send(&numSend, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
202552f7358SJed Brown     ierr = MPI_Send(localVertices, numSend, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
203552f7358SJed Brown     ierr = PetscFree(localVertices);CHKERRQ(ierr);
204552f7358SJed Brown   }
205552f7358SJed Brown   ierr = ISRestoreIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr);
2066b8858a1SBarry Smith   ierr = PetscFPrintf(comm, fp, "CELL_TYPES %D\n", totCells);CHKERRQ(ierr);
207552f7358SJed Brown   if (!rank) {
208552f7358SJed Brown     PetscInt cellType;
209552f7358SJed Brown 
210552f7358SJed Brown     for (c = 0; c < numCells; ++c) {
211de0a4605SMatthew G. Knepley       ierr = DMPlexVTKGetCellType_Internal(dm, dim, corners[c], &cellType);CHKERRQ(ierr);
2126b8858a1SBarry Smith       ierr = PetscFPrintf(comm, fp, "%D\n", cellType);CHKERRQ(ierr);
213552f7358SJed Brown     }
2149852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
215552f7358SJed Brown       MPI_Status status;
216552f7358SJed Brown 
217552f7358SJed Brown       ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
218552f7358SJed Brown       ierr = MPI_Recv(corners, numCells, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
219552f7358SJed Brown       for (c = 0; c < numCells; ++c) {
220de0a4605SMatthew G. Knepley         ierr = DMPlexVTKGetCellType_Internal(dm, dim, corners[c], &cellType);CHKERRQ(ierr);
2216b8858a1SBarry Smith         ierr = PetscFPrintf(comm, fp, "%D\n", cellType);CHKERRQ(ierr);
222552f7358SJed Brown       }
223552f7358SJed Brown     }
224552f7358SJed Brown   } else {
225552f7358SJed Brown     ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
226552f7358SJed Brown     ierr = MPI_Send(corners, numCells, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
227552f7358SJed Brown   }
228552f7358SJed Brown   ierr        = PetscFree(corners);CHKERRQ(ierr);
229552f7358SJed Brown   *totalCells = totCells;
230552f7358SJed Brown   PetscFunctionReturn(0);
231552f7358SJed Brown }
232552f7358SJed Brown 
2337508e736SMatthew G. Knepley static PetscErrorCode DMPlexVTKWritePartition_ASCII(DM dm, FILE *fp)
234e8f6f0f6SMatthew G. Knepley {
235e8f6f0f6SMatthew G. Knepley   MPI_Comm       comm;
236e8f6f0f6SMatthew G. Knepley   PetscInt       numCells = 0, cellHeight;
2372f029394SStefano Zampini   PetscInt       numLabelCells, cStart, cEnd, c;
2389852e123SBarry Smith   PetscMPIInt    size, rank, proc, tag;
239e8f6f0f6SMatthew G. Knepley   PetscBool      hasLabel;
240e8f6f0f6SMatthew G. Knepley   PetscErrorCode ierr;
241e8f6f0f6SMatthew G. Knepley 
242e8f6f0f6SMatthew G. Knepley   PetscFunctionBegin;
243e8f6f0f6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
244e8f6f0f6SMatthew G. Knepley   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
2459852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
246e8f6f0f6SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
247e8f6f0f6SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
248e8f6f0f6SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
249c58f1c22SToby Isaac   ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
250e8f6f0f6SMatthew G. Knepley   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
251e8f6f0f6SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
252e8f6f0f6SMatthew G. Knepley     if (hasLabel) {
253e8f6f0f6SMatthew G. Knepley       PetscInt value;
254e8f6f0f6SMatthew G. Knepley 
255c58f1c22SToby Isaac       ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
256e8f6f0f6SMatthew G. Knepley       if (value != 1) continue;
257e8f6f0f6SMatthew G. Knepley     }
258e8f6f0f6SMatthew G. Knepley     ++numCells;
259e8f6f0f6SMatthew G. Knepley   }
260e8f6f0f6SMatthew G. Knepley   if (!rank) {
261e8f6f0f6SMatthew G. Knepley     for (c = 0; c < numCells; ++c) {ierr = PetscFPrintf(comm, fp, "%d\n", rank);CHKERRQ(ierr);}
2629852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
263e8f6f0f6SMatthew G. Knepley       MPI_Status status;
264e8f6f0f6SMatthew G. Knepley 
265e8f6f0f6SMatthew G. Knepley       ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
266e8f6f0f6SMatthew G. Knepley       for (c = 0; c < numCells; ++c) {ierr = PetscFPrintf(comm, fp, "%d\n", proc);CHKERRQ(ierr);}
267e8f6f0f6SMatthew G. Knepley     }
268e8f6f0f6SMatthew G. Knepley   } else {
269e8f6f0f6SMatthew G. Knepley     ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
270e8f6f0f6SMatthew G. Knepley   }
271e8f6f0f6SMatthew G. Knepley   PetscFunctionReturn(0);
272e8f6f0f6SMatthew G. Knepley }
273e8f6f0f6SMatthew G. Knepley 
2747508e736SMatthew G. Knepley static PetscErrorCode DMPlexVTKWriteSection_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec v, FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale)
2750adebc6cSBarry Smith {
27682f516ccSBarry Smith   MPI_Comm           comm;
277552f7358SJed Brown   const MPI_Datatype mpiType = MPIU_SCALAR;
278552f7358SJed Brown   PetscScalar        *array;
279552f7358SJed Brown   PetscInt           numDof = 0, maxDof;
2802f029394SStefano Zampini   PetscInt           numLabelCells, cellHeight, cStart, cEnd, numLabelVertices, vMax, vStart, vEnd, pStart, pEnd, p;
2819852e123SBarry Smith   PetscMPIInt        size, rank, proc, tag;
282552f7358SJed Brown   PetscBool          hasLabel;
283552f7358SJed Brown   PetscErrorCode     ierr;
284552f7358SJed Brown 
285552f7358SJed Brown   PetscFunctionBegin;
28682f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
287552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
288552f7358SJed Brown   PetscValidHeaderSpecific(v,VEC_CLASSID,4);
289552f7358SJed Brown   if (precision < 0) precision = 6;
290552f7358SJed Brown   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
2919852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
292552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
293552f7358SJed Brown   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
294552f7358SJed Brown   /* VTK only wants the values at cells or vertices */
295552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
296552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
297552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
2982f029394SStefano Zampini   ierr = DMPlexGetHybridBounds(dm, NULL, NULL, NULL, &vMax);CHKERRQ(ierr);
2998865f1eaSKarl Rupp   if (vMax >= 0) vEnd = PetscMin(vEnd, vMax);
300552f7358SJed Brown   pStart   = PetscMax(PetscMin(cStart, vStart), pStart);
301552f7358SJed Brown   pEnd     = PetscMin(PetscMax(cEnd,   vEnd),   pEnd);
302c58f1c22SToby Isaac   ierr     = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
303c58f1c22SToby Isaac   ierr     = DMGetStratumSize(dm, "vtk", 2, &numLabelVertices);CHKERRQ(ierr);
304552f7358SJed Brown   hasLabel = numLabelCells > 0 || numLabelVertices > 0 ? PETSC_TRUE : PETSC_FALSE;
305552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
306552f7358SJed Brown     /* Reject points not either cells or vertices */
307552f7358SJed Brown     if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
308552f7358SJed Brown     if (hasLabel) {
309552f7358SJed Brown       PetscInt value;
310552f7358SJed Brown 
311552f7358SJed Brown       if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
312552f7358SJed Brown           ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
313c58f1c22SToby Isaac         ierr = DMGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
314552f7358SJed Brown         if (value != 1) continue;
315552f7358SJed Brown       }
316552f7358SJed Brown     }
317552f7358SJed Brown     ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr);
3188865f1eaSKarl Rupp     if (numDof) break;
319552f7358SJed Brown   }
320b2566f29SBarry Smith   ierr = MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
321552f7358SJed Brown   enforceDof = PetscMax(enforceDof, maxDof);
322552f7358SJed Brown   ierr = VecGetArray(v, &array);CHKERRQ(ierr);
323552f7358SJed Brown   if (!rank) {
3242dea1156SStefano Zampini #if defined(PETSC_USE_REAL___FLOAT128)
3252dea1156SStefano Zampini     double dval;
3262dea1156SStefano Zampini #endif
327552f7358SJed Brown     char formatString[8];
328552f7358SJed Brown 
329552f7358SJed Brown     ierr = PetscSNPrintf(formatString, 8, "%%.%de", precision);CHKERRQ(ierr);
330552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
331552f7358SJed Brown       /* Here we lose a way to filter points by keeping them out of the Numbering */
332552f7358SJed Brown       PetscInt dof, off, goff, d;
333552f7358SJed Brown 
334552f7358SJed Brown       /* Reject points not either cells or vertices */
335552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
336552f7358SJed Brown       if (hasLabel) {
337552f7358SJed Brown         PetscInt value;
338552f7358SJed Brown 
339552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
340552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
341c58f1c22SToby Isaac           ierr = DMGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
342552f7358SJed Brown           if (value != 1) continue;
343552f7358SJed Brown         }
344552f7358SJed Brown       }
345552f7358SJed Brown       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
346552f7358SJed Brown       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
347552f7358SJed Brown       ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
348552f7358SJed Brown       if (dof && goff >= 0) {
349552f7358SJed Brown         for (d = 0; d < dof; d++) {
350552f7358SJed Brown           if (d > 0) {
351552f7358SJed Brown             ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr);
352552f7358SJed Brown           }
3532dea1156SStefano Zampini #if defined(PETSC_USE_REAL___FLOAT128)
3542dea1156SStefano Zampini           dval = (double)PetscRealPart(array[off+d])*scale;
3552dea1156SStefano Zampini           ierr = PetscFPrintf(comm, fp, formatString, dval);CHKERRQ(ierr);
3562dea1156SStefano Zampini #else
357552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(array[off+d])*scale);CHKERRQ(ierr);
3582dea1156SStefano Zampini #endif
359552f7358SJed Brown         }
360552f7358SJed Brown         for (d = dof; d < enforceDof; d++) {
361552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr);
362552f7358SJed Brown         }
363552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
364552f7358SJed Brown       }
365552f7358SJed Brown     }
3669852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
367552f7358SJed Brown       PetscScalar *remoteValues;
368d892089bSMatthew G. Knepley       PetscInt    size = 0, d;
369552f7358SJed Brown       MPI_Status  status;
370552f7358SJed Brown 
371552f7358SJed Brown       ierr = MPI_Recv(&size, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
372785e854fSJed Brown       ierr = PetscMalloc1(size, &remoteValues);CHKERRQ(ierr);
373552f7358SJed Brown       ierr = MPI_Recv(remoteValues, size, mpiType, proc, tag, comm, &status);CHKERRQ(ierr);
374552f7358SJed Brown       for (p = 0; p < size/maxDof; ++p) {
375552f7358SJed Brown         for (d = 0; d < maxDof; ++d) {
376552f7358SJed Brown           if (d > 0) {
377552f7358SJed Brown             ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr);
378552f7358SJed Brown           }
3792dea1156SStefano Zampini #if defined(PETSC_USE_REAL___FLOAT128)
3802dea1156SStefano Zampini           dval = PetscRealPart(remoteValues[p*maxDof+d])*scale;
3812dea1156SStefano Zampini           ierr = PetscFPrintf(comm, fp, formatString, dval);CHKERRQ(ierr);
3822dea1156SStefano Zampini #else
383552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(remoteValues[p*maxDof+d])*scale);CHKERRQ(ierr);
3842dea1156SStefano Zampini #endif
385552f7358SJed Brown         }
386552f7358SJed Brown         for (d = maxDof; d < enforceDof; ++d) {
387552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr);
388552f7358SJed Brown         }
389552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
390552f7358SJed Brown       }
391552f7358SJed Brown       ierr = PetscFree(remoteValues);CHKERRQ(ierr);
392552f7358SJed Brown     }
393552f7358SJed Brown   } else {
394552f7358SJed Brown     PetscScalar *localValues;
395552f7358SJed Brown     PetscInt    size, k = 0;
396552f7358SJed Brown 
397552f7358SJed Brown     ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
398785e854fSJed Brown     ierr = PetscMalloc1(size, &localValues);CHKERRQ(ierr);
399552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
400552f7358SJed Brown       PetscInt dof, off, goff, d;
401552f7358SJed Brown 
402552f7358SJed Brown       /* Reject points not either cells or vertices */
403552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
404552f7358SJed Brown       if (hasLabel) {
405552f7358SJed Brown         PetscInt value;
406552f7358SJed Brown 
407552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
408552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
409c58f1c22SToby Isaac           ierr = DMGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
410552f7358SJed Brown           if (value != 1) continue;
411552f7358SJed Brown         }
412552f7358SJed Brown       }
413552f7358SJed Brown       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
414552f7358SJed Brown       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
415552f7358SJed Brown       ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
416552f7358SJed Brown       if (goff >= 0) {
417552f7358SJed Brown         for (d = 0; d < dof; ++d) {
418552f7358SJed Brown           localValues[k++] = array[off+d];
419552f7358SJed Brown         }
420552f7358SJed Brown       }
421552f7358SJed Brown     }
422552f7358SJed Brown     ierr = MPI_Send(&k, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
423552f7358SJed Brown     ierr = MPI_Send(localValues, k, mpiType, 0, tag, comm);CHKERRQ(ierr);
424552f7358SJed Brown     ierr = PetscFree(localValues);CHKERRQ(ierr);
425552f7358SJed Brown   }
426552f7358SJed Brown   ierr = VecRestoreArray(v, &array);CHKERRQ(ierr);
427552f7358SJed Brown   PetscFunctionReturn(0);
428552f7358SJed Brown }
429552f7358SJed Brown 
4307508e736SMatthew G. Knepley static PetscErrorCode DMPlexVTKWriteField_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec field, const char name[], FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale)
431552f7358SJed Brown {
43282f516ccSBarry Smith   MPI_Comm       comm;
433552f7358SJed Brown   PetscInt       numDof = 0, maxDof;
434552f7358SJed Brown   PetscInt       pStart, pEnd, p;
435552f7358SJed Brown   PetscErrorCode ierr;
436552f7358SJed Brown 
437552f7358SJed Brown   PetscFunctionBegin;
43882f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
439552f7358SJed Brown   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
440552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
441552f7358SJed Brown     ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr);
4428865f1eaSKarl Rupp     if (numDof) break;
443552f7358SJed Brown   }
444552f7358SJed Brown   numDof = PetscMax(numDof, enforceDof);
445b2566f29SBarry Smith   ierr = MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
446552f7358SJed Brown   if (!name) name = "Unknown";
447552f7358SJed Brown   if (maxDof == 3) {
448552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "VECTORS %s double\n", name);CHKERRQ(ierr);
449552f7358SJed Brown   } else {
4506b8858a1SBarry Smith     ierr = PetscFPrintf(comm, fp, "SCALARS %s double %D\n", name, maxDof);CHKERRQ(ierr);
451552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr);
452552f7358SJed Brown   }
453552f7358SJed Brown   ierr = DMPlexVTKWriteSection_ASCII(dm, section, globalSection, field, fp, enforceDof, precision, scale);CHKERRQ(ierr);
454552f7358SJed Brown   PetscFunctionReturn(0);
455552f7358SJed Brown }
456552f7358SJed Brown 
457552f7358SJed Brown static PetscErrorCode DMPlexVTKWriteAll_ASCII(DM dm, PetscViewer viewer)
458552f7358SJed Brown {
45982f516ccSBarry Smith   MPI_Comm                 comm;
460552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*) viewer->data;
461552f7358SJed Brown   FILE                     *fp;
462552f7358SJed Brown   PetscViewerVTKObjectLink link;
463552f7358SJed Brown   PetscSection             coordSection, globalCoordSection;
464552f7358SJed Brown   PetscLayout              vLayout;
465552f7358SJed Brown   Vec                      coordinates;
466552f7358SJed Brown   PetscReal                lengthScale;
467b6f3c738SMatthew G. Knepley   PetscInt                 vMax, totVertices, totCells = 0;
468d72c04d6SStefano Zampini   PetscBool                hasPoint = PETSC_FALSE, hasCell = PETSC_FALSE, writePartition = PETSC_FALSE, localized;
469552f7358SJed Brown   PetscErrorCode           ierr;
470552f7358SJed Brown 
471552f7358SJed Brown   PetscFunctionBegin;
472d72c04d6SStefano Zampini   ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr);
473d72c04d6SStefano Zampini   if (localized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"VTK output with localized coordinates not yet supported");
47482f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
475552f7358SJed Brown   ierr = PetscFOpen(comm, vtk->filename, "wb", &fp);CHKERRQ(ierr);
476552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "# vtk DataFile Version 2.0\n");CHKERRQ(ierr);
477552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "Simplicial Mesh Example\n");CHKERRQ(ierr);
478552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "ASCII\n");CHKERRQ(ierr);
479552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "DATASET UNSTRUCTURED_GRID\n");CHKERRQ(ierr);
480552f7358SJed Brown   /* Vertices */
481552f7358SJed Brown   ierr = DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);CHKERRQ(ierr);
48269d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
48315b58121SMatthew G. Knepley   ierr = PetscSectionCreateGlobalSection(coordSection, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalCoordSection);CHKERRQ(ierr);
484552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
4850298fd71SBarry Smith   ierr = DMPlexGetHybridBounds(dm, NULL, NULL, NULL, &vMax);CHKERRQ(ierr);
486552f7358SJed Brown   if (vMax >= 0) {
487552f7358SJed Brown     PetscInt pStart, pEnd, p, localSize = 0;
488552f7358SJed Brown 
489552f7358SJed Brown     ierr = PetscSectionGetChart(globalCoordSection, &pStart, &pEnd);CHKERRQ(ierr);
490552f7358SJed Brown     pEnd = PetscMin(pEnd, vMax);
491552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
492552f7358SJed Brown       PetscInt dof;
493552f7358SJed Brown 
494552f7358SJed Brown       ierr = PetscSectionGetDof(globalCoordSection, p, &dof);CHKERRQ(ierr);
4958865f1eaSKarl Rupp       if (dof > 0) ++localSize;
496552f7358SJed Brown     }
49782f516ccSBarry Smith     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &vLayout);CHKERRQ(ierr);
498552f7358SJed Brown     ierr = PetscLayoutSetLocalSize(vLayout, localSize);CHKERRQ(ierr);
499552f7358SJed Brown     ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr);
500552f7358SJed Brown     ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr);
501552f7358SJed Brown   } else {
50282f516ccSBarry Smith     ierr = PetscSectionGetPointLayout(PetscObjectComm((PetscObject)dm), globalCoordSection, &vLayout);CHKERRQ(ierr);
503552f7358SJed Brown   }
504552f7358SJed Brown   ierr = PetscLayoutGetSize(vLayout, &totVertices);CHKERRQ(ierr);
5056b8858a1SBarry Smith   ierr = PetscFPrintf(comm, fp, "POINTS %D double\n", totVertices);CHKERRQ(ierr);
506552f7358SJed Brown   ierr = DMPlexVTKWriteSection_ASCII(dm, coordSection, globalCoordSection, coordinates, fp, 3, PETSC_DETERMINE, lengthScale);CHKERRQ(ierr);
507552f7358SJed Brown   /* Cells */
508552f7358SJed Brown   ierr = DMPlexVTKWriteCells_ASCII(dm, fp, &totCells);CHKERRQ(ierr);
509552f7358SJed Brown   /* Vertex fields */
510552f7358SJed Brown   for (link = vtk->link; link; link = link->next) {
511552f7358SJed Brown     if ((link->ft == PETSC_VTK_POINT_FIELD) || (link->ft == PETSC_VTK_POINT_VECTOR_FIELD)) hasPoint = PETSC_TRUE;
512552f7358SJed Brown     if ((link->ft == PETSC_VTK_CELL_FIELD)  || (link->ft == PETSC_VTK_CELL_VECTOR_FIELD))  hasCell  = PETSC_TRUE;
513552f7358SJed Brown   }
514552f7358SJed Brown   if (hasPoint) {
5156b8858a1SBarry Smith     ierr = PetscFPrintf(comm, fp, "POINT_DATA %D\n", totVertices);CHKERRQ(ierr);
516552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
517552f7358SJed Brown       Vec          X = (Vec) link->vec;
518*e630c359SToby Isaac       PetscSection section = NULL, globalSection, newSection = NULL;
519*e630c359SToby Isaac       char         namebuf[256];
520552f7358SJed Brown       const char   *name;
521552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
522552f7358SJed Brown 
523552f7358SJed Brown       if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
524552f7358SJed Brown       if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) enforceDof = 3;
525552f7358SJed Brown       ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr);
526*e630c359SToby Isaac       ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
527*e630c359SToby Isaac       if (!section) {
528*e630c359SToby Isaac         DM           dmX;
529*e630c359SToby Isaac 
530552f7358SJed Brown         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
531552f7358SJed Brown         if (dmX) {
532088580cfSMatthew G Knepley           DMLabel  subpointMap, subpointMapX;
533839dd189SMatthew G Knepley           PetscInt dim, dimX, pStart, pEnd, qStart, qEnd;
534839dd189SMatthew G Knepley 
53592fd8e1eSJed Brown           ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);
536839dd189SMatthew G Knepley           /* Here is where we check whether dmX is a submesh of dm */
537c73cfb54SMatthew G. Knepley           ierr = DMGetDimension(dm,  &dim);CHKERRQ(ierr);
538c73cfb54SMatthew G. Knepley           ierr = DMGetDimension(dmX, &dimX);CHKERRQ(ierr);
539839dd189SMatthew G Knepley           ierr = DMPlexGetChart(dm,  &pStart, &pEnd);CHKERRQ(ierr);
540839dd189SMatthew G Knepley           ierr = DMPlexGetChart(dmX, &qStart, &qEnd);CHKERRQ(ierr);
541839dd189SMatthew G Knepley           ierr = DMPlexGetSubpointMap(dm,  &subpointMap);CHKERRQ(ierr);
542839dd189SMatthew G Knepley           ierr = DMPlexGetSubpointMap(dmX, &subpointMapX);CHKERRQ(ierr);
543839dd189SMatthew G Knepley           if (((dim != dimX) || ((pEnd-pStart) < (qEnd-qStart))) && subpointMap && !subpointMapX) {
544434e42b5SMatthew G. Knepley             const PetscInt *ind = NULL;
545e6ccafaeSMatthew G Knepley             IS              subpointIS;
546434e42b5SMatthew G. Knepley             PetscInt        n = 0, q;
547839dd189SMatthew G Knepley 
548839dd189SMatthew G Knepley             ierr = PetscSectionGetChart(section, &qStart, &qEnd);CHKERRQ(ierr);
549e6ccafaeSMatthew G Knepley             ierr = DMPlexCreateSubpointIS(dm, &subpointIS);CHKERRQ(ierr);
550434e42b5SMatthew G. Knepley             if (subpointIS) {
551e6ccafaeSMatthew G Knepley               ierr = ISGetLocalSize(subpointIS, &n);CHKERRQ(ierr);
552e6ccafaeSMatthew G Knepley               ierr = ISGetIndices(subpointIS, &ind);CHKERRQ(ierr);
553434e42b5SMatthew G. Knepley             }
554839dd189SMatthew G Knepley             ierr = PetscSectionCreate(comm, &newSection);CHKERRQ(ierr);
555839dd189SMatthew G Knepley             ierr = PetscSectionSetChart(newSection, pStart, pEnd);CHKERRQ(ierr);
556839dd189SMatthew G Knepley             for (q = qStart; q < qEnd; ++q) {
557839dd189SMatthew G Knepley               PetscInt dof, off, p;
558839dd189SMatthew G Knepley 
559839dd189SMatthew G Knepley               ierr = PetscSectionGetDof(section, q, &dof);CHKERRQ(ierr);
560839dd189SMatthew G Knepley               if (dof) {
561839dd189SMatthew G Knepley                 ierr = PetscFindInt(q, n, ind, &p);CHKERRQ(ierr);
562839dd189SMatthew G Knepley                 if (p >= pStart) {
563839dd189SMatthew G Knepley                   ierr = PetscSectionSetDof(newSection, p, dof);CHKERRQ(ierr);
564839dd189SMatthew G Knepley                   ierr = PetscSectionGetOffset(section, q, &off);CHKERRQ(ierr);
565839dd189SMatthew G Knepley                   ierr = PetscSectionSetOffset(newSection, p, off);CHKERRQ(ierr);
566839dd189SMatthew G Knepley                 }
567839dd189SMatthew G Knepley               }
568839dd189SMatthew G Knepley             }
569434e42b5SMatthew G. Knepley             if (subpointIS) {
570e6ccafaeSMatthew G Knepley               ierr = ISRestoreIndices(subpointIS, &ind);CHKERRQ(ierr);
571e6ccafaeSMatthew G Knepley               ierr = ISDestroy(&subpointIS);CHKERRQ(ierr);
572434e42b5SMatthew G. Knepley             }
573839dd189SMatthew G Knepley             /* No need to setup section */
574839dd189SMatthew G Knepley             section = newSection;
575839dd189SMatthew G Knepley           }
576552f7358SJed Brown         }
577*e630c359SToby Isaac       }
578*e630c359SToby Isaac       if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it and could not create one from VecGetDM()", name);
579*e630c359SToby Isaac       if (link->field >= 0) {
580*e630c359SToby Isaac         const char *fieldname;
581*e630c359SToby Isaac 
582*e630c359SToby Isaac         ierr = PetscSectionGetFieldName(section, link->field, &fieldname);CHKERRQ(ierr);
583*e630c359SToby Isaac         ierr = PetscSectionGetField(section, link->field, &section);CHKERRQ(ierr);
584*e630c359SToby Isaac         if (fieldname) {
585*e630c359SToby Isaac           ierr = PetscSNPrintf(namebuf, 255, "%s%s", name, fieldname);CHKERRQ(ierr);
586*e630c359SToby Isaac         } else {
587*e630c359SToby Isaac           ierr = PetscSNPrintf(namebuf, 255, "%s%D", name, link->field);CHKERRQ(ierr);
588*e630c359SToby Isaac         }
589*e630c359SToby Isaac         name = &namebuf[0];
590*e630c359SToby Isaac       }
59115b58121SMatthew G. Knepley       ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection);CHKERRQ(ierr);
592552f7358SJed Brown       ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr);
593552f7358SJed Brown       ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr);
594839dd189SMatthew G Knepley       if (newSection) {ierr = PetscSectionDestroy(&newSection);CHKERRQ(ierr);}
595552f7358SJed Brown     }
596552f7358SJed Brown   }
597552f7358SJed Brown   /* Cell Fields */
598c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_view_partition", &writePartition, NULL);CHKERRQ(ierr);
599e8f6f0f6SMatthew G. Knepley   if (hasCell || writePartition) {
6006b8858a1SBarry Smith     ierr = PetscFPrintf(comm, fp, "CELL_DATA %D\n", totCells);CHKERRQ(ierr);
601552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
602552f7358SJed Brown       Vec          X = (Vec) link->vec;
603*e630c359SToby Isaac       PetscSection section = NULL, globalSection;
604552f7358SJed Brown       const char   *name;
605*e630c359SToby Isaac       char         namebuf[256];
606552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
607552f7358SJed Brown 
608552f7358SJed Brown       if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
609552f7358SJed Brown       if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) enforceDof = 3;
610552f7358SJed Brown       ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr);
611*e630c359SToby Isaac       ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
612*e630c359SToby Isaac       if (!section) {
613*e630c359SToby Isaac         DM           dmX;
614*e630c359SToby Isaac 
615552f7358SJed Brown         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
616552f7358SJed Brown         if (dmX) {
61792fd8e1eSJed Brown           ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);
618552f7358SJed Brown         }
619*e630c359SToby Isaac       }
620*e630c359SToby Isaac       if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it and could not create one from VecGetDM()", name);
621*e630c359SToby Isaac       if (link->field >= 0) {
622*e630c359SToby Isaac         const char *fieldname;
623*e630c359SToby Isaac 
624*e630c359SToby Isaac         ierr = PetscSectionGetFieldName(section, link->field, &fieldname);CHKERRQ(ierr);
625*e630c359SToby Isaac         ierr = PetscSectionGetField(section, link->field, &section);CHKERRQ(ierr);
626*e630c359SToby Isaac         if (fieldname) {
627*e630c359SToby Isaac           ierr = PetscSNPrintf(namebuf, 255, "%s%s", name, fieldname);CHKERRQ(ierr);
628*e630c359SToby Isaac         } else {
629*e630c359SToby Isaac           ierr = PetscSNPrintf(namebuf, 255, "%s%D", name, link->field);CHKERRQ(ierr);
630*e630c359SToby Isaac         }
631*e630c359SToby Isaac         name = &namebuf[0];
632*e630c359SToby Isaac       }
63315b58121SMatthew G. Knepley       ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection);CHKERRQ(ierr);
634552f7358SJed Brown       ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr);
635552f7358SJed Brown       ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr);
636552f7358SJed Brown     }
637e8f6f0f6SMatthew G. Knepley     if (writePartition) {
638e8f6f0f6SMatthew G. Knepley       ierr = PetscFPrintf(comm, fp, "SCALARS partition int 1\n");CHKERRQ(ierr);
639e8f6f0f6SMatthew G. Knepley       ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr);
640e8f6f0f6SMatthew G. Knepley       ierr = DMPlexVTKWritePartition_ASCII(dm, fp);CHKERRQ(ierr);
641e8f6f0f6SMatthew G. Knepley     }
642552f7358SJed Brown   }
643552f7358SJed Brown   /* Cleanup */
644552f7358SJed Brown   ierr = PetscSectionDestroy(&globalCoordSection);CHKERRQ(ierr);
645552f7358SJed Brown   ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr);
646552f7358SJed Brown   ierr = PetscFClose(comm, fp);CHKERRQ(ierr);
647552f7358SJed Brown   PetscFunctionReturn(0);
648552f7358SJed Brown }
649552f7358SJed Brown 
650552f7358SJed Brown /*@C
651552f7358SJed Brown   DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
652552f7358SJed Brown 
653552f7358SJed Brown   Collective
654552f7358SJed Brown 
655552f7358SJed Brown   Input Arguments:
656552f7358SJed Brown + odm - The DMPlex specifying the mesh, passed as a PetscObject
657552f7358SJed Brown - viewer - viewer of type VTK
658552f7358SJed Brown 
659552f7358SJed Brown   Level: developer
660552f7358SJed Brown 
661552f7358SJed Brown   Note:
662552f7358SJed Brown   This function is a callback used by the VTK viewer to actually write the file.
663552f7358SJed 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.
664552f7358SJed Brown   Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
665552f7358SJed Brown 
666552f7358SJed Brown .seealso: PETSCVIEWERVTK
667552f7358SJed Brown @*/
668552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer)
669552f7358SJed Brown {
670552f7358SJed Brown   DM             dm = (DM) odm;
671552f7358SJed Brown   PetscBool      isvtk;
672552f7358SJed Brown   PetscErrorCode ierr;
673552f7358SJed Brown 
674552f7358SJed Brown   PetscFunctionBegin;
675552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
676552f7358SJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
677552f7358SJed Brown   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr);
67882f516ccSBarry Smith   if (!isvtk) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name);
679552f7358SJed Brown   switch (viewer->format) {
680552f7358SJed Brown   case PETSC_VIEWER_ASCII_VTK:
681552f7358SJed Brown     ierr = DMPlexVTKWriteAll_ASCII(dm, viewer);CHKERRQ(ierr);
682552f7358SJed Brown     break;
683552f7358SJed Brown   case PETSC_VIEWER_VTK_VTU:
684552f7358SJed Brown     ierr = DMPlexVTKWriteAll_VTU(dm, viewer);CHKERRQ(ierr);
685552f7358SJed Brown     break;
68682f516ccSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]);
687552f7358SJed Brown   }
688552f7358SJed Brown   PetscFunctionReturn(0);
689552f7358SJed Brown }
690