xref: /petsc/src/dm/impls/plex/plexvtk.c (revision 1dca8a0504492127e77eac64bc165d7372dd6d63)
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 
85552f7358SJed Brown   PetscFunctionBegin;
869566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
879566063dSJacob Faibussowitsch   PetscCall(PetscCommGetNewTag(comm, &tag));
889566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
899566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
909566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
919566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
929566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
939566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
949566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "vtk", &label));
959566063dSJacob Faibussowitsch   PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells));
961c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&numLabelCells, &maxLabelCells, 1, MPIU_INT, MPI_MAX, comm));
97452b7979SMatthew G. Knepley   if (!maxLabelCells) label = NULL;
98552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
990298fd71SBarry Smith     PetscInt *closure = NULL;
10002281ff3SMatthew G. Knepley     PetscInt closureSize, value;
101552f7358SJed Brown 
10202281ff3SMatthew G. Knepley     if (label) {
1039566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(label, c, &value));
104552f7358SJed Brown       if (value != 1) continue;
105552f7358SJed Brown     }
1069566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
107552f7358SJed Brown     for (v = 0; v < closureSize*2; v += 2) {
108552f7358SJed Brown       if ((closure[v] >= vStart) && (closure[v] < vEnd)) ++numCorners;
109552f7358SJed Brown     }
1109566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
111552f7358SJed Brown     ++numCells;
112552f7358SJed Brown   }
113552f7358SJed Brown   maxCells = numCells;
1149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Reduce(&numCells, &totCells, 1, MPIU_INT, MPI_SUM, 0, comm));
1159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Reduce(&numCells, &maxCells, 1, MPIU_INT, MPI_MAX, 0, comm));
1169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Reduce(&numCorners, &totCorners, 1, MPIU_INT, MPI_SUM, 0, comm));
1179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Reduce(&numCorners, &maxCorners, 1, MPIU_INT, MPI_MAX, 0, comm));
1189566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVertexNumbering(dm, &globalVertexNumbers));
1199566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(globalVertexNumbers, &gvertex));
1209566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxCells, &corners));
12163a3b9bcSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "CELLS %" PetscInt_FMT " %" PetscInt_FMT "\n", totCells, totCorners+totCells));
122dd400576SPatrick Sanan   if (rank == 0) {
123a4a685f2SJacob Faibussowitsch     PetscInt *remoteVertices, *vertices;
124552f7358SJed Brown 
1259566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(maxCorners, &vertices));
126552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
1270298fd71SBarry Smith       PetscInt *closure = NULL;
12802281ff3SMatthew G. Knepley       PetscInt closureSize, value, nC = 0;
129552f7358SJed Brown 
13002281ff3SMatthew G. Knepley       if (label) {
1319566063dSJacob Faibussowitsch         PetscCall(DMLabelGetValue(label, c, &value));
132552f7358SJed Brown         if (value != 1) continue;
133552f7358SJed Brown       }
1349566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
135552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
136552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
137552f7358SJed Brown           const PetscInt gv = gvertex[closure[v] - vStart];
1380dcf9c70SMatthew G. Knepley           vertices[nC++] = gv < 0 ? -(gv+1) : gv;
139552f7358SJed Brown         }
140552f7358SJed Brown       }
1419566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
1429566063dSJacob Faibussowitsch       PetscCall(DMPlexReorderCell(dm, c, vertices));
143552f7358SJed Brown       corners[numCells++] = nC;
14463a3b9bcSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "%" PetscInt_FMT " ", nC));
145552f7358SJed Brown       for (v = 0; v < nC; ++v) {
14663a3b9bcSJacob Faibussowitsch         PetscCall(PetscFPrintf(comm, fp, " %" PetscInt_FMT, vertices[v]));
147552f7358SJed Brown       }
1489566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "\n"));
149552f7358SJed Brown     }
1509566063dSJacob Faibussowitsch     if (size > 1) PetscCall(PetscMalloc1(maxCorners+maxCells, &remoteVertices));
1519852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
152552f7358SJed Brown       MPI_Status status;
153552f7358SJed Brown 
1549566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(&numCorners, 1, MPIU_INT, proc, tag, comm, &status));
1559566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(remoteVertices, numCorners, MPIU_INT, proc, tag, comm, &status));
156552f7358SJed Brown       for (c = 0; c < numCorners;) {
157552f7358SJed Brown         PetscInt nC = remoteVertices[c++];
158552f7358SJed Brown 
159552f7358SJed Brown         for (v = 0; v < nC; ++v, ++c) {
1600dcf9c70SMatthew G. Knepley           vertices[v] = remoteVertices[c];
1610dcf9c70SMatthew G. Knepley         }
16263a3b9bcSJacob Faibussowitsch         PetscCall(PetscFPrintf(comm, fp, "%" PetscInt_FMT " ", nC));
1630dcf9c70SMatthew G. Knepley         for (v = 0; v < nC; ++v) {
16463a3b9bcSJacob Faibussowitsch           PetscCall(PetscFPrintf(comm, fp, " %" PetscInt_FMT, vertices[v]));
165552f7358SJed Brown         }
1669566063dSJacob Faibussowitsch         PetscCall(PetscFPrintf(comm, fp, "\n"));
167552f7358SJed Brown       }
168552f7358SJed Brown     }
1699566063dSJacob Faibussowitsch     if (size > 1) PetscCall(PetscFree(remoteVertices));
1709566063dSJacob Faibussowitsch     PetscCall(PetscFree(vertices));
171552f7358SJed Brown   } else {
172552f7358SJed Brown     PetscInt *localVertices, numSend = numCells+numCorners, k = 0;
173552f7358SJed Brown 
1749566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numSend, &localVertices));
175552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
1760298fd71SBarry Smith       PetscInt *closure = NULL;
17702281ff3SMatthew G. Knepley       PetscInt closureSize, value, nC = 0;
178552f7358SJed Brown 
17902281ff3SMatthew G. Knepley       if (label) {
1809566063dSJacob Faibussowitsch         PetscCall(DMLabelGetValue(label, c, &value));
181552f7358SJed Brown         if (value != 1) continue;
182552f7358SJed Brown       }
1839566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
184552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
185552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
186552f7358SJed Brown           const PetscInt gv = gvertex[closure[v] - vStart];
187552f7358SJed Brown           closure[nC++] = gv < 0 ? -(gv+1) : gv;
188552f7358SJed Brown         }
189552f7358SJed Brown       }
190552f7358SJed Brown       corners[numCells++] = nC;
191552f7358SJed Brown       localVertices[k++]  = nC;
192552f7358SJed Brown       for (v = 0; v < nC; ++v, ++k) {
193552f7358SJed Brown         localVertices[k] = closure[v];
194552f7358SJed Brown       }
1959566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
1969566063dSJacob Faibussowitsch       PetscCall(DMPlexReorderCell(dm, c, localVertices+k-nC));
197552f7358SJed Brown     }
19863a3b9bcSJacob Faibussowitsch     PetscCheck(k == numSend,PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid number of vertices to send %" PetscInt_FMT " should be %" PetscInt_FMT, k, numSend);
1999566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&numSend, 1, MPIU_INT, 0, tag, comm));
2009566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(localVertices, numSend, MPIU_INT, 0, tag, comm));
2019566063dSJacob Faibussowitsch     PetscCall(PetscFree(localVertices));
202552f7358SJed Brown   }
2039566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(globalVertexNumbers, &gvertex));
20463a3b9bcSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "CELL_TYPES %" PetscInt_FMT "\n", totCells));
205dd400576SPatrick Sanan   if (rank == 0) {
206552f7358SJed Brown     PetscInt cellType;
207552f7358SJed Brown 
208552f7358SJed Brown     for (c = 0; c < numCells; ++c) {
2099566063dSJacob Faibussowitsch       PetscCall(DMPlexVTKGetCellType_Internal(dm, dim, corners[c], &cellType));
21063a3b9bcSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "%" PetscInt_FMT "\n", cellType));
211552f7358SJed Brown     }
2129852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
213552f7358SJed Brown       MPI_Status status;
214552f7358SJed Brown 
2159566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status));
2169566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(corners, numCells, MPIU_INT, proc, tag, comm, &status));
217552f7358SJed Brown       for (c = 0; c < numCells; ++c) {
2189566063dSJacob Faibussowitsch         PetscCall(DMPlexVTKGetCellType_Internal(dm, dim, corners[c], &cellType));
21963a3b9bcSJacob Faibussowitsch         PetscCall(PetscFPrintf(comm, fp, "%" PetscInt_FMT "\n", cellType));
220552f7358SJed Brown       }
221552f7358SJed Brown     }
222552f7358SJed Brown   } else {
2239566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm));
2249566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(corners, numCells, MPIU_INT, 0, tag, comm));
225552f7358SJed Brown   }
2269566063dSJacob Faibussowitsch   PetscCall(PetscFree(corners));
227552f7358SJed Brown   *totalCells = totCells;
228552f7358SJed Brown   PetscFunctionReturn(0);
229552f7358SJed Brown }
230552f7358SJed Brown 
2317508e736SMatthew G. Knepley static PetscErrorCode DMPlexVTKWritePartition_ASCII(DM dm, FILE *fp)
232e8f6f0f6SMatthew G. Knepley {
233e8f6f0f6SMatthew G. Knepley   MPI_Comm       comm;
234e8f6f0f6SMatthew G. Knepley   PetscInt       numCells = 0, cellHeight;
2352f029394SStefano Zampini   PetscInt       numLabelCells, cStart, cEnd, c;
2369852e123SBarry Smith   PetscMPIInt    size, rank, proc, tag;
237e8f6f0f6SMatthew G. Knepley   PetscBool      hasLabel;
238e8f6f0f6SMatthew G. Knepley 
239e8f6f0f6SMatthew G. Knepley   PetscFunctionBegin;
2409566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
2419566063dSJacob Faibussowitsch   PetscCall(PetscCommGetNewTag(comm, &tag));
2429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2439566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2449566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
2459566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
2469566063dSJacob Faibussowitsch   PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells));
247e8f6f0f6SMatthew G. Knepley   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
248e8f6f0f6SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
249e8f6f0f6SMatthew G. Knepley     if (hasLabel) {
250e8f6f0f6SMatthew G. Knepley       PetscInt value;
251e8f6f0f6SMatthew G. Knepley 
2529566063dSJacob Faibussowitsch       PetscCall(DMGetLabelValue(dm, "vtk", c, &value));
253e8f6f0f6SMatthew G. Knepley       if (value != 1) continue;
254e8f6f0f6SMatthew G. Knepley     }
255e8f6f0f6SMatthew G. Knepley     ++numCells;
256e8f6f0f6SMatthew G. Knepley   }
257dd400576SPatrick Sanan   if (rank == 0) {
2589566063dSJacob Faibussowitsch     for (c = 0; c < numCells; ++c) PetscCall(PetscFPrintf(comm, fp, "%d\n", rank));
2599852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
260e8f6f0f6SMatthew G. Knepley       MPI_Status status;
261e8f6f0f6SMatthew G. Knepley 
2629566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status));
2639566063dSJacob Faibussowitsch       for (c = 0; c < numCells; ++c) PetscCall(PetscFPrintf(comm, fp, "%d\n", proc));
264e8f6f0f6SMatthew G. Knepley     }
265e8f6f0f6SMatthew G. Knepley   } else {
2669566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm));
267e8f6f0f6SMatthew G. Knepley   }
268e8f6f0f6SMatthew G. Knepley   PetscFunctionReturn(0);
269e8f6f0f6SMatthew G. Knepley }
270e8f6f0f6SMatthew G. Knepley 
27176b700caSToby Isaac #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
27276b700caSToby Isaac typedef double PetscVTKReal;
27376b700caSToby Isaac #elif defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL___FP16)
27476b700caSToby Isaac typedef float PetscVTKReal;
27576b700caSToby Isaac #else
27676b700caSToby Isaac typedef PetscReal PetscVTKReal;
27776b700caSToby Isaac #endif
27876b700caSToby Isaac 
27976b700caSToby Isaac static PetscErrorCode DMPlexVTKWriteSection_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec v, FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale, PetscInt imag)
2800adebc6cSBarry Smith {
28182f516ccSBarry Smith   MPI_Comm           comm;
282552f7358SJed Brown   const MPI_Datatype mpiType = MPIU_SCALAR;
283552f7358SJed Brown   PetscScalar        *array;
284552f7358SJed Brown   PetscInt           numDof = 0, maxDof;
285412e9a14SMatthew G. Knepley   PetscInt           numLabelCells, cellHeight, cStart, cEnd, numLabelVertices, vStart, vEnd, pStart, pEnd, p;
2869852e123SBarry Smith   PetscMPIInt        size, rank, proc, tag;
287552f7358SJed Brown   PetscBool          hasLabel;
288552f7358SJed Brown 
289552f7358SJed Brown   PetscFunctionBegin;
2909566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
291552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
292552f7358SJed Brown   PetscValidHeaderSpecific(v,VEC_CLASSID,4);
293552f7358SJed Brown   if (precision < 0) precision = 6;
2949566063dSJacob Faibussowitsch   PetscCall(PetscCommGetNewTag(comm, &tag));
2959566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2969566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2979566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
298552f7358SJed Brown   /* VTK only wants the values at cells or vertices */
2999566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
3009566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
3019566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
302552f7358SJed Brown   pStart   = PetscMax(PetscMin(cStart, vStart), pStart);
303552f7358SJed Brown   pEnd     = PetscMin(PetscMax(cEnd,   vEnd),   pEnd);
3049566063dSJacob Faibussowitsch   PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells));
3059566063dSJacob Faibussowitsch   PetscCall(DMGetStratumSize(dm, "vtk", 2, &numLabelVertices));
306552f7358SJed Brown   hasLabel = numLabelCells > 0 || numLabelVertices > 0 ? PETSC_TRUE : PETSC_FALSE;
307552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
308552f7358SJed Brown     /* Reject points not either cells or vertices */
309552f7358SJed Brown     if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
310552f7358SJed Brown     if (hasLabel) {
311552f7358SJed Brown       PetscInt value;
312552f7358SJed Brown 
313552f7358SJed Brown       if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
314552f7358SJed Brown           ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
3159566063dSJacob Faibussowitsch         PetscCall(DMGetLabelValue(dm, "vtk", p, &value));
316552f7358SJed Brown         if (value != 1) continue;
317552f7358SJed Brown       }
318552f7358SJed Brown     }
3199566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &numDof));
3208865f1eaSKarl Rupp     if (numDof) break;
321552f7358SJed Brown   }
3221c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, comm));
323552f7358SJed Brown   enforceDof = PetscMax(enforceDof, maxDof);
3249566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &array));
325dd400576SPatrick Sanan   if (rank == 0) {
32676b700caSToby Isaac     PetscVTKReal dval;
32776b700caSToby Isaac     PetscScalar  val;
328552f7358SJed Brown     char formatString[8];
329552f7358SJed Brown 
33063a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(formatString, 8, "%%.%" PetscInt_FMT "e", precision));
331552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
332552f7358SJed Brown       /* Here we lose a way to filter points by keeping them out of the Numbering */
333552f7358SJed Brown       PetscInt dof, off, goff, d;
334552f7358SJed Brown 
335552f7358SJed Brown       /* Reject points not either cells or vertices */
336552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
337552f7358SJed Brown       if (hasLabel) {
338552f7358SJed Brown         PetscInt value;
339552f7358SJed Brown 
340552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
341552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
3429566063dSJacob Faibussowitsch           PetscCall(DMGetLabelValue(dm, "vtk", p, &value));
343552f7358SJed Brown           if (value != 1) continue;
344552f7358SJed Brown         }
345552f7358SJed Brown       }
3469566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, p, &dof));
3479566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(section, p, &off));
3489566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(globalSection, p, &goff));
349552f7358SJed Brown       if (dof && goff >= 0) {
350552f7358SJed Brown         for (d = 0; d < dof; d++) {
351552f7358SJed Brown           if (d > 0) {
3529566063dSJacob Faibussowitsch             PetscCall(PetscFPrintf(comm, fp, " "));
353552f7358SJed Brown           }
35476b700caSToby Isaac           val = array[off+d];
35576b700caSToby Isaac           dval = (PetscVTKReal) ((imag ? PetscImaginaryPart(val) : PetscRealPart(val)) * scale);
3569566063dSJacob Faibussowitsch           PetscCall(PetscFPrintf(comm, fp, formatString, dval));
357552f7358SJed Brown         }
358552f7358SJed Brown         for (d = dof; d < enforceDof; d++) {
3599566063dSJacob Faibussowitsch           PetscCall(PetscFPrintf(comm, fp, " 0.0"));
360552f7358SJed Brown         }
3619566063dSJacob Faibussowitsch         PetscCall(PetscFPrintf(comm, fp, "\n"));
362552f7358SJed Brown       }
363552f7358SJed Brown     }
3649852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
365552f7358SJed Brown       PetscScalar *remoteValues;
366d892089bSMatthew G. Knepley       PetscInt    size = 0, d;
367552f7358SJed Brown       MPI_Status  status;
368552f7358SJed Brown 
3699566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(&size, 1, MPIU_INT, proc, tag, comm, &status));
3709566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(size, &remoteValues));
3719566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(remoteValues, size, mpiType, proc, tag, comm, &status));
372552f7358SJed Brown       for (p = 0; p < size/maxDof; ++p) {
373552f7358SJed Brown         for (d = 0; d < maxDof; ++d) {
374552f7358SJed Brown           if (d > 0) {
3759566063dSJacob Faibussowitsch             PetscCall(PetscFPrintf(comm, fp, " "));
376552f7358SJed Brown           }
37776b700caSToby Isaac           val = remoteValues[p*maxDof+d];
37876b700caSToby Isaac           dval = (PetscVTKReal) ((imag ? PetscImaginaryPart(val) : PetscRealPart(val)) * scale);
3799566063dSJacob Faibussowitsch           PetscCall(PetscFPrintf(comm, fp, formatString, dval));
380552f7358SJed Brown         }
381552f7358SJed Brown         for (d = maxDof; d < enforceDof; ++d) {
3829566063dSJacob Faibussowitsch           PetscCall(PetscFPrintf(comm, fp, " 0.0"));
383552f7358SJed Brown         }
3849566063dSJacob Faibussowitsch         PetscCall(PetscFPrintf(comm, fp, "\n"));
385552f7358SJed Brown       }
3869566063dSJacob Faibussowitsch       PetscCall(PetscFree(remoteValues));
387552f7358SJed Brown     }
388552f7358SJed Brown   } else {
389552f7358SJed Brown     PetscScalar *localValues;
390552f7358SJed Brown     PetscInt    size, k = 0;
391552f7358SJed Brown 
3929566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetStorageSize(section, &size));
3939566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &localValues));
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)) {
4049566063dSJacob Faibussowitsch           PetscCall(DMGetLabelValue(dm, "vtk", p, &value));
405552f7358SJed Brown           if (value != 1) continue;
406552f7358SJed Brown         }
407552f7358SJed Brown       }
4089566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, p, &dof));
4099566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(section, p, &off));
4109566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(globalSection, p, &goff));
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     }
4179566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&k, 1, MPIU_INT, 0, tag, comm));
4189566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(localValues, k, mpiType, 0, tag, comm));
4199566063dSJacob Faibussowitsch     PetscCall(PetscFree(localValues));
420552f7358SJed Brown   }
4219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &array));
422552f7358SJed Brown   PetscFunctionReturn(0);
423552f7358SJed Brown }
424552f7358SJed Brown 
42576b700caSToby 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)
426552f7358SJed Brown {
42782f516ccSBarry Smith   MPI_Comm       comm;
428552f7358SJed Brown   PetscInt       numDof = 0, maxDof;
429552f7358SJed Brown   PetscInt       pStart, pEnd, p;
430552f7358SJed Brown 
431552f7358SJed Brown   PetscFunctionBegin;
4329566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
4339566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
434552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
4359566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &numDof));
4368865f1eaSKarl Rupp     if (numDof) break;
437552f7358SJed Brown   }
438552f7358SJed Brown   numDof = PetscMax(numDof, enforceDof);
4391c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
440552f7358SJed Brown   if (!name) name = "Unknown";
441552f7358SJed Brown   if (maxDof == 3) {
44276b700caSToby Isaac     if (nameComplex) {
4439566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "VECTORS %s.%s double\n", name, imag ? "Im" : "Re"));
44476b700caSToby Isaac     } else {
4459566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "VECTORS %s double\n", name));
44676b700caSToby Isaac     }
44776b700caSToby Isaac   } else {
44876b700caSToby Isaac     if (nameComplex) {
44963a3b9bcSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "SCALARS %s.%s double %" PetscInt_FMT "\n", name, imag ? "Im" : "Re", maxDof));
450552f7358SJed Brown     } else {
45163a3b9bcSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "SCALARS %s double %" PetscInt_FMT "\n", name, maxDof));
45276b700caSToby Isaac     }
4539566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n"));
454552f7358SJed Brown   }
4559566063dSJacob Faibussowitsch   PetscCall(DMPlexVTKWriteSection_ASCII(dm, section, globalSection, field, fp, enforceDof, precision, scale, imag));
456552f7358SJed Brown   PetscFunctionReturn(0);
457552f7358SJed Brown }
458552f7358SJed Brown 
459552f7358SJed Brown static PetscErrorCode DMPlexVTKWriteAll_ASCII(DM dm, PetscViewer viewer)
460552f7358SJed Brown {
46182f516ccSBarry Smith   MPI_Comm                 comm;
462552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*) viewer->data;
463552f7358SJed Brown   FILE                     *fp;
464552f7358SJed Brown   PetscViewerVTKObjectLink link;
465412e9a14SMatthew G. Knepley   PetscInt                 totVertices, totCells = 0, loops_per_scalar, l;
46676b700caSToby Isaac   PetscBool                hasPoint = PETSC_FALSE, hasCell = PETSC_FALSE, writePartition = PETSC_FALSE, localized, writeComplex;
46752b05773SJed Brown   const char               *dmname;
468552f7358SJed Brown 
469552f7358SJed Brown   PetscFunctionBegin;
47076b700caSToby Isaac #if defined(PETSC_USE_COMPLEX)
47176b700caSToby Isaac   loops_per_scalar = 2;
47276b700caSToby Isaac   writeComplex = PETSC_TRUE;
47376b700caSToby Isaac #else
47476b700caSToby Isaac   loops_per_scalar = 1;
47576b700caSToby Isaac   writeComplex = PETSC_FALSE;
47676b700caSToby Isaac #endif
4779566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalized(dm,&localized));
4789566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
479*1dca8a05SBarry Smith   PetscCheck(!localized,comm,PETSC_ERR_SUP,"VTK output with localized coordinates not yet supported");
4809566063dSJacob Faibussowitsch   PetscCall(PetscFOpen(comm, vtk->filename, "wb", &fp));
4819566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)dm, &dmname));
4829566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "# vtk DataFile Version 2.0\n"));
4839566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "%s\n", dmname));
4849566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "ASCII\n"));
4859566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "DATASET UNSTRUCTURED_GRID\n"));
486552f7358SJed Brown   /* Vertices */
4879120ba30SVaclav Hapla   {
4889120ba30SVaclav Hapla     PetscSection  section, coordSection, globalCoordSection;
4899120ba30SVaclav Hapla     Vec           coordinates;
4909120ba30SVaclav Hapla     PetscReal     lengthScale;
4919120ba30SVaclav Hapla     DMLabel       label;
4929120ba30SVaclav Hapla     IS            vStratumIS;
4939120ba30SVaclav Hapla     PetscLayout   vLayout;
4949120ba30SVaclav Hapla 
4959566063dSJacob Faibussowitsch     PetscCall(DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale));
4969566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
4979566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthLabel(dm, &label));
4989566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(label, 0, &vStratumIS));
4999566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(dm, &section));                                 /* This section includes all points */
5006b57843cSMatthew G. Knepley     PetscCall(PetscSectionCreateSubdomainSection(section, vStratumIS, &coordSection)); /* This one includes just vertices */
5019566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreateGlobalSection(coordSection, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalCoordSection));
5029566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetPointLayout(comm, globalCoordSection, &vLayout));
5039566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetSize(vLayout, &totVertices));
50463a3b9bcSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "POINTS %" PetscInt_FMT " double\n", totVertices));
5059566063dSJacob Faibussowitsch     PetscCall(DMPlexVTKWriteSection_ASCII(dm, coordSection, globalCoordSection, coordinates, fp, 3, PETSC_DETERMINE, lengthScale, 0));
5069566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&vStratumIS));
5079566063dSJacob Faibussowitsch     PetscCall(PetscLayoutDestroy(&vLayout));
5089566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&coordSection));
5099566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&globalCoordSection));
5109120ba30SVaclav Hapla   }
511552f7358SJed Brown   /* Cells */
5129566063dSJacob Faibussowitsch   PetscCall(DMPlexVTKWriteCells_ASCII(dm, fp, &totCells));
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) {
51963a3b9bcSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "POINT_DATA %" PetscInt_FMT "\n", totVertices));
520552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
521552f7358SJed Brown       Vec          X = (Vec) link->vec;
522e630c359SToby Isaac       PetscSection section = NULL, globalSection, newSection = NULL;
523e630c359SToby Isaac       char         namebuf[256];
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;
5299566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName(link->vec, &name));
5309566063dSJacob Faibussowitsch       PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject*) &section));
531e630c359SToby Isaac       if (!section) {
532e630c359SToby Isaac         DM           dmX;
533e630c359SToby Isaac 
5349566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dmX));
535552f7358SJed Brown         if (dmX) {
536088580cfSMatthew G Knepley           DMLabel  subpointMap, subpointMapX;
537839dd189SMatthew G Knepley           PetscInt dim, dimX, pStart, pEnd, qStart, qEnd;
538839dd189SMatthew G Knepley 
5399566063dSJacob Faibussowitsch           PetscCall(DMGetLocalSection(dmX, &section));
540839dd189SMatthew G Knepley           /* Here is where we check whether dmX is a submesh of dm */
5419566063dSJacob Faibussowitsch           PetscCall(DMGetDimension(dm,  &dim));
5429566063dSJacob Faibussowitsch           PetscCall(DMGetDimension(dmX, &dimX));
5439566063dSJacob Faibussowitsch           PetscCall(DMPlexGetChart(dm,  &pStart, &pEnd));
5449566063dSJacob Faibussowitsch           PetscCall(DMPlexGetChart(dmX, &qStart, &qEnd));
5459566063dSJacob Faibussowitsch           PetscCall(DMPlexGetSubpointMap(dm,  &subpointMap));
5469566063dSJacob Faibussowitsch           PetscCall(DMPlexGetSubpointMap(dmX, &subpointMapX));
547839dd189SMatthew G Knepley           if (((dim != dimX) || ((pEnd-pStart) < (qEnd-qStart))) && subpointMap && !subpointMapX) {
548434e42b5SMatthew G. Knepley             const PetscInt *ind = NULL;
549e6ccafaeSMatthew G Knepley             IS              subpointIS;
550434e42b5SMatthew G. Knepley             PetscInt        n = 0, q;
551839dd189SMatthew G Knepley 
5529566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetChart(section, &qStart, &qEnd));
5539566063dSJacob Faibussowitsch             PetscCall(DMPlexGetSubpointIS(dm, &subpointIS));
554434e42b5SMatthew G. Knepley             if (subpointIS) {
5559566063dSJacob Faibussowitsch               PetscCall(ISGetLocalSize(subpointIS, &n));
5569566063dSJacob Faibussowitsch               PetscCall(ISGetIndices(subpointIS, &ind));
557434e42b5SMatthew G. Knepley             }
5589566063dSJacob Faibussowitsch             PetscCall(PetscSectionCreate(comm, &newSection));
5599566063dSJacob Faibussowitsch             PetscCall(PetscSectionSetChart(newSection, pStart, pEnd));
560839dd189SMatthew G Knepley             for (q = qStart; q < qEnd; ++q) {
561839dd189SMatthew G Knepley               PetscInt dof, off, p;
562839dd189SMatthew G Knepley 
5639566063dSJacob Faibussowitsch               PetscCall(PetscSectionGetDof(section, q, &dof));
564839dd189SMatthew G Knepley               if (dof) {
5659566063dSJacob Faibussowitsch                 PetscCall(PetscFindInt(q, n, ind, &p));
566839dd189SMatthew G Knepley                 if (p >= pStart) {
5679566063dSJacob Faibussowitsch                   PetscCall(PetscSectionSetDof(newSection, p, dof));
5689566063dSJacob Faibussowitsch                   PetscCall(PetscSectionGetOffset(section, q, &off));
5699566063dSJacob Faibussowitsch                   PetscCall(PetscSectionSetOffset(newSection, p, off));
570839dd189SMatthew G Knepley                 }
571839dd189SMatthew G Knepley               }
572839dd189SMatthew G Knepley             }
573434e42b5SMatthew G. Knepley             if (subpointIS) {
5749566063dSJacob Faibussowitsch               PetscCall(ISRestoreIndices(subpointIS, &ind));
575434e42b5SMatthew G. Knepley             }
576839dd189SMatthew G Knepley             /* No need to setup section */
577839dd189SMatthew G Knepley             section = newSection;
578839dd189SMatthew G Knepley           }
579552f7358SJed Brown         }
580e630c359SToby Isaac       }
58128b400f6SJacob Faibussowitsch       PetscCheck(section,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it and could not create one from VecGetDM()", name);
582e630c359SToby Isaac       if (link->field >= 0) {
583e630c359SToby Isaac         const char *fieldname;
584e630c359SToby Isaac 
5859566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetFieldName(section, link->field, &fieldname));
5869566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetField(section, link->field, &section));
587e630c359SToby Isaac         if (fieldname) {
5889566063dSJacob Faibussowitsch           PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s%s", name, fieldname));
589e630c359SToby Isaac         } else {
59063a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s%" PetscInt_FMT, name, link->field));
591e630c359SToby Isaac         }
5922fb742c9SToby Isaac       } else {
5939566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s", name));
594e630c359SToby Isaac       }
5959566063dSJacob Faibussowitsch       PetscCall(PetscViewerVTKSanitizeName_Internal(namebuf, sizeof(namebuf)));
5969566063dSJacob Faibussowitsch       PetscCall(PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection));
59776b700caSToby Isaac       for (l = 0; l < loops_per_scalar; l++) {
5989566063dSJacob Faibussowitsch         PetscCall(DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, namebuf, fp, enforceDof, PETSC_DETERMINE, 1.0, writeComplex, l));
59976b700caSToby Isaac       }
6009566063dSJacob Faibussowitsch       PetscCall(PetscSectionDestroy(&globalSection));
6019566063dSJacob Faibussowitsch       if (newSection) PetscCall(PetscSectionDestroy(&newSection));
602552f7358SJed Brown     }
603552f7358SJed Brown   }
604552f7358SJed Brown   /* Cell Fields */
6059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_view_partition", &writePartition, NULL));
606e8f6f0f6SMatthew G. Knepley   if (hasCell || writePartition) {
60763a3b9bcSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "CELL_DATA %" PetscInt_FMT "\n", totCells));
608552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
609552f7358SJed Brown       Vec          X = (Vec) link->vec;
610e630c359SToby Isaac       PetscSection section = NULL, globalSection;
6112fb742c9SToby Isaac       const char   *name = "";
612e630c359SToby Isaac       char         namebuf[256];
613552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
614552f7358SJed Brown 
615552f7358SJed Brown       if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
616552f7358SJed Brown       if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) enforceDof = 3;
6179566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName(link->vec, &name));
6189566063dSJacob Faibussowitsch       PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject*) &section));
619e630c359SToby Isaac       if (!section) {
620e630c359SToby Isaac         DM           dmX;
621e630c359SToby Isaac 
6229566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dmX));
623552f7358SJed Brown         if (dmX) {
6249566063dSJacob Faibussowitsch           PetscCall(DMGetLocalSection(dmX, &section));
625552f7358SJed Brown         }
626e630c359SToby Isaac       }
62728b400f6SJacob Faibussowitsch       PetscCheck(section,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it and could not create one from VecGetDM()", name);
628e630c359SToby Isaac       if (link->field >= 0) {
629e630c359SToby Isaac         const char *fieldname;
630e630c359SToby Isaac 
6319566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetFieldName(section, link->field, &fieldname));
6329566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetField(section, link->field, &section));
633e630c359SToby Isaac         if (fieldname) {
6349566063dSJacob Faibussowitsch           PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s%s", name, fieldname));
635e630c359SToby Isaac         } else {
63663a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s%" PetscInt_FMT, name, link->field));
637e630c359SToby Isaac         }
6382fb742c9SToby Isaac       } else {
6399566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s", name));
640e630c359SToby Isaac       }
6419566063dSJacob Faibussowitsch       PetscCall(PetscViewerVTKSanitizeName_Internal(namebuf, sizeof(namebuf)));
6429566063dSJacob Faibussowitsch       PetscCall(PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection));
64376b700caSToby Isaac       for (l = 0; l < loops_per_scalar; l++) {
6449566063dSJacob Faibussowitsch         PetscCall(DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, namebuf, fp, enforceDof, PETSC_DETERMINE, 1.0, writeComplex, l));
64576b700caSToby Isaac       }
6469566063dSJacob Faibussowitsch       PetscCall(PetscSectionDestroy(&globalSection));
647552f7358SJed Brown     }
648e8f6f0f6SMatthew G. Knepley     if (writePartition) {
6499566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "SCALARS partition int 1\n"));
6509566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n"));
6519566063dSJacob Faibussowitsch       PetscCall(DMPlexVTKWritePartition_ASCII(dm, fp));
652e8f6f0f6SMatthew G. Knepley     }
653552f7358SJed Brown   }
654552f7358SJed Brown   /* Cleanup */
6559566063dSJacob Faibussowitsch   PetscCall(PetscFClose(comm, fp));
656552f7358SJed Brown   PetscFunctionReturn(0);
657552f7358SJed Brown }
658552f7358SJed Brown 
659552f7358SJed Brown /*@C
660552f7358SJed Brown   DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
661552f7358SJed Brown 
662552f7358SJed Brown   Collective
663552f7358SJed Brown 
6644165533cSJose E. Roman   Input Parameters:
665552f7358SJed Brown + odm - The DMPlex specifying the mesh, passed as a PetscObject
666552f7358SJed Brown - viewer - viewer of type VTK
667552f7358SJed Brown 
668552f7358SJed Brown   Level: developer
669552f7358SJed Brown 
670552f7358SJed Brown   Note:
671552f7358SJed Brown   This function is a callback used by the VTK viewer to actually write the file.
672552f7358SJed 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.
673552f7358SJed Brown   Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
674552f7358SJed Brown 
675552f7358SJed Brown .seealso: PETSCVIEWERVTK
676552f7358SJed Brown @*/
677552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer)
678552f7358SJed Brown {
679552f7358SJed Brown   DM             dm = (DM) odm;
680552f7358SJed Brown   PetscBool      isvtk;
681552f7358SJed Brown 
682552f7358SJed Brown   PetscFunctionBegin;
683552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
684552f7358SJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
6859566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk));
68628b400f6SJacob Faibussowitsch   PetscCheck(isvtk,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name);
687552f7358SJed Brown   switch (viewer->format) {
6888ec8862eSJed Brown   case PETSC_VIEWER_ASCII_VTK_DEPRECATED:
6899566063dSJacob Faibussowitsch     PetscCall(DMPlexVTKWriteAll_ASCII(dm, viewer));
690552f7358SJed Brown     break;
691552f7358SJed Brown   case PETSC_VIEWER_VTK_VTU:
6929566063dSJacob Faibussowitsch     PetscCall(DMPlexVTKWriteAll_VTU(dm, viewer));
693552f7358SJed Brown     break;
69498921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]);
695552f7358SJed Brown   }
696552f7358SJed Brown   PetscFunctionReturn(0);
697552f7358SJed Brown }
698