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 18*94fbd55eSBarry Smith static PetscErrorCode TransferWrite(PetscViewer viewer,FILE *fp,PetscMPIInt srank,PetscMPIInt root,const void *send,void *recv,PetscMPIInt count,MPI_Datatype mpidatatype,PetscMPIInt tag) 19552f7358SJed Brown { 20552f7358SJed Brown PetscMPIInt rank; 21552f7358SJed Brown PetscErrorCode ierr; 22ce94432eSBarry Smith MPI_Comm comm; 23552f7358SJed Brown 24552f7358SJed Brown PetscFunctionBegin; 25ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); 26552f7358SJed Brown ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 27552f7358SJed Brown 28552f7358SJed Brown if (rank == srank && rank != root) { 29552f7358SJed Brown ierr = MPI_Send((void*)send,count,mpidatatype,root,tag,comm);CHKERRQ(ierr); 30552f7358SJed Brown } else if (rank == root) { 31552f7358SJed Brown const void *buffer; 32552f7358SJed Brown if (root == srank) { /* self */ 33552f7358SJed Brown buffer = send; 34552f7358SJed Brown } else { 35552f7358SJed Brown MPI_Status status; 36552f7358SJed Brown PetscMPIInt nrecv; 37552f7358SJed Brown ierr = MPI_Recv(recv,count,mpidatatype,srank,tag,comm,&status);CHKERRQ(ierr); 38552f7358SJed Brown ierr = MPI_Get_count(&status,mpidatatype,&nrecv);CHKERRQ(ierr); 39552f7358SJed Brown if (count != nrecv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch"); 40552f7358SJed Brown buffer = recv; 41552f7358SJed Brown } 42*94fbd55eSBarry Smith ierr = PetscViewerVTKFWrite(viewer,fp,buffer,count,mpidatatype);CHKERRQ(ierr); 43552f7358SJed Brown } 44552f7358SJed Brown PetscFunctionReturn(0); 45552f7358SJed Brown } 46552f7358SJed Brown 472e529ec8SStefano Zampini static PetscErrorCode DMPlexGetVTKConnectivity(DM dm, PetscBool localized, PieceInfo *piece,PetscVTKInt **oconn,PetscVTKInt **ooffsets,PetscVTKType **otypes) 48552f7358SJed Brown { 49552f7358SJed Brown PetscErrorCode ierr; 502e529ec8SStefano Zampini PetscSection coordSection; 51552f7358SJed Brown PetscVTKInt *conn,*offsets; 52552f7358SJed Brown PetscVTKType *types; 53552f7358SJed Brown PetscInt dim,vStart,vEnd,cStart,cEnd,pStart,pEnd,cellHeight,cMax,numLabelCells,hasLabel,c,v,countcell,countconn; 54552f7358SJed Brown 55552f7358SJed Brown PetscFunctionBegin; 56dcca6d9dSJed Brown ierr = PetscMalloc3(piece->nconn,&conn,piece->ncells,&offsets,piece->ncells,&types);CHKERRQ(ierr); 572e529ec8SStefano Zampini ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 58c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 59552f7358SJed Brown ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 60552f7358SJed Brown ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 61552f7358SJed Brown ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 62552f7358SJed Brown ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 630298fd71SBarry Smith ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 648865f1eaSKarl Rupp if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 65c58f1c22SToby Isaac ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 66552f7358SJed Brown hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 67552f7358SJed Brown 68552f7358SJed Brown countcell = 0; 69552f7358SJed Brown countconn = 0; 70552f7358SJed Brown for (c = cStart; c < cEnd; ++c) { 7119003fb5SStefano Zampini PetscInt nverts,dof = 0,celltype,startoffset,nC=0; 72552f7358SJed Brown 73552f7358SJed Brown if (hasLabel) { 74552f7358SJed Brown PetscInt value; 75552f7358SJed Brown 76c58f1c22SToby Isaac ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 77552f7358SJed Brown if (value != 1) continue; 78552f7358SJed Brown } 79552f7358SJed Brown startoffset = countconn; 8019003fb5SStefano Zampini if (localized) { 8119003fb5SStefano Zampini ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr); 8219003fb5SStefano Zampini } 8319003fb5SStefano Zampini if (!dof) { 842e529ec8SStefano Zampini PetscInt *closure = NULL; 852e529ec8SStefano Zampini PetscInt closureSize; 862e529ec8SStefano Zampini 87552f7358SJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 88552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 89552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 9019003fb5SStefano Zampini if (!localized) conn[countconn++] = closure[v] - vStart; 9119003fb5SStefano Zampini else conn[countconn++] = startoffset + nC; 92724b5320SMatthew G. Knepley ++nC; 93552f7358SJed Brown } 94552f7358SJed Brown } 95552f7358SJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 962e529ec8SStefano Zampini } else { 972e529ec8SStefano Zampini for (nC = 0; nC < dof/dim; nC++) conn[countconn++] = startoffset + nC; 982e529ec8SStefano Zampini } 99724b5320SMatthew G. Knepley ierr = DMPlexInvertCell(dim, nC, &conn[countconn-nC]);CHKERRQ(ierr); 1008865f1eaSKarl Rupp 101552f7358SJed Brown offsets[countcell] = countconn; 1028865f1eaSKarl Rupp 103552f7358SJed Brown nverts = countconn - startoffset; 104de0a4605SMatthew G. Knepley ierr = DMPlexVTKGetCellType_Internal(dm,dim,nverts,&celltype);CHKERRQ(ierr); 1058865f1eaSKarl Rupp 106552f7358SJed Brown types[countcell] = celltype; 107552f7358SJed Brown countcell++; 108552f7358SJed Brown } 109552f7358SJed Brown if (countcell != piece->ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count"); 110552f7358SJed Brown if (countconn != piece->nconn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count"); 111552f7358SJed Brown *oconn = conn; 112552f7358SJed Brown *ooffsets = offsets; 113552f7358SJed Brown *otypes = types; 114552f7358SJed Brown PetscFunctionReturn(0); 115552f7358SJed Brown } 116552f7358SJed Brown 117552f7358SJed Brown /* 118552f7358SJed Brown Write all fields that have been provided to the viewer 119552f7358SJed Brown Multi-block XML format with binary appended data. 120552f7358SJed Brown */ 121552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer) 122552f7358SJed Brown { 123ce94432eSBarry Smith MPI_Comm comm; 1242e529ec8SStefano Zampini PetscSection coordSection; 125552f7358SJed Brown PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data; 126552f7358SJed Brown PetscViewerVTKObjectLink link; 127552f7358SJed Brown FILE *fp; 128552f7358SJed Brown PetscMPIInt rank,size,tag; 129552f7358SJed Brown PetscErrorCode ierr; 13019003fb5SStefano Zampini PetscInt dimEmbed,cellHeight,cStart,cEnd,vStart,vEnd,cMax,numLabelCells,hasLabel,c,v,r,i; 1312e529ec8SStefano Zampini PetscBool localized; 1320298fd71SBarry Smith PieceInfo piece,*gpiece = NULL; 1330298fd71SBarry Smith void *buffer = NULL; 134552f7358SJed Brown 135552f7358SJed Brown PetscFunctionBegin; 136ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 137552f7358SJed Brown #if defined(PETSC_USE_COMPLEX) 138ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Complex values not supported"); 139552f7358SJed Brown #endif 140552f7358SJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 141552f7358SJed Brown ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 142552f7358SJed Brown tag = ((PetscObject)viewer)->tag; 143552f7358SJed Brown 144552f7358SJed Brown ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr); 145552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr); 146519f805aSKarl Rupp #if defined(PETSC_WORDS_BIGENDIAN) 147552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr); 148552f7358SJed Brown #else 149552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr); 150552f7358SJed Brown #endif 151552f7358SJed Brown ierr = PetscFPrintf(comm,fp," <UnstructuredGrid>\n");CHKERRQ(ierr); 152552f7358SJed Brown 1533baa072aSToby Isaac ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 154552f7358SJed Brown ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 155552f7358SJed Brown ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 156552f7358SJed Brown ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1570298fd71SBarry Smith ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1588865f1eaSKarl Rupp if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 159c58f1c22SToby Isaac ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 1602e529ec8SStefano Zampini ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr); 1612e529ec8SStefano Zampini ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1628865f1eaSKarl Rupp 163552f7358SJed Brown hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 1642e529ec8SStefano Zampini piece.nvertices = 0; 165552f7358SJed Brown piece.ncells = 0; 166552f7358SJed Brown piece.nconn = 0; 1672e529ec8SStefano Zampini if (!localized) piece.nvertices = vEnd - vStart; 168552f7358SJed Brown for (c = cStart; c < cEnd; ++c) { 16919003fb5SStefano Zampini PetscInt dof = 0; 170552f7358SJed Brown if (hasLabel) { 171552f7358SJed Brown PetscInt value; 172552f7358SJed Brown 173c58f1c22SToby Isaac ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 174552f7358SJed Brown if (value != 1) continue; 175552f7358SJed Brown } 17619003fb5SStefano Zampini if (localized) { 17719003fb5SStefano Zampini ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr); 17819003fb5SStefano Zampini } 17919003fb5SStefano Zampini if (!dof) { 1802e529ec8SStefano Zampini PetscInt *closure = NULL; 1812e529ec8SStefano Zampini PetscInt closureSize; 1822e529ec8SStefano Zampini 183552f7358SJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 184552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 1852e529ec8SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 1862e529ec8SStefano Zampini piece.nconn++; 18719003fb5SStefano Zampini if (localized) piece.nvertices++; 1882e529ec8SStefano Zampini } 189552f7358SJed Brown } 190552f7358SJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1912e529ec8SStefano Zampini } else { 1922e529ec8SStefano Zampini piece.nvertices += dof/dimEmbed; 1932e529ec8SStefano Zampini piece.nconn += dof/dimEmbed; 1942e529ec8SStefano Zampini } 195552f7358SJed Brown piece.ncells++; 196552f7358SJed Brown } 197785e854fSJed Brown if (!rank) {ierr = PetscMalloc1(size,&gpiece);CHKERRQ(ierr);} 198552f7358SJed Brown ierr = MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm);CHKERRQ(ierr); 199552f7358SJed Brown 200552f7358SJed Brown /* 201552f7358SJed Brown * Write file header 202552f7358SJed Brown */ 203552f7358SJed Brown if (!rank) { 204552f7358SJed Brown PetscInt boffset = 0; 205552f7358SJed Brown 206552f7358SJed Brown for (r=0; r<size; r++) { 207552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Piece NumberOfPoints=\"%D\" NumberOfCells=\"%D\">\n",gpiece[r].nvertices,gpiece[r].ncells);CHKERRQ(ierr); 208552f7358SJed Brown /* Coordinate positions */ 209552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Points>\n");CHKERRQ(ierr); 210552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr); 211552f7358SJed Brown boffset += gpiece[r].nvertices*3*sizeof(PetscScalar) + sizeof(int); 212552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Points>\n");CHKERRQ(ierr); 213552f7358SJed Brown /* Cell connectivity */ 214552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <Cells>\n");CHKERRQ(ierr); 215552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 216552f7358SJed Brown boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int); 217552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"offsets\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 218552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int); 219552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"UInt8\" Name=\"types\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 220552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int); 221552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Cells>\n");CHKERRQ(ierr); 222552f7358SJed Brown 223552f7358SJed Brown /* 224552f7358SJed Brown * Cell Data headers 225552f7358SJed Brown */ 226552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <CellData>\n");CHKERRQ(ierr); 227552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",boffset);CHKERRQ(ierr); 228552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(int) + sizeof(int); 229552f7358SJed Brown /* all the vectors */ 230552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 231552f7358SJed Brown Vec X = (Vec)link->vec; 2321cfafdd3SJed Brown PetscInt bs,nfields,field; 233552f7358SJed Brown const char *vecname = ""; 234552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 235552f7358SJed 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. */ 236552f7358SJed Brown ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr); 237552f7358SJed Brown } 238552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr); 2391cfafdd3SJed Brown ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr); 2407ded4cbaSJed Brown for (field=0,i=0; field<(nfields?nfields:1); field++) { 2411cfafdd3SJed Brown PetscInt fbs,j; 242a00cdb45SToby Isaac PetscFV fv = NULL; 243a00cdb45SToby Isaac PetscObject f; 244a00cdb45SToby Isaac PetscClassId fClass; 2450298fd71SBarry Smith const char *fieldname = NULL; 2461cfafdd3SJed Brown char buf[256]; 2477ded4cbaSJed Brown if (nfields) { /* We have user-defined fields/components */ 2481cfafdd3SJed Brown ierr = PetscSectionGetFieldDof(dm->defaultSection,cStart,field,&fbs);CHKERRQ(ierr); 2491cfafdd3SJed Brown ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr); 2507ded4cbaSJed Brown } else fbs = bs; /* Say we have one field with 'bs' components */ 251a00cdb45SToby Isaac ierr = DMGetField(dm,field,&f);CHKERRQ(ierr); 252a00cdb45SToby Isaac ierr = PetscObjectGetClassId(f,&fClass);CHKERRQ(ierr); 253a00cdb45SToby Isaac if (fClass == PETSCFV_CLASSID) { 254a00cdb45SToby Isaac fv = (PetscFV) f; 255a00cdb45SToby Isaac } 256552f7358SJed Brown if (!fieldname) { 2571cfafdd3SJed Brown ierr = PetscSNPrintf(buf,sizeof(buf),"CellField%D",field);CHKERRQ(ierr); 258552f7358SJed Brown fieldname = buf; 259552f7358SJed Brown } 2601cfafdd3SJed Brown for (j=0; j<fbs; j++) { 261a00cdb45SToby Isaac const char *compName = NULL; 262a00cdb45SToby Isaac if (fv) { 263a00cdb45SToby Isaac ierr = PetscFVGetComponentName(fv,j,&compName);CHKERRQ(ierr); 264a00cdb45SToby Isaac } 265a00cdb45SToby Isaac if (compName) { 266a00cdb45SToby 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); 267a00cdb45SToby Isaac } 268a00cdb45SToby Isaac else { 2691cfafdd3SJed 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); 270a00cdb45SToby Isaac } 271552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(PetscScalar) + sizeof(int); 2721cfafdd3SJed Brown i++; 273552f7358SJed Brown } 274552f7358SJed Brown } 2751cfafdd3SJed Brown if (i != bs) SETERRQ2(comm,PETSC_ERR_PLIB,"Total number of field components %D != block size %D",i,bs); 2761cfafdd3SJed Brown } 277552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </CellData>\n");CHKERRQ(ierr); 278552f7358SJed Brown 279552f7358SJed Brown /* 280552f7358SJed Brown * Point Data headers 281552f7358SJed Brown */ 282552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," <PointData>\n");CHKERRQ(ierr); 283552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 284552f7358SJed Brown Vec X = (Vec)link->vec; 2851cfafdd3SJed Brown PetscInt bs,nfields,field; 286552f7358SJed Brown const char *vecname = ""; 287552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 288552f7358SJed 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. */ 289552f7358SJed Brown ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr); 290552f7358SJed Brown } 291552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 2921cfafdd3SJed Brown ierr = PetscSectionGetNumFields(dm->defaultSection,&nfields);CHKERRQ(ierr); 2937ded4cbaSJed Brown for (field=0,i=0; field<(nfields?nfields:1); field++) { 2941cfafdd3SJed Brown PetscInt fbs,j; 2951cfafdd3SJed Brown const char *fieldname = NULL; 2961cfafdd3SJed Brown char buf[256]; 2977ded4cbaSJed Brown if (nfields) { /* We have user-defined fields/components */ 2981cfafdd3SJed Brown ierr = PetscSectionGetFieldDof(dm->defaultSection,vStart,field,&fbs);CHKERRQ(ierr); 2991cfafdd3SJed Brown ierr = PetscSectionGetFieldName(dm->defaultSection,field,&fieldname);CHKERRQ(ierr); 3007ded4cbaSJed Brown } else fbs = bs; /* Say we have one field with 'bs' components */ 3011cfafdd3SJed Brown if (!fieldname) { 3021cfafdd3SJed Brown ierr = PetscSNPrintf(buf,sizeof(buf),"PointField%D",field);CHKERRQ(ierr); 3031cfafdd3SJed Brown fieldname = buf; 3041cfafdd3SJed Brown } 3051cfafdd3SJed Brown for (j=0; j<fbs; j++) { 3061cfafdd3SJed 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); 307552f7358SJed Brown boffset += gpiece[r].nvertices*sizeof(PetscScalar) + sizeof(int); 308552f7358SJed Brown } 309552f7358SJed Brown } 3101cfafdd3SJed Brown } 311552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </PointData>\n");CHKERRQ(ierr); 312552f7358SJed Brown 313552f7358SJed Brown ierr = PetscFPrintf(PETSC_COMM_SELF,fp," </Piece>\n");CHKERRQ(ierr); 314552f7358SJed Brown } 315552f7358SJed Brown } 316552f7358SJed Brown 317552f7358SJed Brown ierr = PetscFPrintf(comm,fp," </UnstructuredGrid>\n");CHKERRQ(ierr); 318552f7358SJed Brown ierr = PetscFPrintf(comm,fp," <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr); 319552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr); 320552f7358SJed Brown 321552f7358SJed Brown if (!rank) { 322552f7358SJed Brown PetscInt maxsize = 0; 323552f7358SJed Brown for (r=0; r<size; r++) { 324552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscScalar))); 325552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*sizeof(PetscScalar))); 326552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt))); 327552f7358SJed Brown } 328552f7358SJed Brown ierr = PetscMalloc(maxsize,&buffer);CHKERRQ(ierr); 329552f7358SJed Brown } 330552f7358SJed Brown for (r=0; r<size; r++) { 331552f7358SJed Brown if (r == rank) { 332552f7358SJed Brown PetscInt nsend; 333552f7358SJed Brown { /* Position */ 334552f7358SJed Brown const PetscScalar *x; 3350298fd71SBarry Smith PetscScalar *y = NULL; 336552f7358SJed Brown Vec coords; 3372e529ec8SStefano Zampini 338552f7358SJed Brown ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr); 339552f7358SJed Brown ierr = VecGetArrayRead(coords,&x);CHKERRQ(ierr); 34019003fb5SStefano Zampini if (dimEmbed != 3 || localized) { 341785e854fSJed Brown ierr = PetscMalloc1(piece.nvertices*3,&y);CHKERRQ(ierr); 34219003fb5SStefano Zampini if (localized) { 34319003fb5SStefano Zampini PetscInt cnt; 34419003fb5SStefano Zampini for (c=cStart,cnt=0; c<cEnd; c++) { 34519003fb5SStefano Zampini PetscInt off, dof; 34619003fb5SStefano Zampini 34719003fb5SStefano Zampini ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr); 34819003fb5SStefano Zampini if (!dof) { 34919003fb5SStefano Zampini PetscInt *closure = NULL; 35019003fb5SStefano Zampini PetscInt closureSize; 35119003fb5SStefano Zampini 35219003fb5SStefano Zampini ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 35319003fb5SStefano Zampini for (v = 0; v < closureSize*2; v += 2) { 35419003fb5SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 35519003fb5SStefano Zampini ierr = PetscSectionGetOffset(coordSection, closure[v], &off);CHKERRQ(ierr); 35619003fb5SStefano Zampini if (dimEmbed != 3) { 35719003fb5SStefano Zampini y[cnt*3+0] = x[off+0]; 35819003fb5SStefano Zampini y[cnt*3+1] = (dimEmbed > 1) ? x[off+1] : 0.0; 35919003fb5SStefano Zampini y[cnt*3+2] = 0.0; 36019003fb5SStefano Zampini } else { 36119003fb5SStefano Zampini y[cnt*3+0] = x[off+0]; 36219003fb5SStefano Zampini y[cnt*3+1] = x[off+1]; 36319003fb5SStefano Zampini y[cnt*3+2] = x[off+2]; 36419003fb5SStefano Zampini } 36519003fb5SStefano Zampini cnt++; 36619003fb5SStefano Zampini } 36719003fb5SStefano Zampini } 36819003fb5SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 36919003fb5SStefano Zampini } else { 37019003fb5SStefano Zampini ierr = PetscSectionGetOffset(coordSection, c, &off);CHKERRQ(ierr); 37119003fb5SStefano Zampini if (dimEmbed != 3) { 37219003fb5SStefano Zampini for (i=0; i<dof/dimEmbed; i++) { 37319003fb5SStefano Zampini y[cnt*3+0] = x[off + i*dimEmbed + 0]; 37419003fb5SStefano Zampini y[cnt*3+1] = (dimEmbed > 1) ? x[off + i*dimEmbed + 1] : 0.0; 37519003fb5SStefano Zampini y[cnt*3+2] = 0.0; 37619003fb5SStefano Zampini cnt++; 37719003fb5SStefano Zampini } 37819003fb5SStefano Zampini } else { 37919003fb5SStefano Zampini for (i=0; i<dof; i ++) { 38019003fb5SStefano Zampini y[cnt*3+i] = x[off + i]; 38119003fb5SStefano Zampini } 38219003fb5SStefano Zampini cnt += dof/dimEmbed; 38319003fb5SStefano Zampini } 38419003fb5SStefano Zampini } 38519003fb5SStefano Zampini } 38619003fb5SStefano Zampini if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 38719003fb5SStefano Zampini } else { 388552f7358SJed Brown for (i=0; i<piece.nvertices; i++) { 3893baa072aSToby Isaac y[i*3+0] = x[i*dimEmbed+0]; 3903baa072aSToby Isaac y[i*3+1] = (dimEmbed > 1) ? x[i*dimEmbed+1] : 0; 39119003fb5SStefano Zampini y[i*3+2] = 0.0; 39219003fb5SStefano Zampini } 393552f7358SJed Brown } 394552f7358SJed Brown } 3952e529ec8SStefano Zampini nsend = piece.nvertices*3; 396*94fbd55eSBarry Smith ierr = TransferWrite(viewer,fp,r,0,y ? y : x,buffer,nsend,MPIU_SCALAR,tag);CHKERRQ(ierr); 397552f7358SJed Brown ierr = PetscFree(y);CHKERRQ(ierr); 398552f7358SJed Brown ierr = VecRestoreArrayRead(coords,&x);CHKERRQ(ierr); 399552f7358SJed Brown } 400552f7358SJed Brown { /* Connectivity, offsets, types */ 4013e869605SMatthew G. Knepley PetscVTKInt *connectivity = NULL, *offsets = NULL; 4023e869605SMatthew G. Knepley PetscVTKType *types = NULL; 4032e529ec8SStefano Zampini ierr = DMPlexGetVTKConnectivity(dm,localized,&piece,&connectivity,&offsets,&types);CHKERRQ(ierr); 404*94fbd55eSBarry Smith ierr = TransferWrite(viewer,fp,r,0,connectivity,buffer,piece.nconn,MPI_INT,tag);CHKERRQ(ierr); 405*94fbd55eSBarry Smith ierr = TransferWrite(viewer,fp,r,0,offsets,buffer,piece.ncells,MPI_INT,tag);CHKERRQ(ierr); 406*94fbd55eSBarry Smith ierr = TransferWrite(viewer,fp,r,0,types,buffer,piece.ncells,MPI_CHAR,tag);CHKERRQ(ierr); 407552f7358SJed Brown ierr = PetscFree3(connectivity,offsets,types);CHKERRQ(ierr); 408552f7358SJed Brown } 409552f7358SJed Brown { /* Owners (cell data) */ 410552f7358SJed Brown PetscVTKInt *owners; 411785e854fSJed Brown ierr = PetscMalloc1(piece.ncells,&owners);CHKERRQ(ierr); 412552f7358SJed Brown for (i=0; i<piece.ncells; i++) owners[i] = rank; 413*94fbd55eSBarry Smith ierr = TransferWrite(viewer,fp,r,0,owners,buffer,piece.ncells,MPI_INT,tag);CHKERRQ(ierr); 414552f7358SJed Brown ierr = PetscFree(owners);CHKERRQ(ierr); 415552f7358SJed Brown } 416552f7358SJed Brown /* Cell data */ 417552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 418552f7358SJed Brown Vec X = (Vec)link->vec; 419552f7358SJed Brown const PetscScalar *x; 420552f7358SJed Brown PetscScalar *y; 421552f7358SJed Brown PetscInt bs; 422552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 423552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr); 424552f7358SJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 425785e854fSJed Brown ierr = PetscMalloc1(piece.ncells,&y);CHKERRQ(ierr); 426552f7358SJed Brown for (i=0; i<bs; i++) { 427552f7358SJed Brown PetscInt cnt; 428552f7358SJed Brown for (c=cStart,cnt=0; c<cEnd; c++) { 429640bce14SSatish Balay PetscScalar *xpoint; 430552f7358SJed Brown if (hasLabel) { /* Ignore some cells */ 431552f7358SJed Brown PetscInt value; 432c58f1c22SToby Isaac ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 433552f7358SJed Brown if (value != 1) continue; 434552f7358SJed Brown } 435552f7358SJed Brown ierr = DMPlexPointLocalRead(dm,c,x,&xpoint);CHKERRQ(ierr); 436552f7358SJed Brown y[cnt++] = xpoint[i]; 437552f7358SJed Brown } 438552f7358SJed Brown if (cnt != piece.ncells) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 439*94fbd55eSBarry Smith ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.ncells,MPIU_SCALAR,tag);CHKERRQ(ierr); 440552f7358SJed Brown } 441552f7358SJed Brown ierr = PetscFree(y);CHKERRQ(ierr); 442552f7358SJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 443552f7358SJed Brown } 444552f7358SJed Brown 445552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 446552f7358SJed Brown Vec X = (Vec)link->vec; 447552f7358SJed Brown const PetscScalar *x; 448552f7358SJed Brown PetscScalar *y; 449552f7358SJed Brown PetscInt bs; 450552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 451552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 452552f7358SJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 453785e854fSJed Brown ierr = PetscMalloc1(piece.nvertices,&y);CHKERRQ(ierr); 454552f7358SJed Brown for (i=0; i<bs; i++) { 455552f7358SJed Brown PetscInt cnt; 4562e529ec8SStefano Zampini if (!localized) { 457552f7358SJed Brown for (v=vStart,cnt=0; v<vEnd; v++) { 458640bce14SSatish Balay PetscScalar *xpoint; 45937045cddSJed Brown ierr = DMPlexPointLocalRead(dm,v,x,&xpoint);CHKERRQ(ierr); 460552f7358SJed Brown y[cnt++] = xpoint[i]; 461552f7358SJed Brown } 4622e529ec8SStefano Zampini } else { 4632e529ec8SStefano Zampini for (c=cStart,cnt=0; c<cEnd; c++) { 4642e529ec8SStefano Zampini PetscInt *closure = NULL; 4652e529ec8SStefano Zampini PetscInt closureSize, off; 4662e529ec8SStefano Zampini 4672e529ec8SStefano Zampini ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 4682e529ec8SStefano Zampini for (v = 0, off = 0; v < closureSize*2; v += 2) { 4692e529ec8SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 4702e529ec8SStefano Zampini PetscScalar *xpoint; 4712e529ec8SStefano Zampini 4722e529ec8SStefano Zampini ierr = DMPlexPointLocalRead(dm,v,x,&xpoint);CHKERRQ(ierr); 4732e529ec8SStefano Zampini y[cnt + off++] = xpoint[i]; 4742e529ec8SStefano Zampini } 4752e529ec8SStefano Zampini } 4762e529ec8SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 4772e529ec8SStefano Zampini } 4782e529ec8SStefano Zampini } 479552f7358SJed Brown if (cnt != piece.nvertices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 480*94fbd55eSBarry Smith ierr = TransferWrite(viewer,fp,r,0,y,buffer,piece.nvertices,MPIU_SCALAR,tag);CHKERRQ(ierr); 481552f7358SJed Brown } 482552f7358SJed Brown ierr = PetscFree(y);CHKERRQ(ierr); 483552f7358SJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 484552f7358SJed Brown } 485552f7358SJed Brown } else if (!rank) { 486*94fbd55eSBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,MPIU_SCALAR,tag);CHKERRQ(ierr); /* positions */ 487*94fbd55eSBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nconn,MPI_INT,tag);CHKERRQ(ierr); /* connectivity */ 488*94fbd55eSBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag);CHKERRQ(ierr); /* offsets */ 489*94fbd55eSBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_CHAR,tag);CHKERRQ(ierr); /* types */ 490*94fbd55eSBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag);CHKERRQ(ierr); /* owner rank (cells) */ 491552f7358SJed Brown /* all cell data */ 492552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 493552f7358SJed Brown PetscInt bs; 494552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 495fd68c46aSJed Brown ierr = PetscSectionGetDof(dm->defaultSection,cStart,&bs);CHKERRQ(ierr); 496552f7358SJed Brown for (i=0; i<bs; i++) { 497*94fbd55eSBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPIU_SCALAR,tag);CHKERRQ(ierr); 498552f7358SJed Brown } 499552f7358SJed Brown } 500552f7358SJed Brown /* all point data */ 501552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 502552f7358SJed Brown PetscInt bs; 503552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 504552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,vStart,&bs);CHKERRQ(ierr); 505552f7358SJed Brown for (i=0; i<bs; i++) { 506*94fbd55eSBarry Smith ierr = TransferWrite(viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices,MPIU_SCALAR,tag);CHKERRQ(ierr); 507552f7358SJed Brown } 508552f7358SJed Brown } 509552f7358SJed Brown } 510552f7358SJed Brown } 511552f7358SJed Brown ierr = PetscFree(gpiece);CHKERRQ(ierr); 512552f7358SJed Brown ierr = PetscFree(buffer);CHKERRQ(ierr); 513552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr); 514552f7358SJed Brown ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr); 5155e32de72SMatthew G. Knepley ierr = PetscFClose(comm,fp);CHKERRQ(ierr); 516552f7358SJed Brown PetscFunctionReturn(0); 517552f7358SJed Brown } 518