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 492e529ec8SStefano Zampini static PetscErrorCode DMPlexGetVTKConnectivity(DM dm, PetscBool localized, PieceInfo *piece,PetscVTKInt **oconn,PetscVTKInt **ooffsets,PetscVTKType **otypes) 50552f7358SJed Brown { 51552f7358SJed Brown PetscErrorCode ierr; 522e529ec8SStefano 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); 592e529ec8SStefano 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*19003fb5SStefano Zampini PetscInt nverts,dof = 0,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*19003fb5SStefano Zampini if (localized) { 83*19003fb5SStefano Zampini ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr); 84*19003fb5SStefano Zampini } 85*19003fb5SStefano Zampini if (!dof) { 862e529ec8SStefano Zampini PetscInt *closure = NULL; 872e529ec8SStefano Zampini PetscInt closureSize; 882e529ec8SStefano Zampini 89552f7358SJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 90552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 91552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 92*19003fb5SStefano Zampini if (!localized) conn[countconn++] = closure[v] - vStart; 93*19003fb5SStefano Zampini else conn[countconn++] = startoffset + nC; 94724b5320SMatthew G. Knepley ++nC; 95552f7358SJed Brown } 96552f7358SJed Brown } 97552f7358SJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 982e529ec8SStefano Zampini } else { 992e529ec8SStefano Zampini for (nC = 0; nC < dof/dim; nC++) conn[countconn++] = startoffset + nC; 1002e529ec8SStefano Zampini } 101724b5320SMatthew G. Knepley ierr = DMPlexInvertCell(dim, nC, &conn[countconn-nC]);CHKERRQ(ierr); 1028865f1eaSKarl Rupp 103552f7358SJed Brown offsets[countcell] = countconn; 1048865f1eaSKarl Rupp 105552f7358SJed Brown nverts = countconn - startoffset; 106de0a4605SMatthew G. Knepley ierr = DMPlexVTKGetCellType_Internal(dm,dim,nverts,&celltype);CHKERRQ(ierr); 1078865f1eaSKarl Rupp 108552f7358SJed Brown types[countcell] = celltype; 109552f7358SJed Brown countcell++; 110552f7358SJed Brown } 111552f7358SJed Brown if (countcell != piece->ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count"); 112552f7358SJed Brown if (countconn != piece->nconn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count"); 113552f7358SJed Brown *oconn = conn; 114552f7358SJed Brown *ooffsets = offsets; 115552f7358SJed Brown *otypes = types; 116552f7358SJed Brown PetscFunctionReturn(0); 117552f7358SJed Brown } 118552f7358SJed Brown 119552f7358SJed Brown /* 120552f7358SJed Brown Write all fields that have been provided to the viewer 121552f7358SJed Brown Multi-block XML format with binary appended data. 122552f7358SJed Brown */ 123552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer) 124552f7358SJed Brown { 125ce94432eSBarry Smith MPI_Comm comm; 1262e529ec8SStefano Zampini PetscSection coordSection; 127552f7358SJed Brown PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data; 128552f7358SJed Brown PetscViewerVTKObjectLink link; 129552f7358SJed Brown FILE *fp; 130552f7358SJed Brown PetscMPIInt rank,size,tag; 131552f7358SJed Brown PetscErrorCode ierr; 132*19003fb5SStefano Zampini PetscInt dimEmbed,cellHeight,cStart,cEnd,vStart,vEnd,cMax,numLabelCells,hasLabel,c,v,r,i; 1332e529ec8SStefano Zampini PetscBool localized; 1340298fd71SBarry Smith PieceInfo piece,*gpiece = NULL; 1350298fd71SBarry Smith void *buffer = NULL; 136552f7358SJed Brown 137552f7358SJed Brown PetscFunctionBegin; 138ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 139552f7358SJed Brown #if defined(PETSC_USE_COMPLEX) 140ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Complex values not supported"); 141552f7358SJed Brown #endif 142552f7358SJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 143552f7358SJed Brown ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 144552f7358SJed Brown tag = ((PetscObject)viewer)->tag; 145552f7358SJed Brown 146552f7358SJed Brown ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr); 147552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr); 148519f805aSKarl Rupp #if defined(PETSC_WORDS_BIGENDIAN) 149552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr); 150552f7358SJed Brown #else 151552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr); 152552f7358SJed Brown #endif 153552f7358SJed Brown ierr = PetscFPrintf(comm,fp," <UnstructuredGrid>\n");CHKERRQ(ierr); 154552f7358SJed Brown 1553baa072aSToby Isaac ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 156552f7358SJed Brown ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 157552f7358SJed Brown ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 158552f7358SJed Brown ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1590298fd71SBarry Smith ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1608865f1eaSKarl Rupp if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 161c58f1c22SToby Isaac ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 1622e529ec8SStefano Zampini ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr); 1632e529ec8SStefano Zampini ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1648865f1eaSKarl Rupp 165552f7358SJed Brown hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 1662e529ec8SStefano Zampini piece.nvertices = 0; 167552f7358SJed Brown piece.ncells = 0; 168552f7358SJed Brown piece.nconn = 0; 1692e529ec8SStefano Zampini if (!localized) piece.nvertices = vEnd - vStart; 170552f7358SJed Brown for (c = cStart; c < cEnd; ++c) { 171*19003fb5SStefano Zampini PetscInt dof = 0; 172552f7358SJed Brown if (hasLabel) { 173552f7358SJed Brown PetscInt value; 174552f7358SJed Brown 175c58f1c22SToby Isaac ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 176552f7358SJed Brown if (value != 1) continue; 177552f7358SJed Brown } 178*19003fb5SStefano Zampini if (localized) { 179*19003fb5SStefano Zampini ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr); 180*19003fb5SStefano Zampini } 181*19003fb5SStefano Zampini if (!dof) { 1822e529ec8SStefano Zampini PetscInt *closure = NULL; 1832e529ec8SStefano Zampini PetscInt closureSize; 1842e529ec8SStefano Zampini 185552f7358SJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 186552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 1872e529ec8SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 1882e529ec8SStefano Zampini piece.nconn++; 189*19003fb5SStefano Zampini if (localized) piece.nvertices++; 1902e529ec8SStefano Zampini } 191552f7358SJed Brown } 192552f7358SJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1932e529ec8SStefano Zampini } else { 1942e529ec8SStefano Zampini piece.nvertices += dof/dimEmbed; 1952e529ec8SStefano Zampini piece.nconn += dof/dimEmbed; 1962e529ec8SStefano Zampini } 197552f7358SJed Brown piece.ncells++; 198552f7358SJed Brown } 199785e854fSJed Brown if (!rank) {ierr = PetscMalloc1(size,&gpiece);CHKERRQ(ierr);} 200552f7358SJed Brown ierr = MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm);CHKERRQ(ierr); 201552f7358SJed Brown 202552f7358SJed Brown /* 203552f7358SJed Brown * Write file header 204552f7358SJed Brown */ 205552f7358SJed Brown if (!rank) { 206552f7358SJed Brown PetscInt boffset = 0; 207552f7358SJed Brown 208552f7358SJed Brown for (r=0; r<size; r++) { 209552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells);CHKERRQ(ierr); 210552f7358SJed Brown /* Coordinate positions */ 211552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Points>\n");CHKERRQ(ierr); 212552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr); 213552f7358SJed Brown boffset += gpiece[r].nvertices*3*sizeof(PetscScalar) + sizeof(int); 214552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Points>\n");CHKERRQ(ierr); 215552f7358SJed Brown /* Cell connectivity */ 216552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Cells>\n");CHKERRQ(ierr); 217552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 218552f7358SJed Brown boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int); 219552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"offsets\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 220552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int); 221552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"UInt8\" Name=\"types\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 222552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int); 223552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Cells>\n");CHKERRQ(ierr); 224552f7358SJed Brown 225552f7358SJed Brown /* 226552f7358SJed Brown * Cell Data headers 227552f7358SJed Brown */ 228552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <CellData>\n");CHKERRQ(ierr); 229552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 230552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(int) + sizeof(int); 231552f7358SJed Brown /* all the vectors */ 232552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 233552f7358SJed Brown Vec X = (Vec)link->vec; 2341cfafdd3SJed Brown PetscInt bs,nfields,field; 235552f7358SJed Brown const char *vecname = ""; 236552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 237552f7358SJed 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. */ 238552f7358SJed Brown ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr); 239552f7358SJed Brown } 240552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr); 2411cfafdd3SJed Brown ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr); 2427ded4cbaSJed Brown for (field=0,i=0; field<(nfields?nfields:1); field++) { 2431cfafdd3SJed Brown PetscInt fbs,j; 244a00cdb45SToby Isaac PetscFV fv = NULL; 245a00cdb45SToby Isaac PetscObject f; 246a00cdb45SToby Isaac PetscClassId fClass; 2470298fd71SBarry Smith const char *fieldname = NULL; 2481cfafdd3SJed Brown char buf[256]; 2497ded4cbaSJed Brown if (nfields) { /* We have user-defined fields/components */ 2501cfafdd3SJed Brown ierr = PetscSectionGetFieldDof(dm->defaultSection,cStart,field,&fbs);CHKERRQ(ierr); 2511cfafdd3SJed Brown ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr); 2527ded4cbaSJed Brown } else fbs = bs; /* Say we have one field with 'bs' components */ 253a00cdb45SToby Isaac ierr = DMGetField(dm,field,&f);CHKERRQ(ierr); 254a00cdb45SToby Isaac ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr); 255a00cdb45SToby Isaac if (fClass == PETSCFV_CLASSID) { 256a00cdb45SToby Isaac fv = (PetscFV) f; 257a00cdb45SToby Isaac } 258552f7358SJed Brown if (!fieldname) { 2591cfafdd3SJed Brown ierr = PetscSNPrintf(buf,sizeof(buf),"CellField%D",field);CHKERRQ(ierr); 260552f7358SJed Brown fieldname = buf; 261552f7358SJed Brown } 2621cfafdd3SJed Brown for (j=0; j<fbs; j++) { 263a00cdb45SToby Isaac const char *compName = NULL; 264a00cdb45SToby Isaac if (fv) { 265a00cdb45SToby Isaac ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr); 266a00cdb45SToby Isaac } 267a00cdb45SToby Isaac if (compName) { 268a00cdb45SToby 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); 269a00cdb45SToby Isaac } 270a00cdb45SToby Isaac else { 2711cfafdd3SJed 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); 272a00cdb45SToby Isaac } 273552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(PetscScalar) + sizeof(int); 2741cfafdd3SJed Brown i++; 275552f7358SJed Brown } 276552f7358SJed Brown } 2771cfafdd3SJed Brown if (i != bs) SETERRQ2(comm,PETSC_ERR_PLIB,"Total number of field components %D != block size %D",i,bs); 2781cfafdd3SJed Brown } 279552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </CellData>\n");CHKERRQ(ierr); 280552f7358SJed Brown 281552f7358SJed Brown /* 282552f7358SJed Brown * Point Data headers 283552f7358SJed Brown */ 284552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <PointData>\n");CHKERRQ(ierr); 285552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 286552f7358SJed Brown Vec X = (Vec)link->vec; 2871cfafdd3SJed Brown PetscInt bs,nfields,field; 288552f7358SJed Brown const char *vecname = ""; 289552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 290552f7358SJed 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. */ 291552f7358SJed Brown ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr); 292552f7358SJed Brown } 293552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 2941cfafdd3SJed Brown ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr); 2957ded4cbaSJed Brown for (field=0,i=0; field<(nfields?nfields:1); field++) { 2961cfafdd3SJed Brown PetscInt fbs,j; 2971cfafdd3SJed Brown const char *fieldname = NULL; 2981cfafdd3SJed Brown char buf[256]; 2997ded4cbaSJed Brown if (nfields) { /* We have user-defined fields/components */ 3001cfafdd3SJed Brown ierr = PetscSectionGetFieldDof(dm->defaultSection,vStart,field,&fbs);CHKERRQ(ierr); 3011cfafdd3SJed Brown ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr); 3027ded4cbaSJed Brown } else fbs = bs; /* Say we have one field with 'bs' components */ 3031cfafdd3SJed Brown if (!fieldname) { 3041cfafdd3SJed Brown ierr = PetscSNPrintf(buf,sizeof(buf),"PointField%D",field);CHKERRQ(ierr); 3051cfafdd3SJed Brown fieldname = buf; 3061cfafdd3SJed Brown } 3071cfafdd3SJed Brown for (j=0; j<fbs; j++) { 3081cfafdd3SJed 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); 309552f7358SJed Brown boffset += gpiece[r].nvertices*sizeof(PetscScalar) + sizeof(int); 310552f7358SJed Brown } 311552f7358SJed Brown } 3121cfafdd3SJed Brown } 313552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </PointData>\n");CHKERRQ(ierr); 314552f7358SJed Brown 315552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Piece>\n");CHKERRQ(ierr); 316552f7358SJed Brown } 317552f7358SJed Brown } 318552f7358SJed Brown 319552f7358SJed Brown ierr = PetscFPrintf(comm,fp," </UnstructuredGrid>\n");CHKERRQ(ierr); 320552f7358SJed Brown ierr = PetscFPrintf(comm,fp," <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr); 321552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr); 322552f7358SJed Brown 323552f7358SJed Brown if (!rank) { 324552f7358SJed Brown PetscInt maxsize = 0; 325552f7358SJed Brown for (r=0; r<size; r++) { 326552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscScalar))); 327552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*sizeof(PetscScalar))); 328552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt))); 329552f7358SJed Brown } 330552f7358SJed Brown ierr = PetscMalloc(maxsize,&buffer);CHKERRQ(ierr); 331552f7358SJed Brown } 332552f7358SJed Brown for (r=0; r<size; r++) { 333552f7358SJed Brown if (r == rank) { 334552f7358SJed Brown PetscInt nsend; 335552f7358SJed Brown { /* Position */ 336552f7358SJed Brown const PetscScalar *x; 3370298fd71SBarry Smith PetscScalar *y = NULL; 338552f7358SJed Brown Vec coords; 3392e529ec8SStefano Zampini 340552f7358SJed Brown ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr); 341552f7358SJed Brown ierr = VecGetArrayRead(coords,&x);CHKERRQ(ierr); 342*19003fb5SStefano Zampini if (dimEmbed != 3 || localized) { 343785e854fSJed Brown ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr); 344*19003fb5SStefano Zampini if (localized) { 345*19003fb5SStefano Zampini PetscInt cnt; 346*19003fb5SStefano Zampini for (c=cStart,cnt=0; c<cEnd; c++) { 347*19003fb5SStefano Zampini PetscInt off, dof; 348*19003fb5SStefano Zampini 349*19003fb5SStefano Zampini ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr); 350*19003fb5SStefano Zampini if (!dof) { 351*19003fb5SStefano Zampini PetscInt *closure = NULL; 352*19003fb5SStefano Zampini PetscInt closureSize; 353*19003fb5SStefano Zampini 354*19003fb5SStefano Zampini ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 355*19003fb5SStefano Zampini for (v = 0; v < closureSize*2; v += 2) { 356*19003fb5SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 357*19003fb5SStefano Zampini ierr = PetscSectionGetOffset(coordSection, closure[v], &off);CHKERRQ(ierr); 358*19003fb5SStefano Zampini if (dimEmbed != 3) { 359*19003fb5SStefano Zampini y[cnt*3+0] = x[off+0]; 360*19003fb5SStefano Zampini y[cnt*3+1] = (dimEmbed > 1) ? x[off+1] : 0.0; 361*19003fb5SStefano Zampini y[cnt*3+2] = 0.0; 362*19003fb5SStefano Zampini } else { 363*19003fb5SStefano Zampini y[cnt*3+0] = x[off+0]; 364*19003fb5SStefano Zampini y[cnt*3+1] = x[off+1]; 365*19003fb5SStefano Zampini y[cnt*3+2] = x[off+2]; 366*19003fb5SStefano Zampini } 367*19003fb5SStefano Zampini cnt++; 368*19003fb5SStefano Zampini } 369*19003fb5SStefano Zampini } 370*19003fb5SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 371*19003fb5SStefano Zampini } else { 372*19003fb5SStefano Zampini ierr = PetscSectionGetOffset(coordSection, c, &off);CHKERRQ(ierr); 373*19003fb5SStefano Zampini if (dimEmbed != 3) { 374*19003fb5SStefano Zampini for (i=0; i<dof/dimEmbed; i++) { 375*19003fb5SStefano Zampini y[cnt*3+0] = x[off + i*dimEmbed + 0]; 376*19003fb5SStefano Zampini y[cnt*3+1] = (dimEmbed > 1) ? x[off + i*dimEmbed + 1] : 0.0; 377*19003fb5SStefano Zampini y[cnt*3+2] = 0.0; 378*19003fb5SStefano Zampini cnt++; 379*19003fb5SStefano Zampini } 380*19003fb5SStefano Zampini } else { 381*19003fb5SStefano Zampini for (i=0; i<dof; i ++) { 382*19003fb5SStefano Zampini y[cnt*3+i] = x[off + i]; 383*19003fb5SStefano Zampini } 384*19003fb5SStefano Zampini cnt += dof/dimEmbed; 385*19003fb5SStefano Zampini } 386*19003fb5SStefano Zampini } 387*19003fb5SStefano Zampini } 388*19003fb5SStefano Zampini if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 389*19003fb5SStefano Zampini } else { 390552f7358SJed Brown for (i=0; i<piece.nvertices; i++) { 3913baa072aSToby Isaac y[i*3+0] = x[i*dimEmbed+0]; 3923baa072aSToby Isaac y[i*3+1] = (dimEmbed > 1) ? x[i*dimEmbed+1] : 0; 393*19003fb5SStefano Zampini y[i*3+2] = 0.0; 394*19003fb5SStefano Zampini } 395552f7358SJed Brown } 396552f7358SJed Brown } 3972e529ec8SStefano Zampini nsend = piece.nvertices*3; 398552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,y ? y : x,buffer,nsend,PETSC_SCALAR,tag);CHKERRQ(ierr); 399552f7358SJed Brown ierr = PetscFree(y);CHKERRQ(ierr); 400552f7358SJed Brown ierr = VecRestoreArrayRead(coords,&x);CHKERRQ(ierr); 401552f7358SJed Brown } 402552f7358SJed Brown { /* Connectivity, offsets, types */ 4033e869605SMatthew G. Knepley PetscVTKInt *connectivity = NULL, *offsets = NULL; 4043e869605SMatthew G. Knepley PetscVTKType *types = NULL; 4052e529ec8SStefano Zampini ierr = DMPlexGetVTKConnectivity(dm,localized,&piece,&connectivity,&offsets,&types);CHKERRQ(ierr); 406552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,connectivity,buffer,piece.nconn,PETSC_INT32,tag);CHKERRQ(ierr); 407552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,offsets,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr); 408552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,types,buffer,piece.ncells,PETSC_UINT8,tag);CHKERRQ(ierr); 409552f7358SJed Brown ierr = PetscFree3(connectivity,offsets,types);CHKERRQ(ierr); 410552f7358SJed Brown } 411552f7358SJed Brown { /* Owners (cell data) */ 412552f7358SJed Brown PetscVTKInt *owners; 413785e854fSJed Brown ierr = PetscMalloc1(piece.ncells,&owners);CHKERRQ(ierr); 414552f7358SJed Brown for (i=0; i<piece.ncells; i++) owners[i] = rank; 415552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,owners,buffer,piece.ncells,PETSC_INT32,tag);CHKERRQ(ierr); 416552f7358SJed Brown ierr = PetscFree(owners);CHKERRQ(ierr); 417552f7358SJed Brown } 418552f7358SJed Brown /* Cell data */ 419552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 420552f7358SJed Brown Vec X = (Vec)link->vec; 421552f7358SJed Brown const PetscScalar *x; 422552f7358SJed Brown PetscScalar *y; 423552f7358SJed Brown PetscInt bs; 424552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 425552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr); 426552f7358SJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 427785e854fSJed Brown ierr = PetscMalloc1(piece.ncells,&y);CHKERRQ(ierr); 428552f7358SJed Brown for (i=0; i<bs; i++) { 429552f7358SJed Brown PetscInt cnt; 430552f7358SJed Brown for (c=cStart,cnt=0; c<cEnd; c++) { 431640bce14SSatish Balay PetscScalar *xpoint; 432552f7358SJed Brown if (hasLabel) { /* Ignore some cells */ 433552f7358SJed Brown PetscInt value; 434c58f1c22SToby Isaac ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 435552f7358SJed Brown if (value != 1) continue; 436552f7358SJed Brown } 437552f7358SJed Brown ierr = DMPlexPointLocalRead(dm,c,x,&xpoint);CHKERRQ(ierr); 438552f7358SJed Brown y[cnt++] = xpoint[i]; 439552f7358SJed Brown } 440552f7358SJed Brown if (cnt != piece.ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 441552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells,PETSC_SCALAR,tag);CHKERRQ(ierr); 442552f7358SJed Brown } 443552f7358SJed Brown ierr = PetscFree(y);CHKERRQ(ierr); 444552f7358SJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 445552f7358SJed Brown } 446552f7358SJed Brown 447552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 448552f7358SJed Brown Vec X = (Vec)link->vec; 449552f7358SJed Brown const PetscScalar *x; 450552f7358SJed Brown PetscScalar *y; 451552f7358SJed Brown PetscInt bs; 452552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 453552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 454552f7358SJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 455785e854fSJed Brown ierr = PetscMalloc1(piece.nvertices,&y);CHKERRQ(ierr); 456552f7358SJed Brown for (i=0; i<bs; i++) { 457552f7358SJed Brown PetscInt cnt; 4582e529ec8SStefano Zampini if (!localized) { 459552f7358SJed Brown for (v=vStart,cnt=0; v<vEnd; v++) { 460640bce14SSatish Balay PetscScalar *xpoint; 46137045cddSJed Brown ierr = DMPlexPointLocalRead(dm,v,x,&xpoint);CHKERRQ(ierr); 462552f7358SJed Brown y[cnt++] = xpoint[i]; 463552f7358SJed Brown } 4642e529ec8SStefano Zampini } else { 4652e529ec8SStefano Zampini for (c=cStart,cnt=0; c<cEnd; c++) { 4662e529ec8SStefano Zampini PetscInt *closure = NULL; 4672e529ec8SStefano Zampini PetscInt closureSize, off; 4682e529ec8SStefano Zampini 4692e529ec8SStefano Zampini ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 4702e529ec8SStefano Zampini for (v = 0, off = 0; v < closureSize*2; v += 2) { 4712e529ec8SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 4722e529ec8SStefano Zampini PetscScalar *xpoint; 4732e529ec8SStefano Zampini 4742e529ec8SStefano Zampini ierr = DMPlexPointLocalRead(dm,v,x,&xpoint);CHKERRQ(ierr); 4752e529ec8SStefano Zampini y[cnt + off++] = xpoint[i]; 4762e529ec8SStefano Zampini } 4772e529ec8SStefano Zampini } 4782e529ec8SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 4792e529ec8SStefano Zampini } 4802e529ec8SStefano Zampini } 481552f7358SJed Brown if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 482552f7358SJed Brown ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr); 483552f7358SJed Brown } 484552f7358SJed Brown ierr = PetscFree(y);CHKERRQ(ierr); 485552f7358SJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 486552f7358SJed Brown } 487552f7358SJed Brown } else if (!rank) { 4880298fd71SBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,PETSC_SCALAR,tag);CHKERRQ(ierr); /* positions */ 4890298fd71SBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nconn,PETSC_INT32,tag);CHKERRQ(ierr); /* connectivity */ 4900298fd71SBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* offsets */ 4910298fd71SBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_UINT8,tag);CHKERRQ(ierr); /* types */ 4920298fd71SBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_INT32,tag);CHKERRQ(ierr); /* owner rank (cells) */ 493552f7358SJed Brown /* all cell data */ 494552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 495552f7358SJed Brown PetscInt bs; 496552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 497fd68c46aSJed Brown ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr); 498552f7358SJed Brown for (i=0; i<bs; i++) { 4990298fd71SBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,PETSC_SCALAR,tag);CHKERRQ(ierr); 500552f7358SJed Brown } 501552f7358SJed Brown } 502552f7358SJed Brown /* all point data */ 503552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 504552f7358SJed Brown PetscInt bs; 505552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 506552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 507552f7358SJed Brown for (i=0; i<bs; i++) { 5080298fd71SBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices,PETSC_SCALAR,tag);CHKERRQ(ierr); 509552f7358SJed Brown } 510552f7358SJed Brown } 511552f7358SJed Brown } 512552f7358SJed Brown } 513552f7358SJed Brown ierr = PetscFree(gpiece);CHKERRQ(ierr); 514552f7358SJed Brown ierr = PetscFree(buffer);CHKERRQ(ierr); 515552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr); 516552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr); 5175e32de72SMatthew G. Knepley ierr = PetscFClose(comm,fp);CHKERRQ(ierr); 518552f7358SJed Brown PetscFunctionReturn(0); 519552f7358SJed Brown } 520