xref: /petsc/src/dm/impls/plex/plexvtu.c (revision 2f029394abf406ed8c3d96e556ffca625f5fd18c)
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;
53*2f029394SStefano 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;
128*2f029394SStefano 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;
132552f7358SJed Brown 
133552f7358SJed Brown   PetscFunctionBegin;
134ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
135552f7358SJed Brown #if defined(PETSC_USE_COMPLEX)
136ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Complex values not supported");
137552f7358SJed Brown #endif
138552f7358SJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
139552f7358SJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
140552f7358SJed Brown   tag  = ((PetscObject)viewer)->tag;
141552f7358SJed Brown 
142552f7358SJed Brown   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
143552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
144519f805aSKarl Rupp #if defined(PETSC_WORDS_BIGENDIAN)
145552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr);
146552f7358SJed Brown #else
147552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr);
148552f7358SJed Brown #endif
149552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <UnstructuredGrid>\n");CHKERRQ(ierr);
150552f7358SJed Brown 
1513baa072aSToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
152552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
153552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
154552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
155c58f1c22SToby Isaac   ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
1562e529ec8SStefano Zampini   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
1572e529ec8SStefano Zampini   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1588865f1eaSKarl Rupp 
159552f7358SJed Brown   hasLabel        = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
1602e529ec8SStefano Zampini   piece.nvertices = 0;
161552f7358SJed Brown   piece.ncells    = 0;
162552f7358SJed Brown   piece.nconn     = 0;
1632e529ec8SStefano Zampini   if (!localized) piece.nvertices = vEnd - vStart;
164552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
16519003fb5SStefano Zampini     PetscInt dof = 0;
166552f7358SJed Brown     if (hasLabel) {
167552f7358SJed Brown       PetscInt value;
168552f7358SJed Brown 
169c58f1c22SToby Isaac       ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
170552f7358SJed Brown       if (value != 1) continue;
171552f7358SJed Brown     }
17219003fb5SStefano Zampini     if (localized) {
17319003fb5SStefano Zampini       ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
17419003fb5SStefano Zampini     }
17519003fb5SStefano Zampini     if (!dof) {
1762e529ec8SStefano Zampini       PetscInt *closure = NULL;
1772e529ec8SStefano Zampini       PetscInt closureSize;
1782e529ec8SStefano Zampini 
179552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
180552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
1812e529ec8SStefano Zampini         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
1822e529ec8SStefano Zampini           piece.nconn++;
18319003fb5SStefano Zampini           if (localized) piece.nvertices++;
1842e529ec8SStefano Zampini         }
185552f7358SJed Brown       }
186552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1872e529ec8SStefano Zampini     } else {
1882e529ec8SStefano Zampini       piece.nvertices += dof/dimEmbed;
1892e529ec8SStefano Zampini       piece.nconn     += dof/dimEmbed;
1902e529ec8SStefano Zampini     }
191552f7358SJed Brown     piece.ncells++;
192552f7358SJed Brown   }
193785e854fSJed Brown   if (!rank) {ierr = PetscMalloc1(size,&gpiece);CHKERRQ(ierr);}
194552f7358SJed Brown   ierr = MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm);CHKERRQ(ierr);
195552f7358SJed Brown 
196552f7358SJed Brown   /*
197552f7358SJed Brown    * Write file header
198552f7358SJed Brown    */
199552f7358SJed Brown   if (!rank) {
200552f7358SJed Brown     PetscInt boffset = 0;
201552f7358SJed Brown 
202552f7358SJed Brown     for (r=0; r<size; r++) {
203552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells);CHKERRQ(ierr);
204552f7358SJed Brown       /* Coordinate positions */
205552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Points>\n");CHKERRQ(ierr);
206552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
207552f7358SJed Brown       boffset += gpiece[r].nvertices*3*sizeof(PetscScalar) + sizeof(int);
208552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Points>\n");CHKERRQ(ierr);
209552f7358SJed Brown       /* Cell connectivity */
210552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Cells>\n");CHKERRQ(ierr);
211552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
212552f7358SJed Brown       boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int);
213552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"offsets\"      NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
214552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int);
215552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"UInt8\" Name=\"types\"        NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
216552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int);
217552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Cells>\n");CHKERRQ(ierr);
218552f7358SJed Brown 
219552f7358SJed Brown       /*
220552f7358SJed Brown        * Cell Data headers
221552f7358SJed Brown        */
222552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <CellData>\n");CHKERRQ(ierr);
223552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
224552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(int) + sizeof(int);
225552f7358SJed Brown       /* all the vectors */
226552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
227552f7358SJed Brown         Vec        X = (Vec)link->vec;
2281cfafdd3SJed Brown         PetscInt   bs,nfields,field;
229552f7358SJed Brown         const char *vecname = "";
230552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
231552f7358SJed 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. */
232552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
233552f7358SJed Brown         }
234552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
2351cfafdd3SJed Brown         ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr);
2367ded4cbaSJed Brown         for (field=0,i=0; field<(nfields?nfields:1); field++) {
2371cfafdd3SJed Brown           PetscInt     fbs,j;
238a00cdb45SToby Isaac           PetscFV      fv = NULL;
239a00cdb45SToby Isaac           PetscObject  f;
240a00cdb45SToby Isaac           PetscClassId fClass;
2410298fd71SBarry Smith           const char *fieldname = NULL;
2421cfafdd3SJed Brown           char       buf[256];
2437ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
2441cfafdd3SJed Brown             ierr = PetscSectionGetFieldDof(dm->defaultSection,cStart,field,&fbs);CHKERRQ(ierr);
2451cfafdd3SJed Brown             ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr);
2467ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
247a00cdb45SToby Isaac           ierr = DMGetField(dm,field,&f);CHKERRQ(ierr);
248a00cdb45SToby Isaac           ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr);
249a00cdb45SToby Isaac           if (fClass == PETSCFV_CLASSID) {
250a00cdb45SToby Isaac             fv = (PetscFV) f;
251a00cdb45SToby Isaac           }
252552f7358SJed Brown           if (!fieldname) {
2531cfafdd3SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"CellField%D",field);CHKERRQ(ierr);
254552f7358SJed Brown             fieldname = buf;
255552f7358SJed Brown           }
2561cfafdd3SJed Brown           for (j=0; j<fbs; j++) {
257a00cdb45SToby Isaac             const char *compName = NULL;
258a00cdb45SToby Isaac             if (fv) {
259a00cdb45SToby Isaac               ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
260a00cdb45SToby Isaac             }
261a00cdb45SToby Isaac             if (compName) {
262a00cdb45SToby 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);
263a00cdb45SToby Isaac             }
264a00cdb45SToby Isaac             else {
2651cfafdd3SJed 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);
266a00cdb45SToby Isaac             }
267552f7358SJed Brown             boffset += gpiece[r].ncells*sizeof(PetscScalar) + sizeof(int);
2681cfafdd3SJed Brown             i++;
269552f7358SJed Brown           }
270552f7358SJed Brown         }
2711cfafdd3SJed Brown         if (i != bs) SETERRQ2(comm,PETSC_ERR_PLIB,"Total number of field components %D != block size %D",i,bs);
2721cfafdd3SJed Brown       }
273552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </CellData>\n");CHKERRQ(ierr);
274552f7358SJed Brown 
275552f7358SJed Brown       /*
276552f7358SJed Brown        * Point Data headers
277552f7358SJed Brown        */
278552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      <PointData>\n");CHKERRQ(ierr);
279552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
280552f7358SJed Brown         Vec        X = (Vec)link->vec;
2811cfafdd3SJed Brown         PetscInt   bs,nfields,field;
282552f7358SJed Brown         const char *vecname = "";
283552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
284552f7358SJed 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. */
285552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
286552f7358SJed Brown         }
287552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
2881cfafdd3SJed Brown         ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr);
2897ded4cbaSJed Brown         for (field=0,i=0; field<(nfields?nfields:1); field++) {
2901cfafdd3SJed Brown           PetscInt   fbs,j;
2911cfafdd3SJed Brown           const char *fieldname = NULL;
2921cfafdd3SJed Brown           char       buf[256];
2937ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
2941cfafdd3SJed Brown             ierr = PetscSectionGetFieldDof(dm->defaultSection,vStart,field,&fbs);CHKERRQ(ierr);
2951cfafdd3SJed Brown             ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr);
2967ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
2971cfafdd3SJed Brown           if (!fieldname) {
2981cfafdd3SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"PointField%D",field);CHKERRQ(ierr);
2991cfafdd3SJed Brown             fieldname = buf;
3001cfafdd3SJed Brown           }
3011cfafdd3SJed Brown           for (j=0; j<fbs; j++) {
3021cfafdd3SJed 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);
303552f7358SJed Brown             boffset += gpiece[r].nvertices*sizeof(PetscScalar) + sizeof(int);
304552f7358SJed Brown           }
305552f7358SJed Brown         }
3061cfafdd3SJed Brown       }
307552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </PointData>\n");CHKERRQ(ierr);
308552f7358SJed Brown 
309552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    </Piece>\n");CHKERRQ(ierr);
310552f7358SJed Brown     }
311552f7358SJed Brown   }
312552f7358SJed Brown 
313552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  </UnstructuredGrid>\n");CHKERRQ(ierr);
314552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
315552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
316552f7358SJed Brown 
317552f7358SJed Brown   if (!rank) {
318552f7358SJed Brown     PetscInt maxsize = 0;
319552f7358SJed Brown     for (r=0; r<size; r++) {
320552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscScalar)));
321552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*sizeof(PetscScalar)));
322552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt)));
323552f7358SJed Brown     }
324552f7358SJed Brown     ierr = PetscMalloc(maxsize,&buffer);CHKERRQ(ierr);
325552f7358SJed Brown   }
326552f7358SJed Brown   for (r=0; r<size; r++) {
327552f7358SJed Brown     if (r == rank) {
328552f7358SJed Brown       PetscInt nsend;
329552f7358SJed Brown       {                         /* Position */
330552f7358SJed Brown         const PetscScalar *x;
3310298fd71SBarry Smith         PetscScalar       *y = NULL;
332552f7358SJed Brown         Vec               coords;
3332e529ec8SStefano Zampini 
334552f7358SJed Brown         ierr  = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
335552f7358SJed Brown         ierr  = VecGetArrayRead(coords,&x);CHKERRQ(ierr);
33619003fb5SStefano Zampini         if (dimEmbed != 3 || localized) {
337785e854fSJed Brown           ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr);
33819003fb5SStefano Zampini           if (localized) {
33919003fb5SStefano Zampini             PetscInt cnt;
34019003fb5SStefano Zampini             for (c=cStart,cnt=0; c<cEnd; c++) {
34119003fb5SStefano Zampini               PetscInt off, dof;
34219003fb5SStefano Zampini 
34319003fb5SStefano Zampini               ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
34419003fb5SStefano Zampini               if (!dof) {
34519003fb5SStefano Zampini                 PetscInt *closure = NULL;
34619003fb5SStefano Zampini                 PetscInt closureSize;
34719003fb5SStefano Zampini 
34819003fb5SStefano Zampini                 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
34919003fb5SStefano Zampini                 for (v = 0; v < closureSize*2; v += 2) {
35019003fb5SStefano Zampini                   if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
35119003fb5SStefano Zampini                     ierr = PetscSectionGetOffset(coordSection, closure[v], &off);CHKERRQ(ierr);
35219003fb5SStefano Zampini                     if (dimEmbed != 3) {
35319003fb5SStefano Zampini                       y[cnt*3+0] = x[off+0];
35419003fb5SStefano Zampini                       y[cnt*3+1] = (dimEmbed > 1) ? x[off+1] : 0.0;
35519003fb5SStefano Zampini                       y[cnt*3+2] = 0.0;
35619003fb5SStefano Zampini                     } else {
35719003fb5SStefano Zampini                       y[cnt*3+0] = x[off+0];
35819003fb5SStefano Zampini                       y[cnt*3+1] = x[off+1];
35919003fb5SStefano Zampini                       y[cnt*3+2] = x[off+2];
36019003fb5SStefano Zampini                     }
36119003fb5SStefano Zampini                     cnt++;
36219003fb5SStefano Zampini                   }
36319003fb5SStefano Zampini                 }
36419003fb5SStefano Zampini                 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
36519003fb5SStefano Zampini               } else {
36619003fb5SStefano Zampini                 ierr = PetscSectionGetOffset(coordSection, c, &off);CHKERRQ(ierr);
36719003fb5SStefano Zampini                 if (dimEmbed != 3) {
36819003fb5SStefano Zampini                   for (i=0; i<dof/dimEmbed; i++) {
36919003fb5SStefano Zampini                     y[cnt*3+0] = x[off + i*dimEmbed + 0];
37019003fb5SStefano Zampini                     y[cnt*3+1] = (dimEmbed > 1) ? x[off + i*dimEmbed + 1] : 0.0;
37119003fb5SStefano Zampini                     y[cnt*3+2] = 0.0;
37219003fb5SStefano Zampini                     cnt++;
37319003fb5SStefano Zampini                   }
37419003fb5SStefano Zampini                 } else {
37519003fb5SStefano Zampini                   for (i=0; i<dof; i ++) {
37619003fb5SStefano Zampini 		    y[cnt*3+i] = x[off + i];
37719003fb5SStefano Zampini                   }
37819003fb5SStefano Zampini                   cnt += dof/dimEmbed;
37919003fb5SStefano Zampini                 }
38019003fb5SStefano Zampini               }
38119003fb5SStefano Zampini             }
38219003fb5SStefano Zampini             if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
38319003fb5SStefano Zampini           } else {
384552f7358SJed Brown             for (i=0; i<piece.nvertices; i++) {
3853baa072aSToby Isaac               y[i*3+0] = x[i*dimEmbed+0];
3863baa072aSToby Isaac               y[i*3+1] = (dimEmbed > 1) ? x[i*dimEmbed+1] : 0;
38719003fb5SStefano Zampini               y[i*3+2] = 0.0;
38819003fb5SStefano Zampini             }
389552f7358SJed Brown           }
390552f7358SJed Brown         }
3912e529ec8SStefano Zampini         nsend = piece.nvertices*3;
39294fbd55eSBarry Smith         ierr  = TransferWrite(viewer,fp,r,0,y ? y : x,buffer,nsend,MPIU_SCALAR,tag);CHKERRQ(ierr);
393552f7358SJed Brown         ierr  = PetscFree(y);CHKERRQ(ierr);
394552f7358SJed Brown         ierr  = VecRestoreArrayRead(coords,&x);CHKERRQ(ierr);
395552f7358SJed Brown       }
396552f7358SJed Brown       {                           /* Connectivity, offsets, types */
3973e869605SMatthew G. Knepley         PetscVTKInt  *connectivity = NULL, *offsets = NULL;
3983e869605SMatthew G. Knepley         PetscVTKType *types = NULL;
3992e529ec8SStefano Zampini         ierr = DMPlexGetVTKConnectivity(dm,localized,&piece,&connectivity,&offsets,&types);CHKERRQ(ierr);
40094fbd55eSBarry Smith         ierr = TransferWrite(viewer,fp,r,0,connectivity,buffer,piece.nconn,MPI_INT,tag);CHKERRQ(ierr);
40194fbd55eSBarry Smith         ierr = TransferWrite(viewer,fp,r,0,offsets,buffer,piece.ncells,MPI_INT,tag);CHKERRQ(ierr);
40294fbd55eSBarry Smith         ierr = TransferWrite(viewer,fp,r,0,types,buffer,piece.ncells,MPI_CHAR,tag);CHKERRQ(ierr);
403552f7358SJed Brown         ierr = PetscFree3(connectivity,offsets,types);CHKERRQ(ierr);
404552f7358SJed Brown       }
405552f7358SJed Brown       {                         /* Owners (cell data) */
406552f7358SJed Brown         PetscVTKInt *owners;
407785e854fSJed Brown         ierr = PetscMalloc1(piece.ncells,&owners);CHKERRQ(ierr);
408552f7358SJed Brown         for (i=0; i<piece.ncells; i++) owners[i] = rank;
40994fbd55eSBarry Smith         ierr = TransferWrite(viewer,fp,r,0,owners,buffer,piece.ncells,MPI_INT,tag);CHKERRQ(ierr);
410552f7358SJed Brown         ierr = PetscFree(owners);CHKERRQ(ierr);
411552f7358SJed Brown       }
412552f7358SJed Brown       /* Cell data */
413552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
414552f7358SJed Brown         Vec               X = (Vec)link->vec;
415552f7358SJed Brown         const PetscScalar *x;
416552f7358SJed Brown         PetscScalar       *y;
417552f7358SJed Brown         PetscInt          bs;
418552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
419552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
420552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
421785e854fSJed Brown         ierr = PetscMalloc1(piece.ncells,&y);CHKERRQ(ierr);
422552f7358SJed Brown         for (i=0; i<bs; i++) {
423552f7358SJed Brown           PetscInt cnt;
424552f7358SJed Brown           for (c=cStart,cnt=0; c<cEnd; c++) {
425640bce14SSatish Balay             PetscScalar *xpoint;
426552f7358SJed Brown             if (hasLabel) {     /* Ignore some cells */
427552f7358SJed Brown               PetscInt value;
428c58f1c22SToby Isaac               ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
429552f7358SJed Brown               if (value != 1) continue;
430552f7358SJed Brown             }
431552f7358SJed Brown             ierr     = DMPlexPointLocalRead(dm,c,x,&xpoint);CHKERRQ(ierr);
432552f7358SJed Brown             y[cnt++] = xpoint[i];
433552f7358SJed Brown           }
434552f7358SJed Brown           if (cnt != piece.ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
43594fbd55eSBarry Smith           ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells,MPIU_SCALAR,tag);CHKERRQ(ierr);
436552f7358SJed Brown         }
437552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
438552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
439552f7358SJed Brown       }
440552f7358SJed Brown 
441552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
442552f7358SJed Brown         Vec               X = (Vec)link->vec;
443552f7358SJed Brown         const PetscScalar *x;
444552f7358SJed Brown         PetscScalar       *y;
445552f7358SJed Brown         PetscInt          bs;
446552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
447552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
448552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
449785e854fSJed Brown         ierr = PetscMalloc1(piece.nvertices,&y);CHKERRQ(ierr);
450552f7358SJed Brown         for (i=0; i<bs; i++) {
451552f7358SJed Brown           PetscInt cnt;
4522e529ec8SStefano Zampini           if (!localized) {
453552f7358SJed Brown             for (v=vStart,cnt=0; v<vEnd; v++) {
454640bce14SSatish Balay               PetscScalar *xpoint;
45537045cddSJed Brown               ierr     = DMPlexPointLocalRead(dm,v,x,&xpoint);CHKERRQ(ierr);
456552f7358SJed Brown               y[cnt++] = xpoint[i];
457552f7358SJed Brown             }
4582e529ec8SStefano Zampini           } else {
4592e529ec8SStefano Zampini             for (c=cStart,cnt=0; c<cEnd; c++) {
4602e529ec8SStefano Zampini               PetscInt *closure = NULL;
4612e529ec8SStefano Zampini               PetscInt  closureSize, off;
4622e529ec8SStefano Zampini 
4632e529ec8SStefano Zampini               ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
4642e529ec8SStefano Zampini               for (v = 0, off = 0; v < closureSize*2; v += 2) {
4652e529ec8SStefano Zampini                 if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
4662e529ec8SStefano Zampini                   PetscScalar *xpoint;
4672e529ec8SStefano Zampini 
4682e529ec8SStefano Zampini                   ierr = DMPlexPointLocalRead(dm,v,x,&xpoint);CHKERRQ(ierr);
4692e529ec8SStefano Zampini                   y[cnt + off++] = xpoint[i];
4702e529ec8SStefano Zampini                 }
4712e529ec8SStefano Zampini               }
4722e529ec8SStefano Zampini               ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
4732e529ec8SStefano Zampini             }
4742e529ec8SStefano Zampini           }
475552f7358SJed Brown           if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
47694fbd55eSBarry Smith           ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices,MPIU_SCALAR,tag);CHKERRQ(ierr);
477552f7358SJed Brown         }
478552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
479552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
480552f7358SJed Brown       }
481552f7358SJed Brown     } else if (!rank) {
48294fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,MPIU_SCALAR,tag);CHKERRQ(ierr); /* positions */
48394fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nconn,MPI_INT,tag);CHKERRQ(ierr); /* connectivity */
48494fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag);CHKERRQ(ierr); /* offsets */
48594fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_CHAR,tag);CHKERRQ(ierr); /* types */
48694fbd55eSBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag);CHKERRQ(ierr); /* owner rank (cells) */
487552f7358SJed Brown       /* all cell data */
488552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
489552f7358SJed Brown         PetscInt bs;
490552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
491fd68c46aSJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
492552f7358SJed Brown         for (i=0; i<bs; i++) {
49394fbd55eSBarry Smith           ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPIU_SCALAR,tag);CHKERRQ(ierr);
494552f7358SJed Brown         }
495552f7358SJed Brown       }
496552f7358SJed Brown       /* all point data */
497552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
498552f7358SJed Brown         PetscInt bs;
499552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
500552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
501552f7358SJed Brown         for (i=0; i<bs; i++) {
50294fbd55eSBarry Smith           ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices,MPIU_SCALAR,tag);CHKERRQ(ierr);
503552f7358SJed Brown         }
504552f7358SJed Brown       }
505552f7358SJed Brown     }
506552f7358SJed Brown   }
507552f7358SJed Brown   ierr = PetscFree(gpiece);CHKERRQ(ierr);
508552f7358SJed Brown   ierr = PetscFree(buffer);CHKERRQ(ierr);
509552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"\n  </AppendedData>\n");CHKERRQ(ierr);
510552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
5115e32de72SMatthew G. Knepley   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
512552f7358SJed Brown   PetscFunctionReturn(0);
513552f7358SJed Brown }
514