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 static PetscErrorCode TransferWrite(PetscViewer viewer,FILE *fp,PetscMPIInt srank,PetscMPIInt root,const void *send,void *recv,PetscMPIInt count,MPI_Datatype mpidatatype,PetscMPIInt tag) 19 { 20 PetscMPIInt rank; 21 PetscErrorCode ierr; 22 MPI_Comm comm; 23 24 PetscFunctionBegin; 25 ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); 26 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 27 28 if (rank == srank && rank != root) { 29 ierr = MPI_Send((void*)send,count,mpidatatype,root,tag,comm);CHKERRQ(ierr); 30 } else if (rank == root) { 31 const void *buffer; 32 if (root == srank) { /* self */ 33 buffer = send; 34 } else { 35 MPI_Status status; 36 PetscMPIInt nrecv; 37 ierr = MPI_Recv(recv,count,mpidatatype,srank,tag,comm,&status);CHKERRQ(ierr); 38 ierr = MPI_Get_count(&status,mpidatatype,&nrecv);CHKERRQ(ierr); 39 if (count != nrecv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch"); 40 buffer = recv; 41 } 42 ierr = PetscViewerVTKFWrite(viewer,fp,buffer,count,mpidatatype);CHKERRQ(ierr); 43 } 44 PetscFunctionReturn(0); 45 } 46 47 static PetscErrorCode DMPlexGetVTKConnectivity(DM dm, PetscBool localized, PieceInfo *piece,PetscVTKInt **oconn,PetscVTKInt **ooffsets,PetscVTKType **otypes) 48 { 49 PetscErrorCode ierr; 50 PetscSection coordSection; 51 PetscVTKInt *conn,*offsets; 52 PetscVTKType *types; 53 PetscInt dim,vStart,vEnd,cStart,cEnd,pStart,pEnd,cellHeight,numLabelCells,hasLabel,c,v,countcell,countconn; 54 55 PetscFunctionBegin; 56 ierr = PetscMalloc3(piece->nconn,&conn,piece->ncells,&offsets,piece->ncells,&types);CHKERRQ(ierr); 57 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 58 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 59 ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 60 ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 61 ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 62 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 63 ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 64 hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 65 66 countcell = 0; 67 countconn = 0; 68 for (c = cStart; c < cEnd; ++c) { 69 PetscInt nverts,dof = 0,celltype,startoffset,nC=0; 70 71 if (hasLabel) { 72 PetscInt value; 73 74 ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 75 if (value != 1) continue; 76 } 77 startoffset = countconn; 78 if (localized) { 79 ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr); 80 } 81 if (!dof) { 82 PetscInt *closure = NULL; 83 PetscInt closureSize; 84 85 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 86 for (v = 0; v < closureSize*2; v += 2) { 87 if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 88 if (!localized) conn[countconn++] = closure[v] - vStart; 89 else conn[countconn++] = startoffset + nC; 90 ++nC; 91 } 92 } 93 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 94 } else { 95 for (nC = 0; nC < dof/dim; nC++) conn[countconn++] = startoffset + nC; 96 } 97 ierr = DMPlexInvertCell(dim, nC, &conn[countconn-nC]);CHKERRQ(ierr); 98 99 offsets[countcell] = countconn; 100 101 nverts = countconn - startoffset; 102 ierr = DMPlexVTKGetCellType_Internal(dm,dim,nverts,&celltype);CHKERRQ(ierr); 103 104 types[countcell] = celltype; 105 countcell++; 106 } 107 if (countcell != piece->ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count"); 108 if (countconn != piece->nconn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count"); 109 *oconn = conn; 110 *ooffsets = offsets; 111 *otypes = types; 112 PetscFunctionReturn(0); 113 } 114 115 /* 116 Write all fields that have been provided to the viewer 117 Multi-block XML format with binary appended data. 118 */ 119 PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer) 120 { 121 MPI_Comm comm; 122 PetscSection coordSection; 123 PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data; 124 PetscViewerVTKObjectLink link; 125 FILE *fp; 126 PetscMPIInt rank,size,tag; 127 PetscErrorCode ierr; 128 PetscInt dimEmbed,cellHeight,cStart,cEnd,vStart,vEnd,numLabelCells,hasLabel,c,v,r,i; 129 PetscBool localized; 130 PieceInfo piece,*gpiece = NULL; 131 void *buffer = NULL; 132 const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian"; 133 134 PetscFunctionBegin; 135 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 136 #if defined(PETSC_USE_COMPLEX) 137 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Complex values not supported"); 138 #endif 139 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 140 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 141 tag = ((PetscObject)viewer)->tag; 142 143 ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr); 144 ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr); 145 ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"%s\">\n", byte_order);CHKERRQ(ierr); 146 ierr = PetscFPrintf(comm,fp," <UnstructuredGrid>\n");CHKERRQ(ierr); 147 148 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 149 ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 150 ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 151 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 152 ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 153 ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr); 154 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 155 156 hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 157 piece.nvertices = 0; 158 piece.ncells = 0; 159 piece.nconn = 0; 160 if (!localized) piece.nvertices = vEnd - vStart; 161 for (c = cStart; c < cEnd; ++c) { 162 PetscInt dof = 0; 163 if (hasLabel) { 164 PetscInt value; 165 166 ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 167 if (value != 1) continue; 168 } 169 if (localized) { 170 ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr); 171 } 172 if (!dof) { 173 PetscInt *closure = NULL; 174 PetscInt closureSize; 175 176 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 177 for (v = 0; v < closureSize*2; v += 2) { 178 if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 179 piece.nconn++; 180 if (localized) piece.nvertices++; 181 } 182 } 183 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 184 } else { 185 piece.nvertices += dof/dimEmbed; 186 piece.nconn += dof/dimEmbed; 187 } 188 piece.ncells++; 189 } 190 if (!rank) {ierr = PetscMalloc1(size,&gpiece);CHKERRQ(ierr);} 191 ierr = MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm);CHKERRQ(ierr); 192 193 /* 194 * Write file header 195 */ 196 if (!rank) { 197 PetscInt boffset = 0; 198 199 for (r=0; r<size; r++) { 200 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells);CHKERRQ(ierr); 201 /* Coordinate positions */ 202 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Points>\n");CHKERRQ(ierr); 203 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr); 204 boffset += gpiece[r].nvertices*3*sizeof(PetscScalar) + sizeof(int); 205 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Points>\n");CHKERRQ(ierr); 206 /* Cell connectivity */ 207 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Cells>\n");CHKERRQ(ierr); 208 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 209 boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int); 210 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"offsets\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 211 boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int); 212 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"UInt8\" Name=\"types\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 213 boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int); 214 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Cells>\n");CHKERRQ(ierr); 215 216 /* 217 * Cell Data headers 218 */ 219 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <CellData>\n");CHKERRQ(ierr); 220 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 221 boffset += gpiece[r].ncells*sizeof(int) + sizeof(int); 222 /* all the vectors */ 223 for (link=vtk->link; link; link=link->next) { 224 Vec X = (Vec)link->vec; 225 PetscInt bs,nfields,field; 226 const char *vecname = ""; 227 if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 228 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. */ 229 ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr); 230 } 231 ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr); 232 ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr); 233 for (field=0,i=0; field<(nfields?nfields:1); field++) { 234 PetscInt fbs,j; 235 PetscFV fv = NULL; 236 PetscObject f; 237 PetscClassId fClass; 238 const char *fieldname = NULL; 239 char buf[256]; 240 if (nfields) { /* We have user-defined fields/components */ 241 ierr = PetscSectionGetFieldDof(dm->defaultSection,cStart,field,&fbs);CHKERRQ(ierr); 242 ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr); 243 } else fbs = bs; /* Say we have one field with 'bs' components */ 244 ierr = DMGetField(dm,field,NULL,&f);CHKERRQ(ierr); 245 ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr); 246 if (fClass == PETSCFV_CLASSID) { 247 fv = (PetscFV) f; 248 } 249 if (!fieldname) { 250 ierr = PetscSNPrintf(buf,sizeof(buf),"CellField%D",field);CHKERRQ(ierr); 251 fieldname = buf; 252 } 253 for (j=0; j<fbs; j++) { 254 const char *compName = NULL; 255 if (fv) { 256 ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr); 257 } 258 if (compName) { 259 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); 260 } 261 else { 262 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); 263 } 264 boffset += gpiece[r].ncells*sizeof(PetscScalar) + sizeof(int); 265 i++; 266 } 267 } 268 if (i != bs) SETERRQ2(comm,PETSC_ERR_PLIB,"Total number of field components %D != block size %D",i,bs); 269 } 270 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </CellData>\n");CHKERRQ(ierr); 271 272 /* 273 * Point Data headers 274 */ 275 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <PointData>\n");CHKERRQ(ierr); 276 for (link=vtk->link; link; link=link->next) { 277 Vec X = (Vec)link->vec; 278 PetscInt bs,nfields,field; 279 const char *vecname = ""; 280 if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 281 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. */ 282 ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr); 283 } 284 ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 285 ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr); 286 for (field=0,i=0; field<(nfields?nfields:1); field++) { 287 PetscInt fbs,j; 288 const char *fieldname = NULL; 289 char buf[256]; 290 if (nfields) { /* We have user-defined fields/components */ 291 ierr = PetscSectionGetFieldDof(dm->defaultSection,vStart,field,&fbs);CHKERRQ(ierr); 292 ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr); 293 } else fbs = bs; /* Say we have one field with 'bs' components */ 294 if (!fieldname) { 295 ierr = PetscSNPrintf(buf,sizeof(buf),"PointField%D",field);CHKERRQ(ierr); 296 fieldname = buf; 297 } 298 for (j=0; j<fbs; j++) { 299 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); 300 boffset += gpiece[r].nvertices*sizeof(PetscScalar) + sizeof(int); 301 } 302 } 303 } 304 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </PointData>\n");CHKERRQ(ierr); 305 306 ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Piece>\n");CHKERRQ(ierr); 307 } 308 } 309 310 ierr = PetscFPrintf(comm,fp," </UnstructuredGrid>\n");CHKERRQ(ierr); 311 ierr = PetscFPrintf(comm,fp," <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr); 312 ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr); 313 314 if (!rank) { 315 PetscInt maxsize = 0; 316 for (r=0; r<size; r++) { 317 maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscScalar))); 318 maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*sizeof(PetscScalar))); 319 maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt))); 320 } 321 ierr = PetscMalloc(maxsize,&buffer);CHKERRQ(ierr); 322 } 323 for (r=0; r<size; r++) { 324 if (r == rank) { 325 PetscInt nsend; 326 { /* Position */ 327 const PetscScalar *x; 328 PetscScalar *y = NULL; 329 Vec coords; 330 331 ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr); 332 ierr = VecGetArrayRead(coords,&x);CHKERRQ(ierr); 333 if (dimEmbed != 3 || localized) { 334 ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr); 335 if (localized) { 336 PetscInt cnt; 337 for (c=cStart,cnt=0; c<cEnd; c++) { 338 PetscInt off, dof; 339 340 ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr); 341 if (!dof) { 342 PetscInt *closure = NULL; 343 PetscInt closureSize; 344 345 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 346 for (v = 0; v < closureSize*2; v += 2) { 347 if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 348 ierr = PetscSectionGetOffset(coordSection, closure[v], &off);CHKERRQ(ierr); 349 if (dimEmbed != 3) { 350 y[cnt*3+0] = x[off+0]; 351 y[cnt*3+1] = (dimEmbed > 1) ? x[off+1] : 0.0; 352 y[cnt*3+2] = 0.0; 353 } else { 354 y[cnt*3+0] = x[off+0]; 355 y[cnt*3+1] = x[off+1]; 356 y[cnt*3+2] = x[off+2]; 357 } 358 cnt++; 359 } 360 } 361 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 362 } else { 363 ierr = PetscSectionGetOffset(coordSection, c, &off);CHKERRQ(ierr); 364 if (dimEmbed != 3) { 365 for (i=0; i<dof/dimEmbed; i++) { 366 y[cnt*3+0] = x[off + i*dimEmbed + 0]; 367 y[cnt*3+1] = (dimEmbed > 1) ? x[off + i*dimEmbed + 1] : 0.0; 368 y[cnt*3+2] = 0.0; 369 cnt++; 370 } 371 } else { 372 for (i=0; i<dof; i ++) { 373 y[cnt*3+i] = x[off + i]; 374 } 375 cnt += dof/dimEmbed; 376 } 377 } 378 } 379 if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 380 } else { 381 for (i=0; i<piece.nvertices; i++) { 382 y[i*3+0] = x[i*dimEmbed+0]; 383 y[i*3+1] = (dimEmbed > 1) ? x[i*dimEmbed+1] : 0; 384 y[i*3+2] = 0.0; 385 } 386 } 387 } 388 nsend = piece.nvertices*3; 389 ierr = TransferWrite(viewer,fp,r,0,y ? y : x,buffer,nsend,MPIU_SCALAR,tag);CHKERRQ(ierr); 390 ierr = PetscFree(y);CHKERRQ(ierr); 391 ierr = VecRestoreArrayRead(coords,&x);CHKERRQ(ierr); 392 } 393 { /* Connectivity, offsets, types */ 394 PetscVTKInt *connectivity = NULL, *offsets = NULL; 395 PetscVTKType *types = NULL; 396 ierr = DMPlexGetVTKConnectivity(dm,localized,&piece,&connectivity,&offsets,&types);CHKERRQ(ierr); 397 ierr = TransferWrite(viewer,fp,r,0,connectivity,buffer,piece.nconn,MPI_INT,tag);CHKERRQ(ierr); 398 ierr = TransferWrite(viewer,fp,r,0,offsets,buffer,piece.ncells,MPI_INT,tag);CHKERRQ(ierr); 399 ierr = TransferWrite(viewer,fp,r,0,types,buffer,piece.ncells,MPI_CHAR,tag);CHKERRQ(ierr); 400 ierr = PetscFree3(connectivity,offsets,types);CHKERRQ(ierr); 401 } 402 { /* Owners (cell data) */ 403 PetscVTKInt *owners; 404 ierr = PetscMalloc1(piece.ncells,&owners);CHKERRQ(ierr); 405 for (i=0; i<piece.ncells; i++) owners[i] = rank; 406 ierr = TransferWrite(viewer,fp,r,0,owners,buffer,piece.ncells,MPI_INT,tag);CHKERRQ(ierr); 407 ierr = PetscFree(owners);CHKERRQ(ierr); 408 } 409 /* Cell data */ 410 for (link=vtk->link; link; link=link->next) { 411 Vec X = (Vec)link->vec; 412 const PetscScalar *x; 413 PetscScalar *y; 414 PetscInt bs; 415 if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 416 ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr); 417 ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 418 ierr = PetscMalloc1(piece.ncells,&y);CHKERRQ(ierr); 419 for (i=0; i<bs; i++) { 420 PetscInt cnt; 421 for (c=cStart,cnt=0; c<cEnd; c++) { 422 PetscScalar *xpoint; 423 if (hasLabel) { /* Ignore some cells */ 424 PetscInt value; 425 ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 426 if (value != 1) continue; 427 } 428 ierr = DMPlexPointLocalRead(dm,c,x,&xpoint);CHKERRQ(ierr); 429 y[cnt++] = xpoint[i]; 430 } 431 if (cnt != piece.ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 432 ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells,MPIU_SCALAR,tag);CHKERRQ(ierr); 433 } 434 ierr = PetscFree(y);CHKERRQ(ierr); 435 ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 436 } 437 438 for (link=vtk->link; link; link=link->next) { 439 Vec X = (Vec)link->vec; 440 const PetscScalar *x; 441 PetscScalar *y; 442 PetscInt bs; 443 if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 444 ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 445 ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 446 ierr = PetscMalloc1(piece.nvertices,&y);CHKERRQ(ierr); 447 for (i=0; i<bs; i++) { 448 PetscInt cnt; 449 if (!localized) { 450 for (v=vStart,cnt=0; v<vEnd; v++) { 451 PetscScalar *xpoint; 452 ierr = DMPlexPointLocalRead(dm,v,x,&xpoint);CHKERRQ(ierr); 453 y[cnt++] = xpoint[i]; 454 } 455 } else { 456 for (c=cStart,cnt=0; c<cEnd; c++) { 457 PetscInt *closure = NULL; 458 PetscInt closureSize, off; 459 460 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 461 for (v = 0, off = 0; v < closureSize*2; v += 2) { 462 if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 463 PetscScalar *xpoint; 464 465 ierr = DMPlexPointLocalRead(dm,closure[v],x,&xpoint);CHKERRQ(ierr); 466 y[cnt + off++] = xpoint[i]; 467 } 468 } 469 cnt += off; 470 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 471 } 472 } 473 if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 474 ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices,MPIU_SCALAR,tag);CHKERRQ(ierr); 475 } 476 ierr = PetscFree(y);CHKERRQ(ierr); 477 ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 478 } 479 } else if (!rank) { 480 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,MPIU_SCALAR,tag);CHKERRQ(ierr); /* positions */ 481 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nconn,MPI_INT,tag);CHKERRQ(ierr); /* connectivity */ 482 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag);CHKERRQ(ierr); /* offsets */ 483 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_CHAR,tag);CHKERRQ(ierr); /* types */ 484 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag);CHKERRQ(ierr); /* owner rank (cells) */ 485 /* all cell data */ 486 for (link=vtk->link; link; link=link->next) { 487 PetscInt bs; 488 if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 489 ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr); 490 for (i=0; i<bs; i++) { 491 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPIU_SCALAR,tag);CHKERRQ(ierr); 492 } 493 } 494 /* all point data */ 495 for (link=vtk->link; link; link=link->next) { 496 PetscInt bs; 497 if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 498 ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 499 for (i=0; i<bs; i++) { 500 ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices,MPIU_SCALAR,tag);CHKERRQ(ierr); 501 } 502 } 503 } 504 } 505 ierr = PetscFree(gpiece);CHKERRQ(ierr); 506 ierr = PetscFree(buffer);CHKERRQ(ierr); 507 ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr); 508 ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr); 509 ierr = PetscFClose(comm,fp);CHKERRQ(ierr); 510 PetscFunctionReturn(0); 511 } 512