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