xref: /petsc/src/dm/impls/plex/plexvtu.c (revision 2e529ec88214be11db4abfc2ebf7b2fbd6e371cc)
1 #include <petsc/private/dmpleximpl.h>
2 #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h>
3 
4 typedef struct {
5   PetscInt nvertices;
6   PetscInt ncells;
7   PetscInt nconn;               /* number of entries in cell->vertex connectivity array */
8 } PieceInfo;
9 
10 #if defined(PETSC_USE_REAL_SINGLE)
11 static const char precision[] = "Float32";
12 #elif defined(PETSC_USE_REAL_DOUBLE)
13 static const char precision[] = "Float64";
14 #else
15 static const char precision[] = "UnknownPrecision";
16 #endif
17 
18 static PetscErrorCode TransferWrite(PetscViewer viewer,FILE *fp,PetscMPIInt srank,PetscMPIInt root,const void *send,void *recv,PetscMPIInt count,PetscDataType datatype,PetscMPIInt tag)
19 {
20   PetscMPIInt    rank;
21   PetscErrorCode ierr;
22   MPI_Comm       comm;
23   MPI_Datatype   mpidatatype;
24 
25   PetscFunctionBegin;
26   ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
27   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
28   ierr = PetscDataTypeToMPIDataType(datatype,&mpidatatype);CHKERRQ(ierr);
29 
30   if (rank == srank && rank != root) {
31     ierr = MPI_Send((void*)send,count,mpidatatype,root,tag,comm);CHKERRQ(ierr);
32   } else if (rank == root) {
33     const void *buffer;
34     if (root == srank) {        /* self */
35       buffer = send;
36     } else {
37       MPI_Status  status;
38       PetscMPIInt nrecv;
39       ierr = MPI_Recv(recv,count,mpidatatype,srank,tag,comm,&status);CHKERRQ(ierr);
40       ierr = MPI_Get_count(&status,mpidatatype,&nrecv);CHKERRQ(ierr);
41       if (count != nrecv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
42       buffer = recv;
43     }
44     ierr = PetscViewerVTKFWrite(viewer,fp,buffer,count,datatype);CHKERRQ(ierr);
45   }
46   PetscFunctionReturn(0);
47 }
48 
49 static PetscErrorCode DMPlexGetVTKConnectivity(DM dm, PetscBool localized, PieceInfo *piece,PetscVTKInt **oconn,PetscVTKInt **ooffsets,PetscVTKType **otypes)
50 {
51   PetscErrorCode ierr;
52   PetscSection   coordSection;
53   PetscVTKInt    *conn,*offsets;
54   PetscVTKType   *types;
55   PetscInt       dim,vStart,vEnd,cStart,cEnd,pStart,pEnd,cellHeight,cMax,numLabelCells,hasLabel,c,v,countcell,countconn;
56 
57   PetscFunctionBegin;
58   ierr = PetscMalloc3(piece->nconn,&conn,piece->ncells,&offsets,piece->ncells,&types);CHKERRQ(ierr);
59   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
60   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
61   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
62   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
63   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
64   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
65   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
66   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
67   ierr     = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
68   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
69 
70   countcell = 0;
71   countconn = 0;
72   for (c = cStart; c < cEnd; ++c) {
73     PetscInt nverts,celltype,startoffset,nC=0;
74 
75     if (hasLabel) {
76       PetscInt value;
77 
78       ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
79       if (value != 1) continue;
80     }
81     startoffset = countconn;
82     if (!localized) {
83       PetscInt *closure = NULL;
84       PetscInt  closureSize;
85 
86       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
87       for (v = 0; v < closureSize*2; v += 2) {
88         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
89           conn[countconn++] = closure[v] - vStart;
90           ++nC;
91         }
92       }
93       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
94     } else {
95       PetscInt dof;
96 
97       ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
98       for (nC = 0; nC < dof/dim; nC++) conn[countconn++] = startoffset + nC;
99     }
100     ierr = DMPlexInvertCell(dim, nC, &conn[countconn-nC]);CHKERRQ(ierr);
101 
102     offsets[countcell] = countconn;
103 
104     nverts = countconn - startoffset;
105     ierr   = DMPlexVTKGetCellType_Internal(dm,dim,nverts,&celltype);CHKERRQ(ierr);
106 
107     types[countcell] = celltype;
108     countcell++;
109   }
110   if (countcell != piece->ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count");
111   if (countconn != piece->nconn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count");
112   *oconn    = conn;
113   *ooffsets = offsets;
114   *otypes   = types;
115   PetscFunctionReturn(0);
116 }
117 
118 /*
119   Write all fields that have been provided to the viewer
120   Multi-block XML format with binary appended data.
121 */
122 PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer)
123 {
124   MPI_Comm                 comm;
125   PetscSection             coordSection;
126   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
127   PetscViewerVTKObjectLink link;
128   FILE                     *fp;
129   PetscMPIInt              rank,size,tag;
130   PetscErrorCode           ierr;
131   PetscInt                 dim, dimEmbed,cellHeight,cStart,cEnd,vStart,vEnd,cMax,numLabelCells,hasLabel,c,v,r,i;
132   PetscBool                localized;
133   PieceInfo                piece,*gpiece = NULL;
134   void                     *buffer = NULL;
135 
136   PetscFunctionBegin;
137   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
138 #if defined(PETSC_USE_COMPLEX)
139   SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Complex values not supported");
140 #endif
141   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
142   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
143   tag  = ((PetscObject)viewer)->tag;
144 
145   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
146   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
147 #if defined(PETSC_WORDS_BIGENDIAN)
148   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr);
149 #else
150   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr);
151 #endif
152   ierr = PetscFPrintf(comm,fp,"  <UnstructuredGrid>\n");CHKERRQ(ierr);
153 
154   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
155   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
156   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
157   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
158   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
159   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
160   ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
161   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
162   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
163 
164   hasLabel        = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
165   piece.nvertices = 0;
166   piece.ncells    = 0;
167   piece.nconn     = 0;
168   if (!localized) piece.nvertices = vEnd - vStart;
169   for (c = cStart; c < cEnd; ++c) {
170     if (hasLabel) {
171       PetscInt value;
172 
173       ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
174       if (value != 1) continue;
175     }
176     if (!localized) {
177       PetscInt *closure = NULL;
178       PetscInt closureSize;
179 
180       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
181       for (v = 0; v < closureSize*2; v += 2) {
182         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
183           piece.nconn++;
184         }
185       }
186       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
187     } else {
188       PetscInt dof;
189 
190       ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
191       piece.nvertices += dof/dimEmbed;
192       piece.nconn     += dof/dimEmbed;
193     }
194     piece.ncells++;
195   }
196   if (!rank) {ierr = PetscMalloc1(size,&gpiece);CHKERRQ(ierr);}
197   ierr = MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm);CHKERRQ(ierr);
198 
199   /*
200    * Write file header
201    */
202   if (!rank) {
203     PetscInt boffset = 0;
204 
205     for (r=0; r<size; r++) {
206       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells);CHKERRQ(ierr);
207       /* Coordinate positions */
208       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Points>\n");CHKERRQ(ierr);
209       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
210       boffset += gpiece[r].nvertices*3*sizeof(PetscScalar) + sizeof(int);
211       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Points>\n");CHKERRQ(ierr);
212       /* Cell connectivity */
213       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Cells>\n");CHKERRQ(ierr);
214       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
215       boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int);
216       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"offsets\"      NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
217       boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int);
218       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"UInt8\" Name=\"types\"        NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
219       boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int);
220       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Cells>\n");CHKERRQ(ierr);
221 
222       /*
223        * Cell Data headers
224        */
225       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <CellData>\n");CHKERRQ(ierr);
226       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
227       boffset += gpiece[r].ncells*sizeof(int) + sizeof(int);
228       /* all the vectors */
229       for (link=vtk->link; link; link=link->next) {
230         Vec        X = (Vec)link->vec;
231         PetscInt   bs,nfields,field;
232         const char *vecname = "";
233         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
234         if (((PetscObject)X)->name || link != vtk->link) { /* If the object is already named, use it. If it is past the first link, name it to disambiguate. */
235           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
236         }
237         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
238         ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr);
239         for (field=0,i=0; field<(nfields?nfields:1); field++) {
240           PetscInt     fbs,j;
241           PetscFV      fv = NULL;
242           PetscObject  f;
243           PetscClassId fClass;
244           const char *fieldname = NULL;
245           char       buf[256];
246           if (nfields) {        /* We have user-defined fields/components */
247             ierr = PetscSectionGetFieldDof(dm->defaultSection,cStart,field,&fbs);CHKERRQ(ierr);
248             ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr);
249           } else fbs = bs;      /* Say we have one field with 'bs' components */
250           ierr = DMGetField(dm,field,&f);CHKERRQ(ierr);
251           ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr);
252           if (fClass == PETSCFV_CLASSID) {
253             fv = (PetscFV) f;
254           }
255           if (!fieldname) {
256             ierr = PetscSNPrintf(buf,sizeof(buf),"CellField%D",field);CHKERRQ(ierr);
257             fieldname = buf;
258           }
259           for (j=0; j<fbs; j++) {
260             const char *compName = NULL;
261             if (fv) {
262               ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
263             }
264             if (compName) {
265               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,compName,boffset);CHKERRQ(ierr);
266             }
267             else {
268               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s.%D\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,j,boffset);CHKERRQ(ierr);
269             }
270             boffset += gpiece[r].ncells*sizeof(PetscScalar) + sizeof(int);
271             i++;
272           }
273         }
274         if (i != bs) SETERRQ2(comm,PETSC_ERR_PLIB,"Total number of field components %D != block size %D",i,bs);
275       }
276       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </CellData>\n");CHKERRQ(ierr);
277 
278       /*
279        * Point Data headers
280        */
281       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      <PointData>\n");CHKERRQ(ierr);
282       for (link=vtk->link; link; link=link->next) {
283         Vec        X = (Vec)link->vec;
284         PetscInt   bs,nfields,field;
285         const char *vecname = "";
286         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
287         if (((PetscObject)X)->name || link != vtk->link) { /* If the object is already named, use it. If it is past the first link, name it to disambiguate. */
288           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
289         }
290         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
291         ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr);
292         for (field=0,i=0; field<(nfields?nfields:1); field++) {
293           PetscInt   fbs,j;
294           const char *fieldname = NULL;
295           char       buf[256];
296           if (nfields) {        /* We have user-defined fields/components */
297             ierr = PetscSectionGetFieldDof(dm->defaultSection,vStart,field,&fbs);CHKERRQ(ierr);
298             ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr);
299           } else fbs = bs;      /* Say we have one field with 'bs' components */
300           if (!fieldname) {
301             ierr = PetscSNPrintf(buf,sizeof(buf),"PointField%D",field);CHKERRQ(ierr);
302             fieldname = buf;
303           }
304           for (j=0; j<fbs; j++) {
305             ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s.%D\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,j,boffset);CHKERRQ(ierr);
306             boffset += gpiece[r].nvertices*sizeof(PetscScalar) + sizeof(int);
307           }
308         }
309       }
310       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </PointData>\n");CHKERRQ(ierr);
311 
312       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    </Piece>\n");CHKERRQ(ierr);
313     }
314   }
315 
316   ierr = PetscFPrintf(comm,fp,"  </UnstructuredGrid>\n");CHKERRQ(ierr);
317   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
318   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
319 
320   if (!rank) {
321     PetscInt maxsize = 0;
322     for (r=0; r<size; r++) {
323       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscScalar)));
324       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*sizeof(PetscScalar)));
325       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt)));
326     }
327     ierr = PetscMalloc(maxsize,&buffer);CHKERRQ(ierr);
328   }
329   for (r=0; r<size; r++) {
330     if (r == rank) {
331       PetscInt nsend;
332       {                         /* Position */
333         const PetscScalar *x;
334         PetscScalar       *y = NULL;
335         Vec               coords;
336 
337         ierr  = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
338         ierr  = VecGetArrayRead(coords,&x);CHKERRQ(ierr);
339         if (dimEmbed != 3) {
340           ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr);
341           for (i=0; i<piece.nvertices; i++) {
342             y[i*3+0] = x[i*dimEmbed+0];
343             y[i*3+1] = (dimEmbed > 1) ? x[i*dimEmbed+1] : 0;
344             y[i*3+2] = 0;
345           }
346         }
347         if (localized) {
348           PetscInt cnt = 0;
349           if (!y) {
350             PetscInt off;
351 
352             ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr);
353             ierr = PetscSectionGetOffset(coordSection, cStart, &off);CHKERRQ(ierr);
354             ierr = PetscMemcpy(y,x+off,piece.nvertices*3*sizeof(PetscScalar));CHKERRQ(ierr);
355           }
356           for (c=cStart; c<cEnd; c++) {
357             PetscInt dof, off;
358 
359             ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
360             ierr = PetscSectionGetOffset(coordSection, c, &off);CHKERRQ(ierr);
361             cnt += dof/dimEmbed;
362           }
363           if (cnt != piece.nvertices) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match %D != %D",cnt,piece.nvertices);
364         }
365         nsend = piece.nvertices*3;
366         ierr  = TransferWrite(viewer,fp,r,0,y ? y : x,buffer,nsend,PETSC_SCALAR,tag);CHKERRQ(ierr);
367         ierr  = PetscFree(y);CHKERRQ(ierr);
368         ierr  = VecRestoreArrayRead(coords,&x);CHKERRQ(ierr);
369       }
370       {                           /* Connectivity, offsets, types */
371         PetscVTKInt  *connectivity = NULL, *offsets = NULL;
372         PetscVTKType *types = NULL;
373         ierr = DMPlexGetVTKConnectivity(dm,localized,&piece,&connectivity,&offsets,&types);CHKERRQ(ierr);
374         ierr = TransferWrite(viewer,fp,r,0,connectivity,buffer,piece.nconn,PETSC_INT32,tag);CHKERRQ(ierr);
375         ierr = TransferWrite(viewer,fp,r,0,offsets,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr);
376         ierr = TransferWrite(viewer,fp,r,0,types,buffer,piece.ncells,PETSC_UINT8,tag);CHKERRQ(ierr);
377         ierr = PetscFree3(connectivity,offsets,types);CHKERRQ(ierr);
378       }
379       {                         /* Owners (cell data) */
380         PetscVTKInt *owners;
381         ierr = PetscMalloc1(piece.ncells,&owners);CHKERRQ(ierr);
382         for (i=0; i<piece.ncells; i++) owners[i] = rank;
383         ierr = TransferWrite(viewer,fp,r,0,owners,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr);
384         ierr = PetscFree(owners);CHKERRQ(ierr);
385       }
386       /* Cell data */
387       for (link=vtk->link; link; link=link->next) {
388         Vec               X = (Vec)link->vec;
389         const PetscScalar *x;
390         PetscScalar       *y;
391         PetscInt          bs;
392         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
393         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
394         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
395         ierr = PetscMalloc1(piece.ncells,&y);CHKERRQ(ierr);
396         for (i=0; i<bs; i++) {
397           PetscInt cnt;
398           for (c=cStart,cnt=0; c<cEnd; c++) {
399             PetscScalar *xpoint;
400             if (hasLabel) {     /* Ignore some cells */
401               PetscInt value;
402               ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
403               if (value != 1) continue;
404             }
405             ierr     = DMPlexPointLocalRead(dm,c,x,&xpoint);CHKERRQ(ierr);
406             y[cnt++] = xpoint[i];
407           }
408           if (cnt != piece.ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
409           ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells,PETSC_SCALAR,tag);CHKERRQ(ierr);
410         }
411         ierr = PetscFree(y);CHKERRQ(ierr);
412         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
413       }
414 
415       for (link=vtk->link; link; link=link->next) {
416         Vec               X = (Vec)link->vec;
417         const PetscScalar *x;
418         PetscScalar       *y;
419         PetscInt          bs;
420         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
421         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
422         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
423         ierr = PetscMalloc1(piece.nvertices,&y);CHKERRQ(ierr);
424         for (i=0; i<bs; i++) {
425           PetscInt cnt;
426           if (!localized) {
427             for (v=vStart,cnt=0; v<vEnd; v++) {
428               PetscScalar *xpoint;
429               ierr     = DMPlexPointLocalRead(dm,v,x,&xpoint);CHKERRQ(ierr);
430               y[cnt++] = xpoint[i];
431             }
432           } else {
433             for (c=cStart,cnt=0; c<cEnd; c++) {
434               PetscInt *closure = NULL;
435               PetscInt  closureSize, off;
436 
437               ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
438               for (v = 0, off = 0; v < closureSize*2; v += 2) {
439                 if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
440                   PetscScalar *xpoint;
441 
442                   ierr = DMPlexPointLocalRead(dm,v,x,&xpoint);CHKERRQ(ierr);
443                   y[cnt + off++] = xpoint[i];
444                 }
445               }
446               ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
447             }
448           }
449           if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
450           ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr);
451         }
452         ierr = PetscFree(y);CHKERRQ(ierr);
453         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
454       }
455     } else if (!rank) {
456       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,PETSC_SCALAR,tag);CHKERRQ(ierr); /* positions */
457       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nconn,PETSC_INT32,tag);CHKERRQ(ierr); /* connectivity */
458       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* offsets */
459       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_UINT8,tag);CHKERRQ(ierr); /* types */
460       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* owner rank (cells) */
461       /* all cell data */
462       for (link=vtk->link; link; link=link->next) {
463         PetscInt bs;
464         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
465         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
466         for (i=0; i<bs; i++) {
467           ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_SCALAR,tag);CHKERRQ(ierr);
468         }
469       }
470       /* all point data */
471       for (link=vtk->link; link; link=link->next) {
472         PetscInt bs;
473         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
474         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
475         for (i=0; i<bs; i++) {
476           ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr);
477         }
478       }
479     }
480   }
481   ierr = PetscFree(gpiece);CHKERRQ(ierr);
482   ierr = PetscFree(buffer);CHKERRQ(ierr);
483   ierr = PetscFPrintf(comm,fp,"\n  </AppendedData>\n");CHKERRQ(ierr);
484   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
485   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
486   PetscFunctionReturn(0);
487 }
488