xref: /petsc/src/dm/impls/plex/plexvtk.c (revision 0dcf9c705d9e3a176e05c9fc6d3f909d0a86dcfd)
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;
77552f7358SJed Brown   IS             globalVertexNumbers;
78552f7358SJed Brown   const PetscInt *gvertex;
79552f7358SJed Brown   PetscInt       dim;
80552f7358SJed Brown   PetscInt       numCorners = 0, totCorners = 0, maxCorners, *corners;
81552f7358SJed Brown   PetscInt       numCells   = 0, totCells   = 0, maxCells, cellHeight;
82552f7358SJed Brown   PetscInt       numLabelCells, cMax, cStart, cEnd, c, vStart, vEnd, v;
83552f7358SJed Brown   PetscMPIInt    numProcs, rank, proc, tag;
84552f7358SJed Brown   PetscBool      hasLabel;
85552f7358SJed Brown   PetscErrorCode ierr;
86552f7358SJed Brown 
87552f7358SJed Brown   PetscFunctionBegin;
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);
98552f7358SJed Brown   ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
998865f1eaSKarl Rupp 
100552f7358SJed Brown   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
101552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
1020298fd71SBarry Smith     PetscInt *closure = NULL;
103552f7358SJed Brown     PetscInt closureSize;
104552f7358SJed Brown 
105552f7358SJed Brown     if (hasLabel) {
106552f7358SJed Brown       PetscInt value;
107552f7358SJed Brown 
108552f7358SJed Brown       ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
109552f7358SJed Brown       if (value != 1) continue;
110552f7358SJed Brown     }
111552f7358SJed Brown     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
112552f7358SJed Brown     for (v = 0; v < closureSize*2; v += 2) {
113552f7358SJed Brown       if ((closure[v] >= vStart) && (closure[v] < vEnd)) ++numCorners;
114552f7358SJed Brown     }
115552f7358SJed Brown     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
116552f7358SJed Brown     ++numCells;
117552f7358SJed Brown   }
118552f7358SJed Brown   maxCells = numCells;
119552f7358SJed Brown   ierr     = MPI_Reduce(&numCells, &totCells, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRQ(ierr);
120552f7358SJed Brown   ierr     = MPI_Reduce(&numCells, &maxCells, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr);
121552f7358SJed Brown   ierr     = MPI_Reduce(&numCorners, &totCorners, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRQ(ierr);
122552f7358SJed Brown   ierr     = MPI_Reduce(&numCorners, &maxCorners, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr);
123552f7358SJed Brown   ierr     = DMPlexGetVertexNumbering(dm, &globalVertexNumbers);CHKERRQ(ierr);
124552f7358SJed Brown   ierr     = ISGetIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr);
125552f7358SJed Brown   ierr     = PetscMalloc(maxCells * sizeof(PetscInt), &corners);CHKERRQ(ierr);
126552f7358SJed Brown   ierr     = PetscFPrintf(comm, fp, "CELLS %d %d\n", totCells, totCorners+totCells);CHKERRQ(ierr);
127552f7358SJed Brown   if (!rank) {
128552f7358SJed Brown     PetscInt *remoteVertices;
129*0dcf9c70SMatthew G. Knepley     int      *vertices;
130552f7358SJed Brown 
131*0dcf9c70SMatthew G. Knepley     ierr = PetscMalloc(maxCorners * sizeof(int), &vertices);CHKERRQ(ierr);
132552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
1330298fd71SBarry Smith       PetscInt *closure = NULL;
134bd09f6e0SJed Brown       PetscInt closureSize, nC = 0;
135552f7358SJed Brown 
136552f7358SJed Brown       if (hasLabel) {
137552f7358SJed Brown         PetscInt value;
138552f7358SJed Brown 
139552f7358SJed Brown         ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
140552f7358SJed Brown         if (value != 1) continue;
141552f7358SJed Brown       }
142552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
143552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
144552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
145552f7358SJed Brown           const PetscInt gv = gvertex[closure[v] - vStart];
146*0dcf9c70SMatthew G. Knepley           vertices[nC++] = gv < 0 ? -(gv+1) : gv;
147552f7358SJed Brown         }
148552f7358SJed Brown       }
149*0dcf9c70SMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
150552f7358SJed Brown       corners[numCells++] = nC;
151552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "%d ", nC);CHKERRQ(ierr);
152*0dcf9c70SMatthew G. Knepley       ierr = DMPlexInvertCell(dim, nC, vertices);CHKERRQ(ierr);
153552f7358SJed Brown       for (v = 0; v < nC; ++v) {
154*0dcf9c70SMatthew G. Knepley         ierr = PetscFPrintf(comm, fp, " %d", vertices[v]);CHKERRQ(ierr);
155552f7358SJed Brown       }
156552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
157552f7358SJed Brown     }
158552f7358SJed Brown     if (numProcs > 1) {ierr = PetscMalloc((maxCorners+maxCells) * sizeof(PetscInt), &remoteVertices);CHKERRQ(ierr);}
159552f7358SJed Brown     for (proc = 1; proc < numProcs; ++proc) {
160552f7358SJed Brown       MPI_Status status;
161552f7358SJed Brown 
162552f7358SJed Brown       ierr = MPI_Recv(&numCorners, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
163552f7358SJed Brown       ierr = MPI_Recv(remoteVertices, numCorners, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
164552f7358SJed Brown       for (c = 0; c < numCorners;) {
165552f7358SJed Brown         PetscInt nC = remoteVertices[c++];
166552f7358SJed Brown 
167552f7358SJed Brown         for (v = 0; v < nC; ++v, ++c) {
168*0dcf9c70SMatthew G. Knepley           vertices[v] = remoteVertices[c];
169*0dcf9c70SMatthew G. Knepley         }
170*0dcf9c70SMatthew G. Knepley         ierr = DMPlexInvertCell(dim, nC, vertices);CHKERRQ(ierr);
171*0dcf9c70SMatthew G. Knepley         ierr = PetscFPrintf(comm, fp, "%d ", nC);CHKERRQ(ierr);
172*0dcf9c70SMatthew G. Knepley         for (v = 0; v < nC; ++v) {
173*0dcf9c70SMatthew G. Knepley           ierr = PetscFPrintf(comm, fp, " %d", vertices[v]);CHKERRQ(ierr);
174552f7358SJed Brown         }
175552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
176552f7358SJed Brown       }
177552f7358SJed Brown     }
178552f7358SJed Brown     if (numProcs > 1) {ierr = PetscFree(remoteVertices);CHKERRQ(ierr);}
179*0dcf9c70SMatthew G. Knepley     ierr = PetscFree(vertices);CHKERRQ(ierr);
180552f7358SJed Brown   } else {
181552f7358SJed Brown     PetscInt *localVertices, numSend = numCells+numCorners, k = 0;
182552f7358SJed Brown 
183552f7358SJed Brown     ierr = PetscMalloc(numSend * sizeof(PetscInt), &localVertices);CHKERRQ(ierr);
184552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
1850298fd71SBarry Smith       PetscInt *closure = NULL;
186552f7358SJed Brown       PetscInt closureSize, nC = 0;
187552f7358SJed Brown 
188552f7358SJed Brown       if (hasLabel) {
189552f7358SJed Brown         PetscInt value;
190552f7358SJed Brown 
191552f7358SJed Brown         ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
192552f7358SJed Brown         if (value != 1) continue;
193552f7358SJed Brown       }
194552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
195552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
196552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
197552f7358SJed Brown           const PetscInt gv = gvertex[closure[v] - vStart];
198552f7358SJed Brown           closure[nC++] = gv < 0 ? -(gv+1) : gv;
199552f7358SJed Brown         }
200552f7358SJed Brown       }
201552f7358SJed Brown       corners[numCells++] = nC;
202552f7358SJed Brown       localVertices[k++]  = nC;
203552f7358SJed Brown       for (v = 0; v < nC; ++v, ++k) {
204552f7358SJed Brown         localVertices[k] = closure[v];
205552f7358SJed Brown       }
206552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
207552f7358SJed Brown     }
208552f7358SJed Brown     if (k != numSend) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid number of vertices to send %d should be %d", k, numSend);
209552f7358SJed Brown     ierr = MPI_Send(&numSend, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
210552f7358SJed Brown     ierr = MPI_Send(localVertices, numSend, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
211552f7358SJed Brown     ierr = PetscFree(localVertices);CHKERRQ(ierr);
212552f7358SJed Brown   }
213552f7358SJed Brown   ierr = ISRestoreIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr);
214552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "CELL_TYPES %d\n", totCells);CHKERRQ(ierr);
215552f7358SJed Brown   if (!rank) {
216552f7358SJed Brown     PetscInt cellType;
217552f7358SJed Brown 
218552f7358SJed Brown     for (c = 0; c < numCells; ++c) {
219552f7358SJed Brown       ierr = DMPlexVTKGetCellType(dm, dim, corners[c], &cellType);CHKERRQ(ierr);
220552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "%d\n", cellType);CHKERRQ(ierr);
221552f7358SJed Brown     }
222552f7358SJed Brown     for (proc = 1; proc < numProcs; ++proc) {
223552f7358SJed Brown       MPI_Status status;
224552f7358SJed Brown 
225552f7358SJed Brown       ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
226552f7358SJed Brown       ierr = MPI_Recv(corners, numCells, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
227552f7358SJed Brown       for (c = 0; c < numCells; ++c) {
228552f7358SJed Brown         ierr = DMPlexVTKGetCellType(dm, dim, corners[c], &cellType);CHKERRQ(ierr);
229552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "%d\n", cellType);CHKERRQ(ierr);
230552f7358SJed Brown       }
231552f7358SJed Brown     }
232552f7358SJed Brown   } else {
233552f7358SJed Brown     ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
234552f7358SJed Brown     ierr = MPI_Send(corners, numCells, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
235552f7358SJed Brown   }
236552f7358SJed Brown   ierr        = PetscFree(corners);CHKERRQ(ierr);
237552f7358SJed Brown   *totalCells = totCells;
238552f7358SJed Brown   PetscFunctionReturn(0);
239552f7358SJed Brown }
240552f7358SJed Brown 
241552f7358SJed Brown #undef __FUNCT__
242e8f6f0f6SMatthew G. Knepley #define __FUNCT__ "DMPlexVTKWritePartition_ASCII"
243e8f6f0f6SMatthew G. Knepley PetscErrorCode DMPlexVTKWritePartition_ASCII(DM dm, FILE *fp)
244e8f6f0f6SMatthew G. Knepley {
245e8f6f0f6SMatthew G. Knepley   MPI_Comm       comm;
246e8f6f0f6SMatthew G. Knepley   PetscInt       numCells = 0, cellHeight;
247e8f6f0f6SMatthew G. Knepley   PetscInt       numLabelCells, cMax, cStart, cEnd, c;
248e8f6f0f6SMatthew G. Knepley   PetscMPIInt    numProcs, rank, proc, tag;
249e8f6f0f6SMatthew G. Knepley   PetscBool      hasLabel;
250e8f6f0f6SMatthew G. Knepley   PetscErrorCode ierr;
251e8f6f0f6SMatthew G. Knepley 
252e8f6f0f6SMatthew G. Knepley   PetscFunctionBegin;
253e8f6f0f6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
254e8f6f0f6SMatthew G. Knepley   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
255e8f6f0f6SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
256e8f6f0f6SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
257e8f6f0f6SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
258e8f6f0f6SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
259e8f6f0f6SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
260e8f6f0f6SMatthew G. Knepley   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
261e8f6f0f6SMatthew G. Knepley   ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
262e8f6f0f6SMatthew G. Knepley   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
263e8f6f0f6SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
264e8f6f0f6SMatthew G. Knepley     if (hasLabel) {
265e8f6f0f6SMatthew G. Knepley       PetscInt value;
266e8f6f0f6SMatthew G. Knepley 
267e8f6f0f6SMatthew G. Knepley       ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
268e8f6f0f6SMatthew G. Knepley       if (value != 1) continue;
269e8f6f0f6SMatthew G. Knepley     }
270e8f6f0f6SMatthew G. Knepley     ++numCells;
271e8f6f0f6SMatthew G. Knepley   }
272e8f6f0f6SMatthew G. Knepley   if (!rank) {
273e8f6f0f6SMatthew G. Knepley     for (c = 0; c < numCells; ++c) {ierr = PetscFPrintf(comm, fp, "%d\n", rank);CHKERRQ(ierr);}
274e8f6f0f6SMatthew G. Knepley     for (proc = 1; proc < numProcs; ++proc) {
275e8f6f0f6SMatthew G. Knepley       MPI_Status status;
276e8f6f0f6SMatthew G. Knepley 
277e8f6f0f6SMatthew G. Knepley       ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
278e8f6f0f6SMatthew G. Knepley       for (c = 0; c < numCells; ++c) {ierr = PetscFPrintf(comm, fp, "%d\n", proc);CHKERRQ(ierr);}
279e8f6f0f6SMatthew G. Knepley     }
280e8f6f0f6SMatthew G. Knepley   } else {
281e8f6f0f6SMatthew G. Knepley     ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
282e8f6f0f6SMatthew G. Knepley   }
283e8f6f0f6SMatthew G. Knepley   PetscFunctionReturn(0);
284e8f6f0f6SMatthew G. Knepley }
285e8f6f0f6SMatthew G. Knepley 
286e8f6f0f6SMatthew G. Knepley #undef __FUNCT__
287552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteSection_ASCII"
2880adebc6cSBarry Smith PetscErrorCode DMPlexVTKWriteSection_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec v, FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale)
2890adebc6cSBarry Smith {
29082f516ccSBarry Smith   MPI_Comm           comm;
291552f7358SJed Brown   const MPI_Datatype mpiType = MPIU_SCALAR;
292552f7358SJed Brown   PetscScalar        *array;
293552f7358SJed Brown   PetscInt           numDof = 0, maxDof;
294552f7358SJed Brown   PetscInt           numLabelCells, cellHeight, cMax, cStart, cEnd, numLabelVertices, vMax, vStart, vEnd, pStart, pEnd, p;
295552f7358SJed Brown   PetscMPIInt        numProcs, rank, proc, tag;
296552f7358SJed Brown   PetscBool          hasLabel;
297552f7358SJed Brown   PetscErrorCode     ierr;
298552f7358SJed Brown 
299552f7358SJed Brown   PetscFunctionBegin;
30082f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
301552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
302552f7358SJed Brown   PetscValidHeaderSpecific(v,VEC_CLASSID,4);
303552f7358SJed Brown   if (precision < 0) precision = 6;
304552f7358SJed Brown   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
305552f7358SJed Brown   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
306552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
307552f7358SJed Brown   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
308552f7358SJed Brown   /* VTK only wants the values at cells or vertices */
309552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
310552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
311552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
3120298fd71SBarry Smith   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, &vMax);CHKERRQ(ierr);
3138865f1eaSKarl Rupp   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
3148865f1eaSKarl Rupp   if (vMax >= 0) vEnd = PetscMin(vEnd, vMax);
315552f7358SJed Brown   pStart   = PetscMax(PetscMin(cStart, vStart), pStart);
316552f7358SJed Brown   pEnd     = PetscMin(PetscMax(cEnd,   vEnd),   pEnd);
317552f7358SJed Brown   ierr     = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
318552f7358SJed Brown   ierr     = DMPlexGetStratumSize(dm, "vtk", 2, &numLabelVertices);CHKERRQ(ierr);
319552f7358SJed Brown   hasLabel = numLabelCells > 0 || numLabelVertices > 0 ? PETSC_TRUE : PETSC_FALSE;
320552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
321552f7358SJed Brown     /* Reject points not either cells or vertices */
322552f7358SJed Brown     if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
323552f7358SJed Brown     if (hasLabel) {
324552f7358SJed Brown       PetscInt value;
325552f7358SJed Brown 
326552f7358SJed Brown       if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
327552f7358SJed Brown           ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
328552f7358SJed Brown         ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
329552f7358SJed Brown         if (value != 1) continue;
330552f7358SJed Brown       }
331552f7358SJed Brown     }
332552f7358SJed Brown     ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr);
3338865f1eaSKarl Rupp     if (numDof) break;
334552f7358SJed Brown   }
335552f7358SJed Brown   ierr = MPI_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
336552f7358SJed Brown   enforceDof = PetscMax(enforceDof, maxDof);
337552f7358SJed Brown   ierr = VecGetArray(v, &array);CHKERRQ(ierr);
338552f7358SJed Brown   if (!rank) {
339552f7358SJed Brown     char formatString[8];
340552f7358SJed Brown 
341552f7358SJed Brown     ierr = PetscSNPrintf(formatString, 8, "%%.%de", precision);CHKERRQ(ierr);
342552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
343552f7358SJed Brown       /* Here we lose a way to filter points by keeping them out of the Numbering */
344552f7358SJed Brown       PetscInt dof, off, goff, d;
345552f7358SJed Brown 
346552f7358SJed Brown       /* Reject points not either cells or vertices */
347552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
348552f7358SJed Brown       if (hasLabel) {
349552f7358SJed Brown         PetscInt value;
350552f7358SJed Brown 
351552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
352552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
353552f7358SJed Brown           ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
354552f7358SJed Brown           if (value != 1) continue;
355552f7358SJed Brown         }
356552f7358SJed Brown       }
357552f7358SJed Brown       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
358552f7358SJed Brown       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
359552f7358SJed Brown       ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
360552f7358SJed Brown       if (dof && goff >= 0) {
361552f7358SJed Brown         for (d = 0; d < dof; d++) {
362552f7358SJed Brown           if (d > 0) {
363552f7358SJed Brown             ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr);
364552f7358SJed Brown           }
365552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(array[off+d])*scale);CHKERRQ(ierr);
366552f7358SJed Brown         }
367552f7358SJed Brown         for (d = dof; d < enforceDof; d++) {
368552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr);
369552f7358SJed Brown         }
370552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
371552f7358SJed Brown       }
372552f7358SJed Brown     }
373552f7358SJed Brown     for (proc = 1; proc < numProcs; ++proc) {
374552f7358SJed Brown       PetscScalar *remoteValues;
375552f7358SJed Brown       PetscInt    size, d;
376552f7358SJed Brown       MPI_Status  status;
377552f7358SJed Brown 
378552f7358SJed Brown       ierr = MPI_Recv(&size, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
379552f7358SJed Brown       ierr = PetscMalloc(size * sizeof(PetscScalar), &remoteValues);CHKERRQ(ierr);
380552f7358SJed Brown       ierr = MPI_Recv(remoteValues, size, mpiType, proc, tag, comm, &status);CHKERRQ(ierr);
381552f7358SJed Brown       for (p = 0; p < size/maxDof; ++p) {
382552f7358SJed Brown         for (d = 0; d < maxDof; ++d) {
383552f7358SJed Brown           if (d > 0) {
384552f7358SJed Brown             ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr);
385552f7358SJed Brown           }
386552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(remoteValues[p*maxDof+d])*scale);CHKERRQ(ierr);
387552f7358SJed Brown         }
388552f7358SJed Brown         for (d = maxDof; d < enforceDof; ++d) {
389552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr);
390552f7358SJed Brown         }
391552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
392552f7358SJed Brown       }
393552f7358SJed Brown       ierr = PetscFree(remoteValues);CHKERRQ(ierr);
394552f7358SJed Brown     }
395552f7358SJed Brown   } else {
396552f7358SJed Brown     PetscScalar *localValues;
397552f7358SJed Brown     PetscInt    size, k = 0;
398552f7358SJed Brown 
399552f7358SJed Brown     ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
400552f7358SJed Brown     ierr = PetscMalloc(size * sizeof(PetscScalar), &localValues);CHKERRQ(ierr);
401552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
402552f7358SJed Brown       PetscInt dof, off, goff, d;
403552f7358SJed Brown 
404552f7358SJed Brown       /* Reject points not either cells or vertices */
405552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
406552f7358SJed Brown       if (hasLabel) {
407552f7358SJed Brown         PetscInt value;
408552f7358SJed Brown 
409552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
410552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
411552f7358SJed Brown           ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
412552f7358SJed Brown           if (value != 1) continue;
413552f7358SJed Brown         }
414552f7358SJed Brown       }
415552f7358SJed Brown       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
416552f7358SJed Brown       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
417552f7358SJed Brown       ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
418552f7358SJed Brown       if (goff >= 0) {
419552f7358SJed Brown         for (d = 0; d < dof; ++d) {
420552f7358SJed Brown           localValues[k++] = array[off+d];
421552f7358SJed Brown         }
422552f7358SJed Brown       }
423552f7358SJed Brown     }
424552f7358SJed Brown     ierr = MPI_Send(&k, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
425552f7358SJed Brown     ierr = MPI_Send(localValues, k, mpiType, 0, tag, comm);CHKERRQ(ierr);
426552f7358SJed Brown     ierr = PetscFree(localValues);CHKERRQ(ierr);
427552f7358SJed Brown   }
428552f7358SJed Brown   ierr = VecRestoreArray(v, &array);CHKERRQ(ierr);
429552f7358SJed Brown   PetscFunctionReturn(0);
430552f7358SJed Brown }
431552f7358SJed Brown 
432552f7358SJed Brown #undef __FUNCT__
433552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteField_ASCII"
434552f7358SJed Brown PetscErrorCode DMPlexVTKWriteField_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec field, const char name[], FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale)
435552f7358SJed Brown {
43682f516ccSBarry Smith   MPI_Comm       comm;
437552f7358SJed Brown   PetscInt       numDof = 0, maxDof;
438552f7358SJed Brown   PetscInt       pStart, pEnd, p;
439552f7358SJed Brown   PetscErrorCode ierr;
440552f7358SJed Brown 
441552f7358SJed Brown   PetscFunctionBegin;
44282f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
443552f7358SJed Brown   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
444552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
445552f7358SJed Brown     ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr);
4468865f1eaSKarl Rupp     if (numDof) break;
447552f7358SJed Brown   }
448552f7358SJed Brown   numDof = PetscMax(numDof, enforceDof);
44982f516ccSBarry Smith   ierr = MPI_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
450552f7358SJed Brown   if (!name) name = "Unknown";
451552f7358SJed Brown   if (maxDof == 3) {
452552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "VECTORS %s double\n", name);CHKERRQ(ierr);
453552f7358SJed Brown   } else {
454552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "SCALARS %s double %d\n", name, maxDof);CHKERRQ(ierr);
455552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr);
456552f7358SJed Brown   }
457552f7358SJed Brown   ierr = DMPlexVTKWriteSection_ASCII(dm, section, globalSection, field, fp, enforceDof, precision, scale);CHKERRQ(ierr);
458552f7358SJed Brown   PetscFunctionReturn(0);
459552f7358SJed Brown }
460552f7358SJed Brown 
461552f7358SJed Brown #undef __FUNCT__
462552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll_ASCII"
463552f7358SJed Brown static PetscErrorCode DMPlexVTKWriteAll_ASCII(DM dm, PetscViewer viewer)
464552f7358SJed Brown {
46582f516ccSBarry Smith   MPI_Comm                 comm;
466552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*) viewer->data;
467552f7358SJed Brown   FILE                     *fp;
468552f7358SJed Brown   PetscViewerVTKObjectLink link;
469552f7358SJed Brown   PetscSection             coordSection, globalCoordSection;
470552f7358SJed Brown   PetscLayout              vLayout;
471552f7358SJed Brown   Vec                      coordinates;
472552f7358SJed Brown   PetscReal                lengthScale;
473552f7358SJed Brown   PetscInt                 vMax, totVertices, totCells;
474e8f6f0f6SMatthew G. Knepley   PetscBool                hasPoint = PETSC_FALSE, hasCell = PETSC_FALSE, writePartition = PETSC_FALSE;
475552f7358SJed Brown   PetscErrorCode           ierr;
476552f7358SJed Brown 
477552f7358SJed Brown   PetscFunctionBegin;
47882f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
479552f7358SJed Brown   ierr = PetscFOpen(comm, vtk->filename, "wb", &fp);CHKERRQ(ierr);
480552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "# vtk DataFile Version 2.0\n");CHKERRQ(ierr);
481552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "Simplicial Mesh Example\n");CHKERRQ(ierr);
482552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "ASCII\n");CHKERRQ(ierr);
483552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "DATASET UNSTRUCTURED_GRID\n");CHKERRQ(ierr);
484552f7358SJed Brown   /* Vertices */
485552f7358SJed Brown   ierr = DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);CHKERRQ(ierr);
486552f7358SJed Brown   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
487552f7358SJed Brown   ierr = PetscSectionCreateGlobalSection(coordSection, dm->sf, PETSC_FALSE, &globalCoordSection);CHKERRQ(ierr);
488552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
4890298fd71SBarry Smith   ierr = DMPlexGetHybridBounds(dm, NULL, NULL, NULL, &vMax);CHKERRQ(ierr);
490552f7358SJed Brown   if (vMax >= 0) {
491552f7358SJed Brown     PetscInt pStart, pEnd, p, localSize = 0;
492552f7358SJed Brown 
493552f7358SJed Brown     ierr = PetscSectionGetChart(globalCoordSection, &pStart, &pEnd);CHKERRQ(ierr);
494552f7358SJed Brown     pEnd = PetscMin(pEnd, vMax);
495552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
496552f7358SJed Brown       PetscInt dof;
497552f7358SJed Brown 
498552f7358SJed Brown       ierr = PetscSectionGetDof(globalCoordSection, p, &dof);CHKERRQ(ierr);
4998865f1eaSKarl Rupp       if (dof > 0) ++localSize;
500552f7358SJed Brown     }
50182f516ccSBarry Smith     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &vLayout);CHKERRQ(ierr);
502552f7358SJed Brown     ierr = PetscLayoutSetLocalSize(vLayout, localSize);CHKERRQ(ierr);
503552f7358SJed Brown     ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr);
504552f7358SJed Brown     ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr);
505552f7358SJed Brown   } else {
50682f516ccSBarry Smith     ierr = PetscSectionGetPointLayout(PetscObjectComm((PetscObject)dm), globalCoordSection, &vLayout);CHKERRQ(ierr);
507552f7358SJed Brown   }
508552f7358SJed Brown   ierr = PetscLayoutGetSize(vLayout, &totVertices);CHKERRQ(ierr);
509552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "POINTS %d double\n", totVertices);CHKERRQ(ierr);
510552f7358SJed Brown   ierr = DMPlexVTKWriteSection_ASCII(dm, coordSection, globalCoordSection, coordinates, fp, 3, PETSC_DETERMINE, lengthScale);CHKERRQ(ierr);
511552f7358SJed Brown   /* Cells */
512552f7358SJed Brown   ierr = DMPlexVTKWriteCells_ASCII(dm, fp, &totCells);CHKERRQ(ierr);
513552f7358SJed Brown   /* Vertex fields */
514552f7358SJed Brown   for (link = vtk->link; link; link = link->next) {
515552f7358SJed Brown     if ((link->ft == PETSC_VTK_POINT_FIELD) || (link->ft == PETSC_VTK_POINT_VECTOR_FIELD)) hasPoint = PETSC_TRUE;
516552f7358SJed Brown     if ((link->ft == PETSC_VTK_CELL_FIELD)  || (link->ft == PETSC_VTK_CELL_VECTOR_FIELD))  hasCell  = PETSC_TRUE;
517552f7358SJed Brown   }
518552f7358SJed Brown   if (hasPoint) {
51922d28d08SBarry Smith     ierr = PetscFPrintf(comm, fp, "POINT_DATA %d\n", totVertices);CHKERRQ(ierr);
520552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
521552f7358SJed Brown       Vec          X = (Vec) link->vec;
522552f7358SJed Brown       DM           dmX;
5230298fd71SBarry Smith       PetscSection section, globalSection, newSection = NULL;
524552f7358SJed Brown       const char   *name;
525552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
526552f7358SJed Brown 
527552f7358SJed Brown       if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
528552f7358SJed Brown       if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) enforceDof = 3;
529552f7358SJed Brown       ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr);
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 
53522d28d08SBarry Smith         ierr = DMGetDefaultSection(dmX, &section);CHKERRQ(ierr);
536839dd189SMatthew G Knepley         /* Here is where we check whether dmX is a submesh of dm */
537839dd189SMatthew G Knepley         ierr = DMPlexGetDimension(dm,  &dim);CHKERRQ(ierr);
538839dd189SMatthew G Knepley         ierr = DMPlexGetDimension(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) {
544839dd189SMatthew G Knepley           const PetscInt *ind;
545e6ccafaeSMatthew G Knepley           IS              subpointIS;
546839dd189SMatthew G Knepley           PetscInt        n, q;
547839dd189SMatthew G Knepley 
548839dd189SMatthew G Knepley           ierr = PetscPrintf(PETSC_COMM_SELF, "Making translation PetscSection\n");CHKERRQ(ierr);
549839dd189SMatthew G Knepley           ierr = PetscSectionGetChart(section, &qStart, &qEnd);CHKERRQ(ierr);
550e6ccafaeSMatthew G Knepley           ierr = DMPlexCreateSubpointIS(dm, &subpointIS);CHKERRQ(ierr);
551e6ccafaeSMatthew G Knepley           ierr = ISGetLocalSize(subpointIS, &n);CHKERRQ(ierr);
552e6ccafaeSMatthew G Knepley           ierr = ISGetIndices(subpointIS, &ind);CHKERRQ(ierr);
553839dd189SMatthew G Knepley           ierr = PetscSectionCreate(comm, &newSection);CHKERRQ(ierr);
554839dd189SMatthew G Knepley           ierr = PetscSectionSetChart(newSection, pStart, pEnd);CHKERRQ(ierr);
555839dd189SMatthew G Knepley           for (q = qStart; q < qEnd; ++q) {
556839dd189SMatthew G Knepley             PetscInt dof, off, p;
557839dd189SMatthew G Knepley 
558839dd189SMatthew G Knepley             ierr = PetscSectionGetDof(section, q, &dof);CHKERRQ(ierr);
559839dd189SMatthew G Knepley             if (dof) {
560839dd189SMatthew G Knepley               ierr = PetscFindInt(q, n, ind, &p);CHKERRQ(ierr);
561839dd189SMatthew G Knepley               if (p >= pStart) {
562839dd189SMatthew G Knepley                 ierr = PetscSectionSetDof(newSection, p, dof);CHKERRQ(ierr);
563839dd189SMatthew G Knepley                 ierr = PetscSectionGetOffset(section, q, &off);CHKERRQ(ierr);
564839dd189SMatthew G Knepley                 ierr = PetscSectionSetOffset(newSection, p, off);CHKERRQ(ierr);
565839dd189SMatthew G Knepley               }
566839dd189SMatthew G Knepley             }
567839dd189SMatthew G Knepley           }
568e6ccafaeSMatthew G Knepley           ierr = ISRestoreIndices(subpointIS, &ind);CHKERRQ(ierr);
569e6ccafaeSMatthew G Knepley           ierr = ISDestroy(&subpointIS);CHKERRQ(ierr);
570839dd189SMatthew G Knepley           /* No need to setup section */
571839dd189SMatthew G Knepley           section = newSection;
572839dd189SMatthew G Knepley         }
573552f7358SJed Brown       } else {
574426ae2f1SMatthew G Knepley         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
575426ae2f1SMatthew G Knepley         if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
576552f7358SJed Brown       }
57782f516ccSBarry Smith       if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
578552f7358SJed Brown       ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, &globalSection);CHKERRQ(ierr);
579552f7358SJed Brown       ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr);
580552f7358SJed Brown       ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr);
581839dd189SMatthew G Knepley       if (newSection) {ierr = PetscSectionDestroy(&newSection);CHKERRQ(ierr);}
582552f7358SJed Brown     }
583552f7358SJed Brown   }
584552f7358SJed Brown   /* Cell Fields */
585e8f6f0f6SMatthew G. Knepley   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_view_partition", &writePartition, NULL);CHKERRQ(ierr);
586e8f6f0f6SMatthew G. Knepley   if (hasCell || writePartition) {
58722d28d08SBarry Smith     ierr = PetscFPrintf(comm, fp, "CELL_DATA %d\n", totCells);CHKERRQ(ierr);
588552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
589552f7358SJed Brown       Vec          X = (Vec) link->vec;
590552f7358SJed Brown       DM           dmX;
591552f7358SJed Brown       PetscSection section, globalSection;
592552f7358SJed Brown       const char   *name;
593552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
594552f7358SJed Brown 
595552f7358SJed Brown       if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
596552f7358SJed Brown       if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) enforceDof = 3;
597552f7358SJed Brown       ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr);
598552f7358SJed Brown       ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
599552f7358SJed Brown       if (dmX) {
60022d28d08SBarry Smith         ierr = DMGetDefaultSection(dmX, &section);CHKERRQ(ierr);
601552f7358SJed Brown       } else {
602552f7358SJed Brown         PetscContainer c;
603552f7358SJed Brown 
604552f7358SJed Brown         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &c);CHKERRQ(ierr);
60582f516ccSBarry Smith         if (!c) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
606552f7358SJed Brown         ierr = PetscContainerGetPointer(c, (void**) &section);CHKERRQ(ierr);
607552f7358SJed Brown       }
60882f516ccSBarry Smith       if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
609552f7358SJed Brown       ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, &globalSection);CHKERRQ(ierr);
610552f7358SJed Brown       ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr);
611552f7358SJed Brown       ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr);
612552f7358SJed Brown     }
613e8f6f0f6SMatthew G. Knepley     if (writePartition) {
614e8f6f0f6SMatthew G. Knepley       ierr = PetscFPrintf(comm, fp, "SCALARS partition int 1\n");CHKERRQ(ierr);
615e8f6f0f6SMatthew G. Knepley       ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr);
616e8f6f0f6SMatthew G. Knepley       ierr = DMPlexVTKWritePartition_ASCII(dm, fp);CHKERRQ(ierr);
617e8f6f0f6SMatthew G. Knepley     }
618552f7358SJed Brown   }
619552f7358SJed Brown   /* Cleanup */
620552f7358SJed Brown   ierr = PetscSectionDestroy(&globalCoordSection);CHKERRQ(ierr);
621552f7358SJed Brown   ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr);
622552f7358SJed Brown   ierr = PetscFClose(comm, fp);CHKERRQ(ierr);
623552f7358SJed Brown   PetscFunctionReturn(0);
624552f7358SJed Brown }
625552f7358SJed Brown 
626552f7358SJed Brown #undef __FUNCT__
627552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll"
628552f7358SJed Brown /*@C
629552f7358SJed Brown   DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
630552f7358SJed Brown 
631552f7358SJed Brown   Collective
632552f7358SJed Brown 
633552f7358SJed Brown   Input Arguments:
634552f7358SJed Brown + odm - The DMPlex specifying the mesh, passed as a PetscObject
635552f7358SJed Brown - viewer - viewer of type VTK
636552f7358SJed Brown 
637552f7358SJed Brown   Level: developer
638552f7358SJed Brown 
639552f7358SJed Brown   Note:
640552f7358SJed Brown   This function is a callback used by the VTK viewer to actually write the file.
641552f7358SJed 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.
642552f7358SJed Brown   Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
643552f7358SJed Brown 
644552f7358SJed Brown .seealso: PETSCVIEWERVTK
645552f7358SJed Brown @*/
646552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer)
647552f7358SJed Brown {
648552f7358SJed Brown   DM             dm = (DM) odm;
649552f7358SJed Brown   PetscBool      isvtk;
650552f7358SJed Brown   PetscErrorCode ierr;
651552f7358SJed Brown 
652552f7358SJed Brown   PetscFunctionBegin;
653552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
654552f7358SJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
655552f7358SJed Brown   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr);
65682f516ccSBarry Smith   if (!isvtk) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name);
657552f7358SJed Brown   switch (viewer->format) {
658552f7358SJed Brown   case PETSC_VIEWER_ASCII_VTK:
659552f7358SJed Brown     ierr = DMPlexVTKWriteAll_ASCII(dm, viewer);CHKERRQ(ierr);
660552f7358SJed Brown     break;
661552f7358SJed Brown   case PETSC_VIEWER_VTK_VTU:
662552f7358SJed Brown     ierr = DMPlexVTKWriteAll_VTU(dm, viewer);CHKERRQ(ierr);
663552f7358SJed Brown     break;
66482f516ccSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]);
665552f7358SJed Brown   }
666552f7358SJed Brown   PetscFunctionReturn(0);
667552f7358SJed Brown }
668