1552f7358SJed Brown #include <petsc-private/pleximpl.h> 2552f7358SJed Brown #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h> 3552f7358SJed Brown 4552f7358SJed Brown typedef struct { 5552f7358SJed Brown PetscInt nvertices; 6552f7358SJed Brown PetscInt ncells; 7552f7358SJed Brown PetscInt nconn; /* number of entries in cell->vertex connectivity array */ 8552f7358SJed Brown } PieceInfo; 9552f7358SJed Brown 10552f7358SJed Brown #if defined(PETSC_USE_REAL_SINGLE) 11552f7358SJed Brown static const char precision[] = "Float32"; 12552f7358SJed Brown #elif defined(PETSC_USE_REAL_DOUBLE) 13552f7358SJed Brown static const char precision[] = "Float64"; 14552f7358SJed Brown #else 15552f7358SJed Brown static const char precision[] = "UnknownPrecision"; 16552f7358SJed Brown #endif 17552f7358SJed Brown 18552f7358SJed Brown #undef __FUNCT__ 19552f7358SJed Brown #define __FUNCT__ "TransferWrite" 20552f7358SJed Brown static PetscErrorCode TransferWrite(PetscViewer viewer,FILE *fp,PetscMPIInt srank,PetscMPIInt root,const void *send,void *recv,PetscMPIInt count,PetscDataType datatype,PetscMPIInt tag) 21552f7358SJed Brown { 22552f7358SJed Brown PetscMPIInt rank; 23552f7358SJed Brown PetscErrorCode ierr; 24552f7358SJed Brown MPI_Comm comm = ((PetscObject)viewer)->comm; 25552f7358SJed Brown MPI_Datatype mpidatatype; 26552f7358SJed Brown 27552f7358SJed Brown PetscFunctionBegin; 28552f7358SJed Brown ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 29552f7358SJed Brown ierr = PetscDataTypeToMPIDataType(datatype,&mpidatatype);CHKERRQ(ierr); 30552f7358SJed Brown 31552f7358SJed Brown if (rank == srank && rank != root) { 32552f7358SJed Brown ierr = MPI_Send((void*)send,count,mpidatatype,root,tag,comm);CHKERRQ(ierr); 33552f7358SJed Brown } else if (rank == root) { 34552f7358SJed Brown const void *buffer; 35552f7358SJed Brown if (root == srank) { /* self */ 36552f7358SJed Brown buffer = send; 37552f7358SJed Brown } else { 38552f7358SJed Brown MPI_Status status; 39552f7358SJed Brown PetscMPIInt nrecv; 40552f7358SJed Brown ierr = MPI_Recv(recv,count,mpidatatype,srank,tag,comm,&status);CHKERRQ(ierr); 41552f7358SJed Brown ierr = MPI_Get_count(&status,mpidatatype,&nrecv);CHKERRQ(ierr); 42552f7358SJed Brown if (count != nrecv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch"); 43552f7358SJed Brown buffer = recv; 44552f7358SJed Brown } 45552f7358SJed Brown ierr = PetscViewerVTKFWrite(viewer,fp,buffer,count,datatype);CHKERRQ(ierr); 46552f7358SJed Brown } 47552f7358SJed Brown PetscFunctionReturn(0); 48552f7358SJed Brown } 49552f7358SJed Brown 50552f7358SJed Brown #undef __FUNCT__ 51552f7358SJed Brown #define __FUNCT__ "DMPlexGetVTKConnectivity" 52552f7358SJed Brown static PetscErrorCode DMPlexGetVTKConnectivity(DM dm,PieceInfo *piece,PetscVTKInt **oconn,PetscVTKInt **ooffsets,PetscVTKType **otypes) 53552f7358SJed Brown { 54552f7358SJed Brown PetscErrorCode ierr; 55552f7358SJed Brown PetscVTKInt *conn,*offsets; 56552f7358SJed Brown PetscVTKType *types; 57552f7358SJed Brown PetscInt dim,vStart,vEnd,cStart,cEnd,pStart,pEnd,cellHeight,cMax,numLabelCells,hasLabel,c,v,countcell,countconn; 58552f7358SJed Brown 59552f7358SJed Brown PetscFunctionBegin; 60552f7358SJed Brown ierr = PetscMalloc3(piece->nconn,PetscVTKInt,&conn,piece->ncells,PetscVTKInt,&offsets,piece->ncells,PetscVTKType,&types);CHKERRQ(ierr); 61552f7358SJed Brown 62552f7358SJed Brown ierr = DMPlexGetDimension(dm,&dim);CHKERRQ(ierr); 63552f7358SJed Brown ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 64552f7358SJed Brown ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 65552f7358SJed Brown ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 66552f7358SJed Brown ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 67770b213bSMatthew G Knepley ierr = DMPlexGetHybridBounds(dm, &cMax, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 68*8865f1eaSKarl Rupp if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 69552f7358SJed Brown ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 70552f7358SJed Brown hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 71552f7358SJed Brown 72552f7358SJed Brown countcell = 0; 73552f7358SJed Brown countconn = 0; 74552f7358SJed Brown for (c = cStart; c < cEnd; ++c) { 75552f7358SJed Brown PetscInt *closure = PETSC_NULL; 76552f7358SJed Brown PetscInt closureSize,nverts,celltype,startoffset; 77552f7358SJed Brown 78552f7358SJed Brown if (hasLabel) { 79552f7358SJed Brown PetscInt value; 80552f7358SJed Brown 81552f7358SJed Brown ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 82552f7358SJed Brown if (value != 1) continue; 83552f7358SJed Brown } 84552f7358SJed Brown startoffset = countconn; 85552f7358SJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 86552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 87552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 88552f7358SJed Brown conn[countconn++] = closure[v] - vStart; 89552f7358SJed Brown } 90552f7358SJed Brown } 91552f7358SJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 92*8865f1eaSKarl Rupp 93552f7358SJed Brown offsets[countcell] = countconn; 94*8865f1eaSKarl Rupp 95552f7358SJed Brown nverts = countconn - startoffset; 96552f7358SJed Brown ierr = DMPlexVTKGetCellType(dm,dim,nverts,&celltype);CHKERRQ(ierr); 97*8865f1eaSKarl Rupp 98552f7358SJed Brown types[countcell] = celltype; 99552f7358SJed Brown countcell++; 100552f7358SJed Brown } 101552f7358SJed Brown if (countcell != piece->ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count"); 102552f7358SJed Brown if (countconn != piece->nconn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count"); 103552f7358SJed Brown *oconn = conn; 104552f7358SJed Brown *ooffsets = offsets; 105552f7358SJed Brown *otypes = types; 106552f7358SJed Brown PetscFunctionReturn(0); 107552f7358SJed Brown } 108552f7358SJed Brown 109552f7358SJed Brown #undef __FUNCT__ 110552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll_VTU" 111552f7358SJed Brown /* 112552f7358SJed Brown Write all fields that have been provided to the viewer 113552f7358SJed Brown Multi-block XML format with binary appended data. 114552f7358SJed Brown */ 115552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer) 116552f7358SJed Brown { 117552f7358SJed Brown MPI_Comm comm = ((PetscObject)dm)->comm; 118552f7358SJed Brown PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data; 119552f7358SJed Brown PetscViewerVTKObjectLink link; 120552f7358SJed Brown FILE *fp; 121552f7358SJed Brown PetscMPIInt rank,size,tag; 122552f7358SJed Brown PetscErrorCode ierr; 123552f7358SJed Brown PetscInt dim,cellHeight,cStart,cEnd,vStart,vEnd,cMax,numLabelCells,hasLabel,c,v,r,i; 124552f7358SJed Brown PieceInfo piece,*gpiece = PETSC_NULL; 125552f7358SJed Brown void *buffer = PETSC_NULL; 126552f7358SJed Brown 127552f7358SJed Brown PetscFunctionBegin; 128552f7358SJed Brown #if defined(PETSC_USE_COMPLEX) 129552f7358SJed Brown SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Complex values not supported"); 130552f7358SJed Brown #endif 131552f7358SJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 132552f7358SJed Brown ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 133552f7358SJed Brown tag = ((PetscObject)viewer)->tag; 134552f7358SJed Brown 135552f7358SJed Brown ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr); 136552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr); 137519f805aSKarl Rupp #if defined(PETSC_WORDS_BIGENDIAN) 138552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr); 139552f7358SJed Brown #else 140552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr); 141552f7358SJed Brown #endif 142552f7358SJed Brown ierr = PetscFPrintf(comm,fp," <UnstructuredGrid>\n");CHKERRQ(ierr); 143552f7358SJed Brown 144552f7358SJed Brown ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 145552f7358SJed Brown ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 146552f7358SJed Brown ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 147552f7358SJed Brown ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 148770b213bSMatthew G Knepley ierr = DMPlexGetHybridBounds(dm, &cMax, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 149*8865f1eaSKarl Rupp if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 150552f7358SJed Brown ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 151*8865f1eaSKarl Rupp 152552f7358SJed Brown hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 153552f7358SJed Brown piece.nvertices = vEnd - vStart; 154552f7358SJed Brown piece.ncells = 0; 155552f7358SJed Brown piece.nconn = 0; 156552f7358SJed Brown for (c = cStart; c < cEnd; ++c) { 157552f7358SJed Brown PetscInt *closure = PETSC_NULL; 158552f7358SJed Brown PetscInt closureSize; 159552f7358SJed Brown 160552f7358SJed Brown if (hasLabel) { 161552f7358SJed Brown PetscInt value; 162552f7358SJed Brown 163552f7358SJed Brown ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 164552f7358SJed Brown if (value != 1) continue; 165552f7358SJed Brown } 166552f7358SJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 167552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 168552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) piece.nconn++; 169552f7358SJed Brown } 170552f7358SJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 171552f7358SJed Brown piece.ncells++; 172552f7358SJed Brown } 173552f7358SJed Brown if (!rank) {ierr = PetscMalloc(size*sizeof(piece),&gpiece);CHKERRQ(ierr);} 174552f7358SJed Brown ierr = MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm);CHKERRQ(ierr); 175552f7358SJed Brown 176552f7358SJed Brown /* 177552f7358SJed Brown * Write file header 178552f7358SJed Brown */ 179552f7358SJed Brown if (!rank) { 180552f7358SJed Brown PetscInt boffset = 0; 181552f7358SJed Brown 182552f7358SJed Brown for (r=0; r<size; r++) { 183552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells);CHKERRQ(ierr); 184552f7358SJed Brown /* Coordinate positions */ 185552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Points>\n");CHKERRQ(ierr); 186552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr); 187552f7358SJed Brown boffset += gpiece[r].nvertices*3*sizeof(PetscScalar) + sizeof(int); 188552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Points>\n");CHKERRQ(ierr); 189552f7358SJed Brown /* Cell connectivity */ 190552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Cells>\n");CHKERRQ(ierr); 191552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 192552f7358SJed Brown boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int); 193552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"offsets\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 194552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int); 195552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"UInt8\" Name=\"types\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 196552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int); 197552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Cells>\n");CHKERRQ(ierr); 198552f7358SJed Brown 199552f7358SJed Brown /* 200552f7358SJed Brown * Cell Data headers 201552f7358SJed Brown */ 202552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <CellData>\n");CHKERRQ(ierr); 203552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 204552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(int) + sizeof(int); 205552f7358SJed Brown /* all the vectors */ 206552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 207552f7358SJed Brown Vec X = (Vec)link->vec; 208552f7358SJed Brown PetscInt bs; 209552f7358SJed Brown const char *vecname = ""; 210552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 211552f7358SJed 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. */ 212552f7358SJed Brown ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr); 213552f7358SJed Brown } 214552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr); 215552f7358SJed Brown for (i=0; i<bs; i++) { 216552f7358SJed Brown char buf[256]; 217552f7358SJed Brown const char *fieldname = PETSC_NULL; 218552f7358SJed Brown /* ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr); */ 219552f7358SJed Brown if (!fieldname) { 220552f7358SJed Brown ierr = PetscSNPrintf(buf,sizeof(buf),"Unnamed%D",i);CHKERRQ(ierr); 221552f7358SJed Brown fieldname = buf; 222552f7358SJed Brown } 223552f7358SJed Brown ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr); 224552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(PetscScalar) + sizeof(int); 225552f7358SJed Brown } 226552f7358SJed Brown } 227552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </CellData>\n");CHKERRQ(ierr); 228552f7358SJed Brown 229552f7358SJed Brown /* 230552f7358SJed Brown * Point Data headers 231552f7358SJed Brown */ 232552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <PointData>\n");CHKERRQ(ierr); 233552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 234552f7358SJed Brown Vec X = (Vec)link->vec; 235552f7358SJed Brown PetscInt bs; 236552f7358SJed Brown const char *vecname = ""; 237552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 238552f7358SJed 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. */ 239552f7358SJed Brown ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr); 240552f7358SJed Brown } 241552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 242552f7358SJed Brown for (i=0; i<bs; i++) { 243552f7358SJed Brown char fieldname[256]; 244552f7358SJed Brown ierr = PetscSNPrintf(fieldname,sizeof(fieldname),"Point%D",i);CHKERRQ(ierr); 245552f7358SJed Brown ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr); 246552f7358SJed Brown boffset += gpiece[r].nvertices*sizeof(PetscScalar) + sizeof(int); 247552f7358SJed Brown } 248552f7358SJed Brown } 249552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </PointData>\n");CHKERRQ(ierr); 250552f7358SJed Brown 251552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Piece>\n");CHKERRQ(ierr); 252552f7358SJed Brown } 253552f7358SJed Brown } 254552f7358SJed Brown 255552f7358SJed Brown ierr = PetscFPrintf(comm,fp," </UnstructuredGrid>\n");CHKERRQ(ierr); 256552f7358SJed Brown ierr = PetscFPrintf(comm,fp," <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr); 257552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr); 258552f7358SJed Brown 259552f7358SJed Brown if (!rank) { 260552f7358SJed Brown PetscInt maxsize = 0; 261552f7358SJed Brown for (r=0; r<size; r++) { 262552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscScalar))); 263552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*sizeof(PetscScalar))); 264552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt))); 265552f7358SJed Brown } 266552f7358SJed Brown ierr = PetscMalloc(maxsize,&buffer);CHKERRQ(ierr); 267552f7358SJed Brown } 268552f7358SJed Brown for (r=0; r<size; r++) { 269552f7358SJed Brown if (r == rank) { 270552f7358SJed Brown PetscInt nsend; 271552f7358SJed Brown { /* Position */ 272552f7358SJed Brown const PetscScalar *x; 273552f7358SJed Brown PetscScalar *y = PETSC_NULL; 274552f7358SJed Brown Vec coords; 275552f7358SJed Brown nsend = piece.nvertices*3; 276552f7358SJed Brown ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr); 277552f7358SJed Brown ierr = VecGetArrayRead(coords,&x);CHKERRQ(ierr); 278552f7358SJed Brown if (dim != 3) { 279552f7358SJed Brown ierr = PetscMalloc(piece.nvertices*3*sizeof(PetscScalar),&y);CHKERRQ(ierr); 280552f7358SJed Brown for (i=0; i<piece.nvertices; i++) { 281552f7358SJed Brown y[i*3+0] = x[i*dim+0]; 282552f7358SJed Brown y[i*3+1] = (dim > 1) ? x[i*dim+1] : 0; 283552f7358SJed Brown y[i*3+2] = 0; 284552f7358SJed Brown } 285552f7358SJed Brown } 286552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,y ? y : x,buffer,nsend,PETSC_SCALAR,tag);CHKERRQ(ierr); 287552f7358SJed Brown ierr = PetscFree(y);CHKERRQ(ierr); 288552f7358SJed Brown ierr = VecRestoreArrayRead(coords,&x);CHKERRQ(ierr); 289552f7358SJed Brown } 290552f7358SJed Brown { /* Connectivity, offsets, types */ 2916a415b8fSMatthew G Knepley PetscVTKInt *connectivity = PETSC_NULL,*offsets; 292552f7358SJed Brown PetscVTKType *types; 293552f7358SJed Brown ierr = DMPlexGetVTKConnectivity(dm,&piece,&connectivity,&offsets,&types);CHKERRQ(ierr); 294552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,connectivity,buffer,piece.nconn,PETSC_INT32,tag);CHKERRQ(ierr); 295552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,offsets,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr); 296552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,types,buffer,piece.ncells,PETSC_UINT8,tag);CHKERRQ(ierr); 297552f7358SJed Brown ierr = PetscFree3(connectivity,offsets,types);CHKERRQ(ierr); 298552f7358SJed Brown } 299552f7358SJed Brown { /* Owners (cell data) */ 300552f7358SJed Brown PetscVTKInt *owners; 301552f7358SJed Brown ierr = PetscMalloc(piece.ncells*sizeof(PetscVTKInt),&owners);CHKERRQ(ierr); 302552f7358SJed Brown for (i=0; i<piece.ncells; i++) owners[i] = rank; 303552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,owners,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr); 304552f7358SJed Brown ierr = PetscFree(owners);CHKERRQ(ierr); 305552f7358SJed Brown } 306552f7358SJed Brown /* Cell data */ 307552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 308552f7358SJed Brown Vec X = (Vec)link->vec; 309552f7358SJed Brown const PetscScalar *x; 310552f7358SJed Brown PetscScalar *y; 311552f7358SJed Brown PetscInt bs; 312552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 313552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr); 314552f7358SJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 315552f7358SJed Brown ierr = PetscMalloc(piece.ncells*sizeof(PetscScalar),&y);CHKERRQ(ierr); 316552f7358SJed Brown for (i=0; i<bs; i++) { 317552f7358SJed Brown PetscInt cnt; 318552f7358SJed Brown for (c=cStart,cnt=0; c<cEnd; c++) { 319552f7358SJed Brown const PetscScalar *xpoint; 320552f7358SJed Brown if (hasLabel) { /* Ignore some cells */ 321552f7358SJed Brown PetscInt value; 322552f7358SJed Brown ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 323552f7358SJed Brown if (value != 1) continue; 324552f7358SJed Brown } 325552f7358SJed Brown ierr = DMPlexPointLocalRead(dm,c,x,&xpoint);CHKERRQ(ierr); 326552f7358SJed Brown y[cnt++] = xpoint[i]; 327552f7358SJed Brown } 328552f7358SJed Brown if (cnt != piece.ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 329552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells,PETSC_SCALAR,tag);CHKERRQ(ierr); 330552f7358SJed Brown } 331552f7358SJed Brown ierr = PetscFree(y);CHKERRQ(ierr); 332552f7358SJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 333552f7358SJed Brown } 334552f7358SJed Brown 335552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 336552f7358SJed Brown Vec X = (Vec)link->vec; 337552f7358SJed Brown const PetscScalar *x; 338552f7358SJed Brown PetscScalar *y; 339552f7358SJed Brown PetscInt bs; 340552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 341552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 342552f7358SJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 343552f7358SJed Brown ierr = PetscMalloc(piece.nvertices*sizeof(PetscScalar),&y);CHKERRQ(ierr); 344552f7358SJed Brown for (i=0; i<bs; i++) { 345552f7358SJed Brown PetscInt cnt; 346552f7358SJed Brown for (v=vStart,cnt=0; v<vEnd; v++) { 347552f7358SJed Brown const PetscScalar *xpoint; 348552f7358SJed Brown ierr = DMPlexPointLocalRead(dm,c,x,&xpoint);CHKERRQ(ierr); 349552f7358SJed Brown y[cnt++] = xpoint[i]; 350552f7358SJed Brown } 351552f7358SJed Brown if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 352552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr); 353552f7358SJed Brown } 354552f7358SJed Brown ierr = PetscFree(y);CHKERRQ(ierr); 355552f7358SJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 356552f7358SJed Brown } 357552f7358SJed Brown } else if (!rank) { 358552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,PETSC_NULL,buffer,gpiece[r].nvertices*3,PETSC_SCALAR,tag);CHKERRQ(ierr); /* positions */ 359552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,PETSC_NULL,buffer,gpiece[r].nconn,PETSC_INT32,tag);CHKERRQ(ierr); /* connectivity */ 360552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,PETSC_NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* offsets */ 361552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,PETSC_NULL,buffer,gpiece[r].ncells,PETSC_UINT8,tag);CHKERRQ(ierr); /* types */ 362552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,PETSC_NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* owner rank (cells) */ 363552f7358SJed Brown /* all cell data */ 364552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 365552f7358SJed Brown PetscInt bs; 366552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 367552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 368552f7358SJed Brown for (i=0; i<bs; i++) { 369552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,PETSC_NULL,buffer,gpiece[r].ncells,PETSC_SCALAR,tag);CHKERRQ(ierr); 370552f7358SJed Brown } 371552f7358SJed Brown } 372552f7358SJed Brown /* all point data */ 373552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 374552f7358SJed Brown PetscInt bs; 375552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 376552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 377552f7358SJed Brown for (i=0; i<bs; i++) { 378552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,PETSC_NULL,buffer,gpiece[r].nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr); 379552f7358SJed Brown } 380552f7358SJed Brown } 381552f7358SJed Brown } 382552f7358SJed Brown } 383552f7358SJed Brown ierr = PetscFree(gpiece);CHKERRQ(ierr); 384552f7358SJed Brown ierr = PetscFree(buffer);CHKERRQ(ierr); 385552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr); 386552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr); 387552f7358SJed Brown PetscFunctionReturn(0); 388552f7358SJed Brown } 389