xref: /petsc/src/dm/impls/plex/plexvtu.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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 
30552f7358SJed Brown   PetscFunctionBegin;
319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
32552f7358SJed Brown   if (rank == srank && rank != root) {
339566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send((void*)send,count,mpidatatype,root,tag,comm));
34552f7358SJed Brown   } else if (rank == root) {
35552f7358SJed Brown     const void *buffer;
36552f7358SJed Brown     if (root == srank) {        /* self */
37552f7358SJed Brown       buffer = send;
38552f7358SJed Brown     } else {
39552f7358SJed Brown       MPI_Status  status;
40552f7358SJed Brown       PetscMPIInt nrecv;
419566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(recv,count,mpidatatype,srank,tag,comm,&status));
429566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Get_count(&status,mpidatatype,&nrecv));
43*08401ef6SPierre Jolivet       PetscCheck(count == nrecv,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
44552f7358SJed Brown       buffer = recv;
45552f7358SJed Brown     }
469566063dSJacob Faibussowitsch     PetscCall(PetscViewerVTKFWrite(viewer,fp,buffer,count,mpidatatype));
47552f7358SJed Brown   }
48552f7358SJed Brown   PetscFunctionReturn(0);
49552f7358SJed Brown }
50552f7358SJed Brown 
512e529ec8SStefano Zampini static PetscErrorCode DMPlexGetVTKConnectivity(DM dm, PetscBool localized, PieceInfo *piece,PetscVTKInt **oconn,PetscVTKInt **ooffsets,PetscVTKType **otypes)
52552f7358SJed Brown {
532e529ec8SStefano Zampini   PetscSection   coordSection;
54552f7358SJed Brown   PetscVTKInt    *conn,*offsets;
55552f7358SJed Brown   PetscVTKType   *types;
562f029394SStefano Zampini   PetscInt       dim,vStart,vEnd,cStart,cEnd,pStart,pEnd,cellHeight,numLabelCells,hasLabel,c,v,countcell,countconn;
57552f7358SJed Brown 
58552f7358SJed Brown   PetscFunctionBegin;
599566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(piece->nconn,&conn,piece->ncells,&offsets,piece->ncells,&types));
609566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dm, &coordSection));
619566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
629566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm,&pStart,&pEnd));
639566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
649566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
659566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
669566063dSJacob Faibussowitsch   PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells));
67552f7358SJed Brown   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
68552f7358SJed Brown 
69552f7358SJed Brown   countcell = 0;
70552f7358SJed Brown   countconn = 0;
71552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
7219003fb5SStefano Zampini     PetscInt nverts,dof = 0,celltype,startoffset,nC=0;
73552f7358SJed Brown 
74552f7358SJed Brown     if (hasLabel) {
75552f7358SJed Brown       PetscInt value;
76552f7358SJed Brown 
779566063dSJacob Faibussowitsch       PetscCall(DMGetLabelValue(dm, "vtk", c, &value));
78552f7358SJed Brown       if (value != 1) continue;
79552f7358SJed Brown     }
80552f7358SJed Brown     startoffset = countconn;
8119003fb5SStefano Zampini     if (localized) {
829566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(coordSection, c, &dof));
8319003fb5SStefano Zampini     }
8419003fb5SStefano Zampini     if (!dof) {
852e529ec8SStefano Zampini       PetscInt *closure = NULL;
862e529ec8SStefano Zampini       PetscInt  closureSize;
872e529ec8SStefano Zampini 
889566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
89552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
90552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
9119003fb5SStefano Zampini           if (!localized) conn[countconn++] = closure[v] - vStart;
9219003fb5SStefano Zampini           else conn[countconn++] = startoffset + nC;
93724b5320SMatthew G. Knepley           ++nC;
94552f7358SJed Brown         }
95552f7358SJed Brown       }
969566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
972e529ec8SStefano Zampini     } else {
982e529ec8SStefano Zampini       for (nC = 0; nC < dof/dim; nC++) conn[countconn++] = startoffset + nC;
992e529ec8SStefano Zampini     }
10096ca5757SLisandro Dalcin 
10196ca5757SLisandro Dalcin     {
10296ca5757SLisandro Dalcin       PetscInt n = PetscMin(nC, 8), s = countconn - nC, i, cone[8];
10396ca5757SLisandro Dalcin       for (i = 0; i < n; ++i) cone[i] = conn[s+i];
1049566063dSJacob Faibussowitsch       PetscCall(DMPlexReorderCell(dm, c, cone));
10596ca5757SLisandro Dalcin       for (i = 0; i < n; ++i) conn[s+i] = (int)cone[i];
10696ca5757SLisandro Dalcin     }
1078865f1eaSKarl Rupp 
108552f7358SJed Brown     offsets[countcell] = countconn;
1098865f1eaSKarl Rupp 
110552f7358SJed Brown     nverts = countconn - startoffset;
1119566063dSJacob Faibussowitsch     PetscCall(DMPlexVTKGetCellType_Internal(dm,dim,nverts,&celltype));
1128865f1eaSKarl Rupp 
113552f7358SJed Brown     types[countcell] = celltype;
114552f7358SJed Brown     countcell++;
115552f7358SJed Brown   }
116*08401ef6SPierre Jolivet   PetscCheck(countcell == piece->ncells,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count");
117*08401ef6SPierre Jolivet   PetscCheck(countconn == piece->nconn,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count");
118552f7358SJed Brown   *oconn    = conn;
119552f7358SJed Brown   *ooffsets = offsets;
120552f7358SJed Brown   *otypes   = types;
121552f7358SJed Brown   PetscFunctionReturn(0);
122552f7358SJed Brown }
123552f7358SJed Brown 
124552f7358SJed Brown /*
125552f7358SJed Brown   Write all fields that have been provided to the viewer
126552f7358SJed Brown   Multi-block XML format with binary appended data.
127552f7358SJed Brown */
128552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer)
129552f7358SJed Brown {
130ce94432eSBarry Smith   MPI_Comm                 comm;
1312e529ec8SStefano Zampini   PetscSection             coordSection;
132552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
133552f7358SJed Brown   PetscViewerVTKObjectLink link;
134552f7358SJed Brown   FILE                     *fp;
135552f7358SJed Brown   PetscMPIInt              rank,size,tag;
1362f029394SStefano Zampini   PetscInt                 dimEmbed,cellHeight,cStart,cEnd,vStart,vEnd,numLabelCells,hasLabel,c,v,r,i;
1372e529ec8SStefano Zampini   PetscBool                localized;
1380298fd71SBarry Smith   PieceInfo                piece,*gpiece = NULL;
1390298fd71SBarry Smith   void                     *buffer = NULL;
14030815ce0SLisandro Dalcin   const char               *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian";
141955d60d1SToby Isaac   PetscInt                 loops_per_scalar;
142552f7358SJed Brown 
143552f7358SJed Brown   PetscFunctionBegin;
1449566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
145552f7358SJed Brown #if defined(PETSC_USE_COMPLEX)
146955d60d1SToby Isaac   loops_per_scalar = 2;
147955d60d1SToby Isaac #else
148955d60d1SToby Isaac   loops_per_scalar = 1;
149552f7358SJed Brown #endif
1509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
1519566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
1529566063dSJacob Faibussowitsch   PetscCall(PetscCommGetNewTag(comm,&tag));
153552f7358SJed Brown 
1549566063dSJacob Faibussowitsch   PetscCall(PetscFOpen(comm,vtk->filename,"wb",&fp));
1559566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n"));
1569566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"%s\">\n", byte_order));
1579566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm,fp,"  <UnstructuredGrid>\n"));
158552f7358SJed Brown 
1599566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
1609566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
1619566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
1629566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1639566063dSJacob Faibussowitsch   PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells));
1649566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalized(dm, &localized));
1659566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dm, &coordSection));
1668865f1eaSKarl Rupp 
167552f7358SJed Brown   hasLabel        = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
1682e529ec8SStefano Zampini   piece.nvertices = 0;
169552f7358SJed Brown   piece.ncells    = 0;
170552f7358SJed Brown   piece.nconn     = 0;
1712e529ec8SStefano Zampini   if (!localized) piece.nvertices = vEnd - vStart;
172552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
17319003fb5SStefano Zampini     PetscInt dof = 0;
174552f7358SJed Brown     if (hasLabel) {
175552f7358SJed Brown       PetscInt value;
176552f7358SJed Brown 
1779566063dSJacob Faibussowitsch       PetscCall(DMGetLabelValue(dm, "vtk", c, &value));
178552f7358SJed Brown       if (value != 1) continue;
179552f7358SJed Brown     }
18019003fb5SStefano Zampini     if (localized) {
1819566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(coordSection, c, &dof));
18219003fb5SStefano Zampini     }
18319003fb5SStefano Zampini     if (!dof) {
1842e529ec8SStefano Zampini       PetscInt *closure = NULL;
1852e529ec8SStefano Zampini       PetscInt closureSize;
1862e529ec8SStefano Zampini 
1879566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
188552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
1892e529ec8SStefano Zampini         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
1902e529ec8SStefano Zampini           piece.nconn++;
19119003fb5SStefano Zampini           if (localized) piece.nvertices++;
1922e529ec8SStefano Zampini         }
193552f7358SJed Brown       }
1949566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
1952e529ec8SStefano Zampini     } else {
1962e529ec8SStefano Zampini       piece.nvertices += dof/dimEmbed;
1972e529ec8SStefano Zampini       piece.nconn     += dof/dimEmbed;
1982e529ec8SStefano Zampini     }
199552f7358SJed Brown     piece.ncells++;
200552f7358SJed Brown   }
2019566063dSJacob Faibussowitsch   if (rank == 0) PetscCall(PetscMalloc1(size,&gpiece));
2029566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm));
203552f7358SJed Brown 
204552f7358SJed Brown   /*
205552f7358SJed Brown    * Write file header
206552f7358SJed Brown    */
207dd400576SPatrick Sanan   if (rank == 0) {
208552f7358SJed Brown     PetscInt boffset = 0;
209552f7358SJed Brown 
210552f7358SJed Brown     for (r=0; r<size; r++) {
2119566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp,"    <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells));
212552f7358SJed Brown       /* Coordinate positions */
2139566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp,"      <Points>\n"));
2149566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset));
215955d60d1SToby Isaac       boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int);
2169566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp,"      </Points>\n"));
217552f7358SJed Brown       /* Cell connectivity */
2189566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp,"      <Cells>\n"));
2199566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset));
220552f7358SJed Brown       boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int);
2219566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"offsets\"      NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset));
222552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int);
2239566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"UInt8\" Name=\"types\"        NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset));
224552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int);
2259566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp,"      </Cells>\n"));
226552f7358SJed Brown 
227552f7358SJed Brown       /*
228552f7358SJed Brown        * Cell Data headers
229552f7358SJed Brown        */
2309566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp,"      <CellData>\n"));
2319566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset));
232552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(int) + sizeof(int);
233552f7358SJed Brown       /* all the vectors */
234552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
235552f7358SJed Brown         Vec          X = (Vec)link->vec;
236e630c359SToby Isaac         DM           dmX = NULL;
2376030a318SStefano Zampini         PetscInt     bs = 1,nfields,field;
238552f7358SJed Brown         const char   *vecname = "";
239e630c359SToby Isaac         PetscSection section;
240552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
241552f7358SJed Brown         if (((PetscObject)X)->name || link != vtk->link) { /* If the object is already named, use it. If it is past the first link, name it to disambiguate. */
2429566063dSJacob Faibussowitsch           PetscCall(PetscObjectGetName((PetscObject)X,&vecname));
243552f7358SJed Brown         }
2449566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dmX));
245e630c359SToby Isaac         if (!dmX) dmX = dm;
2469566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject*) &section));
2479566063dSJacob Faibussowitsch         if (!section) PetscCall(DMGetLocalSection(dmX, &section));
2489566063dSJacob Faibussowitsch         if (cEnd > cStart) PetscCall(PetscSectionGetDof(section,cStart,&bs));
2499566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetNumFields(section,&nfields));
250e630c359SToby Isaac         field = 0;
251e630c359SToby Isaac         if (link->field >= 0) {
252e630c359SToby Isaac           field = link->field;
253e630c359SToby Isaac           nfields = field + 1;
254e630c359SToby Isaac         }
255e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
2561cfafdd3SJed Brown           PetscInt     fbs,j;
257a00cdb45SToby Isaac           PetscFV      fv = NULL;
258a00cdb45SToby Isaac           PetscObject  f;
259a00cdb45SToby Isaac           PetscClassId fClass;
2600298fd71SBarry Smith           const char *fieldname = NULL;
2611cfafdd3SJed Brown           char       buf[256];
262e630c359SToby Isaac           PetscBool    vector;
2637ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
2649566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldDof(section,cStart,field,&fbs));
2659566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldName(section,field,&fieldname));
2667ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
2679566063dSJacob Faibussowitsch           PetscCall(DMGetField(dmX,field,NULL,&f));
2689566063dSJacob Faibussowitsch           PetscCall(PetscObjectGetClassId(f,&fClass));
269a00cdb45SToby Isaac           if (fClass == PETSCFV_CLASSID) {
270a00cdb45SToby Isaac             fv = (PetscFV) f;
271a00cdb45SToby Isaac           }
272e630c359SToby Isaac           if (nfields && !fieldname) {
2739566063dSJacob Faibussowitsch             PetscCall(PetscSNPrintf(buf,sizeof(buf),"CellField%D",field));
274552f7358SJed Brown             fieldname = buf;
275552f7358SJed Brown           }
276e630c359SToby Isaac           vector = PETSC_FALSE;
277e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
278e630c359SToby Isaac             vector = PETSC_TRUE;
279*08401ef6SPierre Jolivet             PetscCheck(fbs <= 3,PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_SIZ,"Cell vector fields can have at most 3 components, %D given", fbs);
280e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
281e630c359SToby Isaac               const char *compName = NULL;
282e630c359SToby Isaac               if (fv) {
2839566063dSJacob Faibussowitsch                 PetscCall(PetscFVGetComponentName(fv,j,&compName));
284e630c359SToby Isaac                 if (compName) break;
285e630c359SToby Isaac               }
286e630c359SToby Isaac             }
287e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
288e630c359SToby Isaac           }
289e630c359SToby Isaac           if (vector) {
290955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
2919566063dSJacob Faibussowitsch             PetscCall(PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s.Re\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset));
292955d60d1SToby Isaac             boffset += gpiece[r].ncells*3*sizeof(PetscVTUReal) + sizeof(int);
2939566063dSJacob Faibussowitsch             PetscCall(PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s.Im\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset));
294955d60d1SToby Isaac             boffset += gpiece[r].ncells*3*sizeof(PetscVTUReal) + sizeof(int);
295955d60d1SToby Isaac #else
2969566063dSJacob Faibussowitsch             PetscCall(PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset));
297955d60d1SToby Isaac             boffset += gpiece[r].ncells*3*sizeof(PetscVTUReal) + sizeof(int);
298955d60d1SToby Isaac #endif
299e630c359SToby Isaac             i+=fbs;
300e630c359SToby Isaac           } else {
3011cfafdd3SJed Brown             for (j=0; j<fbs; j++) {
302a00cdb45SToby Isaac               const char *compName = NULL;
303955d60d1SToby Isaac               char finalname[256];
304a00cdb45SToby Isaac               if (fv) {
3059566063dSJacob Faibussowitsch                 PetscCall(PetscFVGetComponentName(fv,j,&compName));
306a00cdb45SToby Isaac               }
307a00cdb45SToby Isaac               if (compName) {
3089566063dSJacob Faibussowitsch                 PetscCall(PetscSNPrintf(finalname,255,"%s%s.%s",vecname,fieldname,compName));
309a00cdb45SToby Isaac               }
310e630c359SToby Isaac               else if (fbs > 1) {
3119566063dSJacob Faibussowitsch                 PetscCall(PetscSNPrintf(finalname,255,"%s%s.%D",vecname,fieldname,j));
312e630c359SToby Isaac               } else {
3139566063dSJacob Faibussowitsch                 PetscCall(PetscSNPrintf(finalname,255,"%s%s",vecname,fieldname));
314a00cdb45SToby Isaac               }
315955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
3169566063dSJacob Faibussowitsch               PetscCall(PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset));
317955d60d1SToby Isaac               boffset += gpiece[r].ncells*sizeof(PetscVTUReal) + sizeof(int);
3189566063dSJacob Faibussowitsch               PetscCall(PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset));
319955d60d1SToby Isaac               boffset += gpiece[r].ncells*sizeof(PetscVTUReal) + sizeof(int);
320955d60d1SToby Isaac #else
3219566063dSJacob Faibussowitsch               PetscCall(PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset));
322955d60d1SToby Isaac               boffset += gpiece[r].ncells*sizeof(PetscVTUReal) + sizeof(int);
323955d60d1SToby Isaac #endif
3241cfafdd3SJed Brown               i++;
325552f7358SJed Brown             }
326552f7358SJed Brown           }
327e630c359SToby Isaac         }
328*08401ef6SPierre Jolivet         //PetscCheck(i == bs,comm,PETSC_ERR_PLIB,"Total number of field components %D != block size %D",i,bs);
3291cfafdd3SJed Brown       }
3309566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp,"      </CellData>\n"));
331552f7358SJed Brown 
332552f7358SJed Brown       /*
333552f7358SJed Brown        * Point Data headers
334552f7358SJed Brown        */
3359566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp,"      <PointData>\n"));
336552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
337552f7358SJed Brown         Vec          X = (Vec)link->vec;
338e630c359SToby Isaac         DM           dmX;
3396030a318SStefano Zampini         PetscInt     bs = 1,nfields,field;
340552f7358SJed Brown         const char   *vecname = "";
341e630c359SToby Isaac         PetscSection section;
342552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
343552f7358SJed Brown         if (((PetscObject)X)->name || link != vtk->link) { /* If the object is already named, use it. If it is past the first link, name it to disambiguate. */
3449566063dSJacob Faibussowitsch           PetscCall(PetscObjectGetName((PetscObject)X,&vecname));
345552f7358SJed Brown         }
3469566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dmX));
347e630c359SToby Isaac         if (!dmX) dmX = dm;
3489566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject*) &section));
3499566063dSJacob Faibussowitsch         if (!section) PetscCall(DMGetLocalSection(dmX, &section));
3509566063dSJacob Faibussowitsch         if (vEnd > vStart) PetscCall(PetscSectionGetDof(section,vStart,&bs));
3519566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetNumFields(section,&nfields));
352e630c359SToby Isaac         field = 0;
353e630c359SToby Isaac         if (link->field >= 0) {
354e630c359SToby Isaac           field = link->field;
355e630c359SToby Isaac           nfields = field + 1;
356e630c359SToby Isaac         }
357e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
3581cfafdd3SJed Brown           PetscInt   fbs,j;
3591cfafdd3SJed Brown           const char *fieldname = NULL;
3601cfafdd3SJed Brown           char       buf[256];
3617ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
3629566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldDof(section,vStart,field,&fbs));
3639566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldName(section,field,&fieldname));
3647ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
365e630c359SToby Isaac           if (nfields && !fieldname) {
3669566063dSJacob Faibussowitsch             PetscCall(PetscSNPrintf(buf,sizeof(buf),"PointField%D",field));
3671cfafdd3SJed Brown             fieldname = buf;
3681cfafdd3SJed Brown           }
369e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
370*08401ef6SPierre Jolivet             PetscCheck(fbs <= 3,PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_SIZ,"Point vector fields can have at most 3 components, %D given", fbs);
371955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
3729566063dSJacob Faibussowitsch             PetscCall(PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s.Re\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset));
373955d60d1SToby Isaac             boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int);
3749566063dSJacob Faibussowitsch             PetscCall(PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s.Im\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset));
375955d60d1SToby Isaac             boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int);
376955d60d1SToby Isaac #else
3779566063dSJacob Faibussowitsch             PetscCall(PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset));
378955d60d1SToby Isaac             boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int);
379955d60d1SToby Isaac #endif
380e630c359SToby Isaac           } else {
3811cfafdd3SJed Brown             for (j=0; j<fbs; j++) {
382b778fa18SValeria Barra               const char *compName = NULL;
383955d60d1SToby Isaac               char finalname[256];
3849566063dSJacob Faibussowitsch               PetscCall(PetscSectionGetComponentName(section,field,j,&compName));
3859566063dSJacob Faibussowitsch               PetscCall(PetscSNPrintf(finalname,255,"%s%s.%s",vecname,fieldname,compName));
386955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
3879566063dSJacob Faibussowitsch               PetscCall(PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset));
388955d60d1SToby Isaac               boffset += gpiece[r].nvertices*sizeof(PetscVTUReal) + sizeof(int);
3899566063dSJacob Faibussowitsch               PetscCall(PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset));
390955d60d1SToby Isaac               boffset += gpiece[r].nvertices*sizeof(PetscVTUReal) + sizeof(int);
391955d60d1SToby Isaac #else
3929566063dSJacob Faibussowitsch               PetscCall(PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,finalname,boffset));
393955d60d1SToby Isaac               boffset += gpiece[r].nvertices*sizeof(PetscVTUReal) + sizeof(int);
394955d60d1SToby Isaac #endif
395552f7358SJed Brown             }
396552f7358SJed Brown           }
3971cfafdd3SJed Brown         }
398e630c359SToby Isaac       }
3999566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp,"      </PointData>\n"));
4009566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp,"    </Piece>\n"));
401552f7358SJed Brown     }
402552f7358SJed Brown   }
403552f7358SJed Brown 
4049566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm,fp,"  </UnstructuredGrid>\n"));
4059566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n"));
4069566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm,fp,"_"));
407552f7358SJed Brown 
408dd400576SPatrick Sanan   if (rank == 0) {
409552f7358SJed Brown     PetscInt maxsize = 0;
410552f7358SJed Brown     for (r=0; r<size; r++) {
411955d60d1SToby Isaac       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscVTUReal)));
412955d60d1SToby Isaac       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*3*sizeof(PetscVTUReal)));
413552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt)));
414552f7358SJed Brown     }
4159566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(maxsize,&buffer));
416552f7358SJed Brown   }
417552f7358SJed Brown   for (r=0; r<size; r++) {
418552f7358SJed Brown     if (r == rank) {
419552f7358SJed Brown       PetscInt nsend;
420552f7358SJed Brown       {                         /* Position */
421552f7358SJed Brown         const PetscScalar *x;
422955d60d1SToby Isaac         PetscVTUReal      *y = NULL;
423552f7358SJed Brown         Vec               coords;
424955d60d1SToby Isaac         PetscBool         copy;
4252e529ec8SStefano Zampini 
4269566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinatesLocal(dm,&coords));
4279566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(coords,&x));
428955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
429955d60d1SToby Isaac         copy = PETSC_TRUE;
430955d60d1SToby Isaac #else
431955d60d1SToby Isaac         copy = (PetscBool) (dimEmbed != 3 || localized || (sizeof(PetscReal) != sizeof(PetscVTUReal)));
432955d60d1SToby Isaac #endif
433955d60d1SToby Isaac         if (copy) {
4349566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(piece.nvertices*3,&y));
43519003fb5SStefano Zampini           if (localized) {
43619003fb5SStefano Zampini             PetscInt cnt;
43719003fb5SStefano Zampini             for (c=cStart,cnt=0; c<cEnd; c++) {
43819003fb5SStefano Zampini               PetscInt off, dof;
43919003fb5SStefano Zampini 
4409566063dSJacob Faibussowitsch               PetscCall(PetscSectionGetDof(coordSection, c, &dof));
44119003fb5SStefano Zampini               if (!dof) {
44219003fb5SStefano Zampini                 PetscInt *closure = NULL;
44319003fb5SStefano Zampini                 PetscInt closureSize;
44419003fb5SStefano Zampini 
4459566063dSJacob Faibussowitsch                 PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
44619003fb5SStefano Zampini                 for (v = 0; v < closureSize*2; v += 2) {
44719003fb5SStefano Zampini                   if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
4489566063dSJacob Faibussowitsch                     PetscCall(PetscSectionGetOffset(coordSection, closure[v], &off));
44919003fb5SStefano Zampini                     if (dimEmbed != 3) {
450955d60d1SToby Isaac                       y[cnt*3+0] = (PetscVTUReal) PetscRealPart(x[off+0]);
451955d60d1SToby Isaac                       y[cnt*3+1] = (PetscVTUReal) ((dimEmbed > 1) ? PetscRealPart(x[off+1]) : 0.0);
452955d60d1SToby Isaac                       y[cnt*3+2] = (PetscVTUReal) 0.0;
45319003fb5SStefano Zampini                     } else {
454955d60d1SToby Isaac                       y[cnt*3+0] = (PetscVTUReal) PetscRealPart(x[off+0]);
455955d60d1SToby Isaac                       y[cnt*3+1] = (PetscVTUReal) PetscRealPart(x[off+1]);
456955d60d1SToby Isaac                       y[cnt*3+2] = (PetscVTUReal) PetscRealPart(x[off+2]);
45719003fb5SStefano Zampini                     }
45819003fb5SStefano Zampini                     cnt++;
45919003fb5SStefano Zampini                   }
46019003fb5SStefano Zampini                 }
4619566063dSJacob Faibussowitsch                 PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
46219003fb5SStefano Zampini               } else {
4639566063dSJacob Faibussowitsch                 PetscCall(PetscSectionGetOffset(coordSection, c, &off));
46419003fb5SStefano Zampini                 if (dimEmbed != 3) {
46519003fb5SStefano Zampini                   for (i=0; i<dof/dimEmbed; i++) {
466955d60d1SToby Isaac                     y[cnt*3+0] = (PetscVTUReal) PetscRealPart(x[off + i*dimEmbed + 0]);
467955d60d1SToby Isaac                     y[cnt*3+1] = (PetscVTUReal) ((dimEmbed > 1) ? PetscRealPart(x[off + i*dimEmbed + 1]) : 0.0);
468955d60d1SToby Isaac                     y[cnt*3+2] = (PetscVTUReal) 0.0;
46919003fb5SStefano Zampini                     cnt++;
47019003fb5SStefano Zampini                   }
47119003fb5SStefano Zampini                 } else {
47219003fb5SStefano Zampini                   for (i=0; i<dof; i ++) {
473955d60d1SToby Isaac                     y[cnt*3+i] = (PetscVTUReal) PetscRealPart(x[off + i]);
47419003fb5SStefano Zampini                   }
47519003fb5SStefano Zampini                   cnt += dof/dimEmbed;
47619003fb5SStefano Zampini                 }
47719003fb5SStefano Zampini               }
47819003fb5SStefano Zampini             }
479*08401ef6SPierre Jolivet             PetscCheck(cnt == piece.nvertices,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
48019003fb5SStefano Zampini           } else {
481552f7358SJed Brown             for (i=0; i<piece.nvertices; i++) {
482955d60d1SToby Isaac               y[i*3+0] = (PetscVTUReal) PetscRealPart(x[i*dimEmbed+0]);
483955d60d1SToby Isaac               y[i*3+1] = (PetscVTUReal) ((dimEmbed > 1) ? PetscRealPart(x[i*dimEmbed+1]) : 0.);
484955d60d1SToby Isaac               y[i*3+2] = (PetscVTUReal) ((dimEmbed > 2) ? PetscRealPart(x[i*dimEmbed+2]) : 0.);
48519003fb5SStefano Zampini             }
486552f7358SJed Brown           }
487552f7358SJed Brown         }
4882e529ec8SStefano Zampini         nsend = piece.nvertices*3;
4899566063dSJacob Faibussowitsch         PetscCall(TransferWrite(comm,viewer,fp,r,0,copy ? (const void *) y : (const void *) x,buffer,nsend,MPIU_VTUREAL,tag));
4909566063dSJacob Faibussowitsch         PetscCall(PetscFree(y));
4919566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(coords,&x));
492552f7358SJed Brown       }
493552f7358SJed Brown       {                           /* Connectivity, offsets, types */
4943e869605SMatthew G. Knepley         PetscVTKInt  *connectivity = NULL, *offsets = NULL;
4953e869605SMatthew G. Knepley         PetscVTKType *types = NULL;
4969566063dSJacob Faibussowitsch         PetscCall(DMPlexGetVTKConnectivity(dm,localized,&piece,&connectivity,&offsets,&types));
4979566063dSJacob Faibussowitsch         PetscCall(TransferWrite(comm,viewer,fp,r,0,connectivity,buffer,piece.nconn,MPI_INT,tag));
4989566063dSJacob Faibussowitsch         PetscCall(TransferWrite(comm,viewer,fp,r,0,offsets,buffer,piece.ncells,MPI_INT,tag));
4999566063dSJacob Faibussowitsch         PetscCall(TransferWrite(comm,viewer,fp,r,0,types,buffer,piece.ncells,MPI_CHAR,tag));
5009566063dSJacob Faibussowitsch         PetscCall(PetscFree3(connectivity,offsets,types));
501552f7358SJed Brown       }
502552f7358SJed Brown       {                         /* Owners (cell data) */
503552f7358SJed Brown         PetscVTKInt *owners;
5049566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(piece.ncells,&owners));
505552f7358SJed Brown         for (i=0; i<piece.ncells; i++) owners[i] = rank;
5069566063dSJacob Faibussowitsch         PetscCall(TransferWrite(comm,viewer,fp,r,0,owners,buffer,piece.ncells,MPI_INT,tag));
5079566063dSJacob Faibussowitsch         PetscCall(PetscFree(owners));
508552f7358SJed Brown       }
509552f7358SJed Brown       /* Cell data */
510552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
511552f7358SJed Brown         Vec               X = (Vec)link->vec;
512e630c359SToby Isaac         DM                dmX;
513552f7358SJed Brown         const PetscScalar *x;
514955d60d1SToby Isaac         PetscVTUReal      *y;
5156030a318SStefano Zampini         PetscInt          bs = 1, nfields, field;
516e630c359SToby Isaac         PetscSection      section = NULL;
517e630c359SToby Isaac 
518552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
5199566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dmX));
520e630c359SToby Isaac         if (!dmX) dmX = dm;
5219566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject*) &section));
5229566063dSJacob Faibussowitsch         if (!section) PetscCall(DMGetLocalSection(dmX, &section));
5239566063dSJacob Faibussowitsch         if (cEnd > cStart) PetscCall(PetscSectionGetDof(section,cStart,&bs));
5249566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetNumFields(section,&nfields));
525e630c359SToby Isaac         field = 0;
526e630c359SToby Isaac         if (link->field >= 0) {
527e630c359SToby Isaac           field = link->field;
528e630c359SToby Isaac           nfields = field + 1;
529e630c359SToby Isaac         }
5309566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(X,&x));
5319566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(piece.ncells*3,&y));
532e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
533e630c359SToby Isaac           PetscInt     fbs,j;
534e630c359SToby Isaac           PetscFV      fv = NULL;
535e630c359SToby Isaac           PetscObject  f;
536e630c359SToby Isaac           PetscClassId fClass;
537e630c359SToby Isaac           PetscBool    vector;
5386030a318SStefano Zampini           if (nfields && cEnd > cStart) {        /* We have user-defined fields/components */
5399566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldDof(section,cStart,field,&fbs));
540e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
5419566063dSJacob Faibussowitsch           PetscCall(DMGetField(dmX,field,NULL,&f));
5429566063dSJacob Faibussowitsch           PetscCall(PetscObjectGetClassId(f,&fClass));
543e630c359SToby Isaac           if (fClass == PETSCFV_CLASSID) {
544e630c359SToby Isaac             fv = (PetscFV) f;
545e630c359SToby Isaac           }
546e630c359SToby Isaac           vector = PETSC_FALSE;
547e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
548e630c359SToby Isaac             vector = PETSC_TRUE;
549e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
550e630c359SToby Isaac               const char *compName = NULL;
551e630c359SToby Isaac               if (fv) {
5529566063dSJacob Faibussowitsch                 PetscCall(PetscFVGetComponentName(fv,j,&compName));
553e630c359SToby Isaac                 if (compName) break;
554e630c359SToby Isaac               }
555e630c359SToby Isaac             }
556e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
557e630c359SToby Isaac           }
558e630c359SToby Isaac           if (vector) {
559955d60d1SToby Isaac             PetscInt cnt, l;
560955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
561552f7358SJed Brown               for (c=cStart,cnt=0; c<cEnd; c++) {
562e630c359SToby Isaac                 const PetscScalar *xpoint;
563e630c359SToby Isaac                 PetscInt off, j;
564e630c359SToby Isaac 
565552f7358SJed Brown                 if (hasLabel) {     /* Ignore some cells */
566552f7358SJed Brown                   PetscInt value;
5679566063dSJacob Faibussowitsch                   PetscCall(DMGetLabelValue(dmX, "vtk", c, &value));
568552f7358SJed Brown                   if (value != 1) continue;
569552f7358SJed Brown                 }
570e630c359SToby Isaac                 if (nfields) {
5719566063dSJacob Faibussowitsch                   PetscCall(PetscSectionGetFieldOffset(section, c, field, &off));
572e630c359SToby Isaac                 } else {
5739566063dSJacob Faibussowitsch                   PetscCall(PetscSectionGetOffset(section, c, &off));
574e630c359SToby Isaac                 }
575e630c359SToby Isaac                 xpoint = &x[off];
576e630c359SToby Isaac                 for (j = 0; j < fbs; j++) {
577955d60d1SToby Isaac                   y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j]));
578e630c359SToby Isaac                 }
579e630c359SToby Isaac                 for (; j < 3; j++) y[cnt++] = 0.;
580e630c359SToby Isaac               }
581*08401ef6SPierre Jolivet               PetscCheck(cnt == piece.ncells*3,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
5829566063dSJacob Faibussowitsch               PetscCall(TransferWrite(comm,viewer,fp,r,0,y,buffer,piece.ncells*3,MPIU_VTUREAL,tag));
583955d60d1SToby Isaac             }
584e630c359SToby Isaac           } else {
585e630c359SToby Isaac             for (i=0; i<fbs; i++) {
586955d60d1SToby Isaac               PetscInt cnt, l;
587955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
588e630c359SToby Isaac                 for (c=cStart,cnt=0; c<cEnd; c++) {
589e630c359SToby Isaac                   const PetscScalar *xpoint;
590e630c359SToby Isaac                   PetscInt off;
591e630c359SToby Isaac 
592e630c359SToby Isaac                   if (hasLabel) {     /* Ignore some cells */
593e630c359SToby Isaac                     PetscInt value;
5949566063dSJacob Faibussowitsch                     PetscCall(DMGetLabelValue(dmX, "vtk", c, &value));
595e630c359SToby Isaac                     if (value != 1) continue;
596e630c359SToby Isaac                   }
597e630c359SToby Isaac                   if (nfields) {
5989566063dSJacob Faibussowitsch                     PetscCall(PetscSectionGetFieldOffset(section, c, field, &off));
599e630c359SToby Isaac                   } else {
6009566063dSJacob Faibussowitsch                     PetscCall(PetscSectionGetOffset(section, c, &off));
601e630c359SToby Isaac                   }
602e630c359SToby Isaac                   xpoint   = &x[off];
603955d60d1SToby Isaac                   y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i]));
604552f7358SJed Brown                 }
605*08401ef6SPierre Jolivet                 PetscCheck(cnt == piece.ncells,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
6069566063dSJacob Faibussowitsch                 PetscCall(TransferWrite(comm,viewer,fp,r,0,y,buffer,piece.ncells,MPIU_VTUREAL,tag));
607955d60d1SToby Isaac               }
608552f7358SJed Brown             }
609e630c359SToby Isaac           }
610e630c359SToby Isaac         }
6119566063dSJacob Faibussowitsch         PetscCall(PetscFree(y));
6129566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(X,&x));
613552f7358SJed Brown       }
614e630c359SToby Isaac       /* point data */
615552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
616552f7358SJed Brown         Vec               X = (Vec)link->vec;
617e630c359SToby Isaac         DM                dmX;
618552f7358SJed Brown         const PetscScalar *x;
619955d60d1SToby Isaac         PetscVTUReal      *y;
6206030a318SStefano Zampini         PetscInt          bs = 1, nfields, field;
621e630c359SToby Isaac         PetscSection      section = NULL;
622e630c359SToby Isaac 
623552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
6249566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dmX));
625e630c359SToby Isaac         if (!dmX) dmX = dm;
6269566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject*) &section));
6279566063dSJacob Faibussowitsch         if (!section) PetscCall(DMGetLocalSection(dmX, &section));
6289566063dSJacob Faibussowitsch         if (vEnd > vStart) PetscCall(PetscSectionGetDof(section,vStart,&bs));
6299566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetNumFields(section,&nfields));
630e630c359SToby Isaac         field = 0;
631e630c359SToby Isaac         if (link->field >= 0) {
632e630c359SToby Isaac           field = link->field;
633e630c359SToby Isaac           nfields = field + 1;
634e630c359SToby Isaac         }
6359566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(X,&x));
6369566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(piece.nvertices*3,&y));
637e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
638e630c359SToby Isaac           PetscInt fbs,j;
6396030a318SStefano Zampini           if (nfields && vEnd > vStart) {        /* We have user-defined fields/components */
6409566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldDof(section,vStart,field,&fbs));
641e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
642e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
643955d60d1SToby Isaac             PetscInt cnt, l;
644955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
6452e529ec8SStefano Zampini               if (!localized) {
646552f7358SJed Brown                 for (v=vStart,cnt=0; v<vEnd; v++) {
647e630c359SToby Isaac                   PetscInt          off;
648e630c359SToby Isaac                   const PetscScalar *xpoint;
649e630c359SToby Isaac 
650e630c359SToby Isaac                   if (nfields) {
6519566063dSJacob Faibussowitsch                     PetscCall(PetscSectionGetFieldOffset(section,v,field,&off));
652e630c359SToby Isaac                   } else {
6539566063dSJacob Faibussowitsch                     PetscCall(PetscSectionGetOffset(section,v,&off));
654e630c359SToby Isaac                   }
655e630c359SToby Isaac                   xpoint = &x[off];
656e630c359SToby Isaac                   for (j = 0; j < fbs; j++) {
657955d60d1SToby Isaac                     y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j]));
658e630c359SToby Isaac                   }
659e630c359SToby Isaac                   for (; j < 3; j++) y[cnt++] = 0.;
660e630c359SToby Isaac                 }
661e630c359SToby Isaac               } else {
662e630c359SToby Isaac                 for (c=cStart,cnt=0; c<cEnd; c++) {
663e630c359SToby Isaac                   PetscInt *closure = NULL;
664e630c359SToby Isaac                   PetscInt  closureSize, off;
665e630c359SToby Isaac 
6669566063dSJacob Faibussowitsch                   PetscCall(DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure));
667e630c359SToby Isaac                   for (v = 0, off = 0; v < closureSize*2; v += 2) {
668e630c359SToby Isaac                     if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
669e630c359SToby Isaac                       PetscInt    voff;
670e630c359SToby Isaac                       const PetscScalar *xpoint;
671e630c359SToby Isaac 
672e630c359SToby Isaac                       if (nfields) {
6739566063dSJacob Faibussowitsch                         PetscCall(PetscSectionGetFieldOffset(section,closure[v],field,&voff));
674e630c359SToby Isaac                       } else {
6759566063dSJacob Faibussowitsch                         PetscCall(PetscSectionGetOffset(section,closure[v],&voff));
676e630c359SToby Isaac                       }
677e630c359SToby Isaac                       xpoint = &x[voff];
678e630c359SToby Isaac                       for (j = 0; j < fbs; j++) {
679955d60d1SToby Isaac                         y[cnt + off++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j]));
680e630c359SToby Isaac                       }
681e630c359SToby Isaac                       for (; j < 3; j++) y[cnt + off++] = 0.;
682e630c359SToby Isaac                     }
683e630c359SToby Isaac                   }
684e630c359SToby Isaac                   cnt += off;
6859566063dSJacob Faibussowitsch                   PetscCall(DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure));
686e630c359SToby Isaac                 }
687e630c359SToby Isaac               }
688*08401ef6SPierre Jolivet               PetscCheck(cnt == piece.nvertices*3,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
6899566063dSJacob Faibussowitsch               PetscCall(TransferWrite(comm,viewer,fp,r,0,y,buffer,piece.nvertices*3,MPIU_VTUREAL,tag));
690955d60d1SToby Isaac             }
691e630c359SToby Isaac           } else {
692e630c359SToby Isaac             for (i=0; i<fbs; i++) {
693955d60d1SToby Isaac               PetscInt cnt, l;
694955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
695e630c359SToby Isaac                 if (!localized) {
696e630c359SToby Isaac                   for (v=vStart,cnt=0; v<vEnd; v++) {
697e630c359SToby Isaac                     PetscInt          off;
698e630c359SToby Isaac                     const PetscScalar *xpoint;
699e630c359SToby Isaac 
700e630c359SToby Isaac                     if (nfields) {
7019566063dSJacob Faibussowitsch                       PetscCall(PetscSectionGetFieldOffset(section,v,field,&off));
702e630c359SToby Isaac                     } else {
7039566063dSJacob Faibussowitsch                       PetscCall(PetscSectionGetOffset(section,v,&off));
704e630c359SToby Isaac                     }
705e630c359SToby Isaac                     xpoint   = &x[off];
706955d60d1SToby Isaac                     y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i]));
707552f7358SJed Brown                   }
7082e529ec8SStefano Zampini                 } else {
7092e529ec8SStefano Zampini                   for (c=cStart,cnt=0; c<cEnd; c++) {
7102e529ec8SStefano Zampini                     PetscInt *closure = NULL;
7112e529ec8SStefano Zampini                     PetscInt  closureSize, off;
7122e529ec8SStefano Zampini 
7139566063dSJacob Faibussowitsch                     PetscCall(DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure));
7142e529ec8SStefano Zampini                     for (v = 0, off = 0; v < closureSize*2; v += 2) {
7152e529ec8SStefano Zampini                       if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
716e630c359SToby Isaac                         PetscInt    voff;
717e630c359SToby Isaac                         const PetscScalar *xpoint;
7182e529ec8SStefano Zampini 
719e630c359SToby Isaac                         if (nfields) {
7209566063dSJacob Faibussowitsch                           PetscCall(PetscSectionGetFieldOffset(section,closure[v],field,&voff));
721e630c359SToby Isaac                         } else {
7229566063dSJacob Faibussowitsch                           PetscCall(PetscSectionGetOffset(section,closure[v],&voff));
723e630c359SToby Isaac                         }
724e630c359SToby Isaac                         xpoint         = &x[voff];
725955d60d1SToby Isaac                         y[cnt + off++] = (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i]));
7262e529ec8SStefano Zampini                       }
7272e529ec8SStefano Zampini                     }
7289bda96f6SStefano Zampini                     cnt += off;
7299566063dSJacob Faibussowitsch                     PetscCall(DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure));
7302e529ec8SStefano Zampini                   }
7312e529ec8SStefano Zampini                 }
732*08401ef6SPierre Jolivet                 PetscCheck(cnt == piece.nvertices,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
7339566063dSJacob Faibussowitsch                 PetscCall(TransferWrite(comm,viewer,fp,r,0,y,buffer,piece.nvertices,MPIU_VTUREAL,tag));
734955d60d1SToby Isaac               }
735552f7358SJed Brown             }
736e630c359SToby Isaac           }
737e630c359SToby Isaac         }
7389566063dSJacob Faibussowitsch         PetscCall(PetscFree(y));
7399566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(X,&x));
740552f7358SJed Brown       }
741dd400576SPatrick Sanan     } else if (rank == 0) {
742955d60d1SToby Isaac       PetscInt l;
743955d60d1SToby Isaac 
7449566063dSJacob Faibussowitsch       PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,MPIU_VTUREAL,tag)); /* positions */
7459566063dSJacob Faibussowitsch       PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].nconn,MPI_INT,tag)); /* connectivity */
7469566063dSJacob Faibussowitsch       PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag)); /* offsets */
7479566063dSJacob Faibussowitsch       PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_CHAR,tag)); /* types */
7489566063dSJacob Faibussowitsch       PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag)); /* owner rank (cells) */
749552f7358SJed Brown       /* all cell data */
750552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
751e630c359SToby Isaac         Vec          X = (Vec)link->vec;
7526030a318SStefano Zampini         PetscInt     bs = 1, nfields, field;
753e630c359SToby Isaac         DM           dmX;
754e630c359SToby Isaac         PetscSection section = NULL;
755e630c359SToby Isaac 
756552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
7579566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dmX));
758e630c359SToby Isaac         if (!dmX) dmX = dm;
7599566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject*) &section));
7609566063dSJacob Faibussowitsch         if (!section) PetscCall(DMGetLocalSection(dmX, &section));
7619566063dSJacob Faibussowitsch         if (cEnd > cStart) PetscCall(PetscSectionGetDof(section,cStart,&bs));
7629566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetNumFields(section,&nfields));
763e630c359SToby Isaac         field = 0;
764e630c359SToby Isaac         if (link->field >= 0) {
765e630c359SToby Isaac           field = link->field;
766e630c359SToby Isaac           nfields = field + 1;
767e630c359SToby Isaac         }
768e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
769e630c359SToby Isaac           PetscInt     fbs,j;
770e630c359SToby Isaac           PetscFV      fv = NULL;
771e630c359SToby Isaac           PetscObject  f;
772e630c359SToby Isaac           PetscClassId fClass;
773e630c359SToby Isaac           PetscBool    vector;
7746030a318SStefano Zampini           if (nfields && cEnd > cStart) {        /* We have user-defined fields/components */
7759566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldDof(section,cStart,field,&fbs));
776e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
7779566063dSJacob Faibussowitsch           PetscCall(DMGetField(dmX,field,NULL,&f));
7789566063dSJacob Faibussowitsch           PetscCall(PetscObjectGetClassId(f,&fClass));
779e630c359SToby Isaac           if (fClass == PETSCFV_CLASSID) {
780e630c359SToby Isaac             fv = (PetscFV) f;
781e630c359SToby Isaac           }
782e630c359SToby Isaac           vector = PETSC_FALSE;
783e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
784e630c359SToby Isaac             vector = PETSC_TRUE;
785e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
786e630c359SToby Isaac               const char *compName = NULL;
787e630c359SToby Isaac               if (fv) {
7889566063dSJacob Faibussowitsch                 PetscCall(PetscFVGetComponentName(fv,j,&compName));
789e630c359SToby Isaac                 if (compName) break;
790e630c359SToby Isaac               }
791e630c359SToby Isaac             }
792e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
793e630c359SToby Isaac           }
794e630c359SToby Isaac           if (vector) {
795955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
7969566063dSJacob Faibussowitsch               PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].ncells*3,MPIU_VTUREAL,tag));
797955d60d1SToby Isaac             }
798e630c359SToby Isaac           } else {
799e630c359SToby Isaac             for (i=0; i<fbs; i++) {
800955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
8019566063dSJacob Faibussowitsch                 PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPIU_VTUREAL,tag));
802955d60d1SToby Isaac               }
803552f7358SJed Brown             }
804552f7358SJed Brown           }
805e630c359SToby Isaac         }
806e630c359SToby Isaac       }
807552f7358SJed Brown       /* all point data */
808552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
809e630c359SToby Isaac         Vec          X = (Vec)link->vec;
810e630c359SToby Isaac         DM           dmX;
8116030a318SStefano Zampini         PetscInt     bs = 1, nfields, field;
812e630c359SToby Isaac         PetscSection section = NULL;
813e630c359SToby Isaac 
814552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
8159566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dmX));
816e630c359SToby Isaac         if (!dmX) dmX = dm;
8179566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject*) &section));
8189566063dSJacob Faibussowitsch         if (!section) PetscCall(DMGetLocalSection(dmX, &section));
8199566063dSJacob Faibussowitsch         if (vEnd > vStart) PetscCall(PetscSectionGetDof(section,vStart,&bs));
8209566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetNumFields(section,&nfields));
821e630c359SToby Isaac         field = 0;
822e630c359SToby Isaac         if (link->field >= 0) {
823e630c359SToby Isaac           field = link->field;
824e630c359SToby Isaac           nfields = field + 1;
825e630c359SToby Isaac         }
826e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
827e630c359SToby Isaac           PetscInt   fbs;
8286030a318SStefano Zampini           if (nfields && vEnd > vStart) {        /* We have user-defined fields/components */
8299566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldDof(section,vStart,field,&fbs));
830e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
831e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
832955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
8339566063dSJacob Faibussowitsch               PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,MPIU_VTUREAL,tag));
834955d60d1SToby Isaac             }
835e630c359SToby Isaac           } else {
836e630c359SToby Isaac             for (i=0; i<fbs; i++) {
837955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
8389566063dSJacob Faibussowitsch                 PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices,MPIU_VTUREAL,tag));
839955d60d1SToby Isaac               }
840552f7358SJed Brown             }
841552f7358SJed Brown           }
842552f7358SJed Brown         }
843552f7358SJed Brown       }
844e630c359SToby Isaac     }
845e630c359SToby Isaac   }
8469566063dSJacob Faibussowitsch   PetscCall(PetscFree(gpiece));
8479566063dSJacob Faibussowitsch   PetscCall(PetscFree(buffer));
8489566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm,fp,"\n  </AppendedData>\n"));
8499566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm,fp,"</VTKFile>\n"));
8509566063dSJacob Faibussowitsch   PetscCall(PetscFClose(comm,fp));
851552f7358SJed Brown   PetscFunctionReturn(0);
852552f7358SJed Brown }
853