xref: /petsc/src/dm/impls/plex/plexvtk.c (revision a1cb98fac0cdf0eb4d3e8a0c8b58f3fe8f800bc6)
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 
5d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexVTKGetCellType_Internal(DM dm, PetscInt dim, PetscInt corners, PetscInt *cellType)
6d71ae5a4SJacob Faibussowitsch {
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;
15d71ae5a4SJacob Faibussowitsch     default:
16d71ae5a4SJacob Faibussowitsch       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;
27d71ae5a4SJacob Faibussowitsch     default:
28d71ae5a4SJacob Faibussowitsch       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;
45d71ae5a4SJacob Faibussowitsch     default:
46d71ae5a4SJacob Faibussowitsch       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;
54cf4edabeSMatthew G. Knepley     case 5:
55cf4edabeSMatthew G. Knepley       *cellType = 14; /* VTK_PYRAMID */
56cf4edabeSMatthew G. Knepley       break;
572f029394SStefano Zampini     case 6:
582f029394SStefano Zampini       *cellType = 13; /* VTK_WEDGE */
592f029394SStefano Zampini       break;
60552f7358SJed Brown     case 8:
61552f7358SJed Brown       *cellType = 12; /* VTK_HEXAHEDRON */
62552f7358SJed Brown       break;
63552f7358SJed Brown     case 10:
64552f7358SJed Brown       *cellType = 24; /* VTK_QUADRATIC_TETRA */
65552f7358SJed Brown       break;
66552f7358SJed Brown     case 27:
67552f7358SJed Brown       *cellType = 29; /* VTK_QUADRATIC_HEXAHEDRON */
68552f7358SJed Brown       break;
69d71ae5a4SJacob Faibussowitsch     default:
70d71ae5a4SJacob Faibussowitsch       break;
71552f7358SJed Brown     }
72552f7358SJed Brown   }
73552f7358SJed Brown   PetscFunctionReturn(0);
74552f7358SJed Brown }
75552f7358SJed Brown 
76d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexVTKWriteCells_ASCII(DM dm, FILE *fp, PetscInt *totalCells)
77d71ae5a4SJacob Faibussowitsch {
7882f516ccSBarry Smith   MPI_Comm        comm;
7902281ff3SMatthew G. Knepley   DMLabel         label;
8002281ff3SMatthew G. Knepley   IS              globalVertexNumbers = NULL;
81552f7358SJed Brown   const PetscInt *gvertex;
82552f7358SJed Brown   PetscInt        dim;
83552f7358SJed Brown   PetscInt        numCorners = 0, totCorners = 0, maxCorners, *corners;
84552f7358SJed Brown   PetscInt        numCells = 0, totCells = 0, maxCells, cellHeight;
852f029394SStefano Zampini   PetscInt        numLabelCells, maxLabelCells, cStart, cEnd, c, vStart, vEnd, v;
869852e123SBarry Smith   PetscMPIInt     size, rank, proc, tag;
87552f7358SJed Brown 
88552f7358SJed Brown   PetscFunctionBegin;
899566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
909566063dSJacob Faibussowitsch   PetscCall(PetscCommGetNewTag(comm, &tag));
919566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
929566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
939566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
949566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
959566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
969566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
979566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "vtk", &label));
989566063dSJacob Faibussowitsch   PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells));
991c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&numLabelCells, &maxLabelCells, 1, MPIU_INT, MPI_MAX, comm));
100452b7979SMatthew G. Knepley   if (!maxLabelCells) label = NULL;
101552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
1020298fd71SBarry Smith     PetscInt *closure = NULL;
10302281ff3SMatthew G. Knepley     PetscInt  closureSize, value;
104552f7358SJed Brown 
10502281ff3SMatthew G. Knepley     if (label) {
1069566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(label, c, &value));
107552f7358SJed Brown       if (value != 1) continue;
108552f7358SJed Brown     }
1099566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
110552f7358SJed Brown     for (v = 0; v < closureSize * 2; v += 2) {
111552f7358SJed Brown       if ((closure[v] >= vStart) && (closure[v] < vEnd)) ++numCorners;
112552f7358SJed Brown     }
1139566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
114552f7358SJed Brown     ++numCells;
115552f7358SJed Brown   }
116552f7358SJed Brown   maxCells = numCells;
1179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Reduce(&numCells, &totCells, 1, MPIU_INT, MPI_SUM, 0, comm));
1189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Reduce(&numCells, &maxCells, 1, MPIU_INT, MPI_MAX, 0, comm));
1199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Reduce(&numCorners, &totCorners, 1, MPIU_INT, MPI_SUM, 0, comm));
1209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Reduce(&numCorners, &maxCorners, 1, MPIU_INT, MPI_MAX, 0, comm));
1219566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVertexNumbering(dm, &globalVertexNumbers));
1229566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(globalVertexNumbers, &gvertex));
1239566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxCells, &corners));
12463a3b9bcSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "CELLS %" PetscInt_FMT " %" PetscInt_FMT "\n", totCells, totCorners + totCells));
125dd400576SPatrick Sanan   if (rank == 0) {
126a4a685f2SJacob Faibussowitsch     PetscInt *remoteVertices, *vertices;
127552f7358SJed Brown 
1289566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(maxCorners, &vertices));
129552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
1300298fd71SBarry Smith       PetscInt *closure = NULL;
13102281ff3SMatthew G. Knepley       PetscInt  closureSize, value, nC = 0;
132552f7358SJed Brown 
13302281ff3SMatthew G. Knepley       if (label) {
1349566063dSJacob Faibussowitsch         PetscCall(DMLabelGetValue(label, c, &value));
135552f7358SJed Brown         if (value != 1) continue;
136552f7358SJed Brown       }
1379566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
138552f7358SJed Brown       for (v = 0; v < closureSize * 2; v += 2) {
139552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
140552f7358SJed Brown           const PetscInt gv = gvertex[closure[v] - vStart];
1410dcf9c70SMatthew G. Knepley           vertices[nC++]    = gv < 0 ? -(gv + 1) : gv;
142552f7358SJed Brown         }
143552f7358SJed Brown       }
1449566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
1459566063dSJacob Faibussowitsch       PetscCall(DMPlexReorderCell(dm, c, vertices));
146552f7358SJed Brown       corners[numCells++] = nC;
14763a3b9bcSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "%" PetscInt_FMT " ", nC));
14848a46eb9SPierre Jolivet       for (v = 0; v < nC; ++v) PetscCall(PetscFPrintf(comm, fp, " %" PetscInt_FMT, vertices[v]));
1499566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "\n"));
150552f7358SJed Brown     }
1519566063dSJacob Faibussowitsch     if (size > 1) PetscCall(PetscMalloc1(maxCorners + maxCells, &remoteVertices));
1529852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
153552f7358SJed Brown       MPI_Status status;
154552f7358SJed Brown 
1559566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(&numCorners, 1, MPIU_INT, proc, tag, comm, &status));
1569566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(remoteVertices, numCorners, MPIU_INT, proc, tag, comm, &status));
157552f7358SJed Brown       for (c = 0; c < numCorners;) {
158552f7358SJed Brown         PetscInt nC = remoteVertices[c++];
159552f7358SJed Brown 
160ad540459SPierre Jolivet         for (v = 0; v < nC; ++v, ++c) vertices[v] = remoteVertices[c];
16163a3b9bcSJacob Faibussowitsch         PetscCall(PetscFPrintf(comm, fp, "%" PetscInt_FMT " ", nC));
16248a46eb9SPierre Jolivet         for (v = 0; v < nC; ++v) PetscCall(PetscFPrintf(comm, fp, " %" PetscInt_FMT, vertices[v]));
1639566063dSJacob Faibussowitsch         PetscCall(PetscFPrintf(comm, fp, "\n"));
164552f7358SJed Brown       }
165552f7358SJed Brown     }
1669566063dSJacob Faibussowitsch     if (size > 1) PetscCall(PetscFree(remoteVertices));
1679566063dSJacob Faibussowitsch     PetscCall(PetscFree(vertices));
168552f7358SJed Brown   } else {
169552f7358SJed Brown     PetscInt *localVertices, numSend = numCells + numCorners, k = 0;
170552f7358SJed Brown 
1719566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numSend, &localVertices));
172552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
1730298fd71SBarry Smith       PetscInt *closure = NULL;
17402281ff3SMatthew G. Knepley       PetscInt  closureSize, value, nC = 0;
175552f7358SJed Brown 
17602281ff3SMatthew G. Knepley       if (label) {
1779566063dSJacob Faibussowitsch         PetscCall(DMLabelGetValue(label, c, &value));
178552f7358SJed Brown         if (value != 1) continue;
179552f7358SJed Brown       }
1809566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
181552f7358SJed Brown       for (v = 0; v < closureSize * 2; v += 2) {
182552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
183552f7358SJed Brown           const PetscInt gv = gvertex[closure[v] - vStart];
184552f7358SJed Brown           closure[nC++]     = gv < 0 ? -(gv + 1) : gv;
185552f7358SJed Brown         }
186552f7358SJed Brown       }
187552f7358SJed Brown       corners[numCells++] = nC;
188552f7358SJed Brown       localVertices[k++]  = nC;
189ad540459SPierre Jolivet       for (v = 0; v < nC; ++v, ++k) localVertices[k] = closure[v];
1909566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
1919566063dSJacob Faibussowitsch       PetscCall(DMPlexReorderCell(dm, c, localVertices + k - nC));
192552f7358SJed Brown     }
19363a3b9bcSJacob Faibussowitsch     PetscCheck(k == numSend, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of vertices to send %" PetscInt_FMT " should be %" PetscInt_FMT, k, numSend);
1949566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&numSend, 1, MPIU_INT, 0, tag, comm));
1959566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(localVertices, numSend, MPIU_INT, 0, tag, comm));
1969566063dSJacob Faibussowitsch     PetscCall(PetscFree(localVertices));
197552f7358SJed Brown   }
1989566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(globalVertexNumbers, &gvertex));
19963a3b9bcSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "CELL_TYPES %" PetscInt_FMT "\n", totCells));
200dd400576SPatrick Sanan   if (rank == 0) {
201552f7358SJed Brown     PetscInt cellType;
202552f7358SJed Brown 
203552f7358SJed Brown     for (c = 0; c < numCells; ++c) {
2049566063dSJacob Faibussowitsch       PetscCall(DMPlexVTKGetCellType_Internal(dm, dim, corners[c], &cellType));
20563a3b9bcSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "%" PetscInt_FMT "\n", cellType));
206552f7358SJed Brown     }
2079852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
208552f7358SJed Brown       MPI_Status status;
209552f7358SJed Brown 
2109566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status));
2119566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(corners, numCells, MPIU_INT, proc, tag, comm, &status));
212552f7358SJed Brown       for (c = 0; c < numCells; ++c) {
2139566063dSJacob Faibussowitsch         PetscCall(DMPlexVTKGetCellType_Internal(dm, dim, corners[c], &cellType));
21463a3b9bcSJacob Faibussowitsch         PetscCall(PetscFPrintf(comm, fp, "%" PetscInt_FMT "\n", cellType));
215552f7358SJed Brown       }
216552f7358SJed Brown     }
217552f7358SJed Brown   } else {
2189566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm));
2199566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(corners, numCells, MPIU_INT, 0, tag, comm));
220552f7358SJed Brown   }
2219566063dSJacob Faibussowitsch   PetscCall(PetscFree(corners));
222552f7358SJed Brown   *totalCells = totCells;
223552f7358SJed Brown   PetscFunctionReturn(0);
224552f7358SJed Brown }
225552f7358SJed Brown 
226d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexVTKWritePartition_ASCII(DM dm, FILE *fp)
227d71ae5a4SJacob Faibussowitsch {
228e8f6f0f6SMatthew G. Knepley   MPI_Comm    comm;
229e8f6f0f6SMatthew G. Knepley   PetscInt    numCells = 0, cellHeight;
2302f029394SStefano Zampini   PetscInt    numLabelCells, cStart, cEnd, c;
2319852e123SBarry Smith   PetscMPIInt size, rank, proc, tag;
232e8f6f0f6SMatthew G. Knepley   PetscBool   hasLabel;
233e8f6f0f6SMatthew G. Knepley 
234e8f6f0f6SMatthew G. Knepley   PetscFunctionBegin;
2359566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2369566063dSJacob Faibussowitsch   PetscCall(PetscCommGetNewTag(comm, &tag));
2379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2389566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2399566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
2409566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
2419566063dSJacob Faibussowitsch   PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells));
242e8f6f0f6SMatthew G. Knepley   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
243e8f6f0f6SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
244e8f6f0f6SMatthew G. Knepley     if (hasLabel) {
245e8f6f0f6SMatthew G. Knepley       PetscInt value;
246e8f6f0f6SMatthew G. Knepley 
2479566063dSJacob Faibussowitsch       PetscCall(DMGetLabelValue(dm, "vtk", c, &value));
248e8f6f0f6SMatthew G. Knepley       if (value != 1) continue;
249e8f6f0f6SMatthew G. Knepley     }
250e8f6f0f6SMatthew G. Knepley     ++numCells;
251e8f6f0f6SMatthew G. Knepley   }
252dd400576SPatrick Sanan   if (rank == 0) {
2539566063dSJacob Faibussowitsch     for (c = 0; c < numCells; ++c) PetscCall(PetscFPrintf(comm, fp, "%d\n", rank));
2549852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
255e8f6f0f6SMatthew G. Knepley       MPI_Status status;
256e8f6f0f6SMatthew G. Knepley 
2579566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status));
2589566063dSJacob Faibussowitsch       for (c = 0; c < numCells; ++c) PetscCall(PetscFPrintf(comm, fp, "%d\n", proc));
259e8f6f0f6SMatthew G. Knepley     }
260e8f6f0f6SMatthew G. Knepley   } else {
2619566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm));
262e8f6f0f6SMatthew G. Knepley   }
263e8f6f0f6SMatthew G. Knepley   PetscFunctionReturn(0);
264e8f6f0f6SMatthew G. Knepley }
265e8f6f0f6SMatthew G. Knepley 
26676b700caSToby Isaac #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
26776b700caSToby Isaac typedef double PetscVTKReal;
26876b700caSToby Isaac #elif defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL___FP16)
26976b700caSToby Isaac typedef float PetscVTKReal;
27076b700caSToby Isaac #else
27176b700caSToby Isaac typedef PetscReal PetscVTKReal;
27276b700caSToby Isaac #endif
27376b700caSToby Isaac 
274d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexVTKWriteSection_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec v, FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale, PetscInt imag)
275d71ae5a4SJacob Faibussowitsch {
27682f516ccSBarry Smith   MPI_Comm           comm;
277552f7358SJed Brown   const MPI_Datatype mpiType = MPIU_SCALAR;
278552f7358SJed Brown   PetscScalar       *array;
279552f7358SJed Brown   PetscInt           numDof = 0, maxDof;
280412e9a14SMatthew G. Knepley   PetscInt           numLabelCells, cellHeight, cStart, cEnd, numLabelVertices, vStart, vEnd, pStart, pEnd, p;
2819852e123SBarry Smith   PetscMPIInt        size, rank, proc, tag;
282552f7358SJed Brown   PetscBool          hasLabel;
283552f7358SJed Brown 
284552f7358SJed Brown   PetscFunctionBegin;
2859566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
286552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
287552f7358SJed Brown   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
288552f7358SJed Brown   if (precision < 0) precision = 6;
2899566063dSJacob Faibussowitsch   PetscCall(PetscCommGetNewTag(comm, &tag));
2909566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2919566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2929566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
293552f7358SJed Brown   /* VTK only wants the values at cells or vertices */
2949566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
2959566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
2969566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
297552f7358SJed Brown   pStart = PetscMax(PetscMin(cStart, vStart), pStart);
298552f7358SJed Brown   pEnd   = PetscMin(PetscMax(cEnd, vEnd), pEnd);
2999566063dSJacob Faibussowitsch   PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells));
3009566063dSJacob Faibussowitsch   PetscCall(DMGetStratumSize(dm, "vtk", 2, &numLabelVertices));
301552f7358SJed Brown   hasLabel = numLabelCells > 0 || numLabelVertices > 0 ? PETSC_TRUE : PETSC_FALSE;
302552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
303552f7358SJed Brown     /* Reject points not either cells or vertices */
304552f7358SJed Brown     if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
305552f7358SJed Brown     if (hasLabel) {
306552f7358SJed Brown       PetscInt value;
307552f7358SJed Brown 
3089371c9d4SSatish Balay       if (((p >= cStart) && (p < cEnd) && numLabelCells) || ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
3099566063dSJacob Faibussowitsch         PetscCall(DMGetLabelValue(dm, "vtk", p, &value));
310552f7358SJed Brown         if (value != 1) continue;
311552f7358SJed Brown       }
312552f7358SJed Brown     }
3139566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &numDof));
3148865f1eaSKarl Rupp     if (numDof) break;
315552f7358SJed Brown   }
3161c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, comm));
317552f7358SJed Brown   enforceDof = PetscMax(enforceDof, maxDof);
3189566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &array));
319dd400576SPatrick Sanan   if (rank == 0) {
32076b700caSToby Isaac     PetscVTKReal dval;
32176b700caSToby Isaac     PetscScalar  val;
322552f7358SJed Brown     char         formatString[8];
323552f7358SJed Brown 
32463a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(formatString, 8, "%%.%" PetscInt_FMT "e", precision));
325552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
326552f7358SJed Brown       /* Here we lose a way to filter points by keeping them out of the Numbering */
327552f7358SJed Brown       PetscInt dof, off, goff, d;
328552f7358SJed Brown 
329552f7358SJed Brown       /* Reject points not either cells or vertices */
330552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
331552f7358SJed Brown       if (hasLabel) {
332552f7358SJed Brown         PetscInt value;
333552f7358SJed Brown 
3349371c9d4SSatish Balay         if (((p >= cStart) && (p < cEnd) && numLabelCells) || ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
3359566063dSJacob Faibussowitsch           PetscCall(DMGetLabelValue(dm, "vtk", p, &value));
336552f7358SJed Brown           if (value != 1) continue;
337552f7358SJed Brown         }
338552f7358SJed Brown       }
3399566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, p, &dof));
3409566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(section, p, &off));
3419566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(globalSection, p, &goff));
342552f7358SJed Brown       if (dof && goff >= 0) {
343552f7358SJed Brown         for (d = 0; d < dof; d++) {
34448a46eb9SPierre Jolivet           if (d > 0) PetscCall(PetscFPrintf(comm, fp, " "));
34576b700caSToby Isaac           val  = array[off + d];
34676b700caSToby Isaac           dval = (PetscVTKReal)((imag ? PetscImaginaryPart(val) : PetscRealPart(val)) * scale);
3479566063dSJacob Faibussowitsch           PetscCall(PetscFPrintf(comm, fp, formatString, dval));
348552f7358SJed Brown         }
34948a46eb9SPierre Jolivet         for (d = dof; d < enforceDof; d++) PetscCall(PetscFPrintf(comm, fp, " 0.0"));
3509566063dSJacob Faibussowitsch         PetscCall(PetscFPrintf(comm, fp, "\n"));
351552f7358SJed Brown       }
352552f7358SJed Brown     }
3539852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
354552f7358SJed Brown       PetscScalar *remoteValues;
355d892089bSMatthew G. Knepley       PetscInt     size = 0, d;
356552f7358SJed Brown       MPI_Status   status;
357552f7358SJed Brown 
3589566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(&size, 1, MPIU_INT, proc, tag, comm, &status));
3599566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(size, &remoteValues));
3609566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(remoteValues, size, mpiType, proc, tag, comm, &status));
361552f7358SJed Brown       for (p = 0; p < size / maxDof; ++p) {
362552f7358SJed Brown         for (d = 0; d < maxDof; ++d) {
36348a46eb9SPierre Jolivet           if (d > 0) PetscCall(PetscFPrintf(comm, fp, " "));
36476b700caSToby Isaac           val  = remoteValues[p * maxDof + d];
36576b700caSToby Isaac           dval = (PetscVTKReal)((imag ? PetscImaginaryPart(val) : PetscRealPart(val)) * scale);
3669566063dSJacob Faibussowitsch           PetscCall(PetscFPrintf(comm, fp, formatString, dval));
367552f7358SJed Brown         }
36848a46eb9SPierre Jolivet         for (d = maxDof; d < enforceDof; ++d) PetscCall(PetscFPrintf(comm, fp, " 0.0"));
3699566063dSJacob Faibussowitsch         PetscCall(PetscFPrintf(comm, fp, "\n"));
370552f7358SJed Brown       }
3719566063dSJacob Faibussowitsch       PetscCall(PetscFree(remoteValues));
372552f7358SJed Brown     }
373552f7358SJed Brown   } else {
374552f7358SJed Brown     PetscScalar *localValues;
375552f7358SJed Brown     PetscInt     size, k = 0;
376552f7358SJed Brown 
3779566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetStorageSize(section, &size));
3789566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &localValues));
379552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
380552f7358SJed Brown       PetscInt dof, off, goff, d;
381552f7358SJed Brown 
382552f7358SJed Brown       /* Reject points not either cells or vertices */
383552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
384552f7358SJed Brown       if (hasLabel) {
385552f7358SJed Brown         PetscInt value;
386552f7358SJed Brown 
3879371c9d4SSatish Balay         if (((p >= cStart) && (p < cEnd) && numLabelCells) || ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
3889566063dSJacob Faibussowitsch           PetscCall(DMGetLabelValue(dm, "vtk", p, &value));
389552f7358SJed Brown           if (value != 1) continue;
390552f7358SJed Brown         }
391552f7358SJed Brown       }
3929566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, p, &dof));
3939566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(section, p, &off));
3949566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(globalSection, p, &goff));
395552f7358SJed Brown       if (goff >= 0) {
396ad540459SPierre Jolivet         for (d = 0; d < dof; ++d) localValues[k++] = array[off + d];
397552f7358SJed Brown       }
398552f7358SJed Brown     }
3999566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&k, 1, MPIU_INT, 0, tag, comm));
4009566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(localValues, k, mpiType, 0, tag, comm));
4019566063dSJacob Faibussowitsch     PetscCall(PetscFree(localValues));
402552f7358SJed Brown   }
4039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &array));
404552f7358SJed Brown   PetscFunctionReturn(0);
405552f7358SJed Brown }
406552f7358SJed Brown 
407d71ae5a4SJacob Faibussowitsch 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)
408d71ae5a4SJacob Faibussowitsch {
40982f516ccSBarry Smith   MPI_Comm comm;
410552f7358SJed Brown   PetscInt numDof = 0, maxDof;
411552f7358SJed Brown   PetscInt pStart, pEnd, p;
412552f7358SJed Brown 
413552f7358SJed Brown   PetscFunctionBegin;
4149566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
4159566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
416552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
4179566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &numDof));
4188865f1eaSKarl Rupp     if (numDof) break;
419552f7358SJed Brown   }
420552f7358SJed Brown   numDof = PetscMax(numDof, enforceDof);
4211c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
422552f7358SJed Brown   if (!name) name = "Unknown";
423552f7358SJed Brown   if (maxDof == 3) {
42476b700caSToby Isaac     if (nameComplex) {
4259566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "VECTORS %s.%s double\n", name, imag ? "Im" : "Re"));
42676b700caSToby Isaac     } else {
4279566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "VECTORS %s double\n", name));
42876b700caSToby Isaac     }
42976b700caSToby Isaac   } else {
43076b700caSToby Isaac     if (nameComplex) {
43163a3b9bcSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "SCALARS %s.%s double %" PetscInt_FMT "\n", name, imag ? "Im" : "Re", maxDof));
432552f7358SJed Brown     } else {
43363a3b9bcSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "SCALARS %s double %" PetscInt_FMT "\n", name, maxDof));
43476b700caSToby Isaac     }
4359566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n"));
436552f7358SJed Brown   }
4379566063dSJacob Faibussowitsch   PetscCall(DMPlexVTKWriteSection_ASCII(dm, section, globalSection, field, fp, enforceDof, precision, scale, imag));
438552f7358SJed Brown   PetscFunctionReturn(0);
439552f7358SJed Brown }
440552f7358SJed Brown 
441d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexVTKWriteAll_ASCII(DM dm, PetscViewer viewer)
442d71ae5a4SJacob Faibussowitsch {
44382f516ccSBarry Smith   MPI_Comm                 comm;
444552f7358SJed Brown   PetscViewer_VTK         *vtk = (PetscViewer_VTK *)viewer->data;
445552f7358SJed Brown   FILE                    *fp;
446552f7358SJed Brown   PetscViewerVTKObjectLink link;
447412e9a14SMatthew G. Knepley   PetscInt                 totVertices, totCells = 0, loops_per_scalar, l;
44876b700caSToby Isaac   PetscBool                hasPoint = PETSC_FALSE, hasCell = PETSC_FALSE, writePartition = PETSC_FALSE, localized, writeComplex;
44952b05773SJed Brown   const char              *dmname;
450552f7358SJed Brown 
451552f7358SJed Brown   PetscFunctionBegin;
45276b700caSToby Isaac #if defined(PETSC_USE_COMPLEX)
45376b700caSToby Isaac   loops_per_scalar = 2;
45476b700caSToby Isaac   writeComplex     = PETSC_TRUE;
45576b700caSToby Isaac #else
45676b700caSToby Isaac   loops_per_scalar = 1;
45776b700caSToby Isaac   writeComplex     = PETSC_FALSE;
45876b700caSToby Isaac #endif
4599566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalized(dm, &localized));
4609566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
4611dca8a05SBarry Smith   PetscCheck(!localized, comm, PETSC_ERR_SUP, "VTK output with localized coordinates not yet supported");
4629566063dSJacob Faibussowitsch   PetscCall(PetscFOpen(comm, vtk->filename, "wb", &fp));
4639566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)dm, &dmname));
4649566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "# vtk DataFile Version 2.0\n"));
4659566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "%s\n", dmname));
4669566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "ASCII\n"));
4679566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "DATASET UNSTRUCTURED_GRID\n"));
468552f7358SJed Brown   /* Vertices */
4699120ba30SVaclav Hapla   {
4709120ba30SVaclav Hapla     PetscSection section, coordSection, globalCoordSection;
4719120ba30SVaclav Hapla     Vec          coordinates;
4729120ba30SVaclav Hapla     PetscReal    lengthScale;
4739120ba30SVaclav Hapla     DMLabel      label;
4749120ba30SVaclav Hapla     IS           vStratumIS;
4759120ba30SVaclav Hapla     PetscLayout  vLayout;
4769120ba30SVaclav Hapla 
4779566063dSJacob Faibussowitsch     PetscCall(DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale));
4789566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
4799566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthLabel(dm, &label));
4809566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(label, 0, &vStratumIS));
4819566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(dm, &section));                                   /* This section includes all points */
4826b57843cSMatthew G. Knepley     PetscCall(PetscSectionCreateSubdomainSection(section, vStratumIS, &coordSection)); /* This one includes just vertices */
4839566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreateGlobalSection(coordSection, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalCoordSection));
4849566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetPointLayout(comm, globalCoordSection, &vLayout));
4859566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetSize(vLayout, &totVertices));
48663a3b9bcSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "POINTS %" PetscInt_FMT " double\n", totVertices));
4879566063dSJacob Faibussowitsch     PetscCall(DMPlexVTKWriteSection_ASCII(dm, coordSection, globalCoordSection, coordinates, fp, 3, PETSC_DETERMINE, lengthScale, 0));
4889566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&vStratumIS));
4899566063dSJacob Faibussowitsch     PetscCall(PetscLayoutDestroy(&vLayout));
4909566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&coordSection));
4919566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&globalCoordSection));
4929120ba30SVaclav Hapla   }
493552f7358SJed Brown   /* Cells */
4949566063dSJacob Faibussowitsch   PetscCall(DMPlexVTKWriteCells_ASCII(dm, fp, &totCells));
495552f7358SJed Brown   /* Vertex fields */
496552f7358SJed Brown   for (link = vtk->link; link; link = link->next) {
497552f7358SJed Brown     if ((link->ft == PETSC_VTK_POINT_FIELD) || (link->ft == PETSC_VTK_POINT_VECTOR_FIELD)) hasPoint = PETSC_TRUE;
498552f7358SJed Brown     if ((link->ft == PETSC_VTK_CELL_FIELD) || (link->ft == PETSC_VTK_CELL_VECTOR_FIELD)) hasCell = PETSC_TRUE;
499552f7358SJed Brown   }
500552f7358SJed Brown   if (hasPoint) {
50163a3b9bcSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "POINT_DATA %" PetscInt_FMT "\n", totVertices));
502552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
503552f7358SJed Brown       Vec          X       = (Vec)link->vec;
504e630c359SToby Isaac       PetscSection section = NULL, globalSection, newSection = NULL;
505e630c359SToby Isaac       char         namebuf[256];
506552f7358SJed Brown       const char  *name;
507552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
508552f7358SJed Brown 
509552f7358SJed Brown       if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
510552f7358SJed Brown       if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) enforceDof = 3;
5119566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName(link->vec, &name));
5129566063dSJacob Faibussowitsch       PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)&section));
513e630c359SToby Isaac       if (!section) {
514e630c359SToby Isaac         DM dmX;
515e630c359SToby Isaac 
5169566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dmX));
517552f7358SJed Brown         if (dmX) {
518088580cfSMatthew G Knepley           DMLabel  subpointMap, subpointMapX;
519839dd189SMatthew G Knepley           PetscInt dim, dimX, pStart, pEnd, qStart, qEnd;
520839dd189SMatthew G Knepley 
5219566063dSJacob Faibussowitsch           PetscCall(DMGetLocalSection(dmX, &section));
522839dd189SMatthew G Knepley           /* Here is where we check whether dmX is a submesh of dm */
5239566063dSJacob Faibussowitsch           PetscCall(DMGetDimension(dm, &dim));
5249566063dSJacob Faibussowitsch           PetscCall(DMGetDimension(dmX, &dimX));
5259566063dSJacob Faibussowitsch           PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
5269566063dSJacob Faibussowitsch           PetscCall(DMPlexGetChart(dmX, &qStart, &qEnd));
5279566063dSJacob Faibussowitsch           PetscCall(DMPlexGetSubpointMap(dm, &subpointMap));
5289566063dSJacob Faibussowitsch           PetscCall(DMPlexGetSubpointMap(dmX, &subpointMapX));
529839dd189SMatthew G Knepley           if (((dim != dimX) || ((pEnd - pStart) < (qEnd - qStart))) && subpointMap && !subpointMapX) {
530434e42b5SMatthew G. Knepley             const PetscInt *ind = NULL;
531e6ccafaeSMatthew G Knepley             IS              subpointIS;
532434e42b5SMatthew G. Knepley             PetscInt        n = 0, q;
533839dd189SMatthew G Knepley 
5349566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetChart(section, &qStart, &qEnd));
5359566063dSJacob Faibussowitsch             PetscCall(DMPlexGetSubpointIS(dm, &subpointIS));
536434e42b5SMatthew G. Knepley             if (subpointIS) {
5379566063dSJacob Faibussowitsch               PetscCall(ISGetLocalSize(subpointIS, &n));
5389566063dSJacob Faibussowitsch               PetscCall(ISGetIndices(subpointIS, &ind));
539434e42b5SMatthew G. Knepley             }
5409566063dSJacob Faibussowitsch             PetscCall(PetscSectionCreate(comm, &newSection));
5419566063dSJacob Faibussowitsch             PetscCall(PetscSectionSetChart(newSection, pStart, pEnd));
542839dd189SMatthew G Knepley             for (q = qStart; q < qEnd; ++q) {
543839dd189SMatthew G Knepley               PetscInt dof, off, p;
544839dd189SMatthew G Knepley 
5459566063dSJacob Faibussowitsch               PetscCall(PetscSectionGetDof(section, q, &dof));
546839dd189SMatthew G Knepley               if (dof) {
5479566063dSJacob Faibussowitsch                 PetscCall(PetscFindInt(q, n, ind, &p));
548839dd189SMatthew G Knepley                 if (p >= pStart) {
5499566063dSJacob Faibussowitsch                   PetscCall(PetscSectionSetDof(newSection, p, dof));
5509566063dSJacob Faibussowitsch                   PetscCall(PetscSectionGetOffset(section, q, &off));
5519566063dSJacob Faibussowitsch                   PetscCall(PetscSectionSetOffset(newSection, p, off));
552839dd189SMatthew G Knepley                 }
553839dd189SMatthew G Knepley               }
554839dd189SMatthew G Knepley             }
55548a46eb9SPierre Jolivet             if (subpointIS) PetscCall(ISRestoreIndices(subpointIS, &ind));
556839dd189SMatthew G Knepley             /* No need to setup section */
557839dd189SMatthew G Knepley             section = newSection;
558839dd189SMatthew G Knepley           }
559552f7358SJed Brown         }
560e630c359SToby Isaac       }
56128b400f6SJacob 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);
562e630c359SToby Isaac       if (link->field >= 0) {
563e630c359SToby Isaac         const char *fieldname;
564e630c359SToby Isaac 
5659566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetFieldName(section, link->field, &fieldname));
5669566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetField(section, link->field, &section));
567e630c359SToby Isaac         if (fieldname) {
5689566063dSJacob Faibussowitsch           PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s%s", name, fieldname));
569e630c359SToby Isaac         } else {
57063a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s%" PetscInt_FMT, name, link->field));
571e630c359SToby Isaac         }
5722fb742c9SToby Isaac       } else {
5739566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s", name));
574e630c359SToby Isaac       }
5759566063dSJacob Faibussowitsch       PetscCall(PetscViewerVTKSanitizeName_Internal(namebuf, sizeof(namebuf)));
5769566063dSJacob Faibussowitsch       PetscCall(PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection));
57748a46eb9SPierre Jolivet       for (l = 0; l < loops_per_scalar; l++) PetscCall(DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, namebuf, fp, enforceDof, PETSC_DETERMINE, 1.0, writeComplex, l));
5789566063dSJacob Faibussowitsch       PetscCall(PetscSectionDestroy(&globalSection));
5799566063dSJacob Faibussowitsch       if (newSection) PetscCall(PetscSectionDestroy(&newSection));
580552f7358SJed Brown     }
581552f7358SJed Brown   }
582552f7358SJed Brown   /* Cell Fields */
5839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_view_partition", &writePartition, NULL));
584e8f6f0f6SMatthew G. Knepley   if (hasCell || writePartition) {
58563a3b9bcSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "CELL_DATA %" PetscInt_FMT "\n", totCells));
586552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
587552f7358SJed Brown       Vec          X       = (Vec)link->vec;
588e630c359SToby Isaac       PetscSection section = NULL, globalSection;
5892fb742c9SToby Isaac       const char  *name    = "";
590e630c359SToby Isaac       char         namebuf[256];
591552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
592552f7358SJed Brown 
593552f7358SJed Brown       if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
594552f7358SJed Brown       if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) enforceDof = 3;
5959566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName(link->vec, &name));
5969566063dSJacob Faibussowitsch       PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)&section));
597e630c359SToby Isaac       if (!section) {
598e630c359SToby Isaac         DM dmX;
599e630c359SToby Isaac 
6009566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dmX));
60148a46eb9SPierre Jolivet         if (dmX) PetscCall(DMGetLocalSection(dmX, &section));
602e630c359SToby Isaac       }
60328b400f6SJacob 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);
604e630c359SToby Isaac       if (link->field >= 0) {
605e630c359SToby Isaac         const char *fieldname;
606e630c359SToby Isaac 
6079566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetFieldName(section, link->field, &fieldname));
6089566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetField(section, link->field, &section));
609e630c359SToby Isaac         if (fieldname) {
6109566063dSJacob Faibussowitsch           PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s%s", name, fieldname));
611e630c359SToby Isaac         } else {
61263a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s%" PetscInt_FMT, name, link->field));
613e630c359SToby Isaac         }
6142fb742c9SToby Isaac       } else {
6159566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s", name));
616e630c359SToby Isaac       }
6179566063dSJacob Faibussowitsch       PetscCall(PetscViewerVTKSanitizeName_Internal(namebuf, sizeof(namebuf)));
6189566063dSJacob Faibussowitsch       PetscCall(PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection));
61948a46eb9SPierre Jolivet       for (l = 0; l < loops_per_scalar; l++) PetscCall(DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, namebuf, fp, enforceDof, PETSC_DETERMINE, 1.0, writeComplex, l));
6209566063dSJacob Faibussowitsch       PetscCall(PetscSectionDestroy(&globalSection));
621552f7358SJed Brown     }
622e8f6f0f6SMatthew G. Knepley     if (writePartition) {
6239566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "SCALARS partition int 1\n"));
6249566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n"));
6259566063dSJacob Faibussowitsch       PetscCall(DMPlexVTKWritePartition_ASCII(dm, fp));
626e8f6f0f6SMatthew G. Knepley     }
627552f7358SJed Brown   }
628552f7358SJed Brown   /* Cleanup */
6299566063dSJacob Faibussowitsch   PetscCall(PetscFClose(comm, fp));
630552f7358SJed Brown   PetscFunctionReturn(0);
631552f7358SJed Brown }
632552f7358SJed Brown 
633552f7358SJed Brown /*@C
634552f7358SJed Brown   DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
635552f7358SJed Brown 
636552f7358SJed Brown   Collective
637552f7358SJed Brown 
6384165533cSJose E. Roman   Input Parameters:
639*a1cb98faSBarry Smith + odm - The `DMPLEX` specifying the mesh, passed as a `PetscObject`
640*a1cb98faSBarry Smith - viewer - viewer of type `PETSCVIEWERVTK`
641552f7358SJed Brown 
642552f7358SJed Brown   Level: developer
643552f7358SJed Brown 
644552f7358SJed Brown   Note:
645552f7358SJed Brown   This function is a callback used by the VTK viewer to actually write the file.
646552f7358SJed 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.
647552f7358SJed Brown   Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
648552f7358SJed Brown 
649*a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `PETSCVIEWEREXODUSII`, `DMPLEX`, `PETSCVIEWERVTK`
650552f7358SJed Brown @*/
651d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer)
652d71ae5a4SJacob Faibussowitsch {
653552f7358SJed Brown   DM dm = (DM)odm;
654552f7358SJed Brown 
655552f7358SJed Brown   PetscFunctionBegin;
656552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
657*a1cb98faSBarry Smith   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 2, PETSCVIEWERVTK);
658552f7358SJed Brown   switch (viewer->format) {
659d71ae5a4SJacob Faibussowitsch   case PETSC_VIEWER_ASCII_VTK_DEPRECATED:
660d71ae5a4SJacob Faibussowitsch     PetscCall(DMPlexVTKWriteAll_ASCII(dm, viewer));
661d71ae5a4SJacob Faibussowitsch     break;
662d71ae5a4SJacob Faibussowitsch   case PETSC_VIEWER_VTK_VTU:
663d71ae5a4SJacob Faibussowitsch     PetscCall(DMPlexVTKWriteAll_VTU(dm, viewer));
664d71ae5a4SJacob Faibussowitsch     break;
665d71ae5a4SJacob Faibussowitsch   default:
666d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]);
667552f7358SJed Brown   }
668552f7358SJed Brown   PetscFunctionReturn(0);
669552f7358SJed Brown }
670