xref: /petsc/src/dm/impls/plex/plexvtu.c (revision 2e529ec88214be11db4abfc2ebf7b2fbd6e371cc)
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 
49*2e529ec8SStefano Zampini static PetscErrorCode DMPlexGetVTKConnectivity(DM dm, PetscBool localized, PieceInfo *piece,PetscVTKInt **oconn,PetscVTKInt **ooffsets,PetscVTKType **otypes)
50552f7358SJed Brown {
51552f7358SJed Brown   PetscErrorCode ierr;
52*2e529ec8SStefano Zampini   PetscSection   coordSection;
53552f7358SJed Brown   PetscVTKInt    *conn,*offsets;
54552f7358SJed Brown   PetscVTKType   *types;
55552f7358SJed Brown   PetscInt       dim,vStart,vEnd,cStart,cEnd,pStart,pEnd,cellHeight,cMax,numLabelCells,hasLabel,c,v,countcell,countconn;
56552f7358SJed Brown 
57552f7358SJed Brown   PetscFunctionBegin;
58dcca6d9dSJed Brown   ierr = PetscMalloc3(piece->nconn,&conn,piece->ncells,&offsets,piece->ncells,&types);CHKERRQ(ierr);
59*2e529ec8SStefano Zampini   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
60c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
61552f7358SJed Brown   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
62552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
63552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
64552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
650298fd71SBarry Smith   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
668865f1eaSKarl Rupp   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
67c58f1c22SToby Isaac   ierr     = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
68552f7358SJed Brown   hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
69552f7358SJed Brown 
70552f7358SJed Brown   countcell = 0;
71552f7358SJed Brown   countconn = 0;
72552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
73*2e529ec8SStefano Zampini     PetscInt 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;
82*2e529ec8SStefano Zampini     if (!localized) {
83*2e529ec8SStefano Zampini       PetscInt *closure = NULL;
84*2e529ec8SStefano Zampini       PetscInt  closureSize;
85*2e529ec8SStefano Zampini 
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);
94*2e529ec8SStefano Zampini     } else {
95*2e529ec8SStefano Zampini       PetscInt dof;
96*2e529ec8SStefano Zampini 
97*2e529ec8SStefano Zampini       ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
98*2e529ec8SStefano Zampini       for (nC = 0; nC < dof/dim; nC++) conn[countconn++] = startoffset + nC;
99*2e529ec8SStefano Zampini     }
100724b5320SMatthew G. Knepley     ierr = DMPlexInvertCell(dim, nC, &conn[countconn-nC]);CHKERRQ(ierr);
1018865f1eaSKarl Rupp 
102552f7358SJed Brown     offsets[countcell] = countconn;
1038865f1eaSKarl Rupp 
104552f7358SJed Brown     nverts = countconn - startoffset;
105de0a4605SMatthew G. Knepley     ierr   = DMPlexVTKGetCellType_Internal(dm,dim,nverts,&celltype);CHKERRQ(ierr);
1068865f1eaSKarl Rupp 
107552f7358SJed Brown     types[countcell] = celltype;
108552f7358SJed Brown     countcell++;
109552f7358SJed Brown   }
110552f7358SJed Brown   if (countcell != piece->ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count");
111552f7358SJed Brown   if (countconn != piece->nconn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count");
112552f7358SJed Brown   *oconn    = conn;
113552f7358SJed Brown   *ooffsets = offsets;
114552f7358SJed Brown   *otypes   = types;
115552f7358SJed Brown   PetscFunctionReturn(0);
116552f7358SJed Brown }
117552f7358SJed Brown 
118552f7358SJed Brown /*
119552f7358SJed Brown   Write all fields that have been provided to the viewer
120552f7358SJed Brown   Multi-block XML format with binary appended data.
121552f7358SJed Brown */
122552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer)
123552f7358SJed Brown {
124ce94432eSBarry Smith   MPI_Comm                 comm;
125*2e529ec8SStefano Zampini   PetscSection             coordSection;
126552f7358SJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
127552f7358SJed Brown   PetscViewerVTKObjectLink link;
128552f7358SJed Brown   FILE                     *fp;
129552f7358SJed Brown   PetscMPIInt              rank,size,tag;
130552f7358SJed Brown   PetscErrorCode           ierr;
131*2e529ec8SStefano Zampini   PetscInt                 dim, dimEmbed,cellHeight,cStart,cEnd,vStart,vEnd,cMax,numLabelCells,hasLabel,c,v,r,i;
132*2e529ec8SStefano Zampini   PetscBool                localized;
1330298fd71SBarry Smith   PieceInfo                piece,*gpiece = NULL;
1340298fd71SBarry Smith   void                     *buffer = NULL;
135552f7358SJed Brown 
136552f7358SJed Brown   PetscFunctionBegin;
137ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
138552f7358SJed Brown #if defined(PETSC_USE_COMPLEX)
139ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Complex values not supported");
140552f7358SJed Brown #endif
141552f7358SJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
142552f7358SJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
143552f7358SJed Brown   tag  = ((PetscObject)viewer)->tag;
144552f7358SJed Brown 
145552f7358SJed Brown   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
146552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
147519f805aSKarl Rupp #if defined(PETSC_WORDS_BIGENDIAN)
148552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr);
149552f7358SJed Brown #else
150552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr);
151552f7358SJed Brown #endif
152552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <UnstructuredGrid>\n");CHKERRQ(ierr);
153552f7358SJed Brown 
1543baa072aSToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
155552f7358SJed Brown   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
156552f7358SJed Brown   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
157552f7358SJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1580298fd71SBarry Smith   ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr);
1598865f1eaSKarl Rupp   if (cMax >= 0) cEnd = PetscMin(cEnd, cMax);
160c58f1c22SToby Isaac   ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr);
161*2e529ec8SStefano Zampini   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
162*2e529ec8SStefano Zampini   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
1638865f1eaSKarl Rupp 
164552f7358SJed Brown   hasLabel        = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE;
165*2e529ec8SStefano Zampini   piece.nvertices = 0;
166552f7358SJed Brown   piece.ncells    = 0;
167552f7358SJed Brown   piece.nconn     = 0;
168*2e529ec8SStefano Zampini   if (!localized) piece.nvertices = vEnd - vStart;
169552f7358SJed Brown   for (c = cStart; c < cEnd; ++c) {
170552f7358SJed Brown     if (hasLabel) {
171552f7358SJed Brown       PetscInt value;
172552f7358SJed Brown 
173c58f1c22SToby Isaac       ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
174552f7358SJed Brown       if (value != 1) continue;
175552f7358SJed Brown     }
176*2e529ec8SStefano Zampini     if (!localized) {
177*2e529ec8SStefano Zampini       PetscInt *closure = NULL;
178*2e529ec8SStefano Zampini       PetscInt closureSize;
179*2e529ec8SStefano Zampini 
180552f7358SJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
181552f7358SJed Brown       for (v = 0; v < closureSize*2; v += 2) {
182*2e529ec8SStefano Zampini         if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
183*2e529ec8SStefano Zampini           piece.nconn++;
184*2e529ec8SStefano Zampini         }
185552f7358SJed Brown       }
186552f7358SJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
187*2e529ec8SStefano Zampini     } else {
188*2e529ec8SStefano Zampini       PetscInt dof;
189*2e529ec8SStefano Zampini 
190*2e529ec8SStefano Zampini       ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
191*2e529ec8SStefano Zampini       piece.nvertices += dof/dimEmbed;
192*2e529ec8SStefano Zampini       piece.nconn     += dof/dimEmbed;
193*2e529ec8SStefano Zampini     }
194552f7358SJed Brown     piece.ncells++;
195552f7358SJed Brown   }
196785e854fSJed Brown   if (!rank) {ierr = PetscMalloc1(size,&gpiece);CHKERRQ(ierr);}
197552f7358SJed Brown   ierr = MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm);CHKERRQ(ierr);
198552f7358SJed Brown 
199552f7358SJed Brown   /*
200552f7358SJed Brown    * Write file header
201552f7358SJed Brown    */
202552f7358SJed Brown   if (!rank) {
203552f7358SJed Brown     PetscInt boffset = 0;
204552f7358SJed Brown 
205552f7358SJed Brown     for (r=0; r<size; r++) {
206552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells);CHKERRQ(ierr);
207552f7358SJed Brown       /* Coordinate positions */
208552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Points>\n");CHKERRQ(ierr);
209552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
210552f7358SJed Brown       boffset += gpiece[r].nvertices*3*sizeof(PetscScalar) + sizeof(int);
211552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Points>\n");CHKERRQ(ierr);
212552f7358SJed Brown       /* Cell connectivity */
213552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <Cells>\n");CHKERRQ(ierr);
214552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
215552f7358SJed Brown       boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int);
216552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"offsets\"      NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
217552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int);
218552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"UInt8\" Name=\"types\"        NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
219552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int);
220552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      </Cells>\n");CHKERRQ(ierr);
221552f7358SJed Brown 
222552f7358SJed Brown       /*
223552f7358SJed Brown        * Cell Data headers
224552f7358SJed Brown        */
225552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"      <CellData>\n");CHKERRQ(ierr);
226552f7358SJed Brown       ierr     = PetscFPrintf(PETSC_COMM_SELF,fp,"        <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr);
227552f7358SJed Brown       boffset += gpiece[r].ncells*sizeof(int) + sizeof(int);
228552f7358SJed Brown       /* all the vectors */
229552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
230552f7358SJed Brown         Vec        X = (Vec)link->vec;
2311cfafdd3SJed Brown         PetscInt   bs,nfields,field;
232552f7358SJed Brown         const char *vecname = "";
233552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
234552f7358SJed 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. */
235552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
236552f7358SJed Brown         }
237552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
2381cfafdd3SJed Brown         ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr);
2397ded4cbaSJed Brown         for (field=0,i=0; field<(nfields?nfields:1); field++) {
2401cfafdd3SJed Brown           PetscInt     fbs,j;
241a00cdb45SToby Isaac           PetscFV      fv = NULL;
242a00cdb45SToby Isaac           PetscObject  f;
243a00cdb45SToby Isaac           PetscClassId fClass;
2440298fd71SBarry Smith           const char *fieldname = NULL;
2451cfafdd3SJed Brown           char       buf[256];
2467ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
2471cfafdd3SJed Brown             ierr = PetscSectionGetFieldDof(dm->defaultSection,cStart,field,&fbs);CHKERRQ(ierr);
2481cfafdd3SJed Brown             ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr);
2497ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
250a00cdb45SToby Isaac           ierr = DMGetField(dm,field,&f);CHKERRQ(ierr);
251a00cdb45SToby Isaac           ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr);
252a00cdb45SToby Isaac           if (fClass == PETSCFV_CLASSID) {
253a00cdb45SToby Isaac             fv = (PetscFV) f;
254a00cdb45SToby Isaac           }
255552f7358SJed Brown           if (!fieldname) {
2561cfafdd3SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"CellField%D",field);CHKERRQ(ierr);
257552f7358SJed Brown             fieldname = buf;
258552f7358SJed Brown           }
2591cfafdd3SJed Brown           for (j=0; j<fbs; j++) {
260a00cdb45SToby Isaac             const char *compName = NULL;
261a00cdb45SToby Isaac             if (fv) {
262a00cdb45SToby Isaac               ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr);
263a00cdb45SToby Isaac             }
264a00cdb45SToby Isaac             if (compName) {
265a00cdb45SToby 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);
266a00cdb45SToby Isaac             }
267a00cdb45SToby Isaac             else {
2681cfafdd3SJed 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);
269a00cdb45SToby Isaac             }
270552f7358SJed Brown             boffset += gpiece[r].ncells*sizeof(PetscScalar) + sizeof(int);
2711cfafdd3SJed Brown             i++;
272552f7358SJed Brown           }
273552f7358SJed Brown         }
2741cfafdd3SJed Brown         if (i != bs) SETERRQ2(comm,PETSC_ERR_PLIB,"Total number of field components %D != block size %D",i,bs);
2751cfafdd3SJed Brown       }
276552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </CellData>\n");CHKERRQ(ierr);
277552f7358SJed Brown 
278552f7358SJed Brown       /*
279552f7358SJed Brown        * Point Data headers
280552f7358SJed Brown        */
281552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      <PointData>\n");CHKERRQ(ierr);
282552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
283552f7358SJed Brown         Vec        X = (Vec)link->vec;
2841cfafdd3SJed Brown         PetscInt   bs,nfields,field;
285552f7358SJed Brown         const char *vecname = "";
286552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
287552f7358SJed 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. */
288552f7358SJed Brown           ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
289552f7358SJed Brown         }
290552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
2911cfafdd3SJed Brown         ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr);
2927ded4cbaSJed Brown         for (field=0,i=0; field<(nfields?nfields:1); field++) {
2931cfafdd3SJed Brown           PetscInt   fbs,j;
2941cfafdd3SJed Brown           const char *fieldname = NULL;
2951cfafdd3SJed Brown           char       buf[256];
2967ded4cbaSJed Brown           if (nfields) {        /* We have user-defined fields/components */
2971cfafdd3SJed Brown             ierr = PetscSectionGetFieldDof(dm->defaultSection,vStart,field,&fbs);CHKERRQ(ierr);
2981cfafdd3SJed Brown             ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr);
2997ded4cbaSJed Brown           } else fbs = bs;      /* Say we have one field with 'bs' components */
3001cfafdd3SJed Brown           if (!fieldname) {
3011cfafdd3SJed Brown             ierr = PetscSNPrintf(buf,sizeof(buf),"PointField%D",field);CHKERRQ(ierr);
3021cfafdd3SJed Brown             fieldname = buf;
3031cfafdd3SJed Brown           }
3041cfafdd3SJed Brown           for (j=0; j<fbs; j++) {
3051cfafdd3SJed 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);
306552f7358SJed Brown             boffset += gpiece[r].nvertices*sizeof(PetscScalar) + sizeof(int);
307552f7358SJed Brown           }
308552f7358SJed Brown         }
3091cfafdd3SJed Brown       }
310552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"      </PointData>\n");CHKERRQ(ierr);
311552f7358SJed Brown 
312552f7358SJed Brown       ierr = PetscFPrintf(PETSC_COMM_SELF,fp,"    </Piece>\n");CHKERRQ(ierr);
313552f7358SJed Brown     }
314552f7358SJed Brown   }
315552f7358SJed Brown 
316552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  </UnstructuredGrid>\n");CHKERRQ(ierr);
317552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
318552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
319552f7358SJed Brown 
320552f7358SJed Brown   if (!rank) {
321552f7358SJed Brown     PetscInt maxsize = 0;
322552f7358SJed Brown     for (r=0; r<size; r++) {
323552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscScalar)));
324552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*sizeof(PetscScalar)));
325552f7358SJed Brown       maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt)));
326552f7358SJed Brown     }
327552f7358SJed Brown     ierr = PetscMalloc(maxsize,&buffer);CHKERRQ(ierr);
328552f7358SJed Brown   }
329552f7358SJed Brown   for (r=0; r<size; r++) {
330552f7358SJed Brown     if (r == rank) {
331552f7358SJed Brown       PetscInt nsend;
332552f7358SJed Brown       {                         /* Position */
333552f7358SJed Brown         const PetscScalar *x;
3340298fd71SBarry Smith         PetscScalar       *y = NULL;
335552f7358SJed Brown         Vec               coords;
336*2e529ec8SStefano Zampini 
337552f7358SJed Brown         ierr  = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
338552f7358SJed Brown         ierr  = VecGetArrayRead(coords,&x);CHKERRQ(ierr);
3393baa072aSToby Isaac         if (dimEmbed != 3) {
340785e854fSJed Brown           ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr);
341552f7358SJed Brown           for (i=0; i<piece.nvertices; i++) {
3423baa072aSToby Isaac             y[i*3+0] = x[i*dimEmbed+0];
3433baa072aSToby Isaac             y[i*3+1] = (dimEmbed > 1) ? x[i*dimEmbed+1] : 0;
344552f7358SJed Brown             y[i*3+2] = 0;
345552f7358SJed Brown           }
346552f7358SJed Brown         }
347*2e529ec8SStefano Zampini         if (localized) {
348*2e529ec8SStefano Zampini           PetscInt cnt = 0;
349*2e529ec8SStefano Zampini           if (!y) {
350*2e529ec8SStefano Zampini             PetscInt off;
351*2e529ec8SStefano Zampini 
352*2e529ec8SStefano Zampini             ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr);
353*2e529ec8SStefano Zampini             ierr = PetscSectionGetOffset(coordSection, cStart, &off);CHKERRQ(ierr);
354*2e529ec8SStefano Zampini             ierr = PetscMemcpy(y,x+off,piece.nvertices*3*sizeof(PetscScalar));CHKERRQ(ierr);
355*2e529ec8SStefano Zampini           }
356*2e529ec8SStefano Zampini           for (c=cStart; c<cEnd; c++) {
357*2e529ec8SStefano Zampini             PetscInt dof, off;
358*2e529ec8SStefano Zampini 
359*2e529ec8SStefano Zampini             ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr);
360*2e529ec8SStefano Zampini             ierr = PetscSectionGetOffset(coordSection, c, &off);CHKERRQ(ierr);
361*2e529ec8SStefano Zampini             cnt += dof/dimEmbed;
362*2e529ec8SStefano Zampini           }
363*2e529ec8SStefano Zampini           if (cnt != piece.nvertices) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match %D != %D",cnt,piece.nvertices);
364*2e529ec8SStefano Zampini         }
365*2e529ec8SStefano Zampini         nsend = piece.nvertices*3;
366552f7358SJed Brown         ierr  = TransferWrite(viewer,fp,r,0,y ? y : x,buffer,nsend,PETSC_SCALAR,tag);CHKERRQ(ierr);
367552f7358SJed Brown         ierr  = PetscFree(y);CHKERRQ(ierr);
368552f7358SJed Brown         ierr  = VecRestoreArrayRead(coords,&x);CHKERRQ(ierr);
369552f7358SJed Brown       }
370552f7358SJed Brown       {                           /* Connectivity, offsets, types */
3713e869605SMatthew G. Knepley         PetscVTKInt  *connectivity = NULL, *offsets = NULL;
3723e869605SMatthew G. Knepley         PetscVTKType *types = NULL;
373*2e529ec8SStefano Zampini         ierr = DMPlexGetVTKConnectivity(dm,localized,&piece,&connectivity,&offsets,&types);CHKERRQ(ierr);
374552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,connectivity,buffer,piece.nconn,PETSC_INT32,tag);CHKERRQ(ierr);
375552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,offsets,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr);
376552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,types,buffer,piece.ncells,PETSC_UINT8,tag);CHKERRQ(ierr);
377552f7358SJed Brown         ierr = PetscFree3(connectivity,offsets,types);CHKERRQ(ierr);
378552f7358SJed Brown       }
379552f7358SJed Brown       {                         /* Owners (cell data) */
380552f7358SJed Brown         PetscVTKInt *owners;
381785e854fSJed Brown         ierr = PetscMalloc1(piece.ncells,&owners);CHKERRQ(ierr);
382552f7358SJed Brown         for (i=0; i<piece.ncells; i++) owners[i] = rank;
383552f7358SJed Brown         ierr = TransferWrite(viewer,fp,r,0,owners,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr);
384552f7358SJed Brown         ierr = PetscFree(owners);CHKERRQ(ierr);
385552f7358SJed Brown       }
386552f7358SJed Brown       /* Cell data */
387552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
388552f7358SJed Brown         Vec               X = (Vec)link->vec;
389552f7358SJed Brown         const PetscScalar *x;
390552f7358SJed Brown         PetscScalar       *y;
391552f7358SJed Brown         PetscInt          bs;
392552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
393552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
394552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
395785e854fSJed Brown         ierr = PetscMalloc1(piece.ncells,&y);CHKERRQ(ierr);
396552f7358SJed Brown         for (i=0; i<bs; i++) {
397552f7358SJed Brown           PetscInt cnt;
398552f7358SJed Brown           for (c=cStart,cnt=0; c<cEnd; c++) {
399640bce14SSatish Balay             PetscScalar *xpoint;
400552f7358SJed Brown             if (hasLabel) {     /* Ignore some cells */
401552f7358SJed Brown               PetscInt value;
402c58f1c22SToby Isaac               ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr);
403552f7358SJed Brown               if (value != 1) continue;
404552f7358SJed Brown             }
405552f7358SJed Brown             ierr     = DMPlexPointLocalRead(dm,c,x,&xpoint);CHKERRQ(ierr);
406552f7358SJed Brown             y[cnt++] = xpoint[i];
407552f7358SJed Brown           }
408552f7358SJed Brown           if (cnt != piece.ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
409552f7358SJed Brown           ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells,PETSC_SCALAR,tag);CHKERRQ(ierr);
410552f7358SJed Brown         }
411552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
412552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
413552f7358SJed Brown       }
414552f7358SJed Brown 
415552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
416552f7358SJed Brown         Vec               X = (Vec)link->vec;
417552f7358SJed Brown         const PetscScalar *x;
418552f7358SJed Brown         PetscScalar       *y;
419552f7358SJed Brown         PetscInt          bs;
420552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
421552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
422552f7358SJed Brown         ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
423785e854fSJed Brown         ierr = PetscMalloc1(piece.nvertices,&y);CHKERRQ(ierr);
424552f7358SJed Brown         for (i=0; i<bs; i++) {
425552f7358SJed Brown           PetscInt cnt;
426*2e529ec8SStefano Zampini           if (!localized) {
427552f7358SJed Brown             for (v=vStart,cnt=0; v<vEnd; v++) {
428640bce14SSatish Balay               PetscScalar *xpoint;
42937045cddSJed Brown               ierr     = DMPlexPointLocalRead(dm,v,x,&xpoint);CHKERRQ(ierr);
430552f7358SJed Brown               y[cnt++] = xpoint[i];
431552f7358SJed Brown             }
432*2e529ec8SStefano Zampini           } else {
433*2e529ec8SStefano Zampini             for (c=cStart,cnt=0; c<cEnd; c++) {
434*2e529ec8SStefano Zampini               PetscInt *closure = NULL;
435*2e529ec8SStefano Zampini               PetscInt  closureSize, off;
436*2e529ec8SStefano Zampini 
437*2e529ec8SStefano Zampini               ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
438*2e529ec8SStefano Zampini               for (v = 0, off = 0; v < closureSize*2; v += 2) {
439*2e529ec8SStefano Zampini                 if ((closure[v] >= vStart) && (closure[v] < vEnd)) {
440*2e529ec8SStefano Zampini                   PetscScalar *xpoint;
441*2e529ec8SStefano Zampini 
442*2e529ec8SStefano Zampini                   ierr = DMPlexPointLocalRead(dm,v,x,&xpoint);CHKERRQ(ierr);
443*2e529ec8SStefano Zampini                   y[cnt + off++] = xpoint[i];
444*2e529ec8SStefano Zampini                 }
445*2e529ec8SStefano Zampini               }
446*2e529ec8SStefano Zampini               ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
447*2e529ec8SStefano Zampini             }
448*2e529ec8SStefano Zampini           }
449552f7358SJed Brown           if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match");
450552f7358SJed Brown           ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr);
451552f7358SJed Brown         }
452552f7358SJed Brown         ierr = PetscFree(y);CHKERRQ(ierr);
453552f7358SJed Brown         ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
454552f7358SJed Brown       }
455552f7358SJed Brown     } else if (!rank) {
4560298fd71SBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,PETSC_SCALAR,tag);CHKERRQ(ierr); /* positions */
4570298fd71SBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nconn,PETSC_INT32,tag);CHKERRQ(ierr); /* connectivity */
4580298fd71SBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* offsets */
4590298fd71SBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_UINT8,tag);CHKERRQ(ierr); /* types */
4600298fd71SBarry Smith       ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* owner rank (cells) */
461552f7358SJed Brown       /* all cell data */
462552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
463552f7358SJed Brown         PetscInt bs;
464552f7358SJed Brown         if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue;
465fd68c46aSJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr);
466552f7358SJed Brown         for (i=0; i<bs; i++) {
4670298fd71SBarry Smith           ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_SCALAR,tag);CHKERRQ(ierr);
468552f7358SJed Brown         }
469552f7358SJed Brown       }
470552f7358SJed Brown       /* all point data */
471552f7358SJed Brown       for (link=vtk->link; link; link=link->next) {
472552f7358SJed Brown         PetscInt bs;
473552f7358SJed Brown         if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue;
474552f7358SJed Brown         ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr);
475552f7358SJed Brown         for (i=0; i<bs; i++) {
4760298fd71SBarry Smith           ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr);
477552f7358SJed Brown         }
478552f7358SJed Brown       }
479552f7358SJed Brown     }
480552f7358SJed Brown   }
481552f7358SJed Brown   ierr = PetscFree(gpiece);CHKERRQ(ierr);
482552f7358SJed Brown   ierr = PetscFree(buffer);CHKERRQ(ierr);
483552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"\n  </AppendedData>\n");CHKERRQ(ierr);
484552f7358SJed Brown   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
4855e32de72SMatthew G. Knepley   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
486552f7358SJed Brown   PetscFunctionReturn(0);
487552f7358SJed Brown }
488