xref: /petsc/src/dm/impls/plex/plexvtu.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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 
266030a318SStefano Zampini static PetscErrorCode TransferWrite(MPI_Comm comm, 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;
30552f7358SJed Brown 
31552f7358SJed Brown   PetscFunctionBegin;
32ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
33552f7358SJed Brown   if (rank == srank && rank != root) {
34ffc4695bSBarry Smith     ierr = MPI_Send((void*)send,count,mpidatatype,root,tag,comm);CHKERRMPI(ierr);
35552f7358SJed Brown   } else if (rank == root) {
36552f7358SJed Brown     const void *buffer;
37552f7358SJed Brown     if (root == srank) {        /* self */
38552f7358SJed Brown       buffer = send;
39552f7358SJed Brown     } else {
40552f7358SJed Brown       MPI_Status  status;
41552f7358SJed Brown       PetscMPIInt nrecv;
42ffc4695bSBarry Smith       ierr = MPI_Recv(recv,count,mpidatatype,srank,tag,comm,&status);CHKERRMPI(ierr);
43ffc4695bSBarry Smith       ierr = MPI_Get_count(&status,mpidatatype,&nrecv);CHKERRMPI(ierr);
44*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(count != nrecv,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
45552f7358SJed Brown       buffer = recv;
46552f7358SJed Brown     }
4794fbd55eSBarry Smith     ierr = PetscViewerVTKFWrite(viewer,fp,buffer,count,mpidatatype);CHKERRQ(ierr);
48552f7358SJed Brown   }
49552f7358SJed Brown   PetscFunctionReturn(0);
50552f7358SJed Brown }
51552f7358SJed Brown 
522e529ec8SStefano Zampini static PetscErrorCode DMPlexGetVTKConnectivity(DM dm, PetscBool localized, PieceInfo *piece,PetscVTKInt **oconn,PetscVTKInt **ooffsets,PetscVTKType **otypes)
53552f7358SJed Brown {
54552f7358SJed Brown   PetscErrorCode ierr;
552e529ec8SStefano Zampini   PetscSection   coordSection;
56552f7358SJed Brown   PetscVTKInt    *conn,*offsets;
57552f7358SJed Brown   PetscVTKType   *types;
582f029394SStefano Zampini   PetscInt       dim,vStart,vEnd,cStart,cEnd,pStart,pEnd,cellHeight,numLabelCells,hasLabel,c,v,countcell,countconn;
59552f7358SJed Brown 
60552f7358SJed Brown   PetscFunctionBegin;
61dcca6d9dSJed Brown   ierr = PetscMalloc3(piece->nconn,&conn,piece->ncells,&offsets,piece->ncells,&types);CHKERRQ(ierr);
622e529ec8SStefano Zampini   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
63c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
64552f7358SJed Brown   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
65552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
66552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
67552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
68c58f1c22SToby Isaac   ierr     = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
69552f7358SJed Brown   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
70552f7358SJed Brown 
71552f7358SJed Brown   countcell = 0;
72552f7358SJed Brown   countconn = 0;
73552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
7419003fb5SStefano Zampini     PetscInt nverts,dof = 0,celltype,startoffset,nC=0;
75552f7358SJed Brown 
76552f7358SJed Brown     if (hasLabel) {
77552f7358SJed Brown       PetscInt value;
78552f7358SJed Brown 
79c58f1c22SToby Isaac       ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
80552f7358SJed Brown       if (value != 1) continue;
81552f7358SJed Brown     }
82552f7358SJed Brown     startoffset = countconn;
8319003fb5SStefano Zampini     if (localized) {
8419003fb5SStefano Zampini       ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
8519003fb5SStefano Zampini     }
8619003fb5SStefano Zampini     if (!dof) {
872e529ec8SStefano Zampini       PetscInt *closure = NULL;
882e529ec8SStefano Zampini       PetscInt  closureSize;
892e529ec8SStefano Zampini 
90552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
91552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
92552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
9319003fb5SStefano Zampini           if (!localized) conn[countconn++] = closure[v] - vStart;
9419003fb5SStefano Zampini           else conn[countconn++] = startoffset + nC;
95724b5320SMatthew G. Knepley           ++nC;
96552f7358SJed Brown         }
97552f7358SJed Brown       }
98552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
992e529ec8SStefano Zampini     } else {
1002e529ec8SStefano Zampini       for (nC = 0; nC < dof/dim; nC++) conn[countconn++] = startoffset + nC;
1012e529ec8SStefano Zampini     }
10296ca5757SLisandro Dalcin 
10396ca5757SLisandro Dalcin     {
10496ca5757SLisandro Dalcin       PetscInt n = PetscMin(nC, 8), s = countconn - nC, i, cone[8];
10596ca5757SLisandro Dalcin       for (i = 0; i < n; ++i) cone[i] = conn[s+i];
10696ca5757SLisandro Dalcin       ierr = DMPlexReorderCell(dm, c, cone);CHKERRQ(ierr);
10796ca5757SLisandro Dalcin       for (i = 0; i < n; ++i) conn[s+i] = (int)cone[i];
10896ca5757SLisandro Dalcin     }
1098865f1eaSKarl Rupp 
110552f7358SJed Brown     offsets[countcell] = countconn;
1118865f1eaSKarl Rupp 
112552f7358SJed Brown     nverts = countconn - startoffset;
113de0a4605SMatthew G. Knepley     ierr   = DMPlexVTKGetCellType_Internal(dm,dim,nverts,&celltype);CHKERRQ(ierr);
1148865f1eaSKarl Rupp 
115552f7358SJed Brown     types[countcell] = celltype;
116552f7358SJed Brown     countcell++;
117552f7358SJed Brown   }
118*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(countcell != piece->ncells,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count");
119*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(countconn != piece->nconn,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count");
120552f7358SJed Brown   *oconn    = conn;
121552f7358SJed Brown   *ooffsets = offsets;
122552f7358SJed Brown   *otypes   = types;
123552f7358SJed Brown   PetscFunctionReturn(0);
124552f7358SJed Brown }
125552f7358SJed Brown 
126552f7358SJed Brown /*
127552f7358SJed Brown   Write all fields that have been provided to the viewer
128552f7358SJed Brown   Multi-block XML format with binary appended data.
129552f7358SJed Brown */
130552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer)
131552f7358SJed Brown {
132ce94432eSBarry Smith   MPI_Comm                 comm;
1332e529ec8SStefano Zampini   PetscSection             coordSection;
134552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
135552f7358SJed Brown   PetscViewerVTKObjectLink link;
136552f7358SJed Brown   FILE                     *fp;
137552f7358SJed Brown   PetscMPIInt              rank,size,tag;
138552f7358SJed Brown   PetscErrorCode           ierr;
1392f029394SStefano Zampini   PetscInt                 dimEmbed,cellHeight,cStart,cEnd,vStart,vEnd,numLabelCells,hasLabel,c,v,r,i;
1402e529ec8SStefano Zampini   PetscBool                localized;
1410298fd71SBarry Smith   PieceInfo                piece,*gpiece = NULL;
1420298fd71SBarry Smith   void                     *buffer = NULL;
14330815ce0SLisandro Dalcin   const char               *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian";
144955d60d1SToby Isaac   PetscInt                 loops_per_scalar;
145552f7358SJed Brown 
146552f7358SJed Brown   PetscFunctionBegin;
147ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
148552f7358SJed Brown #if defined(PETSC_USE_COMPLEX)
149955d60d1SToby Isaac   loops_per_scalar = 2;
150955d60d1SToby Isaac #else
151955d60d1SToby Isaac   loops_per_scalar = 1;
152552f7358SJed Brown #endif
153ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
154ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
1556030a318SStefano Zampini   ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr);
156552f7358SJed Brown 
157552f7358SJed Brown   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
158552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
15930815ce0SLisandro Dalcin   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"%s\">\n", byte_order);CHKERRQ(ierr);
160552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <UnstructuredGrid>\n");CHKERRQ(ierr);
161552f7358SJed Brown 
1623baa072aSToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
163552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
164552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
165552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
166c58f1c22SToby Isaac   ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
1672e529ec8SStefano Zampini   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
1682e529ec8SStefano Zampini   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1698865f1eaSKarl Rupp 
170552f7358SJed Brown   hasLabel        = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
1712e529ec8SStefano Zampini   piece.nvertices = 0;
172552f7358SJed Brown   piece.ncells    = 0;
173552f7358SJed Brown   piece.nconn     = 0;
1742e529ec8SStefano Zampini   if (!localized) piece.nvertices = vEnd - vStart;
175552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
17619003fb5SStefano Zampini     PetscInt dof = 0;
177552f7358SJed Brown     if (hasLabel) {
178552f7358SJed Brown       PetscInt value;
179552f7358SJed Brown 
180c58f1c22SToby Isaac       ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
181552f7358SJed Brown       if (value != 1) continue;
182552f7358SJed Brown     }
18319003fb5SStefano Zampini     if (localized) {
18419003fb5SStefano Zampini       ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
18519003fb5SStefano Zampini     }
18619003fb5SStefano Zampini     if (!dof) {
1872e529ec8SStefano Zampini       PetscInt *closure = NULL;
1882e529ec8SStefano Zampini       PetscInt closureSize;
1892e529ec8SStefano Zampini 
190552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
191552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
1922e529ec8SStefano Zampini         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
1932e529ec8SStefano Zampini           piece.nconn++;
19419003fb5SStefano Zampini           if (localized) piece.nvertices++;
1952e529ec8SStefano Zampini         }
196552f7358SJed Brown       }
197552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1982e529ec8SStefano Zampini     } else {
1992e529ec8SStefano Zampini       piece.nvertices += dof/dimEmbed;
2002e529ec8SStefano Zampini       piece.nconn     += dof/dimEmbed;
2012e529ec8SStefano Zampini     }
202552f7358SJed Brown     piece.ncells++;
203552f7358SJed Brown   }
204dd400576SPatrick Sanan   if (rank == 0) {ierr = PetscMalloc1(size,&gpiece);CHKERRQ(ierr);}
205ffc4695bSBarry Smith   ierr = MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm);CHKERRMPI(ierr);
206552f7358SJed Brown 
207552f7358SJed Brown   /*
208552f7358SJed Brown    * Write file header
209552f7358SJed Brown    */
210dd400576SPatrick Sanan   if (rank == 0) {
211552f7358SJed Brown     PetscInt boffset = 0;
212552f7358SJed Brown 
213552f7358SJed Brown     for (r=0; r<size; r++) {
214552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells);CHKERRQ(ierr);
215552f7358SJed Brown       /* Coordinate positions */
216552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Points>\n");CHKERRQ(ierr);
217552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
218955d60d1SToby Isaac       boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int);
219552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Points>\n");CHKERRQ(ierr);
220552f7358SJed Brown       /* Cell connectivity */
221552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Cells>\n");CHKERRQ(ierr);
222552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
223552f7358SJed Brown       boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int);
224552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"offsets\"      NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
225552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int);
226552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"UInt8\" Name=\"types\"        NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
227552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int);
228552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Cells>\n");CHKERRQ(ierr);
229552f7358SJed Brown 
230552f7358SJed Brown       /*
231552f7358SJed Brown        * Cell Data headers
232552f7358SJed Brown        */
233552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <CellData>\n");CHKERRQ(ierr);
234552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
235552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(int) + sizeof(int);
236552f7358SJed Brown       /* all the vectors */
237552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
238552f7358SJed Brown         Vec          X = (Vec)link->vec;
239e630c359SToby Isaac         DM           dmX = NULL;
2406030a318SStefano Zampini         PetscInt     bs = 1,nfields,field;
241552f7358SJed Brown         const char   *vecname = "";
242e630c359SToby Isaac         PetscSection section;
243552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
244552f7358SJed 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. */
245552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
246552f7358SJed Brown         }
247e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
248e630c359SToby Isaac         if (!dmX) dmX = dm;
249e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
250e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
2516030a318SStefano Zampini         if (cEnd > cStart) { ierr = PetscSectionGetDof(section,cStart,&bs);CHKERRQ(ierr); }
252e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
253e630c359SToby Isaac         field = 0;
254e630c359SToby Isaac         if (link->field >= 0) {
255e630c359SToby Isaac           field = link->field;
256e630c359SToby Isaac           nfields = field + 1;
257e630c359SToby Isaac         }
258e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
2591cfafdd3SJed Brown           PetscInt     fbs,j;
260a00cdb45SToby Isaac           PetscFV      fv = NULL;
261a00cdb45SToby Isaac           PetscObject  f;
262a00cdb45SToby Isaac           PetscClassId fClass;
2630298fd71SBarry Smith           const char *fieldname = NULL;
2641cfafdd3SJed Brown           char       buf[256];
265e630c359SToby Isaac           PetscBool    vector;
2667ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
267e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,cStart,field,&fbs);CHKERRQ(ierr);
268e630c359SToby Isaac             ierr = PetscSectionGetFieldName(section,field,&fieldname);CHKERRQ(ierr);
2697ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
270e630c359SToby Isaac           ierr = DMGetField(dmX,field,NULL,&f);CHKERRQ(ierr);
271a00cdb45SToby Isaac           ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr);
272a00cdb45SToby Isaac           if (fClass == PETSCFV_CLASSID) {
273a00cdb45SToby Isaac             fv = (PetscFV) f;
274a00cdb45SToby Isaac           }
275e630c359SToby Isaac           if (nfields && !fieldname) {
2761cfafdd3SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"CellField%D",field);CHKERRQ(ierr);
277552f7358SJed Brown             fieldname = buf;
278552f7358SJed Brown           }
279e630c359SToby Isaac           vector = PETSC_FALSE;
280e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
281e630c359SToby Isaac             vector = PETSC_TRUE;
282*2c71b3e2SJacob Faibussowitsch             PetscCheckFalse(fbs > 3,PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_SIZ,"Cell vector fields can have at most 3 components, %D given", fbs);
283e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
284e630c359SToby Isaac               const char *compName = NULL;
285e630c359SToby Isaac               if (fv) {
286e630c359SToby Isaac                 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
287e630c359SToby Isaac                 if (compName) break;
288e630c359SToby Isaac               }
289e630c359SToby Isaac             }
290e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
291e630c359SToby Isaac           }
292e630c359SToby Isaac           if (vector) {
293955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
294955d60d1SToby 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);
295955d60d1SToby Isaac             boffset += gpiece[r].ncells*3*sizeof(PetscVTUReal) + sizeof(int);
296955d60d1SToby 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);
297955d60d1SToby Isaac             boffset += gpiece[r].ncells*3*sizeof(PetscVTUReal) + sizeof(int);
298955d60d1SToby Isaac #else
299e630c359SToby Isaac             ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" 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 #endif
302e630c359SToby Isaac             i+=fbs;
303e630c359SToby Isaac           } else {
3041cfafdd3SJed Brown             for (j=0; j<fbs; j++) {
305a00cdb45SToby Isaac               const char *compName = NULL;
306955d60d1SToby Isaac               char finalname[256];
307a00cdb45SToby Isaac               if (fv) {
308a00cdb45SToby Isaac                 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
309a00cdb45SToby Isaac               }
310a00cdb45SToby Isaac               if (compName) {
311955d60d1SToby Isaac                 ierr = PetscSNPrintf(finalname,255,"%s%s.%s",vecname,fieldname,compName);CHKERRQ(ierr);
312a00cdb45SToby Isaac               }
313e630c359SToby Isaac               else if (fbs > 1) {
314955d60d1SToby Isaac                 ierr = PetscSNPrintf(finalname,255,"%s%s.%D",vecname,fieldname,j);CHKERRQ(ierr);
315e630c359SToby Isaac               } else {
316955d60d1SToby Isaac                 ierr = PetscSNPrintf(finalname,255,"%s%s",vecname,fieldname);CHKERRQ(ierr);
317a00cdb45SToby Isaac               }
318955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
319955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
320955d60d1SToby Isaac               boffset += gpiece[r].ncells*sizeof(PetscVTUReal) + sizeof(int);
321955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
322955d60d1SToby Isaac               boffset += gpiece[r].ncells*sizeof(PetscVTUReal) + sizeof(int);
323955d60d1SToby Isaac #else
324955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
325955d60d1SToby Isaac               boffset += gpiece[r].ncells*sizeof(PetscVTUReal) + sizeof(int);
326955d60d1SToby Isaac #endif
3271cfafdd3SJed Brown               i++;
328552f7358SJed Brown             }
329552f7358SJed Brown           }
330e630c359SToby Isaac         }
331*2c71b3e2SJacob Faibussowitsch         //PetscCheckFalse(i != bs,comm,PETSC_ERR_PLIB,"Total number of field components %D != block size %D",i,bs);
3321cfafdd3SJed Brown       }
333552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </CellData>\n");CHKERRQ(ierr);
334552f7358SJed Brown 
335552f7358SJed Brown       /*
336552f7358SJed Brown        * Point Data headers
337552f7358SJed Brown        */
338552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      <PointData>\n");CHKERRQ(ierr);
339552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
340552f7358SJed Brown         Vec          X = (Vec)link->vec;
341e630c359SToby Isaac         DM           dmX;
3426030a318SStefano Zampini         PetscInt     bs = 1,nfields,field;
343552f7358SJed Brown         const char   *vecname = "";
344e630c359SToby Isaac         PetscSection section;
345552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
346552f7358SJed 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. */
347552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
348552f7358SJed Brown         }
349e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
350e630c359SToby Isaac         if (!dmX) dmX = dm;
351e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
352e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
3536030a318SStefano Zampini         if (vEnd > vStart) { ierr = PetscSectionGetDof(section,vStart,&bs);CHKERRQ(ierr); }
354e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
355e630c359SToby Isaac         field = 0;
356e630c359SToby Isaac         if (link->field >= 0) {
357e630c359SToby Isaac           field = link->field;
358e630c359SToby Isaac           nfields = field + 1;
359e630c359SToby Isaac         }
360e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
3611cfafdd3SJed Brown           PetscInt   fbs,j;
3621cfafdd3SJed Brown           const char *fieldname = NULL;
3631cfafdd3SJed Brown           char       buf[256];
3647ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
365e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,vStart,field,&fbs);CHKERRQ(ierr);
366e630c359SToby Isaac             ierr = PetscSectionGetFieldName(section,field,&fieldname);CHKERRQ(ierr);
3677ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
368e630c359SToby Isaac           if (nfields && !fieldname) {
3691cfafdd3SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"PointField%D",field);CHKERRQ(ierr);
3701cfafdd3SJed Brown             fieldname = buf;
3711cfafdd3SJed Brown           }
372e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
373*2c71b3e2SJacob Faibussowitsch             PetscCheckFalse(fbs > 3,PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_SIZ,"Point vector fields can have at most 3 components, %D given", fbs);
374955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
375955d60d1SToby 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);
376955d60d1SToby Isaac             boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int);
377955d60d1SToby 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);
378955d60d1SToby Isaac             boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int);
379955d60d1SToby Isaac #else
380e630c359SToby Isaac             ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" 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 #endif
383e630c359SToby Isaac           } else {
3841cfafdd3SJed Brown             for (j=0; j<fbs; j++) {
385b778fa18SValeria Barra               const char *compName = NULL;
386955d60d1SToby Isaac               char finalname[256];
387b778fa18SValeria Barra               ierr = PetscSectionGetComponentName(section,field,j,&compName);CHKERRQ(ierr);
388b778fa18SValeria Barra               ierr = PetscSNPrintf(finalname,255,"%s%s.%s",vecname,fieldname,compName);CHKERRQ(ierr);
389955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
390955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
391955d60d1SToby Isaac               boffset += gpiece[r].nvertices*sizeof(PetscVTUReal) + sizeof(int);
392955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
393955d60d1SToby Isaac               boffset += gpiece[r].nvertices*sizeof(PetscVTUReal) + sizeof(int);
394955d60d1SToby Isaac #else
395955d60d1SToby Isaac               ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset);CHKERRQ(ierr);
396955d60d1SToby Isaac               boffset += gpiece[r].nvertices*sizeof(PetscVTUReal) + sizeof(int);
397955d60d1SToby Isaac #endif
398552f7358SJed Brown             }
399552f7358SJed Brown           }
4001cfafdd3SJed Brown         }
401e630c359SToby Isaac       }
402552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </PointData>\n");CHKERRQ(ierr);
403552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    </Piece>\n");CHKERRQ(ierr);
404552f7358SJed Brown     }
405552f7358SJed Brown   }
406552f7358SJed Brown 
407552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  </UnstructuredGrid>\n");CHKERRQ(ierr);
408552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
409552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
410552f7358SJed Brown 
411dd400576SPatrick Sanan   if (rank == 0) {
412552f7358SJed Brown     PetscInt maxsize = 0;
413552f7358SJed Brown     for (r=0; r<size; r++) {
414955d60d1SToby Isaac       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscVTUReal)));
415955d60d1SToby Isaac       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*3*sizeof(PetscVTUReal)));
416552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt)));
417552f7358SJed Brown     }
418552f7358SJed Brown     ierr = PetscMalloc(maxsize,&buffer);CHKERRQ(ierr);
419552f7358SJed Brown   }
420552f7358SJed Brown   for (r=0; r<size; r++) {
421552f7358SJed Brown     if (r == rank) {
422552f7358SJed Brown       PetscInt nsend;
423552f7358SJed Brown       {                         /* Position */
424552f7358SJed Brown         const PetscScalar *x;
425955d60d1SToby Isaac         PetscVTUReal      *y = NULL;
426552f7358SJed Brown         Vec               coords;
427955d60d1SToby Isaac         PetscBool         copy;
4282e529ec8SStefano Zampini 
429552f7358SJed Brown         ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
430552f7358SJed Brown         ierr = VecGetArrayRead(coords,&x);CHKERRQ(ierr);
431955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
432955d60d1SToby Isaac         copy = PETSC_TRUE;
433955d60d1SToby Isaac #else
434955d60d1SToby Isaac         copy = (PetscBool) (dimEmbed != 3 || localized || (sizeof(PetscReal) != sizeof(PetscVTUReal)));
435955d60d1SToby Isaac #endif
436955d60d1SToby Isaac         if (copy) {
437785e854fSJed Brown           ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr);
43819003fb5SStefano Zampini           if (localized) {
43919003fb5SStefano Zampini             PetscInt cnt;
44019003fb5SStefano Zampini             for (c=cStart,cnt=0; c<cEnd; c++) {
44119003fb5SStefano Zampini               PetscInt off, dof;
44219003fb5SStefano Zampini 
44319003fb5SStefano Zampini               ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
44419003fb5SStefano Zampini               if (!dof) {
44519003fb5SStefano Zampini                 PetscInt *closure = NULL;
44619003fb5SStefano Zampini                 PetscInt closureSize;
44719003fb5SStefano Zampini 
44819003fb5SStefano Zampini                 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
44919003fb5SStefano Zampini                 for (v = 0; v < closureSize*2; v += 2) {
45019003fb5SStefano Zampini                   if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
45119003fb5SStefano Zampini                     ierr = PetscSectionGetOffset(coordSection, closure[v], &off);CHKERRQ(ierr);
45219003fb5SStefano Zampini                     if (dimEmbed != 3) {
453955d60d1SToby Isaac                       y[cnt*3+0] = (PetscVTUReal) PetscRealPart(x[off+0]);
454955d60d1SToby Isaac                       y[cnt*3+1] = (PetscVTUReal) ((dimEmbed > 1) ? PetscRealPart(x[off+1]) : 0.0);
455955d60d1SToby Isaac                       y[cnt*3+2] = (PetscVTUReal) 0.0;
45619003fb5SStefano Zampini                     } else {
457955d60d1SToby Isaac                       y[cnt*3+0] = (PetscVTUReal) PetscRealPart(x[off+0]);
458955d60d1SToby Isaac                       y[cnt*3+1] = (PetscVTUReal) PetscRealPart(x[off+1]);
459955d60d1SToby Isaac                       y[cnt*3+2] = (PetscVTUReal) PetscRealPart(x[off+2]);
46019003fb5SStefano Zampini                     }
46119003fb5SStefano Zampini                     cnt++;
46219003fb5SStefano Zampini                   }
46319003fb5SStefano Zampini                 }
46419003fb5SStefano Zampini                 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
46519003fb5SStefano Zampini               } else {
46619003fb5SStefano Zampini                 ierr = PetscSectionGetOffset(coordSection, c, &off);CHKERRQ(ierr);
46719003fb5SStefano Zampini                 if (dimEmbed != 3) {
46819003fb5SStefano Zampini                   for (i=0; i<dof/dimEmbed; i++) {
469955d60d1SToby Isaac                     y[cnt*3+0] = (PetscVTUReal) PetscRealPart(x[off + i*dimEmbed + 0]);
470955d60d1SToby Isaac                     y[cnt*3+1] = (PetscVTUReal) ((dimEmbed > 1) ? PetscRealPart(x[off + i*dimEmbed + 1]) : 0.0);
471955d60d1SToby Isaac                     y[cnt*3+2] = (PetscVTUReal) 0.0;
47219003fb5SStefano Zampini                     cnt++;
47319003fb5SStefano Zampini                   }
47419003fb5SStefano Zampini                 } else {
47519003fb5SStefano Zampini                   for (i=0; i<dof; i ++) {
476955d60d1SToby Isaac                     y[cnt*3+i] = (PetscVTUReal) PetscRealPart(x[off + i]);
47719003fb5SStefano Zampini                   }
47819003fb5SStefano Zampini                   cnt += dof/dimEmbed;
47919003fb5SStefano Zampini                 }
48019003fb5SStefano Zampini               }
48119003fb5SStefano Zampini             }
482*2c71b3e2SJacob Faibussowitsch             PetscCheckFalse(cnt != piece.nvertices,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
48319003fb5SStefano Zampini           } else {
484552f7358SJed Brown             for (i=0; i<piece.nvertices; i++) {
485955d60d1SToby Isaac               y[i*3+0] = (PetscVTUReal) PetscRealPart(x[i*dimEmbed+0]);
486955d60d1SToby Isaac               y[i*3+1] = (PetscVTUReal) ((dimEmbed > 1) ? PetscRealPart(x[i*dimEmbed+1]) : 0.);
487955d60d1SToby Isaac               y[i*3+2] = (PetscVTUReal) ((dimEmbed > 2) ? PetscRealPart(x[i*dimEmbed+2]) : 0.);
48819003fb5SStefano Zampini             }
489552f7358SJed Brown           }
490552f7358SJed Brown         }
4912e529ec8SStefano Zampini         nsend = piece.nvertices*3;
4926030a318SStefano Zampini         ierr  = TransferWrite(comm,viewer,fp,r,0,copy ? (const void *) y : (const void *) x,buffer,nsend,MPIU_VTUREAL,tag);CHKERRQ(ierr);
493552f7358SJed Brown         ierr  = PetscFree(y);CHKERRQ(ierr);
494552f7358SJed Brown         ierr  = VecRestoreArrayRead(coords,&x);CHKERRQ(ierr);
495552f7358SJed Brown       }
496552f7358SJed Brown       {                           /* Connectivity, offsets, types */
4973e869605SMatthew G. Knepley         PetscVTKInt  *connectivity = NULL, *offsets = NULL;
4983e869605SMatthew G. Knepley         PetscVTKType *types = NULL;
4992e529ec8SStefano Zampini         ierr = DMPlexGetVTKConnectivity(dm,localized,&piece,&connectivity,&offsets,&types);CHKERRQ(ierr);
5006030a318SStefano Zampini         ierr = TransferWrite(comm,viewer,fp,r,0,connectivity,buffer,piece.nconn,MPI_INT,tag);CHKERRQ(ierr);
5016030a318SStefano Zampini         ierr = TransferWrite(comm,viewer,fp,r,0,offsets,buffer,piece.ncells,MPI_INT,tag);CHKERRQ(ierr);
5026030a318SStefano Zampini         ierr = TransferWrite(comm,viewer,fp,r,0,types,buffer,piece.ncells,MPI_CHAR,tag);CHKERRQ(ierr);
503552f7358SJed Brown         ierr = PetscFree3(connectivity,offsets,types);CHKERRQ(ierr);
504552f7358SJed Brown       }
505552f7358SJed Brown       {                         /* Owners (cell data) */
506552f7358SJed Brown         PetscVTKInt *owners;
507785e854fSJed Brown         ierr = PetscMalloc1(piece.ncells,&owners);CHKERRQ(ierr);
508552f7358SJed Brown         for (i=0; i<piece.ncells; i++) owners[i] = rank;
5096030a318SStefano Zampini         ierr = TransferWrite(comm,viewer,fp,r,0,owners,buffer,piece.ncells,MPI_INT,tag);CHKERRQ(ierr);
510552f7358SJed Brown         ierr = PetscFree(owners);CHKERRQ(ierr);
511552f7358SJed Brown       }
512552f7358SJed Brown       /* Cell data */
513552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
514552f7358SJed Brown         Vec               X = (Vec)link->vec;
515e630c359SToby Isaac         DM                dmX;
516552f7358SJed Brown         const PetscScalar *x;
517955d60d1SToby Isaac         PetscVTUReal      *y;
5186030a318SStefano Zampini         PetscInt          bs = 1, nfields, field;
519e630c359SToby Isaac         PetscSection      section = NULL;
520e630c359SToby Isaac 
521552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
522e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
523e630c359SToby Isaac         if (!dmX) dmX = dm;
524e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
525e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
5266030a318SStefano Zampini         if (cEnd > cStart) { ierr = PetscSectionGetDof(section,cStart,&bs);CHKERRQ(ierr); }
527e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
528e630c359SToby Isaac         field = 0;
529e630c359SToby Isaac         if (link->field >= 0) {
530e630c359SToby Isaac           field = link->field;
531e630c359SToby Isaac           nfields = field + 1;
532e630c359SToby Isaac         }
533552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
534e630c359SToby Isaac         ierr = PetscMalloc1(piece.ncells*3,&y);CHKERRQ(ierr);
535e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
536e630c359SToby Isaac           PetscInt     fbs,j;
537e630c359SToby Isaac           PetscFV      fv = NULL;
538e630c359SToby Isaac           PetscObject  f;
539e630c359SToby Isaac           PetscClassId fClass;
540e630c359SToby Isaac           PetscBool    vector;
5416030a318SStefano Zampini           if (nfields && cEnd > cStart) {        /* We have user-defined fields/components */
542e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,cStart,field,&fbs);CHKERRQ(ierr);
543e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
544e630c359SToby Isaac           ierr = DMGetField(dmX,field,NULL,&f);CHKERRQ(ierr);
545e630c359SToby Isaac           ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr);
546e630c359SToby Isaac           if (fClass == PETSCFV_CLASSID) {
547e630c359SToby Isaac             fv = (PetscFV) f;
548e630c359SToby Isaac           }
549e630c359SToby Isaac           vector = PETSC_FALSE;
550e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
551e630c359SToby Isaac             vector = PETSC_TRUE;
552e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
553e630c359SToby Isaac               const char *compName = NULL;
554e630c359SToby Isaac               if (fv) {
555e630c359SToby Isaac                 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
556e630c359SToby Isaac                 if (compName) break;
557e630c359SToby Isaac               }
558e630c359SToby Isaac             }
559e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
560e630c359SToby Isaac           }
561e630c359SToby Isaac           if (vector) {
562955d60d1SToby Isaac             PetscInt cnt, l;
563955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
564552f7358SJed Brown               for (c=cStart,cnt=0; c<cEnd; c++) {
565e630c359SToby Isaac                 const PetscScalar *xpoint;
566e630c359SToby Isaac                 PetscInt off, j;
567e630c359SToby Isaac 
568552f7358SJed Brown                 if (hasLabel) {     /* Ignore some cells */
569552f7358SJed Brown                   PetscInt value;
570e630c359SToby Isaac                   ierr = DMGetLabelValue(dmX, "vtk", c, &value);CHKERRQ(ierr);
571552f7358SJed Brown                   if (value != 1) continue;
572552f7358SJed Brown                 }
573e630c359SToby Isaac                 if (nfields) {
574e630c359SToby Isaac                   ierr = PetscSectionGetFieldOffset(section, c, field, &off);CHKERRQ(ierr);
575e630c359SToby Isaac                 } else {
576e630c359SToby Isaac                   ierr = PetscSectionGetOffset(section, c, &off);CHKERRQ(ierr);
577e630c359SToby Isaac                 }
578e630c359SToby Isaac                 xpoint = &x[off];
579e630c359SToby Isaac                 for (j = 0; j < fbs; j++) {
580955d60d1SToby Isaac                   y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j]));
581e630c359SToby Isaac                 }
582e630c359SToby Isaac                 for (; j < 3; j++) y[cnt++] = 0.;
583e630c359SToby Isaac               }
584*2c71b3e2SJacob Faibussowitsch               PetscCheckFalse(cnt != piece.ncells*3,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
5856030a318SStefano Zampini               ierr = TransferWrite(comm,viewer,fp,r,0,y,buffer,piece.ncells*3,MPIU_VTUREAL,tag);CHKERRQ(ierr);
586955d60d1SToby Isaac             }
587e630c359SToby Isaac           } else {
588e630c359SToby Isaac             for (i=0; i<fbs; i++) {
589955d60d1SToby Isaac               PetscInt cnt, l;
590955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
591e630c359SToby Isaac                 for (c=cStart,cnt=0; c<cEnd; c++) {
592e630c359SToby Isaac                   const PetscScalar *xpoint;
593e630c359SToby Isaac                   PetscInt off;
594e630c359SToby Isaac 
595e630c359SToby Isaac                   if (hasLabel) {     /* Ignore some cells */
596e630c359SToby Isaac                     PetscInt value;
597e630c359SToby Isaac                     ierr = DMGetLabelValue(dmX, "vtk", c, &value);CHKERRQ(ierr);
598e630c359SToby Isaac                     if (value != 1) continue;
599e630c359SToby Isaac                   }
600e630c359SToby Isaac                   if (nfields) {
601e630c359SToby Isaac                     ierr = PetscSectionGetFieldOffset(section, c, field, &off);CHKERRQ(ierr);
602e630c359SToby Isaac                   } else {
603e630c359SToby Isaac                     ierr = PetscSectionGetOffset(section, c, &off);CHKERRQ(ierr);
604e630c359SToby Isaac                   }
605e630c359SToby Isaac                   xpoint   = &x[off];
606955d60d1SToby Isaac                   y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i]));
607552f7358SJed Brown                 }
608*2c71b3e2SJacob Faibussowitsch                 PetscCheckFalse(cnt != piece.ncells,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
6096030a318SStefano Zampini                 ierr = TransferWrite(comm,viewer,fp,r,0,y,buffer,piece.ncells,MPIU_VTUREAL,tag);CHKERRQ(ierr);
610955d60d1SToby Isaac               }
611552f7358SJed Brown             }
612e630c359SToby Isaac           }
613e630c359SToby Isaac         }
614552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
615552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
616552f7358SJed Brown       }
617e630c359SToby Isaac       /* point data */
618552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
619552f7358SJed Brown         Vec               X = (Vec)link->vec;
620e630c359SToby Isaac         DM                dmX;
621552f7358SJed Brown         const PetscScalar *x;
622955d60d1SToby Isaac         PetscVTUReal      *y;
6236030a318SStefano Zampini         PetscInt          bs = 1, nfields, field;
624e630c359SToby Isaac         PetscSection      section = NULL;
625e630c359SToby Isaac 
626552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
627e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
628e630c359SToby Isaac         if (!dmX) dmX = dm;
629e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
630e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
6316030a318SStefano Zampini         if (vEnd > vStart) { ierr = PetscSectionGetDof(section,vStart,&bs);CHKERRQ(ierr); }
632e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
633e630c359SToby Isaac         field = 0;
634e630c359SToby Isaac         if (link->field >= 0) {
635e630c359SToby Isaac           field = link->field;
636e630c359SToby Isaac           nfields = field + 1;
637e630c359SToby Isaac         }
638552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
639e630c359SToby Isaac         ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr);
640e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
641e630c359SToby Isaac           PetscInt fbs,j;
6426030a318SStefano Zampini           if (nfields && vEnd > vStart) {        /* We have user-defined fields/components */
643e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,vStart,field,&fbs);CHKERRQ(ierr);
644e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
645e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
646955d60d1SToby Isaac             PetscInt cnt, l;
647955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
6482e529ec8SStefano Zampini               if (!localized) {
649552f7358SJed Brown                 for (v=vStart,cnt=0; v<vEnd; v++) {
650e630c359SToby Isaac                   PetscInt          off;
651e630c359SToby Isaac                   const PetscScalar *xpoint;
652e630c359SToby Isaac 
653e630c359SToby Isaac                   if (nfields) {
654e630c359SToby Isaac                     ierr = PetscSectionGetFieldOffset(section,v,field,&off);CHKERRQ(ierr);
655e630c359SToby Isaac                   } else {
656e630c359SToby Isaac                     ierr = PetscSectionGetOffset(section,v,&off);CHKERRQ(ierr);
657e630c359SToby Isaac                   }
658e630c359SToby Isaac                   xpoint = &x[off];
659e630c359SToby Isaac                   for (j = 0; j < fbs; j++) {
660955d60d1SToby Isaac                     y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j]));
661e630c359SToby Isaac                   }
662e630c359SToby Isaac                   for (; j < 3; j++) y[cnt++] = 0.;
663e630c359SToby Isaac                 }
664e630c359SToby Isaac               } else {
665e630c359SToby Isaac                 for (c=cStart,cnt=0; c<cEnd; c++) {
666e630c359SToby Isaac                   PetscInt *closure = NULL;
667e630c359SToby Isaac                   PetscInt  closureSize, off;
668e630c359SToby Isaac 
669e630c359SToby Isaac                   ierr = DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
670e630c359SToby Isaac                   for (v = 0, off = 0; v < closureSize*2; v += 2) {
671e630c359SToby Isaac                     if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
672e630c359SToby Isaac                       PetscInt    voff;
673e630c359SToby Isaac                       const PetscScalar *xpoint;
674e630c359SToby Isaac 
675e630c359SToby Isaac                       if (nfields) {
676e630c359SToby Isaac                         ierr = PetscSectionGetFieldOffset(section,closure[v],field,&voff);CHKERRQ(ierr);
677e630c359SToby Isaac                       } else {
678e630c359SToby Isaac                         ierr = PetscSectionGetOffset(section,closure[v],&voff);CHKERRQ(ierr);
679e630c359SToby Isaac                       }
680e630c359SToby Isaac                       xpoint = &x[voff];
681e630c359SToby Isaac                       for (j = 0; j < fbs; j++) {
682955d60d1SToby Isaac                         y[cnt + off++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j]));
683e630c359SToby Isaac                       }
684e630c359SToby Isaac                       for (; j < 3; j++) y[cnt + off++] = 0.;
685e630c359SToby Isaac                     }
686e630c359SToby Isaac                   }
687e630c359SToby Isaac                   cnt += off;
688e630c359SToby Isaac                   ierr = DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
689e630c359SToby Isaac                 }
690e630c359SToby Isaac               }
691*2c71b3e2SJacob Faibussowitsch               PetscCheckFalse(cnt != piece.nvertices*3,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
6926030a318SStefano Zampini               ierr = TransferWrite(comm,viewer,fp,r,0,y,buffer,piece.nvertices*3,MPIU_VTUREAL,tag);CHKERRQ(ierr);
693955d60d1SToby Isaac             }
694e630c359SToby Isaac           } else {
695e630c359SToby Isaac             for (i=0; i<fbs; i++) {
696955d60d1SToby Isaac               PetscInt cnt, l;
697955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
698e630c359SToby Isaac                 if (!localized) {
699e630c359SToby Isaac                   for (v=vStart,cnt=0; v<vEnd; v++) {
700e630c359SToby Isaac                     PetscInt          off;
701e630c359SToby Isaac                     const PetscScalar *xpoint;
702e630c359SToby Isaac 
703e630c359SToby Isaac                     if (nfields) {
704e630c359SToby Isaac                       ierr = PetscSectionGetFieldOffset(section,v,field,&off);CHKERRQ(ierr);
705e630c359SToby Isaac                     } else {
706e630c359SToby Isaac                       ierr = PetscSectionGetOffset(section,v,&off);CHKERRQ(ierr);
707e630c359SToby Isaac                     }
708e630c359SToby Isaac                     xpoint   = &x[off];
709955d60d1SToby Isaac                     y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i]));
710552f7358SJed Brown                   }
7112e529ec8SStefano Zampini                 } else {
7122e529ec8SStefano Zampini                   for (c=cStart,cnt=0; c<cEnd; c++) {
7132e529ec8SStefano Zampini                     PetscInt *closure = NULL;
7142e529ec8SStefano Zampini                     PetscInt  closureSize, off;
7152e529ec8SStefano Zampini 
716e630c359SToby Isaac                     ierr = DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
7172e529ec8SStefano Zampini                     for (v = 0, off = 0; v < closureSize*2; v += 2) {
7182e529ec8SStefano Zampini                       if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
719e630c359SToby Isaac                         PetscInt    voff;
720e630c359SToby Isaac                         const PetscScalar *xpoint;
7212e529ec8SStefano Zampini 
722e630c359SToby Isaac                         if (nfields) {
723e630c359SToby Isaac                           ierr = PetscSectionGetFieldOffset(section,closure[v],field,&voff);CHKERRQ(ierr);
724e630c359SToby Isaac                         } else {
725e630c359SToby Isaac                           ierr = PetscSectionGetOffset(section,closure[v],&voff);CHKERRQ(ierr);
726e630c359SToby Isaac                         }
727e630c359SToby Isaac                         xpoint         = &x[voff];
728955d60d1SToby Isaac                         y[cnt + off++] = (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i]));
7292e529ec8SStefano Zampini                       }
7302e529ec8SStefano Zampini                     }
7319bda96f6SStefano Zampini                     cnt += off;
732e630c359SToby Isaac                     ierr = DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
7332e529ec8SStefano Zampini                   }
7342e529ec8SStefano Zampini                 }
735*2c71b3e2SJacob Faibussowitsch                 PetscCheckFalse(cnt != piece.nvertices,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
7366030a318SStefano Zampini                 ierr = TransferWrite(comm,viewer,fp,r,0,y,buffer,piece.nvertices,MPIU_VTUREAL,tag);CHKERRQ(ierr);
737955d60d1SToby Isaac               }
738552f7358SJed Brown             }
739e630c359SToby Isaac           }
740e630c359SToby Isaac         }
741552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
742552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
743552f7358SJed Brown       }
744dd400576SPatrick Sanan     } else if (rank == 0) {
745955d60d1SToby Isaac       PetscInt l;
746955d60d1SToby Isaac 
7476030a318SStefano Zampini       ierr = TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,MPIU_VTUREAL,tag);CHKERRQ(ierr); /* positions */
7486030a318SStefano Zampini       ierr = TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].nconn,MPI_INT,tag);CHKERRQ(ierr); /* connectivity */
7496030a318SStefano Zampini       ierr = TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag);CHKERRQ(ierr); /* offsets */
7506030a318SStefano Zampini       ierr = TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_CHAR,tag);CHKERRQ(ierr); /* types */
7516030a318SStefano Zampini       ierr = TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag);CHKERRQ(ierr); /* owner rank (cells) */
752552f7358SJed Brown       /* all cell data */
753552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
754e630c359SToby Isaac         Vec          X = (Vec)link->vec;
7556030a318SStefano Zampini         PetscInt     bs = 1, nfields, field;
756e630c359SToby Isaac         DM           dmX;
757e630c359SToby Isaac         PetscSection section = NULL;
758e630c359SToby Isaac 
759552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
760e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
761e630c359SToby Isaac         if (!dmX) dmX = dm;
762e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
763e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
7646030a318SStefano Zampini         if (cEnd > cStart) { ierr = PetscSectionGetDof(section,cStart,&bs);CHKERRQ(ierr); }
765e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
766e630c359SToby Isaac         field = 0;
767e630c359SToby Isaac         if (link->field >= 0) {
768e630c359SToby Isaac           field = link->field;
769e630c359SToby Isaac           nfields = field + 1;
770e630c359SToby Isaac         }
771e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
772e630c359SToby Isaac           PetscInt     fbs,j;
773e630c359SToby Isaac           PetscFV      fv = NULL;
774e630c359SToby Isaac           PetscObject  f;
775e630c359SToby Isaac           PetscClassId fClass;
776e630c359SToby Isaac           PetscBool    vector;
7776030a318SStefano Zampini           if (nfields && cEnd > cStart) {        /* We have user-defined fields/components */
778e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,cStart,field,&fbs);CHKERRQ(ierr);
779e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
780e630c359SToby Isaac           ierr = DMGetField(dmX,field,NULL,&f);CHKERRQ(ierr);
781e630c359SToby Isaac           ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr);
782e630c359SToby Isaac           if (fClass == PETSCFV_CLASSID) {
783e630c359SToby Isaac             fv = (PetscFV) f;
784e630c359SToby Isaac           }
785e630c359SToby Isaac           vector = PETSC_FALSE;
786e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
787e630c359SToby Isaac             vector = PETSC_TRUE;
788e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
789e630c359SToby Isaac               const char *compName = NULL;
790e630c359SToby Isaac               if (fv) {
791e630c359SToby Isaac                 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
792e630c359SToby Isaac                 if (compName) break;
793e630c359SToby Isaac               }
794e630c359SToby Isaac             }
795e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
796e630c359SToby Isaac           }
797e630c359SToby Isaac           if (vector) {
798955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
7996030a318SStefano Zampini               ierr = TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].ncells*3,MPIU_VTUREAL,tag);CHKERRQ(ierr);
800955d60d1SToby Isaac             }
801e630c359SToby Isaac           } else {
802e630c359SToby Isaac             for (i=0; i<fbs; i++) {
803955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
8046030a318SStefano Zampini                 ierr = TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPIU_VTUREAL,tag);CHKERRQ(ierr);
805955d60d1SToby Isaac               }
806552f7358SJed Brown             }
807552f7358SJed Brown           }
808e630c359SToby Isaac         }
809e630c359SToby Isaac       }
810552f7358SJed Brown       /* all point data */
811552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
812e630c359SToby Isaac         Vec          X = (Vec)link->vec;
813e630c359SToby Isaac         DM           dmX;
8146030a318SStefano Zampini         PetscInt     bs = 1, nfields, field;
815e630c359SToby Isaac         PetscSection section = NULL;
816e630c359SToby Isaac 
817552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
818e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
819e630c359SToby Isaac         if (!dmX) dmX = dm;
820e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
821e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
8226030a318SStefano Zampini         if (vEnd > vStart) { ierr = PetscSectionGetDof(section,vStart,&bs);CHKERRQ(ierr); }
823e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
824e630c359SToby Isaac         field = 0;
825e630c359SToby Isaac         if (link->field >= 0) {
826e630c359SToby Isaac           field = link->field;
827e630c359SToby Isaac           nfields = field + 1;
828e630c359SToby Isaac         }
829e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
830e630c359SToby Isaac           PetscInt   fbs;
8316030a318SStefano Zampini           if (nfields && vEnd > vStart) {        /* We have user-defined fields/components */
832e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,vStart,field,&fbs);CHKERRQ(ierr);
833e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
834e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
835955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
8366030a318SStefano Zampini               ierr = TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,MPIU_VTUREAL,tag);CHKERRQ(ierr);
837955d60d1SToby Isaac             }
838e630c359SToby Isaac           } else {
839e630c359SToby Isaac             for (i=0; i<fbs; i++) {
840955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
8416030a318SStefano Zampini                 ierr = TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices,MPIU_VTUREAL,tag);CHKERRQ(ierr);
842955d60d1SToby Isaac               }
843552f7358SJed Brown             }
844552f7358SJed Brown           }
845552f7358SJed Brown         }
846552f7358SJed Brown       }
847e630c359SToby Isaac     }
848e630c359SToby Isaac   }
849552f7358SJed Brown   ierr = PetscFree(gpiece);CHKERRQ(ierr);
850552f7358SJed Brown   ierr = PetscFree(buffer);CHKERRQ(ierr);
851552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"\n  </AppendedData>\n");CHKERRQ(ierr);
852552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
8535e32de72SMatthew G. Knepley   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
854552f7358SJed Brown   PetscFunctionReturn(0);
855552f7358SJed Brown }
856