1 #include <petsc-private/dmpleximpl.h> 2 #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h> 3 4 typedef struct { 5 PetscInt nvertices; 6 PetscInt ncells; 7 PetscInt nconn; /* number of entries in cell->vertex connectivity array */ 8 } PieceInfo; 9 10 #if defined(PETSC_USE_REAL_SINGLE) 11 static const char precision[] = "Float32"; 12 #elif defined(PETSC_USE_REAL_DOUBLE) 13 static const char precision[] = "Float64"; 14 #else 15 static const char precision[] = "UnknownPrecision"; 16 #endif 17 18 #undef __FUNCT__ 19 #define __FUNCT__ "TransferWrite" 20 static PetscErrorCode TransferWrite(PetscViewer viewer,FILE *fp,PetscMPIInt srank,PetscMPIInt root,const void *send,void *recv,PetscMPIInt count,PetscDataType datatype,PetscMPIInt tag) 21 { 22 PetscMPIInt rank; 23 PetscErrorCode ierr; 24 MPI_Comm comm; 25 MPI_Datatype mpidatatype; 26 27 PetscFunctionBegin; 28 ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); 29 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 30 ierr = PetscDataTypeToMPIDataType(datatype,&mpidatatype);CHKERRQ(ierr); 31 32 if (rank == srank && rank != root) { 33 ierr = MPI_Send((void*)send,count,mpidatatype,root,tag,comm);CHKERRQ(ierr); 34 } else if (rank == root) { 35 const void *buffer; 36 if (root == srank) { /* self */ 37 buffer = send; 38 } else { 39 MPI_Status status; 40 PetscMPIInt nrecv; 41 ierr = MPI_Recv(recv,count,mpidatatype,srank,tag,comm,&status);CHKERRQ(ierr); 42 ierr = MPI_Get_count(&status,mpidatatype,&nrecv);CHKERRQ(ierr); 43 if (count != nrecv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch"); 44 buffer = recv; 45 } 46 ierr = PetscViewerVTKFWrite(viewer,fp,buffer,count,datatype);CHKERRQ(ierr); 47 } 48 PetscFunctionReturn(0); 49 } 50 51 #undef __FUNCT__ 52 #define __FUNCT__ "DMPlexGetVTKConnectivity" 53 static PetscErrorCode DMPlexGetVTKConnectivity(DM dm,PieceInfo *piece,PetscVTKInt **oconn,PetscVTKInt **ooffsets,PetscVTKType **otypes) 54 { 55 PetscErrorCode ierr; 56 PetscVTKInt *conn,*offsets; 57 PetscVTKType *types; 58 PetscInt dim,vStart,vEnd,cStart,cEnd,pStart,pEnd,cellHeight,cMax,numLabelCells,hasLabel,c,v,countcell,countconn; 59 60 PetscFunctionBegin; 61 ierr = PetscMalloc3(piece->nconn,PetscVTKInt,&conn,piece->ncells,PetscVTKInt,&offsets,piece->ncells,PetscVTKType,&types);CHKERRQ(ierr); 62 63 ierr = DMPlexGetDimension(dm,&dim);CHKERRQ(ierr); 64 ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 65 ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 66 ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 67 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 68 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 69 if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 70 ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 71 hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 72 73 countcell = 0; 74 countconn = 0; 75 for (c = cStart; c < cEnd; ++c) { 76 PetscInt *closure = NULL; 77 PetscInt closureSize,nverts,celltype,startoffset; 78 79 if (hasLabel) { 80 PetscInt value; 81 82 ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 83 if (value != 1) continue; 84 } 85 startoffset = countconn; 86 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 87 for (v = 0; v < closureSize*2; v += 2) { 88 if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 89 conn[countconn++] = closure[v] - vStart; 90 } 91 } 92 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 93 94 offsets[countcell] = countconn; 95 96 nverts = countconn - startoffset; 97 ierr = DMPlexVTKGetCellType(dm,dim,nverts,&celltype);CHKERRQ(ierr); 98 99 types[countcell] = celltype; 100 countcell++; 101 } 102 if (countcell != piece->ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count"); 103 if (countconn != piece->nconn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count"); 104 *oconn = conn; 105 *ooffsets = offsets; 106 *otypes = types; 107 PetscFunctionReturn(0); 108 } 109 110 #undef __FUNCT__ 111 #define __FUNCT__ "DMPlexVTKWriteAll_VTU" 112 /* 113 Write all fields that have been provided to the viewer 114 Multi-block XML format with binary appended data. 115 */ 116 PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer) 117 { 118 MPI_Comm comm; 119 PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data; 120 PetscViewerVTKObjectLink link; 121 FILE *fp; 122 PetscMPIInt rank,size,tag; 123 PetscErrorCode ierr; 124 PetscInt dim,cellHeight,cStart,cEnd,vStart,vEnd,cMax,numLabelCells,hasLabel,c,v,r,i; 125 PieceInfo piece,*gpiece = NULL; 126 void *buffer = NULL; 127 128 PetscFunctionBegin; 129 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 130 #if defined(PETSC_USE_COMPLEX) 131 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Complex values not supported"); 132 #endif 133 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 134 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 135 tag = ((PetscObject)viewer)->tag; 136 137 ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr); 138 ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr); 139 #if defined(PETSC_WORDS_BIGENDIAN) 140 ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr); 141 #else 142 ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr); 143 #endif 144 ierr = PetscFPrintf(comm,fp," <UnstructuredGrid>\n");CHKERRQ(ierr); 145 146 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 147 ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 148 ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 149 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 150 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 151 if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 152 ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 153 154 hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 155 piece.nvertices = vEnd - vStart; 156 piece.ncells = 0; 157 piece.nconn = 0; 158 for (c = cStart; c < cEnd; ++c) { 159 PetscInt *closure = NULL; 160 PetscInt closureSize; 161 162 if (hasLabel) { 163 PetscInt value; 164 165 ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 166 if (value != 1) continue; 167 } 168 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 169 for (v = 0; v < closureSize*2; v += 2) { 170 if ((closure[v] >= vStart) && (closure[v] < vEnd)) piece.nconn++; 171 } 172 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 173 piece.ncells++; 174 } 175 if (!rank) {ierr = PetscMalloc(size*sizeof(piece),&gpiece);CHKERRQ(ierr);} 176 ierr = MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm);CHKERRQ(ierr); 177 178 /* 179 * Write file header 180 */ 181 if (!rank) { 182 PetscInt boffset = 0; 183 184 for (r=0; r<size; r++) { 185 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells);CHKERRQ(ierr); 186 /* Coordinate positions */ 187 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Points>\n");CHKERRQ(ierr); 188 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr); 189 boffset += gpiece[r].nvertices*3*sizeof(PetscScalar) + sizeof(int); 190 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Points>\n");CHKERRQ(ierr); 191 /* Cell connectivity */ 192 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Cells>\n");CHKERRQ(ierr); 193 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 194 boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int); 195 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"offsets\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 196 boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int); 197 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"UInt8\" Name=\"types\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 198 boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int); 199 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Cells>\n");CHKERRQ(ierr); 200 201 /* 202 * Cell Data headers 203 */ 204 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <CellData>\n");CHKERRQ(ierr); 205 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 206 boffset += gpiece[r].ncells*sizeof(int) + sizeof(int); 207 /* all the vectors */ 208 for (link=vtk->link; link; link=link->next) { 209 Vec X = (Vec)link->vec; 210 PetscInt bs,nfields,field; 211 const char *vecname = ""; 212 if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 213 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. */ 214 ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr); 215 } 216 ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr); 217 ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr); 218 for (field=0,i=0; field<nfields; field++) { 219 PetscInt fbs,j; 220 const char *fieldname = NULL; 221 char buf[256]; 222 ierr = PetscSectionGetFieldDof(dm->defaultSection,cStart,field,&fbs);CHKERRQ(ierr); 223 ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr); 224 if (!fieldname) { 225 ierr = PetscSNPrintf(buf,sizeof(buf),"CellField%D",field);CHKERRQ(ierr); 226 fieldname = buf; 227 } 228 for (j=0; j<fbs; j++) { 229 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); 230 boffset += gpiece[r].ncells*sizeof(PetscScalar) + sizeof(int); 231 i++; 232 } 233 } 234 if (i != bs) SETERRQ2(comm,PETSC_ERR_PLIB,"Total number of field components %D != block size %D",i,bs); 235 } 236 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </CellData>\n");CHKERRQ(ierr); 237 238 /* 239 * Point Data headers 240 */ 241 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <PointData>\n");CHKERRQ(ierr); 242 for (link=vtk->link; link; link=link->next) { 243 Vec X = (Vec)link->vec; 244 PetscInt bs,nfields,field; 245 const char *vecname = ""; 246 if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 247 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. */ 248 ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr); 249 } 250 ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 251 ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr); 252 for (field=0,i=0; field<nfields; field++) { 253 PetscInt fbs,j; 254 const char *fieldname = NULL; 255 char buf[256]; 256 ierr = PetscSectionGetFieldDof(dm->defaultSection,vStart,field,&fbs);CHKERRQ(ierr); 257 ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr); 258 if (!fieldname) { 259 ierr = PetscSNPrintf(buf,sizeof(buf),"PointField%D",field);CHKERRQ(ierr); 260 fieldname = buf; 261 } 262 for (j=0; j<fbs; j++) { 263 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); 264 boffset += gpiece[r].nvertices*sizeof(PetscScalar) + sizeof(int); 265 } 266 } 267 } 268 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </PointData>\n");CHKERRQ(ierr); 269 270 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Piece>\n");CHKERRQ(ierr); 271 } 272 } 273 274 ierr = PetscFPrintf(comm,fp," </UnstructuredGrid>\n");CHKERRQ(ierr); 275 ierr = PetscFPrintf(comm,fp," <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr); 276 ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr); 277 278 if (!rank) { 279 PetscInt maxsize = 0; 280 for (r=0; r<size; r++) { 281 maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscScalar))); 282 maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*sizeof(PetscScalar))); 283 maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt))); 284 } 285 ierr = PetscMalloc(maxsize,&buffer);CHKERRQ(ierr); 286 } 287 for (r=0; r<size; r++) { 288 if (r == rank) { 289 PetscInt nsend; 290 { /* Position */ 291 const PetscScalar *x; 292 PetscScalar *y = NULL; 293 Vec coords; 294 nsend = piece.nvertices*3; 295 ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr); 296 ierr = VecGetArrayRead(coords,&x);CHKERRQ(ierr); 297 if (dim != 3) { 298 ierr = PetscMalloc(piece.nvertices*3*sizeof(PetscScalar),&y);CHKERRQ(ierr); 299 for (i=0; i<piece.nvertices; i++) { 300 y[i*3+0] = x[i*dim+0]; 301 y[i*3+1] = (dim > 1) ? x[i*dim+1] : 0; 302 y[i*3+2] = 0; 303 } 304 } 305 ierr = TransferWrite(viewer,fp,r,0,y ? y : x,buffer,nsend,PETSC_SCALAR,tag);CHKERRQ(ierr); 306 ierr = PetscFree(y);CHKERRQ(ierr); 307 ierr = VecRestoreArrayRead(coords,&x);CHKERRQ(ierr); 308 } 309 { /* Connectivity, offsets, types */ 310 PetscVTKInt *connectivity = NULL,*offsets; 311 PetscVTKType *types; 312 ierr = DMPlexGetVTKConnectivity(dm,&piece,&connectivity,&offsets,&types);CHKERRQ(ierr); 313 ierr = TransferWrite(viewer,fp,r,0,connectivity,buffer,piece.nconn,PETSC_INT32,tag);CHKERRQ(ierr); 314 ierr = TransferWrite(viewer,fp,r,0,offsets,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr); 315 ierr = TransferWrite(viewer,fp,r,0,types,buffer,piece.ncells,PETSC_UINT8,tag);CHKERRQ(ierr); 316 ierr = PetscFree3(connectivity,offsets,types);CHKERRQ(ierr); 317 } 318 { /* Owners (cell data) */ 319 PetscVTKInt *owners; 320 ierr = PetscMalloc(piece.ncells*sizeof(PetscVTKInt),&owners);CHKERRQ(ierr); 321 for (i=0; i<piece.ncells; i++) owners[i] = rank; 322 ierr = TransferWrite(viewer,fp,r,0,owners,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr); 323 ierr = PetscFree(owners);CHKERRQ(ierr); 324 } 325 /* Cell data */ 326 for (link=vtk->link; link; link=link->next) { 327 Vec X = (Vec)link->vec; 328 const PetscScalar *x; 329 PetscScalar *y; 330 PetscInt bs; 331 if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 332 ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr); 333 ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 334 ierr = PetscMalloc(piece.ncells*sizeof(PetscScalar),&y);CHKERRQ(ierr); 335 for (i=0; i<bs; i++) { 336 PetscInt cnt; 337 for (c=cStart,cnt=0; c<cEnd; c++) { 338 const PetscScalar *xpoint; 339 if (hasLabel) { /* Ignore some cells */ 340 PetscInt value; 341 ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 342 if (value != 1) continue; 343 } 344 ierr = DMPlexPointLocalRead(dm,c,x,&xpoint);CHKERRQ(ierr); 345 y[cnt++] = xpoint[i]; 346 } 347 if (cnt != piece.ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 348 ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells,PETSC_SCALAR,tag);CHKERRQ(ierr); 349 } 350 ierr = PetscFree(y);CHKERRQ(ierr); 351 ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 352 } 353 354 for (link=vtk->link; link; link=link->next) { 355 Vec X = (Vec)link->vec; 356 const PetscScalar *x; 357 PetscScalar *y; 358 PetscInt bs; 359 if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 360 ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 361 ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 362 ierr = PetscMalloc(piece.nvertices*sizeof(PetscScalar),&y);CHKERRQ(ierr); 363 for (i=0; i<bs; i++) { 364 PetscInt cnt; 365 for (v=vStart,cnt=0; v<vEnd; v++) { 366 const PetscScalar *xpoint; 367 ierr = DMPlexPointLocalRead(dm,v,x,&xpoint);CHKERRQ(ierr); 368 y[cnt++] = xpoint[i]; 369 } 370 if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 371 ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr); 372 } 373 ierr = PetscFree(y);CHKERRQ(ierr); 374 ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 375 } 376 } else if (!rank) { 377 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,PETSC_SCALAR,tag);CHKERRQ(ierr); /* positions */ 378 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nconn,PETSC_INT32,tag);CHKERRQ(ierr); /* connectivity */ 379 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* offsets */ 380 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_UINT8,tag);CHKERRQ(ierr); /* types */ 381 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* owner rank (cells) */ 382 /* all cell data */ 383 for (link=vtk->link; link; link=link->next) { 384 PetscInt bs; 385 if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 386 ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 387 for (i=0; i<bs; i++) { 388 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_SCALAR,tag);CHKERRQ(ierr); 389 } 390 } 391 /* all point data */ 392 for (link=vtk->link; link; link=link->next) { 393 PetscInt bs; 394 if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 395 ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 396 for (i=0; i<bs; i++) { 397 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr); 398 } 399 } 400 } 401 } 402 ierr = PetscFree(gpiece);CHKERRQ(ierr); 403 ierr = PetscFree(buffer);CHKERRQ(ierr); 404 ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr); 405 ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr); 406 PetscFunctionReturn(0); 407 } 408