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 269371c9d4SSatish Balay static PetscErrorCode TransferWrite(MPI_Comm comm, PetscViewer viewer, FILE *fp, PetscMPIInt srank, PetscMPIInt root, const void *send, void *recv, PetscMPIInt count, MPI_Datatype mpidatatype, PetscMPIInt tag) { 27552f7358SJed Brown PetscMPIInt rank; 28552f7358SJed Brown 29552f7358SJed Brown PetscFunctionBegin; 309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 31552f7358SJed Brown if (rank == srank && rank != root) { 329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void *)send, count, mpidatatype, root, tag, comm)); 33552f7358SJed Brown } else if (rank == root) { 34552f7358SJed Brown const void *buffer; 35552f7358SJed Brown if (root == srank) { /* self */ 36552f7358SJed Brown buffer = send; 37552f7358SJed Brown } else { 38552f7358SJed Brown MPI_Status status; 39552f7358SJed Brown PetscMPIInt nrecv; 409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(recv, count, mpidatatype, srank, tag, comm, &status)); 419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status, mpidatatype, &nrecv)); 4208401ef6SPierre Jolivet PetscCheck(count == nrecv, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch"); 43552f7358SJed Brown buffer = recv; 44552f7358SJed Brown } 459566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, buffer, count, mpidatatype)); 46552f7358SJed Brown } 47552f7358SJed Brown PetscFunctionReturn(0); 48552f7358SJed Brown } 49552f7358SJed Brown 509371c9d4SSatish Balay static PetscErrorCode DMPlexGetVTKConnectivity(DM dm, PetscBool localized, PieceInfo *piece, PetscVTKInt **oconn, PetscVTKInt **ooffsets, PetscVTKType **otypes) { 516858538eSMatthew G. Knepley PetscSection coordSection, cellCoordSection; 52552f7358SJed Brown PetscVTKInt *conn, *offsets; 53552f7358SJed Brown PetscVTKType *types; 542f029394SStefano Zampini PetscInt dim, vStart, vEnd, cStart, cEnd, pStart, pEnd, cellHeight, numLabelCells, hasLabel, c, v, countcell, countconn; 55552f7358SJed Brown 56552f7358SJed Brown PetscFunctionBegin; 579566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(piece->nconn, &conn, piece->ncells, &offsets, piece->ncells, &types)); 589566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 596858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinateSection(dm, &cellCoordSection)); 609566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 619566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 629566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight)); 639566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd)); 649566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 659566063dSJacob Faibussowitsch PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells)); 66552f7358SJed Brown hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 67552f7358SJed Brown 68552f7358SJed Brown countcell = 0; 69552f7358SJed Brown countconn = 0; 70552f7358SJed Brown for (c = cStart; c < cEnd; ++c) { 7119003fb5SStefano Zampini PetscInt nverts, dof = 0, celltype, startoffset, nC = 0; 72552f7358SJed Brown 73552f7358SJed Brown if (hasLabel) { 74552f7358SJed Brown PetscInt value; 75552f7358SJed Brown 769566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dm, "vtk", c, &value)); 77552f7358SJed Brown if (value != 1) continue; 78552f7358SJed Brown } 79552f7358SJed Brown startoffset = countconn; 806858538eSMatthew G. Knepley if (localized) PetscCall(PetscSectionGetDof(cellCoordSection, c, &dof)); 8119003fb5SStefano Zampini if (!dof) { 822e529ec8SStefano Zampini PetscInt *closure = NULL; 832e529ec8SStefano Zampini PetscInt closureSize; 842e529ec8SStefano Zampini 859566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 86552f7358SJed Brown for (v = 0; v < closureSize * 2; v += 2) { 87552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 8819003fb5SStefano Zampini if (!localized) conn[countconn++] = closure[v] - vStart; 8919003fb5SStefano Zampini else conn[countconn++] = startoffset + nC; 90724b5320SMatthew G. Knepley ++nC; 91552f7358SJed Brown } 92552f7358SJed Brown } 939566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 942e529ec8SStefano Zampini } else { 952e529ec8SStefano Zampini for (nC = 0; nC < dof / dim; nC++) conn[countconn++] = startoffset + nC; 962e529ec8SStefano Zampini } 9796ca5757SLisandro Dalcin 9896ca5757SLisandro Dalcin { 9996ca5757SLisandro Dalcin PetscInt n = PetscMin(nC, 8), s = countconn - nC, i, cone[8]; 10096ca5757SLisandro Dalcin for (i = 0; i < n; ++i) cone[i] = conn[s + i]; 1019566063dSJacob Faibussowitsch PetscCall(DMPlexReorderCell(dm, c, cone)); 10296ca5757SLisandro Dalcin for (i = 0; i < n; ++i) conn[s + i] = (int)cone[i]; 10396ca5757SLisandro Dalcin } 1048865f1eaSKarl Rupp 105552f7358SJed Brown offsets[countcell] = countconn; 1068865f1eaSKarl Rupp 107552f7358SJed Brown nverts = countconn - startoffset; 1089566063dSJacob Faibussowitsch PetscCall(DMPlexVTKGetCellType_Internal(dm, dim, nverts, &celltype)); 1098865f1eaSKarl Rupp 110552f7358SJed Brown types[countcell] = celltype; 111552f7358SJed Brown countcell++; 112552f7358SJed Brown } 11308401ef6SPierre Jolivet PetscCheck(countcell == piece->ncells, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent cell count"); 11408401ef6SPierre Jolivet PetscCheck(countconn == piece->nconn, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent connectivity count"); 115552f7358SJed Brown *oconn = conn; 116552f7358SJed Brown *ooffsets = offsets; 117552f7358SJed Brown *otypes = types; 118552f7358SJed Brown PetscFunctionReturn(0); 119552f7358SJed Brown } 120552f7358SJed Brown 121552f7358SJed Brown /* 122552f7358SJed Brown Write all fields that have been provided to the viewer 123552f7358SJed Brown Multi-block XML format with binary appended data. 124552f7358SJed Brown */ 1259371c9d4SSatish Balay PetscErrorCode DMPlexVTKWriteAll_VTU(DM dm, PetscViewer viewer) { 126ce94432eSBarry Smith MPI_Comm comm; 1276858538eSMatthew G. Knepley PetscSection coordSection, cellCoordSection; 128552f7358SJed Brown PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data; 129552f7358SJed Brown PetscViewerVTKObjectLink link; 130552f7358SJed Brown FILE *fp; 131552f7358SJed Brown PetscMPIInt rank, size, tag; 1322f029394SStefano Zampini PetscInt dimEmbed, cellHeight, cStart, cEnd, vStart, vEnd, numLabelCells, hasLabel, c, v, r, i; 1332e529ec8SStefano Zampini PetscBool localized; 1340298fd71SBarry Smith PieceInfo piece, *gpiece = NULL; 1350298fd71SBarry Smith void *buffer = NULL; 13630815ce0SLisandro Dalcin const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian"; 137955d60d1SToby Isaac PetscInt loops_per_scalar; 138552f7358SJed Brown 139552f7358SJed Brown PetscFunctionBegin; 1409566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 141552f7358SJed Brown #if defined(PETSC_USE_COMPLEX) 142955d60d1SToby Isaac loops_per_scalar = 2; 143955d60d1SToby Isaac #else 144955d60d1SToby Isaac loops_per_scalar = 1; 145552f7358SJed Brown #endif 1469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 1479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1489566063dSJacob Faibussowitsch PetscCall(PetscCommGetNewTag(comm, &tag)); 149552f7358SJed Brown 1509566063dSJacob Faibussowitsch PetscCall(PetscFOpen(comm, vtk->filename, "wb", &fp)); 1519566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "<?xml version=\"1.0\"?>\n")); 1529566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"%s\">\n", byte_order)); 1539566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <UnstructuredGrid>\n")); 154552f7358SJed Brown 1559566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dimEmbed)); 1569566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight)); 1579566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd)); 1589566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1599566063dSJacob Faibussowitsch PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells)); 1609566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalized(dm, &localized)); 1619566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 1626858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinateSection(dm, &cellCoordSection)); 1638865f1eaSKarl Rupp 164552f7358SJed Brown hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 1652e529ec8SStefano Zampini piece.nvertices = 0; 166552f7358SJed Brown piece.ncells = 0; 167552f7358SJed Brown piece.nconn = 0; 1682e529ec8SStefano Zampini if (!localized) piece.nvertices = vEnd - vStart; 169552f7358SJed Brown for (c = cStart; c < cEnd; ++c) { 17019003fb5SStefano Zampini PetscInt dof = 0; 171552f7358SJed Brown if (hasLabel) { 172552f7358SJed Brown PetscInt value; 173552f7358SJed Brown 1749566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dm, "vtk", c, &value)); 175552f7358SJed Brown if (value != 1) continue; 176552f7358SJed Brown } 1776858538eSMatthew G. Knepley if (localized) PetscCall(PetscSectionGetDof(cellCoordSection, c, &dof)); 17819003fb5SStefano Zampini if (!dof) { 1792e529ec8SStefano Zampini PetscInt *closure = NULL; 1802e529ec8SStefano Zampini PetscInt closureSize; 1812e529ec8SStefano Zampini 1829566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 183552f7358SJed Brown for (v = 0; v < closureSize * 2; v += 2) { 1842e529ec8SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 1852e529ec8SStefano Zampini piece.nconn++; 18619003fb5SStefano Zampini if (localized) piece.nvertices++; 1872e529ec8SStefano Zampini } 188552f7358SJed Brown } 1899566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 1902e529ec8SStefano Zampini } else { 1912e529ec8SStefano Zampini piece.nvertices += dof / dimEmbed; 1922e529ec8SStefano Zampini piece.nconn += dof / dimEmbed; 1932e529ec8SStefano Zampini } 194552f7358SJed Brown piece.ncells++; 195552f7358SJed Brown } 1969566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc1(size, &gpiece)); 1979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gather((PetscInt *)&piece, sizeof(piece) / sizeof(PetscInt), MPIU_INT, (PetscInt *)gpiece, sizeof(piece) / sizeof(PetscInt), MPIU_INT, 0, comm)); 198552f7358SJed Brown 199552f7358SJed Brown /* 200552f7358SJed Brown * Write file header 201552f7358SJed Brown */ 202dd400576SPatrick Sanan if (rank == 0) { 203552f7358SJed Brown PetscInt boffset = 0; 204552f7358SJed Brown 205552f7358SJed Brown for (r = 0; r < size; r++) { 20663a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <Piece NumberOfPoints=\"%" PetscInt_FMT "\" NumberOfCells=\"%" PetscInt_FMT "\">\n", gpiece[r].nvertices, gpiece[r].ncells)); 207552f7358SJed Brown /* Coordinate positions */ 2089566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <Points>\n")); 20963a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, boffset)); 210955d60d1SToby Isaac boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + sizeof(int); 2119566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </Points>\n")); 212552f7358SJed Brown /* Cell connectivity */ 2139566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <Cells>\n")); 21463a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", boffset)); 215552f7358SJed Brown boffset += gpiece[r].nconn * sizeof(PetscInt) + sizeof(int); 21663a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"Int32\" Name=\"offsets\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", boffset)); 217552f7358SJed Brown boffset += gpiece[r].ncells * sizeof(PetscInt) + sizeof(int); 21863a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"UInt8\" Name=\"types\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", boffset)); 219552f7358SJed Brown boffset += gpiece[r].ncells * sizeof(unsigned char) + sizeof(int); 2209566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </Cells>\n")); 221552f7358SJed Brown 222552f7358SJed Brown /* 223552f7358SJed Brown * Cell Data headers 224552f7358SJed Brown */ 2259566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <CellData>\n")); 22663a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <DataArray type=\"Int32\" Name=\"Rank\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", boffset)); 227552f7358SJed Brown boffset += gpiece[r].ncells * sizeof(int) + sizeof(int); 228552f7358SJed Brown /* all the vectors */ 229552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 230552f7358SJed Brown Vec X = (Vec)link->vec; 231e630c359SToby Isaac DM dmX = NULL; 2326030a318SStefano Zampini PetscInt bs = 1, nfields, field; 233552f7358SJed Brown const char *vecname = ""; 234e630c359SToby Isaac PetscSection section; 235552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 236552f7358SJed 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. */ 2379566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)X, &vecname)); 238552f7358SJed Brown } 2399566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 240e630c359SToby Isaac if (!dmX) dmX = dm; 2419566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 2429566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 2439566063dSJacob Faibussowitsch if (cEnd > cStart) PetscCall(PetscSectionGetDof(section, cStart, &bs)); 2449566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 245e630c359SToby Isaac field = 0; 246e630c359SToby Isaac if (link->field >= 0) { 247e630c359SToby Isaac field = link->field; 248e630c359SToby Isaac nfields = field + 1; 249e630c359SToby Isaac } 250e630c359SToby Isaac for (i = 0; field < (nfields ? nfields : 1); field++) { 2511cfafdd3SJed Brown PetscInt fbs, j; 252a00cdb45SToby Isaac PetscFV fv = NULL; 253a00cdb45SToby Isaac PetscObject f; 254a00cdb45SToby Isaac PetscClassId fClass; 2550298fd71SBarry Smith const char *fieldname = NULL; 2561cfafdd3SJed Brown char buf[256]; 257e630c359SToby Isaac PetscBool vector; 2587ded4cbaSJed Brown if (nfields) { /* We have user-defined fields/components */ 2599566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, cStart, field, &fbs)); 2609566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(section, field, &fieldname)); 2617ded4cbaSJed Brown } else fbs = bs; /* Say we have one field with 'bs' components */ 2629566063dSJacob Faibussowitsch PetscCall(DMGetField(dmX, field, NULL, &f)); 2639566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(f, &fClass)); 2649371c9d4SSatish Balay if (fClass == PETSCFV_CLASSID) { fv = (PetscFV)f; } 265e630c359SToby Isaac if (nfields && !fieldname) { 26663a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "CellField%" PetscInt_FMT, field)); 267552f7358SJed Brown fieldname = buf; 268552f7358SJed Brown } 269e630c359SToby Isaac vector = PETSC_FALSE; 270e630c359SToby Isaac if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) { 271e630c359SToby Isaac vector = PETSC_TRUE; 27263a3b9bcSJacob Faibussowitsch PetscCheck(fbs <= 3, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_SIZ, "Cell vector fields can have at most 3 components, %" PetscInt_FMT " given", fbs); 273e630c359SToby Isaac for (j = 0; j < fbs; j++) { 274e630c359SToby Isaac const char *compName = NULL; 275e630c359SToby Isaac if (fv) { 2769566063dSJacob Faibussowitsch PetscCall(PetscFVGetComponentName(fv, j, &compName)); 277e630c359SToby Isaac if (compName) break; 278e630c359SToby Isaac } 279e630c359SToby Isaac } 280e630c359SToby Isaac if (j < fbs) vector = PETSC_FALSE; 281e630c359SToby Isaac } 282e630c359SToby Isaac if (vector) { 283955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 28463a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s.Re\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, vecname, fieldname, boffset)); 285955d60d1SToby Isaac boffset += gpiece[r].ncells * 3 * sizeof(PetscVTUReal) + sizeof(int); 28663a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s.Im\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, vecname, fieldname, boffset)); 287955d60d1SToby Isaac boffset += gpiece[r].ncells * 3 * sizeof(PetscVTUReal) + sizeof(int); 288955d60d1SToby Isaac #else 28963a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, vecname, fieldname, boffset)); 290955d60d1SToby Isaac boffset += gpiece[r].ncells * 3 * sizeof(PetscVTUReal) + sizeof(int); 291955d60d1SToby Isaac #endif 292e630c359SToby Isaac i += fbs; 293e630c359SToby Isaac } else { 2941cfafdd3SJed Brown for (j = 0; j < fbs; j++) { 295a00cdb45SToby Isaac const char *compName = NULL; 296955d60d1SToby Isaac char finalname[256]; 297*48a46eb9SPierre Jolivet if (fv) PetscCall(PetscFVGetComponentName(fv, j, &compName)); 298a00cdb45SToby Isaac if (compName) { 2999566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname, 255, "%s%s.%s", vecname, fieldname, compName)); 3009371c9d4SSatish Balay } else if (fbs > 1) { 30163a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname, 255, "%s%s.%" PetscInt_FMT, vecname, fieldname, j)); 302e630c359SToby Isaac } else { 3039566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname, 255, "%s%s", vecname, fieldname)); 304a00cdb45SToby Isaac } 305955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 30663a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, finalname, boffset)); 307955d60d1SToby Isaac boffset += gpiece[r].ncells * sizeof(PetscVTUReal) + sizeof(int); 30863a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, finalname, boffset)); 309955d60d1SToby Isaac boffset += gpiece[r].ncells * sizeof(PetscVTUReal) + sizeof(int); 310955d60d1SToby Isaac #else 31163a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, finalname, boffset)); 312955d60d1SToby Isaac boffset += gpiece[r].ncells * sizeof(PetscVTUReal) + sizeof(int); 313955d60d1SToby Isaac #endif 3141cfafdd3SJed Brown i++; 315552f7358SJed Brown } 316552f7358SJed Brown } 317e630c359SToby Isaac } 31863a3b9bcSJacob Faibussowitsch //PetscCheck(i == bs,comm,PETSC_ERR_PLIB,"Total number of field components %" PetscInt_FMT " != block size %" PetscInt_FMT,i,bs); 3191cfafdd3SJed Brown } 3209566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </CellData>\n")); 321552f7358SJed Brown 322552f7358SJed Brown /* 323552f7358SJed Brown * Point Data headers 324552f7358SJed Brown */ 3259566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " <PointData>\n")); 326552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 327552f7358SJed Brown Vec X = (Vec)link->vec; 328e630c359SToby Isaac DM dmX; 3296030a318SStefano Zampini PetscInt bs = 1, nfields, field; 330552f7358SJed Brown const char *vecname = ""; 331e630c359SToby Isaac PetscSection section; 332552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 333552f7358SJed 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. */ 3349566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)X, &vecname)); 335552f7358SJed Brown } 3369566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 337e630c359SToby Isaac if (!dmX) dmX = dm; 3389566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 3399566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 3409566063dSJacob Faibussowitsch if (vEnd > vStart) PetscCall(PetscSectionGetDof(section, vStart, &bs)); 3419566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 342e630c359SToby Isaac field = 0; 343e630c359SToby Isaac if (link->field >= 0) { 344e630c359SToby Isaac field = link->field; 345e630c359SToby Isaac nfields = field + 1; 346e630c359SToby Isaac } 347e630c359SToby Isaac for (i = 0; field < (nfields ? nfields : 1); field++) { 3481cfafdd3SJed Brown PetscInt fbs, j; 3491cfafdd3SJed Brown const char *fieldname = NULL; 3501cfafdd3SJed Brown char buf[256]; 3517ded4cbaSJed Brown if (nfields) { /* We have user-defined fields/components */ 3529566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, vStart, field, &fbs)); 3539566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(section, field, &fieldname)); 3547ded4cbaSJed Brown } else fbs = bs; /* Say we have one field with 'bs' components */ 355e630c359SToby Isaac if (nfields && !fieldname) { 35663a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "PointField%" PetscInt_FMT, field)); 3571cfafdd3SJed Brown fieldname = buf; 3581cfafdd3SJed Brown } 359e630c359SToby Isaac if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) { 36063a3b9bcSJacob Faibussowitsch PetscCheck(fbs <= 3, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_SIZ, "Point vector fields can have at most 3 components, %" PetscInt_FMT " given", fbs); 361955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 36263a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s.Re\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, vecname, fieldname, boffset)); 363955d60d1SToby Isaac boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + sizeof(int); 36463a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s.Im\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, vecname, fieldname, boffset)); 365955d60d1SToby Isaac boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + sizeof(int); 366955d60d1SToby Isaac #else 36763a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, vecname, fieldname, boffset)); 368955d60d1SToby Isaac boffset += gpiece[r].nvertices * 3 * sizeof(PetscVTUReal) + sizeof(int); 369955d60d1SToby Isaac #endif 370e630c359SToby Isaac } else { 3711cfafdd3SJed Brown for (j = 0; j < fbs; j++) { 372b778fa18SValeria Barra const char *compName = NULL; 373955d60d1SToby Isaac char finalname[256]; 3749566063dSJacob Faibussowitsch PetscCall(PetscSectionGetComponentName(section, field, j, &compName)); 3759566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(finalname, 255, "%s%s.%s", vecname, fieldname, compName)); 376955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 37763a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.Re\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, finalname, boffset)); 378955d60d1SToby Isaac boffset += gpiece[r].nvertices * sizeof(PetscVTUReal) + sizeof(int); 37963a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.Im\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, finalname, boffset)); 380955d60d1SToby Isaac boffset += gpiece[r].nvertices * sizeof(PetscVTUReal) + sizeof(int); 381955d60d1SToby Isaac #else 38263a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, finalname, boffset)); 383955d60d1SToby Isaac boffset += gpiece[r].nvertices * sizeof(PetscVTUReal) + sizeof(int); 384955d60d1SToby Isaac #endif 385552f7358SJed Brown } 386552f7358SJed Brown } 3871cfafdd3SJed Brown } 388e630c359SToby Isaac } 3899566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </PointData>\n")); 3909566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF, fp, " </Piece>\n")); 391552f7358SJed Brown } 392552f7358SJed Brown } 393552f7358SJed Brown 3949566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </UnstructuredGrid>\n")); 3959566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <AppendedData encoding=\"raw\">\n")); 3969566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "_")); 397552f7358SJed Brown 398dd400576SPatrick Sanan if (rank == 0) { 399552f7358SJed Brown PetscInt maxsize = 0; 400552f7358SJed Brown for (r = 0; r < size; r++) { 401955d60d1SToby Isaac maxsize = PetscMax(maxsize, (PetscInt)(gpiece[r].nvertices * 3 * sizeof(PetscVTUReal))); 402955d60d1SToby Isaac maxsize = PetscMax(maxsize, (PetscInt)(gpiece[r].ncells * 3 * sizeof(PetscVTUReal))); 403552f7358SJed Brown maxsize = PetscMax(maxsize, (PetscInt)(gpiece[r].nconn * sizeof(PetscVTKInt))); 404552f7358SJed Brown } 4059566063dSJacob Faibussowitsch PetscCall(PetscMalloc(maxsize, &buffer)); 406552f7358SJed Brown } 407552f7358SJed Brown for (r = 0; r < size; r++) { 408552f7358SJed Brown if (r == rank) { 409552f7358SJed Brown PetscInt nsend; 410552f7358SJed Brown { /* Position */ 4116858538eSMatthew G. Knepley const PetscScalar *x, *cx = NULL; 412955d60d1SToby Isaac PetscVTUReal *y = NULL; 4136858538eSMatthew G. Knepley Vec coords, cellCoords; 414955d60d1SToby Isaac PetscBool copy; 4152e529ec8SStefano Zampini 4169566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coords)); 4179566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coords, &x)); 4186858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinatesLocal(dm, &cellCoords)); 4196858538eSMatthew G. Knepley if (cellCoords) PetscCall(VecGetArrayRead(cellCoords, &cx)); 420955d60d1SToby Isaac #if defined(PETSC_USE_COMPLEX) 421955d60d1SToby Isaac copy = PETSC_TRUE; 422955d60d1SToby Isaac #else 423955d60d1SToby Isaac copy = (PetscBool)(dimEmbed != 3 || localized || (sizeof(PetscReal) != sizeof(PetscVTUReal))); 424955d60d1SToby Isaac #endif 425955d60d1SToby Isaac if (copy) { 4269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.nvertices * 3, &y)); 42719003fb5SStefano Zampini if (localized) { 42819003fb5SStefano Zampini PetscInt cnt; 42919003fb5SStefano Zampini for (c = cStart, cnt = 0; c < cEnd; c++) { 43019003fb5SStefano Zampini PetscInt off, dof; 43119003fb5SStefano Zampini 4326858538eSMatthew G. Knepley PetscCall(PetscSectionGetDof(cellCoordSection, c, &dof)); 43319003fb5SStefano Zampini if (!dof) { 43419003fb5SStefano Zampini PetscInt *closure = NULL; 43519003fb5SStefano Zampini PetscInt closureSize; 43619003fb5SStefano Zampini 4379566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 43819003fb5SStefano Zampini for (v = 0; v < closureSize * 2; v += 2) { 43919003fb5SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 4409566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, closure[v], &off)); 44119003fb5SStefano Zampini if (dimEmbed != 3) { 442955d60d1SToby Isaac y[cnt * 3 + 0] = (PetscVTUReal)PetscRealPart(x[off + 0]); 443955d60d1SToby Isaac y[cnt * 3 + 1] = (PetscVTUReal)((dimEmbed > 1) ? PetscRealPart(x[off + 1]) : 0.0); 444955d60d1SToby Isaac y[cnt * 3 + 2] = (PetscVTUReal)0.0; 44519003fb5SStefano Zampini } else { 446955d60d1SToby Isaac y[cnt * 3 + 0] = (PetscVTUReal)PetscRealPart(x[off + 0]); 447955d60d1SToby Isaac y[cnt * 3 + 1] = (PetscVTUReal)PetscRealPart(x[off + 1]); 448955d60d1SToby Isaac y[cnt * 3 + 2] = (PetscVTUReal)PetscRealPart(x[off + 2]); 44919003fb5SStefano Zampini } 45019003fb5SStefano Zampini cnt++; 45119003fb5SStefano Zampini } 45219003fb5SStefano Zampini } 4539566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 45419003fb5SStefano Zampini } else { 4556858538eSMatthew G. Knepley PetscCall(PetscSectionGetOffset(cellCoordSection, c, &off)); 45619003fb5SStefano Zampini if (dimEmbed != 3) { 45719003fb5SStefano Zampini for (i = 0; i < dof / dimEmbed; i++) { 4586858538eSMatthew G. Knepley y[cnt * 3 + 0] = (PetscVTUReal)PetscRealPart(cx[off + i * dimEmbed + 0]); 4596858538eSMatthew G. Knepley y[cnt * 3 + 1] = (PetscVTUReal)((dimEmbed > 1) ? PetscRealPart(cx[off + i * dimEmbed + 1]) : 0.0); 460955d60d1SToby Isaac y[cnt * 3 + 2] = (PetscVTUReal)0.0; 46119003fb5SStefano Zampini cnt++; 46219003fb5SStefano Zampini } 46319003fb5SStefano Zampini } else { 4649371c9d4SSatish Balay for (i = 0; i < dof; i++) { y[cnt * 3 + i] = (PetscVTUReal)PetscRealPart(cx[off + i]); } 46519003fb5SStefano Zampini cnt += dof / dimEmbed; 46619003fb5SStefano Zampini } 46719003fb5SStefano Zampini } 46819003fb5SStefano Zampini } 46908401ef6SPierre Jolivet PetscCheck(cnt == piece.nvertices, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 47019003fb5SStefano Zampini } else { 471552f7358SJed Brown for (i = 0; i < piece.nvertices; i++) { 472955d60d1SToby Isaac y[i * 3 + 0] = (PetscVTUReal)PetscRealPart(x[i * dimEmbed + 0]); 473955d60d1SToby Isaac y[i * 3 + 1] = (PetscVTUReal)((dimEmbed > 1) ? PetscRealPart(x[i * dimEmbed + 1]) : 0.); 474955d60d1SToby Isaac y[i * 3 + 2] = (PetscVTUReal)((dimEmbed > 2) ? PetscRealPart(x[i * dimEmbed + 2]) : 0.); 47519003fb5SStefano Zampini } 476552f7358SJed Brown } 477552f7358SJed Brown } 4782e529ec8SStefano Zampini nsend = piece.nvertices * 3; 4799566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, copy ? (const void *)y : (const void *)x, buffer, nsend, MPIU_VTUREAL, tag)); 4809566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 4819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coords, &x)); 4826858538eSMatthew G. Knepley if (cellCoords) PetscCall(VecRestoreArrayRead(cellCoords, &cx)); 483552f7358SJed Brown } 484552f7358SJed Brown { /* Connectivity, offsets, types */ 4853e869605SMatthew G. Knepley PetscVTKInt *connectivity = NULL, *offsets = NULL; 4863e869605SMatthew G. Knepley PetscVTKType *types = NULL; 4879566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKConnectivity(dm, localized, &piece, &connectivity, &offsets, &types)); 4889566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, connectivity, buffer, piece.nconn, MPI_INT, tag)); 4899566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, offsets, buffer, piece.ncells, MPI_INT, tag)); 4909566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, types, buffer, piece.ncells, MPI_CHAR, tag)); 4919566063dSJacob Faibussowitsch PetscCall(PetscFree3(connectivity, offsets, types)); 492552f7358SJed Brown } 493552f7358SJed Brown { /* Owners (cell data) */ 494552f7358SJed Brown PetscVTKInt *owners; 4959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.ncells, &owners)); 496552f7358SJed Brown for (i = 0; i < piece.ncells; i++) owners[i] = rank; 4979566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, owners, buffer, piece.ncells, MPI_INT, tag)); 4989566063dSJacob Faibussowitsch PetscCall(PetscFree(owners)); 499552f7358SJed Brown } 500552f7358SJed Brown /* Cell data */ 501552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 502552f7358SJed Brown Vec X = (Vec)link->vec; 503e630c359SToby Isaac DM dmX; 504552f7358SJed Brown const PetscScalar *x; 505955d60d1SToby Isaac PetscVTUReal *y; 5066030a318SStefano Zampini PetscInt bs = 1, nfields, field; 507e630c359SToby Isaac PetscSection section = NULL; 508e630c359SToby Isaac 509552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 5109566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 511e630c359SToby Isaac if (!dmX) dmX = dm; 5129566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 5139566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 5149566063dSJacob Faibussowitsch if (cEnd > cStart) PetscCall(PetscSectionGetDof(section, cStart, &bs)); 5159566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 516e630c359SToby Isaac field = 0; 517e630c359SToby Isaac if (link->field >= 0) { 518e630c359SToby Isaac field = link->field; 519e630c359SToby Isaac nfields = field + 1; 520e630c359SToby Isaac } 5219566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 5229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.ncells * 3, &y)); 523e630c359SToby Isaac for (i = 0; field < (nfields ? nfields : 1); field++) { 524e630c359SToby Isaac PetscInt fbs, j; 525e630c359SToby Isaac PetscFV fv = NULL; 526e630c359SToby Isaac PetscObject f; 527e630c359SToby Isaac PetscClassId fClass; 528e630c359SToby Isaac PetscBool vector; 5296030a318SStefano Zampini if (nfields && cEnd > cStart) { /* We have user-defined fields/components */ 5309566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, cStart, field, &fbs)); 531e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 5329566063dSJacob Faibussowitsch PetscCall(DMGetField(dmX, field, NULL, &f)); 5339566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(f, &fClass)); 5349371c9d4SSatish Balay if (fClass == PETSCFV_CLASSID) { fv = (PetscFV)f; } 535e630c359SToby Isaac vector = PETSC_FALSE; 536e630c359SToby Isaac if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) { 537e630c359SToby Isaac vector = PETSC_TRUE; 538e630c359SToby Isaac for (j = 0; j < fbs; j++) { 539e630c359SToby Isaac const char *compName = NULL; 540e630c359SToby Isaac if (fv) { 5419566063dSJacob Faibussowitsch PetscCall(PetscFVGetComponentName(fv, j, &compName)); 542e630c359SToby Isaac if (compName) break; 543e630c359SToby Isaac } 544e630c359SToby Isaac } 545e630c359SToby Isaac if (j < fbs) vector = PETSC_FALSE; 546e630c359SToby Isaac } 547e630c359SToby Isaac if (vector) { 548955d60d1SToby Isaac PetscInt cnt, l; 549955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 550552f7358SJed Brown for (c = cStart, cnt = 0; c < cEnd; c++) { 551e630c359SToby Isaac const PetscScalar *xpoint; 552e630c359SToby Isaac PetscInt off, j; 553e630c359SToby Isaac 554552f7358SJed Brown if (hasLabel) { /* Ignore some cells */ 555552f7358SJed Brown PetscInt value; 5569566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dmX, "vtk", c, &value)); 557552f7358SJed Brown if (value != 1) continue; 558552f7358SJed Brown } 559e630c359SToby Isaac if (nfields) { 5609566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, c, field, &off)); 561e630c359SToby Isaac } else { 5629566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, c, &off)); 563e630c359SToby Isaac } 564e630c359SToby Isaac xpoint = &x[off]; 5659371c9d4SSatish Balay for (j = 0; j < fbs; j++) { y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j])); } 566e630c359SToby Isaac for (; j < 3; j++) y[cnt++] = 0.; 567e630c359SToby Isaac } 56808401ef6SPierre Jolivet PetscCheck(cnt == piece.ncells * 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 5699566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.ncells * 3, MPIU_VTUREAL, tag)); 570955d60d1SToby Isaac } 571e630c359SToby Isaac } else { 572e630c359SToby Isaac for (i = 0; i < fbs; i++) { 573955d60d1SToby Isaac PetscInt cnt, l; 574955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 575e630c359SToby Isaac for (c = cStart, cnt = 0; c < cEnd; c++) { 576e630c359SToby Isaac const PetscScalar *xpoint; 577e630c359SToby Isaac PetscInt off; 578e630c359SToby Isaac 579e630c359SToby Isaac if (hasLabel) { /* Ignore some cells */ 580e630c359SToby Isaac PetscInt value; 5819566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dmX, "vtk", c, &value)); 582e630c359SToby Isaac if (value != 1) continue; 583e630c359SToby Isaac } 584e630c359SToby Isaac if (nfields) { 5859566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, c, field, &off)); 586e630c359SToby Isaac } else { 5879566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, c, &off)); 588e630c359SToby Isaac } 589e630c359SToby Isaac xpoint = &x[off]; 590955d60d1SToby Isaac y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i])); 591552f7358SJed Brown } 59208401ef6SPierre Jolivet PetscCheck(cnt == piece.ncells, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 5939566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.ncells, MPIU_VTUREAL, tag)); 594955d60d1SToby Isaac } 595552f7358SJed Brown } 596e630c359SToby Isaac } 597e630c359SToby Isaac } 5989566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 5999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 600552f7358SJed Brown } 601e630c359SToby Isaac /* point data */ 602552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 603552f7358SJed Brown Vec X = (Vec)link->vec; 604e630c359SToby Isaac DM dmX; 605552f7358SJed Brown const PetscScalar *x; 606955d60d1SToby Isaac PetscVTUReal *y; 6076030a318SStefano Zampini PetscInt bs = 1, nfields, field; 608e630c359SToby Isaac PetscSection section = NULL; 609e630c359SToby Isaac 610552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 6119566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 612e630c359SToby Isaac if (!dmX) dmX = dm; 6139566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 6149566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 6159566063dSJacob Faibussowitsch if (vEnd > vStart) PetscCall(PetscSectionGetDof(section, vStart, &bs)); 6169566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 617e630c359SToby Isaac field = 0; 618e630c359SToby Isaac if (link->field >= 0) { 619e630c359SToby Isaac field = link->field; 620e630c359SToby Isaac nfields = field + 1; 621e630c359SToby Isaac } 6229566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 6239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(piece.nvertices * 3, &y)); 624e630c359SToby Isaac for (i = 0; field < (nfields ? nfields : 1); field++) { 625e630c359SToby Isaac PetscInt fbs, j; 6266030a318SStefano Zampini if (nfields && vEnd > vStart) { /* We have user-defined fields/components */ 6279566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, vStart, field, &fbs)); 628e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 629e630c359SToby Isaac if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) { 630955d60d1SToby Isaac PetscInt cnt, l; 631955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 6322e529ec8SStefano Zampini if (!localized) { 633552f7358SJed Brown for (v = vStart, cnt = 0; v < vEnd; v++) { 634e630c359SToby Isaac PetscInt off; 635e630c359SToby Isaac const PetscScalar *xpoint; 636e630c359SToby Isaac 637e630c359SToby Isaac if (nfields) { 6389566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, v, field, &off)); 639e630c359SToby Isaac } else { 6409566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, v, &off)); 641e630c359SToby Isaac } 642e630c359SToby Isaac xpoint = &x[off]; 6439371c9d4SSatish Balay for (j = 0; j < fbs; j++) { y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j])); } 644e630c359SToby Isaac for (; j < 3; j++) y[cnt++] = 0.; 645e630c359SToby Isaac } 646e630c359SToby Isaac } else { 647e630c359SToby Isaac for (c = cStart, cnt = 0; c < cEnd; c++) { 648e630c359SToby Isaac PetscInt *closure = NULL; 649e630c359SToby Isaac PetscInt closureSize, off; 650e630c359SToby Isaac 6519566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 652e630c359SToby Isaac for (v = 0, off = 0; v < closureSize * 2; v += 2) { 653e630c359SToby Isaac if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 654e630c359SToby Isaac PetscInt voff; 655e630c359SToby Isaac const PetscScalar *xpoint; 656e630c359SToby Isaac 657e630c359SToby Isaac if (nfields) { 6589566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, closure[v], field, &voff)); 659e630c359SToby Isaac } else { 6609566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, closure[v], &voff)); 661e630c359SToby Isaac } 662e630c359SToby Isaac xpoint = &x[voff]; 6639371c9d4SSatish Balay for (j = 0; j < fbs; j++) { y[cnt + off++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[j]) : PetscRealPart(xpoint[j])); } 664e630c359SToby Isaac for (; j < 3; j++) y[cnt + off++] = 0.; 665e630c359SToby Isaac } 666e630c359SToby Isaac } 667e630c359SToby Isaac cnt += off; 6689566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 669e630c359SToby Isaac } 670e630c359SToby Isaac } 67108401ef6SPierre Jolivet PetscCheck(cnt == piece.nvertices * 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 6729566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.nvertices * 3, MPIU_VTUREAL, tag)); 673955d60d1SToby Isaac } 674e630c359SToby Isaac } else { 675e630c359SToby Isaac for (i = 0; i < fbs; i++) { 676955d60d1SToby Isaac PetscInt cnt, l; 677955d60d1SToby Isaac for (l = 0; l < loops_per_scalar; l++) { 678e630c359SToby Isaac if (!localized) { 679e630c359SToby Isaac for (v = vStart, cnt = 0; v < vEnd; v++) { 680e630c359SToby Isaac PetscInt off; 681e630c359SToby Isaac const PetscScalar *xpoint; 682e630c359SToby Isaac 683e630c359SToby Isaac if (nfields) { 6849566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, v, field, &off)); 685e630c359SToby Isaac } else { 6869566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, v, &off)); 687e630c359SToby Isaac } 688e630c359SToby Isaac xpoint = &x[off]; 689955d60d1SToby Isaac y[cnt++] = (PetscVTUReal)(l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i])); 690552f7358SJed Brown } 6912e529ec8SStefano Zampini } else { 6922e529ec8SStefano Zampini for (c = cStart, cnt = 0; c < cEnd; c++) { 6932e529ec8SStefano Zampini PetscInt *closure = NULL; 6942e529ec8SStefano Zampini PetscInt closureSize, off; 6952e529ec8SStefano Zampini 6969566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 6972e529ec8SStefano Zampini for (v = 0, off = 0; v < closureSize * 2; v += 2) { 6982e529ec8SStefano Zampini if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 699e630c359SToby Isaac PetscInt voff; 700e630c359SToby Isaac const PetscScalar *xpoint; 7012e529ec8SStefano Zampini 702e630c359SToby Isaac if (nfields) { 7039566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(section, closure[v], field, &voff)); 704e630c359SToby Isaac } else { 7059566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, closure[v], &voff)); 706e630c359SToby Isaac } 707e630c359SToby Isaac xpoint = &x[voff]; 708955d60d1SToby Isaac y[cnt + off++] = (l ? PetscImaginaryPart(xpoint[i]) : PetscRealPart(xpoint[i])); 7092e529ec8SStefano Zampini } 7102e529ec8SStefano Zampini } 7119bda96f6SStefano Zampini cnt += off; 7129566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dmX, c, PETSC_TRUE, &closureSize, &closure)); 7132e529ec8SStefano Zampini } 7142e529ec8SStefano Zampini } 71508401ef6SPierre Jolivet PetscCheck(cnt == piece.nvertices, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count does not match"); 7169566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, y, buffer, piece.nvertices, MPIU_VTUREAL, tag)); 717955d60d1SToby Isaac } 718552f7358SJed Brown } 719e630c359SToby Isaac } 720e630c359SToby Isaac } 7219566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 7229566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 723552f7358SJed Brown } 724dd400576SPatrick Sanan } else if (rank == 0) { 725955d60d1SToby Isaac PetscInt l; 726955d60d1SToby Isaac 7279566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nvertices * 3, MPIU_VTUREAL, tag)); /* positions */ 7289566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nconn, MPI_INT, tag)); /* connectivity */ 7299566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPI_INT, tag)); /* offsets */ 7309566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPI_CHAR, tag)); /* types */ 7319566063dSJacob Faibussowitsch PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPI_INT, tag)); /* owner rank (cells) */ 732552f7358SJed Brown /* all cell data */ 733552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 734e630c359SToby Isaac Vec X = (Vec)link->vec; 7356030a318SStefano Zampini PetscInt bs = 1, nfields, field; 736e630c359SToby Isaac DM dmX; 737e630c359SToby Isaac PetscSection section = NULL; 738e630c359SToby Isaac 739552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 7409566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 741e630c359SToby Isaac if (!dmX) dmX = dm; 7429566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 7439566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 7449566063dSJacob Faibussowitsch if (cEnd > cStart) PetscCall(PetscSectionGetDof(section, cStart, &bs)); 7459566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 746e630c359SToby Isaac field = 0; 747e630c359SToby Isaac if (link->field >= 0) { 748e630c359SToby Isaac field = link->field; 749e630c359SToby Isaac nfields = field + 1; 750e630c359SToby Isaac } 751e630c359SToby Isaac for (i = 0; field < (nfields ? nfields : 1); field++) { 752e630c359SToby Isaac PetscInt fbs, j; 753e630c359SToby Isaac PetscFV fv = NULL; 754e630c359SToby Isaac PetscObject f; 755e630c359SToby Isaac PetscClassId fClass; 756e630c359SToby Isaac PetscBool vector; 7576030a318SStefano Zampini if (nfields && cEnd > cStart) { /* We have user-defined fields/components */ 7589566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, cStart, field, &fbs)); 759e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 7609566063dSJacob Faibussowitsch PetscCall(DMGetField(dmX, field, NULL, &f)); 7619566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(f, &fClass)); 7629371c9d4SSatish Balay if (fClass == PETSCFV_CLASSID) { fv = (PetscFV)f; } 763e630c359SToby Isaac vector = PETSC_FALSE; 764e630c359SToby Isaac if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) { 765e630c359SToby Isaac vector = PETSC_TRUE; 766e630c359SToby Isaac for (j = 0; j < fbs; j++) { 767e630c359SToby Isaac const char *compName = NULL; 768e630c359SToby Isaac if (fv) { 7699566063dSJacob Faibussowitsch PetscCall(PetscFVGetComponentName(fv, j, &compName)); 770e630c359SToby Isaac if (compName) break; 771e630c359SToby Isaac } 772e630c359SToby Isaac } 773e630c359SToby Isaac if (j < fbs) vector = PETSC_FALSE; 774e630c359SToby Isaac } 775e630c359SToby Isaac if (vector) { 776*48a46eb9SPierre 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)); 777e630c359SToby Isaac } else { 778e630c359SToby Isaac for (i = 0; i < fbs; i++) { 779*48a46eb9SPierre Jolivet for (l = 0; l < loops_per_scalar; l++) PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].ncells, MPIU_VTUREAL, tag)); 780552f7358SJed Brown } 781552f7358SJed Brown } 782e630c359SToby Isaac } 783e630c359SToby Isaac } 784552f7358SJed Brown /* all point data */ 785552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 786e630c359SToby Isaac Vec X = (Vec)link->vec; 787e630c359SToby Isaac DM dmX; 7886030a318SStefano Zampini PetscInt bs = 1, nfields, field; 789e630c359SToby Isaac PetscSection section = NULL; 790e630c359SToby Isaac 791552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 7929566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 793e630c359SToby Isaac if (!dmX) dmX = dm; 7949566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject *)§ion)); 7959566063dSJacob Faibussowitsch if (!section) PetscCall(DMGetLocalSection(dmX, §ion)); 7969566063dSJacob Faibussowitsch if (vEnd > vStart) PetscCall(PetscSectionGetDof(section, vStart, &bs)); 7979566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &nfields)); 798e630c359SToby Isaac field = 0; 799e630c359SToby Isaac if (link->field >= 0) { 800e630c359SToby Isaac field = link->field; 801e630c359SToby Isaac nfields = field + 1; 802e630c359SToby Isaac } 803e630c359SToby Isaac for (i = 0; field < (nfields ? nfields : 1); field++) { 804e630c359SToby Isaac PetscInt fbs; 8056030a318SStefano Zampini if (nfields && vEnd > vStart) { /* We have user-defined fields/components */ 8069566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, vStart, field, &fbs)); 807e630c359SToby Isaac } else fbs = bs; /* Say we have one field with 'bs' components */ 808e630c359SToby Isaac if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) { 809*48a46eb9SPierre 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)); 810e630c359SToby Isaac } else { 811e630c359SToby Isaac for (i = 0; i < fbs; i++) { 812*48a46eb9SPierre Jolivet for (l = 0; l < loops_per_scalar; l++) PetscCall(TransferWrite(comm, viewer, fp, r, 0, NULL, buffer, gpiece[r].nvertices, MPIU_VTUREAL, tag)); 813552f7358SJed Brown } 814552f7358SJed Brown } 815552f7358SJed Brown } 816552f7358SJed Brown } 817e630c359SToby Isaac } 818e630c359SToby Isaac } 8199566063dSJacob Faibussowitsch PetscCall(PetscFree(gpiece)); 8209566063dSJacob Faibussowitsch PetscCall(PetscFree(buffer)); 8219566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "\n </AppendedData>\n")); 8229566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "</VTKFile>\n")); 8239566063dSJacob Faibussowitsch PetscCall(PetscFClose(comm, fp)); 824552f7358SJed Brown PetscFunctionReturn(0); 825552f7358SJed Brown } 826