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 10955d60d1SToby Isaac #if defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL___FP16) 11955d60d1SToby Isaac /* output in float if single or half precision in memory */ 12552f7358SJed Brown static const char precision[] = "Float32"; 13955d60d1SToby Isaac typedef float PetscVTUReal; 14955d60d1SToby Isaac #define MPIU_VTUREAL MPI_FLOAT 15955d60d1SToby Isaac #elif defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128) 16955d60d1SToby Isaac /* output in double if double or quad precision in memory */ 17552f7358SJed Brown static const char precision[] = "Float64"; 18955d60d1SToby Isaac typedef double PetscVTUReal; 19955d60d1SToby Isaac #define MPIU_VTUREAL MPI_DOUBLE 20552f7358SJed Brown #else 21552f7358SJed Brown static const char precision[] = "UnknownPrecision"; 22955d60d1SToby Isaac typedef PetscReal PetscVTUReal; 23955d60d1SToby Isaac #define MPIU_VTUREAL MPIU_REAL 24552f7358SJed Brown #endif 25552f7358SJed Brown 266030a318SStefano Zampini static PetscErrorCode TransferWrite(MPI_Comm comm, PetscViewer viewer,FILE *fp,PetscMPIInt srank,PetscMPIInt root,const void *send,void *recv,PetscMPIInt count,MPI_Datatype mpidatatype,PetscMPIInt tag) 27552f7358SJed Brown { 28552f7358SJed Brown PetscMPIInt rank; 29552f7358SJed Brown 30552f7358SJed Brown PetscFunctionBegin; 319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 32552f7358SJed Brown if (rank == srank && rank != root) { 339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void*)send,count,mpidatatype,root,tag,comm)); 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; 419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(recv,count,mpidatatype,srank,tag,comm,&status)); 429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status,mpidatatype,&nrecv)); 4308401ef6SPierre Jolivet PetscCheck(count == nrecv,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch"); 44552f7358SJed Brown buffer = recv; 45552f7358SJed Brown } 469566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer,fp,buffer,count,mpidatatype)); 47552f7358SJed Brown } 48552f7358SJed Brown PetscFunctionReturn(0); 49552f7358SJed Brown } 50552f7358SJed Brown 512e529ec8SStefano Zampini static PetscErrorCode DMPlexGetVTKConnectivity(DM dm, PetscBool localized, PieceInfo *piece,PetscVTKInt **oconn,PetscVTKInt **ooffsets,PetscVTKType **otypes) 52552f7358SJed Brown { 53*6858538eSMatthew G. Knepley PetscSection coordSection, cellCoordSection; 54552f7358SJed Brown PetscVTKInt *conn,*offsets; 55552f7358SJed Brown PetscVTKType *types; 562f029394SStefano Zampini PetscInt dim,vStart,vEnd,cStart,cEnd,pStart,pEnd,cellHeight,numLabelCells,hasLabel,c,v,countcell,countconn; 57552f7358SJed Brown 58552f7358SJed Brown PetscFunctionBegin; 599566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(piece->nconn,&conn,piece->ncells,&offsets,piece->ncells,&types)); 609566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 61*6858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinateSection(dm, &cellCoordSection)); 629566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm,&dim)); 639566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm,&pStart,&pEnd)); 649566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight)); 659566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd)); 669566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 679566063dSJacob Faibussowitsch PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells)); 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) { 7319003fb5SStefano Zampini PetscInt nverts,dof = 0,celltype,startoffset,nC=0; 74552f7358SJed Brown 75552f7358SJed Brown if (hasLabel) { 76552f7358SJed Brown PetscInt value; 77552f7358SJed Brown 789566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dm, "vtk", c, &value)); 79552f7358SJed Brown if (value != 1) continue; 80552f7358SJed Brown } 81552f7358SJed Brown startoffset = countconn; 82*6858538eSMatthew G. Knepley if (localized) PetscCall(PetscSectionGetDof(cellCoordSection, c, &dof)); 8319003fb5SStefano Zampini if (!dof) { 842e529ec8SStefano Zampini PetscInt *closure = NULL; 852e529ec8SStefano Zampini PetscInt closureSize; 862e529ec8SStefano Zampini 879566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 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 } 959566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 962e529ec8SStefano Zampini } else { 972e529ec8SStefano Zampini for (nC = 0; nC < dof/dim; nC++) conn[countconn++] = startoffset + nC; 982e529ec8SStefano Zampini } 9996ca5757SLisandro Dalcin 10096ca5757SLisandro Dalcin { 10196ca5757SLisandro Dalcin PetscInt n = PetscMin(nC, 8), s = countconn - nC, i, cone[8]; 10296ca5757SLisandro Dalcin for (i = 0; i < n; ++i) cone[i] = conn[s+i]; 1039566063dSJacob Faibussowitsch PetscCall(DMPlexReorderCell(dm, c, cone)); 10496ca5757SLisandro Dalcin for (i = 0; i < n; ++i) conn[s+i] = (int)cone[i]; 10596ca5757SLisandro Dalcin } 1068865f1eaSKarl Rupp 107552f7358SJed Brown offsets[countcell] = countconn; 1088865f1eaSKarl Rupp 109552f7358SJed Brown nverts = countconn - startoffset; 1109566063dSJacob Faibussowitsch PetscCall(DMPlexVTKGetCellType_Internal(dm,dim,nverts,&celltype)); 1118865f1eaSKarl Rupp 112552f7358SJed Brown types[countcell] = celltype; 113552f7358SJed Brown countcell++; 114552f7358SJed Brown } 11508401ef6SPierre Jolivet PetscCheck(countcell == piece->ncells,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent cell count"); 11608401ef6SPierre Jolivet PetscCheck(countconn == piece->nconn,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent connectivity count"); 117552f7358SJed Brown *oconn = conn; 118552f7358SJed Brown *ooffsets = offsets; 119552f7358SJed Brown *otypes = types; 120552f7358SJed Brown PetscFunctionReturn(0); 121552f7358SJed Brown } 122552f7358SJed Brown 123552f7358SJed Brown /* 124552f7358SJed Brown Write all fields that have been provided to the viewer 125552f7358SJed Brown Multi-block XML format with binary appended data. 126552f7358SJed Brown */ 127552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm,PetscViewer viewer) 128552f7358SJed Brown { 129ce94432eSBarry Smith MPI_Comm comm; 130*6858538eSMatthew G. Knepley PetscSection coordSection, cellCoordSection; 131552f7358SJed Brown PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data; 132552f7358SJed Brown PetscViewerVTKObjectLink link; 133552f7358SJed Brown FILE *fp; 134552f7358SJed Brown PetscMPIInt rank,size,tag; 1352f029394SStefano Zampini PetscInt dimEmbed,cellHeight,cStart,cEnd,vStart,vEnd,numLabelCells,hasLabel,c,v,r,i; 1362e529ec8SStefano Zampini PetscBool localized; 1370298fd71SBarry Smith PieceInfo piece,*gpiece = NULL; 1380298fd71SBarry Smith void *buffer = NULL; 13930815ce0SLisandro Dalcin const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian"; 140955d60d1SToby Isaac PetscInt loops_per_scalar; 141552f7358SJed Brown 142552f7358SJed Brown PetscFunctionBegin; 1439566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm,&comm)); 144552f7358SJed Brown #if defined(PETSC_USE_COMPLEX) 145955d60d1SToby Isaac loops_per_scalar = 2; 146955d60d1SToby Isaac #else 147955d60d1SToby Isaac loops_per_scalar = 1; 148552f7358SJed Brown #endif 1499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 1509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 1519566063dSJacob Faibussowitsch PetscCall(PetscCommGetNewTag(comm,&tag)); 152552f7358SJed Brown 1539566063dSJacob Faibussowitsch PetscCall(PetscFOpen(comm,vtk->filename,"wb",&fp)); 1549566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n")); 1559566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"%s\">\n", byte_order)); 1569566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <UnstructuredGrid>\n")); 157552f7358SJed Brown 1589566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dimEmbed)); 1599566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight)); 1609566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd)); 1619566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1629566063dSJacob Faibussowitsch PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells)); 1639566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalized(dm, &localized)); 1649566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 165*6858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinateSection(dm, &cellCoordSection)); 1668865f1eaSKarl Rupp 167552f7358SJed Brown hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 1682e529ec8SStefano Zampini piece.nvertices = 0; 169552f7358SJed Brown piece.ncells = 0; 170552f7358SJed Brown piece.nconn = 0; 1712e529ec8SStefano Zampini if (!localized) piece.nvertices = vEnd - vStart; 172552f7358SJed Brown for (c = cStart; c < cEnd; ++c) { 17319003fb5SStefano Zampini PetscInt dof = 0; 174552f7358SJed Brown if (hasLabel) { 175552f7358SJed Brown PetscInt value; 176552f7358SJed Brown 1779566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dm, "vtk", c, &value)); 178552f7358SJed Brown if (value != 1) continue; 179552f7358SJed Brown } 180*6858538eSMatthew G. Knepley if (localized) PetscCall(PetscSectionGetDof(cellCoordSection, c, &dof)); 18119003fb5SStefano Zampini if (!dof) { 1822e529ec8SStefano Zampini PetscInt *closure = NULL; 1832e529ec8SStefano Zampini PetscInt closureSize; 1842e529ec8SStefano Zampini 1859566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 186552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 1872e529ec8SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 1882e529ec8SStefano Zampini piece.nconn++; 18919003fb5SStefano Zampini if (localized) piece.nvertices++; 1902e529ec8SStefano Zampini } 191552f7358SJed Brown } 1929566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 1932e529ec8SStefano Zampini } else { 1942e529ec8SStefano Zampini piece.nvertices += dof/dimEmbed; 1952e529ec8SStefano Zampini piece.nconn += dof/dimEmbed; 1962e529ec8SStefano Zampini } 197552f7358SJed Brown piece.ncells++; 198552f7358SJed Brown } 1999566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc1(size,&gpiece)); 2009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gather((PetscInt*)&piece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,(PetscInt*)gpiece,sizeof(piece)/sizeof(PetscInt),MPIU_INT,0,comm)); 201552f7358SJed Brown 202552f7358SJed Brown /* 203552f7358SJed Brown * Write file header 204552f7358SJed Brown */ 205dd400576SPatrick Sanan if (rank == 0) { 206552f7358SJed Brown PetscInt boffset = 0; 207552f7358SJed Brown 208552f7358SJed Brown for (r=0; r<size; r++) { 20963a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp," <Piece NumberOfPoints=\"%" PetscInt_FMT "\" NumberOfCells=\"%" PetscInt_FMT "\">\n",gpiece[r].nvertices,gpiece[r].ncells)); 210552f7358SJed Brown /* Coordinate positions */ 2119566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp," <Points>\n")); 21263a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,boffset)); 213955d60d1SToby Isaac boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int); 2149566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp," </Points>\n")); 215552f7358SJed Brown /* Cell connectivity */ 2169566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp," <Cells>\n")); 21763a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",boffset)); 218552f7358SJed Brown boffset += gpiece[r].nconn*sizeof(PetscInt) + sizeof(int); 21963a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"offsets\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",boffset)); 220552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(PetscInt) + sizeof(int); 22163a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"UInt8\" Name=\"types\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",boffset)); 222552f7358SJed Brown boffset += gpiece[r].ncells*sizeof(unsigned char) + sizeof(int); 2239566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp," </Cells>\n")); 224552f7358SJed Brown 225552f7358SJed Brown /* 226552f7358SJed Brown * Cell Data headers 227552f7358SJed Brown */ 2289566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp," <CellData>\n")); 22963a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp," <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",boffset)); 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; 234e630c359SToby Isaac DM dmX = NULL; 2356030a318SStefano Zampini PetscInt bs = 1,nfields,field; 236552f7358SJed Brown const char *vecname = ""; 237e630c359SToby Isaac PetscSection section; 238552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 239552f7358SJed 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. */ 2409566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)X,&vecname)); 241552f7358SJed Brown } 2429566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 243e630c359SToby Isaac if (!dmX) dmX = dm; 2449566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject*) §ion)); 2459566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 2469566063dSJacob Faibussowitsch if (cEnd > cStart) PetscCall(PetscSectionGetDof(section,cStart,&bs)); 2479566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section,&nfields)); 248e630c359SToby Isaac field = 0; 249e630c359SToby Isaac if (link->field >= 0) { 250e630c359SToby Isaac field = link->field; 251e630c359SToby Isaac nfields = field + 1; 252e630c359SToby Isaac } 253e630c359SToby Isaac for (i=0; field<(nfields?nfields:1); field++) { 2541cfafdd3SJed Brown PetscInt fbs,j; 255a00cdb45SToby Isaac PetscFV fv = NULL; 256a00cdb45SToby Isaac PetscObject f; 257a00cdb45SToby Isaac PetscClassId fClass; 2580298fd71SBarry Smith const char *fieldname = NULL; 2591cfafdd3SJed Brown char buf[256]; 260e630c359SToby Isaac PetscBool vector; 2617ded4cbaSJed Brown if (nfields) { /* We have user-defined fields/components */ 2629566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section,cStart,field,&fbs)); 2639566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(section,field,&fieldname)); 2647ded4cbaSJed Brown } else fbs = bs; /* Say we have one field with 'bs' components */ 2659566063dSJacob Faibussowitsch PetscCall(DMGetField(dmX,field,NULL,&f)); 2669566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(f,&fClass)); 267a00cdb45SToby Isaac if (fClass == PETSCFV_CLASSID) { 268a00cdb45SToby Isaac fv = (PetscFV) f; 269a00cdb45SToby Isaac } 270e630c359SToby Isaac if (nfields && !fieldname) { 27163a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf,sizeof(buf),"CellField%" PetscInt_FMT,field)); 272552f7358SJed Brown fieldname = buf; 273552f7358SJed Brown } 274e630c359SToby Isaac vector = PETSC_FALSE; 275e630c359SToby Isaac if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) { 276e630c359SToby Isaac vector = PETSC_TRUE; 27763a3b9bcSJacob Faibussowitsch PetscCheck(fbs <= 3,PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_SIZ,"Cell vector fields can have at most 3 components, %" PetscInt_FMT " given", fbs); 278e630c359SToby Isaac for (j = 0; j < fbs; j++) { 279e630c359SToby Isaac const char *compName = NULL; 280e630c359SToby Isaac if (fv) { 2819566063dSJacob Faibussowitsch PetscCall(PetscFVGetComponentName(fv,j,&compName)); 282e630c359SToby Isaac if (compName) break; 283e630c359SToby Isaac } 284e630c359SToby Isaac } 285e630c359SToby Isaac if (j < fbs) vector = PETSC_FALSE; 286e630c359SToby Isaac } 287e630c359SToby Isaac if (vector) { 288955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 28963a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s%s.Re\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,vecname,fieldname,boffset)); 290955d60d1SToby Isaac boffset += gpiece[r].ncells*3*sizeof(PetscVTUReal) + sizeof(int); 29163a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s%s.Im\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,vecname,fieldname,boffset)); 292955d60d1SToby Isaac boffset += gpiece[r].ncells*3*sizeof(PetscVTUReal) + sizeof(int); 293955d60d1SToby Isaac #else 29463a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,vecname,fieldname,boffset)); 295955d60d1SToby Isaac boffset += gpiece[r].ncells*3*sizeof(PetscVTUReal) + sizeof(int); 296955d60d1SToby Isaac #endif 297e630c359SToby Isaac i+=fbs; 298e630c359SToby Isaac } else { 2991cfafdd3SJed Brown for (j=0; j<fbs; j++) { 300a00cdb45SToby Isaac const char *compName = NULL; 301955d60d1SToby Isaac char finalname[256]; 302a00cdb45SToby Isaac if (fv) { 3039566063dSJacob Faibussowitsch PetscCall(PetscFVGetComponentName(fv,j,&compName)); 304a00cdb45SToby Isaac } 305a00cdb45SToby Isaac if (compName) { 3069566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname,255,"%s%s.%s",vecname,fieldname,compName)); 307a00cdb45SToby Isaac } 308e630c359SToby Isaac else if (fbs > 1) { 30963a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname,255,"%s%s.%" PetscInt_FMT,vecname,fieldname,j)); 310e630c359SToby Isaac } else { 3119566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname,255,"%s%s",vecname,fieldname)); 312a00cdb45SToby Isaac } 313955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 31463a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,finalname,boffset)); 315955d60d1SToby Isaac boffset += gpiece[r].ncells*sizeof(PetscVTUReal) + sizeof(int); 31663a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,finalname,boffset)); 317955d60d1SToby Isaac boffset += gpiece[r].ncells*sizeof(PetscVTUReal) + sizeof(int); 318955d60d1SToby Isaac #else 31963a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,finalname,boffset)); 320955d60d1SToby Isaac boffset += gpiece[r].ncells*sizeof(PetscVTUReal) + sizeof(int); 321955d60d1SToby Isaac #endif 3221cfafdd3SJed Brown i++; 323552f7358SJed Brown } 324552f7358SJed Brown } 325e630c359SToby Isaac } 32663a3b9bcSJacob Faibussowitsch //PetscCheck(i == bs,comm,PETSC_ERR_PLIB,"Total number of field components %" PetscInt_FMT " != block size %" PetscInt_FMT,i,bs); 3271cfafdd3SJed Brown } 3289566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp," </CellData>\n")); 329552f7358SJed Brown 330552f7358SJed Brown /* 331552f7358SJed Brown * Point Data headers 332552f7358SJed Brown */ 3339566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp," <PointData>\n")); 334552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 335552f7358SJed Brown Vec X = (Vec)link->vec; 336e630c359SToby Isaac DM dmX; 3376030a318SStefano Zampini PetscInt bs = 1,nfields,field; 338552f7358SJed Brown const char *vecname = ""; 339e630c359SToby Isaac PetscSection section; 340552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 341552f7358SJed 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. */ 3429566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)X,&vecname)); 343552f7358SJed Brown } 3449566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 345e630c359SToby Isaac if (!dmX) dmX = dm; 3469566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject*) §ion)); 3479566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 3489566063dSJacob Faibussowitsch if (vEnd > vStart) PetscCall(PetscSectionGetDof(section,vStart,&bs)); 3499566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section,&nfields)); 350e630c359SToby Isaac field = 0; 351e630c359SToby Isaac if (link->field >= 0) { 352e630c359SToby Isaac field = link->field; 353e630c359SToby Isaac nfields = field + 1; 354e630c359SToby Isaac } 355e630c359SToby Isaac for (i=0; field<(nfields?nfields:1); field++) { 3561cfafdd3SJed Brown PetscInt fbs,j; 3571cfafdd3SJed Brown const char *fieldname = NULL; 3581cfafdd3SJed Brown char buf[256]; 3597ded4cbaSJed Brown if (nfields) { /* We have user-defined fields/components */ 3609566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section,vStart,field,&fbs)); 3619566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(section,field,&fieldname)); 3627ded4cbaSJed Brown } else fbs = bs; /* Say we have one field with 'bs' components */ 363e630c359SToby Isaac if (nfields && !fieldname) { 36463a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf,sizeof(buf),"PointField%" PetscInt_FMT,field)); 3651cfafdd3SJed Brown fieldname = buf; 3661cfafdd3SJed Brown } 367e630c359SToby Isaac if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) { 36863a3b9bcSJacob Faibussowitsch PetscCheck(fbs <= 3,PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_SIZ,"Point vector fields can have at most 3 components, %" PetscInt_FMT " given", fbs); 369955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 37063a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s%s.Re\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,vecname,fieldname,boffset)); 371955d60d1SToby Isaac boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int); 37263a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s%s.Im\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,vecname,fieldname,boffset)); 373955d60d1SToby Isaac boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int); 374955d60d1SToby Isaac #else 37563a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,vecname,fieldname,boffset)); 376955d60d1SToby Isaac boffset += gpiece[r].nvertices*3*sizeof(PetscVTUReal) + sizeof(int); 377955d60d1SToby Isaac #endif 378e630c359SToby Isaac } else { 3791cfafdd3SJed Brown for (j=0; j<fbs; j++) { 380b778fa18SValeria Barra const char *compName = NULL; 381955d60d1SToby Isaac char finalname[256]; 3829566063dSJacob Faibussowitsch PetscCall(PetscSectionGetComponentName(section,field,j,&compName)); 3839566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname,255,"%s%s.%s",vecname,fieldname,compName)); 384955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 38563a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,finalname,boffset)); 386955d60d1SToby Isaac boffset += gpiece[r].nvertices*sizeof(PetscVTUReal) + sizeof(int); 38763a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,finalname,boffset)); 388955d60d1SToby Isaac boffset += gpiece[r].nvertices*sizeof(PetscVTUReal) + sizeof(int); 389955d60d1SToby Isaac #else 39063a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,finalname,boffset)); 391955d60d1SToby Isaac boffset += gpiece[r].nvertices*sizeof(PetscVTUReal) + sizeof(int); 392955d60d1SToby Isaac #endif 393552f7358SJed Brown } 394552f7358SJed Brown } 3951cfafdd3SJed Brown } 396e630c359SToby Isaac } 3979566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp," </PointData>\n")); 3989566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fp," </Piece>\n")); 399552f7358SJed Brown } 400552f7358SJed Brown } 401552f7358SJed Brown 4029566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," </UnstructuredGrid>\n")); 4039566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <AppendedData encoding=\"raw\">\n")); 4049566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"_")); 405552f7358SJed Brown 406dd400576SPatrick Sanan if (rank == 0) { 407552f7358SJed Brown PetscInt maxsize = 0; 408552f7358SJed Brown for (r=0; r<size; r++) { 409955d60d1SToby Isaac maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nvertices*3*sizeof(PetscVTUReal))); 410955d60d1SToby Isaac maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].ncells*3*sizeof(PetscVTUReal))); 411552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt) (gpiece[r].nconn*sizeof(PetscVTKInt))); 412552f7358SJed Brown } 4139566063dSJacob Faibussowitsch PetscCall(PetscMalloc(maxsize,&buffer)); 414552f7358SJed Brown } 415552f7358SJed Brown for (r=0; r<size; r++) { 416552f7358SJed Brown if (r == rank) { 417552f7358SJed Brown PetscInt nsend; 418552f7358SJed Brown { /* Position */ 419*6858538eSMatthew G. Knepley const PetscScalar *x, *cx = NULL; 420955d60d1SToby Isaac PetscVTUReal *y = NULL; 421*6858538eSMatthew G. Knepley Vec coords, cellCoords; 422955d60d1SToby Isaac PetscBool copy; 4232e529ec8SStefano Zampini 4249566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm,&coords)); 4259566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coords,&x)); 426*6858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinatesLocal(dm,&cellCoords)); 427*6858538eSMatthew G. Knepley if (cellCoords) PetscCall(VecGetArrayRead(cellCoords,&cx)); 428955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 429955d60d1SToby Isaac copy = PETSC_TRUE; 430955d60d1SToby Isaac #else 431955d60d1SToby Isaac copy = (PetscBool) (dimEmbed != 3 || localized || (sizeof(PetscReal) != sizeof(PetscVTUReal))); 432955d60d1SToby Isaac #endif 433955d60d1SToby Isaac if (copy) { 4349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.nvertices*3,&y)); 43519003fb5SStefano Zampini if (localized) { 43619003fb5SStefano Zampini PetscInt cnt; 43719003fb5SStefano Zampini for (c=cStart,cnt=0; c<cEnd; c++) { 43819003fb5SStefano Zampini PetscInt off, dof; 43919003fb5SStefano Zampini 440*6858538eSMatthew G. Knepley PetscCall(PetscSectionGetDof(cellCoordSection, c, &dof)); 44119003fb5SStefano Zampini if (!dof) { 44219003fb5SStefano Zampini PetscInt *closure = NULL; 44319003fb5SStefano Zampini PetscInt closureSize; 44419003fb5SStefano Zampini 4459566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 44619003fb5SStefano Zampini for (v = 0; v < closureSize*2; v += 2) { 44719003fb5SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 4489566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, closure[v], &off)); 44919003fb5SStefano Zampini if (dimEmbed != 3) { 450955d60d1SToby Isaac y[cnt*3+0] = (PetscVTUReal) PetscRealPart(x[off+0]); 451955d60d1SToby Isaac y[cnt*3+1] = (PetscVTUReal) ((dimEmbed > 1) ? PetscRealPart(x[off+1]) : 0.0); 452955d60d1SToby Isaac y[cnt*3+2] = (PetscVTUReal) 0.0; 45319003fb5SStefano Zampini } else { 454955d60d1SToby Isaac y[cnt*3+0] = (PetscVTUReal) PetscRealPart(x[off+0]); 455955d60d1SToby Isaac y[cnt*3+1] = (PetscVTUReal) PetscRealPart(x[off+1]); 456955d60d1SToby Isaac y[cnt*3+2] = (PetscVTUReal) PetscRealPart(x[off+2]); 45719003fb5SStefano Zampini } 45819003fb5SStefano Zampini cnt++; 45919003fb5SStefano Zampini } 46019003fb5SStefano Zampini } 4619566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 46219003fb5SStefano Zampini } else { 463*6858538eSMatthew G. Knepley PetscCall(PetscSectionGetOffset(cellCoordSection, c, &off)); 46419003fb5SStefano Zampini if (dimEmbed != 3) { 46519003fb5SStefano Zampini for (i=0; i<dof/dimEmbed; i++) { 466*6858538eSMatthew G. Knepley y[cnt*3+0] = (PetscVTUReal) PetscRealPart(cx[off + i*dimEmbed + 0]); 467*6858538eSMatthew G. Knepley y[cnt*3+1] = (PetscVTUReal) ((dimEmbed > 1) ? PetscRealPart(cx[off + i*dimEmbed + 1]) : 0.0); 468955d60d1SToby Isaac y[cnt*3+2] = (PetscVTUReal) 0.0; 46919003fb5SStefano Zampini cnt++; 47019003fb5SStefano Zampini } 47119003fb5SStefano Zampini } else { 47219003fb5SStefano Zampini for (i=0; i<dof; i ++) { 473*6858538eSMatthew G. Knepley y[cnt*3+i] = (PetscVTUReal) PetscRealPart(cx[off + i]); 47419003fb5SStefano Zampini } 47519003fb5SStefano Zampini cnt += dof/dimEmbed; 47619003fb5SStefano Zampini } 47719003fb5SStefano Zampini } 47819003fb5SStefano Zampini } 47908401ef6SPierre Jolivet PetscCheck(cnt == piece.nvertices,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 48019003fb5SStefano Zampini } else { 481552f7358SJed Brown for (i=0; i<piece.nvertices; i++) { 482955d60d1SToby Isaac y[i*3+0] = (PetscVTUReal) PetscRealPart(x[i*dimEmbed+0]); 483955d60d1SToby Isaac y[i*3+1] = (PetscVTUReal) ((dimEmbed > 1) ? PetscRealPart(x[i*dimEmbed+1]) : 0.); 484955d60d1SToby Isaac y[i*3+2] = (PetscVTUReal) ((dimEmbed > 2) ? PetscRealPart(x[i*dimEmbed+2]) : 0.); 48519003fb5SStefano Zampini } 486552f7358SJed Brown } 487552f7358SJed Brown } 4882e529ec8SStefano Zampini nsend = piece.nvertices*3; 4899566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,copy ? (const void *) y : (const void *) x,buffer,nsend,MPIU_VTUREAL,tag)); 4909566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 4919566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coords,&x)); 492*6858538eSMatthew G. Knepley if (cellCoords) PetscCall(VecRestoreArrayRead(cellCoords,&cx)); 493552f7358SJed Brown } 494552f7358SJed Brown { /* Connectivity, offsets, types */ 4953e869605SMatthew G. Knepley PetscVTKInt *connectivity = NULL, *offsets = NULL; 4963e869605SMatthew G. Knepley PetscVTKType *types = NULL; 4979566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKConnectivity(dm,localized,&piece,&connectivity,&offsets,&types)); 4989566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,connectivity,buffer,piece.nconn,MPI_INT,tag)); 4999566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,offsets,buffer,piece.ncells,MPI_INT,tag)); 5009566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,types,buffer,piece.ncells,MPI_CHAR,tag)); 5019566063dSJacob Faibussowitsch PetscCall(PetscFree3(connectivity,offsets,types)); 502552f7358SJed Brown } 503552f7358SJed Brown { /* Owners (cell data) */ 504552f7358SJed Brown PetscVTKInt *owners; 5059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.ncells,&owners)); 506552f7358SJed Brown for (i=0; i<piece.ncells; i++) owners[i] = rank; 5079566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,owners,buffer,piece.ncells,MPI_INT,tag)); 5089566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 509552f7358SJed Brown } 510552f7358SJed Brown /* Cell data */ 511552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 512552f7358SJed Brown Vec X = (Vec)link->vec; 513e630c359SToby Isaac DM dmX; 514552f7358SJed Brown const PetscScalar *x; 515955d60d1SToby Isaac PetscVTUReal *y; 5166030a318SStefano Zampini PetscInt bs = 1, nfields, field; 517e630c359SToby Isaac PetscSection section = NULL; 518e630c359SToby Isaac 519552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 5209566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 521e630c359SToby Isaac if (!dmX) dmX = dm; 5229566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject*) §ion)); 5239566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 5249566063dSJacob Faibussowitsch if (cEnd > cStart) PetscCall(PetscSectionGetDof(section,cStart,&bs)); 5259566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section,&nfields)); 526e630c359SToby Isaac field = 0; 527e630c359SToby Isaac if (link->field >= 0) { 528e630c359SToby Isaac field = link->field; 529e630c359SToby Isaac nfields = field + 1; 530e630c359SToby Isaac } 5319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 5329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.ncells*3,&y)); 533e630c359SToby Isaac for (i=0; field<(nfields?nfields:1); field++) { 534e630c359SToby Isaac PetscInt fbs,j; 535e630c359SToby Isaac PetscFV fv = NULL; 536e630c359SToby Isaac PetscObject f; 537e630c359SToby Isaac PetscClassId fClass; 538e630c359SToby Isaac PetscBool vector; 5396030a318SStefano Zampini if (nfields && cEnd > cStart) { /* We have user-defined fields/components */ 5409566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section,cStart,field,&fbs)); 541e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 5429566063dSJacob Faibussowitsch PetscCall(DMGetField(dmX,field,NULL,&f)); 5439566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(f,&fClass)); 544e630c359SToby Isaac if (fClass == PETSCFV_CLASSID) { 545e630c359SToby Isaac fv = (PetscFV) f; 546e630c359SToby Isaac } 547e630c359SToby Isaac vector = PETSC_FALSE; 548e630c359SToby Isaac if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) { 549e630c359SToby Isaac vector = PETSC_TRUE; 550e630c359SToby Isaac for (j = 0; j < fbs; j++) { 551e630c359SToby Isaac const char *compName = NULL; 552e630c359SToby Isaac if (fv) { 5539566063dSJacob Faibussowitsch PetscCall(PetscFVGetComponentName(fv,j,&compName)); 554e630c359SToby Isaac if (compName) break; 555e630c359SToby Isaac } 556e630c359SToby Isaac } 557e630c359SToby Isaac if (j < fbs) vector = PETSC_FALSE; 558e630c359SToby Isaac } 559e630c359SToby Isaac if (vector) { 560955d60d1SToby Isaac PetscInt cnt, l; 561955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 562552f7358SJed Brown for (c=cStart,cnt=0; c<cEnd; c++) { 563e630c359SToby Isaac const PetscScalar *xpoint; 564e630c359SToby Isaac PetscInt off, j; 565e630c359SToby Isaac 566552f7358SJed Brown if (hasLabel) { /* Ignore some cells */ 567552f7358SJed Brown PetscInt value; 5689566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dmX, "vtk", c, &value)); 569552f7358SJed Brown if (value != 1) continue; 570552f7358SJed Brown } 571e630c359SToby Isaac if (nfields) { 5729566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, c, field, &off)); 573e630c359SToby Isaac } else { 5749566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, c, &off)); 575e630c359SToby Isaac } 576e630c359SToby Isaac xpoint = &x[off]; 577e630c359SToby Isaac for (j = 0; j < fbs; j++) { 578955d60d1SToby Isaac y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j])); 579e630c359SToby Isaac } 580e630c359SToby Isaac for (; j < 3; j++) y[cnt++] = 0.; 581e630c359SToby Isaac } 58208401ef6SPierre Jolivet PetscCheck(cnt == piece.ncells*3,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 5839566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,y,buffer,piece.ncells*3,MPIU_VTUREAL,tag)); 584955d60d1SToby Isaac } 585e630c359SToby Isaac } else { 586e630c359SToby Isaac for (i=0; i<fbs; i++) { 587955d60d1SToby Isaac PetscInt cnt, l; 588955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 589e630c359SToby Isaac for (c=cStart,cnt=0; c<cEnd; c++) { 590e630c359SToby Isaac const PetscScalar *xpoint; 591e630c359SToby Isaac PetscInt off; 592e630c359SToby Isaac 593e630c359SToby Isaac if (hasLabel) { /* Ignore some cells */ 594e630c359SToby Isaac PetscInt value; 5959566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dmX, "vtk", c, &value)); 596e630c359SToby Isaac if (value != 1) continue; 597e630c359SToby Isaac } 598e630c359SToby Isaac if (nfields) { 5999566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, c, field, &off)); 600e630c359SToby Isaac } else { 6019566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, c, &off)); 602e630c359SToby Isaac } 603e630c359SToby Isaac xpoint = &x[off]; 604955d60d1SToby Isaac y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i])); 605552f7358SJed Brown } 60608401ef6SPierre Jolivet PetscCheck(cnt == piece.ncells,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 6079566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,y,buffer,piece.ncells,MPIU_VTUREAL,tag)); 608955d60d1SToby Isaac } 609552f7358SJed Brown } 610e630c359SToby Isaac } 611e630c359SToby Isaac } 6129566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 6139566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 614552f7358SJed Brown } 615e630c359SToby Isaac /* point data */ 616552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 617552f7358SJed Brown Vec X = (Vec)link->vec; 618e630c359SToby Isaac DM dmX; 619552f7358SJed Brown const PetscScalar *x; 620955d60d1SToby Isaac PetscVTUReal *y; 6216030a318SStefano Zampini PetscInt bs = 1, nfields, field; 622e630c359SToby Isaac PetscSection section = NULL; 623e630c359SToby Isaac 624552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 6259566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 626e630c359SToby Isaac if (!dmX) dmX = dm; 6279566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject*) §ion)); 6289566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 6299566063dSJacob Faibussowitsch if (vEnd > vStart) PetscCall(PetscSectionGetDof(section,vStart,&bs)); 6309566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section,&nfields)); 631e630c359SToby Isaac field = 0; 632e630c359SToby Isaac if (link->field >= 0) { 633e630c359SToby Isaac field = link->field; 634e630c359SToby Isaac nfields = field + 1; 635e630c359SToby Isaac } 6369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 6379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.nvertices*3,&y)); 638e630c359SToby Isaac for (i=0; field<(nfields?nfields:1); field++) { 639e630c359SToby Isaac PetscInt fbs,j; 6406030a318SStefano Zampini if (nfields && vEnd > vStart) { /* We have user-defined fields/components */ 6419566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section,vStart,field,&fbs)); 642e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 643e630c359SToby Isaac if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) { 644955d60d1SToby Isaac PetscInt cnt, l; 645955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 6462e529ec8SStefano Zampini if (!localized) { 647552f7358SJed Brown for (v=vStart,cnt=0; v<vEnd; v++) { 648e630c359SToby Isaac PetscInt off; 649e630c359SToby Isaac const PetscScalar *xpoint; 650e630c359SToby Isaac 651e630c359SToby Isaac if (nfields) { 6529566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section,v,field,&off)); 653e630c359SToby Isaac } else { 6549566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section,v,&off)); 655e630c359SToby Isaac } 656e630c359SToby Isaac xpoint = &x[off]; 657e630c359SToby Isaac for (j = 0; j < fbs; j++) { 658955d60d1SToby Isaac y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j])); 659e630c359SToby Isaac } 660e630c359SToby Isaac for (; j < 3; j++) y[cnt++] = 0.; 661e630c359SToby Isaac } 662e630c359SToby Isaac } else { 663e630c359SToby Isaac for (c=cStart,cnt=0; c<cEnd; c++) { 664e630c359SToby Isaac PetscInt *closure = NULL; 665e630c359SToby Isaac PetscInt closureSize, off; 666e630c359SToby Isaac 6679566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 668e630c359SToby Isaac for (v = 0, off = 0; v < closureSize*2; v += 2) { 669e630c359SToby Isaac if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 670e630c359SToby Isaac PetscInt voff; 671e630c359SToby Isaac const PetscScalar *xpoint; 672e630c359SToby Isaac 673e630c359SToby Isaac if (nfields) { 6749566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section,closure[v],field,&voff)); 675e630c359SToby Isaac } else { 6769566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section,closure[v],&voff)); 677e630c359SToby Isaac } 678e630c359SToby Isaac xpoint = &x[voff]; 679e630c359SToby Isaac for (j = 0; j < fbs; j++) { 680955d60d1SToby Isaac y[cnt + off++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j])); 681e630c359SToby Isaac } 682e630c359SToby Isaac for (; j < 3; j++) y[cnt + off++] = 0.; 683e630c359SToby Isaac } 684e630c359SToby Isaac } 685e630c359SToby Isaac cnt += off; 6869566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 687e630c359SToby Isaac } 688e630c359SToby Isaac } 68908401ef6SPierre Jolivet PetscCheck(cnt == piece.nvertices*3,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 6909566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,y,buffer,piece.nvertices*3,MPIU_VTUREAL,tag)); 691955d60d1SToby Isaac } 692e630c359SToby Isaac } else { 693e630c359SToby Isaac for (i=0; i<fbs; i++) { 694955d60d1SToby Isaac PetscInt cnt, l; 695955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 696e630c359SToby Isaac if (!localized) { 697e630c359SToby Isaac for (v=vStart,cnt=0; v<vEnd; v++) { 698e630c359SToby Isaac PetscInt off; 699e630c359SToby Isaac const PetscScalar *xpoint; 700e630c359SToby Isaac 701e630c359SToby Isaac if (nfields) { 7029566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section,v,field,&off)); 703e630c359SToby Isaac } else { 7049566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section,v,&off)); 705e630c359SToby Isaac } 706e630c359SToby Isaac xpoint = &x[off]; 707955d60d1SToby Isaac y[cnt++] = (PetscVTUReal) (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i])); 708552f7358SJed Brown } 7092e529ec8SStefano Zampini } else { 7102e529ec8SStefano Zampini for (c=cStart,cnt=0; c<cEnd; c++) { 7112e529ec8SStefano Zampini PetscInt *closure = NULL; 7122e529ec8SStefano Zampini PetscInt closureSize, off; 7132e529ec8SStefano Zampini 7149566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 7152e529ec8SStefano Zampini for (v = 0, off = 0; v < closureSize*2; v += 2) { 7162e529ec8SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 717e630c359SToby Isaac PetscInt voff; 718e630c359SToby Isaac const PetscScalar *xpoint; 7192e529ec8SStefano Zampini 720e630c359SToby Isaac if (nfields) { 7219566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section,closure[v],field,&voff)); 722e630c359SToby Isaac } else { 7239566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section,closure[v],&voff)); 724e630c359SToby Isaac } 725e630c359SToby Isaac xpoint = &x[voff]; 726955d60d1SToby Isaac y[cnt + off++] = (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i])); 7272e529ec8SStefano Zampini } 7282e529ec8SStefano Zampini } 7299bda96f6SStefano Zampini cnt += off; 7309566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 7312e529ec8SStefano Zampini } 7322e529ec8SStefano Zampini } 73308401ef6SPierre Jolivet PetscCheck(cnt == piece.nvertices,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not match"); 7349566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,y,buffer,piece.nvertices,MPIU_VTUREAL,tag)); 735955d60d1SToby Isaac } 736552f7358SJed Brown } 737e630c359SToby Isaac } 738e630c359SToby Isaac } 7399566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 7409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 741552f7358SJed Brown } 742dd400576SPatrick Sanan } else if (rank == 0) { 743955d60d1SToby Isaac PetscInt l; 744955d60d1SToby Isaac 7459566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,MPIU_VTUREAL,tag)); /* positions */ 7469566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].nconn,MPI_INT,tag)); /* connectivity */ 7479566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag)); /* offsets */ 7489566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_CHAR,tag)); /* types */ 7499566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPI_INT,tag)); /* owner rank (cells) */ 750552f7358SJed Brown /* all cell data */ 751552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 752e630c359SToby Isaac Vec X = (Vec)link->vec; 7536030a318SStefano Zampini PetscInt bs = 1, nfields, field; 754e630c359SToby Isaac DM dmX; 755e630c359SToby Isaac PetscSection section = NULL; 756e630c359SToby Isaac 757552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 7589566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 759e630c359SToby Isaac if (!dmX) dmX = dm; 7609566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject*) §ion)); 7619566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 7629566063dSJacob Faibussowitsch if (cEnd > cStart) PetscCall(PetscSectionGetDof(section,cStart,&bs)); 7639566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section,&nfields)); 764e630c359SToby Isaac field = 0; 765e630c359SToby Isaac if (link->field >= 0) { 766e630c359SToby Isaac field = link->field; 767e630c359SToby Isaac nfields = field + 1; 768e630c359SToby Isaac } 769e630c359SToby Isaac for (i=0; field<(nfields?nfields:1); field++) { 770e630c359SToby Isaac PetscInt fbs,j; 771e630c359SToby Isaac PetscFV fv = NULL; 772e630c359SToby Isaac PetscObject f; 773e630c359SToby Isaac PetscClassId fClass; 774e630c359SToby Isaac PetscBool vector; 7756030a318SStefano Zampini if (nfields && cEnd > cStart) { /* We have user-defined fields/components */ 7769566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section,cStart,field,&fbs)); 777e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 7789566063dSJacob Faibussowitsch PetscCall(DMGetField(dmX,field,NULL,&f)); 7799566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(f,&fClass)); 780e630c359SToby Isaac if (fClass == PETSCFV_CLASSID) { 781e630c359SToby Isaac fv = (PetscFV) f; 782e630c359SToby Isaac } 783e630c359SToby Isaac vector = PETSC_FALSE; 784e630c359SToby Isaac if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) { 785e630c359SToby Isaac vector = PETSC_TRUE; 786e630c359SToby Isaac for (j = 0; j < fbs; j++) { 787e630c359SToby Isaac const char *compName = NULL; 788e630c359SToby Isaac if (fv) { 7899566063dSJacob Faibussowitsch PetscCall(PetscFVGetComponentName(fv,j,&compName)); 790e630c359SToby Isaac if (compName) break; 791e630c359SToby Isaac } 792e630c359SToby Isaac } 793e630c359SToby Isaac if (j < fbs) vector = PETSC_FALSE; 794e630c359SToby Isaac } 795e630c359SToby Isaac if (vector) { 796955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 7979566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].ncells*3,MPIU_VTUREAL,tag)); 798955d60d1SToby Isaac } 799e630c359SToby Isaac } else { 800e630c359SToby Isaac for (i=0; i<fbs; i++) { 801955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 8029566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].ncells,MPIU_VTUREAL,tag)); 803955d60d1SToby Isaac } 804552f7358SJed Brown } 805552f7358SJed Brown } 806e630c359SToby Isaac } 807e630c359SToby Isaac } 808552f7358SJed Brown /* all point data */ 809552f7358SJed Brown for (link=vtk->link; link; link=link->next) { 810e630c359SToby Isaac Vec X = (Vec)link->vec; 811e630c359SToby Isaac DM dmX; 8126030a318SStefano Zampini PetscInt bs = 1, nfields, field; 813e630c359SToby Isaac PetscSection section = NULL; 814e630c359SToby Isaac 815552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 8169566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 817e630c359SToby Isaac if (!dmX) dmX = dm; 8189566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject*) §ion)); 8199566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 8209566063dSJacob Faibussowitsch if (vEnd > vStart) PetscCall(PetscSectionGetDof(section,vStart,&bs)); 8219566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section,&nfields)); 822e630c359SToby Isaac field = 0; 823e630c359SToby Isaac if (link->field >= 0) { 824e630c359SToby Isaac field = link->field; 825e630c359SToby Isaac nfields = field + 1; 826e630c359SToby Isaac } 827e630c359SToby Isaac for (i=0; field<(nfields?nfields:1); field++) { 828e630c359SToby Isaac PetscInt fbs; 8296030a318SStefano Zampini if (nfields && vEnd > vStart) { /* We have user-defined fields/components */ 8309566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section,vStart,field,&fbs)); 831e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 832e630c359SToby Isaac if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) { 833955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 8349566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices*3,MPIU_VTUREAL,tag)); 835955d60d1SToby Isaac } 836e630c359SToby Isaac } else { 837e630c359SToby Isaac for (i=0; i<fbs; i++) { 838955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 8399566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm,viewer,fp,r,0,NULL,buffer,gpiece[r].nvertices,MPIU_VTUREAL,tag)); 840955d60d1SToby Isaac } 841552f7358SJed Brown } 842552f7358SJed Brown } 843552f7358SJed Brown } 844552f7358SJed Brown } 845e630c359SToby Isaac } 846e630c359SToby Isaac } 8479566063dSJacob Faibussowitsch PetscCall(PetscFree(gpiece)); 8489566063dSJacob Faibussowitsch PetscCall(PetscFree(buffer)); 8499566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"\n </AppendedData>\n")); 8509566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"</VTKFile>\n")); 8519566063dSJacob Faibussowitsch PetscCall(PetscFClose(comm,fp)); 852552f7358SJed Brown PetscFunctionReturn(0); 853552f7358SJed Brown } 854