xref: /petsc/src/dm/impls/plex/plexvtk.c (revision 2fb742c9ec26461b89af95227d23a4873a850bdd)
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 
27476b700caSToby Isaac #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
27576b700caSToby Isaac typedef double PetscVTKReal;
27676b700caSToby Isaac #elif defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL___FP16)
27776b700caSToby Isaac typedef float PetscVTKReal;
27876b700caSToby Isaac #else
27976b700caSToby Isaac typedef PetscReal PetscVTKReal;
28076b700caSToby Isaac #endif
28176b700caSToby Isaac 
28276b700caSToby Isaac static PetscErrorCode DMPlexVTKWriteSection_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec v, FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale, PetscInt imag)
2830adebc6cSBarry Smith {
28482f516ccSBarry Smith   MPI_Comm           comm;
285552f7358SJed Brown   const MPI_Datatype mpiType = MPIU_SCALAR;
286552f7358SJed Brown   PetscScalar        *array;
287552f7358SJed Brown   PetscInt           numDof = 0, maxDof;
2882f029394SStefano Zampini   PetscInt           numLabelCells, cellHeight, cStart, cEnd, numLabelVertices, vMax, vStart, vEnd, pStart, pEnd, p;
2899852e123SBarry Smith   PetscMPIInt        size, rank, proc, tag;
290552f7358SJed Brown   PetscBool          hasLabel;
291552f7358SJed Brown   PetscErrorCode     ierr;
292552f7358SJed Brown 
293552f7358SJed Brown   PetscFunctionBegin;
29482f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
295552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
296552f7358SJed Brown   PetscValidHeaderSpecific(v,VEC_CLASSID,4);
297552f7358SJed Brown   if (precision < 0) precision = 6;
298552f7358SJed Brown   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
2999852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
300552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
301552f7358SJed Brown   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
302552f7358SJed Brown   /* VTK only wants the values at cells or vertices */
303552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
304552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
305552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
3062f029394SStefano Zampini   ierr = DMPlexGetHybridBounds(dm, NULL, NULL, NULL, &vMax);CHKERRQ(ierr);
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);
310c58f1c22SToby Isaac   ierr     = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
311c58f1c22SToby Isaac   ierr     = DMGetStratumSize(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)) {
321c58f1c22SToby Isaac         ierr = DMGetLabelValue(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   }
328b2566f29SBarry Smith   ierr = MPIU_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) {
33276b700caSToby Isaac     PetscVTKReal dval;
33376b700caSToby Isaac     PetscScalar  val;
334552f7358SJed Brown     char formatString[8];
335552f7358SJed Brown 
336552f7358SJed Brown     ierr = PetscSNPrintf(formatString, 8, "%%.%de", precision);CHKERRQ(ierr);
337552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
338552f7358SJed Brown       /* Here we lose a way to filter points by keeping them out of the Numbering */
339552f7358SJed Brown       PetscInt dof, off, goff, d;
340552f7358SJed Brown 
341552f7358SJed Brown       /* Reject points not either cells or vertices */
342552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
343552f7358SJed Brown       if (hasLabel) {
344552f7358SJed Brown         PetscInt value;
345552f7358SJed Brown 
346552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
347552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
348c58f1c22SToby Isaac           ierr = DMGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
349552f7358SJed Brown           if (value != 1) continue;
350552f7358SJed Brown         }
351552f7358SJed Brown       }
352552f7358SJed Brown       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
353552f7358SJed Brown       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
354552f7358SJed Brown       ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
355552f7358SJed Brown       if (dof && goff >= 0) {
356552f7358SJed Brown         for (d = 0; d < dof; d++) {
357552f7358SJed Brown           if (d > 0) {
358552f7358SJed Brown             ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr);
359552f7358SJed Brown           }
36076b700caSToby Isaac           val = array[off+d];
36176b700caSToby Isaac           dval = (PetscVTKReal) ((imag ? PetscImaginaryPart(val) : PetscRealPart(val)) * scale);
3622dea1156SStefano Zampini           ierr = PetscFPrintf(comm, fp, formatString, dval);CHKERRQ(ierr);
363552f7358SJed Brown         }
364552f7358SJed Brown         for (d = dof; d < enforceDof; d++) {
365552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr);
366552f7358SJed Brown         }
367552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
368552f7358SJed Brown       }
369552f7358SJed Brown     }
3709852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
371552f7358SJed Brown       PetscScalar *remoteValues;
372d892089bSMatthew G. Knepley       PetscInt    size = 0, d;
373552f7358SJed Brown       MPI_Status  status;
374552f7358SJed Brown 
375552f7358SJed Brown       ierr = MPI_Recv(&size, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
376785e854fSJed Brown       ierr = PetscMalloc1(size, &remoteValues);CHKERRQ(ierr);
377552f7358SJed Brown       ierr = MPI_Recv(remoteValues, size, mpiType, proc, tag, comm, &status);CHKERRQ(ierr);
378552f7358SJed Brown       for (p = 0; p < size/maxDof; ++p) {
379552f7358SJed Brown         for (d = 0; d < maxDof; ++d) {
380552f7358SJed Brown           if (d > 0) {
381552f7358SJed Brown             ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr);
382552f7358SJed Brown           }
38376b700caSToby Isaac           val = remoteValues[p*maxDof+d];
38476b700caSToby Isaac           dval = (PetscVTKReal) ((imag ? PetscImaginaryPart(val) : PetscRealPart(val)) * scale);
3852dea1156SStefano Zampini           ierr = PetscFPrintf(comm, fp, formatString, dval);CHKERRQ(ierr);
386552f7358SJed Brown         }
387552f7358SJed Brown         for (d = maxDof; d < enforceDof; ++d) {
388552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr);
389552f7358SJed Brown         }
390552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
391552f7358SJed Brown       }
392552f7358SJed Brown       ierr = PetscFree(remoteValues);CHKERRQ(ierr);
393552f7358SJed Brown     }
394552f7358SJed Brown   } else {
395552f7358SJed Brown     PetscScalar *localValues;
396552f7358SJed Brown     PetscInt    size, k = 0;
397552f7358SJed Brown 
398552f7358SJed Brown     ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
399785e854fSJed Brown     ierr = PetscMalloc1(size, &localValues);CHKERRQ(ierr);
400552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
401552f7358SJed Brown       PetscInt dof, off, goff, d;
402552f7358SJed Brown 
403552f7358SJed Brown       /* Reject points not either cells or vertices */
404552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
405552f7358SJed Brown       if (hasLabel) {
406552f7358SJed Brown         PetscInt value;
407552f7358SJed Brown 
408552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
409552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
410c58f1c22SToby Isaac           ierr = DMGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
411552f7358SJed Brown           if (value != 1) continue;
412552f7358SJed Brown         }
413552f7358SJed Brown       }
414552f7358SJed Brown       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
415552f7358SJed Brown       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
416552f7358SJed Brown       ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
417552f7358SJed Brown       if (goff >= 0) {
418552f7358SJed Brown         for (d = 0; d < dof; ++d) {
419552f7358SJed Brown           localValues[k++] = array[off+d];
420552f7358SJed Brown         }
421552f7358SJed Brown       }
422552f7358SJed Brown     }
423552f7358SJed Brown     ierr = MPI_Send(&k, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
424552f7358SJed Brown     ierr = MPI_Send(localValues, k, mpiType, 0, tag, comm);CHKERRQ(ierr);
425552f7358SJed Brown     ierr = PetscFree(localValues);CHKERRQ(ierr);
426552f7358SJed Brown   }
427552f7358SJed Brown   ierr = VecRestoreArray(v, &array);CHKERRQ(ierr);
428552f7358SJed Brown   PetscFunctionReturn(0);
429552f7358SJed Brown }
430552f7358SJed Brown 
43176b700caSToby Isaac static PetscErrorCode DMPlexVTKWriteField_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec field, const char name[], FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale, PetscBool nameComplex, PetscInt imag)
432552f7358SJed Brown {
43382f516ccSBarry Smith   MPI_Comm       comm;
434552f7358SJed Brown   PetscInt       numDof = 0, maxDof;
435552f7358SJed Brown   PetscInt       pStart, pEnd, p;
436552f7358SJed Brown   PetscErrorCode ierr;
437552f7358SJed Brown 
438552f7358SJed Brown   PetscFunctionBegin;
43982f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
440552f7358SJed Brown   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
441552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
442552f7358SJed Brown     ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr);
4438865f1eaSKarl Rupp     if (numDof) break;
444552f7358SJed Brown   }
445552f7358SJed Brown   numDof = PetscMax(numDof, enforceDof);
446b2566f29SBarry Smith   ierr = MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
447552f7358SJed Brown   if (!name) name = "Unknown";
448552f7358SJed Brown   if (maxDof == 3) {
44976b700caSToby Isaac     if (nameComplex) {
45076b700caSToby Isaac       ierr = PetscFPrintf(comm, fp, "VECTORS %s.%s double\n", name, imag ? "Im" : "Re");CHKERRQ(ierr);
45176b700caSToby Isaac     } else {
452552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "VECTORS %s double\n", name);CHKERRQ(ierr);
45376b700caSToby Isaac     }
45476b700caSToby Isaac   } else {
45576b700caSToby Isaac     if (nameComplex) {
45676b700caSToby Isaac       ierr = PetscFPrintf(comm, fp, "SCALARS %s.%s double %D\n", name, imag ? "Im" : "Re", maxDof);CHKERRQ(ierr);
457552f7358SJed Brown     } else {
4586b8858a1SBarry Smith       ierr = PetscFPrintf(comm, fp, "SCALARS %s double %D\n", name, maxDof);CHKERRQ(ierr);
45976b700caSToby Isaac     }
460552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr);
461552f7358SJed Brown   }
46276b700caSToby Isaac   ierr = DMPlexVTKWriteSection_ASCII(dm, section, globalSection, field, fp, enforceDof, precision, scale, imag);CHKERRQ(ierr);
463552f7358SJed Brown   PetscFunctionReturn(0);
464552f7358SJed Brown }
465552f7358SJed Brown 
466552f7358SJed Brown static PetscErrorCode DMPlexVTKWriteAll_ASCII(DM dm, PetscViewer viewer)
467552f7358SJed Brown {
46882f516ccSBarry Smith   MPI_Comm                 comm;
469552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*) viewer->data;
470552f7358SJed Brown   FILE                     *fp;
471552f7358SJed Brown   PetscViewerVTKObjectLink link;
472552f7358SJed Brown   PetscSection             coordSection, globalCoordSection;
473552f7358SJed Brown   PetscLayout              vLayout;
474552f7358SJed Brown   Vec                      coordinates;
475552f7358SJed Brown   PetscReal                lengthScale;
47676b700caSToby Isaac   PetscInt                 vMax, totVertices, totCells = 0, loops_per_scalar, l;
47776b700caSToby Isaac   PetscBool                hasPoint = PETSC_FALSE, hasCell = PETSC_FALSE, writePartition = PETSC_FALSE, localized, writeComplex;
478552f7358SJed Brown   PetscErrorCode           ierr;
479552f7358SJed Brown 
480552f7358SJed Brown   PetscFunctionBegin;
48176b700caSToby Isaac #if defined(PETSC_USE_COMPLEX)
48276b700caSToby Isaac   loops_per_scalar = 2;
48376b700caSToby Isaac   writeComplex = PETSC_TRUE;
48476b700caSToby Isaac #else
48576b700caSToby Isaac   loops_per_scalar = 1;
48676b700caSToby Isaac   writeComplex = PETSC_FALSE;
48776b700caSToby Isaac #endif
488d72c04d6SStefano Zampini   ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr);
489d72c04d6SStefano Zampini   if (localized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"VTK output with localized coordinates not yet supported");
49082f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
491552f7358SJed Brown   ierr = PetscFOpen(comm, vtk->filename, "wb", &fp);CHKERRQ(ierr);
492552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "# vtk DataFile Version 2.0\n");CHKERRQ(ierr);
493552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "Simplicial Mesh Example\n");CHKERRQ(ierr);
494552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "ASCII\n");CHKERRQ(ierr);
495552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "DATASET UNSTRUCTURED_GRID\n");CHKERRQ(ierr);
496552f7358SJed Brown   /* Vertices */
497552f7358SJed Brown   ierr = DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);CHKERRQ(ierr);
49869d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
49915b58121SMatthew G. Knepley   ierr = PetscSectionCreateGlobalSection(coordSection, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalCoordSection);CHKERRQ(ierr);
500552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
5010298fd71SBarry Smith   ierr = DMPlexGetHybridBounds(dm, NULL, NULL, NULL, &vMax);CHKERRQ(ierr);
502552f7358SJed Brown   if (vMax >= 0) {
503552f7358SJed Brown     PetscInt pStart, pEnd, p, localSize = 0;
504552f7358SJed Brown 
505552f7358SJed Brown     ierr = PetscSectionGetChart(globalCoordSection, &pStart, &pEnd);CHKERRQ(ierr);
506552f7358SJed Brown     pEnd = PetscMin(pEnd, vMax);
507552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
508552f7358SJed Brown       PetscInt dof;
509552f7358SJed Brown 
510552f7358SJed Brown       ierr = PetscSectionGetDof(globalCoordSection, p, &dof);CHKERRQ(ierr);
5118865f1eaSKarl Rupp       if (dof > 0) ++localSize;
512552f7358SJed Brown     }
51382f516ccSBarry Smith     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &vLayout);CHKERRQ(ierr);
514552f7358SJed Brown     ierr = PetscLayoutSetLocalSize(vLayout, localSize);CHKERRQ(ierr);
515552f7358SJed Brown     ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr);
516552f7358SJed Brown     ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr);
517552f7358SJed Brown   } else {
51882f516ccSBarry Smith     ierr = PetscSectionGetPointLayout(PetscObjectComm((PetscObject)dm), globalCoordSection, &vLayout);CHKERRQ(ierr);
519552f7358SJed Brown   }
520552f7358SJed Brown   ierr = PetscLayoutGetSize(vLayout, &totVertices);CHKERRQ(ierr);
5216b8858a1SBarry Smith   ierr = PetscFPrintf(comm, fp, "POINTS %D double\n", totVertices);CHKERRQ(ierr);
52276b700caSToby Isaac   ierr = DMPlexVTKWriteSection_ASCII(dm, coordSection, globalCoordSection, coordinates, fp, 3, PETSC_DETERMINE, lengthScale, 0);CHKERRQ(ierr);
523552f7358SJed Brown   /* Cells */
524552f7358SJed Brown   ierr = DMPlexVTKWriteCells_ASCII(dm, fp, &totCells);CHKERRQ(ierr);
525552f7358SJed Brown   /* Vertex fields */
526552f7358SJed Brown   for (link = vtk->link; link; link = link->next) {
527552f7358SJed Brown     if ((link->ft == PETSC_VTK_POINT_FIELD) || (link->ft == PETSC_VTK_POINT_VECTOR_FIELD)) hasPoint = PETSC_TRUE;
528552f7358SJed Brown     if ((link->ft == PETSC_VTK_CELL_FIELD)  || (link->ft == PETSC_VTK_CELL_VECTOR_FIELD))  hasCell  = PETSC_TRUE;
529552f7358SJed Brown   }
530552f7358SJed Brown   if (hasPoint) {
5316b8858a1SBarry Smith     ierr = PetscFPrintf(comm, fp, "POINT_DATA %D\n", totVertices);CHKERRQ(ierr);
532552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
533552f7358SJed Brown       Vec          X = (Vec) link->vec;
534e630c359SToby Isaac       PetscSection section = NULL, globalSection, newSection = NULL;
535e630c359SToby Isaac       char         namebuf[256];
536552f7358SJed Brown       const char   *name;
537552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
538552f7358SJed Brown 
539552f7358SJed Brown       if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
540552f7358SJed Brown       if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) enforceDof = 3;
541552f7358SJed Brown       ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr);
542e630c359SToby Isaac       ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
543e630c359SToby Isaac       if (!section) {
544e630c359SToby Isaac         DM           dmX;
545e630c359SToby Isaac 
546552f7358SJed Brown         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
547552f7358SJed Brown         if (dmX) {
548088580cfSMatthew G Knepley           DMLabel  subpointMap, subpointMapX;
549839dd189SMatthew G Knepley           PetscInt dim, dimX, pStart, pEnd, qStart, qEnd;
550839dd189SMatthew G Knepley 
55192fd8e1eSJed Brown           ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);
552839dd189SMatthew G Knepley           /* Here is where we check whether dmX is a submesh of dm */
553c73cfb54SMatthew G. Knepley           ierr = DMGetDimension(dm,  &dim);CHKERRQ(ierr);
554c73cfb54SMatthew G. Knepley           ierr = DMGetDimension(dmX, &dimX);CHKERRQ(ierr);
555839dd189SMatthew G Knepley           ierr = DMPlexGetChart(dm,  &pStart, &pEnd);CHKERRQ(ierr);
556839dd189SMatthew G Knepley           ierr = DMPlexGetChart(dmX, &qStart, &qEnd);CHKERRQ(ierr);
557839dd189SMatthew G Knepley           ierr = DMPlexGetSubpointMap(dm,  &subpointMap);CHKERRQ(ierr);
558839dd189SMatthew G Knepley           ierr = DMPlexGetSubpointMap(dmX, &subpointMapX);CHKERRQ(ierr);
559839dd189SMatthew G Knepley           if (((dim != dimX) || ((pEnd-pStart) < (qEnd-qStart))) && subpointMap && !subpointMapX) {
560434e42b5SMatthew G. Knepley             const PetscInt *ind = NULL;
561e6ccafaeSMatthew G Knepley             IS              subpointIS;
562434e42b5SMatthew G. Knepley             PetscInt        n = 0, q;
563839dd189SMatthew G Knepley 
564839dd189SMatthew G Knepley             ierr = PetscSectionGetChart(section, &qStart, &qEnd);CHKERRQ(ierr);
565e6ccafaeSMatthew G Knepley             ierr = DMPlexCreateSubpointIS(dm, &subpointIS);CHKERRQ(ierr);
566434e42b5SMatthew G. Knepley             if (subpointIS) {
567e6ccafaeSMatthew G Knepley               ierr = ISGetLocalSize(subpointIS, &n);CHKERRQ(ierr);
568e6ccafaeSMatthew G Knepley               ierr = ISGetIndices(subpointIS, &ind);CHKERRQ(ierr);
569434e42b5SMatthew G. Knepley             }
570839dd189SMatthew G Knepley             ierr = PetscSectionCreate(comm, &newSection);CHKERRQ(ierr);
571839dd189SMatthew G Knepley             ierr = PetscSectionSetChart(newSection, pStart, pEnd);CHKERRQ(ierr);
572839dd189SMatthew G Knepley             for (q = qStart; q < qEnd; ++q) {
573839dd189SMatthew G Knepley               PetscInt dof, off, p;
574839dd189SMatthew G Knepley 
575839dd189SMatthew G Knepley               ierr = PetscSectionGetDof(section, q, &dof);CHKERRQ(ierr);
576839dd189SMatthew G Knepley               if (dof) {
577839dd189SMatthew G Knepley                 ierr = PetscFindInt(q, n, ind, &p);CHKERRQ(ierr);
578839dd189SMatthew G Knepley                 if (p >= pStart) {
579839dd189SMatthew G Knepley                   ierr = PetscSectionSetDof(newSection, p, dof);CHKERRQ(ierr);
580839dd189SMatthew G Knepley                   ierr = PetscSectionGetOffset(section, q, &off);CHKERRQ(ierr);
581839dd189SMatthew G Knepley                   ierr = PetscSectionSetOffset(newSection, p, off);CHKERRQ(ierr);
582839dd189SMatthew G Knepley                 }
583839dd189SMatthew G Knepley               }
584839dd189SMatthew G Knepley             }
585434e42b5SMatthew G. Knepley             if (subpointIS) {
586e6ccafaeSMatthew G Knepley               ierr = ISRestoreIndices(subpointIS, &ind);CHKERRQ(ierr);
587e6ccafaeSMatthew G Knepley               ierr = ISDestroy(&subpointIS);CHKERRQ(ierr);
588434e42b5SMatthew G. Knepley             }
589839dd189SMatthew G Knepley             /* No need to setup section */
590839dd189SMatthew G Knepley             section = newSection;
591839dd189SMatthew G Knepley           }
592552f7358SJed Brown         }
593e630c359SToby Isaac       }
594e630c359SToby 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);
595e630c359SToby Isaac       if (link->field >= 0) {
596e630c359SToby Isaac         const char *fieldname;
597e630c359SToby Isaac 
598e630c359SToby Isaac         ierr = PetscSectionGetFieldName(section, link->field, &fieldname);CHKERRQ(ierr);
599e630c359SToby Isaac         ierr = PetscSectionGetField(section, link->field, &section);CHKERRQ(ierr);
600e630c359SToby Isaac         if (fieldname) {
601*2fb742c9SToby Isaac           ierr = PetscSNPrintf(namebuf, sizeof(namebuf), "%s%s", name, fieldname);CHKERRQ(ierr);
602e630c359SToby Isaac         } else {
603*2fb742c9SToby Isaac           ierr = PetscSNPrintf(namebuf, sizeof(namebuf), "%s%D", name, link->field);CHKERRQ(ierr);
604e630c359SToby Isaac         }
605*2fb742c9SToby Isaac       } else {
606*2fb742c9SToby Isaac         ierr = PetscSNPrintf(namebuf, sizeof(namebuf), "%s", name);CHKERRQ(ierr);
607e630c359SToby Isaac       }
608*2fb742c9SToby Isaac       ierr = PetscViewerVTKSanitizeName_Internal(namebuf, sizeof(namebuf));CHKERRQ(ierr);
60915b58121SMatthew G. Knepley       ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection);CHKERRQ(ierr);
61076b700caSToby Isaac       for (l = 0; l < loops_per_scalar; l++) {
611*2fb742c9SToby Isaac         ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, namebuf, fp, enforceDof, PETSC_DETERMINE, 1.0, writeComplex, l);CHKERRQ(ierr);
61276b700caSToby Isaac       }
613552f7358SJed Brown       ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr);
614839dd189SMatthew G Knepley       if (newSection) {ierr = PetscSectionDestroy(&newSection);CHKERRQ(ierr);}
615552f7358SJed Brown     }
616552f7358SJed Brown   }
617552f7358SJed Brown   /* Cell Fields */
618c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_view_partition", &writePartition, NULL);CHKERRQ(ierr);
619e8f6f0f6SMatthew G. Knepley   if (hasCell || writePartition) {
6206b8858a1SBarry Smith     ierr = PetscFPrintf(comm, fp, "CELL_DATA %D\n", totCells);CHKERRQ(ierr);
621552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
622552f7358SJed Brown       Vec          X = (Vec) link->vec;
623e630c359SToby Isaac       PetscSection section = NULL, globalSection;
624*2fb742c9SToby Isaac       const char   *name = "";
625e630c359SToby Isaac       char         namebuf[256];
626552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
627552f7358SJed Brown 
628552f7358SJed Brown       if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
629552f7358SJed Brown       if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) enforceDof = 3;
630552f7358SJed Brown       ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr);
631e630c359SToby Isaac       ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
632e630c359SToby Isaac       if (!section) {
633e630c359SToby Isaac         DM           dmX;
634e630c359SToby Isaac 
635552f7358SJed Brown         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
636552f7358SJed Brown         if (dmX) {
63792fd8e1eSJed Brown           ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);
638552f7358SJed Brown         }
639e630c359SToby Isaac       }
640e630c359SToby 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);
641e630c359SToby Isaac       if (link->field >= 0) {
642e630c359SToby Isaac         const char *fieldname;
643e630c359SToby Isaac 
644e630c359SToby Isaac         ierr = PetscSectionGetFieldName(section, link->field, &fieldname);CHKERRQ(ierr);
645e630c359SToby Isaac         ierr = PetscSectionGetField(section, link->field, &section);CHKERRQ(ierr);
646e630c359SToby Isaac         if (fieldname) {
647*2fb742c9SToby Isaac           ierr = PetscSNPrintf(namebuf, sizeof(namebuf), "%s%s", name, fieldname);CHKERRQ(ierr);
648e630c359SToby Isaac         } else {
649*2fb742c9SToby Isaac           ierr = PetscSNPrintf(namebuf, sizeof(namebuf), "%s%D", name, link->field);CHKERRQ(ierr);
650e630c359SToby Isaac         }
651*2fb742c9SToby Isaac       } else {
652*2fb742c9SToby Isaac         ierr = PetscSNPrintf(namebuf, sizeof(namebuf), "%s", name);CHKERRQ(ierr);
653e630c359SToby Isaac       }
654*2fb742c9SToby Isaac       ierr = PetscViewerVTKSanitizeName_Internal(namebuf, sizeof(namebuf));CHKERRQ(ierr);
65515b58121SMatthew G. Knepley       ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection);CHKERRQ(ierr);
65676b700caSToby Isaac       for (l = 0; l < loops_per_scalar; l++) {
657*2fb742c9SToby Isaac         ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, namebuf, fp, enforceDof, PETSC_DETERMINE, 1.0, writeComplex, l);CHKERRQ(ierr);
65876b700caSToby Isaac       }
659552f7358SJed Brown       ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr);
660552f7358SJed Brown     }
661e8f6f0f6SMatthew G. Knepley     if (writePartition) {
662e8f6f0f6SMatthew G. Knepley       ierr = PetscFPrintf(comm, fp, "SCALARS partition int 1\n");CHKERRQ(ierr);
663e8f6f0f6SMatthew G. Knepley       ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr);
664e8f6f0f6SMatthew G. Knepley       ierr = DMPlexVTKWritePartition_ASCII(dm, fp);CHKERRQ(ierr);
665e8f6f0f6SMatthew G. Knepley     }
666552f7358SJed Brown   }
667552f7358SJed Brown   /* Cleanup */
668552f7358SJed Brown   ierr = PetscSectionDestroy(&globalCoordSection);CHKERRQ(ierr);
669552f7358SJed Brown   ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr);
670552f7358SJed Brown   ierr = PetscFClose(comm, fp);CHKERRQ(ierr);
671552f7358SJed Brown   PetscFunctionReturn(0);
672552f7358SJed Brown }
673552f7358SJed Brown 
674552f7358SJed Brown /*@C
675552f7358SJed Brown   DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
676552f7358SJed Brown 
677552f7358SJed Brown   Collective
678552f7358SJed Brown 
679552f7358SJed Brown   Input Arguments:
680552f7358SJed Brown + odm - The DMPlex specifying the mesh, passed as a PetscObject
681552f7358SJed Brown - viewer - viewer of type VTK
682552f7358SJed Brown 
683552f7358SJed Brown   Level: developer
684552f7358SJed Brown 
685552f7358SJed Brown   Note:
686552f7358SJed Brown   This function is a callback used by the VTK viewer to actually write the file.
687552f7358SJed 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.
688552f7358SJed Brown   Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
689552f7358SJed Brown 
690552f7358SJed Brown .seealso: PETSCVIEWERVTK
691552f7358SJed Brown @*/
692552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer)
693552f7358SJed Brown {
694552f7358SJed Brown   DM             dm = (DM) odm;
695552f7358SJed Brown   PetscBool      isvtk;
696552f7358SJed Brown   PetscErrorCode ierr;
697552f7358SJed Brown 
698552f7358SJed Brown   PetscFunctionBegin;
699552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
700552f7358SJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
701552f7358SJed Brown   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr);
70282f516ccSBarry Smith   if (!isvtk) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name);
703552f7358SJed Brown   switch (viewer->format) {
704552f7358SJed Brown   case PETSC_VIEWER_ASCII_VTK:
705552f7358SJed Brown     ierr = DMPlexVTKWriteAll_ASCII(dm, viewer);CHKERRQ(ierr);
706552f7358SJed Brown     break;
707552f7358SJed Brown   case PETSC_VIEWER_VTK_VTU:
708552f7358SJed Brown     ierr = DMPlexVTKWriteAll_VTU(dm, viewer);CHKERRQ(ierr);
709552f7358SJed Brown     break;
71082f516ccSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]);
711552f7358SJed Brown   }
712552f7358SJed Brown   PetscFunctionReturn(0);
713552f7358SJed Brown }
714