xref: /petsc/src/dm/impls/plex/plexvtu.c (revision de0a4605751fa9ba53ac89429f80759642e84ce0)
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 
18552f7358SJed Brown static PetscErrorCode TransferWrite(PetscViewer viewer,FILE *fp,PetscMPIInt srank,PetscMPIInt root,const void *send,void *recv,PetscMPIInt count,PetscDataType datatype,PetscMPIInt tag)
19552f7358SJed Brown {
20552f7358SJed Brown   PetscMPIInt    rank;
21552f7358SJed Brown   PetscErrorCode ierr;
22ce94432eSBarry Smith   MPI_Comm       comm;
23552f7358SJed Brown   MPI_Datatype   mpidatatype;
24552f7358SJed Brown 
25552f7358SJed Brown   PetscFunctionBegin;
26ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
27552f7358SJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
28552f7358SJed Brown   ierr = PetscDataTypeToMPIDataType(datatype,&mpidatatype);CHKERRQ(ierr);
29552f7358SJed Brown 
30552f7358SJed Brown   if (rank == srank && rank != root) {
31552f7358SJed Brown     ierr = MPI_Send((void*)send,count,mpidatatype,root,tag,comm);CHKERRQ(ierr);
32552f7358SJed Brown   } else if (rank == root) {
33552f7358SJed Brown     const void *buffer;
34552f7358SJed Brown     if (root == srank) {        /* self */
35552f7358SJed Brown       buffer = send;
36552f7358SJed Brown     } else {
37552f7358SJed Brown       MPI_Status  status;
38552f7358SJed Brown       PetscMPIInt nrecv;
39552f7358SJed Brown       ierr = MPI_Recv(recv,count,mpidatatype,srank,tag,comm,&status);CHKERRQ(ierr);
40552f7358SJed Brown       ierr = MPI_Get_count(&status,mpidatatype,&nrecv);CHKERRQ(ierr);
41552f7358SJed Brown       if (count != nrecv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
42552f7358SJed Brown       buffer = recv;
43552f7358SJed Brown     }
44552f7358SJed Brown     ierr = PetscViewerVTKFWrite(viewer,fp,buffer,count,datatype);CHKERRQ(ierr);
45552f7358SJed Brown   }
46552f7358SJed Brown   PetscFunctionReturn(0);
47552f7358SJed Brown }
48552f7358SJed Brown 
49552f7358SJed Brown static PetscErrorCode DMPlexGetVTKConnectivity(DM dm,PieceInfo *piece,PetscVTKInt **oconn,PetscVTKInt **ooffsets,PetscVTKType **otypes)
50552f7358SJed Brown {
51552f7358SJed Brown   PetscErrorCode ierr;
52552f7358SJed Brown   PetscVTKInt    *conn,*offsets;
53552f7358SJed Brown   PetscVTKType   *types;
54552f7358SJed Brown   PetscInt       dim,vStart,vEnd,cStart,cEnd,pStart,pEnd,cellHeight,cMax,numLabelCells,hasLabel,c,v,countcell,countconn;
55552f7358SJed Brown 
56552f7358SJed Brown   PetscFunctionBegin;
57dcca6d9dSJed Brown   ierr = PetscMalloc3(piece->nconn,&conn,piece->ncells,&offsets,piece->ncells,&types);CHKERRQ(ierr);
58552f7358SJed Brown 
59c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
60552f7358SJed Brown   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
61552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
62552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
63552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
640298fd71SBarry Smith   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
658865f1eaSKarl Rupp   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
66c58f1c22SToby Isaac   ierr     = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
67552f7358SJed Brown   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
68552f7358SJed Brown 
69552f7358SJed Brown   countcell = 0;
70552f7358SJed Brown   countconn = 0;
71552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
720298fd71SBarry Smith     PetscInt *closure = NULL;
73724b5320SMatthew G. Knepley     PetscInt  closureSize,nverts,celltype,startoffset,nC=0;
74552f7358SJed Brown 
75552f7358SJed Brown     if (hasLabel) {
76552f7358SJed Brown       PetscInt value;
77552f7358SJed Brown 
78c58f1c22SToby Isaac       ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
79552f7358SJed Brown       if (value != 1) continue;
80552f7358SJed Brown     }
81552f7358SJed Brown     startoffset = countconn;
82552f7358SJed Brown     ierr        = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
83552f7358SJed Brown     for (v = 0; v < closureSize*2; v += 2) {
84552f7358SJed Brown       if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
85552f7358SJed Brown         conn[countconn++] = closure[v] - vStart;
86724b5320SMatthew G. Knepley         ++nC;
87552f7358SJed Brown       }
88552f7358SJed Brown     }
89552f7358SJed Brown     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
90724b5320SMatthew G. Knepley     ierr = DMPlexInvertCell(dim, nC, &conn[countconn-nC]);CHKERRQ(ierr);
918865f1eaSKarl Rupp 
92552f7358SJed Brown     offsets[countcell] = countconn;
938865f1eaSKarl Rupp 
94552f7358SJed Brown     nverts = countconn - startoffset;
95*de0a4605SMatthew G. Knepley     ierr   = DMPlexVTKGetCellType_Internal(dm,dim,nverts,&celltype);CHKERRQ(ierr);
968865f1eaSKarl Rupp 
97552f7358SJed Brown     types[countcell] = celltype;
98552f7358SJed Brown     countcell++;
99552f7358SJed Brown   }
100552f7358SJed Brown   if (countcell != piece->ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count");
101552f7358SJed Brown   if (countconn != piece->nconn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count");
102552f7358SJed Brown   *oconn    = conn;
103552f7358SJed Brown   *ooffsets = offsets;
104552f7358SJed Brown   *otypes   = types;
105552f7358SJed Brown   PetscFunctionReturn(0);
106552f7358SJed Brown }
107552f7358SJed Brown 
108552f7358SJed Brown /*
109552f7358SJed Brown   Write all fields that have been provided to the viewer
110552f7358SJed Brown   Multi-block XML format with binary appended data.
111552f7358SJed Brown */
112552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer)
113552f7358SJed Brown {
114ce94432eSBarry Smith   MPI_Comm                 comm;
115552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
116552f7358SJed Brown   PetscViewerVTKObjectLink link;
117552f7358SJed Brown   FILE                     *fp;
118552f7358SJed Brown   PetscMPIInt              rank,size,tag;
119552f7358SJed Brown   PetscErrorCode           ierr;
1203baa072aSToby Isaac   PetscInt                 dimEmbed,cellHeight,cStart,cEnd,vStart,vEnd,cMax,numLabelCells,hasLabel,c,v,r,i;
1210298fd71SBarry Smith   PieceInfo                piece,*gpiece = NULL;
1220298fd71SBarry Smith   void                     *buffer = NULL;
123552f7358SJed Brown 
124552f7358SJed Brown   PetscFunctionBegin;
125ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
126552f7358SJed Brown #if defined(PETSC_USE_COMPLEX)
127ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Complex values not supported");
128552f7358SJed Brown #endif
129552f7358SJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
130552f7358SJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
131552f7358SJed Brown   tag  = ((PetscObject)viewer)->tag;
132552f7358SJed Brown 
133552f7358SJed Brown   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
134552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
135519f805aSKarl Rupp #if defined(PETSC_WORDS_BIGENDIAN)
136552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr);
137552f7358SJed Brown #else
138552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr);
139552f7358SJed Brown #endif
140552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <UnstructuredGrid>\n");CHKERRQ(ierr);
141552f7358SJed Brown 
1423baa072aSToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
143552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
144552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
145552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1460298fd71SBarry Smith   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
1478865f1eaSKarl Rupp   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
148c58f1c22SToby Isaac   ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
1498865f1eaSKarl Rupp 
150552f7358SJed Brown   hasLabel        = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
151552f7358SJed Brown   piece.nvertices = vEnd - vStart;
152552f7358SJed Brown   piece.ncells    = 0;
153552f7358SJed Brown   piece.nconn     = 0;
154552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
1550298fd71SBarry Smith     PetscInt *closure = NULL;
156552f7358SJed Brown     PetscInt closureSize;
157552f7358SJed Brown 
158552f7358SJed Brown     if (hasLabel) {
159552f7358SJed Brown       PetscInt value;
160552f7358SJed Brown 
161c58f1c22SToby Isaac       ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
162552f7358SJed Brown       if (value != 1) continue;
163552f7358SJed Brown     }
164552f7358SJed Brown     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
165552f7358SJed Brown     for (v = 0; v < closureSize*2; v += 2) {
166552f7358SJed Brown       if ((closure[v] >= vStart) && (closure[v] < vEnd)) piece.nconn++;
167552f7358SJed Brown     }
168552f7358SJed Brown     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
169552f7358SJed Brown     piece.ncells++;
170552f7358SJed Brown   }
171785e854fSJed Brown   if (!rank) {ierr = PetscMalloc1(size,&gpiece);CHKERRQ(ierr);}
172552f7358SJed Brown   ierr = MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm);CHKERRQ(ierr);
173552f7358SJed Brown 
174552f7358SJed Brown   /*
175552f7358SJed Brown    * Write file header
176552f7358SJed Brown    */
177552f7358SJed Brown   if (!rank) {
178552f7358SJed Brown     PetscInt boffset = 0;
179552f7358SJed Brown 
180552f7358SJed Brown     for (r=0; r<size; r++) {
181552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells);CHKERRQ(ierr);
182552f7358SJed Brown       /* Coordinate positions */
183552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Points>\n");CHKERRQ(ierr);
184552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
185552f7358SJed Brown       boffset += gpiece[r].nvertices*3*sizeof(PetscScalar) + sizeof(int);
186552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Points>\n");CHKERRQ(ierr);
187552f7358SJed Brown       /* Cell connectivity */
188552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Cells>\n");CHKERRQ(ierr);
189552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
190552f7358SJed Brown       boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int);
191552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"offsets\"      NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
192552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int);
193552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"UInt8\" Name=\"types\"        NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
194552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int);
195552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Cells>\n");CHKERRQ(ierr);
196552f7358SJed Brown 
197552f7358SJed Brown       /*
198552f7358SJed Brown        * Cell Data headers
199552f7358SJed Brown        */
200552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <CellData>\n");CHKERRQ(ierr);
201552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
202552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(int) + sizeof(int);
203552f7358SJed Brown       /* all the vectors */
204552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
205552f7358SJed Brown         Vec        X = (Vec)link->vec;
2061cfafdd3SJed Brown         PetscInt   bs,nfields,field;
207552f7358SJed Brown         const char *vecname = "";
208552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
209552f7358SJed 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. */
210552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
211552f7358SJed Brown         }
212552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
2131cfafdd3SJed Brown         ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr);
2147ded4cbaSJed Brown         for (field=0,i=0; field<(nfields?nfields:1); field++) {
2151cfafdd3SJed Brown           PetscInt     fbs,j;
216a00cdb45SToby Isaac           PetscFV      fv = NULL;
217a00cdb45SToby Isaac           PetscObject  f;
218a00cdb45SToby Isaac           PetscClassId fClass;
2190298fd71SBarry Smith           const char *fieldname = NULL;
2201cfafdd3SJed Brown           char       buf[256];
2217ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
2221cfafdd3SJed Brown             ierr = PetscSectionGetFieldDof(dm->defaultSection,cStart,field,&fbs);CHKERRQ(ierr);
2231cfafdd3SJed Brown             ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr);
2247ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
225a00cdb45SToby Isaac           ierr = DMGetField(dm,field,&f);CHKERRQ(ierr);
226a00cdb45SToby Isaac           ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr);
227a00cdb45SToby Isaac           if (fClass == PETSCFV_CLASSID) {
228a00cdb45SToby Isaac             fv = (PetscFV) f;
229a00cdb45SToby Isaac           }
230552f7358SJed Brown           if (!fieldname) {
2311cfafdd3SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"CellField%D",field);CHKERRQ(ierr);
232552f7358SJed Brown             fieldname = buf;
233552f7358SJed Brown           }
2341cfafdd3SJed Brown           for (j=0; j<fbs; j++) {
235a00cdb45SToby Isaac             const char *compName = NULL;
236a00cdb45SToby Isaac             if (fv) {
237a00cdb45SToby Isaac               ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
238a00cdb45SToby Isaac             }
239a00cdb45SToby Isaac             if (compName) {
240a00cdb45SToby 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);
241a00cdb45SToby Isaac             }
242a00cdb45SToby Isaac             else {
2431cfafdd3SJed 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);
244a00cdb45SToby Isaac             }
245552f7358SJed Brown             boffset += gpiece[r].ncells*sizeof(PetscScalar) + sizeof(int);
2461cfafdd3SJed Brown             i++;
247552f7358SJed Brown           }
248552f7358SJed Brown         }
2491cfafdd3SJed Brown         if (i != bs) SETERRQ2(comm,PETSC_ERR_PLIB,"Total number of field components %D != block size %D",i,bs);
2501cfafdd3SJed Brown       }
251552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </CellData>\n");CHKERRQ(ierr);
252552f7358SJed Brown 
253552f7358SJed Brown       /*
254552f7358SJed Brown        * Point Data headers
255552f7358SJed Brown        */
256552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      <PointData>\n");CHKERRQ(ierr);
257552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
258552f7358SJed Brown         Vec        X = (Vec)link->vec;
2591cfafdd3SJed Brown         PetscInt   bs,nfields,field;
260552f7358SJed Brown         const char *vecname = "";
261552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
262552f7358SJed 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. */
263552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
264552f7358SJed Brown         }
265552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
2661cfafdd3SJed Brown         ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr);
2677ded4cbaSJed Brown         for (field=0,i=0; field<(nfields?nfields:1); field++) {
2681cfafdd3SJed Brown           PetscInt   fbs,j;
2691cfafdd3SJed Brown           const char *fieldname = NULL;
2701cfafdd3SJed Brown           char       buf[256];
2717ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
2721cfafdd3SJed Brown             ierr = PetscSectionGetFieldDof(dm->defaultSection,vStart,field,&fbs);CHKERRQ(ierr);
2731cfafdd3SJed Brown             ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr);
2747ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
2751cfafdd3SJed Brown           if (!fieldname) {
2761cfafdd3SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"PointField%D",field);CHKERRQ(ierr);
2771cfafdd3SJed Brown             fieldname = buf;
2781cfafdd3SJed Brown           }
2791cfafdd3SJed Brown           for (j=0; j<fbs; j++) {
2801cfafdd3SJed 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);
281552f7358SJed Brown             boffset += gpiece[r].nvertices*sizeof(PetscScalar) + sizeof(int);
282552f7358SJed Brown           }
283552f7358SJed Brown         }
2841cfafdd3SJed Brown       }
285552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </PointData>\n");CHKERRQ(ierr);
286552f7358SJed Brown 
287552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    </Piece>\n");CHKERRQ(ierr);
288552f7358SJed Brown     }
289552f7358SJed Brown   }
290552f7358SJed Brown 
291552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  </UnstructuredGrid>\n");CHKERRQ(ierr);
292552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
293552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
294552f7358SJed Brown 
295552f7358SJed Brown   if (!rank) {
296552f7358SJed Brown     PetscInt maxsize = 0;
297552f7358SJed Brown     for (r=0; r<size; r++) {
298552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscScalar)));
299552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*sizeof(PetscScalar)));
300552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt)));
301552f7358SJed Brown     }
302552f7358SJed Brown     ierr = PetscMalloc(maxsize,&buffer);CHKERRQ(ierr);
303552f7358SJed Brown   }
304552f7358SJed Brown   for (r=0; r<size; r++) {
305552f7358SJed Brown     if (r == rank) {
306552f7358SJed Brown       PetscInt nsend;
307552f7358SJed Brown       {                         /* Position */
308552f7358SJed Brown         const PetscScalar *x;
3090298fd71SBarry Smith         PetscScalar       *y = NULL;
310552f7358SJed Brown         Vec               coords;
311552f7358SJed Brown         nsend = piece.nvertices*3;
312552f7358SJed Brown         ierr  = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
313552f7358SJed Brown         ierr  = VecGetArrayRead(coords,&x);CHKERRQ(ierr);
3143baa072aSToby Isaac         if (dimEmbed != 3) {
315785e854fSJed Brown           ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr);
316552f7358SJed Brown           for (i=0; i<piece.nvertices; i++) {
3173baa072aSToby Isaac             y[i*3+0] = x[i*dimEmbed+0];
3183baa072aSToby Isaac             y[i*3+1] = (dimEmbed > 1) ? x[i*dimEmbed+1] : 0;
319552f7358SJed Brown             y[i*3+2] = 0;
320552f7358SJed Brown           }
321552f7358SJed Brown         }
322552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,y ? y : x,buffer,nsend,PETSC_SCALAR,tag);CHKERRQ(ierr);
323552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
324552f7358SJed Brown         ierr = VecRestoreArrayRead(coords,&x);CHKERRQ(ierr);
325552f7358SJed Brown       }
326552f7358SJed Brown       {                           /* Connectivity, offsets, types */
3273e869605SMatthew G. Knepley         PetscVTKInt  *connectivity = NULL, *offsets = NULL;
3283e869605SMatthew G. Knepley         PetscVTKType *types = NULL;
329552f7358SJed Brown         ierr = DMPlexGetVTKConnectivity(dm,&piece,&connectivity,&offsets,&types);CHKERRQ(ierr);
330552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,connectivity,buffer,piece.nconn,PETSC_INT32,tag);CHKERRQ(ierr);
331552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,offsets,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr);
332552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,types,buffer,piece.ncells,PETSC_UINT8,tag);CHKERRQ(ierr);
333552f7358SJed Brown         ierr = PetscFree3(connectivity,offsets,types);CHKERRQ(ierr);
334552f7358SJed Brown       }
335552f7358SJed Brown       {                         /* Owners (cell data) */
336552f7358SJed Brown         PetscVTKInt *owners;
337785e854fSJed Brown         ierr = PetscMalloc1(piece.ncells,&owners);CHKERRQ(ierr);
338552f7358SJed Brown         for (i=0; i<piece.ncells; i++) owners[i] = rank;
339552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,owners,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr);
340552f7358SJed Brown         ierr = PetscFree(owners);CHKERRQ(ierr);
341552f7358SJed Brown       }
342552f7358SJed Brown       /* Cell data */
343552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
344552f7358SJed Brown         Vec               X = (Vec)link->vec;
345552f7358SJed Brown         const PetscScalar *x;
346552f7358SJed Brown         PetscScalar       *y;
347552f7358SJed Brown         PetscInt          bs;
348552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
349552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
350552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
351785e854fSJed Brown         ierr = PetscMalloc1(piece.ncells,&y);CHKERRQ(ierr);
352552f7358SJed Brown         for (i=0; i<bs; i++) {
353552f7358SJed Brown           PetscInt cnt;
354552f7358SJed Brown           for (c=cStart,cnt=0; c<cEnd; c++) {
355640bce14SSatish Balay             PetscScalar *xpoint;
356552f7358SJed Brown             if (hasLabel) {     /* Ignore some cells */
357552f7358SJed Brown               PetscInt value;
358c58f1c22SToby Isaac               ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
359552f7358SJed Brown               if (value != 1) continue;
360552f7358SJed Brown             }
361552f7358SJed Brown             ierr     = DMPlexPointLocalRead(dm,c,x,&xpoint);CHKERRQ(ierr);
362552f7358SJed Brown             y[cnt++] = xpoint[i];
363552f7358SJed Brown           }
364552f7358SJed Brown           if (cnt != piece.ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
365552f7358SJed Brown           ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells,PETSC_SCALAR,tag);CHKERRQ(ierr);
366552f7358SJed Brown         }
367552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
368552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
369552f7358SJed Brown       }
370552f7358SJed Brown 
371552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
372552f7358SJed Brown         Vec               X = (Vec)link->vec;
373552f7358SJed Brown         const PetscScalar *x;
374552f7358SJed Brown         PetscScalar       *y;
375552f7358SJed Brown         PetscInt          bs;
376552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
377552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
378552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
379785e854fSJed Brown         ierr = PetscMalloc1(piece.nvertices,&y);CHKERRQ(ierr);
380552f7358SJed Brown         for (i=0; i<bs; i++) {
381552f7358SJed Brown           PetscInt cnt;
382552f7358SJed Brown           for (v=vStart,cnt=0; v<vEnd; v++) {
383640bce14SSatish Balay             PetscScalar *xpoint;
38437045cddSJed Brown             ierr     = DMPlexPointLocalRead(dm,v,x,&xpoint);CHKERRQ(ierr);
385552f7358SJed Brown             y[cnt++] = xpoint[i];
386552f7358SJed Brown           }
387552f7358SJed Brown           if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
388552f7358SJed Brown           ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr);
389552f7358SJed Brown         }
390552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
391552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
392552f7358SJed Brown       }
393552f7358SJed Brown     } else if (!rank) {
3940298fd71SBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,PETSC_SCALAR,tag);CHKERRQ(ierr); /* positions */
3950298fd71SBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nconn,PETSC_INT32,tag);CHKERRQ(ierr); /* connectivity */
3960298fd71SBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* offsets */
3970298fd71SBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_UINT8,tag);CHKERRQ(ierr); /* types */
3980298fd71SBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* owner rank (cells) */
399552f7358SJed Brown       /* all cell data */
400552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
401552f7358SJed Brown         PetscInt bs;
402552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
403fd68c46aSJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
404552f7358SJed Brown         for (i=0; i<bs; i++) {
4050298fd71SBarry Smith           ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_SCALAR,tag);CHKERRQ(ierr);
406552f7358SJed Brown         }
407552f7358SJed Brown       }
408552f7358SJed Brown       /* all point data */
409552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
410552f7358SJed Brown         PetscInt bs;
411552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
412552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
413552f7358SJed Brown         for (i=0; i<bs; i++) {
4140298fd71SBarry Smith           ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr);
415552f7358SJed Brown         }
416552f7358SJed Brown       }
417552f7358SJed Brown     }
418552f7358SJed Brown   }
419552f7358SJed Brown   ierr = PetscFree(gpiece);CHKERRQ(ierr);
420552f7358SJed Brown   ierr = PetscFree(buffer);CHKERRQ(ierr);
421552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"\n  </AppendedData>\n");CHKERRQ(ierr);
422552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
4235e32de72SMatthew G. Knepley   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
424552f7358SJed Brown   PetscFunctionReturn(0);
425552f7358SJed Brown }
426