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