xref: /petsc/src/dm/impls/plex/plexvtu.c (revision c29ce622fa97c279dda1474b2f4e55cdcf02f6de)
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 
26d71ae5a4SJacob Faibussowitsch 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)
27d71ae5a4SJacob Faibussowitsch {
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));
4308401ef6SPierre 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   }
483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
49552f7358SJed Brown }
50552f7358SJed Brown 
51d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexGetVTKConnectivity(DM dm, PetscBool localized, PieceInfo *piece, PetscVTKInt **oconn, PetscVTKInt **ooffsets, PetscVTKType **otypes)
52d71ae5a4SJacob Faibussowitsch {
536858538eSMatthew G. Knepley   PetscSection  coordSection, cellCoordSection;
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));
616858538eSMatthew G. Knepley   PetscCall(DMGetCellCoordinateSection(dm, &cellCoordSection));
629566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
639566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
649566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
659566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
669566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
679566063dSJacob Faibussowitsch   PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells));
68552f7358SJed Brown   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
69552f7358SJed Brown 
70552f7358SJed Brown   countcell = 0;
71552f7358SJed Brown   countconn = 0;
72552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
7319003fb5SStefano Zampini     PetscInt nverts, dof = 0, celltype, startoffset, nC = 0;
74552f7358SJed Brown 
75552f7358SJed Brown     if (hasLabel) {
76552f7358SJed Brown       PetscInt value;
77552f7358SJed Brown 
789566063dSJacob Faibussowitsch       PetscCall(DMGetLabelValue(dm, "vtk", c, &value));
79552f7358SJed Brown       if (value != 1) continue;
80552f7358SJed Brown     }
81552f7358SJed Brown     startoffset = countconn;
826858538eSMatthew G. Knepley     if (localized) PetscCall(PetscSectionGetDof(cellCoordSection, c, &dof));
8319003fb5SStefano Zampini     if (!dof) {
842e529ec8SStefano Zampini       PetscInt *closure = NULL;
852e529ec8SStefano Zampini       PetscInt  closureSize;
862e529ec8SStefano Zampini 
879566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
88552f7358SJed Brown       for (v = 0; v < closureSize * 2; v += 2) {
89552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
9019003fb5SStefano Zampini           if (!localized) conn[countconn++] = closure[v] - vStart;
9119003fb5SStefano Zampini           else conn[countconn++] = startoffset + nC;
92724b5320SMatthew G. Knepley           ++nC;
93552f7358SJed Brown         }
94552f7358SJed Brown       }
959566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
962e529ec8SStefano Zampini     } else {
972e529ec8SStefano Zampini       for (nC = 0; nC < dof / dim; nC++) conn[countconn++] = startoffset + nC;
982e529ec8SStefano Zampini     }
9996ca5757SLisandro Dalcin 
10096ca5757SLisandro Dalcin     {
10196ca5757SLisandro Dalcin       PetscInt n = PetscMin(nC, 8), s = countconn - nC, i, cone[8];
10296ca5757SLisandro Dalcin       for (i = 0; i < n; ++i) cone[i] = conn[s + i];
1039566063dSJacob Faibussowitsch       PetscCall(DMPlexReorderCell(dm, c, cone));
10496ca5757SLisandro Dalcin       for (i = 0; i < n; ++i) conn[s + i] = (int)cone[i];
10596ca5757SLisandro Dalcin     }
1068865f1eaSKarl Rupp 
107552f7358SJed Brown     offsets[countcell] = countconn;
1088865f1eaSKarl Rupp 
109552f7358SJed Brown     nverts = countconn - startoffset;
1109566063dSJacob Faibussowitsch     PetscCall(DMPlexVTKGetCellType_Internal(dm, dim, nverts, &celltype));
1118865f1eaSKarl Rupp 
112552f7358SJed Brown     types[countcell] = celltype;
113552f7358SJed Brown     countcell++;
114552f7358SJed Brown   }
11508401ef6SPierre Jolivet   PetscCheck(countcell == piece->ncells, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent cell count");
11608401ef6SPierre Jolivet   PetscCheck(countconn == piece->nconn, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent connectivity count");
117552f7358SJed Brown   *oconn    = conn;
118552f7358SJed Brown   *ooffsets = offsets;
119552f7358SJed Brown   *otypes   = types;
1203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
121552f7358SJed Brown }
122552f7358SJed Brown 
123*c29ce622SStefano Zampini PETSC_INTERN PetscErrorCode DMPlexGetNonEmptyComm_Private(DM dm, MPI_Comm *comm)
124*c29ce622SStefano Zampini {
125*c29ce622SStefano Zampini   DM_Plex *mesh = (DM_Plex *)dm->data;
126*c29ce622SStefano Zampini 
127*c29ce622SStefano Zampini   PetscFunctionBegin;
128*c29ce622SStefano Zampini   if (mesh->nonempty_comm == MPI_COMM_SELF) { /* Not yet setup */
129*c29ce622SStefano Zampini     PetscInt    cStart, cEnd, cellHeight;
130*c29ce622SStefano Zampini     MPI_Comm    dmcomm = PetscObjectComm((PetscObject)dm);
131*c29ce622SStefano Zampini     PetscMPIInt color, rank;
132*c29ce622SStefano Zampini 
133*c29ce622SStefano Zampini     PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
134*c29ce622SStefano Zampini     PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
135*c29ce622SStefano Zampini     color = (cStart < cEnd) ? 0 : 1;
136*c29ce622SStefano Zampini     PetscCallMPI(MPI_Comm_rank(dmcomm, &rank));
137*c29ce622SStefano Zampini     PetscCallMPI(MPI_Comm_split(dmcomm, color, rank, &mesh->nonempty_comm));
138*c29ce622SStefano Zampini     if (color == 1) {
139*c29ce622SStefano Zampini       PetscCallMPI(MPI_Comm_free(&mesh->nonempty_comm));
140*c29ce622SStefano Zampini       mesh->nonempty_comm = MPI_COMM_NULL;
141*c29ce622SStefano Zampini     }
142*c29ce622SStefano Zampini   }
143*c29ce622SStefano Zampini   *comm = mesh->nonempty_comm;
144*c29ce622SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
145*c29ce622SStefano Zampini }
146*c29ce622SStefano Zampini 
147552f7358SJed Brown /*
148552f7358SJed Brown   Write all fields that have been provided to the viewer
149552f7358SJed Brown   Multi-block XML format with binary appended data.
150552f7358SJed Brown */
151d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm, PetscViewer viewer)
152d71ae5a4SJacob Faibussowitsch {
153ce94432eSBarry Smith   MPI_Comm                 comm;
1546858538eSMatthew G. Knepley   PetscSection             coordSection, cellCoordSection;
155552f7358SJed Brown   PetscViewer_VTK         *vtk = (PetscViewer_VTK *)viewer->data;
156552f7358SJed Brown   PetscViewerVTKObjectLink link;
157552f7358SJed Brown   FILE                    *fp;
158552f7358SJed Brown   PetscMPIInt              rank, size, tag;
1592f029394SStefano Zampini   PetscInt                 dimEmbed, cellHeight, cStart, cEnd, vStart, vEnd, numLabelCells, hasLabel, c, v, r, i;
1602e529ec8SStefano Zampini   PetscBool                localized;
1610298fd71SBarry Smith   PieceInfo                piece, *gpiece = NULL;
1620298fd71SBarry Smith   void                    *buffer     = NULL;
16330815ce0SLisandro Dalcin   const char              *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian";
164955d60d1SToby Isaac   PetscInt                 loops_per_scalar;
165552f7358SJed Brown 
166552f7358SJed Brown   PetscFunctionBegin;
1679566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
1689566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
1699566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
1709566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1719566063dSJacob Faibussowitsch   PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells));
1729566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalized(dm, &localized));
1739566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dm, &coordSection));
1746858538eSMatthew G. Knepley   PetscCall(DMGetCellCoordinateSection(dm, &cellCoordSection));
175*c29ce622SStefano Zampini   PetscCall(PetscCommGetNewTag(PetscObjectComm((PetscObject)dm), &tag));
176*c29ce622SStefano Zampini 
177*c29ce622SStefano Zampini   PetscCall(DMPlexGetNonEmptyComm_Private(dm, &comm));
178*c29ce622SStefano Zampini #if defined(PETSC_USE_COMPLEX)
179*c29ce622SStefano Zampini   loops_per_scalar = 2;
180*c29ce622SStefano Zampini #else
181*c29ce622SStefano Zampini   loops_per_scalar = 1;
182*c29ce622SStefano Zampini #endif
183*c29ce622SStefano Zampini   if (comm == MPI_COMM_NULL) goto finalize;
184*c29ce622SStefano Zampini   PetscCallMPI(MPI_Comm_size(comm, &size));
185*c29ce622SStefano Zampini   PetscCallMPI(MPI_Comm_rank(comm, &rank));
186*c29ce622SStefano Zampini 
187*c29ce622SStefano Zampini   PetscCall(PetscFOpen(comm, vtk->filename, "wb", &fp));
188*c29ce622SStefano Zampini   PetscCall(PetscFPrintf(comm, fp, "<?xml version=\"1.0\"?>\n"));
189*c29ce622SStefano Zampini   PetscCall(PetscFPrintf(comm, fp, "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"%s\" header_type=\"UInt64\">\n", byte_order));
190*c29ce622SStefano Zampini   PetscCall(PetscFPrintf(comm, fp, "  <UnstructuredGrid>\n"));
1918865f1eaSKarl Rupp 
192552f7358SJed Brown   hasLabel        = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
1932e529ec8SStefano Zampini   piece.nvertices = 0;
194552f7358SJed Brown   piece.ncells    = 0;
195552f7358SJed Brown   piece.nconn     = 0;
1962e529ec8SStefano Zampini   if (!localized) piece.nvertices = vEnd - vStart;
197552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
19819003fb5SStefano Zampini     PetscInt dof = 0;
199552f7358SJed Brown     if (hasLabel) {
200552f7358SJed Brown       PetscInt value;
201552f7358SJed Brown 
2029566063dSJacob Faibussowitsch       PetscCall(DMGetLabelValue(dm, "vtk", c, &value));
203552f7358SJed Brown       if (value != 1) continue;
204552f7358SJed Brown     }
2056858538eSMatthew G. Knepley     if (localized) PetscCall(PetscSectionGetDof(cellCoordSection, c, &dof));
20619003fb5SStefano Zampini     if (!dof) {
2072e529ec8SStefano Zampini       PetscInt *closure = NULL;
2082e529ec8SStefano Zampini       PetscInt  closureSize;
2092e529ec8SStefano Zampini 
2109566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
211552f7358SJed Brown       for (v = 0; v < closureSize * 2; v += 2) {
2122e529ec8SStefano Zampini         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
2132e529ec8SStefano Zampini           piece.nconn++;
21419003fb5SStefano Zampini           if (localized) piece.nvertices++;
2152e529ec8SStefano Zampini         }
216552f7358SJed Brown       }
2179566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
2182e529ec8SStefano Zampini     } else {
2192e529ec8SStefano Zampini       piece.nvertices += dof / dimEmbed;
2202e529ec8SStefano Zampini       piece.nconn += dof / dimEmbed;
2212e529ec8SStefano Zampini     }
222552f7358SJed Brown     piece.ncells++;
223552f7358SJed Brown   }
2249566063dSJacob Faibussowitsch   if (rank == 0) PetscCall(PetscMalloc1(size, &gpiece));
2259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Gather((PetscInt *)&piece, sizeof(piece) / sizeof(PetscInt), MPIU_INT, (PetscInt *)gpiece, sizeof(piece) / sizeof(PetscInt), MPIU_INT, 0, comm));
226552f7358SJed Brown 
227552f7358SJed Brown   /*
228552f7358SJed Brown    * Write file header
229552f7358SJed Brown    */
230dd400576SPatrick Sanan   if (rank == 0) {
2310f2609c8SJed Brown     PetscInt64 boffset = 0;
232552f7358SJed Brown 
233552f7358SJed Brown     for (r = 0; r < size; r++) {
23463a3b9bcSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, "    <Piece NumberOfPoints=\"%" PetscInt_FMT "\" NumberOfCells=\"%" PetscInt_FMT "\">\n", gpiece[r].nvertices, gpiece[r].ncells));
235552f7358SJed Brown       /* Coordinate positions */
2369566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, "      <Points>\n"));
2370f2609c8SJed Brown       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, "        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, boffset));
2386f98e5d0SMatthew G. Knepley       boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0);
2399566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, "      </Points>\n"));
240552f7358SJed Brown       /* Cell connectivity */
2419566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, "      <Cells>\n"));
2420f2609c8SJed Brown       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, "        <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", boffset));
2436f98e5d0SMatthew G. Knepley       boffset += gpiece[r].nconn * sizeof(PetscVTKInt) + (gpiece[r].nconn ? sizeof(PetscInt64) : 0);
2440f2609c8SJed Brown       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, "        <DataArray type=\"Int32\" Name=\"offsets\"      NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", boffset));
2456f98e5d0SMatthew G. Knepley       boffset += gpiece[r].ncells * sizeof(PetscVTKInt) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0);
2460f2609c8SJed Brown       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, "        <DataArray type=\"UInt8\" Name=\"types\"        NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", boffset));
2476f98e5d0SMatthew G. Knepley       boffset += gpiece[r].ncells * sizeof(unsigned char) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0);
2489566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, "      </Cells>\n"));
249552f7358SJed Brown 
250552f7358SJed Brown       /*
251552f7358SJed Brown        * Cell Data headers
252552f7358SJed Brown        */
2539566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, "      <CellData>\n"));
2540f2609c8SJed Brown       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, "        <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", boffset));
2556f98e5d0SMatthew G. Knepley       boffset += gpiece[r].ncells * sizeof(PetscVTKInt) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0);
256552f7358SJed Brown       /* all the vectors */
257552f7358SJed Brown       for (link = vtk->link; link; link = link->next) {
258552f7358SJed Brown         Vec          X       = (Vec)link->vec;
259e630c359SToby Isaac         DM           dmX     = NULL;
2606030a318SStefano Zampini         PetscInt     bs      = 1, nfields, field;
261552f7358SJed Brown         const char  *vecname = "";
262e630c359SToby Isaac         PetscSection section;
263552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
264552f7358SJed 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. */
2659566063dSJacob Faibussowitsch           PetscCall(PetscObjectGetName((PetscObject)X, &vecname));
266552f7358SJed Brown         }
2679566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dmX));
268e630c359SToby Isaac         if (!dmX) dmX = dm;
2699566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)&section));
2709566063dSJacob Faibussowitsch         if (!section) PetscCall(DMGetLocalSection(dmX, &section));
2719566063dSJacob Faibussowitsch         if (cEnd > cStart) PetscCall(PetscSectionGetDof(section, cStart, &bs));
2729566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetNumFields(section, &nfields));
273e630c359SToby Isaac         field = 0;
274e630c359SToby Isaac         if (link->field >= 0) {
275e630c359SToby Isaac           field   = link->field;
276e630c359SToby Isaac           nfields = field + 1;
277e630c359SToby Isaac         }
278e630c359SToby Isaac         for (i = 0; field < (nfields ? nfields : 1); field++) {
2791cfafdd3SJed Brown           PetscInt     fbs, j;
280a00cdb45SToby Isaac           PetscFV      fv = NULL;
281a00cdb45SToby Isaac           PetscObject  f;
282a00cdb45SToby Isaac           PetscClassId fClass;
2830298fd71SBarry Smith           const char  *fieldname = NULL;
2841cfafdd3SJed Brown           char         buf[256];
285e630c359SToby Isaac           PetscBool    vector;
2867ded4cbaSJed Brown           if (nfields) { /* We have user-defined fields/components */
2879566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldDof(section, cStart, field, &fbs));
2889566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldName(section, field, &fieldname));
2897ded4cbaSJed Brown           } else fbs = bs; /* Say we have one field with 'bs' components */
2909566063dSJacob Faibussowitsch           PetscCall(DMGetField(dmX, field, NULL, &f));
2919566063dSJacob Faibussowitsch           PetscCall(PetscObjectGetClassId(f, &fClass));
292ad540459SPierre Jolivet           if (fClass == PETSCFV_CLASSID) fv = (PetscFV)f;
293e630c359SToby Isaac           if (nfields && !fieldname) {
29463a3b9bcSJacob Faibussowitsch             PetscCall(PetscSNPrintf(buf, sizeof(buf), "CellField%" PetscInt_FMT, field));
295552f7358SJed Brown             fieldname = buf;
296552f7358SJed Brown           }
297e630c359SToby Isaac           vector = PETSC_FALSE;
298e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
299e630c359SToby Isaac             vector = PETSC_TRUE;
30063a3b9bcSJacob Faibussowitsch             PetscCheck(fbs <= 3, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_SIZ, "Cell vector fields can have at most 3 components, %" PetscInt_FMT " given", fbs);
301e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
302e630c359SToby Isaac               const char *compName = NULL;
303e630c359SToby Isaac               if (fv) {
3049566063dSJacob Faibussowitsch                 PetscCall(PetscFVGetComponentName(fv, j, &compName));
305e630c359SToby Isaac                 if (compName) break;
306e630c359SToby Isaac               }
307e630c359SToby Isaac             }
308e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
309e630c359SToby Isaac           }
310e630c359SToby Isaac           if (vector) {
311955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
3120f2609c8SJed Brown             PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"%s%s.Re\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset));
3136f98e5d0SMatthew G. Knepley             boffset += gpiece[r].ncells * 3 * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0);
3140f2609c8SJed Brown             PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"%s%s.Im\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset));
3156f98e5d0SMatthew G. Knepley             boffset += gpiece[r].ncells * 3 * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0);
316955d60d1SToby Isaac #else
3170f2609c8SJed Brown             PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset));
3186f98e5d0SMatthew G. Knepley             boffset += gpiece[r].ncells * 3 * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0);
319955d60d1SToby Isaac #endif
320e630c359SToby Isaac             i += fbs;
321e630c359SToby Isaac           } else {
3221cfafdd3SJed Brown             for (j = 0; j < fbs; j++) {
323a00cdb45SToby Isaac               const char *compName = NULL;
324955d60d1SToby Isaac               char        finalname[256];
32548a46eb9SPierre Jolivet               if (fv) PetscCall(PetscFVGetComponentName(fv, j, &compName));
326a00cdb45SToby Isaac               if (compName) {
3279566063dSJacob Faibussowitsch                 PetscCall(PetscSNPrintf(finalname, 255, "%s%s.%s", vecname, fieldname, compName));
3289371c9d4SSatish Balay               } else if (fbs > 1) {
32963a3b9bcSJacob Faibussowitsch                 PetscCall(PetscSNPrintf(finalname, 255, "%s%s.%" PetscInt_FMT, vecname, fieldname, j));
330e630c359SToby Isaac               } else {
3319566063dSJacob Faibussowitsch                 PetscCall(PetscSNPrintf(finalname, 255, "%s%s", vecname, fieldname));
332a00cdb45SToby Isaac               }
333955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
3340f2609c8SJed Brown               PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset));
3356f98e5d0SMatthew G. Knepley               boffset += gpiece[r].ncells * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0);
3360f2609c8SJed Brown               PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset));
3376f98e5d0SMatthew G. Knepley               boffset += gpiece[r].ncells * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0);
338955d60d1SToby Isaac #else
3390f2609c8SJed Brown               PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset));
3406f98e5d0SMatthew G. Knepley               boffset += gpiece[r].ncells * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0);
341955d60d1SToby Isaac #endif
3421cfafdd3SJed Brown               i++;
343552f7358SJed Brown             }
344552f7358SJed Brown           }
345e630c359SToby Isaac         }
34663a3b9bcSJacob Faibussowitsch         //PetscCheck(i == bs,comm,PETSC_ERR_PLIB,"Total number of field components %" PetscInt_FMT " != block size %" PetscInt_FMT,i,bs);
3471cfafdd3SJed Brown       }
3489566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, "      </CellData>\n"));
349552f7358SJed Brown 
350552f7358SJed Brown       /*
351552f7358SJed Brown        * Point Data headers
352552f7358SJed Brown        */
3539566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, "      <PointData>\n"));
354552f7358SJed Brown       for (link = vtk->link; link; link = link->next) {
355552f7358SJed Brown         Vec          X = (Vec)link->vec;
356e630c359SToby Isaac         DM           dmX;
3576030a318SStefano Zampini         PetscInt     bs      = 1, nfields, field;
358552f7358SJed Brown         const char  *vecname = "";
359e630c359SToby Isaac         PetscSection section;
360552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
361552f7358SJed 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. */
3629566063dSJacob Faibussowitsch           PetscCall(PetscObjectGetName((PetscObject)X, &vecname));
363552f7358SJed Brown         }
3649566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dmX));
365e630c359SToby Isaac         if (!dmX) dmX = dm;
3669566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)&section));
3679566063dSJacob Faibussowitsch         if (!section) PetscCall(DMGetLocalSection(dmX, &section));
3689566063dSJacob Faibussowitsch         if (vEnd > vStart) PetscCall(PetscSectionGetDof(section, vStart, &bs));
3699566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetNumFields(section, &nfields));
370e630c359SToby Isaac         field = 0;
371e630c359SToby Isaac         if (link->field >= 0) {
372e630c359SToby Isaac           field   = link->field;
373e630c359SToby Isaac           nfields = field + 1;
374e630c359SToby Isaac         }
375e630c359SToby Isaac         for (i = 0; field < (nfields ? nfields : 1); field++) {
3761cfafdd3SJed Brown           PetscInt    fbs, j;
3771cfafdd3SJed Brown           const char *fieldname = NULL;
3781cfafdd3SJed Brown           char        buf[256];
3797ded4cbaSJed Brown           if (nfields) { /* We have user-defined fields/components */
3809566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldDof(section, vStart, field, &fbs));
3819566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldName(section, field, &fieldname));
3827ded4cbaSJed Brown           } else fbs = bs; /* Say we have one field with 'bs' components */
383e630c359SToby Isaac           if (nfields && !fieldname) {
38463a3b9bcSJacob Faibussowitsch             PetscCall(PetscSNPrintf(buf, sizeof(buf), "PointField%" PetscInt_FMT, field));
3851cfafdd3SJed Brown             fieldname = buf;
3861cfafdd3SJed Brown           }
387e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
38863a3b9bcSJacob Faibussowitsch             PetscCheck(fbs <= 3, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_SIZ, "Point vector fields can have at most 3 components, %" PetscInt_FMT " given", fbs);
389955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
3900f2609c8SJed Brown             PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"%s%s.Re\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset));
3916f98e5d0SMatthew G. Knepley             boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0);
3920f2609c8SJed Brown             PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"%s%s.Im\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset));
3936f98e5d0SMatthew G. Knepley             boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0);
394955d60d1SToby Isaac #else
3950f2609c8SJed Brown             PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset));
3966f98e5d0SMatthew G. Knepley             boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0);
397955d60d1SToby Isaac #endif
398e630c359SToby Isaac           } else {
3991cfafdd3SJed Brown             for (j = 0; j < fbs; j++) {
400b778fa18SValeria Barra               const char *compName = NULL;
401955d60d1SToby Isaac               char        finalname[256];
4029566063dSJacob Faibussowitsch               PetscCall(PetscSectionGetComponentName(section, field, j, &compName));
4039566063dSJacob Faibussowitsch               PetscCall(PetscSNPrintf(finalname, 255, "%s%s.%s", vecname, fieldname, compName));
404955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
4050f2609c8SJed Brown               PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset));
4066f98e5d0SMatthew G. Knepley               boffset += gpiece[r].nvertices * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0);
4070f2609c8SJed Brown               PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset));
4086f98e5d0SMatthew G. Knepley               boffset += gpiece[r].nvertices * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0);
409955d60d1SToby Isaac #else
4100f2609c8SJed Brown               PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset));
4116f98e5d0SMatthew G. Knepley               boffset += gpiece[r].nvertices * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0);
412955d60d1SToby Isaac #endif
413552f7358SJed Brown             }
414552f7358SJed Brown           }
4151cfafdd3SJed Brown         }
416e630c359SToby Isaac       }
4179566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, "      </PointData>\n"));
4189566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, "    </Piece>\n"));
419552f7358SJed Brown     }
420552f7358SJed Brown   }
421552f7358SJed Brown 
4229566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "  </UnstructuredGrid>\n"));
4239566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "  <AppendedData encoding=\"raw\">\n"));
4249566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "_"));
425552f7358SJed Brown 
426dd400576SPatrick Sanan   if (rank == 0) {
427552f7358SJed Brown     PetscInt maxsize = 0;
428552f7358SJed Brown     for (r = 0; r < size; r++) {
429955d60d1SToby Isaac       maxsize = PetscMax(maxsize, (PetscInt)(gpiece[r].nvertices * 3 * sizeof(PetscVTUReal)));
430955d60d1SToby Isaac       maxsize = PetscMax(maxsize, (PetscInt)(gpiece[r].ncells * 3 * sizeof(PetscVTUReal)));
431552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt)(gpiece[r].nconn * sizeof(PetscVTKInt)));
432552f7358SJed Brown     }
4339566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(maxsize, &buffer));
434552f7358SJed Brown   }
435552f7358SJed Brown   for (r = 0; r < size; r++) {
436552f7358SJed Brown     if (r == rank) {
437552f7358SJed Brown       PetscInt nsend;
438552f7358SJed Brown       { /* Position */
4396858538eSMatthew G. Knepley         const PetscScalar *x, *cx = NULL;
440955d60d1SToby Isaac         PetscVTUReal      *y = NULL;
4416858538eSMatthew G. Knepley         Vec                coords, cellCoords;
442955d60d1SToby Isaac         PetscBool          copy;
4432e529ec8SStefano Zampini 
4449566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinatesLocal(dm, &coords));
4459566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(coords, &x));
4466858538eSMatthew G. Knepley         PetscCall(DMGetCellCoordinatesLocal(dm, &cellCoords));
4476858538eSMatthew G. Knepley         if (cellCoords) PetscCall(VecGetArrayRead(cellCoords, &cx));
448955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX)
449955d60d1SToby Isaac         copy = PETSC_TRUE;
450955d60d1SToby Isaac #else
451955d60d1SToby Isaac         copy = (PetscBool)(dimEmbed != 3 || localized || (sizeof(PetscReal) != sizeof(PetscVTUReal)));
452955d60d1SToby Isaac #endif
453955d60d1SToby Isaac         if (copy) {
4549566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(piece.nvertices * 3, &y));
45519003fb5SStefano Zampini           if (localized) {
45619003fb5SStefano Zampini             PetscInt cnt;
45719003fb5SStefano Zampini             for (c = cStart, cnt = 0; c < cEnd; c++) {
45819003fb5SStefano Zampini               PetscInt off, dof;
45919003fb5SStefano Zampini 
4606858538eSMatthew G. Knepley               PetscCall(PetscSectionGetDof(cellCoordSection, c, &dof));
46119003fb5SStefano Zampini               if (!dof) {
46219003fb5SStefano Zampini                 PetscInt *closure = NULL;
46319003fb5SStefano Zampini                 PetscInt  closureSize;
46419003fb5SStefano Zampini 
4659566063dSJacob Faibussowitsch                 PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
46619003fb5SStefano Zampini                 for (v = 0; v < closureSize * 2; v += 2) {
46719003fb5SStefano Zampini                   if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
4689566063dSJacob Faibussowitsch                     PetscCall(PetscSectionGetOffset(coordSection, closure[v], &off));
46919003fb5SStefano Zampini                     if (dimEmbed != 3) {
470955d60d1SToby Isaac                       y[cnt * 3 + 0] = (PetscVTUReal)PetscRealPart(x[off + 0]);
471955d60d1SToby Isaac                       y[cnt * 3 + 1] = (PetscVTUReal)((dimEmbed > 1) ? PetscRealPart(x[off + 1]) : 0.0);
472955d60d1SToby Isaac                       y[cnt * 3 + 2] = (PetscVTUReal)0.0;
47319003fb5SStefano Zampini                     } else {
474955d60d1SToby Isaac                       y[cnt * 3 + 0] = (PetscVTUReal)PetscRealPart(x[off + 0]);
475955d60d1SToby Isaac                       y[cnt * 3 + 1] = (PetscVTUReal)PetscRealPart(x[off + 1]);
476955d60d1SToby Isaac                       y[cnt * 3 + 2] = (PetscVTUReal)PetscRealPart(x[off + 2]);
47719003fb5SStefano Zampini                     }
47819003fb5SStefano Zampini                     cnt++;
47919003fb5SStefano Zampini                   }
48019003fb5SStefano Zampini                 }
4819566063dSJacob Faibussowitsch                 PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
48219003fb5SStefano Zampini               } else {
4836858538eSMatthew G. Knepley                 PetscCall(PetscSectionGetOffset(cellCoordSection, c, &off));
48419003fb5SStefano Zampini                 if (dimEmbed != 3) {
48519003fb5SStefano Zampini                   for (i = 0; i < dof / dimEmbed; i++) {
4866858538eSMatthew G. Knepley                     y[cnt * 3 + 0] = (PetscVTUReal)PetscRealPart(cx[off + i * dimEmbed + 0]);
4876858538eSMatthew G. Knepley                     y[cnt * 3 + 1] = (PetscVTUReal)((dimEmbed > 1) ? PetscRealPart(cx[off + i * dimEmbed + 1]) : 0.0);
488955d60d1SToby Isaac                     y[cnt * 3 + 2] = (PetscVTUReal)0.0;
48919003fb5SStefano Zampini                     cnt++;
49019003fb5SStefano Zampini                   }
49119003fb5SStefano Zampini                 } else {
492ad540459SPierre Jolivet                   for (i = 0; i < dof; i++) y[cnt * 3 + i] = (PetscVTUReal)PetscRealPart(cx[off + i]);
49319003fb5SStefano Zampini                   cnt += dof / dimEmbed;
49419003fb5SStefano Zampini                 }
49519003fb5SStefano Zampini               }
49619003fb5SStefano Zampini             }
49708401ef6SPierre Jolivet             PetscCheck(cnt == piece.nvertices, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match");
49819003fb5SStefano Zampini           } else {
499552f7358SJed Brown             for (i = 0; i < piece.nvertices; i++) {
500955d60d1SToby Isaac               y[i * 3 + 0] = (PetscVTUReal)PetscRealPart(x[i * dimEmbed + 0]);
501955d60d1SToby Isaac               y[i * 3 + 1] = (PetscVTUReal)((dimEmbed > 1) ? PetscRealPart(x[i * dimEmbed + 1]) : 0.);
502955d60d1SToby Isaac               y[i * 3 + 2] = (PetscVTUReal)((dimEmbed > 2) ? PetscRealPart(x[i * dimEmbed + 2]) : 0.);
50319003fb5SStefano Zampini             }
504552f7358SJed Brown           }
505552f7358SJed Brown         }
5062e529ec8SStefano Zampini         nsend = piece.nvertices * 3;
5079566063dSJacob Faibussowitsch         PetscCall(TransferWrite(comm, viewer, fp, r, 0, copy ? (const void *)y : (const void *)x, buffer, nsend, MPIU_VTUREAL, tag));
5089566063dSJacob Faibussowitsch         PetscCall(PetscFree(y));
5099566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(coords, &x));
5106858538eSMatthew G. Knepley         if (cellCoords) PetscCall(VecRestoreArrayRead(cellCoords, &cx));
511552f7358SJed Brown       }
512552f7358SJed Brown       { /* Connectivity, offsets, types */
5133e869605SMatthew G. Knepley         PetscVTKInt  *connectivity = NULL, *offsets = NULL;
5143e869605SMatthew G. Knepley         PetscVTKType *types = NULL;
5159566063dSJacob Faibussowitsch         PetscCall(DMPlexGetVTKConnectivity(dm, localized, &piece, &connectivity, &offsets, &types));
5169566063dSJacob Faibussowitsch         PetscCall(TransferWrite(comm, viewer, fp, r, 0, connectivity, buffer, piece.nconn, MPI_INT, tag));
5179566063dSJacob Faibussowitsch         PetscCall(TransferWrite(comm, viewer, fp, r, 0, offsets, buffer, piece.ncells, MPI_INT, tag));
5189566063dSJacob Faibussowitsch         PetscCall(TransferWrite(comm, viewer, fp, r, 0, types, buffer, piece.ncells, MPI_CHAR, tag));
5199566063dSJacob Faibussowitsch         PetscCall(PetscFree3(connectivity, offsets, types));
520552f7358SJed Brown       }
521552f7358SJed Brown       { /* Owners (cell data) */
522552f7358SJed Brown         PetscVTKInt *owners;
523*c29ce622SStefano Zampini         PetscMPIInt  orank;
524*c29ce622SStefano Zampini 
525*c29ce622SStefano Zampini         PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &orank));
5269566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(piece.ncells, &owners));
527*c29ce622SStefano Zampini         for (i = 0; i < piece.ncells; i++) owners[i] = orank;
5289566063dSJacob Faibussowitsch         PetscCall(TransferWrite(comm, viewer, fp, r, 0, owners, buffer, piece.ncells, MPI_INT, tag));
5299566063dSJacob Faibussowitsch         PetscCall(PetscFree(owners));
530552f7358SJed Brown       }
531552f7358SJed Brown       /* Cell data */
532552f7358SJed Brown       for (link = vtk->link; link; link = link->next) {
533552f7358SJed Brown         Vec                X = (Vec)link->vec;
534e630c359SToby Isaac         DM                 dmX;
535552f7358SJed Brown         const PetscScalar *x;
536955d60d1SToby Isaac         PetscVTUReal      *y;
5376030a318SStefano Zampini         PetscInt           bs      = 1, nfields, field;
538e630c359SToby Isaac         PetscSection       section = NULL;
539e630c359SToby Isaac 
540552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
5419566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dmX));
542e630c359SToby Isaac         if (!dmX) dmX = dm;
5439566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)&section));
5449566063dSJacob Faibussowitsch         if (!section) PetscCall(DMGetLocalSection(dmX, &section));
5459566063dSJacob Faibussowitsch         if (cEnd > cStart) PetscCall(PetscSectionGetDof(section, cStart, &bs));
5469566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetNumFields(section, &nfields));
547e630c359SToby Isaac         field = 0;
548e630c359SToby Isaac         if (link->field >= 0) {
549e630c359SToby Isaac           field   = link->field;
550e630c359SToby Isaac           nfields = field + 1;
551e630c359SToby Isaac         }
5529566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(X, &x));
5539566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(piece.ncells * 3, &y));
554e630c359SToby Isaac         for (i = 0; field < (nfields ? nfields : 1); field++) {
555e630c359SToby Isaac           PetscInt     fbs, j;
556e630c359SToby Isaac           PetscFV      fv = NULL;
557e630c359SToby Isaac           PetscObject  f;
558e630c359SToby Isaac           PetscClassId fClass;
559e630c359SToby Isaac           PetscBool    vector;
5606030a318SStefano Zampini           if (nfields && cEnd > cStart) { /* We have user-defined fields/components */
5619566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldDof(section, cStart, field, &fbs));
562e630c359SToby Isaac           } else fbs = bs; /* Say we have one field with 'bs' components */
5639566063dSJacob Faibussowitsch           PetscCall(DMGetField(dmX, field, NULL, &f));
5649566063dSJacob Faibussowitsch           PetscCall(PetscObjectGetClassId(f, &fClass));
565ad540459SPierre Jolivet           if (fClass == PETSCFV_CLASSID) fv = (PetscFV)f;
566e630c359SToby Isaac           vector = PETSC_FALSE;
567e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
568e630c359SToby Isaac             vector = PETSC_TRUE;
569e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
570e630c359SToby Isaac               const char *compName = NULL;
571e630c359SToby Isaac               if (fv) {
5729566063dSJacob Faibussowitsch                 PetscCall(PetscFVGetComponentName(fv, j, &compName));
573e630c359SToby Isaac                 if (compName) break;
574e630c359SToby Isaac               }
575e630c359SToby Isaac             }
576e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
577e630c359SToby Isaac           }
578e630c359SToby Isaac           if (vector) {
579955d60d1SToby Isaac             PetscInt cnt, l;
580955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
581552f7358SJed Brown               for (c = cStart, cnt = 0; c < cEnd; c++) {
582e630c359SToby Isaac                 const PetscScalar *xpoint;
583e630c359SToby Isaac                 PetscInt           off, j;
584e630c359SToby Isaac 
585552f7358SJed Brown                 if (hasLabel) { /* Ignore some cells */
586552f7358SJed Brown                   PetscInt value;
5879566063dSJacob Faibussowitsch                   PetscCall(DMGetLabelValue(dmX, "vtk", c, &value));
588552f7358SJed Brown                   if (value != 1) continue;
589552f7358SJed Brown                 }
590e630c359SToby Isaac                 if (nfields) {
5919566063dSJacob Faibussowitsch                   PetscCall(PetscSectionGetFieldOffset(section, c, field, &off));
592e630c359SToby Isaac                 } else {
5939566063dSJacob Faibussowitsch                   PetscCall(PetscSectionGetOffset(section, c, &off));
594e630c359SToby Isaac                 }
595e630c359SToby Isaac                 xpoint = &x[off];
596ad540459SPierre Jolivet                 for (j = 0; j < fbs; j++) y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j]));
597e630c359SToby Isaac                 for (; j < 3; j++) y[cnt++] = 0.;
598e630c359SToby Isaac               }
59908401ef6SPierre Jolivet               PetscCheck(cnt == piece.ncells * 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match");
6009566063dSJacob Faibussowitsch               PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.ncells * 3, MPIU_VTUREAL, tag));
601955d60d1SToby Isaac             }
602e630c359SToby Isaac           } else {
603e630c359SToby Isaac             for (i = 0; i < fbs; i++) {
604955d60d1SToby Isaac               PetscInt cnt, l;
605955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
606e630c359SToby Isaac                 for (c = cStart, cnt = 0; c < cEnd; c++) {
607e630c359SToby Isaac                   const PetscScalar *xpoint;
608e630c359SToby Isaac                   PetscInt           off;
609e630c359SToby Isaac 
610e630c359SToby Isaac                   if (hasLabel) { /* Ignore some cells */
611e630c359SToby Isaac                     PetscInt value;
6129566063dSJacob Faibussowitsch                     PetscCall(DMGetLabelValue(dmX, "vtk", c, &value));
613e630c359SToby Isaac                     if (value != 1) continue;
614e630c359SToby Isaac                   }
615e630c359SToby Isaac                   if (nfields) {
6169566063dSJacob Faibussowitsch                     PetscCall(PetscSectionGetFieldOffset(section, c, field, &off));
617e630c359SToby Isaac                   } else {
6189566063dSJacob Faibussowitsch                     PetscCall(PetscSectionGetOffset(section, c, &off));
619e630c359SToby Isaac                   }
620e630c359SToby Isaac                   xpoint   = &x[off];
621955d60d1SToby Isaac                   y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i]));
622552f7358SJed Brown                 }
62308401ef6SPierre Jolivet                 PetscCheck(cnt == piece.ncells, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match");
6249566063dSJacob Faibussowitsch                 PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.ncells, MPIU_VTUREAL, tag));
625955d60d1SToby Isaac               }
626552f7358SJed Brown             }
627e630c359SToby Isaac           }
628e630c359SToby Isaac         }
6299566063dSJacob Faibussowitsch         PetscCall(PetscFree(y));
6309566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(X, &x));
631552f7358SJed Brown       }
632e630c359SToby Isaac       /* point data */
633552f7358SJed Brown       for (link = vtk->link; link; link = link->next) {
634552f7358SJed Brown         Vec                X = (Vec)link->vec;
635e630c359SToby Isaac         DM                 dmX;
636552f7358SJed Brown         const PetscScalar *x;
637955d60d1SToby Isaac         PetscVTUReal      *y;
6386030a318SStefano Zampini         PetscInt           bs      = 1, nfields, field;
639e630c359SToby Isaac         PetscSection       section = NULL;
640e630c359SToby Isaac 
641552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
6429566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dmX));
643e630c359SToby Isaac         if (!dmX) dmX = dm;
6449566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)&section));
6459566063dSJacob Faibussowitsch         if (!section) PetscCall(DMGetLocalSection(dmX, &section));
6469566063dSJacob Faibussowitsch         if (vEnd > vStart) PetscCall(PetscSectionGetDof(section, vStart, &bs));
6479566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetNumFields(section, &nfields));
648e630c359SToby Isaac         field = 0;
649e630c359SToby Isaac         if (link->field >= 0) {
650e630c359SToby Isaac           field   = link->field;
651e630c359SToby Isaac           nfields = field + 1;
652e630c359SToby Isaac         }
6539566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(X, &x));
6549566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(piece.nvertices * 3, &y));
655e630c359SToby Isaac         for (i = 0; field < (nfields ? nfields : 1); field++) {
656e630c359SToby Isaac           PetscInt fbs, j;
6576030a318SStefano Zampini           if (nfields && vEnd > vStart) { /* We have user-defined fields/components */
6589566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldDof(section, vStart, field, &fbs));
659e630c359SToby Isaac           } else fbs = bs; /* Say we have one field with 'bs' components */
660e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
661955d60d1SToby Isaac             PetscInt cnt, l;
662955d60d1SToby Isaac             for (l = 0; l < loops_per_scalar; l++) {
6632e529ec8SStefano Zampini               if (!localized) {
664552f7358SJed Brown                 for (v = vStart, cnt = 0; v < vEnd; v++) {
665e630c359SToby Isaac                   PetscInt           off;
666e630c359SToby Isaac                   const PetscScalar *xpoint;
667e630c359SToby Isaac 
668e630c359SToby Isaac                   if (nfields) {
6699566063dSJacob Faibussowitsch                     PetscCall(PetscSectionGetFieldOffset(section, v, field, &off));
670e630c359SToby Isaac                   } else {
6719566063dSJacob Faibussowitsch                     PetscCall(PetscSectionGetOffset(section, v, &off));
672e630c359SToby Isaac                   }
673e630c359SToby Isaac                   xpoint = &x[off];
674ad540459SPierre Jolivet                   for (j = 0; j < fbs; j++) y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j]));
675e630c359SToby Isaac                   for (; j < 3; j++) y[cnt++] = 0.;
676e630c359SToby Isaac                 }
677e630c359SToby Isaac               } else {
678e630c359SToby Isaac                 for (c = cStart, cnt = 0; c < cEnd; c++) {
679e630c359SToby Isaac                   PetscInt *closure = NULL;
680e630c359SToby Isaac                   PetscInt  closureSize, off;
681e630c359SToby Isaac 
6829566063dSJacob Faibussowitsch                   PetscCall(DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure));
683e630c359SToby Isaac                   for (v = 0, off = 0; v < closureSize * 2; v += 2) {
684e630c359SToby Isaac                     if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
685e630c359SToby Isaac                       PetscInt           voff;
686e630c359SToby Isaac                       const PetscScalar *xpoint;
687e630c359SToby Isaac 
688e630c359SToby Isaac                       if (nfields) {
6899566063dSJacob Faibussowitsch                         PetscCall(PetscSectionGetFieldOffset(section, closure[v], field, &voff));
690e630c359SToby Isaac                       } else {
6919566063dSJacob Faibussowitsch                         PetscCall(PetscSectionGetOffset(section, closure[v], &voff));
692e630c359SToby Isaac                       }
693e630c359SToby Isaac                       xpoint = &x[voff];
694ad540459SPierre Jolivet                       for (j = 0; j < fbs; j++) y[cnt + off++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j]));
695e630c359SToby Isaac                       for (; j < 3; j++) y[cnt + off++] = 0.;
696e630c359SToby Isaac                     }
697e630c359SToby Isaac                   }
698e630c359SToby Isaac                   cnt += off;
6999566063dSJacob Faibussowitsch                   PetscCall(DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure));
700e630c359SToby Isaac                 }
701e630c359SToby Isaac               }
70208401ef6SPierre Jolivet               PetscCheck(cnt == piece.nvertices * 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match");
7039566063dSJacob Faibussowitsch               PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.nvertices * 3, MPIU_VTUREAL, tag));
704955d60d1SToby Isaac             }
705e630c359SToby Isaac           } else {
706e630c359SToby Isaac             for (i = 0; i < fbs; i++) {
707955d60d1SToby Isaac               PetscInt cnt, l;
708955d60d1SToby Isaac               for (l = 0; l < loops_per_scalar; l++) {
709e630c359SToby Isaac                 if (!localized) {
710e630c359SToby Isaac                   for (v = vStart, cnt = 0; v < vEnd; v++) {
711e630c359SToby Isaac                     PetscInt           off;
712e630c359SToby Isaac                     const PetscScalar *xpoint;
713e630c359SToby Isaac 
714e630c359SToby Isaac                     if (nfields) {
7159566063dSJacob Faibussowitsch                       PetscCall(PetscSectionGetFieldOffset(section, v, field, &off));
716e630c359SToby Isaac                     } else {
7179566063dSJacob Faibussowitsch                       PetscCall(PetscSectionGetOffset(section, v, &off));
718e630c359SToby Isaac                     }
719e630c359SToby Isaac                     xpoint   = &x[off];
720955d60d1SToby Isaac                     y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i]));
721552f7358SJed Brown                   }
7222e529ec8SStefano Zampini                 } else {
7232e529ec8SStefano Zampini                   for (c = cStart, cnt = 0; c < cEnd; c++) {
7242e529ec8SStefano Zampini                     PetscInt *closure = NULL;
7252e529ec8SStefano Zampini                     PetscInt  closureSize, off;
7262e529ec8SStefano Zampini 
7279566063dSJacob Faibussowitsch                     PetscCall(DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure));
7282e529ec8SStefano Zampini                     for (v = 0, off = 0; v < closureSize * 2; v += 2) {
7292e529ec8SStefano Zampini                       if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
730e630c359SToby Isaac                         PetscInt           voff;
731e630c359SToby Isaac                         const PetscScalar *xpoint;
7322e529ec8SStefano Zampini 
733e630c359SToby Isaac                         if (nfields) {
7349566063dSJacob Faibussowitsch                           PetscCall(PetscSectionGetFieldOffset(section, closure[v], field, &voff));
735e630c359SToby Isaac                         } else {
7369566063dSJacob Faibussowitsch                           PetscCall(PetscSectionGetOffset(section, closure[v], &voff));
737e630c359SToby Isaac                         }
738e630c359SToby Isaac                         xpoint         = &x[voff];
739955d60d1SToby Isaac                         y[cnt + off++] = (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i]));
7402e529ec8SStefano Zampini                       }
7412e529ec8SStefano Zampini                     }
7429bda96f6SStefano Zampini                     cnt += off;
7439566063dSJacob Faibussowitsch                     PetscCall(DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure));
7442e529ec8SStefano Zampini                   }
7452e529ec8SStefano Zampini                 }
74608401ef6SPierre Jolivet                 PetscCheck(cnt == piece.nvertices, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match");
7479566063dSJacob Faibussowitsch                 PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.nvertices, MPIU_VTUREAL, tag));
748955d60d1SToby Isaac               }
749552f7358SJed Brown             }
750e630c359SToby Isaac           }
751e630c359SToby Isaac         }
7529566063dSJacob Faibussowitsch         PetscCall(PetscFree(y));
7539566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(X, &x));
754552f7358SJed Brown       }
755dd400576SPatrick Sanan     } else if (rank == 0) {
756955d60d1SToby Isaac       PetscInt l;
757955d60d1SToby Isaac 
7589566063dSJacob Faibussowitsch       PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nvertices * 3, MPIU_VTUREAL, tag)); /* positions */
7599566063dSJacob Faibussowitsch       PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nconn, MPI_INT, tag));              /* connectivity */
7609566063dSJacob Faibussowitsch       PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPI_INT, tag));             /* offsets */
7619566063dSJacob Faibussowitsch       PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPI_CHAR, tag));            /* types */
7629566063dSJacob Faibussowitsch       PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPI_INT, tag));             /* owner rank (cells) */
763552f7358SJed Brown       /* all cell data */
764552f7358SJed Brown       for (link = vtk->link; link; link = link->next) {
765e630c359SToby Isaac         Vec          X  = (Vec)link->vec;
7666030a318SStefano Zampini         PetscInt     bs = 1, nfields, field;
767e630c359SToby Isaac         DM           dmX;
768e630c359SToby Isaac         PetscSection section = NULL;
769e630c359SToby Isaac 
770552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
7719566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dmX));
772e630c359SToby Isaac         if (!dmX) dmX = dm;
7739566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)&section));
7749566063dSJacob Faibussowitsch         if (!section) PetscCall(DMGetLocalSection(dmX, &section));
7759566063dSJacob Faibussowitsch         if (cEnd > cStart) PetscCall(PetscSectionGetDof(section, cStart, &bs));
7769566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetNumFields(section, &nfields));
777e630c359SToby Isaac         field = 0;
778e630c359SToby Isaac         if (link->field >= 0) {
779e630c359SToby Isaac           field   = link->field;
780e630c359SToby Isaac           nfields = field + 1;
781e630c359SToby Isaac         }
782e630c359SToby Isaac         for (i = 0; field < (nfields ? nfields : 1); field++) {
783e630c359SToby Isaac           PetscInt     fbs, j;
784e630c359SToby Isaac           PetscFV      fv = NULL;
785e630c359SToby Isaac           PetscObject  f;
786e630c359SToby Isaac           PetscClassId fClass;
787e630c359SToby Isaac           PetscBool    vector;
7886030a318SStefano Zampini           if (nfields && cEnd > cStart) { /* We have user-defined fields/components */
7899566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldDof(section, cStart, field, &fbs));
790e630c359SToby Isaac           } else fbs = bs; /* Say we have one field with 'bs' components */
7919566063dSJacob Faibussowitsch           PetscCall(DMGetField(dmX, field, NULL, &f));
7929566063dSJacob Faibussowitsch           PetscCall(PetscObjectGetClassId(f, &fClass));
793ad540459SPierre Jolivet           if (fClass == PETSCFV_CLASSID) fv = (PetscFV)f;
794e630c359SToby Isaac           vector = PETSC_FALSE;
795e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
796e630c359SToby Isaac             vector = PETSC_TRUE;
797e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
798e630c359SToby Isaac               const char *compName = NULL;
799e630c359SToby Isaac               if (fv) {
8009566063dSJacob Faibussowitsch                 PetscCall(PetscFVGetComponentName(fv, j, &compName));
801e630c359SToby Isaac                 if (compName) break;
802e630c359SToby Isaac               }
803e630c359SToby Isaac             }
804e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
805e630c359SToby Isaac           }
806e630c359SToby Isaac           if (vector) {
80748a46eb9SPierre Jolivet             for (l = 0; l < loops_per_scalar; l++) PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells * 3, MPIU_VTUREAL, tag));
808e630c359SToby Isaac           } else {
809e630c359SToby Isaac             for (i = 0; i < fbs; i++) {
81048a46eb9SPierre Jolivet               for (l = 0; l < loops_per_scalar; l++) PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPIU_VTUREAL, tag));
811552f7358SJed Brown             }
812552f7358SJed Brown           }
813e630c359SToby Isaac         }
814e630c359SToby Isaac       }
815552f7358SJed Brown       /* all point data */
816552f7358SJed Brown       for (link = vtk->link; link; link = link->next) {
817e630c359SToby Isaac         Vec          X = (Vec)link->vec;
818e630c359SToby Isaac         DM           dmX;
8196030a318SStefano Zampini         PetscInt     bs      = 1, nfields, field;
820e630c359SToby Isaac         PetscSection section = NULL;
821e630c359SToby Isaac 
822552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
8239566063dSJacob Faibussowitsch         PetscCall(VecGetDM(X, &dmX));
824e630c359SToby Isaac         if (!dmX) dmX = dm;
8259566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)&section));
8269566063dSJacob Faibussowitsch         if (!section) PetscCall(DMGetLocalSection(dmX, &section));
8279566063dSJacob Faibussowitsch         if (vEnd > vStart) PetscCall(PetscSectionGetDof(section, vStart, &bs));
8289566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetNumFields(section, &nfields));
829e630c359SToby Isaac         field = 0;
830e630c359SToby Isaac         if (link->field >= 0) {
831e630c359SToby Isaac           field   = link->field;
832e630c359SToby Isaac           nfields = field + 1;
833e630c359SToby Isaac         }
834e630c359SToby Isaac         for (i = 0; field < (nfields ? nfields : 1); field++) {
835e630c359SToby Isaac           PetscInt fbs;
8366030a318SStefano Zampini           if (nfields && vEnd > vStart) { /* We have user-defined fields/components */
8379566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldDof(section, vStart, field, &fbs));
838e630c359SToby Isaac           } else fbs = bs; /* Say we have one field with 'bs' components */
839e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
84048a46eb9SPierre Jolivet             for (l = 0; l < loops_per_scalar; l++) PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nvertices * 3, MPIU_VTUREAL, tag));
841e630c359SToby Isaac           } else {
842e630c359SToby Isaac             for (i = 0; i < fbs; i++) {
84348a46eb9SPierre Jolivet               for (l = 0; l < loops_per_scalar; l++) PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nvertices, MPIU_VTUREAL, tag));
844552f7358SJed Brown             }
845552f7358SJed Brown           }
846552f7358SJed Brown         }
847552f7358SJed Brown       }
848e630c359SToby Isaac     }
849e630c359SToby Isaac   }
8509566063dSJacob Faibussowitsch   PetscCall(PetscFree(gpiece));
8519566063dSJacob Faibussowitsch   PetscCall(PetscFree(buffer));
8529566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "\n  </AppendedData>\n"));
8539566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "</VTKFile>\n"));
8549566063dSJacob Faibussowitsch   PetscCall(PetscFClose(comm, fp));
855*c29ce622SStefano Zampini finalize:
8563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
857552f7358SJed Brown }
858