xref: /petsc/src/dm/impls/plex/plexvtu.c (revision ffc4695bcb29f4b022f59a5fd6bc99fc280ff6d8)
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 
10955d60d1SToby Isaac #if defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL___FP16)
11955d60d1SToby Isaac /* output in float if single or half precision in memory */
12552f7358SJed Brown static const char precision[] = "Float32";
13955d60d1SToby Isaac typedef float PetscVTUReal;
14955d60d1SToby Isaac #define MPIU_VTUREAL MPI_FLOAT
15955d60d1SToby Isaac #elif defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
16955d60d1SToby Isaac /* output in double if double or quad precision in memory */
17552f7358SJed Brown static const char precision[] = "Float64";
18955d60d1SToby Isaac typedef double PetscVTUReal;
19955d60d1SToby Isaac #define MPIU_VTUREAL MPI_DOUBLE
20552f7358SJed Brown #else
21552f7358SJed Brown static const char precision[] = "UnknownPrecision";
22955d60d1SToby Isaac typedef PetscReal PetscVTUReal;
23955d60d1SToby 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);
34*ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
35552f7358SJed Brown 
36552f7358SJed Brown   if (rank == srank && rank != root) {
37*ffc4695bSBarry Smith     ierr = MPI_Send((void*)send,count,mpidatatype,root,tag,comm);CHKERRMPI(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;
45*ffc4695bSBarry Smith       ierr = MPI_Recv(recv,count,mpidatatype,srank,tag,comm,&status);CHKERRMPI(ierr);
46*ffc4695bSBarry Smith       ierr = MPI_Get_count(&status,mpidatatype,&nrecv);CHKERRMPI(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     }
10596ca5757SLisandro Dalcin 
10696ca5757SLisandro Dalcin     {
10796ca5757SLisandro Dalcin       PetscInt n = PetscMin(nC, 8), s = countconn - nC, i, cone[8];
10896ca5757SLisandro Dalcin       for (i = 0; i < n; ++i) cone[i] = conn[s+i];
10996ca5757SLisandro Dalcin       ierr = DMPlexReorderCell(dm, c, cone);CHKERRQ(ierr);
11096ca5757SLisandro Dalcin       for (i = 0; i < n; ++i) conn[s+i] = (int)cone[i];
11196ca5757SLisandro Dalcin     }
1128865f1eaSKarl Rupp 
113552f7358SJed Brown     offsets[countcell] = countconn;
1148865f1eaSKarl Rupp 
115552f7358SJed Brown     nverts = countconn - startoffset;
116de0a4605SMatthew G. Knepley     ierr   = DMPlexVTKGetCellType_Internal(dm,dim,nverts,&celltype);CHKERRQ(ierr);
1178865f1eaSKarl Rupp 
118552f7358SJed Brown     types[countcell] = celltype;
119552f7358SJed Brown     countcell++;
120552f7358SJed Brown   }
121552f7358SJed Brown   if (countcell != piece->ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count");
122552f7358SJed Brown   if (countconn != piece->nconn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count");
123552f7358SJed Brown   *oconn    = conn;
124552f7358SJed Brown   *ooffsets = offsets;
125552f7358SJed Brown   *otypes   = types;
126552f7358SJed Brown   PetscFunctionReturn(0);
127552f7358SJed Brown }
128552f7358SJed Brown 
129552f7358SJed Brown /*
130552f7358SJed Brown   Write all fields that have been provided to the viewer
131552f7358SJed Brown   Multi-block XML format with binary appended data.
132552f7358SJed Brown */
133552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer)
134552f7358SJed Brown {
135ce94432eSBarry Smith   MPI_Comm                 comm;
1362e529ec8SStefano Zampini   PetscSection             coordSection;
137552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
138552f7358SJed Brown   PetscViewerVTKObjectLink link;
139552f7358SJed Brown   FILE                     *fp;
140552f7358SJed Brown   PetscMPIInt              rank,size,tag;
141552f7358SJed Brown   PetscErrorCode           ierr;
1422f029394SStefano Zampini   PetscInt                 dimEmbed,cellHeight,cStart,cEnd,vStart,vEnd,numLabelCells,hasLabel,c,v,r,i;
1432e529ec8SStefano Zampini   PetscBool                localized;
1440298fd71SBarry Smith   PieceInfo                piece,*gpiece = NULL;
1450298fd71SBarry Smith   void                     *buffer = NULL;
14630815ce0SLisandro Dalcin   const char               *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian";
147955d60d1SToby Isaac   PetscInt                 loops_per_scalar;
148552f7358SJed Brown 
149552f7358SJed Brown   PetscFunctionBegin;
150ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
151552f7358SJed Brown #if defined(PETSC_USE_COMPLEX)
152955d60d1SToby Isaac   loops_per_scalar = 2;
153955d60d1SToby Isaac #else
154955d60d1SToby Isaac   loops_per_scalar = 1;
155552f7358SJed Brown #endif
156*ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
157*ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
158552f7358SJed Brown   tag  = ((PetscObject)viewer)->tag;
159552f7358SJed Brown 
160552f7358SJed Brown   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
161552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
16230815ce0SLisandro Dalcin   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"%s\">\n", byte_order);CHKERRQ(ierr);
163552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <UnstructuredGrid>\n");CHKERRQ(ierr);
164552f7358SJed Brown 
1653baa072aSToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
166552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
167552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
168552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
169c58f1c22SToby Isaac   ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
1702e529ec8SStefano Zampini   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
1712e529ec8SStefano Zampini   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1728865f1eaSKarl Rupp 
173552f7358SJed Brown   hasLabel        = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
1742e529ec8SStefano Zampini   piece.nvertices = 0;
175552f7358SJed Brown   piece.ncells    = 0;
176552f7358SJed Brown   piece.nconn     = 0;
1772e529ec8SStefano Zampini   if (!localized) piece.nvertices = vEnd - vStart;
178552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
17919003fb5SStefano Zampini     PetscInt dof = 0;
180552f7358SJed Brown     if (hasLabel) {
181552f7358SJed Brown       PetscInt value;
182552f7358SJed Brown 
183c58f1c22SToby Isaac       ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
184552f7358SJed Brown       if (value != 1) continue;
185552f7358SJed Brown     }
18619003fb5SStefano Zampini     if (localized) {
18719003fb5SStefano Zampini       ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
18819003fb5SStefano Zampini     }
18919003fb5SStefano Zampini     if (!dof) {
1902e529ec8SStefano Zampini       PetscInt *closure = NULL;
1912e529ec8SStefano Zampini       PetscInt closureSize;
1922e529ec8SStefano Zampini 
193552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
194552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
1952e529ec8SStefano Zampini         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
1962e529ec8SStefano Zampini           piece.nconn++;
19719003fb5SStefano Zampini           if (localized) piece.nvertices++;
1982e529ec8SStefano Zampini         }
199552f7358SJed Brown       }
200552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
2012e529ec8SStefano Zampini     } else {
2022e529ec8SStefano Zampini       piece.nvertices += dof/dimEmbed;
2032e529ec8SStefano Zampini       piece.nconn     += dof/dimEmbed;
2042e529ec8SStefano Zampini     }
205552f7358SJed Brown     piece.ncells++;
206552f7358SJed Brown   }
207785e854fSJed Brown   if (!rank) {ierr = PetscMalloc1(size,&gpiece);CHKERRQ(ierr);}
208*ffc4695bSBarry Smith   ierr = MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm);CHKERRMPI(ierr);
209552f7358SJed Brown 
210552f7358SJed Brown   /*
211552f7358SJed Brown    * Write file header
212552f7358SJed Brown    */
213552f7358SJed Brown   if (!rank) {
214552f7358SJed Brown     PetscInt boffset = 0;
215552f7358SJed Brown 
216552f7358SJed Brown     for (r=0; r<size; r++) {
217552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells);CHKERRQ(ierr);
218552f7358SJed Brown       /* Coordinate positions */
219552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Points>\n");CHKERRQ(ierr);
220552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
221955d60d1SToby Isaac       boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int);
222552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Points>\n");CHKERRQ(ierr);
223552f7358SJed Brown       /* Cell connectivity */
224552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Cells>\n");CHKERRQ(ierr);
225552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
226552f7358SJed Brown       boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int);
227552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"offsets\"      NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
228552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int);
229552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"UInt8\" Name=\"types\"        NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
230552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int);
231552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Cells>\n");CHKERRQ(ierr);
232552f7358SJed Brown 
233552f7358SJed Brown       /*
234552f7358SJed Brown        * Cell Data headers
235552f7358SJed Brown        */
236552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <CellData>\n");CHKERRQ(ierr);
237552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
238552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(int) + sizeof(int);
239552f7358SJed Brown       /* all the vectors */
240552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
241552f7358SJed Brown         Vec        X = (Vec)link->vec;
242e630c359SToby Isaac         DM         dmX = NULL;
2431cfafdd3SJed Brown         PetscInt   bs,nfields,field;
244552f7358SJed Brown         const char *vecname = "";
245e630c359SToby Isaac         PetscSection section;
246552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
247552f7358SJed 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. */
248552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
249552f7358SJed Brown         }
250e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
251e630c359SToby Isaac         if (!dmX) dmX = dm;
252e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
253e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
254e630c359SToby Isaac         ierr = PetscSectionGetDof(section,cStart,&bs);CHKERRQ(ierr);
255e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
256e630c359SToby Isaac         field = 0;
257e630c359SToby Isaac         if (link->field >= 0) {
258e630c359SToby Isaac           field = link->field;
259e630c359SToby Isaac           nfields = field + 1;
260e630c359SToby Isaac         }
261e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
2621cfafdd3SJed Brown           PetscInt     fbs,j;
263a00cdb45SToby Isaac           PetscFV      fv = NULL;
264a00cdb45SToby Isaac           PetscObject  f;
265a00cdb45SToby Isaac           PetscClassId fClass;
2660298fd71SBarry Smith           const char *fieldname = NULL;
2671cfafdd3SJed Brown           char       buf[256];
268e630c359SToby Isaac           PetscBool    vector;
2697ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
270e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,cStart,field,&fbs);CHKERRQ(ierr);
271e630c359SToby Isaac             ierr = PetscSectionGetFieldName(section,field,&fieldname);CHKERRQ(ierr);
2727ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
273e630c359SToby Isaac           ierr = DMGetField(dmX,field,NULL,&f);CHKERRQ(ierr);
274a00cdb45SToby Isaac           ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr);
275a00cdb45SToby Isaac           if (fClass == PETSCFV_CLASSID) {
276a00cdb45SToby Isaac             fv = (PetscFV) f;
277a00cdb45SToby Isaac           }
278e630c359SToby Isaac           if (nfields && !fieldname) {
2791cfafdd3SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"CellField%D",field);CHKERRQ(ierr);
280552f7358SJed Brown             fieldname = buf;
281552f7358SJed Brown           }
282e630c359SToby Isaac           vector = PETSC_FALSE;
283e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
284e630c359SToby Isaac             vector = PETSC_TRUE;
285e630c359SToby 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);
286e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
287e630c359SToby Isaac               const char *compName = NULL;
288e630c359SToby Isaac               if (fv) {
289e630c359SToby Isaac                 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
290e630c359SToby Isaac                 if (compName) break;
291e630c359SToby Isaac               }
292e630c359SToby Isaac             }
293e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
294e630c359SToby Isaac           }
295e630c359SToby Isaac           if (vector) {
296955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
297955d60d1SToby 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);
298955d60d1SToby Isaac             boffset += gpiece[r].ncells*3*sizeof(PetscVTUReal) + sizeof(int);
299955d60d1SToby 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);
300955d60d1SToby Isaac             boffset += gpiece[r].ncells*3*sizeof(PetscVTUReal) + sizeof(int);
301955d60d1SToby Isaac #else
302e630c359SToby Isaac             ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
303955d60d1SToby Isaac             boffset += gpiece[r].ncells*3*sizeof(PetscVTUReal) + sizeof(int);
304955d60d1SToby Isaac #endif
305e630c359SToby Isaac             i+=fbs;
306e630c359SToby Isaac           } else {
3071cfafdd3SJed Brown             for (j=0; j<fbs; j++) {
308a00cdb45SToby Isaac               const char *compName = NULL;
309955d60d1SToby Isaac               char finalname[256];
310a00cdb45SToby Isaac               if (fv) {
311a00cdb45SToby Isaac                 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
312a00cdb45SToby Isaac               }
313a00cdb45SToby Isaac               if (compName) {
314955d60d1SToby Isaac                 ierr = PetscSNPrintf(finalname,255,"%s%s.%s",vecname,fieldname,compName);CHKERRQ(ierr);
315a00cdb45SToby Isaac               }
316e630c359SToby Isaac               else if (fbs > 1) {
317955d60d1SToby Isaac                 ierr = PetscSNPrintf(finalname,255,"%s%s.%D",vecname,fieldname,j);CHKERRQ(ierr);
318e630c359SToby Isaac               } else {
319955d60d1SToby Isaac                 ierr = PetscSNPrintf(finalname,255,"%s%s",vecname,fieldname);CHKERRQ(ierr);
320a00cdb45SToby Isaac               }
321955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
322955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
323955d60d1SToby Isaac               boffset += gpiece[r].ncells*sizeof(PetscVTUReal) + sizeof(int);
324955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
325955d60d1SToby Isaac               boffset += gpiece[r].ncells*sizeof(PetscVTUReal) + sizeof(int);
326955d60d1SToby Isaac #else
327955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
328955d60d1SToby Isaac               boffset += gpiece[r].ncells*sizeof(PetscVTUReal) + sizeof(int);
329955d60d1SToby Isaac #endif
3301cfafdd3SJed Brown               i++;
331552f7358SJed Brown             }
332552f7358SJed Brown           }
333e630c359SToby Isaac         }
3348ec8862eSJed Brown         //if (i != bs) SETERRQ2(comm,PETSC_ERR_PLIB,"Total number of field components %D != block size %D",i,bs);
3351cfafdd3SJed Brown       }
336552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </CellData>\n");CHKERRQ(ierr);
337552f7358SJed Brown 
338552f7358SJed Brown       /*
339552f7358SJed Brown        * Point Data headers
340552f7358SJed Brown        */
341552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      <PointData>\n");CHKERRQ(ierr);
342552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
343552f7358SJed Brown         Vec        X = (Vec)link->vec;
344e630c359SToby Isaac         DM         dmX;
3451cfafdd3SJed Brown         PetscInt   bs,nfields,field;
346552f7358SJed Brown         const char *vecname = "";
347e630c359SToby Isaac         PetscSection section;
348552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
349552f7358SJed 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. */
350552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
351552f7358SJed Brown         }
352e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
353e630c359SToby Isaac         if (!dmX) dmX = dm;
354e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
355e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
356e630c359SToby Isaac         ierr = PetscSectionGetDof(section,vStart,&bs);CHKERRQ(ierr);
357e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
358e630c359SToby Isaac         field = 0;
359e630c359SToby Isaac         if (link->field >= 0) {
360e630c359SToby Isaac           field = link->field;
361e630c359SToby Isaac           nfields = field + 1;
362e630c359SToby Isaac         }
363e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
3641cfafdd3SJed Brown           PetscInt   fbs,j;
3651cfafdd3SJed Brown           const char *fieldname = NULL;
3661cfafdd3SJed Brown           char       buf[256];
3677ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
368e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,vStart,field,&fbs);CHKERRQ(ierr);
369e630c359SToby Isaac             ierr = PetscSectionGetFieldName(section,field,&fieldname);CHKERRQ(ierr);
3707ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
371e630c359SToby Isaac           if (nfields && !fieldname) {
3721cfafdd3SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"PointField%D",field);CHKERRQ(ierr);
3731cfafdd3SJed Brown             fieldname = buf;
3741cfafdd3SJed Brown           }
375e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
376e630c359SToby 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);
377955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
378955d60d1SToby 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);
379955d60d1SToby Isaac             boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int);
380955d60d1SToby 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);
381955d60d1SToby Isaac             boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int);
382955d60d1SToby Isaac #else
383e630c359SToby Isaac             ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
384955d60d1SToby Isaac             boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int);
385955d60d1SToby Isaac #endif
386e630c359SToby Isaac           } else {
3871cfafdd3SJed Brown             for (j=0; j<fbs; j++) {
388b778fa18SValeria Barra               const char *compName = NULL;
389955d60d1SToby Isaac               char finalname[256];
390b778fa18SValeria Barra               ierr = PetscSectionGetComponentName(section,field,j,&compName);CHKERRQ(ierr);
391b778fa18SValeria Barra               ierr = PetscSNPrintf(finalname,255,"%s%s.%s",vecname,fieldname,compName);CHKERRQ(ierr);
392955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
393955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
394955d60d1SToby Isaac               boffset += gpiece[r].nvertices*sizeof(PetscVTUReal) + sizeof(int);
395955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
396955d60d1SToby Isaac               boffset += gpiece[r].nvertices*sizeof(PetscVTUReal) + sizeof(int);
397955d60d1SToby Isaac #else
398955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
399955d60d1SToby Isaac               boffset += gpiece[r].nvertices*sizeof(PetscVTUReal) + sizeof(int);
400955d60d1SToby Isaac #endif
401552f7358SJed Brown             }
402552f7358SJed Brown           }
4031cfafdd3SJed Brown         }
404e630c359SToby Isaac       }
405552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </PointData>\n");CHKERRQ(ierr);
406552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    </Piece>\n");CHKERRQ(ierr);
407552f7358SJed Brown     }
408552f7358SJed Brown   }
409552f7358SJed Brown 
410552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  </UnstructuredGrid>\n");CHKERRQ(ierr);
411552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
412552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
413552f7358SJed Brown 
414552f7358SJed Brown   if (!rank) {
415552f7358SJed Brown     PetscInt maxsize = 0;
416552f7358SJed Brown     for (r=0; r<size; r++) {
417955d60d1SToby Isaac       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscVTUReal)));
418955d60d1SToby Isaac       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*3*sizeof(PetscVTUReal)));
419552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt)));
420552f7358SJed Brown     }
421552f7358SJed Brown     ierr = PetscMalloc(maxsize,&buffer);CHKERRQ(ierr);
422552f7358SJed Brown   }
423552f7358SJed Brown   for (r=0; r<size; r++) {
424552f7358SJed Brown     if (r == rank) {
425552f7358SJed Brown       PetscInt nsend;
426552f7358SJed Brown       {                         /* Position */
427552f7358SJed Brown         const PetscScalar *x;
428955d60d1SToby Isaac         PetscVTUReal      *y = NULL;
429552f7358SJed Brown         Vec               coords;
430955d60d1SToby Isaac         PetscBool         copy;
4312e529ec8SStefano Zampini 
432552f7358SJed Brown         ierr  = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
433552f7358SJed Brown         ierr  = VecGetArrayRead(coords,&x);CHKERRQ(ierr);
434955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
435955d60d1SToby Isaac         copy = PETSC_TRUE;
436955d60d1SToby Isaac #else
437955d60d1SToby Isaac         copy = (PetscBool) (dimEmbed != 3 || localized || (sizeof(PetscReal) != sizeof(PetscVTUReal)));
438955d60d1SToby Isaac #endif
439955d60d1SToby Isaac         if (copy) {
440785e854fSJed Brown           ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr);
44119003fb5SStefano Zampini           if (localized) {
44219003fb5SStefano Zampini             PetscInt cnt;
44319003fb5SStefano Zampini             for (c=cStart,cnt=0; c<cEnd; c++) {
44419003fb5SStefano Zampini               PetscInt off, dof;
44519003fb5SStefano Zampini 
44619003fb5SStefano Zampini               ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
44719003fb5SStefano Zampini               if (!dof) {
44819003fb5SStefano Zampini                 PetscInt *closure = NULL;
44919003fb5SStefano Zampini                 PetscInt closureSize;
45019003fb5SStefano Zampini 
45119003fb5SStefano Zampini                 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
45219003fb5SStefano Zampini                 for (v = 0; v < closureSize*2; v += 2) {
45319003fb5SStefano Zampini                   if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
45419003fb5SStefano Zampini                     ierr = PetscSectionGetOffset(coordSection, closure[v], &off);CHKERRQ(ierr);
45519003fb5SStefano Zampini                     if (dimEmbed != 3) {
456955d60d1SToby Isaac                       y[cnt*3+0] = (PetscVTUReal) PetscRealPart(x[off+0]);
457955d60d1SToby Isaac                       y[cnt*3+1] = (PetscVTUReal) ((dimEmbed > 1) ? PetscRealPart(x[off+1]) : 0.0);
458955d60d1SToby Isaac                       y[cnt*3+2] = (PetscVTUReal) 0.0;
45919003fb5SStefano Zampini                     } else {
460955d60d1SToby Isaac                       y[cnt*3+0] = (PetscVTUReal) PetscRealPart(x[off+0]);
461955d60d1SToby Isaac                       y[cnt*3+1] = (PetscVTUReal) PetscRealPart(x[off+1]);
462955d60d1SToby Isaac                       y[cnt*3+2] = (PetscVTUReal) PetscRealPart(x[off+2]);
46319003fb5SStefano Zampini                     }
46419003fb5SStefano Zampini                     cnt++;
46519003fb5SStefano Zampini                   }
46619003fb5SStefano Zampini                 }
46719003fb5SStefano Zampini                 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
46819003fb5SStefano Zampini               } else {
46919003fb5SStefano Zampini                 ierr = PetscSectionGetOffset(coordSection, c, &off);CHKERRQ(ierr);
47019003fb5SStefano Zampini                 if (dimEmbed != 3) {
47119003fb5SStefano Zampini                   for (i=0; i<dof/dimEmbed; i++) {
472955d60d1SToby Isaac                     y[cnt*3+0] = (PetscVTUReal) PetscRealPart(x[off + i*dimEmbed + 0]);
473955d60d1SToby Isaac                     y[cnt*3+1] = (PetscVTUReal) ((dimEmbed > 1) ? PetscRealPart(x[off + i*dimEmbed + 1]) : 0.0);
474955d60d1SToby Isaac                     y[cnt*3+2] = (PetscVTUReal) 0.0;
47519003fb5SStefano Zampini                     cnt++;
47619003fb5SStefano Zampini                   }
47719003fb5SStefano Zampini                 } else {
47819003fb5SStefano Zampini                   for (i=0; i<dof; i ++) {
479955d60d1SToby Isaac                     y[cnt*3+i] = (PetscVTUReal) PetscRealPart(x[off + i]);
48019003fb5SStefano Zampini                   }
48119003fb5SStefano Zampini                   cnt += dof/dimEmbed;
48219003fb5SStefano Zampini                 }
48319003fb5SStefano Zampini               }
48419003fb5SStefano Zampini             }
48519003fb5SStefano Zampini             if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
48619003fb5SStefano Zampini           } else {
487552f7358SJed Brown             for (i=0; i<piece.nvertices; i++) {
488955d60d1SToby Isaac               y[i*3+0] = (PetscVTUReal) PetscRealPart(x[i*dimEmbed+0]);
489955d60d1SToby Isaac               y[i*3+1] = (PetscVTUReal) ((dimEmbed > 1) ? PetscRealPart(x[i*dimEmbed+1]) : 0.);
490955d60d1SToby Isaac               y[i*3+2] = (PetscVTUReal) ((dimEmbed > 2) ? PetscRealPart(x[i*dimEmbed+2]) : 0.);
49119003fb5SStefano Zampini             }
492552f7358SJed Brown           }
493552f7358SJed Brown         }
4942e529ec8SStefano Zampini         nsend = piece.nvertices*3;
495955d60d1SToby Isaac         ierr  = TransferWrite(viewer,fp,r,0,copy ? (const void *) y : (const void *) x,buffer,nsend,MPIU_VTUREAL,tag);CHKERRQ(ierr);
496552f7358SJed Brown         ierr  = PetscFree(y);CHKERRQ(ierr);
497552f7358SJed Brown         ierr  = VecRestoreArrayRead(coords,&x);CHKERRQ(ierr);
498552f7358SJed Brown       }
499552f7358SJed Brown       {                           /* Connectivity, offsets, types */
5003e869605SMatthew G. Knepley         PetscVTKInt  *connectivity = NULL, *offsets = NULL;
5013e869605SMatthew G. Knepley         PetscVTKType *types = NULL;
5022e529ec8SStefano Zampini         ierr = DMPlexGetVTKConnectivity(dm,localized,&piece,&connectivity,&offsets,&types);CHKERRQ(ierr);
50394fbd55eSBarry Smith         ierr = TransferWrite(viewer,fp,r,0,connectivity,buffer,piece.nconn,MPI_INT,tag);CHKERRQ(ierr);
50494fbd55eSBarry Smith         ierr = TransferWrite(viewer,fp,r,0,offsets,buffer,piece.ncells,MPI_INT,tag);CHKERRQ(ierr);
50594fbd55eSBarry Smith         ierr = TransferWrite(viewer,fp,r,0,types,buffer,piece.ncells,MPI_CHAR,tag);CHKERRQ(ierr);
506552f7358SJed Brown         ierr = PetscFree3(connectivity,offsets,types);CHKERRQ(ierr);
507552f7358SJed Brown       }
508552f7358SJed Brown       {                         /* Owners (cell data) */
509552f7358SJed Brown         PetscVTKInt *owners;
510785e854fSJed Brown         ierr = PetscMalloc1(piece.ncells,&owners);CHKERRQ(ierr);
511552f7358SJed Brown         for (i=0; i<piece.ncells; i++) owners[i] = rank;
51294fbd55eSBarry Smith         ierr = TransferWrite(viewer,fp,r,0,owners,buffer,piece.ncells,MPI_INT,tag);CHKERRQ(ierr);
513552f7358SJed Brown         ierr = PetscFree(owners);CHKERRQ(ierr);
514552f7358SJed Brown       }
515552f7358SJed Brown       /* Cell data */
516552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
517552f7358SJed Brown         Vec               X = (Vec)link->vec;
518e630c359SToby Isaac         DM                dmX;
519552f7358SJed Brown         const PetscScalar *x;
520955d60d1SToby Isaac         PetscVTUReal      *y;
521e630c359SToby Isaac         PetscInt          bs, nfields, field;
522e630c359SToby Isaac         PetscSection      section = NULL;
523e630c359SToby Isaac 
524552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
525e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
526e630c359SToby Isaac         if (!dmX) dmX = dm;
527e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
528e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
529e630c359SToby Isaac         ierr = PetscSectionGetDof(section,cStart,&bs);CHKERRQ(ierr);
530e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
531e630c359SToby Isaac         field = 0;
532e630c359SToby Isaac         if (link->field >= 0) {
533e630c359SToby Isaac           field = link->field;
534e630c359SToby Isaac           nfields = field + 1;
535e630c359SToby Isaac         }
536552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
537e630c359SToby Isaac         ierr = PetscMalloc1(piece.ncells*3,&y);CHKERRQ(ierr);
538e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
539e630c359SToby Isaac           PetscInt     fbs,j;
540e630c359SToby Isaac           PetscFV      fv = NULL;
541e630c359SToby Isaac           PetscObject  f;
542e630c359SToby Isaac           PetscClassId fClass;
543e630c359SToby Isaac           PetscBool    vector;
544e630c359SToby Isaac           if (nfields) {        /* We have user-defined fields/components */
545e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,cStart,field,&fbs);CHKERRQ(ierr);
546e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
547e630c359SToby Isaac           ierr = DMGetField(dmX,field,NULL,&f);CHKERRQ(ierr);
548e630c359SToby Isaac           ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr);
549e630c359SToby Isaac           if (fClass == PETSCFV_CLASSID) {
550e630c359SToby Isaac             fv = (PetscFV) f;
551e630c359SToby Isaac           }
552e630c359SToby Isaac           vector = PETSC_FALSE;
553e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
554e630c359SToby Isaac             vector = PETSC_TRUE;
555e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
556e630c359SToby Isaac               const char *compName = NULL;
557e630c359SToby Isaac               if (fv) {
558e630c359SToby Isaac                 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
559e630c359SToby Isaac                 if (compName) break;
560e630c359SToby Isaac               }
561e630c359SToby Isaac             }
562e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
563e630c359SToby Isaac           }
564e630c359SToby Isaac           if (vector) {
565955d60d1SToby Isaac             PetscInt cnt, l;
566955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
567552f7358SJed Brown               for (c=cStart,cnt=0; c<cEnd; c++) {
568e630c359SToby Isaac                 const PetscScalar *xpoint;
569e630c359SToby Isaac                 PetscInt off, j;
570e630c359SToby Isaac 
571552f7358SJed Brown                 if (hasLabel) {     /* Ignore some cells */
572552f7358SJed Brown                   PetscInt value;
573e630c359SToby Isaac                   ierr = DMGetLabelValue(dmX, "vtk", c, &value);CHKERRQ(ierr);
574552f7358SJed Brown                   if (value != 1) continue;
575552f7358SJed Brown                 }
576e630c359SToby Isaac                 if (nfields) {
577e630c359SToby Isaac                   ierr = PetscSectionGetFieldOffset(section, c, field, &off);CHKERRQ(ierr);
578e630c359SToby Isaac                 } else {
579e630c359SToby Isaac                   ierr = PetscSectionGetOffset(section, c, &off);CHKERRQ(ierr);
580e630c359SToby Isaac                 }
581e630c359SToby Isaac                 xpoint = &x[off];
582e630c359SToby Isaac                 for (j = 0; j < fbs; j++) {
583955d60d1SToby Isaac                   y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j]));
584e630c359SToby Isaac                 }
585e630c359SToby Isaac                 for (; j < 3; j++) y[cnt++] = 0.;
586e630c359SToby Isaac               }
587e630c359SToby Isaac               if (cnt != piece.ncells*3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
588955d60d1SToby Isaac               ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells*3,MPIU_VTUREAL,tag);CHKERRQ(ierr);
589955d60d1SToby Isaac             }
590e630c359SToby Isaac           } else {
591e630c359SToby Isaac             for (i=0; i<fbs; i++) {
592955d60d1SToby Isaac               PetscInt cnt, l;
593955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
594e630c359SToby Isaac                 for (c=cStart,cnt=0; c<cEnd; c++) {
595e630c359SToby Isaac                   const PetscScalar *xpoint;
596e630c359SToby Isaac                   PetscInt off;
597e630c359SToby Isaac 
598e630c359SToby Isaac                   if (hasLabel) {     /* Ignore some cells */
599e630c359SToby Isaac                     PetscInt value;
600e630c359SToby Isaac                     ierr = DMGetLabelValue(dmX, "vtk", c, &value);CHKERRQ(ierr);
601e630c359SToby Isaac                     if (value != 1) continue;
602e630c359SToby Isaac                   }
603e630c359SToby Isaac                   if (nfields) {
604e630c359SToby Isaac                     ierr = PetscSectionGetFieldOffset(section, c, field, &off);CHKERRQ(ierr);
605e630c359SToby Isaac                   } else {
606e630c359SToby Isaac                     ierr = PetscSectionGetOffset(section, c, &off);CHKERRQ(ierr);
607e630c359SToby Isaac                   }
608e630c359SToby Isaac                   xpoint   = &x[off];
609955d60d1SToby Isaac                   y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i]));
610552f7358SJed Brown                 }
611552f7358SJed Brown                 if (cnt != piece.ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
612955d60d1SToby Isaac                 ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells,MPIU_VTUREAL,tag);CHKERRQ(ierr);
613955d60d1SToby Isaac               }
614552f7358SJed Brown             }
615e630c359SToby Isaac           }
616e630c359SToby Isaac         }
617552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
618552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
619552f7358SJed Brown       }
620e630c359SToby Isaac       /* point data */
621552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
622552f7358SJed Brown         Vec               X = (Vec)link->vec;
623e630c359SToby Isaac         DM                dmX;
624552f7358SJed Brown         const PetscScalar *x;
625955d60d1SToby Isaac         PetscVTUReal      *y;
626e630c359SToby Isaac         PetscInt          bs, nfields, field;
627e630c359SToby Isaac         PetscSection      section = NULL;
628e630c359SToby Isaac 
629552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
630e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
631e630c359SToby Isaac         if (!dmX) dmX = dm;
632e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
633e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
634e630c359SToby Isaac         ierr = PetscSectionGetDof(section,vStart,&bs);CHKERRQ(ierr);
635e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
636e630c359SToby Isaac         field = 0;
637e630c359SToby Isaac         if (link->field >= 0) {
638e630c359SToby Isaac           field = link->field;
639e630c359SToby Isaac           nfields = field + 1;
640e630c359SToby Isaac         }
641552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
642e630c359SToby Isaac         ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr);
643e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
644e630c359SToby Isaac           PetscInt   fbs,j;
645e630c359SToby Isaac           if (nfields) {        /* We have user-defined fields/components */
646e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,vStart,field,&fbs);CHKERRQ(ierr);
647e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
648e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
649955d60d1SToby Isaac             PetscInt cnt, l;
650955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
6512e529ec8SStefano Zampini               if (!localized) {
652552f7358SJed Brown                 for (v=vStart,cnt=0; v<vEnd; v++) {
653e630c359SToby Isaac                   PetscInt    off;
654e630c359SToby Isaac                   const PetscScalar *xpoint;
655e630c359SToby Isaac 
656e630c359SToby Isaac                   if (nfields) {
657e630c359SToby Isaac                     ierr     = PetscSectionGetFieldOffset(section,v,field,&off);CHKERRQ(ierr);
658e630c359SToby Isaac                   } else {
659e630c359SToby Isaac                     ierr     = PetscSectionGetOffset(section,v,&off);CHKERRQ(ierr);
660e630c359SToby Isaac                   }
661e630c359SToby Isaac                   xpoint   = &x[off];
662e630c359SToby Isaac                   for (j = 0; j < fbs; j++) {
663955d60d1SToby Isaac                     y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j]));
664e630c359SToby Isaac                   }
665e630c359SToby Isaac                   for (; j < 3; j++) y[cnt++] = 0.;
666e630c359SToby Isaac                 }
667e630c359SToby Isaac               } else {
668e630c359SToby Isaac                 for (c=cStart,cnt=0; c<cEnd; c++) {
669e630c359SToby Isaac                   PetscInt *closure = NULL;
670e630c359SToby Isaac                   PetscInt  closureSize, off;
671e630c359SToby Isaac 
672e630c359SToby Isaac                   ierr = DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
673e630c359SToby Isaac                   for (v = 0, off = 0; v < closureSize*2; v += 2) {
674e630c359SToby Isaac                     if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
675e630c359SToby Isaac                       PetscInt    voff;
676e630c359SToby Isaac                       const PetscScalar *xpoint;
677e630c359SToby Isaac 
678e630c359SToby Isaac                       if (nfields) {
679e630c359SToby Isaac                         ierr = PetscSectionGetFieldOffset(section,closure[v],field,&voff);CHKERRQ(ierr);
680e630c359SToby Isaac                       } else {
681e630c359SToby Isaac                         ierr = PetscSectionGetOffset(section,closure[v],&voff);CHKERRQ(ierr);
682e630c359SToby Isaac                       }
683e630c359SToby Isaac                       xpoint         = &x[voff];
684e630c359SToby Isaac                       for (j = 0; j < fbs; j++) {
685955d60d1SToby Isaac                         y[cnt + off++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j]));
686e630c359SToby Isaac                       }
687e630c359SToby Isaac                       for (; j < 3; j++) y[cnt + off++] = 0.;
688e630c359SToby Isaac                     }
689e630c359SToby Isaac                   }
690e630c359SToby Isaac                   cnt += off;
691e630c359SToby Isaac                   ierr = DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
692e630c359SToby Isaac                 }
693e630c359SToby Isaac               }
694e630c359SToby Isaac               if (cnt != piece.nvertices*3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
695955d60d1SToby Isaac               ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices*3,MPIU_VTUREAL,tag);CHKERRQ(ierr);
696955d60d1SToby Isaac             }
697e630c359SToby Isaac           } else {
698e630c359SToby Isaac             for (i=0; i<fbs; i++) {
699955d60d1SToby Isaac               PetscInt cnt, l;
700955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
701e630c359SToby Isaac                 if (!localized) {
702e630c359SToby Isaac                   for (v=vStart,cnt=0; v<vEnd; v++) {
703e630c359SToby Isaac                     PetscInt    off;
704e630c359SToby Isaac                     const PetscScalar *xpoint;
705e630c359SToby Isaac 
706e630c359SToby Isaac                     if (nfields) {
707e630c359SToby Isaac                       ierr     = PetscSectionGetFieldOffset(section,v,field,&off);CHKERRQ(ierr);
708e630c359SToby Isaac                     } else {
709e630c359SToby Isaac                       ierr     = PetscSectionGetOffset(section,v,&off);CHKERRQ(ierr);
710e630c359SToby Isaac                     }
711e630c359SToby Isaac                     xpoint   = &x[off];
712955d60d1SToby Isaac                     y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i]));
713552f7358SJed Brown                   }
7142e529ec8SStefano Zampini                 } else {
7152e529ec8SStefano Zampini                   for (c=cStart,cnt=0; c<cEnd; c++) {
7162e529ec8SStefano Zampini                     PetscInt *closure = NULL;
7172e529ec8SStefano Zampini                     PetscInt  closureSize, off;
7182e529ec8SStefano Zampini 
719e630c359SToby Isaac                     ierr = DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
7202e529ec8SStefano Zampini                     for (v = 0, off = 0; v < closureSize*2; v += 2) {
7212e529ec8SStefano Zampini                       if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
722e630c359SToby Isaac                         PetscInt    voff;
723e630c359SToby Isaac                         const PetscScalar *xpoint;
7242e529ec8SStefano Zampini 
725e630c359SToby Isaac                         if (nfields) {
726e630c359SToby Isaac                           ierr           = PetscSectionGetFieldOffset(section,closure[v],field,&voff);CHKERRQ(ierr);
727e630c359SToby Isaac                         } else {
728e630c359SToby Isaac                           ierr           = PetscSectionGetOffset(section,closure[v],&voff);CHKERRQ(ierr);
729e630c359SToby Isaac                         }
730e630c359SToby Isaac                         xpoint         = &x[voff];
731955d60d1SToby Isaac                         y[cnt + off++] = (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i]));
7322e529ec8SStefano Zampini                       }
7332e529ec8SStefano Zampini                     }
7349bda96f6SStefano Zampini                     cnt += off;
735e630c359SToby Isaac                     ierr = DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
7362e529ec8SStefano Zampini                   }
7372e529ec8SStefano Zampini                 }
738552f7358SJed Brown                 if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
739955d60d1SToby Isaac                 ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices,MPIU_VTUREAL,tag);CHKERRQ(ierr);
740955d60d1SToby Isaac               }
741552f7358SJed Brown             }
742e630c359SToby Isaac           }
743e630c359SToby Isaac         }
744552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
745552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
746552f7358SJed Brown       }
747552f7358SJed Brown     } else if (!rank) {
748955d60d1SToby Isaac       PetscInt l;
749955d60d1SToby Isaac 
750955d60d1SToby Isaac       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,MPIU_VTUREAL,tag);CHKERRQ(ierr); /* positions */
75194fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nconn,MPI_INT,tag);CHKERRQ(ierr); /* connectivity */
75294fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag);CHKERRQ(ierr); /* offsets */
75394fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_CHAR,tag);CHKERRQ(ierr); /* types */
75494fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag);CHKERRQ(ierr); /* owner rank (cells) */
755552f7358SJed Brown       /* all cell data */
756552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
757e630c359SToby Isaac         Vec               X = (Vec)link->vec;
758e630c359SToby Isaac         PetscInt bs, nfields, field;
759e630c359SToby Isaac         DM           dmX;
760e630c359SToby Isaac         PetscSection section = NULL;
761e630c359SToby Isaac 
762552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
763e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
764e630c359SToby Isaac         if (!dmX) dmX = dm;
765e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
766e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
767e630c359SToby Isaac         ierr = PetscSectionGetDof(section,cStart,&bs);CHKERRQ(ierr);
768e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
769e630c359SToby Isaac         field = 0;
770e630c359SToby Isaac         if (link->field >= 0) {
771e630c359SToby Isaac           field = link->field;
772e630c359SToby Isaac           nfields = field + 1;
773e630c359SToby Isaac         }
774e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
775e630c359SToby Isaac           PetscInt     fbs,j;
776e630c359SToby Isaac           PetscFV      fv = NULL;
777e630c359SToby Isaac           PetscObject  f;
778e630c359SToby Isaac           PetscClassId fClass;
779e630c359SToby Isaac           PetscBool    vector;
780e630c359SToby Isaac           if (nfields) {        /* We have user-defined fields/components */
781e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,cStart,field,&fbs);CHKERRQ(ierr);
782e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
783e630c359SToby Isaac           ierr = DMGetField(dmX,field,NULL,&f);CHKERRQ(ierr);
784e630c359SToby Isaac           ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr);
785e630c359SToby Isaac           if (fClass == PETSCFV_CLASSID) {
786e630c359SToby Isaac             fv = (PetscFV) f;
787e630c359SToby Isaac           }
788e630c359SToby Isaac           vector = PETSC_FALSE;
789e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
790e630c359SToby Isaac             vector = PETSC_TRUE;
791e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
792e630c359SToby Isaac               const char *compName = NULL;
793e630c359SToby Isaac               if (fv) {
794e630c359SToby Isaac                 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
795e630c359SToby Isaac                 if (compName) break;
796e630c359SToby Isaac               }
797e630c359SToby Isaac             }
798e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
799e630c359SToby Isaac           }
800e630c359SToby Isaac           if (vector) {
801955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
802955d60d1SToby Isaac               ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells*3,MPIU_VTUREAL,tag);CHKERRQ(ierr);
803955d60d1SToby Isaac             }
804e630c359SToby Isaac           } else {
805e630c359SToby Isaac             for (i=0; i<fbs; i++) {
806955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
807955d60d1SToby Isaac                 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPIU_VTUREAL,tag);CHKERRQ(ierr);
808955d60d1SToby Isaac               }
809552f7358SJed Brown             }
810552f7358SJed Brown           }
811e630c359SToby Isaac         }
812e630c359SToby Isaac       }
813552f7358SJed Brown       /* all point data */
814552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
815e630c359SToby Isaac         Vec               X = (Vec)link->vec;
816e630c359SToby Isaac         DM                dmX;
817e630c359SToby Isaac         PetscInt bs, nfields, field;
818e630c359SToby Isaac         PetscSection section = NULL;
819e630c359SToby Isaac 
820552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
821e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
822e630c359SToby Isaac         if (!dmX) dmX = dm;
823e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
824e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
825e630c359SToby Isaac         ierr = PetscSectionGetDof(section,vStart,&bs);CHKERRQ(ierr);
826e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
827e630c359SToby Isaac         field = 0;
828e630c359SToby Isaac         if (link->field >= 0) {
829e630c359SToby Isaac           field = link->field;
830e630c359SToby Isaac           nfields = field + 1;
831e630c359SToby Isaac         }
832e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
833e630c359SToby Isaac           PetscInt   fbs;
834e630c359SToby Isaac           if (nfields) {        /* We have user-defined fields/components */
835e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,vStart,field,&fbs);CHKERRQ(ierr);
836e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
837e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
838955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
839955d60d1SToby Isaac               ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,MPIU_VTUREAL,tag);CHKERRQ(ierr);
840955d60d1SToby Isaac             }
841e630c359SToby Isaac           } else {
842e630c359SToby Isaac             for (i=0; i<fbs; i++) {
843955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
844955d60d1SToby Isaac                 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices,MPIU_VTUREAL,tag);CHKERRQ(ierr);
845955d60d1SToby Isaac               }
846552f7358SJed Brown             }
847552f7358SJed Brown           }
848552f7358SJed Brown         }
849552f7358SJed Brown       }
850e630c359SToby Isaac     }
851e630c359SToby Isaac   }
852552f7358SJed Brown   ierr = PetscFree(gpiece);CHKERRQ(ierr);
853552f7358SJed Brown   ierr = PetscFree(buffer);CHKERRQ(ierr);
854552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"\n  </AppendedData>\n");CHKERRQ(ierr);
855552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
8565e32de72SMatthew G. Knepley   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
857552f7358SJed Brown   PetscFunctionReturn(0);
858552f7358SJed Brown }
859