xref: /petsc/src/dm/impls/plex/plexvtk.c (revision 02281ff37adb6085fb26862da565c8baf3221123)
1552f7358SJed Brown #define PETSCDM_DLL
234541f0dSBarry 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;
77*02281ff3SMatthew G. Knepley   DMLabel        label;
78*02281ff3SMatthew 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;
83552f7358SJed Brown   PetscInt       numLabelCells, 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);
92552f7358SJed Brown   ierr = DMPlexGetDimension(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);
98*02281ff3SMatthew G. Knepley   ierr = DMPlexGetLabel(dm, "vtk", &label);CHKERRQ(ierr);
99552f7358SJed Brown   ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
100552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
1010298fd71SBarry Smith     PetscInt *closure = NULL;
102*02281ff3SMatthew G. Knepley     PetscInt closureSize, value;
103552f7358SJed Brown 
104*02281ff3SMatthew G. Knepley     if (label) {
105*02281ff3SMatthew G. Knepley       ierr = DMLabelGetValue(label, c, &value);CHKERRQ(ierr);
106552f7358SJed Brown       if (value != 1) continue;
107552f7358SJed Brown     }
108552f7358SJed Brown     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
109552f7358SJed Brown     for (v = 0; v < closureSize*2; v += 2) {
110552f7358SJed Brown       if ((closure[v] >= vStart) && (closure[v] < vEnd)) ++numCorners;
111552f7358SJed Brown     }
112552f7358SJed Brown     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
113552f7358SJed Brown     ++numCells;
114552f7358SJed Brown   }
115552f7358SJed Brown   maxCells = numCells;
116552f7358SJed Brown   ierr     = MPI_Reduce(&numCells, &totCells, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRQ(ierr);
117552f7358SJed Brown   ierr     = MPI_Reduce(&numCells, &maxCells, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr);
118552f7358SJed Brown   ierr     = MPI_Reduce(&numCorners, &totCorners, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRQ(ierr);
119552f7358SJed Brown   ierr     = MPI_Reduce(&numCorners, &maxCorners, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr);
120552f7358SJed Brown   ierr     = DMPlexGetVertexNumbering(dm, &globalVertexNumbers);CHKERRQ(ierr);
121552f7358SJed Brown   ierr     = ISGetIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr);
122552f7358SJed Brown   ierr     = PetscMalloc(maxCells * sizeof(PetscInt), &corners);CHKERRQ(ierr);
123552f7358SJed Brown   ierr     = PetscFPrintf(comm, fp, "CELLS %d %d\n", totCells, totCorners+totCells);CHKERRQ(ierr);
124552f7358SJed Brown   if (!rank) {
125552f7358SJed Brown     PetscInt *remoteVertices;
1260dcf9c70SMatthew G. Knepley     int      *vertices;
127552f7358SJed Brown 
1280dcf9c70SMatthew G. Knepley     ierr = PetscMalloc(maxCorners * sizeof(int), &vertices);CHKERRQ(ierr);
129552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
1300298fd71SBarry Smith       PetscInt *closure = NULL;
131*02281ff3SMatthew G. Knepley       PetscInt closureSize, value, nC = 0;
132552f7358SJed Brown 
133*02281ff3SMatthew G. Knepley       if (label) {
134*02281ff3SMatthew G. Knepley         ierr = DMLabelGetValue(label, c, &value);CHKERRQ(ierr);
135552f7358SJed Brown         if (value != 1) continue;
136552f7358SJed Brown       }
137552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
138552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
139552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
140552f7358SJed Brown           const PetscInt gv = gvertex[closure[v] - vStart];
1410dcf9c70SMatthew G. Knepley           vertices[nC++] = gv < 0 ? -(gv+1) : gv;
142552f7358SJed Brown         }
143552f7358SJed Brown       }
1440dcf9c70SMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
145552f7358SJed Brown       corners[numCells++] = nC;
146552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "%d ", nC);CHKERRQ(ierr);
1470dcf9c70SMatthew G. Knepley       ierr = DMPlexInvertCell(dim, nC, vertices);CHKERRQ(ierr);
148552f7358SJed Brown       for (v = 0; v < nC; ++v) {
1490dcf9c70SMatthew G. Knepley         ierr = PetscFPrintf(comm, fp, " %d", vertices[v]);CHKERRQ(ierr);
150552f7358SJed Brown       }
151552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
152552f7358SJed Brown     }
153552f7358SJed Brown     if (numProcs > 1) {ierr = PetscMalloc((maxCorners+maxCells) * sizeof(PetscInt), &remoteVertices);CHKERRQ(ierr);}
154552f7358SJed Brown     for (proc = 1; proc < numProcs; ++proc) {
155552f7358SJed Brown       MPI_Status status;
156552f7358SJed Brown 
157552f7358SJed Brown       ierr = MPI_Recv(&numCorners, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
158552f7358SJed Brown       ierr = MPI_Recv(remoteVertices, numCorners, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
159552f7358SJed Brown       for (c = 0; c < numCorners;) {
160552f7358SJed Brown         PetscInt nC = remoteVertices[c++];
161552f7358SJed Brown 
162552f7358SJed Brown         for (v = 0; v < nC; ++v, ++c) {
1630dcf9c70SMatthew G. Knepley           vertices[v] = remoteVertices[c];
1640dcf9c70SMatthew G. Knepley         }
1650dcf9c70SMatthew G. Knepley         ierr = DMPlexInvertCell(dim, nC, vertices);CHKERRQ(ierr);
1660dcf9c70SMatthew G. Knepley         ierr = PetscFPrintf(comm, fp, "%d ", nC);CHKERRQ(ierr);
1670dcf9c70SMatthew G. Knepley         for (v = 0; v < nC; ++v) {
1680dcf9c70SMatthew G. Knepley           ierr = PetscFPrintf(comm, fp, " %d", vertices[v]);CHKERRQ(ierr);
169552f7358SJed Brown         }
170552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
171552f7358SJed Brown       }
172552f7358SJed Brown     }
173552f7358SJed Brown     if (numProcs > 1) {ierr = PetscFree(remoteVertices);CHKERRQ(ierr);}
1740dcf9c70SMatthew G. Knepley     ierr = PetscFree(vertices);CHKERRQ(ierr);
175552f7358SJed Brown   } else {
176552f7358SJed Brown     PetscInt *localVertices, numSend = numCells+numCorners, k = 0;
177552f7358SJed Brown 
178552f7358SJed Brown     ierr = PetscMalloc(numSend * sizeof(PetscInt), &localVertices);CHKERRQ(ierr);
179552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
1800298fd71SBarry Smith       PetscInt *closure = NULL;
181*02281ff3SMatthew G. Knepley       PetscInt closureSize, value, nC = 0;
182552f7358SJed Brown 
183*02281ff3SMatthew G. Knepley       if (label) {
184*02281ff3SMatthew G. Knepley         ierr = DMLabelGetValue(label, c, &value);CHKERRQ(ierr);
185552f7358SJed Brown         if (value != 1) continue;
186552f7358SJed Brown       }
187552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
188552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
189552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
190552f7358SJed Brown           const PetscInt gv = gvertex[closure[v] - vStart];
191552f7358SJed Brown           closure[nC++] = gv < 0 ? -(gv+1) : gv;
192552f7358SJed Brown         }
193552f7358SJed Brown       }
194552f7358SJed Brown       corners[numCells++] = nC;
195552f7358SJed Brown       localVertices[k++]  = nC;
196552f7358SJed Brown       for (v = 0; v < nC; ++v, ++k) {
197552f7358SJed Brown         localVertices[k] = closure[v];
198552f7358SJed Brown       }
199552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
200552f7358SJed Brown     }
201552f7358SJed Brown     if (k != numSend) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid number of vertices to send %d should be %d", k, numSend);
202552f7358SJed Brown     ierr = MPI_Send(&numSend, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
203552f7358SJed Brown     ierr = MPI_Send(localVertices, numSend, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
204552f7358SJed Brown     ierr = PetscFree(localVertices);CHKERRQ(ierr);
205552f7358SJed Brown   }
206552f7358SJed Brown   ierr = ISRestoreIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr);
207552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "CELL_TYPES %d\n", totCells);CHKERRQ(ierr);
208552f7358SJed Brown   if (!rank) {
209552f7358SJed Brown     PetscInt cellType;
210552f7358SJed Brown 
211552f7358SJed Brown     for (c = 0; c < numCells; ++c) {
212552f7358SJed Brown       ierr = DMPlexVTKGetCellType(dm, dim, corners[c], &cellType);CHKERRQ(ierr);
213552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "%d\n", cellType);CHKERRQ(ierr);
214552f7358SJed Brown     }
215552f7358SJed Brown     for (proc = 1; proc < numProcs; ++proc) {
216552f7358SJed Brown       MPI_Status status;
217552f7358SJed Brown 
218552f7358SJed Brown       ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
219552f7358SJed Brown       ierr = MPI_Recv(corners, numCells, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
220552f7358SJed Brown       for (c = 0; c < numCells; ++c) {
221552f7358SJed Brown         ierr = DMPlexVTKGetCellType(dm, dim, corners[c], &cellType);CHKERRQ(ierr);
222552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "%d\n", cellType);CHKERRQ(ierr);
223552f7358SJed Brown       }
224552f7358SJed Brown     }
225552f7358SJed Brown   } else {
226552f7358SJed Brown     ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
227552f7358SJed Brown     ierr = MPI_Send(corners, numCells, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
228552f7358SJed Brown   }
229552f7358SJed Brown   ierr        = PetscFree(corners);CHKERRQ(ierr);
230552f7358SJed Brown   *totalCells = totCells;
231552f7358SJed Brown   PetscFunctionReturn(0);
232552f7358SJed Brown }
233552f7358SJed Brown 
234552f7358SJed Brown #undef __FUNCT__
235e8f6f0f6SMatthew G. Knepley #define __FUNCT__ "DMPlexVTKWritePartition_ASCII"
236e8f6f0f6SMatthew G. Knepley PetscErrorCode DMPlexVTKWritePartition_ASCII(DM dm, FILE *fp)
237e8f6f0f6SMatthew G. Knepley {
238e8f6f0f6SMatthew G. Knepley   MPI_Comm       comm;
239e8f6f0f6SMatthew G. Knepley   PetscInt       numCells = 0, cellHeight;
240e8f6f0f6SMatthew G. Knepley   PetscInt       numLabelCells, cMax, cStart, cEnd, c;
241e8f6f0f6SMatthew G. Knepley   PetscMPIInt    numProcs, rank, proc, tag;
242e8f6f0f6SMatthew G. Knepley   PetscBool      hasLabel;
243e8f6f0f6SMatthew G. Knepley   PetscErrorCode ierr;
244e8f6f0f6SMatthew G. Knepley 
245e8f6f0f6SMatthew G. Knepley   PetscFunctionBegin;
246e8f6f0f6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
247e8f6f0f6SMatthew G. Knepley   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
248e8f6f0f6SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
249e8f6f0f6SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
250e8f6f0f6SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
251e8f6f0f6SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
252e8f6f0f6SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
253e8f6f0f6SMatthew G. Knepley   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
254e8f6f0f6SMatthew G. Knepley   ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
255e8f6f0f6SMatthew G. Knepley   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
256e8f6f0f6SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
257e8f6f0f6SMatthew G. Knepley     if (hasLabel) {
258e8f6f0f6SMatthew G. Knepley       PetscInt value;
259e8f6f0f6SMatthew G. Knepley 
260e8f6f0f6SMatthew G. Knepley       ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
261e8f6f0f6SMatthew G. Knepley       if (value != 1) continue;
262e8f6f0f6SMatthew G. Knepley     }
263e8f6f0f6SMatthew G. Knepley     ++numCells;
264e8f6f0f6SMatthew G. Knepley   }
265e8f6f0f6SMatthew G. Knepley   if (!rank) {
266e8f6f0f6SMatthew G. Knepley     for (c = 0; c < numCells; ++c) {ierr = PetscFPrintf(comm, fp, "%d\n", rank);CHKERRQ(ierr);}
267e8f6f0f6SMatthew G. Knepley     for (proc = 1; proc < numProcs; ++proc) {
268e8f6f0f6SMatthew G. Knepley       MPI_Status status;
269e8f6f0f6SMatthew G. Knepley 
270e8f6f0f6SMatthew G. Knepley       ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
271e8f6f0f6SMatthew G. Knepley       for (c = 0; c < numCells; ++c) {ierr = PetscFPrintf(comm, fp, "%d\n", proc);CHKERRQ(ierr);}
272e8f6f0f6SMatthew G. Knepley     }
273e8f6f0f6SMatthew G. Knepley   } else {
274e8f6f0f6SMatthew G. Knepley     ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
275e8f6f0f6SMatthew G. Knepley   }
276e8f6f0f6SMatthew G. Knepley   PetscFunctionReturn(0);
277e8f6f0f6SMatthew G. Knepley }
278e8f6f0f6SMatthew G. Knepley 
279e8f6f0f6SMatthew G. Knepley #undef __FUNCT__
280552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteSection_ASCII"
2810adebc6cSBarry Smith PetscErrorCode DMPlexVTKWriteSection_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec v, FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale)
2820adebc6cSBarry Smith {
28382f516ccSBarry Smith   MPI_Comm           comm;
284552f7358SJed Brown   const MPI_Datatype mpiType = MPIU_SCALAR;
285552f7358SJed Brown   PetscScalar        *array;
286552f7358SJed Brown   PetscInt           numDof = 0, maxDof;
287552f7358SJed Brown   PetscInt           numLabelCells, cellHeight, cMax, cStart, cEnd, numLabelVertices, vMax, vStart, vEnd, pStart, pEnd, p;
288552f7358SJed Brown   PetscMPIInt        numProcs, rank, proc, tag;
289552f7358SJed Brown   PetscBool          hasLabel;
290552f7358SJed Brown   PetscErrorCode     ierr;
291552f7358SJed Brown 
292552f7358SJed Brown   PetscFunctionBegin;
29382f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
294552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
295552f7358SJed Brown   PetscValidHeaderSpecific(v,VEC_CLASSID,4);
296552f7358SJed Brown   if (precision < 0) precision = 6;
297552f7358SJed Brown   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
298552f7358SJed Brown   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
299552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
300552f7358SJed Brown   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
301552f7358SJed Brown   /* VTK only wants the values at cells or vertices */
302552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
303552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
304552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
3050298fd71SBarry Smith   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, &vMax);CHKERRQ(ierr);
3068865f1eaSKarl Rupp   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
3078865f1eaSKarl Rupp   if (vMax >= 0) vEnd = PetscMin(vEnd, vMax);
308552f7358SJed Brown   pStart   = PetscMax(PetscMin(cStart, vStart), pStart);
309552f7358SJed Brown   pEnd     = PetscMin(PetscMax(cEnd,   vEnd),   pEnd);
310552f7358SJed Brown   ierr     = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
311552f7358SJed Brown   ierr     = DMPlexGetStratumSize(dm, "vtk", 2, &numLabelVertices);CHKERRQ(ierr);
312552f7358SJed Brown   hasLabel = numLabelCells > 0 || numLabelVertices > 0 ? PETSC_TRUE : PETSC_FALSE;
313552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
314552f7358SJed Brown     /* Reject points not either cells or vertices */
315552f7358SJed Brown     if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
316552f7358SJed Brown     if (hasLabel) {
317552f7358SJed Brown       PetscInt value;
318552f7358SJed Brown 
319552f7358SJed Brown       if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
320552f7358SJed Brown           ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
321552f7358SJed Brown         ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
322552f7358SJed Brown         if (value != 1) continue;
323552f7358SJed Brown       }
324552f7358SJed Brown     }
325552f7358SJed Brown     ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr);
3268865f1eaSKarl Rupp     if (numDof) break;
327552f7358SJed Brown   }
328552f7358SJed Brown   ierr = MPI_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
329552f7358SJed Brown   enforceDof = PetscMax(enforceDof, maxDof);
330552f7358SJed Brown   ierr = VecGetArray(v, &array);CHKERRQ(ierr);
331552f7358SJed Brown   if (!rank) {
332552f7358SJed Brown     char formatString[8];
333552f7358SJed Brown 
334552f7358SJed Brown     ierr = PetscSNPrintf(formatString, 8, "%%.%de", precision);CHKERRQ(ierr);
335552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
336552f7358SJed Brown       /* Here we lose a way to filter points by keeping them out of the Numbering */
337552f7358SJed Brown       PetscInt dof, off, goff, d;
338552f7358SJed Brown 
339552f7358SJed Brown       /* Reject points not either cells or vertices */
340552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
341552f7358SJed Brown       if (hasLabel) {
342552f7358SJed Brown         PetscInt value;
343552f7358SJed Brown 
344552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
345552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
346552f7358SJed Brown           ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
347552f7358SJed Brown           if (value != 1) continue;
348552f7358SJed Brown         }
349552f7358SJed Brown       }
350552f7358SJed Brown       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
351552f7358SJed Brown       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
352552f7358SJed Brown       ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
353552f7358SJed Brown       if (dof && goff >= 0) {
354552f7358SJed Brown         for (d = 0; d < dof; d++) {
355552f7358SJed Brown           if (d > 0) {
356552f7358SJed Brown             ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr);
357552f7358SJed Brown           }
358552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(array[off+d])*scale);CHKERRQ(ierr);
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     }
366552f7358SJed Brown     for (proc = 1; proc < numProcs; ++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);
372552f7358SJed Brown       ierr = PetscMalloc(size * sizeof(PetscScalar), &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           }
379552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(remoteValues[p*maxDof+d])*scale);CHKERRQ(ierr);
380552f7358SJed Brown         }
381552f7358SJed Brown         for (d = maxDof; d < enforceDof; ++d) {
382552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr);
383552f7358SJed Brown         }
384552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
385552f7358SJed Brown       }
386552f7358SJed Brown       ierr = PetscFree(remoteValues);CHKERRQ(ierr);
387552f7358SJed Brown     }
388552f7358SJed Brown   } else {
389552f7358SJed Brown     PetscScalar *localValues;
390552f7358SJed Brown     PetscInt    size, k = 0;
391552f7358SJed Brown 
392552f7358SJed Brown     ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
393552f7358SJed Brown     ierr = PetscMalloc(size * sizeof(PetscScalar), &localValues);CHKERRQ(ierr);
394552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
395552f7358SJed Brown       PetscInt dof, off, goff, d;
396552f7358SJed Brown 
397552f7358SJed Brown       /* Reject points not either cells or vertices */
398552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
399552f7358SJed Brown       if (hasLabel) {
400552f7358SJed Brown         PetscInt value;
401552f7358SJed Brown 
402552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
403552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
404552f7358SJed Brown           ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
405552f7358SJed Brown           if (value != 1) continue;
406552f7358SJed Brown         }
407552f7358SJed Brown       }
408552f7358SJed Brown       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
409552f7358SJed Brown       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
410552f7358SJed Brown       ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
411552f7358SJed Brown       if (goff >= 0) {
412552f7358SJed Brown         for (d = 0; d < dof; ++d) {
413552f7358SJed Brown           localValues[k++] = array[off+d];
414552f7358SJed Brown         }
415552f7358SJed Brown       }
416552f7358SJed Brown     }
417552f7358SJed Brown     ierr = MPI_Send(&k, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
418552f7358SJed Brown     ierr = MPI_Send(localValues, k, mpiType, 0, tag, comm);CHKERRQ(ierr);
419552f7358SJed Brown     ierr = PetscFree(localValues);CHKERRQ(ierr);
420552f7358SJed Brown   }
421552f7358SJed Brown   ierr = VecRestoreArray(v, &array);CHKERRQ(ierr);
422552f7358SJed Brown   PetscFunctionReturn(0);
423552f7358SJed Brown }
424552f7358SJed Brown 
425552f7358SJed Brown #undef __FUNCT__
426552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteField_ASCII"
427552f7358SJed Brown PetscErrorCode DMPlexVTKWriteField_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec field, const char name[], FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale)
428552f7358SJed Brown {
42982f516ccSBarry Smith   MPI_Comm       comm;
430552f7358SJed Brown   PetscInt       numDof = 0, maxDof;
431552f7358SJed Brown   PetscInt       pStart, pEnd, p;
432552f7358SJed Brown   PetscErrorCode ierr;
433552f7358SJed Brown 
434552f7358SJed Brown   PetscFunctionBegin;
43582f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
436552f7358SJed Brown   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
437552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
438552f7358SJed Brown     ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr);
4398865f1eaSKarl Rupp     if (numDof) break;
440552f7358SJed Brown   }
441552f7358SJed Brown   numDof = PetscMax(numDof, enforceDof);
44282f516ccSBarry Smith   ierr = MPI_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
443552f7358SJed Brown   if (!name) name = "Unknown";
444552f7358SJed Brown   if (maxDof == 3) {
445552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "VECTORS %s double\n", name);CHKERRQ(ierr);
446552f7358SJed Brown   } else {
447552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "SCALARS %s double %d\n", name, maxDof);CHKERRQ(ierr);
448552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr);
449552f7358SJed Brown   }
450552f7358SJed Brown   ierr = DMPlexVTKWriteSection_ASCII(dm, section, globalSection, field, fp, enforceDof, precision, scale);CHKERRQ(ierr);
451552f7358SJed Brown   PetscFunctionReturn(0);
452552f7358SJed Brown }
453552f7358SJed Brown 
454552f7358SJed Brown #undef __FUNCT__
455552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll_ASCII"
456552f7358SJed Brown static PetscErrorCode DMPlexVTKWriteAll_ASCII(DM dm, PetscViewer viewer)
457552f7358SJed Brown {
45882f516ccSBarry Smith   MPI_Comm                 comm;
459552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*) viewer->data;
460552f7358SJed Brown   FILE                     *fp;
461552f7358SJed Brown   PetscViewerVTKObjectLink link;
462552f7358SJed Brown   PetscSection             coordSection, globalCoordSection;
463552f7358SJed Brown   PetscLayout              vLayout;
464552f7358SJed Brown   Vec                      coordinates;
465552f7358SJed Brown   PetscReal                lengthScale;
466552f7358SJed Brown   PetscInt                 vMax, totVertices, totCells;
467e8f6f0f6SMatthew G. Knepley   PetscBool                hasPoint = PETSC_FALSE, hasCell = PETSC_FALSE, writePartition = PETSC_FALSE;
468552f7358SJed Brown   PetscErrorCode           ierr;
469552f7358SJed Brown 
470552f7358SJed Brown   PetscFunctionBegin;
47182f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
472552f7358SJed Brown   ierr = PetscFOpen(comm, vtk->filename, "wb", &fp);CHKERRQ(ierr);
473552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "# vtk DataFile Version 2.0\n");CHKERRQ(ierr);
474552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "Simplicial Mesh Example\n");CHKERRQ(ierr);
475552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "ASCII\n");CHKERRQ(ierr);
476552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "DATASET UNSTRUCTURED_GRID\n");CHKERRQ(ierr);
477552f7358SJed Brown   /* Vertices */
478552f7358SJed Brown   ierr = DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);CHKERRQ(ierr);
479552f7358SJed Brown   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
480552f7358SJed Brown   ierr = PetscSectionCreateGlobalSection(coordSection, dm->sf, PETSC_FALSE, &globalCoordSection);CHKERRQ(ierr);
481552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
4820298fd71SBarry Smith   ierr = DMPlexGetHybridBounds(dm, NULL, NULL, NULL, &vMax);CHKERRQ(ierr);
483552f7358SJed Brown   if (vMax >= 0) {
484552f7358SJed Brown     PetscInt pStart, pEnd, p, localSize = 0;
485552f7358SJed Brown 
486552f7358SJed Brown     ierr = PetscSectionGetChart(globalCoordSection, &pStart, &pEnd);CHKERRQ(ierr);
487552f7358SJed Brown     pEnd = PetscMin(pEnd, vMax);
488552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
489552f7358SJed Brown       PetscInt dof;
490552f7358SJed Brown 
491552f7358SJed Brown       ierr = PetscSectionGetDof(globalCoordSection, p, &dof);CHKERRQ(ierr);
4928865f1eaSKarl Rupp       if (dof > 0) ++localSize;
493552f7358SJed Brown     }
49482f516ccSBarry Smith     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &vLayout);CHKERRQ(ierr);
495552f7358SJed Brown     ierr = PetscLayoutSetLocalSize(vLayout, localSize);CHKERRQ(ierr);
496552f7358SJed Brown     ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr);
497552f7358SJed Brown     ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr);
498552f7358SJed Brown   } else {
49982f516ccSBarry Smith     ierr = PetscSectionGetPointLayout(PetscObjectComm((PetscObject)dm), globalCoordSection, &vLayout);CHKERRQ(ierr);
500552f7358SJed Brown   }
501552f7358SJed Brown   ierr = PetscLayoutGetSize(vLayout, &totVertices);CHKERRQ(ierr);
502552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "POINTS %d double\n", totVertices);CHKERRQ(ierr);
503552f7358SJed Brown   ierr = DMPlexVTKWriteSection_ASCII(dm, coordSection, globalCoordSection, coordinates, fp, 3, PETSC_DETERMINE, lengthScale);CHKERRQ(ierr);
504552f7358SJed Brown   /* Cells */
505552f7358SJed Brown   ierr = DMPlexVTKWriteCells_ASCII(dm, fp, &totCells);CHKERRQ(ierr);
506552f7358SJed Brown   /* Vertex fields */
507552f7358SJed Brown   for (link = vtk->link; link; link = link->next) {
508552f7358SJed Brown     if ((link->ft == PETSC_VTK_POINT_FIELD) || (link->ft == PETSC_VTK_POINT_VECTOR_FIELD)) hasPoint = PETSC_TRUE;
509552f7358SJed Brown     if ((link->ft == PETSC_VTK_CELL_FIELD)  || (link->ft == PETSC_VTK_CELL_VECTOR_FIELD))  hasCell  = PETSC_TRUE;
510552f7358SJed Brown   }
511552f7358SJed Brown   if (hasPoint) {
51222d28d08SBarry Smith     ierr = PetscFPrintf(comm, fp, "POINT_DATA %d\n", totVertices);CHKERRQ(ierr);
513552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
514552f7358SJed Brown       Vec          X = (Vec) link->vec;
515552f7358SJed Brown       DM           dmX;
5160298fd71SBarry Smith       PetscSection section, globalSection, newSection = NULL;
517552f7358SJed Brown       const char   *name;
518552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
519552f7358SJed Brown 
520552f7358SJed Brown       if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
521552f7358SJed Brown       if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) enforceDof = 3;
522552f7358SJed Brown       ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr);
523552f7358SJed Brown       ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
524552f7358SJed Brown       if (dmX) {
525088580cfSMatthew G Knepley         DMLabel  subpointMap, subpointMapX;
526839dd189SMatthew G Knepley         PetscInt dim, dimX, pStart, pEnd, qStart, qEnd;
527839dd189SMatthew G Knepley 
52822d28d08SBarry Smith         ierr = DMGetDefaultSection(dmX, &section);CHKERRQ(ierr);
529839dd189SMatthew G Knepley         /* Here is where we check whether dmX is a submesh of dm */
530839dd189SMatthew G Knepley         ierr = DMPlexGetDimension(dm,  &dim);CHKERRQ(ierr);
531839dd189SMatthew G Knepley         ierr = DMPlexGetDimension(dmX, &dimX);CHKERRQ(ierr);
532839dd189SMatthew G Knepley         ierr = DMPlexGetChart(dm,  &pStart, &pEnd);CHKERRQ(ierr);
533839dd189SMatthew G Knepley         ierr = DMPlexGetChart(dmX, &qStart, &qEnd);CHKERRQ(ierr);
534839dd189SMatthew G Knepley         ierr = DMPlexGetSubpointMap(dm,  &subpointMap);CHKERRQ(ierr);
535839dd189SMatthew G Knepley         ierr = DMPlexGetSubpointMap(dmX, &subpointMapX);CHKERRQ(ierr);
536839dd189SMatthew G Knepley         if (((dim != dimX) || ((pEnd-pStart) < (qEnd-qStart))) && subpointMap && !subpointMapX) {
537434e42b5SMatthew G. Knepley           const PetscInt *ind = NULL;
538e6ccafaeSMatthew G Knepley           IS              subpointIS;
539434e42b5SMatthew G. Knepley           PetscInt        n = 0, q;
540839dd189SMatthew G Knepley 
541839dd189SMatthew G Knepley           ierr = PetscSectionGetChart(section, &qStart, &qEnd);CHKERRQ(ierr);
542e6ccafaeSMatthew G Knepley           ierr = DMPlexCreateSubpointIS(dm, &subpointIS);CHKERRQ(ierr);
543434e42b5SMatthew G. Knepley           if (subpointIS) {
544e6ccafaeSMatthew G Knepley             ierr = ISGetLocalSize(subpointIS, &n);CHKERRQ(ierr);
545e6ccafaeSMatthew G Knepley             ierr = ISGetIndices(subpointIS, &ind);CHKERRQ(ierr);
546434e42b5SMatthew G. Knepley           }
547839dd189SMatthew G Knepley           ierr = PetscSectionCreate(comm, &newSection);CHKERRQ(ierr);
548839dd189SMatthew G Knepley           ierr = PetscSectionSetChart(newSection, pStart, pEnd);CHKERRQ(ierr);
549839dd189SMatthew G Knepley           for (q = qStart; q < qEnd; ++q) {
550839dd189SMatthew G Knepley             PetscInt dof, off, p;
551839dd189SMatthew G Knepley 
552839dd189SMatthew G Knepley             ierr = PetscSectionGetDof(section, q, &dof);CHKERRQ(ierr);
553839dd189SMatthew G Knepley             if (dof) {
554839dd189SMatthew G Knepley               ierr = PetscFindInt(q, n, ind, &p);CHKERRQ(ierr);
555839dd189SMatthew G Knepley               if (p >= pStart) {
556839dd189SMatthew G Knepley                 ierr = PetscSectionSetDof(newSection, p, dof);CHKERRQ(ierr);
557839dd189SMatthew G Knepley                 ierr = PetscSectionGetOffset(section, q, &off);CHKERRQ(ierr);
558839dd189SMatthew G Knepley                 ierr = PetscSectionSetOffset(newSection, p, off);CHKERRQ(ierr);
559839dd189SMatthew G Knepley               }
560839dd189SMatthew G Knepley             }
561839dd189SMatthew G Knepley           }
562434e42b5SMatthew G. Knepley           if (subpointIS) {
563e6ccafaeSMatthew G Knepley             ierr = ISRestoreIndices(subpointIS, &ind);CHKERRQ(ierr);
564e6ccafaeSMatthew G Knepley             ierr = ISDestroy(&subpointIS);CHKERRQ(ierr);
565434e42b5SMatthew G. Knepley           }
566839dd189SMatthew G Knepley           /* No need to setup section */
567839dd189SMatthew G Knepley           section = newSection;
568839dd189SMatthew G Knepley         }
569552f7358SJed Brown       } else {
570426ae2f1SMatthew G Knepley         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
571426ae2f1SMatthew G Knepley         if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
572552f7358SJed Brown       }
57382f516ccSBarry Smith       if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
574552f7358SJed Brown       ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, &globalSection);CHKERRQ(ierr);
575552f7358SJed Brown       ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr);
576552f7358SJed Brown       ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr);
577839dd189SMatthew G Knepley       if (newSection) {ierr = PetscSectionDestroy(&newSection);CHKERRQ(ierr);}
578552f7358SJed Brown     }
579552f7358SJed Brown   }
580552f7358SJed Brown   /* Cell Fields */
581e8f6f0f6SMatthew G. Knepley   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_view_partition", &writePartition, NULL);CHKERRQ(ierr);
582e8f6f0f6SMatthew G. Knepley   if (hasCell || writePartition) {
58322d28d08SBarry Smith     ierr = PetscFPrintf(comm, fp, "CELL_DATA %d\n", totCells);CHKERRQ(ierr);
584552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
585552f7358SJed Brown       Vec          X = (Vec) link->vec;
586552f7358SJed Brown       DM           dmX;
587552f7358SJed Brown       PetscSection section, globalSection;
588552f7358SJed Brown       const char   *name;
589552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
590552f7358SJed Brown 
591552f7358SJed Brown       if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
592552f7358SJed Brown       if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) enforceDof = 3;
593552f7358SJed Brown       ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr);
594552f7358SJed Brown       ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
595552f7358SJed Brown       if (dmX) {
59622d28d08SBarry Smith         ierr = DMGetDefaultSection(dmX, &section);CHKERRQ(ierr);
597552f7358SJed Brown       } else {
598552f7358SJed Brown         PetscContainer c;
599552f7358SJed Brown 
600552f7358SJed Brown         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &c);CHKERRQ(ierr);
60182f516ccSBarry Smith         if (!c) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
602552f7358SJed Brown         ierr = PetscContainerGetPointer(c, (void**) &section);CHKERRQ(ierr);
603552f7358SJed Brown       }
60482f516ccSBarry Smith       if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
605552f7358SJed Brown       ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, &globalSection);CHKERRQ(ierr);
606552f7358SJed Brown       ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr);
607552f7358SJed Brown       ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr);
608552f7358SJed Brown     }
609e8f6f0f6SMatthew G. Knepley     if (writePartition) {
610e8f6f0f6SMatthew G. Knepley       ierr = PetscFPrintf(comm, fp, "SCALARS partition int 1\n");CHKERRQ(ierr);
611e8f6f0f6SMatthew G. Knepley       ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr);
612e8f6f0f6SMatthew G. Knepley       ierr = DMPlexVTKWritePartition_ASCII(dm, fp);CHKERRQ(ierr);
613e8f6f0f6SMatthew G. Knepley     }
614552f7358SJed Brown   }
615552f7358SJed Brown   /* Cleanup */
616552f7358SJed Brown   ierr = PetscSectionDestroy(&globalCoordSection);CHKERRQ(ierr);
617552f7358SJed Brown   ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr);
618552f7358SJed Brown   ierr = PetscFClose(comm, fp);CHKERRQ(ierr);
619552f7358SJed Brown   PetscFunctionReturn(0);
620552f7358SJed Brown }
621552f7358SJed Brown 
622552f7358SJed Brown #undef __FUNCT__
623552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll"
624552f7358SJed Brown /*@C
625552f7358SJed Brown   DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
626552f7358SJed Brown 
627552f7358SJed Brown   Collective
628552f7358SJed Brown 
629552f7358SJed Brown   Input Arguments:
630552f7358SJed Brown + odm - The DMPlex specifying the mesh, passed as a PetscObject
631552f7358SJed Brown - viewer - viewer of type VTK
632552f7358SJed Brown 
633552f7358SJed Brown   Level: developer
634552f7358SJed Brown 
635552f7358SJed Brown   Note:
636552f7358SJed Brown   This function is a callback used by the VTK viewer to actually write the file.
637552f7358SJed 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.
638552f7358SJed Brown   Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
639552f7358SJed Brown 
640552f7358SJed Brown .seealso: PETSCVIEWERVTK
641552f7358SJed Brown @*/
642552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer)
643552f7358SJed Brown {
644552f7358SJed Brown   DM             dm = (DM) odm;
645552f7358SJed Brown   PetscBool      isvtk;
646552f7358SJed Brown   PetscErrorCode ierr;
647552f7358SJed Brown 
648552f7358SJed Brown   PetscFunctionBegin;
649552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
650552f7358SJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
651552f7358SJed Brown   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr);
65282f516ccSBarry Smith   if (!isvtk) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name);
653552f7358SJed Brown   switch (viewer->format) {
654552f7358SJed Brown   case PETSC_VIEWER_ASCII_VTK:
655552f7358SJed Brown     ierr = DMPlexVTKWriteAll_ASCII(dm, viewer);CHKERRQ(ierr);
656552f7358SJed Brown     break;
657552f7358SJed Brown   case PETSC_VIEWER_VTK_VTU:
658552f7358SJed Brown     ierr = DMPlexVTKWriteAll_VTU(dm, viewer);CHKERRQ(ierr);
659552f7358SJed Brown     break;
66082f516ccSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]);
661552f7358SJed Brown   }
662552f7358SJed Brown   PetscFunctionReturn(0);
663552f7358SJed Brown }
664