xref: /petsc/src/dm/impls/plex/plexvtu.c (revision 955d60d1063793c51bfe60f1f11b88d482c00bf3)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.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 
10*955d60d1SToby Isaac #if defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL___FP16)
11*955d60d1SToby Isaac /* output in float if single or half precision in memory */
12552f7358SJed Brown static const char precision[] = "Float32";
13*955d60d1SToby Isaac typedef float PetscVTUReal;
14*955d60d1SToby Isaac #define MPIU_VTUREAL MPI_FLOAT
15*955d60d1SToby Isaac #elif defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
16*955d60d1SToby Isaac /* output in double if double or quad precision in memory */
17552f7358SJed Brown static const char precision[] = "Float64";
18*955d60d1SToby Isaac typedef double PetscVTUReal;
19*955d60d1SToby Isaac #define MPIU_VTUREAL MPI_DOUBLE
20552f7358SJed Brown #else
21552f7358SJed Brown static const char precision[] = "UnknownPrecision";
22*955d60d1SToby Isaac typedef PetscReal PetscVTUReal;
23*955d60d1SToby Isaac #define MPIU_VTUREAL MPIU_REAL
24552f7358SJed Brown #endif
25552f7358SJed Brown 
2694fbd55eSBarry Smith static PetscErrorCode TransferWrite(PetscViewer viewer,FILE *fp,PetscMPIInt srank,PetscMPIInt root,const void *send,void *recv,PetscMPIInt count,MPI_Datatype mpidatatype,PetscMPIInt tag)
27552f7358SJed Brown {
28552f7358SJed Brown   PetscMPIInt    rank;
29552f7358SJed Brown   PetscErrorCode ierr;
30ce94432eSBarry Smith   MPI_Comm       comm;
31552f7358SJed Brown 
32552f7358SJed Brown   PetscFunctionBegin;
33ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
34552f7358SJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
35552f7358SJed Brown 
36552f7358SJed Brown   if (rank == srank && rank != root) {
37552f7358SJed Brown     ierr = MPI_Send((void*)send,count,mpidatatype,root,tag,comm);CHKERRQ(ierr);
38552f7358SJed Brown   } else if (rank == root) {
39552f7358SJed Brown     const void *buffer;
40552f7358SJed Brown     if (root == srank) {        /* self */
41552f7358SJed Brown       buffer = send;
42552f7358SJed Brown     } else {
43552f7358SJed Brown       MPI_Status  status;
44552f7358SJed Brown       PetscMPIInt nrecv;
45552f7358SJed Brown       ierr = MPI_Recv(recv,count,mpidatatype,srank,tag,comm,&status);CHKERRQ(ierr);
46552f7358SJed Brown       ierr = MPI_Get_count(&status,mpidatatype,&nrecv);CHKERRQ(ierr);
47552f7358SJed Brown       if (count != nrecv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
48552f7358SJed Brown       buffer = recv;
49552f7358SJed Brown     }
5094fbd55eSBarry Smith     ierr = PetscViewerVTKFWrite(viewer,fp,buffer,count,mpidatatype);CHKERRQ(ierr);
51552f7358SJed Brown   }
52552f7358SJed Brown   PetscFunctionReturn(0);
53552f7358SJed Brown }
54552f7358SJed Brown 
552e529ec8SStefano Zampini static PetscErrorCode DMPlexGetVTKConnectivity(DM dm, PetscBool localized, PieceInfo *piece,PetscVTKInt **oconn,PetscVTKInt **ooffsets,PetscVTKType **otypes)
56552f7358SJed Brown {
57552f7358SJed Brown   PetscErrorCode ierr;
582e529ec8SStefano Zampini   PetscSection   coordSection;
59552f7358SJed Brown   PetscVTKInt    *conn,*offsets;
60552f7358SJed Brown   PetscVTKType   *types;
612f029394SStefano Zampini   PetscInt       dim,vStart,vEnd,cStart,cEnd,pStart,pEnd,cellHeight,numLabelCells,hasLabel,c,v,countcell,countconn;
62552f7358SJed Brown 
63552f7358SJed Brown   PetscFunctionBegin;
64dcca6d9dSJed Brown   ierr = PetscMalloc3(piece->nconn,&conn,piece->ncells,&offsets,piece->ncells,&types);CHKERRQ(ierr);
652e529ec8SStefano Zampini   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
66c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
67552f7358SJed Brown   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
68552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
69552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
70552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
71c58f1c22SToby Isaac   ierr     = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
72552f7358SJed Brown   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
73552f7358SJed Brown 
74552f7358SJed Brown   countcell = 0;
75552f7358SJed Brown   countconn = 0;
76552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
7719003fb5SStefano Zampini     PetscInt nverts,dof = 0,celltype,startoffset,nC=0;
78552f7358SJed Brown 
79552f7358SJed Brown     if (hasLabel) {
80552f7358SJed Brown       PetscInt value;
81552f7358SJed Brown 
82c58f1c22SToby Isaac       ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
83552f7358SJed Brown       if (value != 1) continue;
84552f7358SJed Brown     }
85552f7358SJed Brown     startoffset = countconn;
8619003fb5SStefano Zampini     if (localized) {
8719003fb5SStefano Zampini       ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
8819003fb5SStefano Zampini     }
8919003fb5SStefano Zampini     if (!dof) {
902e529ec8SStefano Zampini       PetscInt *closure = NULL;
912e529ec8SStefano Zampini       PetscInt  closureSize;
922e529ec8SStefano Zampini 
93552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
94552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
95552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
9619003fb5SStefano Zampini           if (!localized) conn[countconn++] = closure[v] - vStart;
9719003fb5SStefano Zampini           else conn[countconn++] = startoffset + nC;
98724b5320SMatthew G. Knepley           ++nC;
99552f7358SJed Brown         }
100552f7358SJed Brown       }
101552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1022e529ec8SStefano Zampini     } else {
1032e529ec8SStefano Zampini       for (nC = 0; nC < dof/dim; nC++) conn[countconn++] = startoffset + nC;
1042e529ec8SStefano Zampini     }
105724b5320SMatthew G. Knepley     ierr = DMPlexInvertCell(dim, nC, &conn[countconn-nC]);CHKERRQ(ierr);
1068865f1eaSKarl Rupp 
107552f7358SJed Brown     offsets[countcell] = countconn;
1088865f1eaSKarl Rupp 
109552f7358SJed Brown     nverts = countconn - startoffset;
110de0a4605SMatthew G. Knepley     ierr   = DMPlexVTKGetCellType_Internal(dm,dim,nverts,&celltype);CHKERRQ(ierr);
1118865f1eaSKarl Rupp 
112552f7358SJed Brown     types[countcell] = celltype;
113552f7358SJed Brown     countcell++;
114552f7358SJed Brown   }
115552f7358SJed Brown   if (countcell != piece->ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count");
116552f7358SJed Brown   if (countconn != piece->nconn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count");
117552f7358SJed Brown   *oconn    = conn;
118552f7358SJed Brown   *ooffsets = offsets;
119552f7358SJed Brown   *otypes   = types;
120552f7358SJed Brown   PetscFunctionReturn(0);
121552f7358SJed Brown }
122552f7358SJed Brown 
123552f7358SJed Brown /*
124552f7358SJed Brown   Write all fields that have been provided to the viewer
125552f7358SJed Brown   Multi-block XML format with binary appended data.
126552f7358SJed Brown */
127552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer)
128552f7358SJed Brown {
129ce94432eSBarry Smith   MPI_Comm                 comm;
1302e529ec8SStefano Zampini   PetscSection             coordSection;
131552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
132552f7358SJed Brown   PetscViewerVTKObjectLink link;
133552f7358SJed Brown   FILE                     *fp;
134552f7358SJed Brown   PetscMPIInt              rank,size,tag;
135552f7358SJed Brown   PetscErrorCode           ierr;
1362f029394SStefano Zampini   PetscInt                 dimEmbed,cellHeight,cStart,cEnd,vStart,vEnd,numLabelCells,hasLabel,c,v,r,i;
1372e529ec8SStefano Zampini   PetscBool                localized;
1380298fd71SBarry Smith   PieceInfo                piece,*gpiece = NULL;
1390298fd71SBarry Smith   void                     *buffer = NULL;
14030815ce0SLisandro Dalcin   const char               *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian";
141*955d60d1SToby Isaac   PetscInt                 loops_per_scalar;
142552f7358SJed Brown 
143552f7358SJed Brown   PetscFunctionBegin;
144ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
145552f7358SJed Brown #if defined(PETSC_USE_COMPLEX)
146*955d60d1SToby Isaac   loops_per_scalar = 2;
147*955d60d1SToby Isaac #else
148*955d60d1SToby Isaac   loops_per_scalar = 1;
149552f7358SJed Brown #endif
150552f7358SJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
151552f7358SJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
152552f7358SJed Brown   tag  = ((PetscObject)viewer)->tag;
153552f7358SJed Brown 
154552f7358SJed Brown   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
155552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
15630815ce0SLisandro Dalcin   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"%s\">\n", byte_order);CHKERRQ(ierr);
157552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <UnstructuredGrid>\n");CHKERRQ(ierr);
158552f7358SJed Brown 
1593baa072aSToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
160552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
161552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
162552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
163c58f1c22SToby Isaac   ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
1642e529ec8SStefano Zampini   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
1652e529ec8SStefano Zampini   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1668865f1eaSKarl Rupp 
167552f7358SJed Brown   hasLabel        = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
1682e529ec8SStefano Zampini   piece.nvertices = 0;
169552f7358SJed Brown   piece.ncells    = 0;
170552f7358SJed Brown   piece.nconn     = 0;
1712e529ec8SStefano Zampini   if (!localized) piece.nvertices = vEnd - vStart;
172552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
17319003fb5SStefano Zampini     PetscInt dof = 0;
174552f7358SJed Brown     if (hasLabel) {
175552f7358SJed Brown       PetscInt value;
176552f7358SJed Brown 
177c58f1c22SToby Isaac       ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
178552f7358SJed Brown       if (value != 1) continue;
179552f7358SJed Brown     }
18019003fb5SStefano Zampini     if (localized) {
18119003fb5SStefano Zampini       ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
18219003fb5SStefano Zampini     }
18319003fb5SStefano Zampini     if (!dof) {
1842e529ec8SStefano Zampini       PetscInt *closure = NULL;
1852e529ec8SStefano Zampini       PetscInt closureSize;
1862e529ec8SStefano Zampini 
187552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
188552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
1892e529ec8SStefano Zampini         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
1902e529ec8SStefano Zampini           piece.nconn++;
19119003fb5SStefano Zampini           if (localized) piece.nvertices++;
1922e529ec8SStefano Zampini         }
193552f7358SJed Brown       }
194552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1952e529ec8SStefano Zampini     } else {
1962e529ec8SStefano Zampini       piece.nvertices += dof/dimEmbed;
1972e529ec8SStefano Zampini       piece.nconn     += dof/dimEmbed;
1982e529ec8SStefano Zampini     }
199552f7358SJed Brown     piece.ncells++;
200552f7358SJed Brown   }
201785e854fSJed Brown   if (!rank) {ierr = PetscMalloc1(size,&gpiece);CHKERRQ(ierr);}
202552f7358SJed Brown   ierr = MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm);CHKERRQ(ierr);
203552f7358SJed Brown 
204552f7358SJed Brown   /*
205552f7358SJed Brown    * Write file header
206552f7358SJed Brown    */
207552f7358SJed Brown   if (!rank) {
208552f7358SJed Brown     PetscInt boffset = 0;
209552f7358SJed Brown 
210552f7358SJed Brown     for (r=0; r<size; r++) {
211552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells);CHKERRQ(ierr);
212552f7358SJed Brown       /* Coordinate positions */
213552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Points>\n");CHKERRQ(ierr);
214552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
215*955d60d1SToby Isaac       boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int);
216552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Points>\n");CHKERRQ(ierr);
217552f7358SJed Brown       /* Cell connectivity */
218552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Cells>\n");CHKERRQ(ierr);
219552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
220552f7358SJed Brown       boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int);
221552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"offsets\"      NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
222552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int);
223552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"UInt8\" Name=\"types\"        NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
224552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int);
225552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Cells>\n");CHKERRQ(ierr);
226552f7358SJed Brown 
227552f7358SJed Brown       /*
228552f7358SJed Brown        * Cell Data headers
229552f7358SJed Brown        */
230552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <CellData>\n");CHKERRQ(ierr);
231552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
232552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(int) + sizeof(int);
233552f7358SJed Brown       /* all the vectors */
234552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
235552f7358SJed Brown         Vec        X = (Vec)link->vec;
236e630c359SToby Isaac         DM         dmX = NULL;
2371cfafdd3SJed Brown         PetscInt   bs,nfields,field;
238552f7358SJed Brown         const char *vecname = "";
239e630c359SToby Isaac         PetscSection section;
240552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
241552f7358SJed 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. */
242552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
243552f7358SJed Brown         }
244e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
245e630c359SToby Isaac         if (!dmX) dmX = dm;
246e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
247e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
248e630c359SToby Isaac         ierr = PetscSectionGetDof(section,cStart,&bs);CHKERRQ(ierr);
249e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
250e630c359SToby Isaac         field = 0;
251e630c359SToby Isaac         if (link->field >= 0) {
252e630c359SToby Isaac           field = link->field;
253e630c359SToby Isaac           nfields = field + 1;
254e630c359SToby Isaac         }
255e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
2561cfafdd3SJed Brown           PetscInt     fbs,j;
257a00cdb45SToby Isaac           PetscFV      fv = NULL;
258a00cdb45SToby Isaac           PetscObject  f;
259a00cdb45SToby Isaac           PetscClassId fClass;
2600298fd71SBarry Smith           const char *fieldname = NULL;
2611cfafdd3SJed Brown           char       buf[256];
262e630c359SToby Isaac           PetscBool    vector;
2637ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
264e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,cStart,field,&fbs);CHKERRQ(ierr);
265e630c359SToby Isaac             ierr = PetscSectionGetFieldName(section,field,&fieldname);CHKERRQ(ierr);
2667ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
267e630c359SToby Isaac           ierr = DMGetField(dmX,field,NULL,&f);CHKERRQ(ierr);
268a00cdb45SToby Isaac           ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr);
269a00cdb45SToby Isaac           if (fClass == PETSCFV_CLASSID) {
270a00cdb45SToby Isaac             fv = (PetscFV) f;
271a00cdb45SToby Isaac           }
272e630c359SToby Isaac           if (nfields && !fieldname) {
2731cfafdd3SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"CellField%D",field);CHKERRQ(ierr);
274552f7358SJed Brown             fieldname = buf;
275552f7358SJed Brown           }
276e630c359SToby Isaac           vector = PETSC_FALSE;
277e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
278e630c359SToby Isaac             vector = PETSC_TRUE;
279e630c359SToby Isaac             if (fbs > 3) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_SIZ,"Cell vector fields can have at most 3 components, %D given\n", fbs);
280e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
281e630c359SToby Isaac               const char *compName = NULL;
282e630c359SToby Isaac               if (fv) {
283e630c359SToby Isaac                 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
284e630c359SToby Isaac                 if (compName) break;
285e630c359SToby Isaac               }
286e630c359SToby Isaac             }
287e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
288e630c359SToby Isaac           }
289e630c359SToby Isaac           if (vector) {
290*955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
291*955d60d1SToby Isaac             ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s.Re\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
292*955d60d1SToby Isaac             boffset += gpiece[r].ncells*3*sizeof(PetscVTUReal) + sizeof(int);
293*955d60d1SToby Isaac             ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s.Im\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
294*955d60d1SToby Isaac             boffset += gpiece[r].ncells*3*sizeof(PetscVTUReal) + sizeof(int);
295*955d60d1SToby Isaac #else
296e630c359SToby Isaac             ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
297*955d60d1SToby Isaac             boffset += gpiece[r].ncells*3*sizeof(PetscVTUReal) + sizeof(int);
298*955d60d1SToby Isaac #endif
299e630c359SToby Isaac             i+=fbs;
300e630c359SToby Isaac           } else {
3011cfafdd3SJed Brown             for (j=0; j<fbs; j++) {
302a00cdb45SToby Isaac               const char *compName = NULL;
303*955d60d1SToby Isaac               char finalname[256];
304a00cdb45SToby Isaac               if (fv) {
305a00cdb45SToby Isaac                 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
306a00cdb45SToby Isaac               }
307a00cdb45SToby Isaac               if (compName) {
308*955d60d1SToby Isaac                 ierr = PetscSNPrintf(finalname,255,"%s%s.%s",vecname,fieldname,compName);CHKERRQ(ierr);
309a00cdb45SToby Isaac               }
310e630c359SToby Isaac               else if (fbs > 1) {
311*955d60d1SToby Isaac                 ierr = PetscSNPrintf(finalname,255,"%s%s.%D",vecname,fieldname,j);CHKERRQ(ierr);
312e630c359SToby Isaac               } else {
313*955d60d1SToby Isaac                 ierr = PetscSNPrintf(finalname,255,"%s%s",vecname,fieldname);CHKERRQ(ierr);
314a00cdb45SToby Isaac               }
315*955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
316*955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
317*955d60d1SToby Isaac               boffset += gpiece[r].ncells*sizeof(PetscVTUReal) + sizeof(int);
318*955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
319*955d60d1SToby Isaac               boffset += gpiece[r].ncells*sizeof(PetscVTUReal) + sizeof(int);
320*955d60d1SToby Isaac #else
321*955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
322*955d60d1SToby Isaac               boffset += gpiece[r].ncells*sizeof(PetscVTUReal) + sizeof(int);
323*955d60d1SToby Isaac #endif
3241cfafdd3SJed Brown               i++;
325552f7358SJed Brown             }
326552f7358SJed Brown           }
327e630c359SToby Isaac         }
3281cfafdd3SJed Brown         if (i != bs) SETERRQ2(comm,PETSC_ERR_PLIB,"Total number of field components %D != block size %D",i,bs);
3291cfafdd3SJed Brown       }
330552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </CellData>\n");CHKERRQ(ierr);
331552f7358SJed Brown 
332552f7358SJed Brown       /*
333552f7358SJed Brown        * Point Data headers
334552f7358SJed Brown        */
335552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      <PointData>\n");CHKERRQ(ierr);
336552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
337552f7358SJed Brown         Vec        X = (Vec)link->vec;
338e630c359SToby Isaac         DM         dmX;
3391cfafdd3SJed Brown         PetscInt   bs,nfields,field;
340552f7358SJed Brown         const char *vecname = "";
341e630c359SToby Isaac         PetscSection section;
342552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
343552f7358SJed 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. */
344552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
345552f7358SJed Brown         }
346e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
347e630c359SToby Isaac         if (!dmX) dmX = dm;
348e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
349e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
350e630c359SToby Isaac         ierr = PetscSectionGetDof(section,vStart,&bs);CHKERRQ(ierr);
351e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
352e630c359SToby Isaac         field = 0;
353e630c359SToby Isaac         if (link->field >= 0) {
354e630c359SToby Isaac           field = link->field;
355e630c359SToby Isaac           nfields = field + 1;
356e630c359SToby Isaac         }
357e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
3581cfafdd3SJed Brown           PetscInt   fbs,j;
3591cfafdd3SJed Brown           const char *fieldname = NULL;
3601cfafdd3SJed Brown           char       buf[256];
3617ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
362e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,vStart,field,&fbs);CHKERRQ(ierr);
363e630c359SToby Isaac             ierr = PetscSectionGetFieldName(section,field,&fieldname);CHKERRQ(ierr);
3647ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
365e630c359SToby Isaac           if (nfields && !fieldname) {
3661cfafdd3SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"PointField%D",field);CHKERRQ(ierr);
3671cfafdd3SJed Brown             fieldname = buf;
3681cfafdd3SJed Brown           }
369e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
370e630c359SToby Isaac             if (fbs > 3) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_SIZ,"Point vector fields can have at most 3 components, %D given\n", fbs);
371*955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
372*955d60d1SToby Isaac             ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s.Re\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
373*955d60d1SToby Isaac             boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int);
374*955d60d1SToby Isaac             ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s.Im\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
375*955d60d1SToby Isaac             boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int);
376*955d60d1SToby Isaac #else
377e630c359SToby Isaac             ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
378*955d60d1SToby Isaac             boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int);
379*955d60d1SToby Isaac #endif
380e630c359SToby Isaac           } else {
3811cfafdd3SJed Brown             for (j=0; j<fbs; j++) {
382*955d60d1SToby Isaac               char finalname[256];
383e630c359SToby Isaac               if (fbs > 1) {
384*955d60d1SToby Isaac                 ierr = PetscSNPrintf(finalname,255,"%s%s.%D",vecname,fieldname,j);CHKERRQ(ierr);
385e630c359SToby Isaac               } else {
386*955d60d1SToby Isaac                 ierr = PetscSNPrintf(finalname,255,"%s%s",vecname,fieldname);CHKERRQ(ierr);
387e630c359SToby Isaac               }
388*955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
389*955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
390*955d60d1SToby Isaac               boffset += gpiece[r].nvertices*sizeof(PetscVTUReal) + sizeof(int);
391*955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
392*955d60d1SToby Isaac               boffset += gpiece[r].nvertices*sizeof(PetscVTUReal) + sizeof(int);
393*955d60d1SToby Isaac #else
394*955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
395*955d60d1SToby Isaac               boffset += gpiece[r].nvertices*sizeof(PetscVTUReal) + sizeof(int);
396*955d60d1SToby Isaac #endif
397552f7358SJed Brown             }
398552f7358SJed Brown           }
3991cfafdd3SJed Brown         }
400e630c359SToby Isaac       }
401552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </PointData>\n");CHKERRQ(ierr);
402552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    </Piece>\n");CHKERRQ(ierr);
403552f7358SJed Brown     }
404552f7358SJed Brown   }
405552f7358SJed Brown 
406552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  </UnstructuredGrid>\n");CHKERRQ(ierr);
407552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
408552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
409552f7358SJed Brown 
410552f7358SJed Brown   if (!rank) {
411552f7358SJed Brown     PetscInt maxsize = 0;
412552f7358SJed Brown     for (r=0; r<size; r++) {
413*955d60d1SToby Isaac       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscVTUReal)));
414*955d60d1SToby Isaac       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*3*sizeof(PetscVTUReal)));
415552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt)));
416552f7358SJed Brown     }
417552f7358SJed Brown     ierr = PetscMalloc(maxsize,&buffer);CHKERRQ(ierr);
418552f7358SJed Brown   }
419552f7358SJed Brown   for (r=0; r<size; r++) {
420552f7358SJed Brown     if (r == rank) {
421552f7358SJed Brown       PetscInt nsend;
422552f7358SJed Brown       {                         /* Position */
423552f7358SJed Brown         const PetscScalar *x;
424*955d60d1SToby Isaac         PetscVTUReal      *y = NULL;
425552f7358SJed Brown         Vec               coords;
426*955d60d1SToby Isaac         PetscBool         copy;
4272e529ec8SStefano Zampini 
428552f7358SJed Brown         ierr  = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
429552f7358SJed Brown         ierr  = VecGetArrayRead(coords,&x);CHKERRQ(ierr);
430*955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
431*955d60d1SToby Isaac         copy = PETSC_TRUE;
432*955d60d1SToby Isaac #else
433*955d60d1SToby Isaac         copy = (PetscBool) (dimEmbed != 3 || localized || (sizeof(PetscReal) != sizeof(PetscVTUReal)));
434*955d60d1SToby Isaac #endif
435*955d60d1SToby Isaac         if (copy) {
436785e854fSJed Brown           ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr);
43719003fb5SStefano Zampini           if (localized) {
43819003fb5SStefano Zampini             PetscInt cnt;
43919003fb5SStefano Zampini             for (c=cStart,cnt=0; c<cEnd; c++) {
44019003fb5SStefano Zampini               PetscInt off, dof;
44119003fb5SStefano Zampini 
44219003fb5SStefano Zampini               ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
44319003fb5SStefano Zampini               if (!dof) {
44419003fb5SStefano Zampini                 PetscInt *closure = NULL;
44519003fb5SStefano Zampini                 PetscInt closureSize;
44619003fb5SStefano Zampini 
44719003fb5SStefano Zampini                 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
44819003fb5SStefano Zampini                 for (v = 0; v < closureSize*2; v += 2) {
44919003fb5SStefano Zampini                   if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
45019003fb5SStefano Zampini                     ierr = PetscSectionGetOffset(coordSection, closure[v], &off);CHKERRQ(ierr);
45119003fb5SStefano Zampini                     if (dimEmbed != 3) {
452*955d60d1SToby Isaac                       y[cnt*3+0] = (PetscVTUReal) PetscRealPart(x[off+0]);
453*955d60d1SToby Isaac                       y[cnt*3+1] = (PetscVTUReal) ((dimEmbed > 1) ? PetscRealPart(x[off+1]) : 0.0);
454*955d60d1SToby Isaac                       y[cnt*3+2] = (PetscVTUReal) 0.0;
45519003fb5SStefano Zampini                     } else {
456*955d60d1SToby Isaac                       y[cnt*3+0] = (PetscVTUReal) PetscRealPart(x[off+0]);
457*955d60d1SToby Isaac                       y[cnt*3+1] = (PetscVTUReal) PetscRealPart(x[off+1]);
458*955d60d1SToby Isaac                       y[cnt*3+2] = (PetscVTUReal) PetscRealPart(x[off+2]);
45919003fb5SStefano Zampini                     }
46019003fb5SStefano Zampini                     cnt++;
46119003fb5SStefano Zampini                   }
46219003fb5SStefano Zampini                 }
46319003fb5SStefano Zampini                 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
46419003fb5SStefano Zampini               } else {
46519003fb5SStefano Zampini                 ierr = PetscSectionGetOffset(coordSection, c, &off);CHKERRQ(ierr);
46619003fb5SStefano Zampini                 if (dimEmbed != 3) {
46719003fb5SStefano Zampini                   for (i=0; i<dof/dimEmbed; i++) {
468*955d60d1SToby Isaac                     y[cnt*3+0] = (PetscVTUReal) PetscRealPart(x[off + i*dimEmbed + 0]);
469*955d60d1SToby Isaac                     y[cnt*3+1] = (PetscVTUReal) ((dimEmbed > 1) ? PetscRealPart(x[off + i*dimEmbed + 1]) : 0.0);
470*955d60d1SToby Isaac                     y[cnt*3+2] = (PetscVTUReal) 0.0;
47119003fb5SStefano Zampini                     cnt++;
47219003fb5SStefano Zampini                   }
47319003fb5SStefano Zampini                 } else {
47419003fb5SStefano Zampini                   for (i=0; i<dof; i ++) {
475*955d60d1SToby Isaac                     y[cnt*3+i] = (PetscVTUReal) PetscRealPart(x[off + i]);
47619003fb5SStefano Zampini                   }
47719003fb5SStefano Zampini                   cnt += dof/dimEmbed;
47819003fb5SStefano Zampini                 }
47919003fb5SStefano Zampini               }
48019003fb5SStefano Zampini             }
48119003fb5SStefano Zampini             if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
48219003fb5SStefano Zampini           } else {
483552f7358SJed Brown             for (i=0; i<piece.nvertices; i++) {
484*955d60d1SToby Isaac               y[i*3+0] = (PetscVTUReal) PetscRealPart(x[i*dimEmbed+0]);
485*955d60d1SToby Isaac               y[i*3+1] = (PetscVTUReal) ((dimEmbed > 1) ? PetscRealPart(x[i*dimEmbed+1]) : 0.);
486*955d60d1SToby Isaac               y[i*3+2] = (PetscVTUReal) ((dimEmbed > 2) ? PetscRealPart(x[i*dimEmbed+2]) : 0.);
48719003fb5SStefano Zampini             }
488552f7358SJed Brown           }
489552f7358SJed Brown         }
4902e529ec8SStefano Zampini         nsend = piece.nvertices*3;
491*955d60d1SToby Isaac         ierr  = TransferWrite(viewer,fp,r,0,copy ? (const void *) y : (const void *) x,buffer,nsend,MPIU_VTUREAL,tag);CHKERRQ(ierr);
492552f7358SJed Brown         ierr  = PetscFree(y);CHKERRQ(ierr);
493552f7358SJed Brown         ierr  = VecRestoreArrayRead(coords,&x);CHKERRQ(ierr);
494552f7358SJed Brown       }
495552f7358SJed Brown       {                           /* Connectivity, offsets, types */
4963e869605SMatthew G. Knepley         PetscVTKInt  *connectivity = NULL, *offsets = NULL;
4973e869605SMatthew G. Knepley         PetscVTKType *types = NULL;
4982e529ec8SStefano Zampini         ierr = DMPlexGetVTKConnectivity(dm,localized,&piece,&connectivity,&offsets,&types);CHKERRQ(ierr);
49994fbd55eSBarry Smith         ierr = TransferWrite(viewer,fp,r,0,connectivity,buffer,piece.nconn,MPI_INT,tag);CHKERRQ(ierr);
50094fbd55eSBarry Smith         ierr = TransferWrite(viewer,fp,r,0,offsets,buffer,piece.ncells,MPI_INT,tag);CHKERRQ(ierr);
50194fbd55eSBarry Smith         ierr = TransferWrite(viewer,fp,r,0,types,buffer,piece.ncells,MPI_CHAR,tag);CHKERRQ(ierr);
502552f7358SJed Brown         ierr = PetscFree3(connectivity,offsets,types);CHKERRQ(ierr);
503552f7358SJed Brown       }
504552f7358SJed Brown       {                         /* Owners (cell data) */
505552f7358SJed Brown         PetscVTKInt *owners;
506785e854fSJed Brown         ierr = PetscMalloc1(piece.ncells,&owners);CHKERRQ(ierr);
507552f7358SJed Brown         for (i=0; i<piece.ncells; i++) owners[i] = rank;
50894fbd55eSBarry Smith         ierr = TransferWrite(viewer,fp,r,0,owners,buffer,piece.ncells,MPI_INT,tag);CHKERRQ(ierr);
509552f7358SJed Brown         ierr = PetscFree(owners);CHKERRQ(ierr);
510552f7358SJed Brown       }
511552f7358SJed Brown       /* Cell data */
512552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
513552f7358SJed Brown         Vec               X = (Vec)link->vec;
514e630c359SToby Isaac         DM                dmX;
515552f7358SJed Brown         const PetscScalar *x;
516*955d60d1SToby Isaac         PetscVTUReal      *y;
517e630c359SToby Isaac         PetscInt          bs, nfields, field;
518e630c359SToby Isaac         PetscSection      section = NULL;
519e630c359SToby Isaac 
520552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
521e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
522e630c359SToby Isaac         if (!dmX) dmX = dm;
523e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
524e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
525e630c359SToby Isaac         ierr = PetscSectionGetDof(section,cStart,&bs);CHKERRQ(ierr);
526e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
527e630c359SToby Isaac         field = 0;
528e630c359SToby Isaac         if (link->field >= 0) {
529e630c359SToby Isaac           field = link->field;
530e630c359SToby Isaac           nfields = field + 1;
531e630c359SToby Isaac         }
532552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
533e630c359SToby Isaac         ierr = PetscMalloc1(piece.ncells*3,&y);CHKERRQ(ierr);
534e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
535e630c359SToby Isaac           PetscInt     fbs,j;
536e630c359SToby Isaac           PetscFV      fv = NULL;
537e630c359SToby Isaac           PetscObject  f;
538e630c359SToby Isaac           PetscClassId fClass;
539e630c359SToby Isaac           PetscBool    vector;
540e630c359SToby Isaac           if (nfields) {        /* We have user-defined fields/components */
541e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,cStart,field,&fbs);CHKERRQ(ierr);
542e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
543e630c359SToby Isaac           ierr = DMGetField(dmX,field,NULL,&f);CHKERRQ(ierr);
544e630c359SToby Isaac           ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr);
545e630c359SToby Isaac           if (fClass == PETSCFV_CLASSID) {
546e630c359SToby Isaac             fv = (PetscFV) f;
547e630c359SToby Isaac           }
548e630c359SToby Isaac           vector = PETSC_FALSE;
549e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
550e630c359SToby Isaac             vector = PETSC_TRUE;
551e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
552e630c359SToby Isaac               const char *compName = NULL;
553e630c359SToby Isaac               if (fv) {
554e630c359SToby Isaac                 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
555e630c359SToby Isaac                 if (compName) break;
556e630c359SToby Isaac               }
557e630c359SToby Isaac             }
558e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
559e630c359SToby Isaac           }
560e630c359SToby Isaac           if (vector) {
561*955d60d1SToby Isaac             PetscInt cnt, l;
562*955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
563552f7358SJed Brown               for (c=cStart,cnt=0; c<cEnd; c++) {
564e630c359SToby Isaac                 const PetscScalar *xpoint;
565e630c359SToby Isaac                 PetscInt off, j;
566e630c359SToby Isaac 
567552f7358SJed Brown                 if (hasLabel) {     /* Ignore some cells */
568552f7358SJed Brown                   PetscInt value;
569e630c359SToby Isaac                   ierr = DMGetLabelValue(dmX, "vtk", c, &value);CHKERRQ(ierr);
570552f7358SJed Brown                   if (value != 1) continue;
571552f7358SJed Brown                 }
572e630c359SToby Isaac                 if (nfields) {
573e630c359SToby Isaac                   ierr = PetscSectionGetFieldOffset(section, c, field, &off);CHKERRQ(ierr);
574e630c359SToby Isaac                 } else {
575e630c359SToby Isaac                   ierr = PetscSectionGetOffset(section, c, &off);CHKERRQ(ierr);
576e630c359SToby Isaac                 }
577e630c359SToby Isaac                 xpoint = &x[off];
578e630c359SToby Isaac                 for (j = 0; j < fbs; j++) {
579*955d60d1SToby Isaac                   y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j]));
580e630c359SToby Isaac                 }
581e630c359SToby Isaac                 for (; j < 3; j++) y[cnt++] = 0.;
582e630c359SToby Isaac               }
583e630c359SToby Isaac               if (cnt != piece.ncells*3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
584*955d60d1SToby Isaac               ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells*3,MPIU_VTUREAL,tag);CHKERRQ(ierr);
585*955d60d1SToby Isaac             }
586e630c359SToby Isaac           } else {
587e630c359SToby Isaac             for (i=0; i<fbs; i++) {
588*955d60d1SToby Isaac               PetscInt cnt, l;
589*955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
590e630c359SToby Isaac                 for (c=cStart,cnt=0; c<cEnd; c++) {
591e630c359SToby Isaac                   const PetscScalar *xpoint;
592e630c359SToby Isaac                   PetscInt off;
593e630c359SToby Isaac 
594e630c359SToby Isaac                   if (hasLabel) {     /* Ignore some cells */
595e630c359SToby Isaac                     PetscInt value;
596e630c359SToby Isaac                     ierr = DMGetLabelValue(dmX, "vtk", c, &value);CHKERRQ(ierr);
597e630c359SToby Isaac                     if (value != 1) continue;
598e630c359SToby Isaac                   }
599e630c359SToby Isaac                   if (nfields) {
600e630c359SToby Isaac                     ierr = PetscSectionGetFieldOffset(section, c, field, &off);CHKERRQ(ierr);
601e630c359SToby Isaac                   } else {
602e630c359SToby Isaac                     ierr = PetscSectionGetOffset(section, c, &off);CHKERRQ(ierr);
603e630c359SToby Isaac                   }
604e630c359SToby Isaac                   xpoint   = &x[off];
605*955d60d1SToby Isaac                   y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i]));
606552f7358SJed Brown                 }
607552f7358SJed Brown                 if (cnt != piece.ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
608*955d60d1SToby Isaac                 ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells,MPIU_VTUREAL,tag);CHKERRQ(ierr);
609*955d60d1SToby Isaac               }
610552f7358SJed Brown             }
611e630c359SToby Isaac           }
612e630c359SToby Isaac         }
613552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
614552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
615552f7358SJed Brown       }
616e630c359SToby Isaac       /* point data */
617552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
618552f7358SJed Brown         Vec               X = (Vec)link->vec;
619e630c359SToby Isaac         DM                dmX;
620552f7358SJed Brown         const PetscScalar *x;
621*955d60d1SToby Isaac         PetscVTUReal      *y;
622e630c359SToby Isaac         PetscInt          bs, nfields, field;
623e630c359SToby Isaac         PetscSection      section = NULL;
624e630c359SToby Isaac 
625552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
626e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
627e630c359SToby Isaac         if (!dmX) dmX = dm;
628e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
629e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
630e630c359SToby Isaac         ierr = PetscSectionGetDof(section,vStart,&bs);CHKERRQ(ierr);
631e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
632e630c359SToby Isaac         field = 0;
633e630c359SToby Isaac         if (link->field >= 0) {
634e630c359SToby Isaac           field = link->field;
635e630c359SToby Isaac           nfields = field + 1;
636e630c359SToby Isaac         }
637552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
638e630c359SToby Isaac         ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr);
639e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
640e630c359SToby Isaac           PetscInt   fbs,j;
641e630c359SToby Isaac           if (nfields) {        /* We have user-defined fields/components */
642e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,vStart,field,&fbs);CHKERRQ(ierr);
643e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
644e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
645*955d60d1SToby Isaac             PetscInt cnt, l;
646*955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
6472e529ec8SStefano Zampini               if (!localized) {
648552f7358SJed Brown                 for (v=vStart,cnt=0; v<vEnd; v++) {
649e630c359SToby Isaac                   PetscInt    off;
650e630c359SToby Isaac                   const PetscScalar *xpoint;
651e630c359SToby Isaac 
652e630c359SToby Isaac                   if (nfields) {
653e630c359SToby Isaac                     ierr     = PetscSectionGetFieldOffset(section,v,field,&off);CHKERRQ(ierr);
654e630c359SToby Isaac                   } else {
655e630c359SToby Isaac                     ierr     = PetscSectionGetOffset(section,v,&off);CHKERRQ(ierr);
656e630c359SToby Isaac                   }
657e630c359SToby Isaac                   xpoint   = &x[off];
658e630c359SToby Isaac                   for (j = 0; j < fbs; j++) {
659*955d60d1SToby Isaac                     y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j]));
660e630c359SToby Isaac                   }
661e630c359SToby Isaac                   for (; j < 3; j++) y[cnt++] = 0.;
662e630c359SToby Isaac                 }
663e630c359SToby Isaac               } else {
664e630c359SToby Isaac                 for (c=cStart,cnt=0; c<cEnd; c++) {
665e630c359SToby Isaac                   PetscInt *closure = NULL;
666e630c359SToby Isaac                   PetscInt  closureSize, off;
667e630c359SToby Isaac 
668e630c359SToby Isaac                   ierr = DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
669e630c359SToby Isaac                   for (v = 0, off = 0; v < closureSize*2; v += 2) {
670e630c359SToby Isaac                     if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
671e630c359SToby Isaac                       PetscInt    voff;
672e630c359SToby Isaac                       const PetscScalar *xpoint;
673e630c359SToby Isaac 
674e630c359SToby Isaac                       if (nfields) {
675e630c359SToby Isaac                         ierr = PetscSectionGetFieldOffset(section,closure[v],field,&voff);CHKERRQ(ierr);
676e630c359SToby Isaac                       } else {
677e630c359SToby Isaac                         ierr = PetscSectionGetOffset(section,closure[v],&voff);CHKERRQ(ierr);
678e630c359SToby Isaac                       }
679e630c359SToby Isaac                       xpoint         = &x[voff];
680e630c359SToby Isaac                       for (j = 0; j < fbs; j++) {
681*955d60d1SToby Isaac                         y[cnt + off++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j]));
682e630c359SToby Isaac                       }
683e630c359SToby Isaac                       for (; j < 3; j++) y[cnt + off++] = 0.;
684e630c359SToby Isaac                     }
685e630c359SToby Isaac                   }
686e630c359SToby Isaac                   cnt += off;
687e630c359SToby Isaac                   ierr = DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
688e630c359SToby Isaac                 }
689e630c359SToby Isaac               }
690e630c359SToby Isaac               if (cnt != piece.nvertices*3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
691*955d60d1SToby Isaac               ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices*3,MPIU_VTUREAL,tag);CHKERRQ(ierr);
692*955d60d1SToby Isaac             }
693e630c359SToby Isaac           } else {
694e630c359SToby Isaac             for (i=0; i<fbs; i++) {
695*955d60d1SToby Isaac               PetscInt cnt, l;
696*955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
697e630c359SToby Isaac                 if (!localized) {
698e630c359SToby Isaac                   for (v=vStart,cnt=0; v<vEnd; v++) {
699e630c359SToby Isaac                     PetscInt    off;
700e630c359SToby Isaac                     const PetscScalar *xpoint;
701e630c359SToby Isaac 
702e630c359SToby Isaac                     if (nfields) {
703e630c359SToby Isaac                       ierr     = PetscSectionGetFieldOffset(section,v,field,&off);CHKERRQ(ierr);
704e630c359SToby Isaac                     } else {
705e630c359SToby Isaac                       ierr     = PetscSectionGetOffset(section,v,&off);CHKERRQ(ierr);
706e630c359SToby Isaac                     }
707e630c359SToby Isaac                     xpoint   = &x[off];
708*955d60d1SToby Isaac                     y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i]));
709552f7358SJed Brown                   }
7102e529ec8SStefano Zampini                 } else {
7112e529ec8SStefano Zampini                   for (c=cStart,cnt=0; c<cEnd; c++) {
7122e529ec8SStefano Zampini                     PetscInt *closure = NULL;
7132e529ec8SStefano Zampini                     PetscInt  closureSize, off;
7142e529ec8SStefano Zampini 
715e630c359SToby Isaac                     ierr = DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
7162e529ec8SStefano Zampini                     for (v = 0, off = 0; v < closureSize*2; v += 2) {
7172e529ec8SStefano Zampini                       if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
718e630c359SToby Isaac                         PetscInt    voff;
719e630c359SToby Isaac                         const PetscScalar *xpoint;
7202e529ec8SStefano Zampini 
721e630c359SToby Isaac                         if (nfields) {
722e630c359SToby Isaac                           ierr           = PetscSectionGetFieldOffset(section,closure[v],field,&voff);CHKERRQ(ierr);
723e630c359SToby Isaac                         } else {
724e630c359SToby Isaac                           ierr           = PetscSectionGetOffset(section,closure[v],&voff);CHKERRQ(ierr);
725e630c359SToby Isaac                         }
726e630c359SToby Isaac                         xpoint         = &x[voff];
727*955d60d1SToby Isaac                         y[cnt + off++] = (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i]));
7282e529ec8SStefano Zampini                       }
7292e529ec8SStefano Zampini                     }
7309bda96f6SStefano Zampini                     cnt += off;
731e630c359SToby Isaac                     ierr = DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
7322e529ec8SStefano Zampini                   }
7332e529ec8SStefano Zampini                 }
734552f7358SJed Brown                 if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
735*955d60d1SToby Isaac                 ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices,MPIU_VTUREAL,tag);CHKERRQ(ierr);
736*955d60d1SToby Isaac               }
737552f7358SJed Brown             }
738e630c359SToby Isaac           }
739e630c359SToby Isaac         }
740552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
741552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
742552f7358SJed Brown       }
743552f7358SJed Brown     } else if (!rank) {
744*955d60d1SToby Isaac       PetscInt l;
745*955d60d1SToby Isaac 
746*955d60d1SToby Isaac       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,MPIU_VTUREAL,tag);CHKERRQ(ierr); /* positions */
74794fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nconn,MPI_INT,tag);CHKERRQ(ierr); /* connectivity */
74894fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag);CHKERRQ(ierr); /* offsets */
74994fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_CHAR,tag);CHKERRQ(ierr); /* types */
75094fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag);CHKERRQ(ierr); /* owner rank (cells) */
751552f7358SJed Brown       /* all cell data */
752552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
753e630c359SToby Isaac         Vec               X = (Vec)link->vec;
754e630c359SToby Isaac         PetscInt bs, nfields, field;
755e630c359SToby Isaac         DM           dmX;
756e630c359SToby Isaac         PetscSection section = NULL;
757e630c359SToby Isaac 
758552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
759e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
760e630c359SToby Isaac         if (!dmX) dmX = dm;
761e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
762e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
763e630c359SToby Isaac         ierr = PetscSectionGetDof(section,cStart,&bs);CHKERRQ(ierr);
764e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
765e630c359SToby Isaac         field = 0;
766e630c359SToby Isaac         if (link->field >= 0) {
767e630c359SToby Isaac           field = link->field;
768e630c359SToby Isaac           nfields = field + 1;
769e630c359SToby Isaac         }
770e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
771e630c359SToby Isaac           PetscInt     fbs,j;
772e630c359SToby Isaac           PetscFV      fv = NULL;
773e630c359SToby Isaac           PetscObject  f;
774e630c359SToby Isaac           PetscClassId fClass;
775e630c359SToby Isaac           PetscBool    vector;
776e630c359SToby Isaac           if (nfields) {        /* We have user-defined fields/components */
777e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,cStart,field,&fbs);CHKERRQ(ierr);
778e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
779e630c359SToby Isaac           ierr = DMGetField(dmX,field,NULL,&f);CHKERRQ(ierr);
780e630c359SToby Isaac           ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr);
781e630c359SToby Isaac           if (fClass == PETSCFV_CLASSID) {
782e630c359SToby Isaac             fv = (PetscFV) f;
783e630c359SToby Isaac           }
784e630c359SToby Isaac           vector = PETSC_FALSE;
785e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
786e630c359SToby Isaac             vector = PETSC_TRUE;
787e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
788e630c359SToby Isaac               const char *compName = NULL;
789e630c359SToby Isaac               if (fv) {
790e630c359SToby Isaac                 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
791e630c359SToby Isaac                 if (compName) break;
792e630c359SToby Isaac               }
793e630c359SToby Isaac             }
794e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
795e630c359SToby Isaac           }
796e630c359SToby Isaac           if (vector) {
797*955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
798*955d60d1SToby Isaac               ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells*3,MPIU_VTUREAL,tag);CHKERRQ(ierr);
799*955d60d1SToby Isaac             }
800e630c359SToby Isaac           } else {
801e630c359SToby Isaac             for (i=0; i<fbs; i++) {
802*955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
803*955d60d1SToby Isaac                 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPIU_VTUREAL,tag);CHKERRQ(ierr);
804*955d60d1SToby Isaac               }
805552f7358SJed Brown             }
806552f7358SJed Brown           }
807e630c359SToby Isaac         }
808e630c359SToby Isaac       }
809552f7358SJed Brown       /* all point data */
810552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
811e630c359SToby Isaac         Vec               X = (Vec)link->vec;
812e630c359SToby Isaac         DM                dmX;
813e630c359SToby Isaac         PetscInt bs, nfields, field;
814e630c359SToby Isaac         PetscSection section = NULL;
815e630c359SToby Isaac 
816552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
817e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
818e630c359SToby Isaac         if (!dmX) dmX = dm;
819e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
820e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
821e630c359SToby Isaac         ierr = PetscSectionGetDof(section,vStart,&bs);CHKERRQ(ierr);
822e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
823e630c359SToby Isaac         field = 0;
824e630c359SToby Isaac         if (link->field >= 0) {
825e630c359SToby Isaac           field = link->field;
826e630c359SToby Isaac           nfields = field + 1;
827e630c359SToby Isaac         }
828e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
829e630c359SToby Isaac           PetscInt   fbs;
830e630c359SToby Isaac           if (nfields) {        /* We have user-defined fields/components */
831e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,vStart,field,&fbs);CHKERRQ(ierr);
832e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
833e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
834*955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
835*955d60d1SToby Isaac               ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,MPIU_VTUREAL,tag);CHKERRQ(ierr);
836*955d60d1SToby Isaac             }
837e630c359SToby Isaac           } else {
838e630c359SToby Isaac             for (i=0; i<fbs; i++) {
839*955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
840*955d60d1SToby Isaac                 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices,MPIU_VTUREAL,tag);CHKERRQ(ierr);
841*955d60d1SToby Isaac               }
842552f7358SJed Brown             }
843552f7358SJed Brown           }
844552f7358SJed Brown         }
845552f7358SJed Brown       }
846e630c359SToby Isaac     }
847e630c359SToby Isaac   }
848552f7358SJed Brown   ierr = PetscFree(gpiece);CHKERRQ(ierr);
849552f7358SJed Brown   ierr = PetscFree(buffer);CHKERRQ(ierr);
850552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"\n  </AppendedData>\n");CHKERRQ(ierr);
851552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
8525e32de72SMatthew G. Knepley   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
853552f7358SJed Brown   PetscFunctionReturn(0);
854552f7358SJed Brown }
855