xref: /petsc/src/dm/impls/plex/plexvtu.c (revision e630c359805a7909577b2c5515bf56b45bba1d70)
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 
10552f7358SJed Brown #if defined(PETSC_USE_REAL_SINGLE)
11552f7358SJed Brown static const char precision[] = "Float32";
12552f7358SJed Brown #elif defined(PETSC_USE_REAL_DOUBLE)
13552f7358SJed Brown static const char precision[] = "Float64";
14552f7358SJed Brown #else
15552f7358SJed Brown static const char precision[] = "UnknownPrecision";
16552f7358SJed Brown #endif
17552f7358SJed Brown 
1894fbd55eSBarry Smith static PetscErrorCode TransferWrite(PetscViewer viewer,FILE *fp,PetscMPIInt srank,PetscMPIInt root,const void *send,void *recv,PetscMPIInt count,MPI_Datatype mpidatatype,PetscMPIInt tag)
19552f7358SJed Brown {
20552f7358SJed Brown   PetscMPIInt    rank;
21552f7358SJed Brown   PetscErrorCode ierr;
22ce94432eSBarry Smith   MPI_Comm       comm;
23552f7358SJed Brown 
24552f7358SJed Brown   PetscFunctionBegin;
25ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
26552f7358SJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
27552f7358SJed Brown 
28552f7358SJed Brown   if (rank == srank && rank != root) {
29552f7358SJed Brown     ierr = MPI_Send((void*)send,count,mpidatatype,root,tag,comm);CHKERRQ(ierr);
30552f7358SJed Brown   } else if (rank == root) {
31552f7358SJed Brown     const void *buffer;
32552f7358SJed Brown     if (root == srank) {        /* self */
33552f7358SJed Brown       buffer = send;
34552f7358SJed Brown     } else {
35552f7358SJed Brown       MPI_Status  status;
36552f7358SJed Brown       PetscMPIInt nrecv;
37552f7358SJed Brown       ierr = MPI_Recv(recv,count,mpidatatype,srank,tag,comm,&status);CHKERRQ(ierr);
38552f7358SJed Brown       ierr = MPI_Get_count(&status,mpidatatype,&nrecv);CHKERRQ(ierr);
39552f7358SJed Brown       if (count != nrecv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
40552f7358SJed Brown       buffer = recv;
41552f7358SJed Brown     }
4294fbd55eSBarry Smith     ierr = PetscViewerVTKFWrite(viewer,fp,buffer,count,mpidatatype);CHKERRQ(ierr);
43552f7358SJed Brown   }
44552f7358SJed Brown   PetscFunctionReturn(0);
45552f7358SJed Brown }
46552f7358SJed Brown 
472e529ec8SStefano Zampini static PetscErrorCode DMPlexGetVTKConnectivity(DM dm, PetscBool localized, PieceInfo *piece,PetscVTKInt **oconn,PetscVTKInt **ooffsets,PetscVTKType **otypes)
48552f7358SJed Brown {
49552f7358SJed Brown   PetscErrorCode ierr;
502e529ec8SStefano Zampini   PetscSection   coordSection;
51552f7358SJed Brown   PetscVTKInt    *conn,*offsets;
52552f7358SJed Brown   PetscVTKType   *types;
532f029394SStefano Zampini   PetscInt       dim,vStart,vEnd,cStart,cEnd,pStart,pEnd,cellHeight,numLabelCells,hasLabel,c,v,countcell,countconn;
54552f7358SJed Brown 
55552f7358SJed Brown   PetscFunctionBegin;
56dcca6d9dSJed Brown   ierr = PetscMalloc3(piece->nconn,&conn,piece->ncells,&offsets,piece->ncells,&types);CHKERRQ(ierr);
572e529ec8SStefano Zampini   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
58c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
59552f7358SJed Brown   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
60552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
61552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
62552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
63c58f1c22SToby Isaac   ierr     = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
64552f7358SJed Brown   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
65552f7358SJed Brown 
66552f7358SJed Brown   countcell = 0;
67552f7358SJed Brown   countconn = 0;
68552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
6919003fb5SStefano Zampini     PetscInt nverts,dof = 0,celltype,startoffset,nC=0;
70552f7358SJed Brown 
71552f7358SJed Brown     if (hasLabel) {
72552f7358SJed Brown       PetscInt value;
73552f7358SJed Brown 
74c58f1c22SToby Isaac       ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
75552f7358SJed Brown       if (value != 1) continue;
76552f7358SJed Brown     }
77552f7358SJed Brown     startoffset = countconn;
7819003fb5SStefano Zampini     if (localized) {
7919003fb5SStefano Zampini       ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
8019003fb5SStefano Zampini     }
8119003fb5SStefano Zampini     if (!dof) {
822e529ec8SStefano Zampini       PetscInt *closure = NULL;
832e529ec8SStefano Zampini       PetscInt  closureSize;
842e529ec8SStefano Zampini 
85552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
86552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
87552f7358SJed Brown         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
8819003fb5SStefano Zampini           if (!localized) conn[countconn++] = closure[v] - vStart;
8919003fb5SStefano Zampini           else conn[countconn++] = startoffset + nC;
90724b5320SMatthew G. Knepley           ++nC;
91552f7358SJed Brown         }
92552f7358SJed Brown       }
93552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
942e529ec8SStefano Zampini     } else {
952e529ec8SStefano Zampini       for (nC = 0; nC < dof/dim; nC++) conn[countconn++] = startoffset + nC;
962e529ec8SStefano Zampini     }
97724b5320SMatthew G. Knepley     ierr = DMPlexInvertCell(dim, nC, &conn[countconn-nC]);CHKERRQ(ierr);
988865f1eaSKarl Rupp 
99552f7358SJed Brown     offsets[countcell] = countconn;
1008865f1eaSKarl Rupp 
101552f7358SJed Brown     nverts = countconn - startoffset;
102de0a4605SMatthew G. Knepley     ierr   = DMPlexVTKGetCellType_Internal(dm,dim,nverts,&celltype);CHKERRQ(ierr);
1038865f1eaSKarl Rupp 
104552f7358SJed Brown     types[countcell] = celltype;
105552f7358SJed Brown     countcell++;
106552f7358SJed Brown   }
107552f7358SJed Brown   if (countcell != piece->ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count");
108552f7358SJed Brown   if (countconn != piece->nconn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count");
109552f7358SJed Brown   *oconn    = conn;
110552f7358SJed Brown   *ooffsets = offsets;
111552f7358SJed Brown   *otypes   = types;
112552f7358SJed Brown   PetscFunctionReturn(0);
113552f7358SJed Brown }
114552f7358SJed Brown 
115552f7358SJed Brown /*
116552f7358SJed Brown   Write all fields that have been provided to the viewer
117552f7358SJed Brown   Multi-block XML format with binary appended data.
118552f7358SJed Brown */
119552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer)
120552f7358SJed Brown {
121ce94432eSBarry Smith   MPI_Comm                 comm;
1222e529ec8SStefano Zampini   PetscSection             coordSection;
123552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
124552f7358SJed Brown   PetscViewerVTKObjectLink link;
125552f7358SJed Brown   FILE                     *fp;
126552f7358SJed Brown   PetscMPIInt              rank,size,tag;
127552f7358SJed Brown   PetscErrorCode           ierr;
1282f029394SStefano Zampini   PetscInt                 dimEmbed,cellHeight,cStart,cEnd,vStart,vEnd,numLabelCells,hasLabel,c,v,r,i;
1292e529ec8SStefano Zampini   PetscBool                localized;
1300298fd71SBarry Smith   PieceInfo                piece,*gpiece = NULL;
1310298fd71SBarry Smith   void                     *buffer = NULL;
13230815ce0SLisandro Dalcin   const char               *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian";
133552f7358SJed Brown 
134552f7358SJed Brown   PetscFunctionBegin;
135ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
136552f7358SJed Brown #if defined(PETSC_USE_COMPLEX)
137ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Complex values not supported");
138552f7358SJed Brown #endif
139552f7358SJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
140552f7358SJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
141552f7358SJed Brown   tag  = ((PetscObject)viewer)->tag;
142552f7358SJed Brown 
143552f7358SJed Brown   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
144552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
14530815ce0SLisandro Dalcin   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"%s\">\n", byte_order);CHKERRQ(ierr);
146552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <UnstructuredGrid>\n");CHKERRQ(ierr);
147552f7358SJed Brown 
1483baa072aSToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
149552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
150552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
151552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
152c58f1c22SToby Isaac   ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
1532e529ec8SStefano Zampini   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
1542e529ec8SStefano Zampini   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1558865f1eaSKarl Rupp 
156552f7358SJed Brown   hasLabel        = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
1572e529ec8SStefano Zampini   piece.nvertices = 0;
158552f7358SJed Brown   piece.ncells    = 0;
159552f7358SJed Brown   piece.nconn     = 0;
1602e529ec8SStefano Zampini   if (!localized) piece.nvertices = vEnd - vStart;
161552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
16219003fb5SStefano Zampini     PetscInt dof = 0;
163552f7358SJed Brown     if (hasLabel) {
164552f7358SJed Brown       PetscInt value;
165552f7358SJed Brown 
166c58f1c22SToby Isaac       ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
167552f7358SJed Brown       if (value != 1) continue;
168552f7358SJed Brown     }
16919003fb5SStefano Zampini     if (localized) {
17019003fb5SStefano Zampini       ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
17119003fb5SStefano Zampini     }
17219003fb5SStefano Zampini     if (!dof) {
1732e529ec8SStefano Zampini       PetscInt *closure = NULL;
1742e529ec8SStefano Zampini       PetscInt closureSize;
1752e529ec8SStefano Zampini 
176552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
177552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
1782e529ec8SStefano Zampini         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
1792e529ec8SStefano Zampini           piece.nconn++;
18019003fb5SStefano Zampini           if (localized) piece.nvertices++;
1812e529ec8SStefano Zampini         }
182552f7358SJed Brown       }
183552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1842e529ec8SStefano Zampini     } else {
1852e529ec8SStefano Zampini       piece.nvertices += dof/dimEmbed;
1862e529ec8SStefano Zampini       piece.nconn     += dof/dimEmbed;
1872e529ec8SStefano Zampini     }
188552f7358SJed Brown     piece.ncells++;
189552f7358SJed Brown   }
190785e854fSJed Brown   if (!rank) {ierr = PetscMalloc1(size,&gpiece);CHKERRQ(ierr);}
191552f7358SJed Brown   ierr = MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm);CHKERRQ(ierr);
192552f7358SJed Brown 
193552f7358SJed Brown   /*
194552f7358SJed Brown    * Write file header
195552f7358SJed Brown    */
196552f7358SJed Brown   if (!rank) {
197552f7358SJed Brown     PetscInt boffset = 0;
198552f7358SJed Brown 
199552f7358SJed Brown     for (r=0; r<size; r++) {
200552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells);CHKERRQ(ierr);
201552f7358SJed Brown       /* Coordinate positions */
202552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Points>\n");CHKERRQ(ierr);
203552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
204552f7358SJed Brown       boffset += gpiece[r].nvertices*3*sizeof(PetscScalar) + sizeof(int);
205552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Points>\n");CHKERRQ(ierr);
206552f7358SJed Brown       /* Cell connectivity */
207552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Cells>\n");CHKERRQ(ierr);
208552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
209552f7358SJed Brown       boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int);
210552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"offsets\"      NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
211552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int);
212552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"UInt8\" Name=\"types\"        NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
213552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int);
214552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Cells>\n");CHKERRQ(ierr);
215552f7358SJed Brown 
216552f7358SJed Brown       /*
217552f7358SJed Brown        * Cell Data headers
218552f7358SJed Brown        */
219552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <CellData>\n");CHKERRQ(ierr);
220552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
221552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(int) + sizeof(int);
222552f7358SJed Brown       /* all the vectors */
223552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
224552f7358SJed Brown         Vec        X = (Vec)link->vec;
225*e630c359SToby Isaac         DM         dmX = NULL;
2261cfafdd3SJed Brown         PetscInt   bs,nfields,field;
227552f7358SJed Brown         const char *vecname = "";
228*e630c359SToby Isaac         PetscSection section;
229552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
230552f7358SJed 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. */
231552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
232552f7358SJed Brown         }
233*e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
234*e630c359SToby Isaac         if (!dmX) dmX = dm;
235*e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
236*e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
237*e630c359SToby Isaac         ierr = PetscSectionGetDof(section,cStart,&bs);CHKERRQ(ierr);
238*e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
239*e630c359SToby Isaac         field = 0;
240*e630c359SToby Isaac         if (link->field >= 0) {
241*e630c359SToby Isaac           field = link->field;
242*e630c359SToby Isaac           nfields = field + 1;
243*e630c359SToby Isaac         }
244*e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
2451cfafdd3SJed Brown           PetscInt     fbs,j;
246a00cdb45SToby Isaac           PetscFV      fv = NULL;
247a00cdb45SToby Isaac           PetscObject  f;
248a00cdb45SToby Isaac           PetscClassId fClass;
2490298fd71SBarry Smith           const char *fieldname = NULL;
2501cfafdd3SJed Brown           char       buf[256];
251*e630c359SToby Isaac           PetscBool    vector;
2527ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
253*e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,cStart,field,&fbs);CHKERRQ(ierr);
254*e630c359SToby Isaac             ierr = PetscSectionGetFieldName(section,field,&fieldname);CHKERRQ(ierr);
2557ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
256*e630c359SToby Isaac           ierr = DMGetField(dmX,field,NULL,&f);CHKERRQ(ierr);
257a00cdb45SToby Isaac           ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr);
258a00cdb45SToby Isaac           if (fClass == PETSCFV_CLASSID) {
259a00cdb45SToby Isaac             fv = (PetscFV) f;
260a00cdb45SToby Isaac           }
261*e630c359SToby Isaac           if (nfields && !fieldname) {
2621cfafdd3SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"CellField%D",field);CHKERRQ(ierr);
263552f7358SJed Brown             fieldname = buf;
264552f7358SJed Brown           }
265*e630c359SToby Isaac           vector = PETSC_FALSE;
266*e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
267*e630c359SToby Isaac             vector = PETSC_TRUE;
268*e630c359SToby Isaac             if (fbs > 3) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_SIZ,"Cell vector fields can have at most 3 components, %D given\n", fbs);
269*e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
270*e630c359SToby Isaac               const char *compName = NULL;
271*e630c359SToby Isaac               if (fv) {
272*e630c359SToby Isaac                 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
273*e630c359SToby Isaac                 if (compName) break;
274*e630c359SToby Isaac               }
275*e630c359SToby Isaac             }
276*e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
277*e630c359SToby Isaac           }
278*e630c359SToby Isaac           if (vector) {
279*e630c359SToby Isaac             ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
280*e630c359SToby Isaac             boffset += gpiece[r].ncells*3*sizeof(PetscScalar) + sizeof(int);
281*e630c359SToby Isaac             i+=fbs;
282*e630c359SToby Isaac           } else {
2831cfafdd3SJed Brown             for (j=0; j<fbs; j++) {
284a00cdb45SToby Isaac               const char *compName = NULL;
285a00cdb45SToby Isaac               if (fv) {
286a00cdb45SToby Isaac                 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
287a00cdb45SToby Isaac               }
288a00cdb45SToby Isaac               if (compName) {
289a00cdb45SToby Isaac                 ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,compName,boffset);CHKERRQ(ierr);
290a00cdb45SToby Isaac               }
291*e630c359SToby Isaac               else if (fbs > 1) {
2921cfafdd3SJed Brown                 ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s.%D\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,j,boffset);CHKERRQ(ierr);
293*e630c359SToby Isaac               } else {
294*e630c359SToby Isaac                 ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
295a00cdb45SToby Isaac               }
296552f7358SJed Brown               boffset += gpiece[r].ncells*sizeof(PetscScalar) + sizeof(int);
2971cfafdd3SJed Brown               i++;
298552f7358SJed Brown             }
299552f7358SJed Brown           }
300*e630c359SToby Isaac         }
3011cfafdd3SJed Brown         if (i != bs) SETERRQ2(comm,PETSC_ERR_PLIB,"Total number of field components %D != block size %D",i,bs);
3021cfafdd3SJed Brown       }
303552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </CellData>\n");CHKERRQ(ierr);
304552f7358SJed Brown 
305552f7358SJed Brown       /*
306552f7358SJed Brown        * Point Data headers
307552f7358SJed Brown        */
308552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      <PointData>\n");CHKERRQ(ierr);
309552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
310552f7358SJed Brown         Vec        X = (Vec)link->vec;
311*e630c359SToby Isaac         DM         dmX;
3121cfafdd3SJed Brown         PetscInt   bs,nfields,field;
313552f7358SJed Brown         const char *vecname = "";
314*e630c359SToby Isaac         PetscSection section;
315552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
316552f7358SJed 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. */
317552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
318552f7358SJed Brown         }
319*e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
320*e630c359SToby Isaac         if (!dmX) dmX = dm;
321*e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
322*e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
323*e630c359SToby Isaac         ierr = PetscSectionGetDof(section,vStart,&bs);CHKERRQ(ierr);
324*e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
325*e630c359SToby Isaac         field = 0;
326*e630c359SToby Isaac         if (link->field >= 0) {
327*e630c359SToby Isaac           field = link->field;
328*e630c359SToby Isaac           nfields = field + 1;
329*e630c359SToby Isaac         }
330*e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
3311cfafdd3SJed Brown           PetscInt   fbs,j;
3321cfafdd3SJed Brown           const char *fieldname = NULL;
3331cfafdd3SJed Brown           char       buf[256];
3347ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
335*e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,vStart,field,&fbs);CHKERRQ(ierr);
336*e630c359SToby Isaac             ierr = PetscSectionGetFieldName(section,field,&fieldname);CHKERRQ(ierr);
3377ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
338*e630c359SToby Isaac           if (nfields && !fieldname) {
3391cfafdd3SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"PointField%D",field);CHKERRQ(ierr);
3401cfafdd3SJed Brown             fieldname = buf;
3411cfafdd3SJed Brown           }
342*e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
343*e630c359SToby Isaac             if (fbs > 3) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_SIZ,"Point vector fields can have at most 3 components, %D given\n", fbs);
344*e630c359SToby Isaac             ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
345*e630c359SToby Isaac             boffset += gpiece[r].nvertices*3*sizeof(PetscScalar) + sizeof(int);
346*e630c359SToby Isaac           } else {
3471cfafdd3SJed Brown             for (j=0; j<fbs; j++) {
348*e630c359SToby Isaac               if (fbs > 1) {
3491cfafdd3SJed Brown                 ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s.%D\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,j,boffset);CHKERRQ(ierr);
350*e630c359SToby Isaac               } else {
351*e630c359SToby Isaac                 ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
352*e630c359SToby Isaac               }
353552f7358SJed Brown               boffset += gpiece[r].nvertices*sizeof(PetscScalar) + sizeof(int);
354552f7358SJed Brown             }
355552f7358SJed Brown           }
3561cfafdd3SJed Brown         }
357*e630c359SToby Isaac       }
358552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </PointData>\n");CHKERRQ(ierr);
359552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    </Piece>\n");CHKERRQ(ierr);
360552f7358SJed Brown     }
361552f7358SJed Brown   }
362552f7358SJed Brown 
363552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  </UnstructuredGrid>\n");CHKERRQ(ierr);
364552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
365552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
366552f7358SJed Brown 
367552f7358SJed Brown   if (!rank) {
368552f7358SJed Brown     PetscInt maxsize = 0;
369552f7358SJed Brown     for (r=0; r<size; r++) {
370552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscScalar)));
371*e630c359SToby Isaac       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*3*sizeof(PetscScalar)));
372552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt)));
373552f7358SJed Brown     }
374552f7358SJed Brown     ierr = PetscMalloc(maxsize,&buffer);CHKERRQ(ierr);
375552f7358SJed Brown   }
376552f7358SJed Brown   for (r=0; r<size; r++) {
377552f7358SJed Brown     if (r == rank) {
378552f7358SJed Brown       PetscInt nsend;
379552f7358SJed Brown       {                         /* Position */
380552f7358SJed Brown         const PetscScalar *x;
3810298fd71SBarry Smith         PetscScalar       *y = NULL;
382552f7358SJed Brown         Vec               coords;
3832e529ec8SStefano Zampini 
384552f7358SJed Brown         ierr  = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
385552f7358SJed Brown         ierr  = VecGetArrayRead(coords,&x);CHKERRQ(ierr);
38619003fb5SStefano Zampini         if (dimEmbed != 3 || localized) {
387785e854fSJed Brown           ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr);
38819003fb5SStefano Zampini           if (localized) {
38919003fb5SStefano Zampini             PetscInt cnt;
39019003fb5SStefano Zampini             for (c=cStart,cnt=0; c<cEnd; c++) {
39119003fb5SStefano Zampini               PetscInt off, dof;
39219003fb5SStefano Zampini 
39319003fb5SStefano Zampini               ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
39419003fb5SStefano Zampini               if (!dof) {
39519003fb5SStefano Zampini                 PetscInt *closure = NULL;
39619003fb5SStefano Zampini                 PetscInt closureSize;
39719003fb5SStefano Zampini 
39819003fb5SStefano Zampini                 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
39919003fb5SStefano Zampini                 for (v = 0; v < closureSize*2; v += 2) {
40019003fb5SStefano Zampini                   if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
40119003fb5SStefano Zampini                     ierr = PetscSectionGetOffset(coordSection, closure[v], &off);CHKERRQ(ierr);
40219003fb5SStefano Zampini                     if (dimEmbed != 3) {
40319003fb5SStefano Zampini                       y[cnt*3+0] = x[off+0];
40419003fb5SStefano Zampini                       y[cnt*3+1] = (dimEmbed > 1) ? x[off+1] : 0.0;
40519003fb5SStefano Zampini                       y[cnt*3+2] = 0.0;
40619003fb5SStefano Zampini                     } else {
40719003fb5SStefano Zampini                       y[cnt*3+0] = x[off+0];
40819003fb5SStefano Zampini                       y[cnt*3+1] = x[off+1];
40919003fb5SStefano Zampini                       y[cnt*3+2] = x[off+2];
41019003fb5SStefano Zampini                     }
41119003fb5SStefano Zampini                     cnt++;
41219003fb5SStefano Zampini                   }
41319003fb5SStefano Zampini                 }
41419003fb5SStefano Zampini                 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
41519003fb5SStefano Zampini               } else {
41619003fb5SStefano Zampini                 ierr = PetscSectionGetOffset(coordSection, c, &off);CHKERRQ(ierr);
41719003fb5SStefano Zampini                 if (dimEmbed != 3) {
41819003fb5SStefano Zampini                   for (i=0; i<dof/dimEmbed; i++) {
41919003fb5SStefano Zampini                     y[cnt*3+0] = x[off + i*dimEmbed + 0];
42019003fb5SStefano Zampini                     y[cnt*3+1] = (dimEmbed > 1) ? x[off + i*dimEmbed + 1] : 0.0;
42119003fb5SStefano Zampini                     y[cnt*3+2] = 0.0;
42219003fb5SStefano Zampini                     cnt++;
42319003fb5SStefano Zampini                   }
42419003fb5SStefano Zampini                 } else {
42519003fb5SStefano Zampini                   for (i=0; i<dof; i ++) {
42619003fb5SStefano Zampini                     y[cnt*3+i] = x[off + i];
42719003fb5SStefano Zampini                   }
42819003fb5SStefano Zampini                   cnt += dof/dimEmbed;
42919003fb5SStefano Zampini                 }
43019003fb5SStefano Zampini               }
43119003fb5SStefano Zampini             }
43219003fb5SStefano Zampini             if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
43319003fb5SStefano Zampini           } else {
434552f7358SJed Brown             for (i=0; i<piece.nvertices; i++) {
4353baa072aSToby Isaac               y[i*3+0] = x[i*dimEmbed+0];
4363baa072aSToby Isaac               y[i*3+1] = (dimEmbed > 1) ? x[i*dimEmbed+1] : 0;
43719003fb5SStefano Zampini               y[i*3+2] = 0.0;
43819003fb5SStefano Zampini             }
439552f7358SJed Brown           }
440552f7358SJed Brown         }
4412e529ec8SStefano Zampini         nsend = piece.nvertices*3;
44294fbd55eSBarry Smith         ierr  = TransferWrite(viewer,fp,r,0,y ? y : x,buffer,nsend,MPIU_SCALAR,tag);CHKERRQ(ierr);
443552f7358SJed Brown         ierr  = PetscFree(y);CHKERRQ(ierr);
444552f7358SJed Brown         ierr  = VecRestoreArrayRead(coords,&x);CHKERRQ(ierr);
445552f7358SJed Brown       }
446552f7358SJed Brown       {                           /* Connectivity, offsets, types */
4473e869605SMatthew G. Knepley         PetscVTKInt  *connectivity = NULL, *offsets = NULL;
4483e869605SMatthew G. Knepley         PetscVTKType *types = NULL;
4492e529ec8SStefano Zampini         ierr = DMPlexGetVTKConnectivity(dm,localized,&piece,&connectivity,&offsets,&types);CHKERRQ(ierr);
45094fbd55eSBarry Smith         ierr = TransferWrite(viewer,fp,r,0,connectivity,buffer,piece.nconn,MPI_INT,tag);CHKERRQ(ierr);
45194fbd55eSBarry Smith         ierr = TransferWrite(viewer,fp,r,0,offsets,buffer,piece.ncells,MPI_INT,tag);CHKERRQ(ierr);
45294fbd55eSBarry Smith         ierr = TransferWrite(viewer,fp,r,0,types,buffer,piece.ncells,MPI_CHAR,tag);CHKERRQ(ierr);
453552f7358SJed Brown         ierr = PetscFree3(connectivity,offsets,types);CHKERRQ(ierr);
454552f7358SJed Brown       }
455552f7358SJed Brown       {                         /* Owners (cell data) */
456552f7358SJed Brown         PetscVTKInt *owners;
457785e854fSJed Brown         ierr = PetscMalloc1(piece.ncells,&owners);CHKERRQ(ierr);
458552f7358SJed Brown         for (i=0; i<piece.ncells; i++) owners[i] = rank;
45994fbd55eSBarry Smith         ierr = TransferWrite(viewer,fp,r,0,owners,buffer,piece.ncells,MPI_INT,tag);CHKERRQ(ierr);
460552f7358SJed Brown         ierr = PetscFree(owners);CHKERRQ(ierr);
461552f7358SJed Brown       }
462552f7358SJed Brown       /* Cell data */
463552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
464552f7358SJed Brown         Vec               X = (Vec)link->vec;
465*e630c359SToby Isaac         DM                dmX;
466552f7358SJed Brown         const PetscScalar *x;
467552f7358SJed Brown         PetscScalar       *y;
468*e630c359SToby Isaac         PetscInt          bs, nfields, field;
469*e630c359SToby Isaac         PetscSection      section = NULL;
470*e630c359SToby Isaac 
471552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
472*e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
473*e630c359SToby Isaac         if (!dmX) dmX = dm;
474*e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
475*e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
476*e630c359SToby Isaac         ierr = PetscSectionGetDof(section,cStart,&bs);CHKERRQ(ierr);
477*e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
478*e630c359SToby Isaac         field = 0;
479*e630c359SToby Isaac         if (link->field >= 0) {
480*e630c359SToby Isaac           field = link->field;
481*e630c359SToby Isaac           nfields = field + 1;
482*e630c359SToby Isaac         }
483552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
484*e630c359SToby Isaac         ierr = PetscMalloc1(piece.ncells*3,&y);CHKERRQ(ierr);
485*e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
486*e630c359SToby Isaac           PetscInt     fbs,j;
487*e630c359SToby Isaac           PetscFV      fv = NULL;
488*e630c359SToby Isaac           PetscObject  f;
489*e630c359SToby Isaac           PetscClassId fClass;
490*e630c359SToby Isaac           PetscBool    vector;
491*e630c359SToby Isaac           if (nfields) {        /* We have user-defined fields/components */
492*e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,cStart,field,&fbs);CHKERRQ(ierr);
493*e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
494*e630c359SToby Isaac           ierr = DMGetField(dmX,field,NULL,&f);CHKERRQ(ierr);
495*e630c359SToby Isaac           ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr);
496*e630c359SToby Isaac           if (fClass == PETSCFV_CLASSID) {
497*e630c359SToby Isaac             fv = (PetscFV) f;
498*e630c359SToby Isaac           }
499*e630c359SToby Isaac           vector = PETSC_FALSE;
500*e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
501*e630c359SToby Isaac             vector = PETSC_TRUE;
502*e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
503*e630c359SToby Isaac               const char *compName = NULL;
504*e630c359SToby Isaac               if (fv) {
505*e630c359SToby Isaac                 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
506*e630c359SToby Isaac                 if (compName) break;
507*e630c359SToby Isaac               }
508*e630c359SToby Isaac             }
509*e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
510*e630c359SToby Isaac           }
511*e630c359SToby Isaac           if (vector) {
512552f7358SJed Brown             PetscInt cnt;
513552f7358SJed Brown             for (c=cStart,cnt=0; c<cEnd; c++) {
514*e630c359SToby Isaac               const PetscScalar *xpoint;
515*e630c359SToby Isaac               PetscInt off, j;
516*e630c359SToby Isaac 
517552f7358SJed Brown               if (hasLabel) {     /* Ignore some cells */
518552f7358SJed Brown                 PetscInt value;
519*e630c359SToby Isaac                 ierr = DMGetLabelValue(dmX, "vtk", c, &value);CHKERRQ(ierr);
520552f7358SJed Brown                 if (value != 1) continue;
521552f7358SJed Brown               }
522*e630c359SToby Isaac               if (nfields) {
523*e630c359SToby Isaac                 ierr = PetscSectionGetFieldOffset(section, c, field, &off);CHKERRQ(ierr);
524*e630c359SToby Isaac               } else {
525*e630c359SToby Isaac                 ierr = PetscSectionGetOffset(section, c, &off);CHKERRQ(ierr);
526*e630c359SToby Isaac               }
527*e630c359SToby Isaac               xpoint = &x[off];
528*e630c359SToby Isaac               for (j = 0; j < fbs; j++) {
529*e630c359SToby Isaac                 y[cnt++] = xpoint[j];
530*e630c359SToby Isaac               }
531*e630c359SToby Isaac               for (; j < 3; j++) y[cnt++] = 0.;
532*e630c359SToby Isaac             }
533*e630c359SToby Isaac             if (cnt != piece.ncells*3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
534*e630c359SToby Isaac             ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells*3,MPIU_SCALAR,tag);CHKERRQ(ierr);
535*e630c359SToby Isaac           } else {
536*e630c359SToby Isaac             for (i=0; i<fbs; i++) {
537*e630c359SToby Isaac               PetscInt cnt;
538*e630c359SToby Isaac               for (c=cStart,cnt=0; c<cEnd; c++) {
539*e630c359SToby Isaac                 const PetscScalar *xpoint;
540*e630c359SToby Isaac                 PetscInt off;
541*e630c359SToby Isaac 
542*e630c359SToby Isaac                 if (hasLabel) {     /* Ignore some cells */
543*e630c359SToby Isaac                   PetscInt value;
544*e630c359SToby Isaac                   ierr = DMGetLabelValue(dmX, "vtk", c, &value);CHKERRQ(ierr);
545*e630c359SToby Isaac                   if (value != 1) continue;
546*e630c359SToby Isaac                 }
547*e630c359SToby Isaac                 if (nfields) {
548*e630c359SToby Isaac                   ierr = PetscSectionGetFieldOffset(section, c, field, &off);CHKERRQ(ierr);
549*e630c359SToby Isaac                 } else {
550*e630c359SToby Isaac                   ierr = PetscSectionGetOffset(section, c, &off);CHKERRQ(ierr);
551*e630c359SToby Isaac                 }
552*e630c359SToby Isaac                 xpoint   = &x[off];
553552f7358SJed Brown                 y[cnt++] = xpoint[i];
554552f7358SJed Brown               }
555552f7358SJed Brown               if (cnt != piece.ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
55694fbd55eSBarry Smith               ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells,MPIU_SCALAR,tag);CHKERRQ(ierr);
557552f7358SJed Brown             }
558*e630c359SToby Isaac           }
559*e630c359SToby Isaac         }
560552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
561552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
562552f7358SJed Brown       }
563*e630c359SToby Isaac       /* point data */
564552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
565552f7358SJed Brown         Vec               X = (Vec)link->vec;
566*e630c359SToby Isaac         DM                dmX;
567552f7358SJed Brown         const PetscScalar *x;
568552f7358SJed Brown         PetscScalar       *y;
569*e630c359SToby Isaac         PetscInt          bs, nfields, field;
570*e630c359SToby Isaac         PetscSection      section = NULL;
571*e630c359SToby Isaac 
572552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
573*e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
574*e630c359SToby Isaac         if (!dmX) dmX = dm;
575*e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
576*e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
577*e630c359SToby Isaac         ierr = PetscSectionGetDof(section,vStart,&bs);CHKERRQ(ierr);
578*e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
579*e630c359SToby Isaac         field = 0;
580*e630c359SToby Isaac         if (link->field >= 0) {
581*e630c359SToby Isaac           field = link->field;
582*e630c359SToby Isaac           nfields = field + 1;
583*e630c359SToby Isaac         }
584552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
585*e630c359SToby Isaac         ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr);
586*e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
587*e630c359SToby Isaac           PetscInt   fbs,j;
588*e630c359SToby Isaac           if (nfields) {        /* We have user-defined fields/components */
589*e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,vStart,field,&fbs);CHKERRQ(ierr);
590*e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
591*e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
592552f7358SJed Brown             PetscInt cnt;
5932e529ec8SStefano Zampini             if (!localized) {
594552f7358SJed Brown               for (v=vStart,cnt=0; v<vEnd; v++) {
595*e630c359SToby Isaac                 PetscInt    off;
596*e630c359SToby Isaac                 const PetscScalar *xpoint;
597*e630c359SToby Isaac 
598*e630c359SToby Isaac                 if (nfields) {
599*e630c359SToby Isaac                   ierr     = PetscSectionGetFieldOffset(section,v,field,&off);CHKERRQ(ierr);
600*e630c359SToby Isaac                 } else {
601*e630c359SToby Isaac                   ierr     = PetscSectionGetOffset(section,v,&off);CHKERRQ(ierr);
602*e630c359SToby Isaac                 }
603*e630c359SToby Isaac                 xpoint   = &x[off];
604*e630c359SToby Isaac                 for (j = 0; j < fbs; j++) {
605*e630c359SToby Isaac                   y[cnt++] = xpoint[j];
606*e630c359SToby Isaac                 }
607*e630c359SToby Isaac                 for (; j < 3; j++) y[cnt++] = 0.;
608*e630c359SToby Isaac               }
609*e630c359SToby Isaac             } else {
610*e630c359SToby Isaac               for (c=cStart,cnt=0; c<cEnd; c++) {
611*e630c359SToby Isaac                 PetscInt *closure = NULL;
612*e630c359SToby Isaac                 PetscInt  closureSize, off;
613*e630c359SToby Isaac 
614*e630c359SToby Isaac                 ierr = DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
615*e630c359SToby Isaac                 for (v = 0, off = 0; v < closureSize*2; v += 2) {
616*e630c359SToby Isaac                   if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
617*e630c359SToby Isaac                     PetscInt    voff;
618*e630c359SToby Isaac                     const PetscScalar *xpoint;
619*e630c359SToby Isaac 
620*e630c359SToby Isaac                     if (nfields) {
621*e630c359SToby Isaac                       ierr = PetscSectionGetFieldOffset(section,closure[v],field,&voff);CHKERRQ(ierr);
622*e630c359SToby Isaac                     } else {
623*e630c359SToby Isaac                       ierr = PetscSectionGetOffset(section,closure[v],&voff);CHKERRQ(ierr);
624*e630c359SToby Isaac                     }
625*e630c359SToby Isaac                     xpoint         = &x[voff];
626*e630c359SToby Isaac                     for (j = 0; j < fbs; j++) {
627*e630c359SToby Isaac                       y[cnt + off++] = xpoint[i];
628*e630c359SToby Isaac                     }
629*e630c359SToby Isaac                     for (; j < 3; j++) y[cnt + off++] = 0.;
630*e630c359SToby Isaac                   }
631*e630c359SToby Isaac                 }
632*e630c359SToby Isaac                 cnt += off;
633*e630c359SToby Isaac                 ierr = DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
634*e630c359SToby Isaac               }
635*e630c359SToby Isaac             }
636*e630c359SToby Isaac             if (cnt != piece.nvertices*3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
637*e630c359SToby Isaac             ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices*3,MPIU_SCALAR,tag);CHKERRQ(ierr);
638*e630c359SToby Isaac           } else {
639*e630c359SToby Isaac             for (i=0; i<fbs; i++) {
640*e630c359SToby Isaac               PetscInt cnt;
641*e630c359SToby Isaac               if (!localized) {
642*e630c359SToby Isaac                 for (v=vStart,cnt=0; v<vEnd; v++) {
643*e630c359SToby Isaac                   PetscInt    off;
644*e630c359SToby Isaac                   const PetscScalar *xpoint;
645*e630c359SToby Isaac 
646*e630c359SToby Isaac                   if (nfields) {
647*e630c359SToby Isaac                     ierr     = PetscSectionGetFieldOffset(section,v,field,&off);CHKERRQ(ierr);
648*e630c359SToby Isaac                   } else {
649*e630c359SToby Isaac                     ierr     = PetscSectionGetOffset(section,v,&off);CHKERRQ(ierr);
650*e630c359SToby Isaac                   }
651*e630c359SToby Isaac                   xpoint   = &x[off];
652552f7358SJed Brown                   y[cnt++] = xpoint[i];
653552f7358SJed Brown                 }
6542e529ec8SStefano Zampini               } else {
6552e529ec8SStefano Zampini                 for (c=cStart,cnt=0; c<cEnd; c++) {
6562e529ec8SStefano Zampini                   PetscInt *closure = NULL;
6572e529ec8SStefano Zampini                   PetscInt  closureSize, off;
6582e529ec8SStefano Zampini 
659*e630c359SToby Isaac                   ierr = DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
6602e529ec8SStefano Zampini                   for (v = 0, off = 0; v < closureSize*2; v += 2) {
6612e529ec8SStefano Zampini                     if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
662*e630c359SToby Isaac                       PetscInt    voff;
663*e630c359SToby Isaac                       const PetscScalar *xpoint;
6642e529ec8SStefano Zampini 
665*e630c359SToby Isaac                       if (nfields) {
666*e630c359SToby Isaac                         ierr           = PetscSectionGetFieldOffset(section,closure[v],field,&voff);CHKERRQ(ierr);
667*e630c359SToby Isaac                       } else {
668*e630c359SToby Isaac                         ierr           = PetscSectionGetOffset(section,closure[v],&voff);CHKERRQ(ierr);
669*e630c359SToby Isaac                       }
670*e630c359SToby Isaac                       xpoint         = &x[voff];
6712e529ec8SStefano Zampini                       y[cnt + off++] = xpoint[i];
6722e529ec8SStefano Zampini                     }
6732e529ec8SStefano Zampini                   }
6749bda96f6SStefano Zampini                   cnt += off;
675*e630c359SToby Isaac                   ierr = DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
6762e529ec8SStefano Zampini                 }
6772e529ec8SStefano Zampini               }
678552f7358SJed Brown               if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
67994fbd55eSBarry Smith               ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices,MPIU_SCALAR,tag);CHKERRQ(ierr);
680552f7358SJed Brown             }
681*e630c359SToby Isaac           }
682*e630c359SToby Isaac         }
683552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
684552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
685552f7358SJed Brown       }
686552f7358SJed Brown     } else if (!rank) {
68794fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,MPIU_SCALAR,tag);CHKERRQ(ierr); /* positions */
68894fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nconn,MPI_INT,tag);CHKERRQ(ierr); /* connectivity */
68994fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag);CHKERRQ(ierr); /* offsets */
69094fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_CHAR,tag);CHKERRQ(ierr); /* types */
69194fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag);CHKERRQ(ierr); /* owner rank (cells) */
692552f7358SJed Brown       /* all cell data */
693552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
694*e630c359SToby Isaac         Vec               X = (Vec)link->vec;
695*e630c359SToby Isaac         PetscInt bs, nfields, field;
696*e630c359SToby Isaac         DM           dmX;
697*e630c359SToby Isaac         PetscSection section = NULL;
698*e630c359SToby Isaac 
699552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
700*e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
701*e630c359SToby Isaac         if (!dmX) dmX = dm;
702*e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
703*e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
704*e630c359SToby Isaac         ierr = PetscSectionGetDof(section,cStart,&bs);CHKERRQ(ierr);
705*e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
706*e630c359SToby Isaac         field = 0;
707*e630c359SToby Isaac         if (link->field >= 0) {
708*e630c359SToby Isaac           field = link->field;
709*e630c359SToby Isaac           nfields = field + 1;
710*e630c359SToby Isaac         }
711*e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
712*e630c359SToby Isaac           PetscInt     fbs,j;
713*e630c359SToby Isaac           PetscFV      fv = NULL;
714*e630c359SToby Isaac           PetscObject  f;
715*e630c359SToby Isaac           PetscClassId fClass;
716*e630c359SToby Isaac           PetscBool    vector;
717*e630c359SToby Isaac           if (nfields) {        /* We have user-defined fields/components */
718*e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,cStart,field,&fbs);CHKERRQ(ierr);
719*e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
720*e630c359SToby Isaac           ierr = DMGetField(dmX,field,NULL,&f);CHKERRQ(ierr);
721*e630c359SToby Isaac           ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr);
722*e630c359SToby Isaac           if (fClass == PETSCFV_CLASSID) {
723*e630c359SToby Isaac             fv = (PetscFV) f;
724*e630c359SToby Isaac           }
725*e630c359SToby Isaac           vector = PETSC_FALSE;
726*e630c359SToby Isaac           if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) {
727*e630c359SToby Isaac             vector = PETSC_TRUE;
728*e630c359SToby Isaac             for (j = 0; j < fbs; j++) {
729*e630c359SToby Isaac               const char *compName = NULL;
730*e630c359SToby Isaac               if (fv) {
731*e630c359SToby Isaac                 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
732*e630c359SToby Isaac                 if (compName) break;
733*e630c359SToby Isaac               }
734*e630c359SToby Isaac             }
735*e630c359SToby Isaac             if (j < fbs) vector = PETSC_FALSE;
736*e630c359SToby Isaac           }
737*e630c359SToby Isaac           if (vector) {
738*e630c359SToby Isaac             ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells*3,MPIU_SCALAR,tag);CHKERRQ(ierr);
739*e630c359SToby Isaac           } else {
740*e630c359SToby Isaac             for (i=0; i<fbs; i++) {
74194fbd55eSBarry Smith               ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPIU_SCALAR,tag);CHKERRQ(ierr);
742552f7358SJed Brown             }
743552f7358SJed Brown           }
744*e630c359SToby Isaac         }
745*e630c359SToby Isaac       }
746552f7358SJed Brown       /* all point data */
747552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
748*e630c359SToby Isaac         Vec               X = (Vec)link->vec;
749*e630c359SToby Isaac         DM                dmX;
750*e630c359SToby Isaac         PetscInt bs, nfields, field;
751*e630c359SToby Isaac         PetscSection section = NULL;
752*e630c359SToby Isaac 
753552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
754*e630c359SToby Isaac         ierr = VecGetDM(X, &dmX);CHKERRQ(ierr);
755*e630c359SToby Isaac         if (!dmX) dmX = dm;
756*e630c359SToby Isaac         ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &section);CHKERRQ(ierr);
757*e630c359SToby Isaac         if (!section) {ierr = DMGetLocalSection(dmX, &section);CHKERRQ(ierr);}
758*e630c359SToby Isaac         ierr = PetscSectionGetDof(section,vStart,&bs);CHKERRQ(ierr);
759*e630c359SToby Isaac         ierr = PetscSectionGetNumFields(section,&nfields);CHKERRQ(ierr);
760*e630c359SToby Isaac         field = 0;
761*e630c359SToby Isaac         if (link->field >= 0) {
762*e630c359SToby Isaac           field = link->field;
763*e630c359SToby Isaac           nfields = field + 1;
764*e630c359SToby Isaac         }
765*e630c359SToby Isaac         for (i=0; field<(nfields?nfields:1); field++) {
766*e630c359SToby Isaac           PetscInt   fbs;
767*e630c359SToby Isaac           if (nfields) {        /* We have user-defined fields/components */
768*e630c359SToby Isaac             ierr = PetscSectionGetFieldDof(section,vStart,field,&fbs);CHKERRQ(ierr);
769*e630c359SToby Isaac           } else fbs = bs;      /* Say we have one field with 'bs' components */
770*e630c359SToby Isaac           if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) {
771*e630c359SToby Isaac             ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,MPIU_SCALAR,tag);CHKERRQ(ierr);
772*e630c359SToby Isaac           } else {
773*e630c359SToby Isaac             for (i=0; i<fbs; i++) {
77494fbd55eSBarry Smith               ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices,MPIU_SCALAR,tag);CHKERRQ(ierr);
775552f7358SJed Brown             }
776552f7358SJed Brown           }
777552f7358SJed Brown         }
778552f7358SJed Brown       }
779*e630c359SToby Isaac     }
780*e630c359SToby Isaac   }
781552f7358SJed Brown   ierr = PetscFree(gpiece);CHKERRQ(ierr);
782552f7358SJed Brown   ierr = PetscFree(buffer);CHKERRQ(ierr);
783552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"\n  </AppendedData>\n");CHKERRQ(ierr);
784552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
7855e32de72SMatthew G. Knepley   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
786552f7358SJed Brown   PetscFunctionReturn(0);
787552f7358SJed Brown }
788