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 #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; 24ce94432eSBarry Smith MPI_Comm comm; 25552f7358SJed Brown MPI_Datatype mpidatatype; 26552f7358SJed Brown 27552f7358SJed Brown PetscFunctionBegin; 28ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); 29552f7358SJed Brown ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 30552f7358SJed Brown ierr = PetscDataTypeToMPIDataType(datatype,&mpidatatype);CHKERRQ(ierr); 31552f7358SJed Brown 32552f7358SJed Brown if (rank == srank && rank != root) { 33552f7358SJed Brown ierr = MPI_Send((void*)send,count,mpidatatype,root,tag,comm);CHKERRQ(ierr); 34552f7358SJed Brown } else if (rank == root) { 35552f7358SJed Brown const void *buffer; 36552f7358SJed Brown if (root == srank) { /* self */ 37552f7358SJed Brown buffer = send; 38552f7358SJed Brown } else { 39552f7358SJed Brown MPI_Status status; 40552f7358SJed Brown PetscMPIInt nrecv; 41552f7358SJed Brown ierr = MPI_Recv(recv,count,mpidatatype,srank,tag,comm,&status);CHKERRQ(ierr); 42552f7358SJed Brown ierr = MPI_Get_count(&status,mpidatatype,&nrecv);CHKERRQ(ierr); 43552f7358SJed Brown if (count != nrecv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch"); 44552f7358SJed Brown buffer = recv; 45552f7358SJed Brown } 46552f7358SJed Brown ierr = PetscViewerVTKFWrite(viewer,fp,buffer,count,datatype);CHKERRQ(ierr); 47552f7358SJed Brown } 48552f7358SJed Brown PetscFunctionReturn(0); 49552f7358SJed Brown } 50552f7358SJed Brown 51552f7358SJed Brown #undef __FUNCT__ 52552f7358SJed Brown #define __FUNCT__ "DMPlexGetVTKConnectivity" 53552f7358SJed Brown static PetscErrorCode DMPlexGetVTKConnectivity(DM dm,PieceInfo *piece,PetscVTKInt **oconn,PetscVTKInt **ooffsets,PetscVTKType **otypes) 54552f7358SJed Brown { 55552f7358SJed Brown PetscErrorCode ierr; 56552f7358SJed Brown PetscVTKInt *conn,*offsets; 57552f7358SJed Brown PetscVTKType *types; 58552f7358SJed Brown PetscInt dim,vStart,vEnd,cStart,cEnd,pStart,pEnd,cellHeight,cMax,numLabelCells,hasLabel,c,v,countcell,countconn; 59552f7358SJed Brown 60552f7358SJed Brown PetscFunctionBegin; 61dcca6d9dSJed Brown ierr = PetscMalloc3(piece->nconn,&conn,piece->ncells,&offsets,piece->ncells,&types);CHKERRQ(ierr); 62552f7358SJed Brown 63c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 64552f7358SJed Brown ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 65552f7358SJed Brown ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 66552f7358SJed Brown ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 67552f7358SJed Brown ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 680298fd71SBarry Smith ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 698865f1eaSKarl Rupp if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 70c58f1c22SToby Isaac ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 71552f7358SJed Brown hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 72552f7358SJed Brown 73552f7358SJed Brown countcell = 0; 74552f7358SJed Brown countconn = 0; 75552f7358SJed Brown for (c = cStart; c < cEnd; ++c) { 760298fd71SBarry Smith PetscInt *closure = NULL; 77724b5320SMatthew G. Knepley PetscInt closureSize,nverts,celltype,startoffset,nC=0; 78552f7358SJed Brown 79552f7358SJed Brown if (hasLabel) { 80552f7358SJed Brown PetscInt value; 81552f7358SJed Brown 82c58f1c22SToby Isaac ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 83552f7358SJed Brown if (value != 1) continue; 84552f7358SJed Brown } 85552f7358SJed Brown startoffset = countconn; 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); 94724b5320SMatthew G. Knepley ierr = DMPlexInvertCell(dim, nC, &conn[countconn-nC]);CHKERRQ(ierr); 958865f1eaSKarl Rupp 96552f7358SJed Brown offsets[countcell] = countconn; 978865f1eaSKarl Rupp 98552f7358SJed Brown nverts = countconn - startoffset; 99552f7358SJed Brown ierr = DMPlexVTKGetCellType(dm,dim,nverts,&celltype);CHKERRQ(ierr); 1008865f1eaSKarl Rupp 101552f7358SJed Brown types[countcell] = celltype; 102552f7358SJed Brown countcell++; 103552f7358SJed Brown } 104552f7358SJed Brown if (countcell != piece->ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count"); 105552f7358SJed Brown if (countconn != piece->nconn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count"); 106552f7358SJed Brown *oconn = conn; 107552f7358SJed Brown *ooffsets = offsets; 108552f7358SJed Brown *otypes = types; 109552f7358SJed Brown PetscFunctionReturn(0); 110552f7358SJed Brown } 111552f7358SJed Brown 112552f7358SJed Brown #undef __FUNCT__ 113552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll_VTU" 114552f7358SJed Brown /* 115552f7358SJed Brown Write all fields that have been provided to the viewer 116552f7358SJed Brown Multi-block XML format with binary appended data. 117552f7358SJed Brown */ 118552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer) 119552f7358SJed Brown { 120ce94432eSBarry Smith MPI_Comm comm; 121552f7358SJed Brown PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data; 122552f7358SJed Brown PetscViewerVTKObjectLink link; 123552f7358SJed Brown FILE *fp; 124552f7358SJed Brown PetscMPIInt rank,size,tag; 125552f7358SJed Brown PetscErrorCode ierr; 126552f7358SJed Brown PetscInt dim,cellHeight,cStart,cEnd,vStart,vEnd,cMax,numLabelCells,hasLabel,c,v,r,i; 1270298fd71SBarry Smith PieceInfo piece,*gpiece = NULL; 1280298fd71SBarry Smith void *buffer = NULL; 129552f7358SJed Brown 130552f7358SJed Brown PetscFunctionBegin; 131ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 132552f7358SJed Brown #if defined(PETSC_USE_COMPLEX) 133ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Complex values not supported"); 134552f7358SJed Brown #endif 135552f7358SJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 136552f7358SJed Brown ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 137552f7358SJed Brown tag = ((PetscObject)viewer)->tag; 138552f7358SJed Brown 139552f7358SJed Brown ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr); 140552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr); 141519f805aSKarl Rupp #if defined(PETSC_WORDS_BIGENDIAN) 142552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr); 143552f7358SJed Brown #else 144552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr); 145552f7358SJed Brown #endif 146552f7358SJed Brown ierr = PetscFPrintf(comm,fp," <UnstructuredGrid>\n");CHKERRQ(ierr); 147552f7358SJed Brown 148c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 149552f7358SJed Brown ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 150552f7358SJed Brown ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 151552f7358SJed Brown ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1520298fd71SBarry Smith ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1538865f1eaSKarl Rupp if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 154c58f1c22SToby Isaac ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 1558865f1eaSKarl Rupp 156552f7358SJed Brown hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 157552f7358SJed Brown piece.nvertices = vEnd - vStart; 158552f7358SJed Brown piece.ncells = 0; 159552f7358SJed Brown piece.nconn = 0; 160552f7358SJed Brown for (c = cStart; c < cEnd; ++c) { 1610298fd71SBarry Smith PetscInt *closure = NULL; 162552f7358SJed Brown PetscInt closureSize; 163552f7358SJed Brown 164552f7358SJed Brown if (hasLabel) { 165552f7358SJed Brown PetscInt value; 166552f7358SJed Brown 167c58f1c22SToby Isaac ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 168552f7358SJed Brown if (value != 1) continue; 169552f7358SJed Brown } 170552f7358SJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 171552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 172552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) piece.nconn++; 173552f7358SJed Brown } 174552f7358SJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 175552f7358SJed Brown piece.ncells++; 176552f7358SJed Brown } 177785e854fSJed Brown if (!rank) {ierr = PetscMalloc1(size,&gpiece);CHKERRQ(ierr);} 178552f7358SJed Brown ierr = MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm);CHKERRQ(ierr); 179552f7358SJed Brown 180552f7358SJed Brown /* 181552f7358SJed Brown * Write file header 182552f7358SJed Brown */ 183552f7358SJed Brown if (!rank) { 184552f7358SJed Brown PetscInt boffset = 0; 185552f7358SJed Brown 186552f7358SJed Brown for (r=0; r<size; r++) { 187552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells);CHKERRQ(ierr); 188552f7358SJed Brown /* Coordinate positions */ 189552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Points>\n");CHKERRQ(ierr); 190552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr); 191552f7358SJed Brown boffset += gpiece[r].nvertices*3*sizeof(PetscScalar) + sizeof(int); 192552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Points>\n");CHKERRQ(ierr); 193552f7358SJed Brown /* Cell connectivity */ 194552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Cells>\n");CHKERRQ(ierr); 195552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 196552f7358SJed Brown boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int); 197552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"offsets\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 198552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int); 199552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"UInt8\" Name=\"types\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 200552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int); 201552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Cells>\n");CHKERRQ(ierr); 202552f7358SJed Brown 203552f7358SJed Brown /* 204552f7358SJed Brown * Cell Data headers 205552f7358SJed Brown */ 206552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <CellData>\n");CHKERRQ(ierr); 207552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 208552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(int) + sizeof(int); 209552f7358SJed Brown /* all the vectors */ 210552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 211552f7358SJed Brown Vec X = (Vec)link->vec; 2121cfafdd3SJed Brown PetscInt bs,nfields,field; 213552f7358SJed Brown const char *vecname = ""; 214552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 215552f7358SJed 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. */ 216552f7358SJed Brown ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr); 217552f7358SJed Brown } 218552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr); 2191cfafdd3SJed Brown ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr); 2207ded4cbaSJed Brown for (field=0,i=0; field<(nfields?nfields:1); field++) { 2211cfafdd3SJed Brown PetscInt fbs,j; 222*a00cdb45SToby Isaac PetscFV fv = NULL; 223*a00cdb45SToby Isaac PetscObject f; 224*a00cdb45SToby Isaac PetscClassId fClass; 2250298fd71SBarry Smith const char *fieldname = NULL; 2261cfafdd3SJed Brown char buf[256]; 2277ded4cbaSJed Brown if (nfields) { /* We have user-defined fields/components */ 2281cfafdd3SJed Brown ierr = PetscSectionGetFieldDof(dm->defaultSection,cStart,field,&fbs);CHKERRQ(ierr); 2291cfafdd3SJed Brown ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr); 2307ded4cbaSJed Brown } else fbs = bs; /* Say we have one field with 'bs' components */ 231*a00cdb45SToby Isaac ierr = DMGetField(dm,field,&f);CHKERRQ(ierr); 232*a00cdb45SToby Isaac ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr); 233*a00cdb45SToby Isaac if (fClass == PETSCFV_CLASSID) { 234*a00cdb45SToby Isaac fv = (PetscFV) f; 235*a00cdb45SToby Isaac } 236552f7358SJed Brown if (!fieldname) { 2371cfafdd3SJed Brown ierr = PetscSNPrintf(buf,sizeof(buf),"CellField%D",field);CHKERRQ(ierr); 238552f7358SJed Brown fieldname = buf; 239552f7358SJed Brown } 2401cfafdd3SJed Brown for (j=0; j<fbs; j++) { 241*a00cdb45SToby Isaac const char *compName = NULL; 242*a00cdb45SToby Isaac if (fv) { 243*a00cdb45SToby Isaac ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr); 244*a00cdb45SToby Isaac } 245*a00cdb45SToby Isaac if (compName) { 246*a00cdb45SToby 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); 247*a00cdb45SToby Isaac } 248*a00cdb45SToby Isaac else { 2491cfafdd3SJed 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); 250*a00cdb45SToby Isaac } 251552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(PetscScalar) + sizeof(int); 2521cfafdd3SJed Brown i++; 253552f7358SJed Brown } 254552f7358SJed Brown } 2551cfafdd3SJed Brown if (i != bs) SETERRQ2(comm,PETSC_ERR_PLIB,"Total number of field components %D != block size %D",i,bs); 2561cfafdd3SJed Brown } 257552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </CellData>\n");CHKERRQ(ierr); 258552f7358SJed Brown 259552f7358SJed Brown /* 260552f7358SJed Brown * Point Data headers 261552f7358SJed Brown */ 262552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <PointData>\n");CHKERRQ(ierr); 263552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 264552f7358SJed Brown Vec X = (Vec)link->vec; 2651cfafdd3SJed Brown PetscInt bs,nfields,field; 266552f7358SJed Brown const char *vecname = ""; 267552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 268552f7358SJed 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. */ 269552f7358SJed Brown ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr); 270552f7358SJed Brown } 271552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 2721cfafdd3SJed Brown ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr); 2737ded4cbaSJed Brown for (field=0,i=0; field<(nfields?nfields:1); field++) { 2741cfafdd3SJed Brown PetscInt fbs,j; 2751cfafdd3SJed Brown const char *fieldname = NULL; 2761cfafdd3SJed Brown char buf[256]; 2777ded4cbaSJed Brown if (nfields) { /* We have user-defined fields/components */ 2781cfafdd3SJed Brown ierr = PetscSectionGetFieldDof(dm->defaultSection,vStart,field,&fbs);CHKERRQ(ierr); 2791cfafdd3SJed Brown ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr); 2807ded4cbaSJed Brown } else fbs = bs; /* Say we have one field with 'bs' components */ 2811cfafdd3SJed Brown if (!fieldname) { 2821cfafdd3SJed Brown ierr = PetscSNPrintf(buf,sizeof(buf),"PointField%D",field);CHKERRQ(ierr); 2831cfafdd3SJed Brown fieldname = buf; 2841cfafdd3SJed Brown } 2851cfafdd3SJed Brown for (j=0; j<fbs; j++) { 2861cfafdd3SJed 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); 287552f7358SJed Brown boffset += gpiece[r].nvertices*sizeof(PetscScalar) + sizeof(int); 288552f7358SJed Brown } 289552f7358SJed Brown } 2901cfafdd3SJed Brown } 291552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </PointData>\n");CHKERRQ(ierr); 292552f7358SJed Brown 293552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Piece>\n");CHKERRQ(ierr); 294552f7358SJed Brown } 295552f7358SJed Brown } 296552f7358SJed Brown 297552f7358SJed Brown ierr = PetscFPrintf(comm,fp," </UnstructuredGrid>\n");CHKERRQ(ierr); 298552f7358SJed Brown ierr = PetscFPrintf(comm,fp," <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr); 299552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr); 300552f7358SJed Brown 301552f7358SJed Brown if (!rank) { 302552f7358SJed Brown PetscInt maxsize = 0; 303552f7358SJed Brown for (r=0; r<size; r++) { 304552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscScalar))); 305552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*sizeof(PetscScalar))); 306552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt))); 307552f7358SJed Brown } 308552f7358SJed Brown ierr = PetscMalloc(maxsize,&buffer);CHKERRQ(ierr); 309552f7358SJed Brown } 310552f7358SJed Brown for (r=0; r<size; r++) { 311552f7358SJed Brown if (r == rank) { 312552f7358SJed Brown PetscInt nsend; 313552f7358SJed Brown { /* Position */ 314552f7358SJed Brown const PetscScalar *x; 3150298fd71SBarry Smith PetscScalar *y = NULL; 316552f7358SJed Brown Vec coords; 317552f7358SJed Brown nsend = piece.nvertices*3; 318552f7358SJed Brown ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr); 319552f7358SJed Brown ierr = VecGetArrayRead(coords,&x);CHKERRQ(ierr); 320552f7358SJed Brown if (dim != 3) { 321785e854fSJed Brown ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr); 322552f7358SJed Brown for (i=0; i<piece.nvertices; i++) { 323552f7358SJed Brown y[i*3+0] = x[i*dim+0]; 324552f7358SJed Brown y[i*3+1] = (dim > 1) ? x[i*dim+1] : 0; 325552f7358SJed Brown y[i*3+2] = 0; 326552f7358SJed Brown } 327552f7358SJed Brown } 328552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,y ? y : x,buffer,nsend,PETSC_SCALAR,tag);CHKERRQ(ierr); 329552f7358SJed Brown ierr = PetscFree(y);CHKERRQ(ierr); 330552f7358SJed Brown ierr = VecRestoreArrayRead(coords,&x);CHKERRQ(ierr); 331552f7358SJed Brown } 332552f7358SJed Brown { /* Connectivity, offsets, types */ 3333e869605SMatthew G. Knepley PetscVTKInt *connectivity = NULL, *offsets = NULL; 3343e869605SMatthew G. Knepley PetscVTKType *types = NULL; 335552f7358SJed Brown ierr = DMPlexGetVTKConnectivity(dm,&piece,&connectivity,&offsets,&types);CHKERRQ(ierr); 336552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,connectivity,buffer,piece.nconn,PETSC_INT32,tag);CHKERRQ(ierr); 337552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,offsets,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr); 338552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,types,buffer,piece.ncells,PETSC_UINT8,tag);CHKERRQ(ierr); 339552f7358SJed Brown ierr = PetscFree3(connectivity,offsets,types);CHKERRQ(ierr); 340552f7358SJed Brown } 341552f7358SJed Brown { /* Owners (cell data) */ 342552f7358SJed Brown PetscVTKInt *owners; 343785e854fSJed Brown ierr = PetscMalloc1(piece.ncells,&owners);CHKERRQ(ierr); 344552f7358SJed Brown for (i=0; i<piece.ncells; i++) owners[i] = rank; 345552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,owners,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr); 346552f7358SJed Brown ierr = PetscFree(owners);CHKERRQ(ierr); 347552f7358SJed Brown } 348552f7358SJed Brown /* Cell data */ 349552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 350552f7358SJed Brown Vec X = (Vec)link->vec; 351552f7358SJed Brown const PetscScalar *x; 352552f7358SJed Brown PetscScalar *y; 353552f7358SJed Brown PetscInt bs; 354552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 355552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr); 356552f7358SJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 357785e854fSJed Brown ierr = PetscMalloc1(piece.ncells,&y);CHKERRQ(ierr); 358552f7358SJed Brown for (i=0; i<bs; i++) { 359552f7358SJed Brown PetscInt cnt; 360552f7358SJed Brown for (c=cStart,cnt=0; c<cEnd; c++) { 361552f7358SJed Brown const PetscScalar *xpoint; 362552f7358SJed Brown if (hasLabel) { /* Ignore some cells */ 363552f7358SJed Brown PetscInt value; 364c58f1c22SToby Isaac ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 365552f7358SJed Brown if (value != 1) continue; 366552f7358SJed Brown } 367552f7358SJed Brown ierr = DMPlexPointLocalRead(dm,c,x,&xpoint);CHKERRQ(ierr); 368552f7358SJed Brown y[cnt++] = xpoint[i]; 369552f7358SJed Brown } 370552f7358SJed Brown if (cnt != piece.ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 371552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells,PETSC_SCALAR,tag);CHKERRQ(ierr); 372552f7358SJed Brown } 373552f7358SJed Brown ierr = PetscFree(y);CHKERRQ(ierr); 374552f7358SJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 375552f7358SJed Brown } 376552f7358SJed Brown 377552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 378552f7358SJed Brown Vec X = (Vec)link->vec; 379552f7358SJed Brown const PetscScalar *x; 380552f7358SJed Brown PetscScalar *y; 381552f7358SJed Brown PetscInt bs; 382552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 383552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 384552f7358SJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 385785e854fSJed Brown ierr = PetscMalloc1(piece.nvertices,&y);CHKERRQ(ierr); 386552f7358SJed Brown for (i=0; i<bs; i++) { 387552f7358SJed Brown PetscInt cnt; 388552f7358SJed Brown for (v=vStart,cnt=0; v<vEnd; v++) { 389552f7358SJed Brown const PetscScalar *xpoint; 39037045cddSJed Brown ierr = DMPlexPointLocalRead(dm,v,x,&xpoint);CHKERRQ(ierr); 391552f7358SJed Brown y[cnt++] = xpoint[i]; 392552f7358SJed Brown } 393552f7358SJed Brown if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 394552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr); 395552f7358SJed Brown } 396552f7358SJed Brown ierr = PetscFree(y);CHKERRQ(ierr); 397552f7358SJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 398552f7358SJed Brown } 399552f7358SJed Brown } else if (!rank) { 4000298fd71SBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,PETSC_SCALAR,tag);CHKERRQ(ierr); /* positions */ 4010298fd71SBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nconn,PETSC_INT32,tag);CHKERRQ(ierr); /* connectivity */ 4020298fd71SBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* offsets */ 4030298fd71SBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_UINT8,tag);CHKERRQ(ierr); /* types */ 4040298fd71SBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* owner rank (cells) */ 405552f7358SJed Brown /* all cell data */ 406552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 407552f7358SJed Brown PetscInt bs; 408552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 409fd68c46aSJed Brown ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr); 410552f7358SJed Brown for (i=0; i<bs; i++) { 4110298fd71SBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_SCALAR,tag);CHKERRQ(ierr); 412552f7358SJed Brown } 413552f7358SJed Brown } 414552f7358SJed Brown /* all point data */ 415552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 416552f7358SJed Brown PetscInt bs; 417552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 418552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 419552f7358SJed Brown for (i=0; i<bs; i++) { 4200298fd71SBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr); 421552f7358SJed Brown } 422552f7358SJed Brown } 423552f7358SJed Brown } 424552f7358SJed Brown } 425552f7358SJed Brown ierr = PetscFree(gpiece);CHKERRQ(ierr); 426552f7358SJed Brown ierr = PetscFree(buffer);CHKERRQ(ierr); 427552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr); 428552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr); 4295e32de72SMatthew G. Knepley ierr = PetscFClose(comm,fp);CHKERRQ(ierr); 430552f7358SJed Brown PetscFunctionReturn(0); 431552f7358SJed Brown } 432