xref: /petsc/src/dm/impls/plex/plexvtk.c (revision 770b213bc1b675a346d8e3c224f34bd2ee8b3f57)
1552f7358SJed Brown #define PETSCDM_DLL
2552f7358SJed Brown #include <petsc-private/pleximpl.h>    /*I   "petscdmplex.h"   I*/
3552f7358SJed Brown #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h>
4552f7358SJed Brown 
5552f7358SJed Brown #undef __FUNCT__
6552f7358SJed Brown #define __FUNCT__ "DMPlexVTKGetCellType"
7552f7358SJed Brown PetscErrorCode DMPlexVTKGetCellType(DM dm, PetscInt dim, PetscInt corners, PetscInt *cellType) {
8552f7358SJed Brown   PetscFunctionBegin;
9552f7358SJed Brown   *cellType = -1;
10552f7358SJed Brown   switch(dim) {
11552f7358SJed Brown   case 0:
12552f7358SJed Brown     switch(corners) {
13552f7358SJed Brown     case 1:
14552f7358SJed Brown       *cellType = 1; /* VTK_VERTEX */
15552f7358SJed Brown       break;
16552f7358SJed Brown     default:
17552f7358SJed Brown       break;
18552f7358SJed Brown     }
19552f7358SJed Brown     break;
20552f7358SJed Brown   case 1:
21552f7358SJed Brown     switch(corners) {
22552f7358SJed Brown     case 2:
23552f7358SJed Brown       *cellType = 3; /* VTK_LINE */
24552f7358SJed Brown       break;
25552f7358SJed Brown     case 3:
26552f7358SJed Brown       *cellType = 21; /* VTK_QUADRATIC_EDGE */
27552f7358SJed Brown       break;
28552f7358SJed Brown     default:
29552f7358SJed Brown       break;
30552f7358SJed Brown     }
31552f7358SJed Brown     break;
32552f7358SJed Brown   case 2:
33552f7358SJed Brown     switch(corners) {
34552f7358SJed Brown     case 3:
35552f7358SJed Brown       *cellType = 5; /* VTK_TRIANGLE */
36552f7358SJed Brown       break;
37552f7358SJed Brown     case 4:
38552f7358SJed Brown       *cellType = 9; /* VTK_QUAD */
39552f7358SJed Brown       break;
40552f7358SJed Brown     case 6:
41552f7358SJed Brown       *cellType = 22; /* VTK_QUADRATIC_TRIANGLE */
42552f7358SJed Brown       break;
43552f7358SJed Brown     case 9:
44552f7358SJed Brown       *cellType = 23; /* VTK_QUADRATIC_QUAD */
45552f7358SJed Brown       break;
46552f7358SJed Brown     default:
47552f7358SJed Brown       break;
48552f7358SJed Brown     }
49552f7358SJed Brown     break;
50552f7358SJed Brown   case 3:
51552f7358SJed Brown     switch(corners) {
52552f7358SJed Brown     case 4:
53552f7358SJed Brown       *cellType = 10; /* VTK_TETRA */
54552f7358SJed Brown       break;
55552f7358SJed Brown     case 8:
56552f7358SJed Brown       *cellType = 12; /* VTK_HEXAHEDRON */
57552f7358SJed Brown       break;
58552f7358SJed Brown     case 10:
59552f7358SJed Brown       *cellType = 24; /* VTK_QUADRATIC_TETRA */
60552f7358SJed Brown       break;
61552f7358SJed Brown     case 27:
62552f7358SJed Brown       *cellType = 29; /* VTK_QUADRATIC_HEXAHEDRON */
63552f7358SJed Brown       break;
64552f7358SJed Brown     default:
65552f7358SJed Brown       break;
66552f7358SJed Brown     }
67552f7358SJed Brown   }
68552f7358SJed Brown   PetscFunctionReturn(0);
69552f7358SJed Brown }
70552f7358SJed Brown 
71552f7358SJed Brown #undef __FUNCT__
72552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteCells_ASCII"
73552f7358SJed Brown PetscErrorCode DMPlexVTKWriteCells_ASCII(DM dm, FILE *fp, PetscInt *totalCells)
74552f7358SJed Brown {
75552f7358SJed Brown   MPI_Comm        comm = ((PetscObject) dm)->comm;
76552f7358SJed Brown   IS              globalVertexNumbers;
77552f7358SJed Brown   const PetscInt *gvertex;
78552f7358SJed Brown   PetscInt        dim;
79552f7358SJed Brown   PetscInt        numCorners = 0, totCorners = 0, maxCorners, *corners;
80552f7358SJed Brown   PetscInt        numCells   = 0, totCells   = 0, maxCells, cellHeight;
81552f7358SJed Brown   PetscInt        numLabelCells, cMax, cStart, cEnd, c, vStart, vEnd, v;
82552f7358SJed Brown   PetscMPIInt     numProcs, rank, proc, tag;
83552f7358SJed Brown   PetscBool       hasLabel;
84552f7358SJed Brown   PetscErrorCode  ierr;
85552f7358SJed Brown 
86552f7358SJed Brown   PetscFunctionBegin;
87552f7358SJed Brown   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
88552f7358SJed Brown   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
89552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
90552f7358SJed Brown   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
91552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
92552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
93552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
94*770b213bSMatthew G Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
95552f7358SJed Brown   if (cMax >= 0) {cEnd = PetscMin(cEnd, cMax);}
96552f7358SJed Brown   ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
97552f7358SJed Brown   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
98552f7358SJed Brown   for(c = cStart; c < cEnd; ++c) {
99552f7358SJed Brown     PetscInt *closure = PETSC_NULL;
100552f7358SJed Brown     PetscInt  closureSize;
101552f7358SJed Brown 
102552f7358SJed Brown     if (hasLabel) {
103552f7358SJed Brown       PetscInt value;
104552f7358SJed Brown 
105552f7358SJed Brown       ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
106552f7358SJed Brown       if (value != 1) continue;
107552f7358SJed Brown     }
108552f7358SJed Brown     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
109552f7358SJed Brown     for (v = 0; v < closureSize*2; v += 2) {
110552f7358SJed Brown       if ((closure[v] >= vStart) && (closure[v] < vEnd)) ++numCorners;
111552f7358SJed Brown     }
112552f7358SJed Brown     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
113552f7358SJed Brown     ++numCells;
114552f7358SJed Brown   }
115552f7358SJed Brown   maxCells = numCells;
116552f7358SJed Brown   ierr = MPI_Reduce(&numCells, &totCells, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRQ(ierr);
117552f7358SJed Brown   ierr = MPI_Reduce(&numCells, &maxCells, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr);
118552f7358SJed Brown   ierr = MPI_Reduce(&numCorners, &totCorners, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRQ(ierr);
119552f7358SJed Brown   ierr = MPI_Reduce(&numCorners, &maxCorners, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr);
120552f7358SJed Brown   ierr = DMPlexGetVertexNumbering(dm, &globalVertexNumbers);CHKERRQ(ierr);
121552f7358SJed Brown   ierr = ISGetIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr);
122552f7358SJed Brown   ierr = PetscMalloc(maxCells * sizeof(PetscInt), &corners);CHKERRQ(ierr);
123552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "CELLS %d %d\n", totCells, totCorners+totCells);CHKERRQ(ierr);
124552f7358SJed Brown   if (!rank) {
125552f7358SJed Brown     PetscInt *remoteVertices;
126552f7358SJed Brown 
127552f7358SJed Brown     for(c = cStart, numCells = 0; c < cEnd; ++c) {
128552f7358SJed Brown       PetscInt *closure = PETSC_NULL;
129552f7358SJed Brown       PetscInt  closureSize, nC = 0;
130552f7358SJed Brown 
131552f7358SJed Brown       if (hasLabel) {
132552f7358SJed Brown         PetscInt value;
133552f7358SJed Brown 
134552f7358SJed Brown         ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
135552f7358SJed Brown         if (value != 1) continue;
136552f7358SJed Brown       }
137552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
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];
141552f7358SJed Brown           closure[nC++] = gv < 0 ? -(gv+1) : gv;
142552f7358SJed Brown         }
143552f7358SJed Brown       }
144552f7358SJed Brown       corners[numCells++] = nC;
145552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "%d ", nC);CHKERRQ(ierr);
146552f7358SJed Brown       for(v = 0; v < nC; ++v) {
147552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, " %d", closure[v]);CHKERRQ(ierr);
148552f7358SJed Brown       }
149552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
150552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
151552f7358SJed Brown     }
152552f7358SJed Brown     if (numProcs > 1) {ierr = PetscMalloc((maxCorners+maxCells) * sizeof(PetscInt), &remoteVertices);CHKERRQ(ierr);}
153552f7358SJed Brown     for(proc = 1; proc < numProcs; ++proc) {
154552f7358SJed Brown       MPI_Status status;
155552f7358SJed Brown 
156552f7358SJed Brown       ierr = MPI_Recv(&numCorners, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
157552f7358SJed Brown       ierr = MPI_Recv(remoteVertices, numCorners, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
158552f7358SJed Brown       for(c = 0; c < numCorners;) {
159552f7358SJed Brown         PetscInt nC = remoteVertices[c++];
160552f7358SJed Brown 
161552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "%d ", nC);CHKERRQ(ierr);
162552f7358SJed Brown         for(v = 0; v < nC; ++v, ++c) {
163552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " %d", remoteVertices[c]);CHKERRQ(ierr);
164552f7358SJed Brown         }
165552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
166552f7358SJed Brown       }
167552f7358SJed Brown     }
168552f7358SJed Brown     if (numProcs > 1) {ierr = PetscFree(remoteVertices);CHKERRQ(ierr);}
169552f7358SJed Brown   } else {
170552f7358SJed Brown     PetscInt *localVertices, numSend = numCells+numCorners, k = 0;
171552f7358SJed Brown 
172552f7358SJed Brown     ierr = PetscMalloc(numSend * sizeof(PetscInt), &localVertices);CHKERRQ(ierr);
173552f7358SJed Brown     for (c = cStart, numCells = 0; c < cEnd; ++c) {
174552f7358SJed Brown       PetscInt *closure = PETSC_NULL;
175552f7358SJed Brown       PetscInt  closureSize, nC = 0;
176552f7358SJed Brown 
177552f7358SJed Brown       if (hasLabel) {
178552f7358SJed Brown         PetscInt value;
179552f7358SJed Brown 
180552f7358SJed Brown         ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
181552f7358SJed Brown         if (value != 1) continue;
182552f7358SJed Brown       }
183552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
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       }
195552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
196552f7358SJed Brown     }
197552f7358SJed Brown     if (k != numSend) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid number of vertices to send %d should be %d", k, numSend);
198552f7358SJed Brown     ierr = MPI_Send(&numSend, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
199552f7358SJed Brown     ierr = MPI_Send(localVertices, numSend, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
200552f7358SJed Brown     ierr = PetscFree(localVertices);CHKERRQ(ierr);
201552f7358SJed Brown   }
202552f7358SJed Brown   ierr = ISRestoreIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr);
203552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "CELL_TYPES %d\n", totCells);CHKERRQ(ierr);
204552f7358SJed Brown   if (!rank) {
205552f7358SJed Brown     PetscInt cellType;
206552f7358SJed Brown 
207552f7358SJed Brown     for (c = 0; c < numCells; ++c) {
208552f7358SJed Brown       ierr = DMPlexVTKGetCellType(dm, dim, corners[c], &cellType);CHKERRQ(ierr);
209552f7358SJed Brown       ierr = PetscFPrintf(comm, fp, "%d\n", cellType);CHKERRQ(ierr);
210552f7358SJed Brown     }
211552f7358SJed Brown     for(proc = 1; proc < numProcs; ++proc) {
212552f7358SJed Brown       MPI_Status status;
213552f7358SJed Brown 
214552f7358SJed Brown       ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
215552f7358SJed Brown       ierr = MPI_Recv(corners, numCells, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
216552f7358SJed Brown       for(c = 0; c < numCells; ++c) {
217552f7358SJed Brown         ierr = DMPlexVTKGetCellType(dm, dim, corners[c], &cellType);CHKERRQ(ierr);
218552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "%d\n", cellType);CHKERRQ(ierr);
219552f7358SJed Brown       }
220552f7358SJed Brown     }
221552f7358SJed Brown   } else {
222552f7358SJed Brown     ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
223552f7358SJed Brown     ierr = MPI_Send(corners, numCells, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
224552f7358SJed Brown   }
225552f7358SJed Brown   ierr = PetscFree(corners);CHKERRQ(ierr);
226552f7358SJed Brown   *totalCells = totCells;
227552f7358SJed Brown   PetscFunctionReturn(0);
228552f7358SJed Brown }
229552f7358SJed Brown 
230552f7358SJed Brown #undef __FUNCT__
231552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteSection_ASCII"
232552f7358SJed Brown PetscErrorCode DMPlexVTKWriteSection_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec v, FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale) {
233552f7358SJed Brown   MPI_Comm           comm    = ((PetscObject) dm)->comm;
234552f7358SJed Brown   const MPI_Datatype mpiType = MPIU_SCALAR;
235552f7358SJed Brown   PetscScalar       *array;
236552f7358SJed Brown   PetscInt           numDof = 0, maxDof;
237552f7358SJed Brown   PetscInt           numLabelCells, cellHeight, cMax, cStart, cEnd, numLabelVertices, vMax, vStart, vEnd, pStart, pEnd, p;
238552f7358SJed Brown   PetscMPIInt        numProcs, rank, proc, tag;
239552f7358SJed Brown   PetscBool          hasLabel;
240552f7358SJed Brown   PetscErrorCode     ierr;
241552f7358SJed Brown 
242552f7358SJed Brown   PetscFunctionBegin;
243552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
244552f7358SJed Brown   PetscValidHeaderSpecific(v,VEC_CLASSID,4);
245552f7358SJed Brown   if (precision < 0) precision = 6;
246552f7358SJed Brown   ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr);
247552f7358SJed Brown   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
248552f7358SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
249552f7358SJed Brown   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
250552f7358SJed Brown   /* VTK only wants the values at cells or vertices */
251552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
252552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
253552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
254*770b213bSMatthew G Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, PETSC_NULL, PETSC_NULL, &vMax);CHKERRQ(ierr);
255552f7358SJed Brown   if (cMax >= 0) {cEnd = PetscMin(cEnd, cMax);}
256552f7358SJed Brown   if (vMax >= 0) {vEnd = PetscMin(vEnd, vMax);}
257552f7358SJed Brown   pStart = PetscMax(PetscMin(cStart, vStart), pStart);
258552f7358SJed Brown   pEnd   = PetscMin(PetscMax(cEnd,   vEnd),   pEnd);
259552f7358SJed Brown   ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
260552f7358SJed Brown   ierr = DMPlexGetStratumSize(dm, "vtk", 2, &numLabelVertices);CHKERRQ(ierr);
261552f7358SJed Brown   hasLabel = numLabelCells > 0 || numLabelVertices > 0 ? PETSC_TRUE : PETSC_FALSE;
262552f7358SJed Brown   for(p = pStart; p < pEnd; ++p) {
263552f7358SJed Brown     /* Reject points not either cells or vertices */
264552f7358SJed Brown     if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
265552f7358SJed Brown     if (hasLabel) {
266552f7358SJed Brown       PetscInt value;
267552f7358SJed Brown 
268552f7358SJed Brown       if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
269552f7358SJed Brown           ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
270552f7358SJed Brown         ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
271552f7358SJed Brown         if (value != 1) continue;
272552f7358SJed Brown       }
273552f7358SJed Brown     }
274552f7358SJed Brown     ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr);
275552f7358SJed Brown     if (numDof) {break;}
276552f7358SJed Brown   }
277552f7358SJed Brown   ierr = MPI_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
278552f7358SJed Brown   enforceDof = PetscMax(enforceDof, maxDof);
279552f7358SJed Brown   ierr = VecGetArray(v, &array);CHKERRQ(ierr);
280552f7358SJed Brown   if (!rank) {
281552f7358SJed Brown     char formatString[8];
282552f7358SJed Brown 
283552f7358SJed Brown     ierr = PetscSNPrintf(formatString, 8, "%%.%de", precision);CHKERRQ(ierr);
284552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
285552f7358SJed Brown       /* Here we lose a way to filter points by keeping them out of the Numbering */
286552f7358SJed Brown       PetscInt dof, off, goff, d;
287552f7358SJed Brown 
288552f7358SJed Brown       /* Reject points not either cells or vertices */
289552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
290552f7358SJed Brown       if (hasLabel) {
291552f7358SJed Brown         PetscInt value;
292552f7358SJed Brown 
293552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
294552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
295552f7358SJed Brown           ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
296552f7358SJed Brown           if (value != 1) continue;
297552f7358SJed Brown         }
298552f7358SJed Brown       }
299552f7358SJed Brown       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
300552f7358SJed Brown       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
301552f7358SJed Brown       ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
302552f7358SJed Brown       if (dof && goff >= 0) {
303552f7358SJed Brown         for (d = 0; d < dof; d++) {
304552f7358SJed Brown           if (d > 0) {
305552f7358SJed Brown             ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr);
306552f7358SJed Brown           }
307552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(array[off+d])*scale);CHKERRQ(ierr);
308552f7358SJed Brown         }
309552f7358SJed Brown         for (d = dof; d < enforceDof; d++) {
310552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr);
311552f7358SJed Brown         }
312552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
313552f7358SJed Brown       }
314552f7358SJed Brown     }
315552f7358SJed Brown     for (proc = 1; proc < numProcs; ++proc) {
316552f7358SJed Brown       PetscScalar *remoteValues;
317552f7358SJed Brown       PetscInt     size, d;
318552f7358SJed Brown       MPI_Status   status;
319552f7358SJed Brown 
320552f7358SJed Brown       ierr = MPI_Recv(&size, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr);
321552f7358SJed Brown       ierr = PetscMalloc(size * sizeof(PetscScalar), &remoteValues);CHKERRQ(ierr);
322552f7358SJed Brown       ierr = MPI_Recv(remoteValues, size, mpiType, proc, tag, comm, &status);CHKERRQ(ierr);
323552f7358SJed Brown       for (p = 0; p < size/maxDof; ++p) {
324552f7358SJed Brown         for (d = 0; d < maxDof; ++d) {
325552f7358SJed Brown           if (d > 0) {
326552f7358SJed Brown             ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr);
327552f7358SJed Brown           }
328552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(remoteValues[p*maxDof+d])*scale);CHKERRQ(ierr);
329552f7358SJed Brown         }
330552f7358SJed Brown         for (d = maxDof; d < enforceDof; ++d) {
331552f7358SJed Brown           ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr);
332552f7358SJed Brown         }
333552f7358SJed Brown         ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr);
334552f7358SJed Brown       }
335552f7358SJed Brown       ierr = PetscFree(remoteValues);CHKERRQ(ierr);
336552f7358SJed Brown     }
337552f7358SJed Brown   } else {
338552f7358SJed Brown     PetscScalar *localValues;
339552f7358SJed Brown     PetscInt     size, k = 0;
340552f7358SJed Brown 
341552f7358SJed Brown     ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
342552f7358SJed Brown     ierr = PetscMalloc(size * sizeof(PetscScalar), &localValues);CHKERRQ(ierr);
343552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
344552f7358SJed Brown       PetscInt dof, off, goff, d;
345552f7358SJed Brown 
346552f7358SJed Brown       /* Reject points not either cells or vertices */
347552f7358SJed Brown       if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue;
348552f7358SJed Brown       if (hasLabel) {
349552f7358SJed Brown         PetscInt value;
350552f7358SJed Brown 
351552f7358SJed Brown         if (((p >= cStart) && (p < cEnd) && numLabelCells) ||
352552f7358SJed Brown             ((p >= vStart) && (p < vEnd) && numLabelVertices)) {
353552f7358SJed Brown           ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr);
354552f7358SJed Brown           if (value != 1) continue;
355552f7358SJed Brown         }
356552f7358SJed Brown       }
357552f7358SJed Brown       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
358552f7358SJed Brown       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
359552f7358SJed Brown       ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
360552f7358SJed Brown       if (goff >= 0) {
361552f7358SJed Brown         for (d = 0; d < dof; ++d) {
362552f7358SJed Brown           localValues[k++] = array[off+d];
363552f7358SJed Brown         }
364552f7358SJed Brown       }
365552f7358SJed Brown     }
366552f7358SJed Brown     ierr = MPI_Send(&k, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr);
367552f7358SJed Brown     ierr = MPI_Send(localValues, k, mpiType, 0, tag, comm);CHKERRQ(ierr);
368552f7358SJed Brown     ierr = PetscFree(localValues);CHKERRQ(ierr);
369552f7358SJed Brown   }
370552f7358SJed Brown   ierr = VecRestoreArray(v, &array);CHKERRQ(ierr);
371552f7358SJed Brown   PetscFunctionReturn(0);
372552f7358SJed Brown }
373552f7358SJed Brown 
374552f7358SJed Brown #undef __FUNCT__
375552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteField_ASCII"
376552f7358SJed Brown PetscErrorCode DMPlexVTKWriteField_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec field, const char name[], FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale)
377552f7358SJed Brown {
378552f7358SJed Brown   MPI_Comm       comm = ((PetscObject) dm)->comm;
379552f7358SJed Brown   PetscInt       numDof = 0, maxDof;
380552f7358SJed Brown   PetscInt       pStart, pEnd, p;
381552f7358SJed Brown   PetscErrorCode ierr;
382552f7358SJed Brown 
383552f7358SJed Brown   PetscFunctionBegin;
384552f7358SJed Brown   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
385552f7358SJed Brown   for (p = pStart; p < pEnd; ++p) {
386552f7358SJed Brown     ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr);
387552f7358SJed Brown     if (numDof) {break;}
388552f7358SJed Brown   }
389552f7358SJed Brown   numDof = PetscMax(numDof, enforceDof);
390552f7358SJed Brown   ierr = MPI_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, ((PetscObject) dm)->comm);CHKERRQ(ierr);
391552f7358SJed Brown   if (!name) name = "Unknown";
392552f7358SJed Brown   if (maxDof == 3) {
393552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "VECTORS %s double\n", name);CHKERRQ(ierr);
394552f7358SJed Brown   } else {
395552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "SCALARS %s double %d\n", name, maxDof);CHKERRQ(ierr);
396552f7358SJed Brown     ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr);
397552f7358SJed Brown   }
398552f7358SJed Brown   ierr = DMPlexVTKWriteSection_ASCII(dm, section, globalSection, field, fp, enforceDof, precision, scale);CHKERRQ(ierr);
399552f7358SJed Brown   PetscFunctionReturn(0);
400552f7358SJed Brown }
401552f7358SJed Brown 
402552f7358SJed Brown #undef __FUNCT__
403552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll_ASCII"
404552f7358SJed Brown static PetscErrorCode DMPlexVTKWriteAll_ASCII(DM dm, PetscViewer viewer)
405552f7358SJed Brown {
406552f7358SJed Brown   MPI_Comm                 comm = ((PetscObject) dm)->comm;
407552f7358SJed Brown   PetscViewer_VTK         *vtk  = (PetscViewer_VTK *) viewer->data;
408552f7358SJed Brown   FILE                    *fp;
409552f7358SJed Brown   PetscViewerVTKObjectLink link;
410552f7358SJed Brown   PetscSection             coordSection, globalCoordSection;
411552f7358SJed Brown   PetscLayout              vLayout;
412552f7358SJed Brown   Vec                      coordinates;
413552f7358SJed Brown   PetscReal                lengthScale;
414552f7358SJed Brown   PetscInt                 vMax, totVertices, totCells;
415552f7358SJed Brown   PetscBool                hasPoint = PETSC_FALSE, hasCell = PETSC_FALSE;
416552f7358SJed Brown   PetscErrorCode           ierr;
417552f7358SJed Brown 
418552f7358SJed Brown   PetscFunctionBegin;
419552f7358SJed Brown   ierr = PetscFOpen(comm, vtk->filename, "wb", &fp);CHKERRQ(ierr);
420552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "# vtk DataFile Version 2.0\n");CHKERRQ(ierr);
421552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "Simplicial Mesh Example\n");CHKERRQ(ierr);
422552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "ASCII\n");CHKERRQ(ierr);
423552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "DATASET UNSTRUCTURED_GRID\n");CHKERRQ(ierr);
424552f7358SJed Brown   /* Vertices */
425552f7358SJed Brown   ierr = DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);CHKERRQ(ierr);
426552f7358SJed Brown   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
427552f7358SJed Brown   ierr = PetscSectionCreateGlobalSection(coordSection, dm->sf, PETSC_FALSE, &globalCoordSection);CHKERRQ(ierr);
428552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
429*770b213bSMatthew G Knepley   ierr = DMPlexGetHybridBounds(dm, PETSC_NULL, PETSC_NULL, PETSC_NULL, &vMax);CHKERRQ(ierr);
430552f7358SJed Brown   if (vMax >= 0) {
431552f7358SJed Brown     PetscInt pStart, pEnd, p, localSize = 0;
432552f7358SJed Brown 
433552f7358SJed Brown     ierr = PetscSectionGetChart(globalCoordSection, &pStart, &pEnd);CHKERRQ(ierr);
434552f7358SJed Brown     pEnd = PetscMin(pEnd, vMax);
435552f7358SJed Brown     for (p = pStart; p < pEnd; ++p) {
436552f7358SJed Brown       PetscInt dof;
437552f7358SJed Brown 
438552f7358SJed Brown       ierr = PetscSectionGetDof(globalCoordSection, p, &dof);CHKERRQ(ierr);
439552f7358SJed Brown       if (dof > 0) {++localSize;}
440552f7358SJed Brown     }
441552f7358SJed Brown     ierr = PetscLayoutCreate(((PetscObject) dm)->comm, &vLayout);CHKERRQ(ierr);
442552f7358SJed Brown     ierr = PetscLayoutSetLocalSize(vLayout, localSize);CHKERRQ(ierr);
443552f7358SJed Brown     ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr);
444552f7358SJed Brown     ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr);
445552f7358SJed Brown   } else {
446552f7358SJed Brown     ierr = PetscSectionGetPointLayout(((PetscObject) dm)->comm, globalCoordSection, &vLayout);CHKERRQ(ierr);
447552f7358SJed Brown   }
448552f7358SJed Brown   ierr = PetscLayoutGetSize(vLayout, &totVertices);CHKERRQ(ierr);
449552f7358SJed Brown   ierr = PetscFPrintf(comm, fp, "POINTS %d double\n", totVertices);CHKERRQ(ierr);
450552f7358SJed Brown   ierr = DMPlexVTKWriteSection_ASCII(dm, coordSection, globalCoordSection, coordinates, fp, 3, PETSC_DETERMINE, lengthScale);CHKERRQ(ierr);
451552f7358SJed Brown   /* Cells */
452552f7358SJed Brown   ierr = DMPlexVTKWriteCells_ASCII(dm, fp, &totCells);CHKERRQ(ierr);
453552f7358SJed Brown   /* Vertex fields */
454552f7358SJed Brown   for (link = vtk->link; link; link = link->next) {
455552f7358SJed Brown     if ((link->ft == PETSC_VTK_POINT_FIELD) || (link->ft == PETSC_VTK_POINT_VECTOR_FIELD)) hasPoint = PETSC_TRUE;
456552f7358SJed Brown     if ((link->ft == PETSC_VTK_CELL_FIELD)  || (link->ft == PETSC_VTK_CELL_VECTOR_FIELD))  hasCell  = PETSC_TRUE;
457552f7358SJed Brown   }
458552f7358SJed Brown   if (hasPoint) {
45922d28d08SBarry Smith     ierr = PetscFPrintf(comm, fp, "POINT_DATA %d\n", totVertices);CHKERRQ(ierr);
460552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
461552f7358SJed Brown       Vec          X = (Vec) link->vec;
462552f7358SJed Brown       DM           dmX;
463552f7358SJed Brown       PetscSection section, globalSection;
464552f7358SJed Brown       const char  *name;
465552f7358SJed Brown       PetscInt     enforceDof = PETSC_DETERMINE;
466552f7358SJed Brown 
467552f7358SJed Brown       if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
468552f7358SJed Brown       if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) enforceDof = 3;
469552f7358SJed Brown       ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr);
470552f7358SJed Brown       ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
471552f7358SJed Brown       if (dmX) {
47222d28d08SBarry Smith         ierr = DMGetDefaultSection(dmX, &section);CHKERRQ(ierr);
473552f7358SJed Brown       } else {
474552f7358SJed Brown         PetscContainer c;
475552f7358SJed Brown 
476552f7358SJed Brown         ierr = PetscObjectQuery(link->vec, "section", (PetscObject *) &c);CHKERRQ(ierr);
477552f7358SJed Brown         if (!c) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
478552f7358SJed Brown         ierr = PetscContainerGetPointer(c, (void **) &section);CHKERRQ(ierr);
479552f7358SJed Brown       }
480552f7358SJed Brown       if (!section) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
481552f7358SJed Brown       ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, &globalSection);CHKERRQ(ierr);
482552f7358SJed Brown       ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr);
483552f7358SJed Brown       ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr);
484552f7358SJed Brown     }
485552f7358SJed Brown   }
486552f7358SJed Brown   /* Cell Fields */
487552f7358SJed Brown   if (hasCell) {
48822d28d08SBarry Smith     ierr = PetscFPrintf(comm, fp, "CELL_DATA %d\n", totCells);CHKERRQ(ierr);
489552f7358SJed Brown     for (link = vtk->link; link; link = link->next) {
490552f7358SJed Brown       Vec            X = (Vec) link->vec;
491552f7358SJed Brown       DM             dmX;
492552f7358SJed Brown       PetscSection   section, globalSection;
493552f7358SJed Brown       const char    *name;
494552f7358SJed Brown       PetscInt       enforceDof = PETSC_DETERMINE;
495552f7358SJed Brown 
496552f7358SJed Brown       if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
497552f7358SJed Brown       if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) enforceDof = 3;
498552f7358SJed Brown       ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr);
499552f7358SJed Brown       ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
500552f7358SJed Brown       if (dmX) {
50122d28d08SBarry Smith         ierr = DMGetDefaultSection(dmX, &section);CHKERRQ(ierr);
502552f7358SJed Brown       } else {
503552f7358SJed Brown         PetscContainer c;
504552f7358SJed Brown 
505552f7358SJed Brown         ierr = PetscObjectQuery(link->vec, "section", (PetscObject *) &c);CHKERRQ(ierr);
506552f7358SJed Brown         if (!c) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
507552f7358SJed Brown         ierr = PetscContainerGetPointer(c, (void **) &section);CHKERRQ(ierr);
508552f7358SJed Brown       }
509552f7358SJed Brown       if (!section) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name);
510552f7358SJed Brown       ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, &globalSection);CHKERRQ(ierr);
511552f7358SJed Brown       ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr);
512552f7358SJed Brown       ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr);
513552f7358SJed Brown     }
514552f7358SJed Brown   }
515552f7358SJed Brown   /* Cleanup */
516552f7358SJed Brown   ierr = PetscSectionDestroy(&globalCoordSection);CHKERRQ(ierr);
517552f7358SJed Brown   ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr);
518552f7358SJed Brown   ierr = PetscFClose(comm, fp);CHKERRQ(ierr);
519552f7358SJed Brown   PetscFunctionReturn(0);
520552f7358SJed Brown }
521552f7358SJed Brown 
522552f7358SJed Brown #undef __FUNCT__
523552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll"
524552f7358SJed Brown /*@C
525552f7358SJed Brown   DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
526552f7358SJed Brown 
527552f7358SJed Brown   Collective
528552f7358SJed Brown 
529552f7358SJed Brown   Input Arguments:
530552f7358SJed Brown + odm - The DMPlex specifying the mesh, passed as a PetscObject
531552f7358SJed Brown - viewer - viewer of type VTK
532552f7358SJed Brown 
533552f7358SJed Brown   Level: developer
534552f7358SJed Brown 
535552f7358SJed Brown   Note:
536552f7358SJed Brown   This function is a callback used by the VTK viewer to actually write the file.
537552f7358SJed 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.
538552f7358SJed Brown   Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
539552f7358SJed Brown 
540552f7358SJed Brown .seealso: PETSCVIEWERVTK
541552f7358SJed Brown @*/
542552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer)
543552f7358SJed Brown {
544552f7358SJed Brown   DM              dm   = (DM) odm;
545552f7358SJed Brown   PetscBool       isvtk;
546552f7358SJed Brown   PetscErrorCode  ierr;
547552f7358SJed Brown 
548552f7358SJed Brown   PetscFunctionBegin;
549552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
550552f7358SJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
551552f7358SJed Brown   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr);
552552f7358SJed Brown   if (!isvtk) SETERRQ1(((PetscObject) viewer)->comm, PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name);
553552f7358SJed Brown   switch (viewer->format) {
554552f7358SJed Brown   case PETSC_VIEWER_ASCII_VTK:
555552f7358SJed Brown     ierr = DMPlexVTKWriteAll_ASCII(dm, viewer);CHKERRQ(ierr);
556552f7358SJed Brown     break;
557552f7358SJed Brown   case PETSC_VIEWER_VTK_VTU:
558552f7358SJed Brown     ierr = DMPlexVTKWriteAll_VTU(dm, viewer);CHKERRQ(ierr);
559552f7358SJed Brown     break;
560552f7358SJed Brown   default: SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]);
561552f7358SJed Brown   }
562552f7358SJed Brown   PetscFunctionReturn(0);
563552f7358SJed Brown }
564