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 146552f7358SJed Brown /* 147552f7358SJed Brown Write all fields that have been provided to the viewer 148552f7358SJed Brown Multi-block XML format with binary appended data. 149552f7358SJed Brown */ 150d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm, PetscViewer viewer) 151d71ae5a4SJacob Faibussowitsch { 152ce94432eSBarry Smith MPI_Comm comm; 1536858538eSMatthew G. Knepley PetscSection coordSection, cellCoordSection; 154552f7358SJed Brown PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data; 155552f7358SJed Brown PetscViewerVTKObjectLink link; 156552f7358SJed Brown FILE *fp; 157552f7358SJed Brown PetscMPIInt rank, size, tag; 1586497c311SBarry Smith PetscInt dimEmbed, cellHeight, cStart, cEnd, vStart, vEnd, numLabelCells, hasLabel, c, v, i; 1592e529ec8SStefano Zampini PetscBool localized; 1600298fd71SBarry Smith PieceInfo piece, *gpiece = NULL; 1610298fd71SBarry Smith void *buffer = NULL; 16230815ce0SLisandro Dalcin const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian"; 163955d60d1SToby Isaac PetscInt loops_per_scalar; 164552f7358SJed Brown 165552f7358SJed Brown PetscFunctionBegin; 1669566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dimEmbed)); 1679566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight)); 1689566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd)); 1699566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1709566063dSJacob Faibussowitsch PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells)); 1719566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalized(dm, &localized)); 1729566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 1736858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinateSection(dm, &cellCoordSection)); 174c29ce622SStefano Zampini PetscCall(PetscCommGetNewTag(PetscObjectComm((PetscObject)dm), &tag)); 175c29ce622SStefano Zampini 176c29ce622SStefano Zampini PetscCall(DMPlexGetNonEmptyComm_Private(dm, &comm)); 177c29ce622SStefano Zampini #if defined(PETSC_USE_COMPLEX) 178c29ce622SStefano Zampini loops_per_scalar = 2; 179c29ce622SStefano Zampini #else 180c29ce622SStefano Zampini loops_per_scalar = 1; 181c29ce622SStefano Zampini #endif 182c29ce622SStefano Zampini if (comm == MPI_COMM_NULL) goto finalize; 183c29ce622SStefano Zampini PetscCallMPI(MPI_Comm_size(comm, &size)); 184c29ce622SStefano Zampini PetscCallMPI(MPI_Comm_rank(comm, &rank)); 185c29ce622SStefano Zampini 186c29ce622SStefano Zampini PetscCall(PetscFOpen(comm, vtk->filename, "wb", &fp)); 187c29ce622SStefano Zampini PetscCall(PetscFPrintf(comm, fp, "<?xml version=\"1.0\"?>\n")); 188c29ce622SStefano Zampini PetscCall(PetscFPrintf(comm, fp, "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"%s\" header_type=\"UInt64\">\n", byte_order)); 189c29ce622SStefano Zampini PetscCall(PetscFPrintf(comm, fp, " <UnstructuredGrid>\n")); 1908865f1eaSKarl Rupp 191552f7358SJed Brown hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 1922e529ec8SStefano Zampini piece.nvertices = 0; 193552f7358SJed Brown piece.ncells = 0; 194552f7358SJed Brown piece.nconn = 0; 1952e529ec8SStefano Zampini if (!localized) piece.nvertices = vEnd - vStart; 196552f7358SJed Brown for (c = cStart; c < cEnd; ++c) { 19719003fb5SStefano Zampini PetscInt dof = 0; 198552f7358SJed Brown if (hasLabel) { 199552f7358SJed Brown PetscInt value; 200552f7358SJed Brown 2019566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dm, "vtk", c, &value)); 202552f7358SJed Brown if (value != 1) continue; 203552f7358SJed Brown } 2046858538eSMatthew G. Knepley if (localized) PetscCall(PetscSectionGetDof(cellCoordSection, c, &dof)); 20519003fb5SStefano Zampini if (!dof) { 2062e529ec8SStefano Zampini PetscInt *closure = NULL; 2072e529ec8SStefano Zampini PetscInt closureSize; 2082e529ec8SStefano Zampini 2099566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 210552f7358SJed Brown for (v = 0; v < closureSize * 2; v += 2) { 2112e529ec8SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 2122e529ec8SStefano Zampini piece.nconn++; 21319003fb5SStefano Zampini if (localized) piece.nvertices++; 2142e529ec8SStefano Zampini } 215552f7358SJed Brown } 2169566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 2172e529ec8SStefano Zampini } else { 2182e529ec8SStefano Zampini piece.nvertices += dof / dimEmbed; 2192e529ec8SStefano Zampini piece.nconn += dof / dimEmbed; 2202e529ec8SStefano Zampini } 221552f7358SJed Brown piece.ncells++; 222552f7358SJed Brown } 2239566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc1(size, &gpiece)); 2249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gather((PetscInt *)&piece, sizeof(piece) / sizeof(PetscInt), MPIU_INT, (PetscInt *)gpiece, sizeof(piece) / sizeof(PetscInt), MPIU_INT, 0, comm)); 225552f7358SJed Brown 226552f7358SJed Brown /* 227552f7358SJed Brown * Write file header 228552f7358SJed Brown */ 229dd400576SPatrick Sanan if (rank == 0) { 2300f2609c8SJed Brown PetscInt64 boffset = 0; 231552f7358SJed Brown 2326497c311SBarry Smith for (PetscMPIInt r = 0; r < size; r++) { 23363a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <Piece NumberOfPoints=\"%" PetscInt_FMT "\" NumberOfCells=\"%" PetscInt_FMT "\">\n", gpiece[r].nvertices, gpiece[r].ncells)); 234552f7358SJed Brown /* Coordinate positions */ 2359566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <Points>\n")); 2360f2609c8SJed Brown PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, boffset)); 2376f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 2389566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </Points>\n")); 239552f7358SJed Brown /* Cell connectivity */ 2409566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <Cells>\n")); 2410f2609c8SJed Brown PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", boffset)); 2426f98e5d0SMatthew G. Knepley boffset += gpiece[r].nconn * sizeof(PetscVTKInt) + (gpiece[r].nconn ? sizeof(PetscInt64) : 0); 2430f2609c8SJed Brown PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"Int32\" Name=\"offsets\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", boffset)); 2446f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(PetscVTKInt) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 2450f2609c8SJed Brown PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"UInt8\" Name=\"types\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", boffset)); 2466f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(unsigned char) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 2479566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </Cells>\n")); 248552f7358SJed Brown 249552f7358SJed Brown /* 250552f7358SJed Brown * Cell Data headers 251552f7358SJed Brown */ 2529566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <CellData>\n")); 2530f2609c8SJed Brown PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", boffset)); 2546f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(PetscVTKInt) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 255552f7358SJed Brown /* all the vectors */ 256552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 257552f7358SJed Brown Vec X = (Vec)link->vec; 258e630c359SToby Isaac DM dmX = NULL; 2596030a318SStefano Zampini PetscInt bs = 1, nfields, field; 260552f7358SJed Brown const char *vecname = ""; 261e630c359SToby Isaac PetscSection section; 262552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 263552f7358SJed 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. */ 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)X, &vecname)); 265552f7358SJed Brown } 2669566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 267e630c359SToby Isaac if (!dmX) dmX = dm; 2689566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 2699566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 2709566063dSJacob Faibussowitsch if (cEnd > cStart) PetscCall(PetscSectionGetDof(section, cStart, &bs)); 2719566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 272e630c359SToby Isaac field = 0; 273e630c359SToby Isaac if (link->field >= 0) { 274e630c359SToby Isaac field = link->field; 275e630c359SToby Isaac nfields = field + 1; 276e630c359SToby Isaac } 277e630c359SToby Isaac for (i = 0; field < (nfields ? nfields : 1); field++) { 2781cfafdd3SJed Brown PetscInt fbs, j; 279a00cdb45SToby Isaac PetscFV fv = NULL; 280a00cdb45SToby Isaac PetscObject f; 281a00cdb45SToby Isaac PetscClassId fClass; 2820298fd71SBarry Smith const char *fieldname = NULL; 2831cfafdd3SJed Brown char buf[256]; 284e630c359SToby Isaac PetscBool vector; 2857ded4cbaSJed Brown if (nfields) { /* We have user-defined fields/components */ 2869566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, cStart, field, &fbs)); 2879566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(section, field, &fieldname)); 2887ded4cbaSJed Brown } else fbs = bs; /* Say we have one field with 'bs' components */ 2899566063dSJacob Faibussowitsch PetscCall(DMGetField(dmX, field, NULL, &f)); 2909566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(f, &fClass)); 291ad540459SPierre Jolivet if (fClass == PETSCFV_CLASSID) fv = (PetscFV)f; 292e630c359SToby Isaac if (nfields && !fieldname) { 29363a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "CellField%" PetscInt_FMT, field)); 294552f7358SJed Brown fieldname = buf; 295552f7358SJed Brown } 296e630c359SToby Isaac vector = PETSC_FALSE; 297e630c359SToby Isaac if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) { 298e630c359SToby Isaac vector = PETSC_TRUE; 29963a3b9bcSJacob Faibussowitsch PetscCheck(fbs <= 3, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_SIZ, "Cell vector fields can have at most 3 components, %" PetscInt_FMT " given", fbs); 300e630c359SToby Isaac for (j = 0; j < fbs; j++) { 301e630c359SToby Isaac const char *compName = NULL; 302e630c359SToby Isaac if (fv) { 3039566063dSJacob Faibussowitsch PetscCall(PetscFVGetComponentName(fv, j, &compName)); 304e630c359SToby Isaac if (compName) break; 305e630c359SToby Isaac } 306e630c359SToby Isaac } 307e630c359SToby Isaac if (j < fbs) vector = PETSC_FALSE; 308e630c359SToby Isaac } 309e630c359SToby Isaac if (vector) { 310955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 3110f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s.Re\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 3126f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * 3 * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 3130f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s.Im\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 3146f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * 3 * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 315955d60d1SToby Isaac #else 3160f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 3176f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * 3 * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 318955d60d1SToby Isaac #endif 319e630c359SToby Isaac i += fbs; 320e630c359SToby Isaac } else { 3211cfafdd3SJed Brown for (j = 0; j < fbs; j++) { 322a00cdb45SToby Isaac const char *compName = NULL; 323955d60d1SToby Isaac char finalname[256]; 32448a46eb9SPierre Jolivet if (fv) PetscCall(PetscFVGetComponentName(fv, j, &compName)); 325a00cdb45SToby Isaac if (compName) { 3269566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname, 255, "%s%s.%s", vecname, fieldname, compName)); 3279371c9d4SSatish Balay } else if (fbs > 1) { 32863a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname, 255, "%s%s.%" PetscInt_FMT, vecname, fieldname, j)); 329e630c359SToby Isaac } else { 3309566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname, 255, "%s%s", vecname, fieldname)); 331a00cdb45SToby Isaac } 332955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 3330f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset)); 3346f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 3350f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset)); 3366f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 337955d60d1SToby Isaac #else 3380f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset)); 3396f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 340955d60d1SToby Isaac #endif 3411cfafdd3SJed Brown i++; 342552f7358SJed Brown } 343552f7358SJed Brown } 344e630c359SToby Isaac } 34563a3b9bcSJacob Faibussowitsch //PetscCheck(i == bs,comm,PETSC_ERR_PLIB,"Total number of field components %" PetscInt_FMT " != block size %" PetscInt_FMT,i,bs); 3461cfafdd3SJed Brown } 3479566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </CellData>\n")); 348552f7358SJed Brown 349552f7358SJed Brown /* 350552f7358SJed Brown * Point Data headers 351552f7358SJed Brown */ 3529566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <PointData>\n")); 353552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 354552f7358SJed Brown Vec X = (Vec)link->vec; 355e630c359SToby Isaac DM dmX; 3566030a318SStefano Zampini PetscInt bs = 1, nfields, field; 357552f7358SJed Brown const char *vecname = ""; 358e630c359SToby Isaac PetscSection section; 359552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 360552f7358SJed 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. */ 3619566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)X, &vecname)); 362552f7358SJed Brown } 3639566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 364e630c359SToby Isaac if (!dmX) dmX = dm; 3659566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 3669566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 3679566063dSJacob Faibussowitsch if (vEnd > vStart) PetscCall(PetscSectionGetDof(section, vStart, &bs)); 3689566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 369e630c359SToby Isaac field = 0; 370e630c359SToby Isaac if (link->field >= 0) { 371e630c359SToby Isaac field = link->field; 372e630c359SToby Isaac nfields = field + 1; 373e630c359SToby Isaac } 374*63bfac88SBarry Smith for (; field < (nfields ? nfields : 1); field++) { 3751cfafdd3SJed Brown PetscInt fbs, j; 3761cfafdd3SJed Brown const char *fieldname = NULL; 3771cfafdd3SJed Brown char buf[256]; 3787ded4cbaSJed Brown if (nfields) { /* We have user-defined fields/components */ 3799566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, vStart, field, &fbs)); 3809566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(section, field, &fieldname)); 3817ded4cbaSJed Brown } else fbs = bs; /* Say we have one field with 'bs' components */ 382e630c359SToby Isaac if (nfields && !fieldname) { 38363a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "PointField%" PetscInt_FMT, field)); 3841cfafdd3SJed Brown fieldname = buf; 3851cfafdd3SJed Brown } 386e630c359SToby Isaac if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) { 38763a3b9bcSJacob Faibussowitsch PetscCheck(fbs <= 3, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_SIZ, "Point vector fields can have at most 3 components, %" PetscInt_FMT " given", fbs); 388955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 3890f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s.Re\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 3906f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 3910f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s.Im\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 3926f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 393955d60d1SToby Isaac #else 3940f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 3956f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 396955d60d1SToby Isaac #endif 397e630c359SToby Isaac } else { 3981cfafdd3SJed Brown for (j = 0; j < fbs; j++) { 399b778fa18SValeria Barra const char *compName = NULL; 400955d60d1SToby Isaac char finalname[256]; 4019566063dSJacob Faibussowitsch PetscCall(PetscSectionGetComponentName(section, field, j, &compName)); 4029566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname, 255, "%s%s.%s", vecname, fieldname, compName)); 403955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 4040f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset)); 4056f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 4060f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset)); 4076f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 408955d60d1SToby Isaac #else 4090f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset)); 4106f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 411955d60d1SToby Isaac #endif 412552f7358SJed Brown } 413552f7358SJed Brown } 4141cfafdd3SJed Brown } 415e630c359SToby Isaac } 4169566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </PointData>\n")); 4179566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </Piece>\n")); 418552f7358SJed Brown } 419552f7358SJed Brown } 420552f7358SJed Brown 4219566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </UnstructuredGrid>\n")); 4229566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <AppendedData encoding=\"raw\">\n")); 4239566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "_")); 424552f7358SJed Brown 425dd400576SPatrick Sanan if (rank == 0) { 426552f7358SJed Brown PetscInt maxsize = 0; 4276497c311SBarry Smith for (PetscMPIInt r = 0; r < size; r++) { 428955d60d1SToby Isaac maxsize = PetscMax(maxsize, (PetscInt)(gpiece[r].nvertices * 3 * sizeof(PetscVTUReal))); 429955d60d1SToby Isaac maxsize = PetscMax(maxsize, (PetscInt)(gpiece[r].ncells * 3 * sizeof(PetscVTUReal))); 430552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt)(gpiece[r].nconn * sizeof(PetscVTKInt))); 431552f7358SJed Brown } 4329566063dSJacob Faibussowitsch PetscCall(PetscMalloc(maxsize, &buffer)); 433552f7358SJed Brown } 4346497c311SBarry Smith for (PetscMPIInt r = 0; r < size; r++) { 435552f7358SJed Brown if (r == rank) { 436552f7358SJed Brown PetscInt nsend; 437552f7358SJed Brown { /* Position */ 4386858538eSMatthew G. Knepley const PetscScalar *x, *cx = NULL; 439955d60d1SToby Isaac PetscVTUReal *y = NULL; 4406858538eSMatthew G. Knepley Vec coords, cellCoords; 441955d60d1SToby Isaac PetscBool copy; 4422e529ec8SStefano Zampini 4439566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coords)); 4449566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coords, &x)); 4456858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinatesLocal(dm, &cellCoords)); 4466858538eSMatthew G. Knepley if (cellCoords) PetscCall(VecGetArrayRead(cellCoords, &cx)); 447955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 448955d60d1SToby Isaac copy = PETSC_TRUE; 449955d60d1SToby Isaac #else 450955d60d1SToby Isaac copy = (PetscBool)(dimEmbed != 3 || localized || (sizeof(PetscReal) != sizeof(PetscVTUReal))); 451955d60d1SToby Isaac #endif 452955d60d1SToby Isaac if (copy) { 4539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.nvertices * 3, &y)); 45419003fb5SStefano Zampini if (localized) { 45519003fb5SStefano Zampini PetscInt cnt; 45619003fb5SStefano Zampini for (c = cStart, cnt = 0; c < cEnd; c++) { 45719003fb5SStefano Zampini PetscInt off, dof; 45819003fb5SStefano Zampini 4596858538eSMatthew G. Knepley PetscCall(PetscSectionGetDof(cellCoordSection, c, &dof)); 46019003fb5SStefano Zampini if (!dof) { 46119003fb5SStefano Zampini PetscInt *closure = NULL; 46219003fb5SStefano Zampini PetscInt closureSize; 46319003fb5SStefano Zampini 4649566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 46519003fb5SStefano Zampini for (v = 0; v < closureSize * 2; v += 2) { 46619003fb5SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 4679566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, closure[v], &off)); 46819003fb5SStefano Zampini if (dimEmbed != 3) { 469955d60d1SToby Isaac y[cnt * 3 + 0] = (PetscVTUReal)PetscRealPart(x[off + 0]); 470955d60d1SToby Isaac y[cnt * 3 + 1] = (PetscVTUReal)((dimEmbed > 1) ? PetscRealPart(x[off + 1]) : 0.0); 471955d60d1SToby Isaac y[cnt * 3 + 2] = (PetscVTUReal)0.0; 47219003fb5SStefano Zampini } else { 473955d60d1SToby Isaac y[cnt * 3 + 0] = (PetscVTUReal)PetscRealPart(x[off + 0]); 474955d60d1SToby Isaac y[cnt * 3 + 1] = (PetscVTUReal)PetscRealPart(x[off + 1]); 475955d60d1SToby Isaac y[cnt * 3 + 2] = (PetscVTUReal)PetscRealPart(x[off + 2]); 47619003fb5SStefano Zampini } 47719003fb5SStefano Zampini cnt++; 47819003fb5SStefano Zampini } 47919003fb5SStefano Zampini } 4809566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 48119003fb5SStefano Zampini } else { 4826858538eSMatthew G. Knepley PetscCall(PetscSectionGetOffset(cellCoordSection, c, &off)); 48319003fb5SStefano Zampini if (dimEmbed != 3) { 48419003fb5SStefano Zampini for (i = 0; i < dof / dimEmbed; i++) { 4856858538eSMatthew G. Knepley y[cnt * 3 + 0] = (PetscVTUReal)PetscRealPart(cx[off + i * dimEmbed + 0]); 4866858538eSMatthew G. Knepley y[cnt * 3 + 1] = (PetscVTUReal)((dimEmbed > 1) ? PetscRealPart(cx[off + i * dimEmbed + 1]) : 0.0); 487955d60d1SToby Isaac y[cnt * 3 + 2] = (PetscVTUReal)0.0; 48819003fb5SStefano Zampini cnt++; 48919003fb5SStefano Zampini } 49019003fb5SStefano Zampini } else { 491ad540459SPierre Jolivet for (i = 0; i < dof; i++) y[cnt * 3 + i] = (PetscVTUReal)PetscRealPart(cx[off + i]); 49219003fb5SStefano Zampini cnt += dof / dimEmbed; 49319003fb5SStefano Zampini } 49419003fb5SStefano Zampini } 49519003fb5SStefano Zampini } 49608401ef6SPierre Jolivet PetscCheck(cnt == piece.nvertices, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 49719003fb5SStefano Zampini } else { 498552f7358SJed Brown for (i = 0; i < piece.nvertices; i++) { 499955d60d1SToby Isaac y[i * 3 + 0] = (PetscVTUReal)PetscRealPart(x[i * dimEmbed + 0]); 500955d60d1SToby Isaac y[i * 3 + 1] = (PetscVTUReal)((dimEmbed > 1) ? PetscRealPart(x[i * dimEmbed + 1]) : 0.); 501955d60d1SToby Isaac y[i * 3 + 2] = (PetscVTUReal)((dimEmbed > 2) ? PetscRealPart(x[i * dimEmbed + 2]) : 0.); 50219003fb5SStefano Zampini } 503552f7358SJed Brown } 504552f7358SJed Brown } 5052e529ec8SStefano Zampini nsend = piece.nvertices * 3; 5069566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, copy ? (const void *)y : (const void *)x, buffer, nsend, MPIU_VTUREAL, tag)); 5079566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 5089566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coords, &x)); 5096858538eSMatthew G. Knepley if (cellCoords) PetscCall(VecRestoreArrayRead(cellCoords, &cx)); 510552f7358SJed Brown } 511552f7358SJed Brown { /* Connectivity, offsets, types */ 5123e869605SMatthew G. Knepley PetscVTKInt *connectivity = NULL, *offsets = NULL; 5133e869605SMatthew G. Knepley PetscVTKType *types = NULL; 5149566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKConnectivity(dm, localized, &piece, &connectivity, &offsets, &types)); 5159566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, connectivity, buffer, piece.nconn, MPI_INT, tag)); 5169566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, offsets, buffer, piece.ncells, MPI_INT, tag)); 5179566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, types, buffer, piece.ncells, MPI_CHAR, tag)); 5189566063dSJacob Faibussowitsch PetscCall(PetscFree3(connectivity, offsets, types)); 519552f7358SJed Brown } 520552f7358SJed Brown { /* Owners (cell data) */ 521552f7358SJed Brown PetscVTKInt *owners; 522c29ce622SStefano Zampini PetscMPIInt orank; 523c29ce622SStefano Zampini 524c29ce622SStefano Zampini PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &orank)); 5259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.ncells, &owners)); 526c29ce622SStefano Zampini for (i = 0; i < piece.ncells; i++) owners[i] = orank; 5279566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, owners, buffer, piece.ncells, MPI_INT, tag)); 5289566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 529552f7358SJed Brown } 530552f7358SJed Brown /* Cell data */ 531552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 532552f7358SJed Brown Vec X = (Vec)link->vec; 533e630c359SToby Isaac DM dmX; 534552f7358SJed Brown const PetscScalar *x; 535955d60d1SToby Isaac PetscVTUReal *y; 5366030a318SStefano Zampini PetscInt bs = 1, nfields, field; 537e630c359SToby Isaac PetscSection section = NULL; 538e630c359SToby Isaac 539552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 5409566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 541e630c359SToby Isaac if (!dmX) dmX = dm; 5429566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 5439566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 5449566063dSJacob Faibussowitsch if (cEnd > cStart) PetscCall(PetscSectionGetDof(section, cStart, &bs)); 5459566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 546e630c359SToby Isaac field = 0; 547e630c359SToby Isaac if (link->field >= 0) { 548e630c359SToby Isaac field = link->field; 549e630c359SToby Isaac nfields = field + 1; 550e630c359SToby Isaac } 5519566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 5529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.ncells * 3, &y)); 553*63bfac88SBarry Smith for (; field < (nfields ? nfields : 1); field++) { 554e630c359SToby Isaac PetscInt fbs, j; 555e630c359SToby Isaac PetscFV fv = NULL; 556e630c359SToby Isaac PetscObject f; 557e630c359SToby Isaac PetscClassId fClass; 558e630c359SToby Isaac PetscBool vector; 5596030a318SStefano Zampini if (nfields && cEnd > cStart) { /* We have user-defined fields/components */ 5609566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, cStart, field, &fbs)); 561e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 5629566063dSJacob Faibussowitsch PetscCall(DMGetField(dmX, field, NULL, &f)); 5639566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(f, &fClass)); 564ad540459SPierre Jolivet if (fClass == PETSCFV_CLASSID) fv = (PetscFV)f; 565e630c359SToby Isaac vector = PETSC_FALSE; 566e630c359SToby Isaac if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) { 567e630c359SToby Isaac vector = PETSC_TRUE; 568e630c359SToby Isaac for (j = 0; j < fbs; j++) { 569e630c359SToby Isaac const char *compName = NULL; 570e630c359SToby Isaac if (fv) { 5719566063dSJacob Faibussowitsch PetscCall(PetscFVGetComponentName(fv, j, &compName)); 572e630c359SToby Isaac if (compName) break; 573e630c359SToby Isaac } 574e630c359SToby Isaac } 575e630c359SToby Isaac if (j < fbs) vector = PETSC_FALSE; 576e630c359SToby Isaac } 577e630c359SToby Isaac if (vector) { 578955d60d1SToby Isaac PetscInt cnt, l; 579955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 580552f7358SJed Brown for (c = cStart, cnt = 0; c < cEnd; c++) { 581e630c359SToby Isaac const PetscScalar *xpoint; 582e630c359SToby Isaac PetscInt off, j; 583e630c359SToby Isaac 584552f7358SJed Brown if (hasLabel) { /* Ignore some cells */ 585552f7358SJed Brown PetscInt value; 5869566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dmX, "vtk", c, &value)); 587552f7358SJed Brown if (value != 1) continue; 588552f7358SJed Brown } 589e630c359SToby Isaac if (nfields) { 5909566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, c, field, &off)); 591e630c359SToby Isaac } else { 5929566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, c, &off)); 593e630c359SToby Isaac } 594e630c359SToby Isaac xpoint = &x[off]; 595ad540459SPierre Jolivet for (j = 0; j < fbs; j++) y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j])); 596e630c359SToby Isaac for (; j < 3; j++) y[cnt++] = 0.; 597e630c359SToby Isaac } 59808401ef6SPierre Jolivet PetscCheck(cnt == piece.ncells * 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 5999566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.ncells * 3, MPIU_VTUREAL, tag)); 600955d60d1SToby Isaac } 601e630c359SToby Isaac } else { 602e630c359SToby Isaac for (i = 0; i < fbs; i++) { 603955d60d1SToby Isaac PetscInt cnt, l; 604955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 605e630c359SToby Isaac for (c = cStart, cnt = 0; c < cEnd; c++) { 606e630c359SToby Isaac const PetscScalar *xpoint; 607e630c359SToby Isaac PetscInt off; 608e630c359SToby Isaac 609e630c359SToby Isaac if (hasLabel) { /* Ignore some cells */ 610e630c359SToby Isaac PetscInt value; 6119566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dmX, "vtk", c, &value)); 612e630c359SToby Isaac if (value != 1) continue; 613e630c359SToby Isaac } 614e630c359SToby Isaac if (nfields) { 6159566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, c, field, &off)); 616e630c359SToby Isaac } else { 6179566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, c, &off)); 618e630c359SToby Isaac } 619e630c359SToby Isaac xpoint = &x[off]; 620955d60d1SToby Isaac y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i])); 621552f7358SJed Brown } 62208401ef6SPierre Jolivet PetscCheck(cnt == piece.ncells, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 6239566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.ncells, MPIU_VTUREAL, tag)); 624955d60d1SToby Isaac } 625552f7358SJed Brown } 626e630c359SToby Isaac } 627e630c359SToby Isaac } 6289566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 6299566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 630552f7358SJed Brown } 631e630c359SToby Isaac /* point data */ 632552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 633552f7358SJed Brown Vec X = (Vec)link->vec; 634e630c359SToby Isaac DM dmX; 635552f7358SJed Brown const PetscScalar *x; 636955d60d1SToby Isaac PetscVTUReal *y; 6376030a318SStefano Zampini PetscInt bs = 1, nfields, field; 638e630c359SToby Isaac PetscSection section = NULL; 639e630c359SToby Isaac 640552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 6419566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 642e630c359SToby Isaac if (!dmX) dmX = dm; 6439566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 6449566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 6459566063dSJacob Faibussowitsch if (vEnd > vStart) PetscCall(PetscSectionGetDof(section, vStart, &bs)); 6469566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 647e630c359SToby Isaac field = 0; 648e630c359SToby Isaac if (link->field >= 0) { 649e630c359SToby Isaac field = link->field; 650e630c359SToby Isaac nfields = field + 1; 651e630c359SToby Isaac } 6529566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 6539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.nvertices * 3, &y)); 654*63bfac88SBarry Smith for (; field < (nfields ? nfields : 1); field++) { 655e630c359SToby Isaac PetscInt fbs, j; 6566030a318SStefano Zampini if (nfields && vEnd > vStart) { /* We have user-defined fields/components */ 6579566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, vStart, field, &fbs)); 658e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 659e630c359SToby Isaac if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) { 660955d60d1SToby Isaac PetscInt cnt, l; 661955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 6622e529ec8SStefano Zampini if (!localized) { 663552f7358SJed Brown for (v = vStart, cnt = 0; v < vEnd; v++) { 664e630c359SToby Isaac PetscInt off; 665e630c359SToby Isaac const PetscScalar *xpoint; 666e630c359SToby Isaac 667e630c359SToby Isaac if (nfields) { 6689566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, v, field, &off)); 669e630c359SToby Isaac } else { 6709566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, v, &off)); 671e630c359SToby Isaac } 672e630c359SToby Isaac xpoint = &x[off]; 673ad540459SPierre Jolivet for (j = 0; j < fbs; j++) y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j])); 674e630c359SToby Isaac for (; j < 3; j++) y[cnt++] = 0.; 675e630c359SToby Isaac } 676e630c359SToby Isaac } else { 677e630c359SToby Isaac for (c = cStart, cnt = 0; c < cEnd; c++) { 678e630c359SToby Isaac PetscInt *closure = NULL; 679e630c359SToby Isaac PetscInt closureSize, off; 680e630c359SToby Isaac 6819566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 682e630c359SToby Isaac for (v = 0, off = 0; v < closureSize * 2; v += 2) { 683e630c359SToby Isaac if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 684e630c359SToby Isaac PetscInt voff; 685e630c359SToby Isaac const PetscScalar *xpoint; 686e630c359SToby Isaac 687e630c359SToby Isaac if (nfields) { 6889566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, closure[v], field, &voff)); 689e630c359SToby Isaac } else { 6909566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, closure[v], &voff)); 691e630c359SToby Isaac } 692e630c359SToby Isaac xpoint = &x[voff]; 693ad540459SPierre Jolivet for (j = 0; j < fbs; j++) y[cnt + off++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j])); 694e630c359SToby Isaac for (; j < 3; j++) y[cnt + off++] = 0.; 695e630c359SToby Isaac } 696e630c359SToby Isaac } 697e630c359SToby Isaac cnt += off; 6989566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 699e630c359SToby Isaac } 700e630c359SToby Isaac } 70108401ef6SPierre Jolivet PetscCheck(cnt == piece.nvertices * 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 7029566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.nvertices * 3, MPIU_VTUREAL, tag)); 703955d60d1SToby Isaac } 704e630c359SToby Isaac } else { 705e630c359SToby Isaac for (i = 0; i < fbs; i++) { 706955d60d1SToby Isaac PetscInt cnt, l; 707955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 708e630c359SToby Isaac if (!localized) { 709e630c359SToby Isaac for (v = vStart, cnt = 0; v < vEnd; v++) { 710e630c359SToby Isaac PetscInt off; 711e630c359SToby Isaac const PetscScalar *xpoint; 712e630c359SToby Isaac 713e630c359SToby Isaac if (nfields) { 7149566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, v, field, &off)); 715e630c359SToby Isaac } else { 7169566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, v, &off)); 717e630c359SToby Isaac } 718e630c359SToby Isaac xpoint = &x[off]; 719955d60d1SToby Isaac y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i])); 720552f7358SJed Brown } 7212e529ec8SStefano Zampini } else { 7222e529ec8SStefano Zampini for (c = cStart, cnt = 0; c < cEnd; c++) { 7232e529ec8SStefano Zampini PetscInt *closure = NULL; 7242e529ec8SStefano Zampini PetscInt closureSize, off; 7252e529ec8SStefano Zampini 7269566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 7272e529ec8SStefano Zampini for (v = 0, off = 0; v < closureSize * 2; v += 2) { 7282e529ec8SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 729e630c359SToby Isaac PetscInt voff; 730e630c359SToby Isaac const PetscScalar *xpoint; 7312e529ec8SStefano Zampini 732e630c359SToby Isaac if (nfields) { 7339566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, closure[v], field, &voff)); 734e630c359SToby Isaac } else { 7359566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, closure[v], &voff)); 736e630c359SToby Isaac } 737e630c359SToby Isaac xpoint = &x[voff]; 738955d60d1SToby Isaac y[cnt + off++] = (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i])); 7392e529ec8SStefano Zampini } 7402e529ec8SStefano Zampini } 7419bda96f6SStefano Zampini cnt += off; 7429566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 7432e529ec8SStefano Zampini } 7442e529ec8SStefano Zampini } 74508401ef6SPierre Jolivet PetscCheck(cnt == piece.nvertices, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 7469566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.nvertices, MPIU_VTUREAL, tag)); 747955d60d1SToby Isaac } 748552f7358SJed Brown } 749e630c359SToby Isaac } 750e630c359SToby Isaac } 7519566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 7529566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 753552f7358SJed Brown } 754dd400576SPatrick Sanan } else if (rank == 0) { 755955d60d1SToby Isaac PetscInt l; 756955d60d1SToby Isaac 7579566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nvertices * 3, MPIU_VTUREAL, tag)); /* positions */ 7589566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nconn, MPI_INT, tag)); /* connectivity */ 7599566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPI_INT, tag)); /* offsets */ 7609566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPI_CHAR, tag)); /* types */ 7619566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPI_INT, tag)); /* owner rank (cells) */ 762552f7358SJed Brown /* all cell data */ 763552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 764e630c359SToby Isaac Vec X = (Vec)link->vec; 7656030a318SStefano Zampini PetscInt bs = 1, nfields, field; 766e630c359SToby Isaac DM dmX; 767e630c359SToby Isaac PetscSection section = NULL; 768e630c359SToby Isaac 769552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 7709566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 771e630c359SToby Isaac if (!dmX) dmX = dm; 7729566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 7739566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 7749566063dSJacob Faibussowitsch if (cEnd > cStart) PetscCall(PetscSectionGetDof(section, cStart, &bs)); 7759566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 776e630c359SToby Isaac field = 0; 777e630c359SToby Isaac if (link->field >= 0) { 778e630c359SToby Isaac field = link->field; 779e630c359SToby Isaac nfields = field + 1; 780e630c359SToby Isaac } 781*63bfac88SBarry Smith for (; field < (nfields ? nfields : 1); field++) { 782e630c359SToby Isaac PetscInt fbs, j; 783e630c359SToby Isaac PetscFV fv = NULL; 784e630c359SToby Isaac PetscObject f; 785e630c359SToby Isaac PetscClassId fClass; 786e630c359SToby Isaac PetscBool vector; 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 */ 7909566063dSJacob Faibussowitsch PetscCall(DMGetField(dmX, field, NULL, &f)); 7919566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(f, &fClass)); 792ad540459SPierre Jolivet if (fClass == PETSCFV_CLASSID) fv = (PetscFV)f; 793e630c359SToby Isaac vector = PETSC_FALSE; 794e630c359SToby Isaac if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) { 795e630c359SToby Isaac vector = PETSC_TRUE; 796e630c359SToby Isaac for (j = 0; j < fbs; j++) { 797e630c359SToby Isaac const char *compName = NULL; 798e630c359SToby Isaac if (fv) { 7999566063dSJacob Faibussowitsch PetscCall(PetscFVGetComponentName(fv, j, &compName)); 800e630c359SToby Isaac if (compName) break; 801e630c359SToby Isaac } 802e630c359SToby Isaac } 803e630c359SToby Isaac if (j < fbs) vector = PETSC_FALSE; 804e630c359SToby Isaac } 805e630c359SToby Isaac if (vector) { 80648a46eb9SPierre 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)); 807e630c359SToby Isaac } else { 808e630c359SToby Isaac for (i = 0; i < fbs; i++) { 80948a46eb9SPierre Jolivet for (l = 0; l < loops_per_scalar; l++) PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPIU_VTUREAL, tag)); 810552f7358SJed Brown } 811552f7358SJed Brown } 812e630c359SToby Isaac } 813e630c359SToby Isaac } 814552f7358SJed Brown /* all point data */ 815552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 816e630c359SToby Isaac Vec X = (Vec)link->vec; 817e630c359SToby Isaac DM dmX; 8186030a318SStefano Zampini PetscInt bs = 1, nfields, field; 819e630c359SToby Isaac PetscSection section = NULL; 820e630c359SToby Isaac 821552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 8229566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 823e630c359SToby Isaac if (!dmX) dmX = dm; 8249566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 8259566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 8269566063dSJacob Faibussowitsch if (vEnd > vStart) PetscCall(PetscSectionGetDof(section, vStart, &bs)); 8279566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 828e630c359SToby Isaac field = 0; 829e630c359SToby Isaac if (link->field >= 0) { 830e630c359SToby Isaac field = link->field; 831e630c359SToby Isaac nfields = field + 1; 832e630c359SToby Isaac } 833*63bfac88SBarry Smith for (; field < (nfields ? nfields : 1); field++) { 834e630c359SToby Isaac PetscInt fbs; 8356030a318SStefano Zampini if (nfields && vEnd > vStart) { /* We have user-defined fields/components */ 8369566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, vStart, field, &fbs)); 837e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 838e630c359SToby Isaac if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) { 83948a46eb9SPierre 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)); 840e630c359SToby Isaac } else { 841e630c359SToby Isaac for (i = 0; i < fbs; i++) { 84248a46eb9SPierre Jolivet for (l = 0; l < loops_per_scalar; l++) PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nvertices, MPIU_VTUREAL, tag)); 843552f7358SJed Brown } 844552f7358SJed Brown } 845552f7358SJed Brown } 846552f7358SJed Brown } 847e630c359SToby Isaac } 848e630c359SToby Isaac } 8499566063dSJacob Faibussowitsch PetscCall(PetscFree(gpiece)); 8509566063dSJacob Faibussowitsch PetscCall(PetscFree(buffer)); 8519566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "\n </AppendedData>\n")); 8529566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "</VTKFile>\n")); 8539566063dSJacob Faibussowitsch PetscCall(PetscFClose(comm, fp)); 854c29ce622SStefano Zampini finalize: 85517bfc685SStefano Zampini /* this code sends to rank 0 that writes. 85617bfc685SStefano Zampini It may lead to very unbalanced log_view timings 85717bfc685SStefano Zampini of the next PETSc function logged. 85817bfc685SStefano Zampini Since this call is not performance critical, we 85917bfc685SStefano Zampini issue a barrier here to synchronize the processes */ 86017bfc685SStefano Zampini PetscCall(PetscBarrier((PetscObject)viewer)); 8613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 862552f7358SJed Brown } 863