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 266497c311SBarry Smith static PetscErrorCode TransferWrite(MPI_Comm comm, PetscViewer viewer, FILE *fp, PetscMPIInt srank, PetscMPIInt root, const void *send, void *recv, PetscCount count, MPI_Datatype mpidatatype, PetscMPIInt tag) 27d71ae5a4SJacob Faibussowitsch { 28552f7358SJed Brown PetscMPIInt rank; 29552f7358SJed Brown 30552f7358SJed Brown PetscFunctionBegin; 319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 32552f7358SJed Brown if (rank == srank && rank != root) { 336497c311SBarry Smith PetscCallMPI(MPIU_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; 416497c311SBarry Smith PetscCallMPI(MPIU_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 } 483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49552f7358SJed Brown } 50552f7358SJed Brown 51d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexGetVTKConnectivity(DM dm, PetscBool localized, PieceInfo *piece, PetscVTKInt **oconn, PetscVTKInt **ooffsets, PetscVTKType **otypes) 52d71ae5a4SJacob Faibussowitsch { 536858538eSMatthew 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)); 616858538eSMatthew 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; 826858538eSMatthew 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)) { 906497c311SBarry Smith if (!localized) PetscCall(PetscVTKIntCast(closure[v] - vStart, &conn[countconn++])); 916497c311SBarry Smith else PetscCall(PetscVTKIntCast(startoffset + nC, &conn[countconn++])); 92724b5320SMatthew G. Knepley ++nC; 93552f7358SJed Brown } 94552f7358SJed Brown } 959566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 962e529ec8SStefano Zampini } else { 976497c311SBarry Smith for (nC = 0; nC < dof / dim; nC++) PetscCall(PetscVTKIntCast(startoffset + nC, &conn[countconn++])); 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)); 1046497c311SBarry Smith for (i = 0; i < n; ++i) PetscCall(PetscVTKIntCast(cone[i], &conn[s + i])); 10596ca5757SLisandro Dalcin } 1066497c311SBarry Smith PetscCall(PetscVTKIntCast(countconn, &offsets[countcell])); 1078865f1eaSKarl Rupp 108552f7358SJed Brown nverts = countconn - startoffset; 1099566063dSJacob Faibussowitsch PetscCall(DMPlexVTKGetCellType_Internal(dm, dim, nverts, &celltype)); 1108865f1eaSKarl Rupp 1116497c311SBarry Smith types[countcell] = (PetscVTKType)celltype; 112552f7358SJed Brown countcell++; 113552f7358SJed Brown } 11408401ef6SPierre Jolivet PetscCheck(countcell == piece->ncells, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent cell count"); 11508401ef6SPierre Jolivet PetscCheck(countconn == piece->nconn, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent connectivity count"); 116552f7358SJed Brown *oconn = conn; 117552f7358SJed Brown *ooffsets = offsets; 118552f7358SJed Brown *otypes = types; 1193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 120552f7358SJed Brown } 121552f7358SJed Brown 122c29ce622SStefano Zampini PETSC_INTERN PetscErrorCode DMPlexGetNonEmptyComm_Private(DM dm, MPI_Comm *comm) 123c29ce622SStefano Zampini { 124c29ce622SStefano Zampini DM_Plex *mesh = (DM_Plex *)dm->data; 125c29ce622SStefano Zampini 126c29ce622SStefano Zampini PetscFunctionBegin; 127c29ce622SStefano Zampini if (mesh->nonempty_comm == MPI_COMM_SELF) { /* Not yet setup */ 128c29ce622SStefano Zampini PetscInt cStart, cEnd, cellHeight; 129c29ce622SStefano Zampini MPI_Comm dmcomm = PetscObjectComm((PetscObject)dm); 130c29ce622SStefano Zampini PetscMPIInt color, rank; 131c29ce622SStefano Zampini 132c29ce622SStefano Zampini PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight)); 133c29ce622SStefano Zampini PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd)); 134c29ce622SStefano Zampini color = (cStart < cEnd) ? 0 : 1; 135c29ce622SStefano Zampini PetscCallMPI(MPI_Comm_rank(dmcomm, &rank)); 136c29ce622SStefano Zampini PetscCallMPI(MPI_Comm_split(dmcomm, color, rank, &mesh->nonempty_comm)); 137c29ce622SStefano Zampini if (color == 1) { 138c29ce622SStefano Zampini PetscCallMPI(MPI_Comm_free(&mesh->nonempty_comm)); 139c29ce622SStefano Zampini mesh->nonempty_comm = MPI_COMM_NULL; 140c29ce622SStefano Zampini } 141c29ce622SStefano Zampini } 142c29ce622SStefano Zampini *comm = mesh->nonempty_comm; 143c29ce622SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 144c29ce622SStefano Zampini } 145c29ce622SStefano Zampini 1466b8451b4SStefano Zampini static PetscErrorCode DMGetFieldIfFV_Private(DM dm, PetscInt field, PetscFV *fv) 1476b8451b4SStefano Zampini { 1486b8451b4SStefano Zampini PetscObject f = NULL; 1496b8451b4SStefano Zampini PetscClassId fClass = PETSC_SMALLEST_CLASSID; 1506b8451b4SStefano Zampini PetscInt nf; 1516b8451b4SStefano Zampini 1526b8451b4SStefano Zampini PetscFunctionBegin; 1536b8451b4SStefano Zampini *fv = NULL; 1546b8451b4SStefano Zampini PetscCall(DMGetNumFields(dm, &nf)); 1556b8451b4SStefano Zampini if (nf > 0) { 1566b8451b4SStefano Zampini PetscCall(DMGetField(dm, field, NULL, &f)); 1576b8451b4SStefano Zampini PetscCall(PetscObjectGetClassId(f, &fClass)); 1586b8451b4SStefano Zampini if (fClass == PETSCFV_CLASSID) *fv = (PetscFV)f; 1596b8451b4SStefano Zampini } 1606b8451b4SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 1616b8451b4SStefano Zampini } 1626b8451b4SStefano Zampini 163552f7358SJed Brown /* 164552f7358SJed Brown Write all fields that have been provided to the viewer 165552f7358SJed Brown Multi-block XML format with binary appended data. 166552f7358SJed Brown */ 167d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm, PetscViewer viewer) 168d71ae5a4SJacob Faibussowitsch { 169ce94432eSBarry Smith MPI_Comm comm; 1706858538eSMatthew G. Knepley PetscSection coordSection, cellCoordSection; 171552f7358SJed Brown PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data; 172552f7358SJed Brown PetscViewerVTKObjectLink link; 173552f7358SJed Brown FILE *fp; 174552f7358SJed Brown PetscMPIInt rank, size, tag; 1756497c311SBarry Smith PetscInt dimEmbed, cellHeight, cStart, cEnd, vStart, vEnd, numLabelCells, hasLabel, c, v, i; 1762e529ec8SStefano Zampini PetscBool localized; 1770298fd71SBarry Smith PieceInfo piece, *gpiece = NULL; 1780298fd71SBarry Smith void *buffer = NULL; 17930815ce0SLisandro Dalcin const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian"; 180*fc2fb351SPierre Jolivet PetscInt loops_per_scalar = PetscDefined(USE_COMPLEX) ? 2 : 1; 181552f7358SJed Brown 182552f7358SJed Brown PetscFunctionBegin; 1839566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dimEmbed)); 1849566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight)); 1859566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd)); 1869566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1879566063dSJacob Faibussowitsch PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells)); 1889566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalized(dm, &localized)); 1899566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 1906858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinateSection(dm, &cellCoordSection)); 191c29ce622SStefano Zampini PetscCall(PetscCommGetNewTag(PetscObjectComm((PetscObject)dm), &tag)); 192c29ce622SStefano Zampini 193c29ce622SStefano Zampini PetscCall(DMPlexGetNonEmptyComm_Private(dm, &comm)); 194c29ce622SStefano Zampini if (comm == MPI_COMM_NULL) goto finalize; 195c29ce622SStefano Zampini PetscCallMPI(MPI_Comm_size(comm, &size)); 196c29ce622SStefano Zampini PetscCallMPI(MPI_Comm_rank(comm, &rank)); 197c29ce622SStefano Zampini 198c29ce622SStefano Zampini PetscCall(PetscFOpen(comm, vtk->filename, "wb", &fp)); 199c29ce622SStefano Zampini PetscCall(PetscFPrintf(comm, fp, "<?xml version=\"1.0\"?>\n")); 200c29ce622SStefano Zampini PetscCall(PetscFPrintf(comm, fp, "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"%s\" header_type=\"UInt64\">\n", byte_order)); 201c29ce622SStefano Zampini PetscCall(PetscFPrintf(comm, fp, " <UnstructuredGrid>\n")); 2028865f1eaSKarl Rupp 203552f7358SJed Brown hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 2042e529ec8SStefano Zampini piece.nvertices = 0; 205552f7358SJed Brown piece.ncells = 0; 206552f7358SJed Brown piece.nconn = 0; 2072e529ec8SStefano Zampini if (!localized) piece.nvertices = vEnd - vStart; 208552f7358SJed Brown for (c = cStart; c < cEnd; ++c) { 20919003fb5SStefano Zampini PetscInt dof = 0; 210552f7358SJed Brown if (hasLabel) { 211552f7358SJed Brown PetscInt value; 212552f7358SJed Brown 2139566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dm, "vtk", c, &value)); 214552f7358SJed Brown if (value != 1) continue; 215552f7358SJed Brown } 2166858538eSMatthew G. Knepley if (localized) PetscCall(PetscSectionGetDof(cellCoordSection, c, &dof)); 21719003fb5SStefano Zampini if (!dof) { 2182e529ec8SStefano Zampini PetscInt *closure = NULL; 2192e529ec8SStefano Zampini PetscInt closureSize; 2202e529ec8SStefano Zampini 2219566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 222552f7358SJed Brown for (v = 0; v < closureSize * 2; v += 2) { 2232e529ec8SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 2242e529ec8SStefano Zampini piece.nconn++; 22519003fb5SStefano Zampini if (localized) piece.nvertices++; 2262e529ec8SStefano Zampini } 227552f7358SJed Brown } 2289566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 2292e529ec8SStefano Zampini } else { 2302e529ec8SStefano Zampini piece.nvertices += dof / dimEmbed; 2312e529ec8SStefano Zampini piece.nconn += dof / dimEmbed; 2322e529ec8SStefano Zampini } 233552f7358SJed Brown piece.ncells++; 234552f7358SJed Brown } 2359566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc1(size, &gpiece)); 2369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gather((PetscInt *)&piece, sizeof(piece) / sizeof(PetscInt), MPIU_INT, (PetscInt *)gpiece, sizeof(piece) / sizeof(PetscInt), MPIU_INT, 0, comm)); 237552f7358SJed Brown 238552f7358SJed Brown /* 239552f7358SJed Brown * Write file header 240552f7358SJed Brown */ 241dd400576SPatrick Sanan if (rank == 0) { 2420f2609c8SJed Brown PetscInt64 boffset = 0; 243552f7358SJed Brown 2446497c311SBarry Smith for (PetscMPIInt r = 0; r < size; r++) { 24563a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <Piece NumberOfPoints=\"%" PetscInt_FMT "\" NumberOfCells=\"%" PetscInt_FMT "\">\n", gpiece[r].nvertices, gpiece[r].ncells)); 246552f7358SJed Brown /* Coordinate positions */ 2479566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <Points>\n")); 2480f2609c8SJed Brown PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, boffset)); 2496f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 2509566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </Points>\n")); 251552f7358SJed Brown /* Cell connectivity */ 2529566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <Cells>\n")); 2530f2609c8SJed Brown PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", boffset)); 2546f98e5d0SMatthew G. Knepley boffset += gpiece[r].nconn * sizeof(PetscVTKInt) + (gpiece[r].nconn ? sizeof(PetscInt64) : 0); 2550f2609c8SJed Brown PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"Int32\" Name=\"offsets\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", boffset)); 2566f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(PetscVTKInt) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 2570f2609c8SJed Brown PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"UInt8\" Name=\"types\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", boffset)); 2586f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(unsigned char) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 2599566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </Cells>\n")); 260552f7358SJed Brown 261552f7358SJed Brown /* 262552f7358SJed Brown * Cell Data headers 263552f7358SJed Brown */ 2649566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <CellData>\n")); 2650f2609c8SJed Brown PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", boffset)); 2666f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(PetscVTKInt) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 267552f7358SJed Brown /* all the vectors */ 268552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 269552f7358SJed Brown Vec X = (Vec)link->vec; 270e630c359SToby Isaac DM dmX = NULL; 2716030a318SStefano Zampini PetscInt bs = 1, nfields, field; 272552f7358SJed Brown const char *vecname = ""; 273e630c359SToby Isaac PetscSection section; 274552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 275552f7358SJed 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. */ 2769566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)X, &vecname)); 277552f7358SJed Brown } 2789566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 279e630c359SToby Isaac if (!dmX) dmX = dm; 2809566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 2819566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 2829566063dSJacob Faibussowitsch if (cEnd > cStart) PetscCall(PetscSectionGetDof(section, cStart, &bs)); 2839566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 284e630c359SToby Isaac field = 0; 285e630c359SToby Isaac if (link->field >= 0) { 286e630c359SToby Isaac field = link->field; 287e630c359SToby Isaac nfields = field + 1; 288e630c359SToby Isaac } 289e630c359SToby Isaac for (i = 0; field < (nfields ? nfields : 1); field++) { 2901cfafdd3SJed Brown PetscInt fbs, j; 291a00cdb45SToby Isaac PetscFV fv = NULL; 2920298fd71SBarry Smith const char *fieldname = NULL; 2931cfafdd3SJed Brown char buf[256]; 294e630c359SToby Isaac PetscBool vector; 2956b8451b4SStefano Zampini 2967ded4cbaSJed Brown if (nfields) { /* We have user-defined fields/components */ 2979566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, cStart, field, &fbs)); 2989566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(section, field, &fieldname)); 2997ded4cbaSJed Brown } else fbs = bs; /* Say we have one field with 'bs' components */ 3006b8451b4SStefano Zampini PetscCall(DMGetFieldIfFV_Private(dmX, field, &fv)); 301e630c359SToby Isaac if (nfields && !fieldname) { 30263a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "CellField%" PetscInt_FMT, field)); 303552f7358SJed Brown fieldname = buf; 304552f7358SJed Brown } 305e630c359SToby Isaac vector = PETSC_FALSE; 306e630c359SToby Isaac if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) { 307e630c359SToby Isaac vector = PETSC_TRUE; 30863a3b9bcSJacob Faibussowitsch PetscCheck(fbs <= 3, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_SIZ, "Cell vector fields can have at most 3 components, %" PetscInt_FMT " given", fbs); 309e630c359SToby Isaac for (j = 0; j < fbs; j++) { 310e630c359SToby Isaac const char *compName = NULL; 311e630c359SToby Isaac if (fv) { 3129566063dSJacob Faibussowitsch PetscCall(PetscFVGetComponentName(fv, j, &compName)); 313e630c359SToby Isaac if (compName) break; 314e630c359SToby Isaac } 315e630c359SToby Isaac } 316e630c359SToby Isaac if (j < fbs) vector = PETSC_FALSE; 317e630c359SToby Isaac } 318e630c359SToby Isaac if (vector) { 319955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 3200f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s.Re\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 3216f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * 3 * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 3220f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s.Im\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 3236f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * 3 * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 324955d60d1SToby Isaac #else 3250f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 3266f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * 3 * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 327955d60d1SToby Isaac #endif 328e630c359SToby Isaac i += fbs; 329e630c359SToby Isaac } else { 3301cfafdd3SJed Brown for (j = 0; j < fbs; j++) { 331a00cdb45SToby Isaac const char *compName = NULL; 332955d60d1SToby Isaac char finalname[256]; 33348a46eb9SPierre Jolivet if (fv) PetscCall(PetscFVGetComponentName(fv, j, &compName)); 334a00cdb45SToby Isaac if (compName) { 3359566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname, 255, "%s%s.%s", vecname, fieldname, compName)); 3369371c9d4SSatish Balay } else if (fbs > 1) { 33763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname, 255, "%s%s.%" PetscInt_FMT, vecname, fieldname, j)); 338e630c359SToby Isaac } else { 3399566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname, 255, "%s%s", vecname, fieldname)); 340a00cdb45SToby Isaac } 341955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 3420f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset)); 3436f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 3440f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset)); 3456f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 346955d60d1SToby Isaac #else 3470f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset)); 3486f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 349955d60d1SToby Isaac #endif 3501cfafdd3SJed Brown i++; 351552f7358SJed Brown } 352552f7358SJed Brown } 353e630c359SToby Isaac } 35463a3b9bcSJacob Faibussowitsch //PetscCheck(i == bs,comm,PETSC_ERR_PLIB,"Total number of field components %" PetscInt_FMT " != block size %" PetscInt_FMT,i,bs); 3551cfafdd3SJed Brown } 3569566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </CellData>\n")); 357552f7358SJed Brown 358552f7358SJed Brown /* 359552f7358SJed Brown * Point Data headers 360552f7358SJed Brown */ 3619566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <PointData>\n")); 362552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 363552f7358SJed Brown Vec X = (Vec)link->vec; 364e630c359SToby Isaac DM dmX; 3656030a318SStefano Zampini PetscInt bs = 1, nfields, field; 366552f7358SJed Brown const char *vecname = ""; 367e630c359SToby Isaac PetscSection section; 368552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 369552f7358SJed 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. */ 3709566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)X, &vecname)); 371552f7358SJed Brown } 3729566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 373e630c359SToby Isaac if (!dmX) dmX = dm; 3749566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 3759566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 3769566063dSJacob Faibussowitsch if (vEnd > vStart) PetscCall(PetscSectionGetDof(section, vStart, &bs)); 3779566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 378e630c359SToby Isaac field = 0; 379e630c359SToby Isaac if (link->field >= 0) { 380e630c359SToby Isaac field = link->field; 381e630c359SToby Isaac nfields = field + 1; 382e630c359SToby Isaac } 38363bfac88SBarry Smith for (; field < (nfields ? nfields : 1); field++) { 3841cfafdd3SJed Brown PetscInt fbs, j; 3851cfafdd3SJed Brown const char *fieldname = NULL; 3861cfafdd3SJed Brown char buf[256]; 3877ded4cbaSJed Brown if (nfields) { /* We have user-defined fields/components */ 3889566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, vStart, field, &fbs)); 3899566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(section, field, &fieldname)); 3907ded4cbaSJed Brown } else fbs = bs; /* Say we have one field with 'bs' components */ 391e630c359SToby Isaac if (nfields && !fieldname) { 39263a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "PointField%" PetscInt_FMT, field)); 3931cfafdd3SJed Brown fieldname = buf; 3941cfafdd3SJed Brown } 395e630c359SToby Isaac if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) { 39663a3b9bcSJacob Faibussowitsch PetscCheck(fbs <= 3, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_SIZ, "Point vector fields can have at most 3 components, %" PetscInt_FMT " given", fbs); 397955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 3980f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s.Re\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 3996f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 4000f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s.Im\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 4016f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 402955d60d1SToby Isaac #else 4030f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 4046f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 405955d60d1SToby Isaac #endif 406e630c359SToby Isaac } else { 4071cfafdd3SJed Brown for (j = 0; j < fbs; j++) { 408b778fa18SValeria Barra const char *compName = NULL; 409955d60d1SToby Isaac char finalname[256]; 4109566063dSJacob Faibussowitsch PetscCall(PetscSectionGetComponentName(section, field, j, &compName)); 4119566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname, 255, "%s%s.%s", vecname, fieldname, compName)); 412955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 4130f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset)); 4146f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 4150f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset)); 4166f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 417955d60d1SToby Isaac #else 4180f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset)); 4196f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 420955d60d1SToby Isaac #endif 421552f7358SJed Brown } 422552f7358SJed Brown } 4231cfafdd3SJed Brown } 424e630c359SToby Isaac } 4259566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </PointData>\n")); 4269566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </Piece>\n")); 427552f7358SJed Brown } 428552f7358SJed Brown } 429552f7358SJed Brown 4309566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </UnstructuredGrid>\n")); 4319566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <AppendedData encoding=\"raw\">\n")); 4329566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "_")); 433552f7358SJed Brown 434dd400576SPatrick Sanan if (rank == 0) { 435552f7358SJed Brown PetscInt maxsize = 0; 4366497c311SBarry Smith for (PetscMPIInt r = 0; r < size; r++) { 437955d60d1SToby Isaac maxsize = PetscMax(maxsize, (PetscInt)(gpiece[r].nvertices * 3 * sizeof(PetscVTUReal))); 438955d60d1SToby Isaac maxsize = PetscMax(maxsize, (PetscInt)(gpiece[r].ncells * 3 * sizeof(PetscVTUReal))); 439552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt)(gpiece[r].nconn * sizeof(PetscVTKInt))); 440552f7358SJed Brown } 4419566063dSJacob Faibussowitsch PetscCall(PetscMalloc(maxsize, &buffer)); 442552f7358SJed Brown } 4436497c311SBarry Smith for (PetscMPIInt r = 0; r < size; r++) { 444552f7358SJed Brown if (r == rank) { 445552f7358SJed Brown PetscInt nsend; 446552f7358SJed Brown { /* Position */ 4476858538eSMatthew G. Knepley const PetscScalar *x, *cx = NULL; 448955d60d1SToby Isaac PetscVTUReal *y = NULL; 4496858538eSMatthew G. Knepley Vec coords, cellCoords; 450*fc2fb351SPierre Jolivet const PetscBool copy = PetscDefined(USE_COMPLEX) ? PETSC_TRUE : (PetscBool)(dimEmbed != 3 || localized || (sizeof(PetscReal) != sizeof(PetscVTUReal))); 4512e529ec8SStefano Zampini 4529566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coords)); 4539566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coords, &x)); 4546858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinatesLocal(dm, &cellCoords)); 4556858538eSMatthew G. Knepley if (cellCoords) PetscCall(VecGetArrayRead(cellCoords, &cx)); 456955d60d1SToby Isaac if (copy) { 4579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.nvertices * 3, &y)); 45819003fb5SStefano Zampini if (localized) { 45919003fb5SStefano Zampini PetscInt cnt; 46019003fb5SStefano Zampini for (c = cStart, cnt = 0; c < cEnd; c++) { 46119003fb5SStefano Zampini PetscInt off, dof; 46219003fb5SStefano Zampini 4636858538eSMatthew G. Knepley PetscCall(PetscSectionGetDof(cellCoordSection, c, &dof)); 46419003fb5SStefano Zampini if (!dof) { 46519003fb5SStefano Zampini PetscInt *closure = NULL; 46619003fb5SStefano Zampini PetscInt closureSize; 46719003fb5SStefano Zampini 4689566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 46919003fb5SStefano Zampini for (v = 0; v < closureSize * 2; v += 2) { 47019003fb5SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 4719566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, closure[v], &off)); 47219003fb5SStefano Zampini if (dimEmbed != 3) { 473955d60d1SToby Isaac y[cnt * 3 + 0] = (PetscVTUReal)PetscRealPart(x[off + 0]); 474955d60d1SToby Isaac y[cnt * 3 + 1] = (PetscVTUReal)((dimEmbed > 1) ? PetscRealPart(x[off + 1]) : 0.0); 475955d60d1SToby Isaac y[cnt * 3 + 2] = (PetscVTUReal)0.0; 47619003fb5SStefano Zampini } else { 477955d60d1SToby Isaac y[cnt * 3 + 0] = (PetscVTUReal)PetscRealPart(x[off + 0]); 478955d60d1SToby Isaac y[cnt * 3 + 1] = (PetscVTUReal)PetscRealPart(x[off + 1]); 479955d60d1SToby Isaac y[cnt * 3 + 2] = (PetscVTUReal)PetscRealPart(x[off + 2]); 48019003fb5SStefano Zampini } 48119003fb5SStefano Zampini cnt++; 48219003fb5SStefano Zampini } 48319003fb5SStefano Zampini } 4849566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 48519003fb5SStefano Zampini } else { 4866858538eSMatthew G. Knepley PetscCall(PetscSectionGetOffset(cellCoordSection, c, &off)); 48719003fb5SStefano Zampini if (dimEmbed != 3) { 48819003fb5SStefano Zampini for (i = 0; i < dof / dimEmbed; i++) { 4896858538eSMatthew G. Knepley y[cnt * 3 + 0] = (PetscVTUReal)PetscRealPart(cx[off + i * dimEmbed + 0]); 4906858538eSMatthew G. Knepley y[cnt * 3 + 1] = (PetscVTUReal)((dimEmbed > 1) ? PetscRealPart(cx[off + i * dimEmbed + 1]) : 0.0); 491955d60d1SToby Isaac y[cnt * 3 + 2] = (PetscVTUReal)0.0; 49219003fb5SStefano Zampini cnt++; 49319003fb5SStefano Zampini } 49419003fb5SStefano Zampini } else { 495ad540459SPierre Jolivet for (i = 0; i < dof; i++) y[cnt * 3 + i] = (PetscVTUReal)PetscRealPart(cx[off + i]); 49619003fb5SStefano Zampini cnt += dof / dimEmbed; 49719003fb5SStefano Zampini } 49819003fb5SStefano Zampini } 49919003fb5SStefano Zampini } 50008401ef6SPierre Jolivet PetscCheck(cnt == piece.nvertices, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 50119003fb5SStefano Zampini } else { 502552f7358SJed Brown for (i = 0; i < piece.nvertices; i++) { 503955d60d1SToby Isaac y[i * 3 + 0] = (PetscVTUReal)PetscRealPart(x[i * dimEmbed + 0]); 504955d60d1SToby Isaac y[i * 3 + 1] = (PetscVTUReal)((dimEmbed > 1) ? PetscRealPart(x[i * dimEmbed + 1]) : 0.); 505955d60d1SToby Isaac y[i * 3 + 2] = (PetscVTUReal)((dimEmbed > 2) ? PetscRealPart(x[i * dimEmbed + 2]) : 0.); 50619003fb5SStefano Zampini } 507552f7358SJed Brown } 508552f7358SJed Brown } 5092e529ec8SStefano Zampini nsend = piece.nvertices * 3; 5109566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, copy ? (const void *)y : (const void *)x, buffer, nsend, MPIU_VTUREAL, tag)); 5119566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 5129566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coords, &x)); 5136858538eSMatthew G. Knepley if (cellCoords) PetscCall(VecRestoreArrayRead(cellCoords, &cx)); 514552f7358SJed Brown } 515552f7358SJed Brown { /* Connectivity, offsets, types */ 5163e869605SMatthew G. Knepley PetscVTKInt *connectivity = NULL, *offsets = NULL; 5173e869605SMatthew G. Knepley PetscVTKType *types = NULL; 5189566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKConnectivity(dm, localized, &piece, &connectivity, &offsets, &types)); 5199566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, connectivity, buffer, piece.nconn, MPI_INT, tag)); 5209566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, offsets, buffer, piece.ncells, MPI_INT, tag)); 5219566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, types, buffer, piece.ncells, MPI_CHAR, tag)); 5229566063dSJacob Faibussowitsch PetscCall(PetscFree3(connectivity, offsets, types)); 523552f7358SJed Brown } 524552f7358SJed Brown { /* Owners (cell data) */ 525552f7358SJed Brown PetscVTKInt *owners; 526c29ce622SStefano Zampini PetscMPIInt orank; 527c29ce622SStefano Zampini 528c29ce622SStefano Zampini PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &orank)); 5299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.ncells, &owners)); 530c29ce622SStefano Zampini for (i = 0; i < piece.ncells; i++) owners[i] = orank; 5319566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, owners, buffer, piece.ncells, MPI_INT, tag)); 5329566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 533552f7358SJed Brown } 534552f7358SJed Brown /* Cell data */ 535552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 536552f7358SJed Brown Vec X = (Vec)link->vec; 537e630c359SToby Isaac DM dmX; 538552f7358SJed Brown const PetscScalar *x; 539955d60d1SToby Isaac PetscVTUReal *y; 5406030a318SStefano Zampini PetscInt bs = 1, nfields, field; 541e630c359SToby Isaac PetscSection section = NULL; 542e630c359SToby Isaac 543552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 5449566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 545e630c359SToby Isaac if (!dmX) dmX = dm; 5469566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 5479566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 5489566063dSJacob Faibussowitsch if (cEnd > cStart) PetscCall(PetscSectionGetDof(section, cStart, &bs)); 5499566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 550e630c359SToby Isaac field = 0; 551e630c359SToby Isaac if (link->field >= 0) { 552e630c359SToby Isaac field = link->field; 553e630c359SToby Isaac nfields = field + 1; 554e630c359SToby Isaac } 5559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 5569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.ncells * 3, &y)); 55763bfac88SBarry Smith for (; field < (nfields ? nfields : 1); field++) { 558e630c359SToby Isaac PetscInt fbs, j; 559e630c359SToby Isaac PetscFV fv = NULL; 560e630c359SToby Isaac PetscBool vector; 5616b8451b4SStefano Zampini 5626030a318SStefano Zampini if (nfields && cEnd > cStart) { /* We have user-defined fields/components */ 5639566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, cStart, field, &fbs)); 564e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 5656b8451b4SStefano Zampini PetscCall(DMGetFieldIfFV_Private(dmX, field, &fv)); 566e630c359SToby Isaac vector = PETSC_FALSE; 567e630c359SToby Isaac if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) { 568e630c359SToby Isaac vector = PETSC_TRUE; 569e630c359SToby Isaac for (j = 0; j < fbs; j++) { 570e630c359SToby Isaac const char *compName = NULL; 571e630c359SToby Isaac if (fv) { 5729566063dSJacob Faibussowitsch PetscCall(PetscFVGetComponentName(fv, j, &compName)); 573e630c359SToby Isaac if (compName) break; 574e630c359SToby Isaac } 575e630c359SToby Isaac } 576e630c359SToby Isaac if (j < fbs) vector = PETSC_FALSE; 577e630c359SToby Isaac } 578e630c359SToby Isaac if (vector) { 579955d60d1SToby Isaac PetscInt cnt, l; 580955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 581552f7358SJed Brown for (c = cStart, cnt = 0; c < cEnd; c++) { 582e630c359SToby Isaac const PetscScalar *xpoint; 583e630c359SToby Isaac PetscInt off, j; 584e630c359SToby Isaac 585552f7358SJed Brown if (hasLabel) { /* Ignore some cells */ 586552f7358SJed Brown PetscInt value; 5879566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dmX, "vtk", c, &value)); 588552f7358SJed Brown if (value != 1) continue; 589552f7358SJed Brown } 590e630c359SToby Isaac if (nfields) { 5919566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, c, field, &off)); 592e630c359SToby Isaac } else { 5939566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, c, &off)); 594e630c359SToby Isaac } 595e630c359SToby Isaac xpoint = &x[off]; 596ad540459SPierre Jolivet for (j = 0; j < fbs; j++) y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j])); 597e630c359SToby Isaac for (; j < 3; j++) y[cnt++] = 0.; 598e630c359SToby Isaac } 59908401ef6SPierre Jolivet PetscCheck(cnt == piece.ncells * 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 6009566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.ncells * 3, MPIU_VTUREAL, tag)); 601955d60d1SToby Isaac } 602e630c359SToby Isaac } else { 603e630c359SToby Isaac for (i = 0; i < fbs; i++) { 604955d60d1SToby Isaac PetscInt cnt, l; 605955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 606e630c359SToby Isaac for (c = cStart, cnt = 0; c < cEnd; c++) { 607e630c359SToby Isaac const PetscScalar *xpoint; 608e630c359SToby Isaac PetscInt off; 609e630c359SToby Isaac 610e630c359SToby Isaac if (hasLabel) { /* Ignore some cells */ 611e630c359SToby Isaac PetscInt value; 6129566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dmX, "vtk", c, &value)); 613e630c359SToby Isaac if (value != 1) continue; 614e630c359SToby Isaac } 615e630c359SToby Isaac if (nfields) { 6169566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, c, field, &off)); 617e630c359SToby Isaac } else { 6189566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, c, &off)); 619e630c359SToby Isaac } 620e630c359SToby Isaac xpoint = &x[off]; 621955d60d1SToby Isaac y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i])); 622552f7358SJed Brown } 62308401ef6SPierre Jolivet PetscCheck(cnt == piece.ncells, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 6249566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.ncells, MPIU_VTUREAL, tag)); 625955d60d1SToby Isaac } 626552f7358SJed Brown } 627e630c359SToby Isaac } 628e630c359SToby Isaac } 6299566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 6309566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 631552f7358SJed Brown } 632e630c359SToby Isaac /* point data */ 633552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 634552f7358SJed Brown Vec X = (Vec)link->vec; 635e630c359SToby Isaac DM dmX; 636552f7358SJed Brown const PetscScalar *x; 637955d60d1SToby Isaac PetscVTUReal *y; 6386030a318SStefano Zampini PetscInt bs = 1, nfields, field; 639e630c359SToby Isaac PetscSection section = NULL; 640e630c359SToby Isaac 641552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 6429566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 643e630c359SToby Isaac if (!dmX) dmX = dm; 6449566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 6459566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 6469566063dSJacob Faibussowitsch if (vEnd > vStart) PetscCall(PetscSectionGetDof(section, vStart, &bs)); 6479566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 648e630c359SToby Isaac field = 0; 649e630c359SToby Isaac if (link->field >= 0) { 650e630c359SToby Isaac field = link->field; 651e630c359SToby Isaac nfields = field + 1; 652e630c359SToby Isaac } 6539566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 6549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.nvertices * 3, &y)); 65563bfac88SBarry Smith for (; field < (nfields ? nfields : 1); field++) { 656e630c359SToby Isaac PetscInt fbs, j; 6576030a318SStefano Zampini if (nfields && vEnd > vStart) { /* We have user-defined fields/components */ 6589566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, vStart, field, &fbs)); 659e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 660e630c359SToby Isaac if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) { 661955d60d1SToby Isaac PetscInt cnt, l; 662955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 6632e529ec8SStefano Zampini if (!localized) { 664552f7358SJed Brown for (v = vStart, cnt = 0; v < vEnd; v++) { 665e630c359SToby Isaac PetscInt off; 666e630c359SToby Isaac const PetscScalar *xpoint; 667e630c359SToby Isaac 668e630c359SToby Isaac if (nfields) { 6699566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, v, field, &off)); 670e630c359SToby Isaac } else { 6719566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, v, &off)); 672e630c359SToby Isaac } 673e630c359SToby Isaac xpoint = &x[off]; 674ad540459SPierre Jolivet for (j = 0; j < fbs; j++) y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j])); 675e630c359SToby Isaac for (; j < 3; j++) y[cnt++] = 0.; 676e630c359SToby Isaac } 677e630c359SToby Isaac } else { 678e630c359SToby Isaac for (c = cStart, cnt = 0; c < cEnd; c++) { 679e630c359SToby Isaac PetscInt *closure = NULL; 680e630c359SToby Isaac PetscInt closureSize, off; 681e630c359SToby Isaac 6829566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 683e630c359SToby Isaac for (v = 0, off = 0; v < closureSize * 2; v += 2) { 684e630c359SToby Isaac if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 685e630c359SToby Isaac PetscInt voff; 686e630c359SToby Isaac const PetscScalar *xpoint; 687e630c359SToby Isaac 688e630c359SToby Isaac if (nfields) { 6899566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, closure[v], field, &voff)); 690e630c359SToby Isaac } else { 6919566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, closure[v], &voff)); 692e630c359SToby Isaac } 693e630c359SToby Isaac xpoint = &x[voff]; 694ad540459SPierre Jolivet for (j = 0; j < fbs; j++) y[cnt + off++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j])); 695e630c359SToby Isaac for (; j < 3; j++) y[cnt + off++] = 0.; 696e630c359SToby Isaac } 697e630c359SToby Isaac } 698e630c359SToby Isaac cnt += off; 6999566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 700e630c359SToby Isaac } 701e630c359SToby Isaac } 70208401ef6SPierre Jolivet PetscCheck(cnt == piece.nvertices * 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 7039566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.nvertices * 3, MPIU_VTUREAL, tag)); 704955d60d1SToby Isaac } 705e630c359SToby Isaac } else { 706e630c359SToby Isaac for (i = 0; i < fbs; i++) { 707955d60d1SToby Isaac PetscInt cnt, l; 708955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 709e630c359SToby Isaac if (!localized) { 710e630c359SToby Isaac for (v = vStart, cnt = 0; v < vEnd; v++) { 711e630c359SToby Isaac PetscInt off; 712e630c359SToby Isaac const PetscScalar *xpoint; 713e630c359SToby Isaac 714e630c359SToby Isaac if (nfields) { 7159566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, v, field, &off)); 716e630c359SToby Isaac } else { 7179566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, v, &off)); 718e630c359SToby Isaac } 719e630c359SToby Isaac xpoint = &x[off]; 720955d60d1SToby Isaac y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i])); 721552f7358SJed Brown } 7222e529ec8SStefano Zampini } else { 7232e529ec8SStefano Zampini for (c = cStart, cnt = 0; c < cEnd; c++) { 7242e529ec8SStefano Zampini PetscInt *closure = NULL; 7252e529ec8SStefano Zampini PetscInt closureSize, off; 7262e529ec8SStefano Zampini 7279566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 7282e529ec8SStefano Zampini for (v = 0, off = 0; v < closureSize * 2; v += 2) { 7292e529ec8SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 730e630c359SToby Isaac PetscInt voff; 731e630c359SToby Isaac const PetscScalar *xpoint; 7322e529ec8SStefano Zampini 733e630c359SToby Isaac if (nfields) { 7349566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, closure[v], field, &voff)); 735e630c359SToby Isaac } else { 7369566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, closure[v], &voff)); 737e630c359SToby Isaac } 738e630c359SToby Isaac xpoint = &x[voff]; 739955d60d1SToby Isaac y[cnt + off++] = (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i])); 7402e529ec8SStefano Zampini } 7412e529ec8SStefano Zampini } 7429bda96f6SStefano Zampini cnt += off; 7439566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 7442e529ec8SStefano Zampini } 7452e529ec8SStefano Zampini } 74608401ef6SPierre Jolivet PetscCheck(cnt == piece.nvertices, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 7479566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.nvertices, MPIU_VTUREAL, tag)); 748955d60d1SToby Isaac } 749552f7358SJed Brown } 750e630c359SToby Isaac } 751e630c359SToby Isaac } 7529566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 7539566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 754552f7358SJed Brown } 755dd400576SPatrick Sanan } else if (rank == 0) { 756955d60d1SToby Isaac PetscInt l; 757955d60d1SToby Isaac 7589566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nvertices * 3, MPIU_VTUREAL, tag)); /* positions */ 7599566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nconn, MPI_INT, tag)); /* connectivity */ 7609566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPI_INT, tag)); /* offsets */ 7619566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPI_CHAR, tag)); /* types */ 7629566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPI_INT, tag)); /* owner rank (cells) */ 763552f7358SJed Brown /* all cell data */ 764552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 765e630c359SToby Isaac Vec X = (Vec)link->vec; 7666030a318SStefano Zampini PetscInt bs = 1, nfields, field; 767e630c359SToby Isaac DM dmX; 768e630c359SToby Isaac PetscSection section = NULL; 769e630c359SToby Isaac 770552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 7719566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 772e630c359SToby Isaac if (!dmX) dmX = dm; 7739566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 7749566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 7759566063dSJacob Faibussowitsch if (cEnd > cStart) PetscCall(PetscSectionGetDof(section, cStart, &bs)); 7769566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 777e630c359SToby Isaac field = 0; 778e630c359SToby Isaac if (link->field >= 0) { 779e630c359SToby Isaac field = link->field; 780e630c359SToby Isaac nfields = field + 1; 781e630c359SToby Isaac } 78263bfac88SBarry Smith for (; field < (nfields ? nfields : 1); field++) { 783e630c359SToby Isaac PetscInt fbs, j; 784e630c359SToby Isaac PetscFV fv = NULL; 785e630c359SToby Isaac PetscBool vector; 7866b8451b4SStefano Zampini 7876030a318SStefano Zampini if (nfields && cEnd > cStart) { /* We have user-defined fields/components */ 7889566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, cStart, field, &fbs)); 789e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 7906b8451b4SStefano Zampini PetscCall(DMGetFieldIfFV_Private(dmX, field, &fv)); 791e630c359SToby Isaac vector = PETSC_FALSE; 792e630c359SToby Isaac if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) { 793e630c359SToby Isaac vector = PETSC_TRUE; 794e630c359SToby Isaac for (j = 0; j < fbs; j++) { 795e630c359SToby Isaac const char *compName = NULL; 796e630c359SToby Isaac if (fv) { 7979566063dSJacob Faibussowitsch PetscCall(PetscFVGetComponentName(fv, j, &compName)); 798e630c359SToby Isaac if (compName) break; 799e630c359SToby Isaac } 800e630c359SToby Isaac } 801e630c359SToby Isaac if (j < fbs) vector = PETSC_FALSE; 802e630c359SToby Isaac } 803e630c359SToby Isaac if (vector) { 80448a46eb9SPierre Jolivet for (l = 0; l < loops_per_scalar; l++) PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells * 3, MPIU_VTUREAL, tag)); 805e630c359SToby Isaac } else { 806e630c359SToby Isaac for (i = 0; i < fbs; i++) { 80748a46eb9SPierre Jolivet for (l = 0; l < loops_per_scalar; l++) PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPIU_VTUREAL, tag)); 808552f7358SJed Brown } 809552f7358SJed Brown } 810e630c359SToby Isaac } 811e630c359SToby Isaac } 812552f7358SJed Brown /* all point data */ 813552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 814e630c359SToby Isaac Vec X = (Vec)link->vec; 815e630c359SToby Isaac DM dmX; 8166030a318SStefano Zampini PetscInt bs = 1, nfields, field; 817e630c359SToby Isaac PetscSection section = NULL; 818e630c359SToby Isaac 819552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 8209566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 821e630c359SToby Isaac if (!dmX) dmX = dm; 8229566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 8239566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 8249566063dSJacob Faibussowitsch if (vEnd > vStart) PetscCall(PetscSectionGetDof(section, vStart, &bs)); 8259566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 826e630c359SToby Isaac field = 0; 827e630c359SToby Isaac if (link->field >= 0) { 828e630c359SToby Isaac field = link->field; 829e630c359SToby Isaac nfields = field + 1; 830e630c359SToby Isaac } 83163bfac88SBarry Smith for (; field < (nfields ? nfields : 1); field++) { 832e630c359SToby Isaac PetscInt fbs; 8336030a318SStefano Zampini if (nfields && vEnd > vStart) { /* We have user-defined fields/components */ 8349566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, vStart, field, &fbs)); 835e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 836e630c359SToby Isaac if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) { 83748a46eb9SPierre Jolivet for (l = 0; l < loops_per_scalar; l++) PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nvertices * 3, MPIU_VTUREAL, tag)); 838e630c359SToby Isaac } else { 839e630c359SToby Isaac for (i = 0; i < fbs; i++) { 84048a46eb9SPierre Jolivet for (l = 0; l < loops_per_scalar; l++) PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nvertices, MPIU_VTUREAL, tag)); 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)); 852c29ce622SStefano Zampini finalize: 85317bfc685SStefano Zampini /* this code sends to rank 0 that writes. 85417bfc685SStefano Zampini It may lead to very unbalanced log_view timings 85517bfc685SStefano Zampini of the next PETSc function logged. 85617bfc685SStefano Zampini Since this call is not performance critical, we 85717bfc685SStefano Zampini issue a barrier here to synchronize the processes */ 85817bfc685SStefano Zampini PetscCall(PetscBarrier((PetscObject)viewer)); 8593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 860552f7358SJed Brown } 861