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 146*6b8451b4SStefano Zampini static PetscErrorCode DMGetFieldIfFV_Private(DM dm, PetscInt field, PetscFV *fv) 147*6b8451b4SStefano Zampini { 148*6b8451b4SStefano Zampini PetscObject f = NULL; 149*6b8451b4SStefano Zampini PetscClassId fClass = PETSC_SMALLEST_CLASSID; 150*6b8451b4SStefano Zampini PetscInt nf; 151*6b8451b4SStefano Zampini 152*6b8451b4SStefano Zampini PetscFunctionBegin; 153*6b8451b4SStefano Zampini *fv = NULL; 154*6b8451b4SStefano Zampini PetscCall(DMGetNumFields(dm, &nf)); 155*6b8451b4SStefano Zampini if (nf > 0) { 156*6b8451b4SStefano Zampini PetscCall(DMGetField(dm, field, NULL, &f)); 157*6b8451b4SStefano Zampini PetscCall(PetscObjectGetClassId(f, &fClass)); 158*6b8451b4SStefano Zampini if (fClass == PETSCFV_CLASSID) *fv = (PetscFV)f; 159*6b8451b4SStefano Zampini } 160*6b8451b4SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 161*6b8451b4SStefano Zampini } 162*6b8451b4SStefano 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"; 180955d60d1SToby Isaac PetscInt loops_per_scalar; 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 defined(PETSC_USE_COMPLEX) 195c29ce622SStefano Zampini loops_per_scalar = 2; 196c29ce622SStefano Zampini #else 197c29ce622SStefano Zampini loops_per_scalar = 1; 198c29ce622SStefano Zampini #endif 199c29ce622SStefano Zampini if (comm == MPI_COMM_NULL) goto finalize; 200c29ce622SStefano Zampini PetscCallMPI(MPI_Comm_size(comm, &size)); 201c29ce622SStefano Zampini PetscCallMPI(MPI_Comm_rank(comm, &rank)); 202c29ce622SStefano Zampini 203c29ce622SStefano Zampini PetscCall(PetscFOpen(comm, vtk->filename, "wb", &fp)); 204c29ce622SStefano Zampini PetscCall(PetscFPrintf(comm, fp, "<?xml version=\"1.0\"?>\n")); 205c29ce622SStefano Zampini PetscCall(PetscFPrintf(comm, fp, "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"%s\" header_type=\"UInt64\">\n", byte_order)); 206c29ce622SStefano Zampini PetscCall(PetscFPrintf(comm, fp, " <UnstructuredGrid>\n")); 2078865f1eaSKarl Rupp 208552f7358SJed Brown hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 2092e529ec8SStefano Zampini piece.nvertices = 0; 210552f7358SJed Brown piece.ncells = 0; 211552f7358SJed Brown piece.nconn = 0; 2122e529ec8SStefano Zampini if (!localized) piece.nvertices = vEnd - vStart; 213552f7358SJed Brown for (c = cStart; c < cEnd; ++c) { 21419003fb5SStefano Zampini PetscInt dof = 0; 215552f7358SJed Brown if (hasLabel) { 216552f7358SJed Brown PetscInt value; 217552f7358SJed Brown 2189566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dm, "vtk", c, &value)); 219552f7358SJed Brown if (value != 1) continue; 220552f7358SJed Brown } 2216858538eSMatthew G. Knepley if (localized) PetscCall(PetscSectionGetDof(cellCoordSection, c, &dof)); 22219003fb5SStefano Zampini if (!dof) { 2232e529ec8SStefano Zampini PetscInt *closure = NULL; 2242e529ec8SStefano Zampini PetscInt closureSize; 2252e529ec8SStefano Zampini 2269566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 227552f7358SJed Brown for (v = 0; v < closureSize * 2; v += 2) { 2282e529ec8SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 2292e529ec8SStefano Zampini piece.nconn++; 23019003fb5SStefano Zampini if (localized) piece.nvertices++; 2312e529ec8SStefano Zampini } 232552f7358SJed Brown } 2339566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 2342e529ec8SStefano Zampini } else { 2352e529ec8SStefano Zampini piece.nvertices += dof / dimEmbed; 2362e529ec8SStefano Zampini piece.nconn += dof / dimEmbed; 2372e529ec8SStefano Zampini } 238552f7358SJed Brown piece.ncells++; 239552f7358SJed Brown } 2409566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc1(size, &gpiece)); 2419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gather((PetscInt *)&piece, sizeof(piece) / sizeof(PetscInt), MPIU_INT, (PetscInt *)gpiece, sizeof(piece) / sizeof(PetscInt), MPIU_INT, 0, comm)); 242552f7358SJed Brown 243552f7358SJed Brown /* 244552f7358SJed Brown * Write file header 245552f7358SJed Brown */ 246dd400576SPatrick Sanan if (rank == 0) { 2470f2609c8SJed Brown PetscInt64 boffset = 0; 248552f7358SJed Brown 2496497c311SBarry Smith for (PetscMPIInt r = 0; r < size; r++) { 25063a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <Piece NumberOfPoints=\"%" PetscInt_FMT "\" NumberOfCells=\"%" PetscInt_FMT "\">\n", gpiece[r].nvertices, gpiece[r].ncells)); 251552f7358SJed Brown /* Coordinate positions */ 2529566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <Points>\n")); 2530f2609c8SJed Brown PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, boffset)); 2546f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 2559566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </Points>\n")); 256552f7358SJed Brown /* Cell connectivity */ 2579566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <Cells>\n")); 2580f2609c8SJed Brown PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", boffset)); 2596f98e5d0SMatthew G. Knepley boffset += gpiece[r].nconn * sizeof(PetscVTKInt) + (gpiece[r].nconn ? sizeof(PetscInt64) : 0); 2600f2609c8SJed Brown PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"Int32\" Name=\"offsets\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", boffset)); 2616f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(PetscVTKInt) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 2620f2609c8SJed Brown PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"UInt8\" Name=\"types\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", boffset)); 2636f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(unsigned char) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 2649566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </Cells>\n")); 265552f7358SJed Brown 266552f7358SJed Brown /* 267552f7358SJed Brown * Cell Data headers 268552f7358SJed Brown */ 2699566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <CellData>\n")); 2700f2609c8SJed Brown PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", boffset)); 2716f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(PetscVTKInt) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 272552f7358SJed Brown /* all the vectors */ 273552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 274552f7358SJed Brown Vec X = (Vec)link->vec; 275e630c359SToby Isaac DM dmX = NULL; 2766030a318SStefano Zampini PetscInt bs = 1, nfields, field; 277552f7358SJed Brown const char *vecname = ""; 278e630c359SToby Isaac PetscSection section; 279552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 280552f7358SJed 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. */ 2819566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)X, &vecname)); 282552f7358SJed Brown } 2839566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 284e630c359SToby Isaac if (!dmX) dmX = dm; 2859566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 2869566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 2879566063dSJacob Faibussowitsch if (cEnd > cStart) PetscCall(PetscSectionGetDof(section, cStart, &bs)); 2889566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 289e630c359SToby Isaac field = 0; 290e630c359SToby Isaac if (link->field >= 0) { 291e630c359SToby Isaac field = link->field; 292e630c359SToby Isaac nfields = field + 1; 293e630c359SToby Isaac } 294e630c359SToby Isaac for (i = 0; field < (nfields ? nfields : 1); field++) { 2951cfafdd3SJed Brown PetscInt fbs, j; 296a00cdb45SToby Isaac PetscFV fv = NULL; 2970298fd71SBarry Smith const char *fieldname = NULL; 2981cfafdd3SJed Brown char buf[256]; 299e630c359SToby Isaac PetscBool vector; 300*6b8451b4SStefano Zampini 3017ded4cbaSJed Brown if (nfields) { /* We have user-defined fields/components */ 3029566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, cStart, field, &fbs)); 3039566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(section, field, &fieldname)); 3047ded4cbaSJed Brown } else fbs = bs; /* Say we have one field with 'bs' components */ 305*6b8451b4SStefano Zampini PetscCall(DMGetFieldIfFV_Private(dmX, field, &fv)); 306e630c359SToby Isaac if (nfields && !fieldname) { 30763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "CellField%" PetscInt_FMT, field)); 308552f7358SJed Brown fieldname = buf; 309552f7358SJed Brown } 310e630c359SToby Isaac vector = PETSC_FALSE; 311e630c359SToby Isaac if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) { 312e630c359SToby Isaac vector = PETSC_TRUE; 31363a3b9bcSJacob Faibussowitsch PetscCheck(fbs <= 3, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_SIZ, "Cell vector fields can have at most 3 components, %" PetscInt_FMT " given", fbs); 314e630c359SToby Isaac for (j = 0; j < fbs; j++) { 315e630c359SToby Isaac const char *compName = NULL; 316e630c359SToby Isaac if (fv) { 3179566063dSJacob Faibussowitsch PetscCall(PetscFVGetComponentName(fv, j, &compName)); 318e630c359SToby Isaac if (compName) break; 319e630c359SToby Isaac } 320e630c359SToby Isaac } 321e630c359SToby Isaac if (j < fbs) vector = PETSC_FALSE; 322e630c359SToby Isaac } 323e630c359SToby Isaac if (vector) { 324955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 3250f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s.Re\" 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); 3270f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s.Im\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 3286f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * 3 * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 329955d60d1SToby Isaac #else 3300f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 3316f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * 3 * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 332955d60d1SToby Isaac #endif 333e630c359SToby Isaac i += fbs; 334e630c359SToby Isaac } else { 3351cfafdd3SJed Brown for (j = 0; j < fbs; j++) { 336a00cdb45SToby Isaac const char *compName = NULL; 337955d60d1SToby Isaac char finalname[256]; 33848a46eb9SPierre Jolivet if (fv) PetscCall(PetscFVGetComponentName(fv, j, &compName)); 339a00cdb45SToby Isaac if (compName) { 3409566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname, 255, "%s%s.%s", vecname, fieldname, compName)); 3419371c9d4SSatish Balay } else if (fbs > 1) { 34263a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname, 255, "%s%s.%" PetscInt_FMT, vecname, fieldname, j)); 343e630c359SToby Isaac } else { 3449566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname, 255, "%s%s", vecname, fieldname)); 345a00cdb45SToby Isaac } 346955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 3470f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.Re\" 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); 3490f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset)); 3506f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 351955d60d1SToby Isaac #else 3520f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset)); 3536f98e5d0SMatthew G. Knepley boffset += gpiece[r].ncells * sizeof(PetscVTUReal) + (gpiece[r].ncells ? sizeof(PetscInt64) : 0); 354955d60d1SToby Isaac #endif 3551cfafdd3SJed Brown i++; 356552f7358SJed Brown } 357552f7358SJed Brown } 358e630c359SToby Isaac } 35963a3b9bcSJacob Faibussowitsch //PetscCheck(i == bs,comm,PETSC_ERR_PLIB,"Total number of field components %" PetscInt_FMT " != block size %" PetscInt_FMT,i,bs); 3601cfafdd3SJed Brown } 3619566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </CellData>\n")); 362552f7358SJed Brown 363552f7358SJed Brown /* 364552f7358SJed Brown * Point Data headers 365552f7358SJed Brown */ 3669566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <PointData>\n")); 367552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 368552f7358SJed Brown Vec X = (Vec)link->vec; 369e630c359SToby Isaac DM dmX; 3706030a318SStefano Zampini PetscInt bs = 1, nfields, field; 371552f7358SJed Brown const char *vecname = ""; 372e630c359SToby Isaac PetscSection section; 373552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 374552f7358SJed 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. */ 3759566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)X, &vecname)); 376552f7358SJed Brown } 3779566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 378e630c359SToby Isaac if (!dmX) dmX = dm; 3799566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 3809566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 3819566063dSJacob Faibussowitsch if (vEnd > vStart) PetscCall(PetscSectionGetDof(section, vStart, &bs)); 3829566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 383e630c359SToby Isaac field = 0; 384e630c359SToby Isaac if (link->field >= 0) { 385e630c359SToby Isaac field = link->field; 386e630c359SToby Isaac nfields = field + 1; 387e630c359SToby Isaac } 38863bfac88SBarry Smith for (; field < (nfields ? nfields : 1); field++) { 3891cfafdd3SJed Brown PetscInt fbs, j; 3901cfafdd3SJed Brown const char *fieldname = NULL; 3911cfafdd3SJed Brown char buf[256]; 3927ded4cbaSJed Brown if (nfields) { /* We have user-defined fields/components */ 3939566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, vStart, field, &fbs)); 3949566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(section, field, &fieldname)); 3957ded4cbaSJed Brown } else fbs = bs; /* Say we have one field with 'bs' components */ 396e630c359SToby Isaac if (nfields && !fieldname) { 39763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "PointField%" PetscInt_FMT, field)); 3981cfafdd3SJed Brown fieldname = buf; 3991cfafdd3SJed Brown } 400e630c359SToby Isaac if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) { 40163a3b9bcSJacob Faibussowitsch PetscCheck(fbs <= 3, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_SIZ, "Point vector fields can have at most 3 components, %" PetscInt_FMT " given", fbs); 402955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 4030f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s.Re\" 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); 4050f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s.Im\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 4066f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 407955d60d1SToby Isaac #else 4080f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 4096f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 410955d60d1SToby Isaac #endif 411e630c359SToby Isaac } else { 4121cfafdd3SJed Brown for (j = 0; j < fbs; j++) { 413b778fa18SValeria Barra const char *compName = NULL; 414955d60d1SToby Isaac char finalname[256]; 4159566063dSJacob Faibussowitsch PetscCall(PetscSectionGetComponentName(section, field, j, &compName)); 4169566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname, 255, "%s%s.%s", vecname, fieldname, compName)); 417955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 4180f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.Re\" 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); 4200f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset)); 4216f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 422955d60d1SToby Isaac #else 4230f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, finalname, boffset)); 4246f98e5d0SMatthew G. Knepley boffset += gpiece[r].nvertices * sizeof(PetscVTUReal) + (gpiece[r].nvertices ? sizeof(PetscInt64) : 0); 425955d60d1SToby Isaac #endif 426552f7358SJed Brown } 427552f7358SJed Brown } 4281cfafdd3SJed Brown } 429e630c359SToby Isaac } 4309566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </PointData>\n")); 4319566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </Piece>\n")); 432552f7358SJed Brown } 433552f7358SJed Brown } 434552f7358SJed Brown 4359566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </UnstructuredGrid>\n")); 4369566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <AppendedData encoding=\"raw\">\n")); 4379566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "_")); 438552f7358SJed Brown 439dd400576SPatrick Sanan if (rank == 0) { 440552f7358SJed Brown PetscInt maxsize = 0; 4416497c311SBarry Smith for (PetscMPIInt r = 0; r < size; r++) { 442955d60d1SToby Isaac maxsize = PetscMax(maxsize, (PetscInt)(gpiece[r].nvertices * 3 * sizeof(PetscVTUReal))); 443955d60d1SToby Isaac maxsize = PetscMax(maxsize, (PetscInt)(gpiece[r].ncells * 3 * sizeof(PetscVTUReal))); 444552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt)(gpiece[r].nconn * sizeof(PetscVTKInt))); 445552f7358SJed Brown } 4469566063dSJacob Faibussowitsch PetscCall(PetscMalloc(maxsize, &buffer)); 447552f7358SJed Brown } 4486497c311SBarry Smith for (PetscMPIInt r = 0; r < size; r++) { 449552f7358SJed Brown if (r == rank) { 450552f7358SJed Brown PetscInt nsend; 451552f7358SJed Brown { /* Position */ 4526858538eSMatthew G. Knepley const PetscScalar *x, *cx = NULL; 453955d60d1SToby Isaac PetscVTUReal *y = NULL; 4546858538eSMatthew G. Knepley Vec coords, cellCoords; 455955d60d1SToby Isaac PetscBool copy; 4562e529ec8SStefano Zampini 4579566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coords)); 4589566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coords, &x)); 4596858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinatesLocal(dm, &cellCoords)); 4606858538eSMatthew G. Knepley if (cellCoords) PetscCall(VecGetArrayRead(cellCoords, &cx)); 461955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 462955d60d1SToby Isaac copy = PETSC_TRUE; 463955d60d1SToby Isaac #else 464955d60d1SToby Isaac copy = (PetscBool)(dimEmbed != 3 || localized || (sizeof(PetscReal) != sizeof(PetscVTUReal))); 465955d60d1SToby Isaac #endif 466955d60d1SToby Isaac if (copy) { 4679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.nvertices * 3, &y)); 46819003fb5SStefano Zampini if (localized) { 46919003fb5SStefano Zampini PetscInt cnt; 47019003fb5SStefano Zampini for (c = cStart, cnt = 0; c < cEnd; c++) { 47119003fb5SStefano Zampini PetscInt off, dof; 47219003fb5SStefano Zampini 4736858538eSMatthew G. Knepley PetscCall(PetscSectionGetDof(cellCoordSection, c, &dof)); 47419003fb5SStefano Zampini if (!dof) { 47519003fb5SStefano Zampini PetscInt *closure = NULL; 47619003fb5SStefano Zampini PetscInt closureSize; 47719003fb5SStefano Zampini 4789566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 47919003fb5SStefano Zampini for (v = 0; v < closureSize * 2; v += 2) { 48019003fb5SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 4819566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, closure[v], &off)); 48219003fb5SStefano Zampini if (dimEmbed != 3) { 483955d60d1SToby Isaac y[cnt * 3 + 0] = (PetscVTUReal)PetscRealPart(x[off + 0]); 484955d60d1SToby Isaac y[cnt * 3 + 1] = (PetscVTUReal)((dimEmbed > 1) ? PetscRealPart(x[off + 1]) : 0.0); 485955d60d1SToby Isaac y[cnt * 3 + 2] = (PetscVTUReal)0.0; 48619003fb5SStefano Zampini } else { 487955d60d1SToby Isaac y[cnt * 3 + 0] = (PetscVTUReal)PetscRealPart(x[off + 0]); 488955d60d1SToby Isaac y[cnt * 3 + 1] = (PetscVTUReal)PetscRealPart(x[off + 1]); 489955d60d1SToby Isaac y[cnt * 3 + 2] = (PetscVTUReal)PetscRealPart(x[off + 2]); 49019003fb5SStefano Zampini } 49119003fb5SStefano Zampini cnt++; 49219003fb5SStefano Zampini } 49319003fb5SStefano Zampini } 4949566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 49519003fb5SStefano Zampini } else { 4966858538eSMatthew G. Knepley PetscCall(PetscSectionGetOffset(cellCoordSection, c, &off)); 49719003fb5SStefano Zampini if (dimEmbed != 3) { 49819003fb5SStefano Zampini for (i = 0; i < dof / dimEmbed; i++) { 4996858538eSMatthew G. Knepley y[cnt * 3 + 0] = (PetscVTUReal)PetscRealPart(cx[off + i * dimEmbed + 0]); 5006858538eSMatthew G. Knepley y[cnt * 3 + 1] = (PetscVTUReal)((dimEmbed > 1) ? PetscRealPart(cx[off + i * dimEmbed + 1]) : 0.0); 501955d60d1SToby Isaac y[cnt * 3 + 2] = (PetscVTUReal)0.0; 50219003fb5SStefano Zampini cnt++; 50319003fb5SStefano Zampini } 50419003fb5SStefano Zampini } else { 505ad540459SPierre Jolivet for (i = 0; i < dof; i++) y[cnt * 3 + i] = (PetscVTUReal)PetscRealPart(cx[off + i]); 50619003fb5SStefano Zampini cnt += dof / dimEmbed; 50719003fb5SStefano Zampini } 50819003fb5SStefano Zampini } 50919003fb5SStefano Zampini } 51008401ef6SPierre Jolivet PetscCheck(cnt == piece.nvertices, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 51119003fb5SStefano Zampini } else { 512552f7358SJed Brown for (i = 0; i < piece.nvertices; i++) { 513955d60d1SToby Isaac y[i * 3 + 0] = (PetscVTUReal)PetscRealPart(x[i * dimEmbed + 0]); 514955d60d1SToby Isaac y[i * 3 + 1] = (PetscVTUReal)((dimEmbed > 1) ? PetscRealPart(x[i * dimEmbed + 1]) : 0.); 515955d60d1SToby Isaac y[i * 3 + 2] = (PetscVTUReal)((dimEmbed > 2) ? PetscRealPart(x[i * dimEmbed + 2]) : 0.); 51619003fb5SStefano Zampini } 517552f7358SJed Brown } 518552f7358SJed Brown } 5192e529ec8SStefano Zampini nsend = piece.nvertices * 3; 5209566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, copy ? (const void *)y : (const void *)x, buffer, nsend, MPIU_VTUREAL, tag)); 5219566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 5229566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coords, &x)); 5236858538eSMatthew G. Knepley if (cellCoords) PetscCall(VecRestoreArrayRead(cellCoords, &cx)); 524552f7358SJed Brown } 525552f7358SJed Brown { /* Connectivity, offsets, types */ 5263e869605SMatthew G. Knepley PetscVTKInt *connectivity = NULL, *offsets = NULL; 5273e869605SMatthew G. Knepley PetscVTKType *types = NULL; 5289566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKConnectivity(dm, localized, &piece, &connectivity, &offsets, &types)); 5299566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, connectivity, buffer, piece.nconn, MPI_INT, tag)); 5309566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, offsets, buffer, piece.ncells, MPI_INT, tag)); 5319566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, types, buffer, piece.ncells, MPI_CHAR, tag)); 5329566063dSJacob Faibussowitsch PetscCall(PetscFree3(connectivity, offsets, types)); 533552f7358SJed Brown } 534552f7358SJed Brown { /* Owners (cell data) */ 535552f7358SJed Brown PetscVTKInt *owners; 536c29ce622SStefano Zampini PetscMPIInt orank; 537c29ce622SStefano Zampini 538c29ce622SStefano Zampini PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &orank)); 5399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.ncells, &owners)); 540c29ce622SStefano Zampini for (i = 0; i < piece.ncells; i++) owners[i] = orank; 5419566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, owners, buffer, piece.ncells, MPI_INT, tag)); 5429566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 543552f7358SJed Brown } 544552f7358SJed Brown /* Cell data */ 545552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 546552f7358SJed Brown Vec X = (Vec)link->vec; 547e630c359SToby Isaac DM dmX; 548552f7358SJed Brown const PetscScalar *x; 549955d60d1SToby Isaac PetscVTUReal *y; 5506030a318SStefano Zampini PetscInt bs = 1, nfields, field; 551e630c359SToby Isaac PetscSection section = NULL; 552e630c359SToby Isaac 553552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 5549566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 555e630c359SToby Isaac if (!dmX) dmX = dm; 5569566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 5579566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 5589566063dSJacob Faibussowitsch if (cEnd > cStart) PetscCall(PetscSectionGetDof(section, cStart, &bs)); 5599566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 560e630c359SToby Isaac field = 0; 561e630c359SToby Isaac if (link->field >= 0) { 562e630c359SToby Isaac field = link->field; 563e630c359SToby Isaac nfields = field + 1; 564e630c359SToby Isaac } 5659566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 5669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.ncells * 3, &y)); 56763bfac88SBarry Smith for (; field < (nfields ? nfields : 1); field++) { 568e630c359SToby Isaac PetscInt fbs, j; 569e630c359SToby Isaac PetscFV fv = NULL; 570e630c359SToby Isaac PetscBool vector; 571*6b8451b4SStefano Zampini 5726030a318SStefano Zampini if (nfields && cEnd > cStart) { /* We have user-defined fields/components */ 5739566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, cStart, field, &fbs)); 574e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 575*6b8451b4SStefano Zampini PetscCall(DMGetFieldIfFV_Private(dmX, field, &fv)); 576e630c359SToby Isaac vector = PETSC_FALSE; 577e630c359SToby Isaac if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) { 578e630c359SToby Isaac vector = PETSC_TRUE; 579e630c359SToby Isaac for (j = 0; j < fbs; j++) { 580e630c359SToby Isaac const char *compName = NULL; 581e630c359SToby Isaac if (fv) { 5829566063dSJacob Faibussowitsch PetscCall(PetscFVGetComponentName(fv, j, &compName)); 583e630c359SToby Isaac if (compName) break; 584e630c359SToby Isaac } 585e630c359SToby Isaac } 586e630c359SToby Isaac if (j < fbs) vector = PETSC_FALSE; 587e630c359SToby Isaac } 588e630c359SToby Isaac if (vector) { 589955d60d1SToby Isaac PetscInt cnt, l; 590955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 591552f7358SJed Brown for (c = cStart, cnt = 0; c < cEnd; c++) { 592e630c359SToby Isaac const PetscScalar *xpoint; 593e630c359SToby Isaac PetscInt off, j; 594e630c359SToby Isaac 595552f7358SJed Brown if (hasLabel) { /* Ignore some cells */ 596552f7358SJed Brown PetscInt value; 5979566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dmX, "vtk", c, &value)); 598552f7358SJed Brown if (value != 1) continue; 599552f7358SJed Brown } 600e630c359SToby Isaac if (nfields) { 6019566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, c, field, &off)); 602e630c359SToby Isaac } else { 6039566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, c, &off)); 604e630c359SToby Isaac } 605e630c359SToby Isaac xpoint = &x[off]; 606ad540459SPierre Jolivet for (j = 0; j < fbs; j++) y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j])); 607e630c359SToby Isaac for (; j < 3; j++) y[cnt++] = 0.; 608e630c359SToby Isaac } 60908401ef6SPierre Jolivet PetscCheck(cnt == piece.ncells * 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 6109566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.ncells * 3, MPIU_VTUREAL, tag)); 611955d60d1SToby Isaac } 612e630c359SToby Isaac } else { 613e630c359SToby Isaac for (i = 0; i < fbs; i++) { 614955d60d1SToby Isaac PetscInt cnt, l; 615955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 616e630c359SToby Isaac for (c = cStart, cnt = 0; c < cEnd; c++) { 617e630c359SToby Isaac const PetscScalar *xpoint; 618e630c359SToby Isaac PetscInt off; 619e630c359SToby Isaac 620e630c359SToby Isaac if (hasLabel) { /* Ignore some cells */ 621e630c359SToby Isaac PetscInt value; 6229566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dmX, "vtk", c, &value)); 623e630c359SToby Isaac if (value != 1) continue; 624e630c359SToby Isaac } 625e630c359SToby Isaac if (nfields) { 6269566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, c, field, &off)); 627e630c359SToby Isaac } else { 6289566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, c, &off)); 629e630c359SToby Isaac } 630e630c359SToby Isaac xpoint = &x[off]; 631955d60d1SToby Isaac y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i])); 632552f7358SJed Brown } 63308401ef6SPierre Jolivet PetscCheck(cnt == piece.ncells, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 6349566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.ncells, MPIU_VTUREAL, tag)); 635955d60d1SToby Isaac } 636552f7358SJed Brown } 637e630c359SToby Isaac } 638e630c359SToby Isaac } 6399566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 6409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 641552f7358SJed Brown } 642e630c359SToby Isaac /* point data */ 643552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 644552f7358SJed Brown Vec X = (Vec)link->vec; 645e630c359SToby Isaac DM dmX; 646552f7358SJed Brown const PetscScalar *x; 647955d60d1SToby Isaac PetscVTUReal *y; 6486030a318SStefano Zampini PetscInt bs = 1, nfields, field; 649e630c359SToby Isaac PetscSection section = NULL; 650e630c359SToby Isaac 651552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 6529566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 653e630c359SToby Isaac if (!dmX) dmX = dm; 6549566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 6559566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 6569566063dSJacob Faibussowitsch if (vEnd > vStart) PetscCall(PetscSectionGetDof(section, vStart, &bs)); 6579566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 658e630c359SToby Isaac field = 0; 659e630c359SToby Isaac if (link->field >= 0) { 660e630c359SToby Isaac field = link->field; 661e630c359SToby Isaac nfields = field + 1; 662e630c359SToby Isaac } 6639566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 6649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.nvertices * 3, &y)); 66563bfac88SBarry Smith for (; field < (nfields ? nfields : 1); field++) { 666e630c359SToby Isaac PetscInt fbs, j; 6676030a318SStefano Zampini if (nfields && vEnd > vStart) { /* We have user-defined fields/components */ 6689566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, vStart, field, &fbs)); 669e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 670e630c359SToby Isaac if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) { 671955d60d1SToby Isaac PetscInt cnt, l; 672955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 6732e529ec8SStefano Zampini if (!localized) { 674552f7358SJed Brown for (v = vStart, cnt = 0; v < vEnd; v++) { 675e630c359SToby Isaac PetscInt off; 676e630c359SToby Isaac const PetscScalar *xpoint; 677e630c359SToby Isaac 678e630c359SToby Isaac if (nfields) { 6799566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, v, field, &off)); 680e630c359SToby Isaac } else { 6819566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, v, &off)); 682e630c359SToby Isaac } 683e630c359SToby Isaac xpoint = &x[off]; 684ad540459SPierre Jolivet for (j = 0; j < fbs; j++) y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j])); 685e630c359SToby Isaac for (; j < 3; j++) y[cnt++] = 0.; 686e630c359SToby Isaac } 687e630c359SToby Isaac } else { 688e630c359SToby Isaac for (c = cStart, cnt = 0; c < cEnd; c++) { 689e630c359SToby Isaac PetscInt *closure = NULL; 690e630c359SToby Isaac PetscInt closureSize, off; 691e630c359SToby Isaac 6929566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 693e630c359SToby Isaac for (v = 0, off = 0; v < closureSize * 2; v += 2) { 694e630c359SToby Isaac if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 695e630c359SToby Isaac PetscInt voff; 696e630c359SToby Isaac const PetscScalar *xpoint; 697e630c359SToby Isaac 698e630c359SToby Isaac if (nfields) { 6999566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, closure[v], field, &voff)); 700e630c359SToby Isaac } else { 7019566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, closure[v], &voff)); 702e630c359SToby Isaac } 703e630c359SToby Isaac xpoint = &x[voff]; 704ad540459SPierre Jolivet for (j = 0; j < fbs; j++) y[cnt + off++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j])); 705e630c359SToby Isaac for (; j < 3; j++) y[cnt + off++] = 0.; 706e630c359SToby Isaac } 707e630c359SToby Isaac } 708e630c359SToby Isaac cnt += off; 7099566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 710e630c359SToby Isaac } 711e630c359SToby Isaac } 71208401ef6SPierre Jolivet PetscCheck(cnt == piece.nvertices * 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 7139566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.nvertices * 3, MPIU_VTUREAL, tag)); 714955d60d1SToby Isaac } 715e630c359SToby Isaac } else { 716e630c359SToby Isaac for (i = 0; i < fbs; i++) { 717955d60d1SToby Isaac PetscInt cnt, l; 718955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 719e630c359SToby Isaac if (!localized) { 720e630c359SToby Isaac for (v = vStart, cnt = 0; v < vEnd; v++) { 721e630c359SToby Isaac PetscInt off; 722e630c359SToby Isaac const PetscScalar *xpoint; 723e630c359SToby Isaac 724e630c359SToby Isaac if (nfields) { 7259566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, v, field, &off)); 726e630c359SToby Isaac } else { 7279566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, v, &off)); 728e630c359SToby Isaac } 729e630c359SToby Isaac xpoint = &x[off]; 730955d60d1SToby Isaac y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i])); 731552f7358SJed Brown } 7322e529ec8SStefano Zampini } else { 7332e529ec8SStefano Zampini for (c = cStart, cnt = 0; c < cEnd; c++) { 7342e529ec8SStefano Zampini PetscInt *closure = NULL; 7352e529ec8SStefano Zampini PetscInt closureSize, off; 7362e529ec8SStefano Zampini 7379566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 7382e529ec8SStefano Zampini for (v = 0, off = 0; v < closureSize * 2; v += 2) { 7392e529ec8SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 740e630c359SToby Isaac PetscInt voff; 741e630c359SToby Isaac const PetscScalar *xpoint; 7422e529ec8SStefano Zampini 743e630c359SToby Isaac if (nfields) { 7449566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, closure[v], field, &voff)); 745e630c359SToby Isaac } else { 7469566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, closure[v], &voff)); 747e630c359SToby Isaac } 748e630c359SToby Isaac xpoint = &x[voff]; 749955d60d1SToby Isaac y[cnt + off++] = (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i])); 7502e529ec8SStefano Zampini } 7512e529ec8SStefano Zampini } 7529bda96f6SStefano Zampini cnt += off; 7539566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 7542e529ec8SStefano Zampini } 7552e529ec8SStefano Zampini } 75608401ef6SPierre Jolivet PetscCheck(cnt == piece.nvertices, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 7579566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.nvertices, MPIU_VTUREAL, tag)); 758955d60d1SToby Isaac } 759552f7358SJed Brown } 760e630c359SToby Isaac } 761e630c359SToby Isaac } 7629566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 7639566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 764552f7358SJed Brown } 765dd400576SPatrick Sanan } else if (rank == 0) { 766955d60d1SToby Isaac PetscInt l; 767955d60d1SToby Isaac 7689566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nvertices * 3, MPIU_VTUREAL, tag)); /* positions */ 7699566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nconn, MPI_INT, tag)); /* connectivity */ 7709566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPI_INT, tag)); /* offsets */ 7719566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPI_CHAR, tag)); /* types */ 7729566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPI_INT, tag)); /* owner rank (cells) */ 773552f7358SJed Brown /* all cell data */ 774552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 775e630c359SToby Isaac Vec X = (Vec)link->vec; 7766030a318SStefano Zampini PetscInt bs = 1, nfields, field; 777e630c359SToby Isaac DM dmX; 778e630c359SToby Isaac PetscSection section = NULL; 779e630c359SToby Isaac 780552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 7819566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 782e630c359SToby Isaac if (!dmX) dmX = dm; 7839566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 7849566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 7859566063dSJacob Faibussowitsch if (cEnd > cStart) PetscCall(PetscSectionGetDof(section, cStart, &bs)); 7869566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 787e630c359SToby Isaac field = 0; 788e630c359SToby Isaac if (link->field >= 0) { 789e630c359SToby Isaac field = link->field; 790e630c359SToby Isaac nfields = field + 1; 791e630c359SToby Isaac } 79263bfac88SBarry Smith for (; field < (nfields ? nfields : 1); field++) { 793e630c359SToby Isaac PetscInt fbs, j; 794e630c359SToby Isaac PetscFV fv = NULL; 795e630c359SToby Isaac PetscBool vector; 796*6b8451b4SStefano Zampini 7976030a318SStefano Zampini if (nfields && cEnd > cStart) { /* We have user-defined fields/components */ 7989566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, cStart, field, &fbs)); 799e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 800*6b8451b4SStefano Zampini PetscCall(DMGetFieldIfFV_Private(dmX, field, &fv)); 801e630c359SToby Isaac vector = PETSC_FALSE; 802e630c359SToby Isaac if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) { 803e630c359SToby Isaac vector = PETSC_TRUE; 804e630c359SToby Isaac for (j = 0; j < fbs; j++) { 805e630c359SToby Isaac const char *compName = NULL; 806e630c359SToby Isaac if (fv) { 8079566063dSJacob Faibussowitsch PetscCall(PetscFVGetComponentName(fv, j, &compName)); 808e630c359SToby Isaac if (compName) break; 809e630c359SToby Isaac } 810e630c359SToby Isaac } 811e630c359SToby Isaac if (j < fbs) vector = PETSC_FALSE; 812e630c359SToby Isaac } 813e630c359SToby Isaac if (vector) { 81448a46eb9SPierre 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)); 815e630c359SToby Isaac } else { 816e630c359SToby Isaac for (i = 0; i < fbs; i++) { 81748a46eb9SPierre Jolivet for (l = 0; l < loops_per_scalar; l++) PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPIU_VTUREAL, tag)); 818552f7358SJed Brown } 819552f7358SJed Brown } 820e630c359SToby Isaac } 821e630c359SToby Isaac } 822552f7358SJed Brown /* all point data */ 823552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 824e630c359SToby Isaac Vec X = (Vec)link->vec; 825e630c359SToby Isaac DM dmX; 8266030a318SStefano Zampini PetscInt bs = 1, nfields, field; 827e630c359SToby Isaac PetscSection section = NULL; 828e630c359SToby Isaac 829552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 8309566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 831e630c359SToby Isaac if (!dmX) dmX = dm; 8329566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 8339566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 8349566063dSJacob Faibussowitsch if (vEnd > vStart) PetscCall(PetscSectionGetDof(section, vStart, &bs)); 8359566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 836e630c359SToby Isaac field = 0; 837e630c359SToby Isaac if (link->field >= 0) { 838e630c359SToby Isaac field = link->field; 839e630c359SToby Isaac nfields = field + 1; 840e630c359SToby Isaac } 84163bfac88SBarry Smith for (; field < (nfields ? nfields : 1); field++) { 842e630c359SToby Isaac PetscInt fbs; 8436030a318SStefano Zampini if (nfields && vEnd > vStart) { /* We have user-defined fields/components */ 8449566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, vStart, field, &fbs)); 845e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 846e630c359SToby Isaac if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) { 84748a46eb9SPierre 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)); 848e630c359SToby Isaac } else { 849e630c359SToby Isaac for (i = 0; i < fbs; i++) { 85048a46eb9SPierre Jolivet for (l = 0; l < loops_per_scalar; l++) PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nvertices, MPIU_VTUREAL, tag)); 851552f7358SJed Brown } 852552f7358SJed Brown } 853552f7358SJed Brown } 854552f7358SJed Brown } 855e630c359SToby Isaac } 856e630c359SToby Isaac } 8579566063dSJacob Faibussowitsch PetscCall(PetscFree(gpiece)); 8589566063dSJacob Faibussowitsch PetscCall(PetscFree(buffer)); 8599566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "\n </AppendedData>\n")); 8609566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "</VTKFile>\n")); 8619566063dSJacob Faibussowitsch PetscCall(PetscFClose(comm, fp)); 862c29ce622SStefano Zampini finalize: 86317bfc685SStefano Zampini /* this code sends to rank 0 that writes. 86417bfc685SStefano Zampini It may lead to very unbalanced log_view timings 86517bfc685SStefano Zampini of the next PETSc function logged. 86617bfc685SStefano Zampini Since this call is not performance critical, we 86717bfc685SStefano Zampini issue a barrier here to synchronize the processes */ 86817bfc685SStefano Zampini PetscCall(PetscBarrier((PetscObject)viewer)); 8693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 870552f7358SJed Brown } 871