xref: /petsc/src/dm/impls/plex/plexvtk.c (revision 98921bda46e76d7aaed9e0138c5ff9d0ce93f355)
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);
89ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
90ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(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);
97820f2d46SBarry Smith   ierr = MPIU_Allreduce(&numLabelCells, &maxLabelCells, 1, MPIU_INT, MPI_MAX, comm);CHKERRMPI(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;
11555b25c41SPierre Jolivet   ierr     = MPI_Reduce(&numCells, &totCells, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRMPI(ierr);
11655b25c41SPierre Jolivet   ierr     = MPI_Reduce(&numCells, &maxCells, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRMPI(ierr);
11755b25c41SPierre Jolivet   ierr     = MPI_Reduce(&numCorners, &totCorners, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRMPI(ierr);
11855b25c41SPierre Jolivet   ierr     = MPI_Reduce(&numCorners, &maxCorners, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRMPI(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);
123dd400576SPatrick Sanan   if (rank == 0) {
124a4a685f2SJacob Faibussowitsch     PetscInt *remoteVertices, *vertices;
125552f7358SJed Brown 
126785e854fSJed Brown     ierr = PetscMalloc1(maxCorners, &vertices);CHKERRQ(ierr);
127552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
1280298fd71SBarry Smith       PetscInt *closure = NULL;
12902281ff3SMatthew G. Knepley       PetscInt closureSize, value, nC = 0;
130552f7358SJed Brown 
13102281ff3SMatthew G. Knepley       if (label) {
13202281ff3SMatthew G. Knepley         ierr = DMLabelGetValue(label, c, &value);CHKERRQ(ierr);
133552f7358SJed Brown         if (value != 1) continue;
134552f7358SJed Brown       }
135552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
136552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
137552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
138552f7358SJed Brown           const PetscInt gv = gvertex[closure[v] - vStart];
1390dcf9c70SMatthew G. Knepley           vertices[nC++] = gv < 0 ? -(gv+1) : gv;
140552f7358SJed Brown         }
141552f7358SJed Brown       }
1420dcf9c70SMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
14396ca5757SLisandro Dalcin       ierr = DMPlexReorderCell(dm, c, vertices);CHKERRQ(ierr);
144552f7358SJed Brown       corners[numCells++] = nC;
1456b8858a1SBarry Smith       ierr = PetscFPrintf(comm, fp, "%D ", nC);CHKERRQ(ierr);
146552f7358SJed Brown       for (v = 0; v < nC; ++v) {
14796ca5757SLisandro Dalcin         ierr = PetscFPrintf(comm, fp, " %D", vertices[v]);CHKERRQ(ierr);
148552f7358SJed Brown       }
149552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
150552f7358SJed Brown     }
1519852e123SBarry Smith     if (size > 1) {ierr = PetscMalloc1(maxCorners+maxCells, &remoteVertices);CHKERRQ(ierr);}
1529852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
153552f7358SJed Brown       MPI_Status status;
154552f7358SJed Brown 
155ffc4695bSBarry Smith       ierr = MPI_Recv(&numCorners, 1, MPIU_INT, proc, tag, comm, &status);CHKERRMPI(ierr);
156ffc4695bSBarry Smith       ierr = MPI_Recv(remoteVertices, numCorners, MPIU_INT, proc, tag, comm, &status);CHKERRMPI(ierr);
157552f7358SJed Brown       for (c = 0; c < numCorners;) {
158552f7358SJed Brown         PetscInt nC = remoteVertices[c++];
159552f7358SJed Brown 
160552f7358SJed Brown         for (v = 0; v < nC; ++v, ++c) {
1610dcf9c70SMatthew G. Knepley           vertices[v] = remoteVertices[c];
1620dcf9c70SMatthew G. Knepley         }
1636b8858a1SBarry Smith         ierr = PetscFPrintf(comm, fp, "%D ", nC);CHKERRQ(ierr);
1640dcf9c70SMatthew G. Knepley         for (v = 0; v < nC; ++v) {
16596ca5757SLisandro Dalcin           ierr = PetscFPrintf(comm, fp, " %D", vertices[v]);CHKERRQ(ierr);
166552f7358SJed Brown         }
167552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
168552f7358SJed Brown       }
169552f7358SJed Brown     }
1709852e123SBarry Smith     if (size > 1) {ierr = PetscFree(remoteVertices);CHKERRQ(ierr);}
1710dcf9c70SMatthew G. Knepley     ierr = PetscFree(vertices);CHKERRQ(ierr);
172552f7358SJed Brown   } else {
173552f7358SJed Brown     PetscInt *localVertices, numSend = numCells+numCorners, k = 0;
174552f7358SJed Brown 
175785e854fSJed Brown     ierr = PetscMalloc1(numSend, &localVertices);CHKERRQ(ierr);
176552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
1770298fd71SBarry Smith       PetscInt *closure = NULL;
17802281ff3SMatthew G. Knepley       PetscInt closureSize, value, nC = 0;
179552f7358SJed Brown 
18002281ff3SMatthew G. Knepley       if (label) {
18102281ff3SMatthew G. Knepley         ierr = DMLabelGetValue(label, c, &value);CHKERRQ(ierr);
182552f7358SJed Brown         if (value != 1) continue;
183552f7358SJed Brown       }
184552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
185552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
186552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
187552f7358SJed Brown           const PetscInt gv = gvertex[closure[v] - vStart];
188552f7358SJed Brown           closure[nC++] = gv < 0 ? -(gv+1) : gv;
189552f7358SJed Brown         }
190552f7358SJed Brown       }
191552f7358SJed Brown       corners[numCells++] = nC;
192552f7358SJed Brown       localVertices[k++]  = nC;
193552f7358SJed Brown       for (v = 0; v < nC; ++v, ++k) {
194552f7358SJed Brown         localVertices[k] = closure[v];
195552f7358SJed Brown       }
196552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
19796ca5757SLisandro Dalcin       ierr = DMPlexReorderCell(dm, c, localVertices+k-nC);CHKERRQ(ierr);
198552f7358SJed Brown     }
199*98921bdaSJacob Faibussowitsch     if (k != numSend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid number of vertices to send %D should be %D", k, numSend);
200ffc4695bSBarry Smith     ierr = MPI_Send(&numSend, 1, MPIU_INT, 0, tag, comm);CHKERRMPI(ierr);
201ffc4695bSBarry Smith     ierr = MPI_Send(localVertices, numSend, MPIU_INT, 0, tag, comm);CHKERRMPI(ierr);
202552f7358SJed Brown     ierr = PetscFree(localVertices);CHKERRQ(ierr);
203552f7358SJed Brown   }
204552f7358SJed Brown   ierr = ISRestoreIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr);
2056b8858a1SBarry Smith   ierr = PetscFPrintf(comm, fp, "CELL_TYPES %D\n", totCells);CHKERRQ(ierr);
206dd400576SPatrick Sanan   if (rank == 0) {
207552f7358SJed Brown     PetscInt cellType;
208552f7358SJed Brown 
209552f7358SJed Brown     for (c = 0; c < numCells; ++c) {
210de0a4605SMatthew G. Knepley       ierr = DMPlexVTKGetCellType_Internal(dm, dim, corners[c], &cellType);CHKERRQ(ierr);
2116b8858a1SBarry Smith       ierr = PetscFPrintf(comm, fp, "%D\n", cellType);CHKERRQ(ierr);
212552f7358SJed Brown     }
2139852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
214552f7358SJed Brown       MPI_Status status;
215552f7358SJed Brown 
216ffc4695bSBarry Smith       ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRMPI(ierr);
217ffc4695bSBarry Smith       ierr = MPI_Recv(corners, numCells, MPIU_INT, proc, tag, comm, &status);CHKERRMPI(ierr);
218552f7358SJed Brown       for (c = 0; c < numCells; ++c) {
219de0a4605SMatthew G. Knepley         ierr = DMPlexVTKGetCellType_Internal(dm, dim, corners[c], &cellType);CHKERRQ(ierr);
2206b8858a1SBarry Smith         ierr = PetscFPrintf(comm, fp, "%D\n", cellType);CHKERRQ(ierr);
221552f7358SJed Brown       }
222552f7358SJed Brown     }
223552f7358SJed Brown   } else {
224ffc4695bSBarry Smith     ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRMPI(ierr);
225ffc4695bSBarry Smith     ierr = MPI_Send(corners, numCells, MPIU_INT, 0, tag, comm);CHKERRMPI(ierr);
226552f7358SJed Brown   }
227552f7358SJed Brown   ierr        = PetscFree(corners);CHKERRQ(ierr);
228552f7358SJed Brown   *totalCells = totCells;
229552f7358SJed Brown   PetscFunctionReturn(0);
230552f7358SJed Brown }
231552f7358SJed Brown 
2327508e736SMatthew G. Knepley static PetscErrorCode DMPlexVTKWritePartition_ASCII(DM dm, FILE *fp)
233e8f6f0f6SMatthew G. Knepley {
234e8f6f0f6SMatthew G. Knepley   MPI_Comm       comm;
235e8f6f0f6SMatthew G. Knepley   PetscInt       numCells = 0, cellHeight;
2362f029394SStefano Zampini   PetscInt       numLabelCells, cStart, cEnd, c;
2379852e123SBarry Smith   PetscMPIInt    size, rank, proc, tag;
238e8f6f0f6SMatthew G. Knepley   PetscBool      hasLabel;
239e8f6f0f6SMatthew G. Knepley   PetscErrorCode ierr;
240e8f6f0f6SMatthew G. Knepley 
241e8f6f0f6SMatthew G. Knepley   PetscFunctionBegin;
242e8f6f0f6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
243e8f6f0f6SMatthew G. Knepley   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
244ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
245ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
246e8f6f0f6SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
247e8f6f0f6SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
248c58f1c22SToby Isaac   ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
249e8f6f0f6SMatthew G. Knepley   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
250e8f6f0f6SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
251e8f6f0f6SMatthew G. Knepley     if (hasLabel) {
252e8f6f0f6SMatthew G. Knepley       PetscInt value;
253e8f6f0f6SMatthew G. Knepley 
254c58f1c22SToby Isaac       ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
255e8f6f0f6SMatthew G. Knepley       if (value != 1) continue;
256e8f6f0f6SMatthew G. Knepley     }
257e8f6f0f6SMatthew G. Knepley     ++numCells;
258e8f6f0f6SMatthew G. Knepley   }
259dd400576SPatrick Sanan   if (rank == 0) {
260e8f6f0f6SMatthew G. Knepley     for (c = 0; c < numCells; ++c) {ierr = PetscFPrintf(comm, fp, "%d\n", rank);CHKERRQ(ierr);}
2619852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
262e8f6f0f6SMatthew G. Knepley       MPI_Status status;
263e8f6f0f6SMatthew G. Knepley 
264ffc4695bSBarry Smith       ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRMPI(ierr);
265e8f6f0f6SMatthew G. Knepley       for (c = 0; c < numCells; ++c) {ierr = PetscFPrintf(comm, fp, "%d\n", proc);CHKERRQ(ierr);}
266e8f6f0f6SMatthew G. Knepley     }
267e8f6f0f6SMatthew G. Knepley   } else {
268ffc4695bSBarry Smith     ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRMPI(ierr);
269e8f6f0f6SMatthew G. Knepley   }
270e8f6f0f6SMatthew G. Knepley   PetscFunctionReturn(0);
271e8f6f0f6SMatthew G. Knepley }
272e8f6f0f6SMatthew G. Knepley 
27376b700caSToby Isaac #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
27476b700caSToby Isaac typedef double PetscVTKReal;
27576b700caSToby Isaac #elif defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL___FP16)
27676b700caSToby Isaac typedef float PetscVTKReal;
27776b700caSToby Isaac #else
27876b700caSToby Isaac typedef PetscReal PetscVTKReal;
27976b700caSToby Isaac #endif
28076b700caSToby Isaac 
28176b700caSToby Isaac static PetscErrorCode DMPlexVTKWriteSection_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec v, FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale, PetscInt imag)
2820adebc6cSBarry Smith {
28382f516ccSBarry Smith   MPI_Comm           comm;
284552f7358SJed Brown   const MPI_Datatype mpiType = MPIU_SCALAR;
285552f7358SJed Brown   PetscScalar        *array;
286552f7358SJed Brown   PetscInt           numDof = 0, maxDof;
287412e9a14SMatthew G. Knepley   PetscInt           numLabelCells, cellHeight, cStart, cEnd, numLabelVertices, vStart, vEnd, pStart, pEnd, p;
2889852e123SBarry Smith   PetscMPIInt        size, rank, proc, tag;
289552f7358SJed Brown   PetscBool          hasLabel;
290552f7358SJed Brown   PetscErrorCode     ierr;
291552f7358SJed Brown 
292552f7358SJed Brown   PetscFunctionBegin;
29382f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
294552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
295552f7358SJed Brown   PetscValidHeaderSpecific(v,VEC_CLASSID,4);
296552f7358SJed Brown   if (precision < 0) precision = 6;
297552f7358SJed Brown   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
298ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
299ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
300552f7358SJed Brown   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
301552f7358SJed Brown   /* VTK only wants the values at cells or vertices */
302552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
303552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
304552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
305552f7358SJed Brown   pStart   = PetscMax(PetscMin(cStart, vStart), pStart);
306552f7358SJed Brown   pEnd     = PetscMin(PetscMax(cEnd,   vEnd),   pEnd);
307c58f1c22SToby Isaac   ierr     = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
308c58f1c22SToby Isaac   ierr     = DMGetStratumSize(dm, "vtk", 2, &numLabelVertices);CHKERRQ(ierr);
309552f7358SJed Brown   hasLabel = numLabelCells > 0 || numLabelVertices > 0 ? PETSC_TRUE : PETSC_FALSE;
310552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
311552f7358SJed Brown     /* Reject points not either cells or vertices */
312552f7358SJed Brown     if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
313552f7358SJed Brown     if (hasLabel) {
314552f7358SJed Brown       PetscInt value;
315552f7358SJed Brown 
316552f7358SJed Brown       if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
317552f7358SJed Brown           ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
318c58f1c22SToby Isaac         ierr = DMGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
319552f7358SJed Brown         if (value != 1) continue;
320552f7358SJed Brown       }
321552f7358SJed Brown     }
322552f7358SJed Brown     ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr);
3238865f1eaSKarl Rupp     if (numDof) break;
324552f7358SJed Brown   }
325820f2d46SBarry Smith   ierr = MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, comm);CHKERRMPI(ierr);
326552f7358SJed Brown   enforceDof = PetscMax(enforceDof, maxDof);
327552f7358SJed Brown   ierr = VecGetArray(v, &array);CHKERRQ(ierr);
328dd400576SPatrick Sanan   if (rank == 0) {
32976b700caSToby Isaac     PetscVTKReal dval;
33076b700caSToby Isaac     PetscScalar  val;
331552f7358SJed Brown     char formatString[8];
332552f7358SJed Brown 
333552f7358SJed Brown     ierr = PetscSNPrintf(formatString, 8, "%%.%de", precision);CHKERRQ(ierr);
334552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
335552f7358SJed Brown       /* Here we lose a way to filter points by keeping them out of the Numbering */
336552f7358SJed Brown       PetscInt dof, off, goff, d;
337552f7358SJed Brown 
338552f7358SJed Brown       /* Reject points not either cells or vertices */
339552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
340552f7358SJed Brown       if (hasLabel) {
341552f7358SJed Brown         PetscInt value;
342552f7358SJed Brown 
343552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
344552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
345c58f1c22SToby Isaac           ierr = DMGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
346552f7358SJed Brown           if (value != 1) continue;
347552f7358SJed Brown         }
348552f7358SJed Brown       }
349552f7358SJed Brown       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
350552f7358SJed Brown       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
351552f7358SJed Brown       ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
352552f7358SJed Brown       if (dof && goff >= 0) {
353552f7358SJed Brown         for (d = 0; d < dof; d++) {
354552f7358SJed Brown           if (d > 0) {
355552f7358SJed Brown             ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr);
356552f7358SJed Brown           }
35776b700caSToby Isaac           val = array[off+d];
35876b700caSToby Isaac           dval = (PetscVTKReal) ((imag ? PetscImaginaryPart(val) : PetscRealPart(val)) * scale);
3592dea1156SStefano Zampini           ierr = PetscFPrintf(comm, fp, formatString, dval);CHKERRQ(ierr);
360552f7358SJed Brown         }
361552f7358SJed Brown         for (d = dof; d < enforceDof; d++) {
362552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr);
363552f7358SJed Brown         }
364552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
365552f7358SJed Brown       }
366552f7358SJed Brown     }
3679852e123SBarry Smith     for (proc = 1; proc < size; ++proc) {
368552f7358SJed Brown       PetscScalar *remoteValues;
369d892089bSMatthew G. Knepley       PetscInt    size = 0, d;
370552f7358SJed Brown       MPI_Status  status;
371552f7358SJed Brown 
372ffc4695bSBarry Smith       ierr = MPI_Recv(&size, 1, MPIU_INT, proc, tag, comm, &status);CHKERRMPI(ierr);
373785e854fSJed Brown       ierr = PetscMalloc1(size, &remoteValues);CHKERRQ(ierr);
374ffc4695bSBarry Smith       ierr = MPI_Recv(remoteValues, size, mpiType, proc, tag, comm, &status);CHKERRMPI(ierr);
375552f7358SJed Brown       for (p = 0; p < size/maxDof; ++p) {
376552f7358SJed Brown         for (d = 0; d < maxDof; ++d) {
377552f7358SJed Brown           if (d > 0) {
378552f7358SJed Brown             ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr);
379552f7358SJed Brown           }
38076b700caSToby Isaac           val = remoteValues[p*maxDof+d];
38176b700caSToby Isaac           dval = (PetscVTKReal) ((imag ? PetscImaginaryPart(val) : PetscRealPart(val)) * scale);
3822dea1156SStefano Zampini           ierr = PetscFPrintf(comm, fp, formatString, dval);CHKERRQ(ierr);
383552f7358SJed Brown         }
384552f7358SJed Brown         for (d = maxDof; d < enforceDof; ++d) {
385552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr);
386552f7358SJed Brown         }
387552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
388552f7358SJed Brown       }
389552f7358SJed Brown       ierr = PetscFree(remoteValues);CHKERRQ(ierr);
390552f7358SJed Brown     }
391552f7358SJed Brown   } else {
392552f7358SJed Brown     PetscScalar *localValues;
393552f7358SJed Brown     PetscInt    size, k = 0;
394552f7358SJed Brown 
395552f7358SJed Brown     ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
396785e854fSJed Brown     ierr = PetscMalloc1(size, &localValues);CHKERRQ(ierr);
397552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
398552f7358SJed Brown       PetscInt dof, off, goff, d;
399552f7358SJed Brown 
400552f7358SJed Brown       /* Reject points not either cells or vertices */
401552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
402552f7358SJed Brown       if (hasLabel) {
403552f7358SJed Brown         PetscInt value;
404552f7358SJed Brown 
405552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
406552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
407c58f1c22SToby Isaac           ierr = DMGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
408552f7358SJed Brown           if (value != 1) continue;
409552f7358SJed Brown         }
410552f7358SJed Brown       }
411552f7358SJed Brown       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
412552f7358SJed Brown       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
413552f7358SJed Brown       ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
414552f7358SJed Brown       if (goff >= 0) {
415552f7358SJed Brown         for (d = 0; d < dof; ++d) {
416552f7358SJed Brown           localValues[k++] = array[off+d];
417552f7358SJed Brown         }
418552f7358SJed Brown       }
419552f7358SJed Brown     }
420ffc4695bSBarry Smith     ierr = MPI_Send(&k, 1, MPIU_INT, 0, tag, comm);CHKERRMPI(ierr);
421ffc4695bSBarry Smith     ierr = MPI_Send(localValues, k, mpiType, 0, tag, comm);CHKERRMPI(ierr);
422552f7358SJed Brown     ierr = PetscFree(localValues);CHKERRQ(ierr);
423552f7358SJed Brown   }
424552f7358SJed Brown   ierr = VecRestoreArray(v, &array);CHKERRQ(ierr);
425552f7358SJed Brown   PetscFunctionReturn(0);
426552f7358SJed Brown }
427552f7358SJed Brown 
42876b700caSToby 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)
429552f7358SJed Brown {
43082f516ccSBarry Smith   MPI_Comm       comm;
431552f7358SJed Brown   PetscInt       numDof = 0, maxDof;
432552f7358SJed Brown   PetscInt       pStart, pEnd, p;
433552f7358SJed Brown   PetscErrorCode ierr;
434552f7358SJed Brown 
435552f7358SJed Brown   PetscFunctionBegin;
43682f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
437552f7358SJed Brown   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
438552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
439552f7358SJed Brown     ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr);
4408865f1eaSKarl Rupp     if (numDof) break;
441552f7358SJed Brown   }
442552f7358SJed Brown   numDof = PetscMax(numDof, enforceDof);
443820f2d46SBarry Smith   ierr = MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm));CHKERRMPI(ierr);
444552f7358SJed Brown   if (!name) name = "Unknown";
445552f7358SJed Brown   if (maxDof == 3) {
44676b700caSToby Isaac     if (nameComplex) {
44776b700caSToby Isaac       ierr = PetscFPrintf(comm, fp, "VECTORS %s.%s double\n", name, imag ? "Im" : "Re");CHKERRQ(ierr);
44876b700caSToby Isaac     } else {
449552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "VECTORS %s double\n", name);CHKERRQ(ierr);
45076b700caSToby Isaac     }
45176b700caSToby Isaac   } else {
45276b700caSToby Isaac     if (nameComplex) {
45376b700caSToby Isaac       ierr = PetscFPrintf(comm, fp, "SCALARS %s.%s double %D\n", name, imag ? "Im" : "Re", maxDof);CHKERRQ(ierr);
454552f7358SJed Brown     } else {
4556b8858a1SBarry Smith       ierr = PetscFPrintf(comm, fp, "SCALARS %s double %D\n", name, maxDof);CHKERRQ(ierr);
45676b700caSToby Isaac     }
457552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr);
458552f7358SJed Brown   }
45976b700caSToby Isaac   ierr = DMPlexVTKWriteSection_ASCII(dm, section, globalSection, field, fp, enforceDof, precision, scale, imag);CHKERRQ(ierr);
460552f7358SJed Brown   PetscFunctionReturn(0);
461552f7358SJed Brown }
462552f7358SJed Brown 
463552f7358SJed Brown static PetscErrorCode DMPlexVTKWriteAll_ASCII(DM dm, PetscViewer viewer)
464552f7358SJed Brown {
46582f516ccSBarry Smith   MPI_Comm                 comm;
466552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*) viewer->data;
467552f7358SJed Brown   FILE                     *fp;
468552f7358SJed Brown   PetscViewerVTKObjectLink link;
469552f7358SJed Brown   PetscSection             coordSection, globalCoordSection;
470552f7358SJed Brown   PetscLayout              vLayout;
471552f7358SJed Brown   Vec                      coordinates;
472552f7358SJed Brown   PetscReal                lengthScale;
473412e9a14SMatthew G. Knepley   PetscInt                 totVertices, totCells = 0, loops_per_scalar, l;
47476b700caSToby Isaac   PetscBool                hasPoint = PETSC_FALSE, hasCell = PETSC_FALSE, writePartition = PETSC_FALSE, localized, writeComplex;
47552b05773SJed Brown   const char               *dmname;
476552f7358SJed Brown   PetscErrorCode           ierr;
477552f7358SJed Brown 
478552f7358SJed Brown   PetscFunctionBegin;
47976b700caSToby Isaac #if defined(PETSC_USE_COMPLEX)
48076b700caSToby Isaac   loops_per_scalar = 2;
48176b700caSToby Isaac   writeComplex = PETSC_TRUE;
48276b700caSToby Isaac #else
48376b700caSToby Isaac   loops_per_scalar = 1;
48476b700caSToby Isaac   writeComplex = PETSC_FALSE;
48576b700caSToby Isaac #endif
486d72c04d6SStefano Zampini   ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr);
487d72c04d6SStefano Zampini   if (localized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"VTK output with localized coordinates not yet supported");
48882f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
489552f7358SJed Brown   ierr = PetscFOpen(comm, vtk->filename, "wb", &fp);CHKERRQ(ierr);
49052b05773SJed Brown   ierr = PetscObjectGetName((PetscObject)dm, &dmname);CHKERRQ(ierr);
491552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "# vtk DataFile Version 2.0\n");CHKERRQ(ierr);
49252b05773SJed Brown   ierr = PetscFPrintf(comm, fp, "%s\n", dmname);CHKERRQ(ierr);
493552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "ASCII\n");CHKERRQ(ierr);
494552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "DATASET UNSTRUCTURED_GRID\n");CHKERRQ(ierr);
495552f7358SJed Brown   /* Vertices */
496552f7358SJed Brown   ierr = DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);CHKERRQ(ierr);
49769d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
49815b58121SMatthew G. Knepley   ierr = PetscSectionCreateGlobalSection(coordSection, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalCoordSection);CHKERRQ(ierr);
499552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
50082f516ccSBarry Smith   ierr = PetscSectionGetPointLayout(PetscObjectComm((PetscObject)dm), globalCoordSection, &vLayout);CHKERRQ(ierr);
501552f7358SJed Brown   ierr = PetscLayoutGetSize(vLayout, &totVertices);CHKERRQ(ierr);
5026b8858a1SBarry Smith   ierr = PetscFPrintf(comm, fp, "POINTS %D double\n", totVertices);CHKERRQ(ierr);
50376b700caSToby Isaac   ierr = DMPlexVTKWriteSection_ASCII(dm, coordSection, globalCoordSection, coordinates, fp, 3, PETSC_DETERMINE, lengthScale, 0);CHKERRQ(ierr);
504552f7358SJed Brown   /* Cells */
505552f7358SJed Brown   ierr = DMPlexVTKWriteCells_ASCII(dm, fp, &totCells);CHKERRQ(ierr);
506552f7358SJed Brown   /* Vertex fields */
507552f7358SJed Brown   for (link = vtk->link; link; link = link->next) {
508552f7358SJed Brown     if ((link->ft == PETSC_VTK_POINT_FIELD) || (link->ft == PETSC_VTK_POINT_VECTOR_FIELD)) hasPoint = PETSC_TRUE;
509552f7358SJed Brown     if ((link->ft == PETSC_VTK_CELL_FIELD)  || (link->ft == PETSC_VTK_CELL_VECTOR_FIELD))  hasCell  = PETSC_TRUE;
510552f7358SJed Brown   }
511552f7358SJed Brown   if (hasPoint) {
5126b8858a1SBarry Smith     ierr = PetscFPrintf(comm, fp, "POINT_DATA %D\n", totVertices);CHKERRQ(ierr);
513552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
514552f7358SJed Brown       Vec          X = (Vec) link->vec;
515e630c359SToby Isaac       PetscSection section = NULL, globalSection, newSection = NULL;
516e630c359SToby Isaac       char         namebuf[256];
517552f7358SJed Brown       const char   *name;
518552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
519552f7358SJed Brown 
520552f7358SJed Brown       if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
521552f7358SJed Brown       if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) enforceDof = 3;
522552f7358SJed Brown       ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr);
523e630c359SToby Isaac       ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
524e630c359SToby Isaac       if (!section) {
525e630c359SToby Isaac         DM           dmX;
526e630c359SToby Isaac 
527552f7358SJed Brown         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
528552f7358SJed Brown         if (dmX) {
529088580cfSMatthew G Knepley           DMLabel  subpointMap, subpointMapX;
530839dd189SMatthew G Knepley           PetscInt dim, dimX, pStart, pEnd, qStart, qEnd;
531839dd189SMatthew G Knepley 
53292fd8e1eSJed Brown           ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);
533839dd189SMatthew G Knepley           /* Here is where we check whether dmX is a submesh of dm */
534c73cfb54SMatthew G. Knepley           ierr = DMGetDimension(dm,  &dim);CHKERRQ(ierr);
535c73cfb54SMatthew G. Knepley           ierr = DMGetDimension(dmX, &dimX);CHKERRQ(ierr);
536839dd189SMatthew G Knepley           ierr = DMPlexGetChart(dm,  &pStart, &pEnd);CHKERRQ(ierr);
537839dd189SMatthew G Knepley           ierr = DMPlexGetChart(dmX, &qStart, &qEnd);CHKERRQ(ierr);
538839dd189SMatthew G Knepley           ierr = DMPlexGetSubpointMap(dm,  &subpointMap);CHKERRQ(ierr);
539839dd189SMatthew G Knepley           ierr = DMPlexGetSubpointMap(dmX, &subpointMapX);CHKERRQ(ierr);
540839dd189SMatthew G Knepley           if (((dim != dimX) || ((pEnd-pStart) < (qEnd-qStart))) && subpointMap && !subpointMapX) {
541434e42b5SMatthew G. Knepley             const PetscInt *ind = NULL;
542e6ccafaeSMatthew G Knepley             IS              subpointIS;
543434e42b5SMatthew G. Knepley             PetscInt        n = 0, q;
544839dd189SMatthew G Knepley 
545839dd189SMatthew G Knepley             ierr = PetscSectionGetChart(section, &qStart, &qEnd);CHKERRQ(ierr);
54697d8846cSMatthew Knepley             ierr = DMPlexGetSubpointIS(dm, &subpointIS);CHKERRQ(ierr);
547434e42b5SMatthew G. Knepley             if (subpointIS) {
548e6ccafaeSMatthew G Knepley               ierr = ISGetLocalSize(subpointIS, &n);CHKERRQ(ierr);
549e6ccafaeSMatthew G Knepley               ierr = ISGetIndices(subpointIS, &ind);CHKERRQ(ierr);
550434e42b5SMatthew G. Knepley             }
551839dd189SMatthew G Knepley             ierr = PetscSectionCreate(comm, &newSection);CHKERRQ(ierr);
552839dd189SMatthew G Knepley             ierr = PetscSectionSetChart(newSection, pStart, pEnd);CHKERRQ(ierr);
553839dd189SMatthew G Knepley             for (q = qStart; q < qEnd; ++q) {
554839dd189SMatthew G Knepley               PetscInt dof, off, p;
555839dd189SMatthew G Knepley 
556839dd189SMatthew G Knepley               ierr = PetscSectionGetDof(section, q, &dof);CHKERRQ(ierr);
557839dd189SMatthew G Knepley               if (dof) {
558839dd189SMatthew G Knepley                 ierr = PetscFindInt(q, n, ind, &p);CHKERRQ(ierr);
559839dd189SMatthew G Knepley                 if (p >= pStart) {
560839dd189SMatthew G Knepley                   ierr = PetscSectionSetDof(newSection, p, dof);CHKERRQ(ierr);
561839dd189SMatthew G Knepley                   ierr = PetscSectionGetOffset(section, q, &off);CHKERRQ(ierr);
562839dd189SMatthew G Knepley                   ierr = PetscSectionSetOffset(newSection, p, off);CHKERRQ(ierr);
563839dd189SMatthew G Knepley                 }
564839dd189SMatthew G Knepley               }
565839dd189SMatthew G Knepley             }
566434e42b5SMatthew G. Knepley             if (subpointIS) {
567e6ccafaeSMatthew G Knepley               ierr = ISRestoreIndices(subpointIS, &ind);CHKERRQ(ierr);
568434e42b5SMatthew G. Knepley             }
569839dd189SMatthew G Knepley             /* No need to setup section */
570839dd189SMatthew G Knepley             section = newSection;
571839dd189SMatthew G Knepley           }
572552f7358SJed Brown         }
573e630c359SToby Isaac       }
574*98921bdaSJacob Faibussowitsch       if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it and could not create one from VecGetDM()", name);
575e630c359SToby Isaac       if (link->field >= 0) {
576e630c359SToby Isaac         const char *fieldname;
577e630c359SToby Isaac 
578e630c359SToby Isaac         ierr = PetscSectionGetFieldName(section, link->field, &fieldname);CHKERRQ(ierr);
579e630c359SToby Isaac         ierr = PetscSectionGetField(section, link->field, &section);CHKERRQ(ierr);
580e630c359SToby Isaac         if (fieldname) {
5812fb742c9SToby Isaac           ierr = PetscSNPrintf(namebuf, sizeof(namebuf), "%s%s", name, fieldname);CHKERRQ(ierr);
582e630c359SToby Isaac         } else {
5832fb742c9SToby Isaac           ierr = PetscSNPrintf(namebuf, sizeof(namebuf), "%s%D", name, link->field);CHKERRQ(ierr);
584e630c359SToby Isaac         }
5852fb742c9SToby Isaac       } else {
5862fb742c9SToby Isaac         ierr = PetscSNPrintf(namebuf, sizeof(namebuf), "%s", name);CHKERRQ(ierr);
587e630c359SToby Isaac       }
5882fb742c9SToby Isaac       ierr = PetscViewerVTKSanitizeName_Internal(namebuf, sizeof(namebuf));CHKERRQ(ierr);
58915b58121SMatthew G. Knepley       ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection);CHKERRQ(ierr);
59076b700caSToby Isaac       for (l = 0; l < loops_per_scalar; l++) {
5912fb742c9SToby Isaac         ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, namebuf, fp, enforceDof, PETSC_DETERMINE, 1.0, writeComplex, l);CHKERRQ(ierr);
59276b700caSToby Isaac       }
593552f7358SJed Brown       ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr);
594839dd189SMatthew G Knepley       if (newSection) {ierr = PetscSectionDestroy(&newSection);CHKERRQ(ierr);}
595552f7358SJed Brown     }
596552f7358SJed Brown   }
597552f7358SJed Brown   /* Cell Fields */
598c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_view_partition", &writePartition, NULL);CHKERRQ(ierr);
599e8f6f0f6SMatthew G. Knepley   if (hasCell || writePartition) {
6006b8858a1SBarry Smith     ierr = PetscFPrintf(comm, fp, "CELL_DATA %D\n", totCells);CHKERRQ(ierr);
601552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
602552f7358SJed Brown       Vec          X = (Vec) link->vec;
603e630c359SToby Isaac       PetscSection section = NULL, globalSection;
6042fb742c9SToby Isaac       const char   *name = "";
605e630c359SToby Isaac       char         namebuf[256];
606552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
607552f7358SJed Brown 
608552f7358SJed Brown       if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
609552f7358SJed Brown       if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) enforceDof = 3;
610552f7358SJed Brown       ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr);
611e630c359SToby Isaac       ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
612e630c359SToby Isaac       if (!section) {
613e630c359SToby Isaac         DM           dmX;
614e630c359SToby Isaac 
615552f7358SJed Brown         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
616552f7358SJed Brown         if (dmX) {
61792fd8e1eSJed Brown           ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);
618552f7358SJed Brown         }
619e630c359SToby Isaac       }
620*98921bdaSJacob Faibussowitsch       if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it and could not create one from VecGetDM()", name);
621e630c359SToby Isaac       if (link->field >= 0) {
622e630c359SToby Isaac         const char *fieldname;
623e630c359SToby Isaac 
624e630c359SToby Isaac         ierr = PetscSectionGetFieldName(section, link->field, &fieldname);CHKERRQ(ierr);
625e630c359SToby Isaac         ierr = PetscSectionGetField(section, link->field, &section);CHKERRQ(ierr);
626e630c359SToby Isaac         if (fieldname) {
6272fb742c9SToby Isaac           ierr = PetscSNPrintf(namebuf, sizeof(namebuf), "%s%s", name, fieldname);CHKERRQ(ierr);
628e630c359SToby Isaac         } else {
6292fb742c9SToby Isaac           ierr = PetscSNPrintf(namebuf, sizeof(namebuf), "%s%D", name, link->field);CHKERRQ(ierr);
630e630c359SToby Isaac         }
6312fb742c9SToby Isaac       } else {
6322fb742c9SToby Isaac         ierr = PetscSNPrintf(namebuf, sizeof(namebuf), "%s", name);CHKERRQ(ierr);
633e630c359SToby Isaac       }
6342fb742c9SToby Isaac       ierr = PetscViewerVTKSanitizeName_Internal(namebuf, sizeof(namebuf));CHKERRQ(ierr);
63515b58121SMatthew G. Knepley       ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection);CHKERRQ(ierr);
63676b700caSToby Isaac       for (l = 0; l < loops_per_scalar; l++) {
6372fb742c9SToby Isaac         ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, namebuf, fp, enforceDof, PETSC_DETERMINE, 1.0, writeComplex, l);CHKERRQ(ierr);
63876b700caSToby Isaac       }
639552f7358SJed Brown       ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr);
640552f7358SJed Brown     }
641e8f6f0f6SMatthew G. Knepley     if (writePartition) {
642e8f6f0f6SMatthew G. Knepley       ierr = PetscFPrintf(comm, fp, "SCALARS partition int 1\n");CHKERRQ(ierr);
643e8f6f0f6SMatthew G. Knepley       ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr);
644e8f6f0f6SMatthew G. Knepley       ierr = DMPlexVTKWritePartition_ASCII(dm, fp);CHKERRQ(ierr);
645e8f6f0f6SMatthew G. Knepley     }
646552f7358SJed Brown   }
647552f7358SJed Brown   /* Cleanup */
648552f7358SJed Brown   ierr = PetscSectionDestroy(&globalCoordSection);CHKERRQ(ierr);
649552f7358SJed Brown   ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr);
650552f7358SJed Brown   ierr = PetscFClose(comm, fp);CHKERRQ(ierr);
651552f7358SJed Brown   PetscFunctionReturn(0);
652552f7358SJed Brown }
653552f7358SJed Brown 
654552f7358SJed Brown /*@C
655552f7358SJed Brown   DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
656552f7358SJed Brown 
657552f7358SJed Brown   Collective
658552f7358SJed Brown 
6594165533cSJose E. Roman   Input Parameters:
660552f7358SJed Brown + odm - The DMPlex specifying the mesh, passed as a PetscObject
661552f7358SJed Brown - viewer - viewer of type VTK
662552f7358SJed Brown 
663552f7358SJed Brown   Level: developer
664552f7358SJed Brown 
665552f7358SJed Brown   Note:
666552f7358SJed Brown   This function is a callback used by the VTK viewer to actually write the file.
667552f7358SJed 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.
668552f7358SJed Brown   Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
669552f7358SJed Brown 
670552f7358SJed Brown .seealso: PETSCVIEWERVTK
671552f7358SJed Brown @*/
672552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer)
673552f7358SJed Brown {
674552f7358SJed Brown   DM             dm = (DM) odm;
675552f7358SJed Brown   PetscBool      isvtk;
676552f7358SJed Brown   PetscErrorCode ierr;
677552f7358SJed Brown 
678552f7358SJed Brown   PetscFunctionBegin;
679552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
680552f7358SJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
681552f7358SJed Brown   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr);
682*98921bdaSJacob Faibussowitsch   if (!isvtk) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name);
683552f7358SJed Brown   switch (viewer->format) {
6848ec8862eSJed Brown   case PETSC_VIEWER_ASCII_VTK_DEPRECATED:
685552f7358SJed Brown     ierr = DMPlexVTKWriteAll_ASCII(dm, viewer);CHKERRQ(ierr);
686552f7358SJed Brown     break;
687552f7358SJed Brown   case PETSC_VIEWER_VTK_VTU:
688552f7358SJed Brown     ierr = DMPlexVTKWriteAll_VTU(dm, viewer);CHKERRQ(ierr);
689552f7358SJed Brown     break;
690*98921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]);
691552f7358SJed Brown   }
692552f7358SJed Brown   PetscFunctionReturn(0);
693552f7358SJed Brown }
694