xref: /petsc/src/dm/impls/plex/plexvtu.c (revision 519f805a543c2a7f195bcba21173bb2cdfadb956)
1552f7358SJed Brown #include <petsc-private/pleximpl.h>
2552f7358SJed Brown #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h>
3552f7358SJed Brown 
4552f7358SJed Brown typedef struct {
5552f7358SJed Brown   PetscInt nvertices;
6552f7358SJed Brown   PetscInt ncells;
7552f7358SJed Brown   PetscInt nconn;               /* number of entries in cell->vertex connectivity array */
8552f7358SJed Brown } PieceInfo;
9552f7358SJed Brown 
10552f7358SJed Brown #if defined(PETSC_USE_REAL_SINGLE)
11552f7358SJed Brown   static const char precision[]  = "Float32";
12552f7358SJed Brown #elif defined(PETSC_USE_REAL_DOUBLE)
13552f7358SJed Brown   static const char precision[]  = "Float64";
14552f7358SJed Brown #else
15552f7358SJed Brown   static const char precision[]  = "UnknownPrecision";
16552f7358SJed Brown #endif
17552f7358SJed Brown 
18552f7358SJed Brown #undef __FUNCT__
19552f7358SJed Brown #define __FUNCT__ "TransferWrite"
20552f7358SJed Brown static PetscErrorCode TransferWrite(PetscViewer viewer,FILE *fp,PetscMPIInt srank,PetscMPIInt root,const void *send,void *recv,PetscMPIInt count,PetscDataType datatype,PetscMPIInt tag)
21552f7358SJed Brown {
22552f7358SJed Brown   PetscMPIInt rank;
23552f7358SJed Brown   PetscErrorCode ierr;
24552f7358SJed Brown   MPI_Comm comm = ((PetscObject)viewer)->comm;
25552f7358SJed Brown   MPI_Datatype mpidatatype;
26552f7358SJed Brown 
27552f7358SJed Brown   PetscFunctionBegin;
28552f7358SJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
29552f7358SJed Brown   ierr = PetscDataTypeToMPIDataType(datatype,&mpidatatype);CHKERRQ(ierr);
30552f7358SJed Brown 
31552f7358SJed Brown   if (rank == srank && rank != root) {
32552f7358SJed Brown     ierr = MPI_Send((void*)send,count,mpidatatype,root,tag,comm);CHKERRQ(ierr);
33552f7358SJed Brown   } else if (rank == root) {
34552f7358SJed Brown     const void *buffer;
35552f7358SJed Brown     if (root == srank) {        /* self */
36552f7358SJed Brown       buffer = send;
37552f7358SJed Brown     } else {
38552f7358SJed Brown       MPI_Status status;
39552f7358SJed Brown       PetscMPIInt nrecv;
40552f7358SJed Brown       ierr = MPI_Recv(recv,count,mpidatatype,srank,tag,comm,&status);CHKERRQ(ierr);
41552f7358SJed Brown       ierr = MPI_Get_count(&status,mpidatatype,&nrecv);CHKERRQ(ierr);
42552f7358SJed Brown       if (count != nrecv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
43552f7358SJed Brown       buffer = recv;
44552f7358SJed Brown     }
45552f7358SJed Brown     ierr = PetscViewerVTKFWrite(viewer,fp,buffer,count,datatype);CHKERRQ(ierr);
46552f7358SJed Brown   }
47552f7358SJed Brown   PetscFunctionReturn(0);
48552f7358SJed Brown }
49552f7358SJed Brown 
50552f7358SJed Brown #undef __FUNCT__
51552f7358SJed Brown #define __FUNCT__ "DMPlexGetVTKConnectivity"
52552f7358SJed Brown static PetscErrorCode DMPlexGetVTKConnectivity(DM dm,PieceInfo *piece,PetscVTKInt **oconn,PetscVTKInt **ooffsets,PetscVTKType **otypes)
53552f7358SJed Brown {
54552f7358SJed Brown   PetscErrorCode ierr;
55552f7358SJed Brown   PetscVTKInt *conn,*offsets;
56552f7358SJed Brown   PetscVTKType *types;
57552f7358SJed Brown   PetscInt dim,vStart,vEnd,cStart,cEnd,pStart,pEnd,cellHeight,cMax,numLabelCells,hasLabel,c,v,countcell,countconn;
58552f7358SJed Brown 
59552f7358SJed Brown   PetscFunctionBegin;
60552f7358SJed Brown   ierr = PetscMalloc3(piece->nconn,PetscVTKInt,&conn,piece->ncells,PetscVTKInt,&offsets,piece->ncells,PetscVTKType,&types);CHKERRQ(ierr);
61552f7358SJed Brown 
62552f7358SJed Brown   ierr = DMPlexGetDimension(dm,&dim);CHKERRQ(ierr);
63552f7358SJed Brown   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
64552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
65552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
66552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
67770b213bSMatthew G Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
68552f7358SJed Brown   if (cMax >= 0) {cEnd = PetscMin(cEnd, cMax);}
69552f7358SJed Brown   ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
70552f7358SJed Brown   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
71552f7358SJed Brown 
72552f7358SJed Brown   countcell = 0;
73552f7358SJed Brown   countconn = 0;
74552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
75552f7358SJed Brown     PetscInt *closure = PETSC_NULL;
76552f7358SJed Brown     PetscInt  closureSize,nverts,celltype,startoffset;
77552f7358SJed Brown 
78552f7358SJed Brown     if (hasLabel) {
79552f7358SJed Brown       PetscInt value;
80552f7358SJed Brown 
81552f7358SJed Brown       ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
82552f7358SJed Brown       if (value != 1) continue;
83552f7358SJed Brown     }
84552f7358SJed Brown     startoffset = countconn;
85552f7358SJed Brown     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
86552f7358SJed Brown     for (v = 0; v < closureSize*2; v += 2) {
87552f7358SJed Brown       if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
88552f7358SJed Brown         conn[countconn++] = closure[v] - vStart;
89552f7358SJed Brown       }
90552f7358SJed Brown     }
91552f7358SJed Brown     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
92552f7358SJed Brown     offsets[countcell] = countconn;
93552f7358SJed Brown     nverts = countconn - startoffset;
94552f7358SJed Brown     ierr = DMPlexVTKGetCellType(dm,dim,nverts,&celltype);CHKERRQ(ierr);
95552f7358SJed Brown     types[countcell] = celltype;
96552f7358SJed Brown     countcell++;
97552f7358SJed Brown   }
98552f7358SJed Brown   if (countcell != piece->ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count");
99552f7358SJed Brown   if (countconn != piece->nconn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count");
100552f7358SJed Brown   *oconn = conn;
101552f7358SJed Brown   *ooffsets = offsets;
102552f7358SJed Brown   *otypes = types;
103552f7358SJed Brown   PetscFunctionReturn(0);
104552f7358SJed Brown }
105552f7358SJed Brown 
106552f7358SJed Brown #undef __FUNCT__
107552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll_VTU"
108552f7358SJed Brown /*
109552f7358SJed Brown   Write all fields that have been provided to the viewer
110552f7358SJed Brown   Multi-block XML format with binary appended data.
111552f7358SJed Brown */
112552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer)
113552f7358SJed Brown {
114552f7358SJed Brown   MPI_Comm comm = ((PetscObject)dm)->comm;
115552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
116552f7358SJed Brown   PetscViewerVTKObjectLink link;
117552f7358SJed Brown   FILE                     *fp;
118552f7358SJed Brown   PetscMPIInt              rank,size,tag;
119552f7358SJed Brown   PetscErrorCode ierr;
120552f7358SJed Brown   PetscInt                 dim,cellHeight,cStart,cEnd,vStart,vEnd,cMax,numLabelCells,hasLabel,c,v,r,i;
121552f7358SJed Brown   PieceInfo                piece,*gpiece = PETSC_NULL;
122552f7358SJed Brown   void                     *buffer = PETSC_NULL;
123552f7358SJed Brown 
124552f7358SJed Brown   PetscFunctionBegin;
125552f7358SJed Brown #if defined(PETSC_USE_COMPLEX)
126552f7358SJed Brown   SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Complex values not supported");
127552f7358SJed Brown #endif
128552f7358SJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
129552f7358SJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
130552f7358SJed Brown   tag = ((PetscObject)viewer)->tag;
131552f7358SJed Brown 
132552f7358SJed Brown   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
133552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
134*519f805aSKarl Rupp #if defined(PETSC_WORDS_BIGENDIAN)
135552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr);
136552f7358SJed Brown #else
137552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr);
138552f7358SJed Brown #endif
139552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <UnstructuredGrid>\n");CHKERRQ(ierr);
140552f7358SJed Brown 
141552f7358SJed Brown   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
142552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
143552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
144552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
145770b213bSMatthew G Knepley   ierr = DMPlexGetHybridBounds(dm, &cMax, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
146552f7358SJed Brown   if (cMax >= 0) {cEnd = PetscMin(cEnd, cMax);}
147552f7358SJed Brown   ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
148552f7358SJed Brown   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
149552f7358SJed Brown   piece.nvertices = vEnd - vStart;
150552f7358SJed Brown   piece.ncells = 0;
151552f7358SJed Brown   piece.nconn = 0;
152552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
153552f7358SJed Brown     PetscInt *closure = PETSC_NULL;
154552f7358SJed Brown     PetscInt  closureSize;
155552f7358SJed Brown 
156552f7358SJed Brown     if (hasLabel) {
157552f7358SJed Brown       PetscInt value;
158552f7358SJed Brown 
159552f7358SJed Brown       ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
160552f7358SJed Brown       if (value != 1) continue;
161552f7358SJed Brown     }
162552f7358SJed Brown     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
163552f7358SJed Brown     for (v = 0; v < closureSize*2; v += 2) {
164552f7358SJed Brown       if ((closure[v] >= vStart) && (closure[v] < vEnd)) piece.nconn++;
165552f7358SJed Brown     }
166552f7358SJed Brown     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
167552f7358SJed Brown     piece.ncells++;
168552f7358SJed Brown   }
169552f7358SJed Brown   if (!rank) {ierr = PetscMalloc(size*sizeof(piece),&gpiece);CHKERRQ(ierr);}
170552f7358SJed Brown   ierr = MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm);CHKERRQ(ierr);
171552f7358SJed Brown 
172552f7358SJed Brown   /*
173552f7358SJed Brown    * Write file header
174552f7358SJed Brown    */
175552f7358SJed Brown   if (!rank) {
176552f7358SJed Brown     PetscInt boffset = 0;
177552f7358SJed Brown 
178552f7358SJed Brown     for (r=0; r<size; r++) {
179552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells);CHKERRQ(ierr);
180552f7358SJed Brown       /* Coordinate positions */
181552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Points>\n");CHKERRQ(ierr);
182552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
183552f7358SJed Brown       boffset += gpiece[r].nvertices*3*sizeof(PetscScalar) + sizeof(int);
184552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Points>\n");CHKERRQ(ierr);
185552f7358SJed Brown       /* Cell connectivity */
186552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Cells>\n");CHKERRQ(ierr);
187552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
188552f7358SJed Brown       boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int);
189552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"offsets\"      NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
190552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int);
191552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"UInt8\" Name=\"types\"        NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
192552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int);
193552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Cells>\n");CHKERRQ(ierr);
194552f7358SJed Brown 
195552f7358SJed Brown       /*
196552f7358SJed Brown        * Cell Data headers
197552f7358SJed Brown        */
198552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      <CellData>\n");CHKERRQ(ierr);
199552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
200552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(int) + sizeof(int);
201552f7358SJed Brown       /* all the vectors */
202552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
203552f7358SJed Brown         Vec X = (Vec)link->vec;
204552f7358SJed Brown         PetscInt bs;
205552f7358SJed Brown         const char *vecname = "";
206552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
207552f7358SJed Brown         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. */
208552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
209552f7358SJed Brown         }
210552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
211552f7358SJed Brown         for (i=0; i<bs; i++) {
212552f7358SJed Brown           char buf[256];
213552f7358SJed Brown           const char *fieldname = PETSC_NULL;
214552f7358SJed Brown           /* ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr); */
215552f7358SJed Brown           if (!fieldname) {
216552f7358SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"Unnamed%D",i);CHKERRQ(ierr);
217552f7358SJed Brown             fieldname = buf;
218552f7358SJed Brown           }
219552f7358SJed Brown           ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
220552f7358SJed Brown           boffset += gpiece[r].ncells*sizeof(PetscScalar) + sizeof(int);
221552f7358SJed Brown         }
222552f7358SJed Brown       }
223552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </CellData>\n");CHKERRQ(ierr);
224552f7358SJed Brown 
225552f7358SJed Brown       /*
226552f7358SJed Brown        * Point Data headers
227552f7358SJed Brown        */
228552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      <PointData>\n");CHKERRQ(ierr);
229552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
230552f7358SJed Brown         Vec X = (Vec)link->vec;
231552f7358SJed Brown         PetscInt bs;
232552f7358SJed Brown         const char *vecname = "";
233552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
234552f7358SJed Brown         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. */
235552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
236552f7358SJed Brown         }
237552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
238552f7358SJed Brown         for (i=0; i<bs; i++) {
239552f7358SJed Brown           char fieldname[256];
240552f7358SJed Brown           ierr = PetscSNPrintf(fieldname,sizeof(fieldname),"Point%D",i);CHKERRQ(ierr);
241552f7358SJed Brown           ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
242552f7358SJed Brown           boffset += gpiece[r].nvertices*sizeof(PetscScalar) + sizeof(int);
243552f7358SJed Brown         }
244552f7358SJed Brown       }
245552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </PointData>\n");CHKERRQ(ierr);
246552f7358SJed Brown 
247552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    </Piece>\n");CHKERRQ(ierr);
248552f7358SJed Brown     }
249552f7358SJed Brown   }
250552f7358SJed Brown 
251552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  </UnstructuredGrid>\n");CHKERRQ(ierr);
252552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
253552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
254552f7358SJed Brown 
255552f7358SJed Brown   if (!rank) {
256552f7358SJed Brown     PetscInt maxsize = 0;
257552f7358SJed Brown     for (r=0; r<size; r++) {
258552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscScalar)));
259552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*sizeof(PetscScalar)));
260552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt)));
261552f7358SJed Brown     }
262552f7358SJed Brown     ierr = PetscMalloc(maxsize,&buffer);CHKERRQ(ierr);
263552f7358SJed Brown   }
264552f7358SJed Brown   for (r=0; r<size; r++) {
265552f7358SJed Brown     if (r == rank) {
266552f7358SJed Brown       PetscInt nsend;
267552f7358SJed Brown       {                         /* Position */
268552f7358SJed Brown         const PetscScalar *x;
269552f7358SJed Brown         PetscScalar *y = PETSC_NULL;
270552f7358SJed Brown         Vec coords;
271552f7358SJed Brown         nsend = piece.nvertices*3;
272552f7358SJed Brown         ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
273552f7358SJed Brown         ierr = VecGetArrayRead(coords,&x);CHKERRQ(ierr);
274552f7358SJed Brown         if (dim != 3) {
275552f7358SJed Brown           ierr = PetscMalloc(piece.nvertices*3*sizeof(PetscScalar),&y);CHKERRQ(ierr);
276552f7358SJed Brown           for (i=0; i<piece.nvertices; i++) {
277552f7358SJed Brown             y[i*3+0] = x[i*dim+0];
278552f7358SJed Brown             y[i*3+1] = (dim > 1) ? x[i*dim+1] : 0;
279552f7358SJed Brown             y[i*3+2] = 0;
280552f7358SJed Brown           }
281552f7358SJed Brown         }
282552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,y?y:x,buffer,nsend,PETSC_SCALAR,tag);CHKERRQ(ierr);
283552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
284552f7358SJed Brown         ierr = VecRestoreArrayRead(coords,&x);CHKERRQ(ierr);
285552f7358SJed Brown       }
286552f7358SJed Brown       {                           /* Connectivity, offsets, types */
2876a415b8fSMatthew G Knepley         PetscVTKInt *connectivity = PETSC_NULL,*offsets;
288552f7358SJed Brown         PetscVTKType *types;
289552f7358SJed Brown         ierr = DMPlexGetVTKConnectivity(dm,&piece,&connectivity,&offsets,&types);CHKERRQ(ierr);
290552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,connectivity,buffer,piece.nconn,PETSC_INT32,tag);CHKERRQ(ierr);
291552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,offsets,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr);
292552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,types,buffer,piece.ncells,PETSC_UINT8,tag);CHKERRQ(ierr);
293552f7358SJed Brown         ierr = PetscFree3(connectivity,offsets,types);CHKERRQ(ierr);
294552f7358SJed Brown       }
295552f7358SJed Brown       {                         /* Owners (cell data) */
296552f7358SJed Brown         PetscVTKInt *owners;
297552f7358SJed Brown         ierr = PetscMalloc(piece.ncells*sizeof(PetscVTKInt),&owners);CHKERRQ(ierr);
298552f7358SJed Brown         for (i=0; i<piece.ncells; i++) owners[i] = rank;
299552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,owners,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr);
300552f7358SJed Brown         ierr = PetscFree(owners);CHKERRQ(ierr);
301552f7358SJed Brown       }
302552f7358SJed Brown                                 /* Cell data */
303552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
304552f7358SJed Brown         Vec X = (Vec)link->vec;
305552f7358SJed Brown         const PetscScalar *x;
306552f7358SJed Brown         PetscScalar *y;
307552f7358SJed Brown         PetscInt bs;
308552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
309552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
310552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
311552f7358SJed Brown         ierr = PetscMalloc(piece.ncells*sizeof(PetscScalar),&y);CHKERRQ(ierr);
312552f7358SJed Brown         for (i=0; i<bs; i++) {
313552f7358SJed Brown           PetscInt cnt;
314552f7358SJed Brown           for (c=cStart,cnt=0; c<cEnd; c++) {
315552f7358SJed Brown             const PetscScalar *xpoint;
316552f7358SJed Brown             if (hasLabel) {     /* Ignore some cells */
317552f7358SJed Brown               PetscInt value;
318552f7358SJed Brown               ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
319552f7358SJed Brown               if (value != 1) continue;
320552f7358SJed Brown             }
321552f7358SJed Brown             ierr = DMPlexPointLocalRead(dm,c,x,&xpoint);CHKERRQ(ierr);
322552f7358SJed Brown             y[cnt++] = xpoint[i];
323552f7358SJed Brown           }
324552f7358SJed Brown           if (cnt != piece.ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
325552f7358SJed Brown           ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells,PETSC_SCALAR,tag);CHKERRQ(ierr);
326552f7358SJed Brown         }
327552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
328552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
329552f7358SJed Brown       }
330552f7358SJed Brown 
331552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
332552f7358SJed Brown         Vec X = (Vec)link->vec;
333552f7358SJed Brown         const PetscScalar *x;
334552f7358SJed Brown         PetscScalar *y;
335552f7358SJed Brown         PetscInt bs;
336552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
337552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
338552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
339552f7358SJed Brown         ierr = PetscMalloc(piece.nvertices*sizeof(PetscScalar),&y);CHKERRQ(ierr);
340552f7358SJed Brown         for (i=0; i<bs; i++) {
341552f7358SJed Brown           PetscInt cnt;
342552f7358SJed Brown           for (v=vStart,cnt=0; v<vEnd; v++) {
343552f7358SJed Brown             const PetscScalar *xpoint;
344552f7358SJed Brown             ierr = DMPlexPointLocalRead(dm,c,x,&xpoint);CHKERRQ(ierr);
345552f7358SJed Brown             y[cnt++] = xpoint[i];
346552f7358SJed Brown           }
347552f7358SJed Brown           if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
348552f7358SJed Brown           ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr);
349552f7358SJed Brown         }
350552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
351552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
352552f7358SJed Brown       }
353552f7358SJed Brown     } else if (!rank) {
354552f7358SJed Brown       ierr = TransferWrite(viewer,fp,r,0,PETSC_NULL,buffer,gpiece[r].nvertices*3,PETSC_SCALAR,tag);CHKERRQ(ierr); /* positions */
355552f7358SJed Brown       ierr = TransferWrite(viewer,fp,r,0,PETSC_NULL,buffer,gpiece[r].nconn,PETSC_INT32,tag);CHKERRQ(ierr); /* connectivity */
356552f7358SJed Brown       ierr = TransferWrite(viewer,fp,r,0,PETSC_NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* offsets */
357552f7358SJed Brown       ierr = TransferWrite(viewer,fp,r,0,PETSC_NULL,buffer,gpiece[r].ncells,PETSC_UINT8,tag);CHKERRQ(ierr); /* types */
358552f7358SJed Brown       ierr = TransferWrite(viewer,fp,r,0,PETSC_NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* owner rank (cells) */
359552f7358SJed Brown       /* all cell data */
360552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
361552f7358SJed Brown         PetscInt bs;
362552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
363552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
364552f7358SJed Brown         for (i=0; i<bs; i++) {
365552f7358SJed Brown           ierr = TransferWrite(viewer,fp,r,0,PETSC_NULL,buffer,gpiece[r].ncells,PETSC_SCALAR,tag);CHKERRQ(ierr);
366552f7358SJed Brown         }
367552f7358SJed Brown       }
368552f7358SJed Brown       /* all point data */
369552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
370552f7358SJed Brown         PetscInt bs;
371552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
372552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
373552f7358SJed Brown         for (i=0; i<bs; i++) {
374552f7358SJed Brown           ierr = TransferWrite(viewer,fp,r,0,PETSC_NULL,buffer,gpiece[r].nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr);
375552f7358SJed Brown         }
376552f7358SJed Brown       }
377552f7358SJed Brown     }
378552f7358SJed Brown   }
379552f7358SJed Brown   ierr = PetscFree(gpiece);CHKERRQ(ierr);
380552f7358SJed Brown   ierr = PetscFree(buffer);CHKERRQ(ierr);
381552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"\n  </AppendedData>\n");CHKERRQ(ierr);
382552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
383552f7358SJed Brown   PetscFunctionReturn(0);
384552f7358SJed Brown }
385