xref: /petsc/src/dm/impls/plex/plexvtu.c (revision dcca6d9d80ebd869fe6029bd05a3aa9faafef49e)
134541f0dSBarry 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 
18552f7358SJed Brown #undef __FUNCT__
19552f7358SJed Brown #define __FUNCT__ "TransferWrite"
20552f7358SJed Brown static PetscErrorCode TransferWrite(PetscViewer viewer,FILE *fp,PetscMPIInt srank,PetscMPIInt root,const void *send,void *recv,PetscMPIInt count,PetscDataType datatype,PetscMPIInt tag)
21552f7358SJed Brown {
22552f7358SJed Brown   PetscMPIInt    rank;
23552f7358SJed Brown   PetscErrorCode ierr;
24ce94432eSBarry Smith   MPI_Comm       comm;
25552f7358SJed Brown   MPI_Datatype   mpidatatype;
26552f7358SJed Brown 
27552f7358SJed Brown   PetscFunctionBegin;
28ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
29552f7358SJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
30552f7358SJed Brown   ierr = PetscDataTypeToMPIDataType(datatype,&mpidatatype);CHKERRQ(ierr);
31552f7358SJed Brown 
32552f7358SJed Brown   if (rank == srank && rank != root) {
33552f7358SJed Brown     ierr = MPI_Send((void*)send,count,mpidatatype,root,tag,comm);CHKERRQ(ierr);
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;
41552f7358SJed Brown       ierr = MPI_Recv(recv,count,mpidatatype,srank,tag,comm,&status);CHKERRQ(ierr);
42552f7358SJed Brown       ierr = MPI_Get_count(&status,mpidatatype,&nrecv);CHKERRQ(ierr);
43552f7358SJed Brown       if (count != nrecv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
44552f7358SJed Brown       buffer = recv;
45552f7358SJed Brown     }
46552f7358SJed Brown     ierr = PetscViewerVTKFWrite(viewer,fp,buffer,count,datatype);CHKERRQ(ierr);
47552f7358SJed Brown   }
48552f7358SJed Brown   PetscFunctionReturn(0);
49552f7358SJed Brown }
50552f7358SJed Brown 
51552f7358SJed Brown #undef __FUNCT__
52552f7358SJed Brown #define __FUNCT__ "DMPlexGetVTKConnectivity"
53552f7358SJed Brown static PetscErrorCode DMPlexGetVTKConnectivity(DM dm,PieceInfo *piece,PetscVTKInt **oconn,PetscVTKInt **ooffsets,PetscVTKType **otypes)
54552f7358SJed Brown {
55552f7358SJed Brown   PetscErrorCode ierr;
56552f7358SJed Brown   PetscVTKInt    *conn,*offsets;
57552f7358SJed Brown   PetscVTKType   *types;
58552f7358SJed Brown   PetscInt       dim,vStart,vEnd,cStart,cEnd,pStart,pEnd,cellHeight,cMax,numLabelCells,hasLabel,c,v,countcell,countconn;
59552f7358SJed Brown 
60552f7358SJed Brown   PetscFunctionBegin;
61*dcca6d9dSJed Brown   ierr = PetscMalloc3(piece->nconn,&conn,piece->ncells,&offsets,piece->ncells,&types);CHKERRQ(ierr);
62552f7358SJed Brown 
63552f7358SJed Brown   ierr = DMPlexGetDimension(dm,&dim);CHKERRQ(ierr);
64552f7358SJed Brown   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
65552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
66552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
67552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
680298fd71SBarry Smith   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
698865f1eaSKarl Rupp   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
70552f7358SJed Brown   ierr     = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
71552f7358SJed Brown   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
72552f7358SJed Brown 
73552f7358SJed Brown   countcell = 0;
74552f7358SJed Brown   countconn = 0;
75552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
760298fd71SBarry Smith     PetscInt *closure = NULL;
77724b5320SMatthew G. Knepley     PetscInt  closureSize,nverts,celltype,startoffset,nC=0;
78552f7358SJed Brown 
79552f7358SJed Brown     if (hasLabel) {
80552f7358SJed Brown       PetscInt value;
81552f7358SJed Brown 
82552f7358SJed Brown       ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
83552f7358SJed Brown       if (value != 1) continue;
84552f7358SJed Brown     }
85552f7358SJed Brown     startoffset = countconn;
86552f7358SJed Brown     ierr        = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
87552f7358SJed Brown     for (v = 0; v < closureSize*2; v += 2) {
88552f7358SJed Brown       if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
89552f7358SJed Brown         conn[countconn++] = closure[v] - vStart;
90724b5320SMatthew G. Knepley         ++nC;
91552f7358SJed Brown       }
92552f7358SJed Brown     }
93552f7358SJed Brown     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
94724b5320SMatthew G. Knepley     ierr = DMPlexInvertCell(dim, nC, &conn[countconn-nC]);CHKERRQ(ierr);
958865f1eaSKarl Rupp 
96552f7358SJed Brown     offsets[countcell] = countconn;
978865f1eaSKarl Rupp 
98552f7358SJed Brown     nverts = countconn - startoffset;
99552f7358SJed Brown     ierr   = DMPlexVTKGetCellType(dm,dim,nverts,&celltype);CHKERRQ(ierr);
1008865f1eaSKarl Rupp 
101552f7358SJed Brown     types[countcell] = celltype;
102552f7358SJed Brown     countcell++;
103552f7358SJed Brown   }
104552f7358SJed Brown   if (countcell != piece->ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count");
105552f7358SJed Brown   if (countconn != piece->nconn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count");
106552f7358SJed Brown   *oconn    = conn;
107552f7358SJed Brown   *ooffsets = offsets;
108552f7358SJed Brown   *otypes   = types;
109552f7358SJed Brown   PetscFunctionReturn(0);
110552f7358SJed Brown }
111552f7358SJed Brown 
112552f7358SJed Brown #undef __FUNCT__
113552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll_VTU"
114552f7358SJed Brown /*
115552f7358SJed Brown   Write all fields that have been provided to the viewer
116552f7358SJed Brown   Multi-block XML format with binary appended data.
117552f7358SJed Brown */
118552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer)
119552f7358SJed Brown {
120ce94432eSBarry Smith   MPI_Comm                 comm;
121552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
122552f7358SJed Brown   PetscViewerVTKObjectLink link;
123552f7358SJed Brown   FILE                     *fp;
124552f7358SJed Brown   PetscMPIInt              rank,size,tag;
125552f7358SJed Brown   PetscErrorCode           ierr;
126552f7358SJed Brown   PetscInt                 dim,cellHeight,cStart,cEnd,vStart,vEnd,cMax,numLabelCells,hasLabel,c,v,r,i;
1270298fd71SBarry Smith   PieceInfo                piece,*gpiece = NULL;
1280298fd71SBarry Smith   void                     *buffer = NULL;
129552f7358SJed Brown 
130552f7358SJed Brown   PetscFunctionBegin;
131ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
132552f7358SJed Brown #if defined(PETSC_USE_COMPLEX)
133ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Complex values not supported");
134552f7358SJed Brown #endif
135552f7358SJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
136552f7358SJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
137552f7358SJed Brown   tag  = ((PetscObject)viewer)->tag;
138552f7358SJed Brown 
139552f7358SJed Brown   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
140552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
141519f805aSKarl Rupp #if defined(PETSC_WORDS_BIGENDIAN)
142552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr);
143552f7358SJed Brown #else
144552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr);
145552f7358SJed Brown #endif
146552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <UnstructuredGrid>\n");CHKERRQ(ierr);
147552f7358SJed Brown 
148552f7358SJed Brown   ierr = DMPlexGetDimension(dm, &dim);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);
1520298fd71SBarry Smith   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
1538865f1eaSKarl Rupp   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
154552f7358SJed Brown   ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
1558865f1eaSKarl Rupp 
156552f7358SJed Brown   hasLabel        = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
157552f7358SJed Brown   piece.nvertices = vEnd - vStart;
158552f7358SJed Brown   piece.ncells    = 0;
159552f7358SJed Brown   piece.nconn     = 0;
160552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
1610298fd71SBarry Smith     PetscInt *closure = NULL;
162552f7358SJed Brown     PetscInt closureSize;
163552f7358SJed Brown 
164552f7358SJed Brown     if (hasLabel) {
165552f7358SJed Brown       PetscInt value;
166552f7358SJed Brown 
167552f7358SJed Brown       ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
168552f7358SJed Brown       if (value != 1) continue;
169552f7358SJed Brown     }
170552f7358SJed Brown     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
171552f7358SJed Brown     for (v = 0; v < closureSize*2; v += 2) {
172552f7358SJed Brown       if ((closure[v] >= vStart) && (closure[v] < vEnd)) piece.nconn++;
173552f7358SJed Brown     }
174552f7358SJed Brown     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
175552f7358SJed Brown     piece.ncells++;
176552f7358SJed Brown   }
177552f7358SJed Brown   if (!rank) {ierr = PetscMalloc(size*sizeof(piece),&gpiece);CHKERRQ(ierr);}
178552f7358SJed Brown   ierr = MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm);CHKERRQ(ierr);
179552f7358SJed Brown 
180552f7358SJed Brown   /*
181552f7358SJed Brown    * Write file header
182552f7358SJed Brown    */
183552f7358SJed Brown   if (!rank) {
184552f7358SJed Brown     PetscInt boffset = 0;
185552f7358SJed Brown 
186552f7358SJed Brown     for (r=0; r<size; r++) {
187552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells);CHKERRQ(ierr);
188552f7358SJed Brown       /* Coordinate positions */
189552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Points>\n");CHKERRQ(ierr);
190552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
191552f7358SJed Brown       boffset += gpiece[r].nvertices*3*sizeof(PetscScalar) + sizeof(int);
192552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Points>\n");CHKERRQ(ierr);
193552f7358SJed Brown       /* Cell connectivity */
194552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Cells>\n");CHKERRQ(ierr);
195552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
196552f7358SJed Brown       boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int);
197552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"offsets\"      NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
198552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int);
199552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"UInt8\" Name=\"types\"        NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
200552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int);
201552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Cells>\n");CHKERRQ(ierr);
202552f7358SJed Brown 
203552f7358SJed Brown       /*
204552f7358SJed Brown        * Cell Data headers
205552f7358SJed Brown        */
206552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <CellData>\n");CHKERRQ(ierr);
207552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
208552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(int) + sizeof(int);
209552f7358SJed Brown       /* all the vectors */
210552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
211552f7358SJed Brown         Vec        X = (Vec)link->vec;
2121cfafdd3SJed Brown         PetscInt   bs,nfields,field;
213552f7358SJed Brown         const char *vecname = "";
214552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
215552f7358SJed 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. */
216552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
217552f7358SJed Brown         }
218552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
2191cfafdd3SJed Brown         ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr);
2207ded4cbaSJed Brown         for (field=0,i=0; field<(nfields?nfields:1); field++) {
2211cfafdd3SJed Brown           PetscInt   fbs,j;
2220298fd71SBarry Smith           const char *fieldname = NULL;
2231cfafdd3SJed Brown           char       buf[256];
2247ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
2251cfafdd3SJed Brown             ierr = PetscSectionGetFieldDof(dm->defaultSection,cStart,field,&fbs);CHKERRQ(ierr);
2261cfafdd3SJed Brown             ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr);
2277ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
228552f7358SJed Brown           if (!fieldname) {
2291cfafdd3SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"CellField%D",field);CHKERRQ(ierr);
230552f7358SJed Brown             fieldname = buf;
231552f7358SJed Brown           }
2321cfafdd3SJed Brown           for (j=0; j<fbs; j++) {
2331cfafdd3SJed 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);
234552f7358SJed Brown             boffset += gpiece[r].ncells*sizeof(PetscScalar) + sizeof(int);
2351cfafdd3SJed Brown             i++;
236552f7358SJed Brown           }
237552f7358SJed Brown         }
2381cfafdd3SJed Brown         if (i != bs) SETERRQ2(comm,PETSC_ERR_PLIB,"Total number of field components %D != block size %D",i,bs);
2391cfafdd3SJed Brown       }
240552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </CellData>\n");CHKERRQ(ierr);
241552f7358SJed Brown 
242552f7358SJed Brown       /*
243552f7358SJed Brown        * Point Data headers
244552f7358SJed Brown        */
245552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      <PointData>\n");CHKERRQ(ierr);
246552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
247552f7358SJed Brown         Vec        X = (Vec)link->vec;
2481cfafdd3SJed Brown         PetscInt   bs,nfields,field;
249552f7358SJed Brown         const char *vecname = "";
250552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
251552f7358SJed 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. */
252552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
253552f7358SJed Brown         }
254552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
2551cfafdd3SJed Brown         ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr);
2567ded4cbaSJed Brown         for (field=0,i=0; field<(nfields?nfields:1); field++) {
2571cfafdd3SJed Brown           PetscInt   fbs,j;
2581cfafdd3SJed Brown           const char *fieldname = NULL;
2591cfafdd3SJed Brown           char       buf[256];
2607ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
2611cfafdd3SJed Brown             ierr = PetscSectionGetFieldDof(dm->defaultSection,vStart,field,&fbs);CHKERRQ(ierr);
2621cfafdd3SJed Brown             ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr);
2637ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
2641cfafdd3SJed Brown           if (!fieldname) {
2651cfafdd3SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"PointField%D",field);CHKERRQ(ierr);
2661cfafdd3SJed Brown             fieldname = buf;
2671cfafdd3SJed Brown           }
2681cfafdd3SJed Brown           for (j=0; j<fbs; j++) {
2691cfafdd3SJed 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);
270552f7358SJed Brown             boffset += gpiece[r].nvertices*sizeof(PetscScalar) + sizeof(int);
271552f7358SJed Brown           }
272552f7358SJed Brown         }
2731cfafdd3SJed Brown       }
274552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </PointData>\n");CHKERRQ(ierr);
275552f7358SJed Brown 
276552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    </Piece>\n");CHKERRQ(ierr);
277552f7358SJed Brown     }
278552f7358SJed Brown   }
279552f7358SJed Brown 
280552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  </UnstructuredGrid>\n");CHKERRQ(ierr);
281552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
282552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
283552f7358SJed Brown 
284552f7358SJed Brown   if (!rank) {
285552f7358SJed Brown     PetscInt maxsize = 0;
286552f7358SJed Brown     for (r=0; r<size; r++) {
287552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscScalar)));
288552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*sizeof(PetscScalar)));
289552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt)));
290552f7358SJed Brown     }
291552f7358SJed Brown     ierr = PetscMalloc(maxsize,&buffer);CHKERRQ(ierr);
292552f7358SJed Brown   }
293552f7358SJed Brown   for (r=0; r<size; r++) {
294552f7358SJed Brown     if (r == rank) {
295552f7358SJed Brown       PetscInt nsend;
296552f7358SJed Brown       {                         /* Position */
297552f7358SJed Brown         const PetscScalar *x;
2980298fd71SBarry Smith         PetscScalar       *y = NULL;
299552f7358SJed Brown         Vec               coords;
300552f7358SJed Brown         nsend = piece.nvertices*3;
301552f7358SJed Brown         ierr  = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
302552f7358SJed Brown         ierr  = VecGetArrayRead(coords,&x);CHKERRQ(ierr);
303552f7358SJed Brown         if (dim != 3) {
304552f7358SJed Brown           ierr = PetscMalloc(piece.nvertices*3*sizeof(PetscScalar),&y);CHKERRQ(ierr);
305552f7358SJed Brown           for (i=0; i<piece.nvertices; i++) {
306552f7358SJed Brown             y[i*3+0] = x[i*dim+0];
307552f7358SJed Brown             y[i*3+1] = (dim > 1) ? x[i*dim+1] : 0;
308552f7358SJed Brown             y[i*3+2] = 0;
309552f7358SJed Brown           }
310552f7358SJed Brown         }
311552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,y ? y : x,buffer,nsend,PETSC_SCALAR,tag);CHKERRQ(ierr);
312552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
313552f7358SJed Brown         ierr = VecRestoreArrayRead(coords,&x);CHKERRQ(ierr);
314552f7358SJed Brown       }
315552f7358SJed Brown       {                           /* Connectivity, offsets, types */
3160298fd71SBarry Smith         PetscVTKInt  *connectivity = NULL,*offsets;
317552f7358SJed Brown         PetscVTKType *types;
318552f7358SJed Brown         ierr = DMPlexGetVTKConnectivity(dm,&piece,&connectivity,&offsets,&types);CHKERRQ(ierr);
319552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,connectivity,buffer,piece.nconn,PETSC_INT32,tag);CHKERRQ(ierr);
320552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,offsets,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr);
321552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,types,buffer,piece.ncells,PETSC_UINT8,tag);CHKERRQ(ierr);
322552f7358SJed Brown         ierr = PetscFree3(connectivity,offsets,types);CHKERRQ(ierr);
323552f7358SJed Brown       }
324552f7358SJed Brown       {                         /* Owners (cell data) */
325552f7358SJed Brown         PetscVTKInt *owners;
326552f7358SJed Brown         ierr = PetscMalloc(piece.ncells*sizeof(PetscVTKInt),&owners);CHKERRQ(ierr);
327552f7358SJed Brown         for (i=0; i<piece.ncells; i++) owners[i] = rank;
328552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,owners,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr);
329552f7358SJed Brown         ierr = PetscFree(owners);CHKERRQ(ierr);
330552f7358SJed Brown       }
331552f7358SJed Brown       /* Cell data */
332552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
333552f7358SJed Brown         Vec               X = (Vec)link->vec;
334552f7358SJed Brown         const PetscScalar *x;
335552f7358SJed Brown         PetscScalar       *y;
336552f7358SJed Brown         PetscInt          bs;
337552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
338552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
339552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
340552f7358SJed Brown         ierr = PetscMalloc(piece.ncells*sizeof(PetscScalar),&y);CHKERRQ(ierr);
341552f7358SJed Brown         for (i=0; i<bs; i++) {
342552f7358SJed Brown           PetscInt cnt;
343552f7358SJed Brown           for (c=cStart,cnt=0; c<cEnd; c++) {
344552f7358SJed Brown             const PetscScalar *xpoint;
345552f7358SJed Brown             if (hasLabel) {     /* Ignore some cells */
346552f7358SJed Brown               PetscInt value;
347552f7358SJed Brown               ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
348552f7358SJed Brown               if (value != 1) continue;
349552f7358SJed Brown             }
350552f7358SJed Brown             ierr     = DMPlexPointLocalRead(dm,c,x,&xpoint);CHKERRQ(ierr);
351552f7358SJed Brown             y[cnt++] = xpoint[i];
352552f7358SJed Brown           }
353552f7358SJed Brown           if (cnt != piece.ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
354552f7358SJed Brown           ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells,PETSC_SCALAR,tag);CHKERRQ(ierr);
355552f7358SJed Brown         }
356552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
357552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
358552f7358SJed Brown       }
359552f7358SJed Brown 
360552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
361552f7358SJed Brown         Vec               X = (Vec)link->vec;
362552f7358SJed Brown         const PetscScalar *x;
363552f7358SJed Brown         PetscScalar       *y;
364552f7358SJed Brown         PetscInt          bs;
365552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
366552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
367552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
368552f7358SJed Brown         ierr = PetscMalloc(piece.nvertices*sizeof(PetscScalar),&y);CHKERRQ(ierr);
369552f7358SJed Brown         for (i=0; i<bs; i++) {
370552f7358SJed Brown           PetscInt cnt;
371552f7358SJed Brown           for (v=vStart,cnt=0; v<vEnd; v++) {
372552f7358SJed Brown             const PetscScalar *xpoint;
37337045cddSJed Brown             ierr     = DMPlexPointLocalRead(dm,v,x,&xpoint);CHKERRQ(ierr);
374552f7358SJed Brown             y[cnt++] = xpoint[i];
375552f7358SJed Brown           }
376552f7358SJed Brown           if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
377552f7358SJed Brown           ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr);
378552f7358SJed Brown         }
379552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
380552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
381552f7358SJed Brown       }
382552f7358SJed Brown     } else if (!rank) {
3830298fd71SBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,PETSC_SCALAR,tag);CHKERRQ(ierr); /* positions */
3840298fd71SBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nconn,PETSC_INT32,tag);CHKERRQ(ierr); /* connectivity */
3850298fd71SBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* offsets */
3860298fd71SBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_UINT8,tag);CHKERRQ(ierr); /* types */
3870298fd71SBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* owner rank (cells) */
388552f7358SJed Brown       /* all cell data */
389552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
390552f7358SJed Brown         PetscInt bs;
391552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
392fd68c46aSJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
393552f7358SJed Brown         for (i=0; i<bs; i++) {
3940298fd71SBarry Smith           ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_SCALAR,tag);CHKERRQ(ierr);
395552f7358SJed Brown         }
396552f7358SJed Brown       }
397552f7358SJed Brown       /* all point data */
398552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
399552f7358SJed Brown         PetscInt bs;
400552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
401552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
402552f7358SJed Brown         for (i=0; i<bs; i++) {
4030298fd71SBarry Smith           ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr);
404552f7358SJed Brown         }
405552f7358SJed Brown       }
406552f7358SJed Brown     }
407552f7358SJed Brown   }
408552f7358SJed Brown   ierr = PetscFree(gpiece);CHKERRQ(ierr);
409552f7358SJed Brown   ierr = PetscFree(buffer);CHKERRQ(ierr);
410552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"\n  </AppendedData>\n");CHKERRQ(ierr);
411552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
412552f7358SJed Brown   PetscFunctionReturn(0);
413552f7358SJed Brown }
414