1552f7358SJed Brown #define PETSCDM_DLL 2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3552f7358SJed Brown #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h> 4552f7358SJed Brown 5de0a4605SMatthew G. Knepley PetscErrorCode DMPlexVTKGetCellType_Internal(DM dm, PetscInt dim, PetscInt corners, PetscInt *cellType) 60adebc6cSBarry Smith { 7552f7358SJed Brown PetscFunctionBegin; 8552f7358SJed Brown *cellType = -1; 9552f7358SJed Brown switch (dim) { 10552f7358SJed Brown case 0: 11552f7358SJed Brown switch (corners) { 12552f7358SJed Brown case 1: 13552f7358SJed Brown *cellType = 1; /* VTK_VERTEX */ 14552f7358SJed Brown break; 15552f7358SJed Brown default: 16552f7358SJed Brown break; 17552f7358SJed Brown } 18552f7358SJed Brown break; 19552f7358SJed Brown case 1: 20552f7358SJed Brown switch (corners) { 21552f7358SJed Brown case 2: 22552f7358SJed Brown *cellType = 3; /* VTK_LINE */ 23552f7358SJed Brown break; 24552f7358SJed Brown case 3: 25552f7358SJed Brown *cellType = 21; /* VTK_QUADRATIC_EDGE */ 26552f7358SJed Brown break; 27552f7358SJed Brown default: 28552f7358SJed Brown break; 29552f7358SJed Brown } 30552f7358SJed Brown break; 31552f7358SJed Brown case 2: 32552f7358SJed Brown switch (corners) { 33552f7358SJed Brown case 3: 34552f7358SJed Brown *cellType = 5; /* VTK_TRIANGLE */ 35552f7358SJed Brown break; 36552f7358SJed Brown case 4: 37552f7358SJed Brown *cellType = 9; /* VTK_QUAD */ 38552f7358SJed Brown break; 39552f7358SJed Brown case 6: 40552f7358SJed Brown *cellType = 22; /* VTK_QUADRATIC_TRIANGLE */ 41552f7358SJed Brown break; 42552f7358SJed Brown case 9: 43552f7358SJed Brown *cellType = 23; /* VTK_QUADRATIC_QUAD */ 44552f7358SJed Brown break; 45552f7358SJed Brown default: 46552f7358SJed Brown break; 47552f7358SJed Brown } 48552f7358SJed Brown break; 49552f7358SJed Brown case 3: 50552f7358SJed Brown switch (corners) { 51552f7358SJed Brown case 4: 52552f7358SJed Brown *cellType = 10; /* VTK_TETRA */ 53552f7358SJed Brown break; 54*cf4edabeSMatthew G. Knepley case 5: 55*cf4edabeSMatthew G. Knepley *cellType = 14; /* VTK_PYRAMID */ 56*cf4edabeSMatthew G. Knepley break; 572f029394SStefano Zampini case 6: 582f029394SStefano Zampini *cellType = 13; /* VTK_WEDGE */ 592f029394SStefano Zampini break; 60552f7358SJed Brown case 8: 61552f7358SJed Brown *cellType = 12; /* VTK_HEXAHEDRON */ 62552f7358SJed Brown break; 63552f7358SJed Brown case 10: 64552f7358SJed Brown *cellType = 24; /* VTK_QUADRATIC_TETRA */ 65552f7358SJed Brown break; 66552f7358SJed Brown case 27: 67552f7358SJed Brown *cellType = 29; /* VTK_QUADRATIC_HEXAHEDRON */ 68552f7358SJed Brown break; 69552f7358SJed Brown default: 70552f7358SJed Brown break; 71552f7358SJed Brown } 72552f7358SJed Brown } 73552f7358SJed Brown PetscFunctionReturn(0); 74552f7358SJed Brown } 75552f7358SJed Brown 767508e736SMatthew G. Knepley static PetscErrorCode DMPlexVTKWriteCells_ASCII(DM dm, FILE *fp, PetscInt *totalCells) 77552f7358SJed Brown { 7882f516ccSBarry Smith MPI_Comm comm; 7902281ff3SMatthew G. Knepley DMLabel label; 8002281ff3SMatthew G. Knepley IS globalVertexNumbers = NULL; 81552f7358SJed Brown const PetscInt *gvertex; 82552f7358SJed Brown PetscInt dim; 83552f7358SJed Brown PetscInt numCorners = 0, totCorners = 0, maxCorners, *corners; 84552f7358SJed Brown PetscInt numCells = 0, totCells = 0, maxCells, cellHeight; 852f029394SStefano Zampini PetscInt numLabelCells, maxLabelCells, cStart, cEnd, c, vStart, vEnd, v; 869852e123SBarry Smith PetscMPIInt size, rank, proc, tag; 87552f7358SJed Brown 88552f7358SJed Brown PetscFunctionBegin; 899566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm,&comm)); 909566063dSJacob Faibussowitsch PetscCall(PetscCommGetNewTag(comm, &tag)); 919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 939566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 949566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight)); 959566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd)); 969566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 979566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "vtk", &label)); 989566063dSJacob Faibussowitsch PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells)); 991c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&numLabelCells, &maxLabelCells, 1, MPIU_INT, MPI_MAX, comm)); 100452b7979SMatthew G. Knepley if (!maxLabelCells) label = NULL; 101552f7358SJed Brown for (c = cStart; c < cEnd; ++c) { 1020298fd71SBarry Smith PetscInt *closure = NULL; 10302281ff3SMatthew G. Knepley PetscInt closureSize, value; 104552f7358SJed Brown 10502281ff3SMatthew G. Knepley if (label) { 1069566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, c, &value)); 107552f7358SJed Brown if (value != 1) continue; 108552f7358SJed Brown } 1099566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 110552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 111552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) ++numCorners; 112552f7358SJed Brown } 1139566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 114552f7358SJed Brown ++numCells; 115552f7358SJed Brown } 116552f7358SJed Brown maxCells = numCells; 1179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce(&numCells, &totCells, 1, MPIU_INT, MPI_SUM, 0, comm)); 1189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce(&numCells, &maxCells, 1, MPIU_INT, MPI_MAX, 0, comm)); 1199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce(&numCorners, &totCorners, 1, MPIU_INT, MPI_SUM, 0, comm)); 1209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce(&numCorners, &maxCorners, 1, MPIU_INT, MPI_MAX, 0, comm)); 1219566063dSJacob Faibussowitsch PetscCall(DMPlexGetVertexNumbering(dm, &globalVertexNumbers)); 1229566063dSJacob Faibussowitsch PetscCall(ISGetIndices(globalVertexNumbers, &gvertex)); 1239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxCells, &corners)); 12463a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "CELLS %" PetscInt_FMT " %" PetscInt_FMT "\n", totCells, totCorners+totCells)); 125dd400576SPatrick Sanan if (rank == 0) { 126a4a685f2SJacob Faibussowitsch PetscInt *remoteVertices, *vertices; 127552f7358SJed Brown 1289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxCorners, &vertices)); 129552f7358SJed Brown for (c = cStart, numCells = 0; c < cEnd; ++c) { 1300298fd71SBarry Smith PetscInt *closure = NULL; 13102281ff3SMatthew G. Knepley PetscInt closureSize, value, nC = 0; 132552f7358SJed Brown 13302281ff3SMatthew G. Knepley if (label) { 1349566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, c, &value)); 135552f7358SJed Brown if (value != 1) continue; 136552f7358SJed Brown } 1379566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 138552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 139552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 140552f7358SJed Brown const PetscInt gv = gvertex[closure[v] - vStart]; 1410dcf9c70SMatthew G. Knepley vertices[nC++] = gv < 0 ? -(gv+1) : gv; 142552f7358SJed Brown } 143552f7358SJed Brown } 1449566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 1459566063dSJacob Faibussowitsch PetscCall(DMPlexReorderCell(dm, c, vertices)); 146552f7358SJed Brown corners[numCells++] = nC; 14763a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "%" PetscInt_FMT " ", nC)); 148552f7358SJed Brown for (v = 0; v < nC; ++v) { 14963a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " %" PetscInt_FMT, vertices[v])); 150552f7358SJed Brown } 1519566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "\n")); 152552f7358SJed Brown } 1539566063dSJacob Faibussowitsch if (size > 1) PetscCall(PetscMalloc1(maxCorners+maxCells, &remoteVertices)); 1549852e123SBarry Smith for (proc = 1; proc < size; ++proc) { 155552f7358SJed Brown MPI_Status status; 156552f7358SJed Brown 1579566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(&numCorners, 1, MPIU_INT, proc, tag, comm, &status)); 1589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(remoteVertices, numCorners, MPIU_INT, proc, tag, comm, &status)); 159552f7358SJed Brown for (c = 0; c < numCorners;) { 160552f7358SJed Brown PetscInt nC = remoteVertices[c++]; 161552f7358SJed Brown 162552f7358SJed Brown for (v = 0; v < nC; ++v, ++c) { 1630dcf9c70SMatthew G. Knepley vertices[v] = remoteVertices[c]; 1640dcf9c70SMatthew G. Knepley } 16563a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "%" PetscInt_FMT " ", nC)); 1660dcf9c70SMatthew G. Knepley for (v = 0; v < nC; ++v) { 16763a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " %" PetscInt_FMT, vertices[v])); 168552f7358SJed Brown } 1699566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "\n")); 170552f7358SJed Brown } 171552f7358SJed Brown } 1729566063dSJacob Faibussowitsch if (size > 1) PetscCall(PetscFree(remoteVertices)); 1739566063dSJacob Faibussowitsch PetscCall(PetscFree(vertices)); 174552f7358SJed Brown } else { 175552f7358SJed Brown PetscInt *localVertices, numSend = numCells+numCorners, k = 0; 176552f7358SJed Brown 1779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numSend, &localVertices)); 178552f7358SJed Brown for (c = cStart, numCells = 0; c < cEnd; ++c) { 1790298fd71SBarry Smith PetscInt *closure = NULL; 18002281ff3SMatthew G. Knepley PetscInt closureSize, value, nC = 0; 181552f7358SJed Brown 18202281ff3SMatthew G. Knepley if (label) { 1839566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, c, &value)); 184552f7358SJed Brown if (value != 1) continue; 185552f7358SJed Brown } 1869566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 187552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 188552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 189552f7358SJed Brown const PetscInt gv = gvertex[closure[v] - vStart]; 190552f7358SJed Brown closure[nC++] = gv < 0 ? -(gv+1) : gv; 191552f7358SJed Brown } 192552f7358SJed Brown } 193552f7358SJed Brown corners[numCells++] = nC; 194552f7358SJed Brown localVertices[k++] = nC; 195552f7358SJed Brown for (v = 0; v < nC; ++v, ++k) { 196552f7358SJed Brown localVertices[k] = closure[v]; 197552f7358SJed Brown } 1989566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 1999566063dSJacob Faibussowitsch PetscCall(DMPlexReorderCell(dm, c, localVertices+k-nC)); 200552f7358SJed Brown } 20163a3b9bcSJacob Faibussowitsch PetscCheck(k == numSend,PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid number of vertices to send %" PetscInt_FMT " should be %" PetscInt_FMT, k, numSend); 2029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(&numSend, 1, MPIU_INT, 0, tag, comm)); 2039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(localVertices, numSend, MPIU_INT, 0, tag, comm)); 2049566063dSJacob Faibussowitsch PetscCall(PetscFree(localVertices)); 205552f7358SJed Brown } 2069566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(globalVertexNumbers, &gvertex)); 20763a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "CELL_TYPES %" PetscInt_FMT "\n", totCells)); 208dd400576SPatrick Sanan if (rank == 0) { 209552f7358SJed Brown PetscInt cellType; 210552f7358SJed Brown 211552f7358SJed Brown for (c = 0; c < numCells; ++c) { 2129566063dSJacob Faibussowitsch PetscCall(DMPlexVTKGetCellType_Internal(dm, dim, corners[c], &cellType)); 21363a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "%" PetscInt_FMT "\n", cellType)); 214552f7358SJed Brown } 2159852e123SBarry Smith for (proc = 1; proc < size; ++proc) { 216552f7358SJed Brown MPI_Status status; 217552f7358SJed Brown 2189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status)); 2199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(corners, numCells, MPIU_INT, proc, tag, comm, &status)); 220552f7358SJed Brown for (c = 0; c < numCells; ++c) { 2219566063dSJacob Faibussowitsch PetscCall(DMPlexVTKGetCellType_Internal(dm, dim, corners[c], &cellType)); 22263a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "%" PetscInt_FMT "\n", cellType)); 223552f7358SJed Brown } 224552f7358SJed Brown } 225552f7358SJed Brown } else { 2269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm)); 2279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(corners, numCells, MPIU_INT, 0, tag, comm)); 228552f7358SJed Brown } 2299566063dSJacob Faibussowitsch PetscCall(PetscFree(corners)); 230552f7358SJed Brown *totalCells = totCells; 231552f7358SJed Brown PetscFunctionReturn(0); 232552f7358SJed Brown } 233552f7358SJed Brown 2347508e736SMatthew G. Knepley static PetscErrorCode DMPlexVTKWritePartition_ASCII(DM dm, FILE *fp) 235e8f6f0f6SMatthew G. Knepley { 236e8f6f0f6SMatthew G. Knepley MPI_Comm comm; 237e8f6f0f6SMatthew G. Knepley PetscInt numCells = 0, cellHeight; 2382f029394SStefano Zampini PetscInt numLabelCells, cStart, cEnd, c; 2399852e123SBarry Smith PetscMPIInt size, rank, proc, tag; 240e8f6f0f6SMatthew G. Knepley PetscBool hasLabel; 241e8f6f0f6SMatthew G. Knepley 242e8f6f0f6SMatthew G. Knepley PetscFunctionBegin; 2439566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm,&comm)); 2449566063dSJacob Faibussowitsch PetscCall(PetscCommGetNewTag(comm, &tag)); 2459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2479566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight)); 2489566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd)); 2499566063dSJacob Faibussowitsch PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells)); 250e8f6f0f6SMatthew G. Knepley hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 251e8f6f0f6SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 252e8f6f0f6SMatthew G. Knepley if (hasLabel) { 253e8f6f0f6SMatthew G. Knepley PetscInt value; 254e8f6f0f6SMatthew G. Knepley 2559566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dm, "vtk", c, &value)); 256e8f6f0f6SMatthew G. Knepley if (value != 1) continue; 257e8f6f0f6SMatthew G. Knepley } 258e8f6f0f6SMatthew G. Knepley ++numCells; 259e8f6f0f6SMatthew G. Knepley } 260dd400576SPatrick Sanan if (rank == 0) { 2619566063dSJacob Faibussowitsch for (c = 0; c < numCells; ++c) PetscCall(PetscFPrintf(comm, fp, "%d\n", rank)); 2629852e123SBarry Smith for (proc = 1; proc < size; ++proc) { 263e8f6f0f6SMatthew G. Knepley MPI_Status status; 264e8f6f0f6SMatthew G. Knepley 2659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status)); 2669566063dSJacob Faibussowitsch for (c = 0; c < numCells; ++c) PetscCall(PetscFPrintf(comm, fp, "%d\n", proc)); 267e8f6f0f6SMatthew G. Knepley } 268e8f6f0f6SMatthew G. Knepley } else { 2699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm)); 270e8f6f0f6SMatthew G. Knepley } 271e8f6f0f6SMatthew G. Knepley PetscFunctionReturn(0); 272e8f6f0f6SMatthew G. Knepley } 273e8f6f0f6SMatthew G. Knepley 27476b700caSToby Isaac #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128) 27576b700caSToby Isaac typedef double PetscVTKReal; 27676b700caSToby Isaac #elif defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL___FP16) 27776b700caSToby Isaac typedef float PetscVTKReal; 27876b700caSToby Isaac #else 27976b700caSToby Isaac typedef PetscReal PetscVTKReal; 28076b700caSToby Isaac #endif 28176b700caSToby Isaac 28276b700caSToby Isaac static PetscErrorCode DMPlexVTKWriteSection_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec v, FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale, PetscInt imag) 2830adebc6cSBarry Smith { 28482f516ccSBarry Smith MPI_Comm comm; 285552f7358SJed Brown const MPI_Datatype mpiType = MPIU_SCALAR; 286552f7358SJed Brown PetscScalar *array; 287552f7358SJed Brown PetscInt numDof = 0, maxDof; 288412e9a14SMatthew G. Knepley PetscInt numLabelCells, cellHeight, cStart, cEnd, numLabelVertices, vStart, vEnd, pStart, pEnd, p; 2899852e123SBarry Smith PetscMPIInt size, rank, proc, tag; 290552f7358SJed Brown PetscBool hasLabel; 291552f7358SJed Brown 292552f7358SJed Brown PetscFunctionBegin; 2939566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm,&comm)); 294552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 295552f7358SJed Brown PetscValidHeaderSpecific(v,VEC_CLASSID,4); 296552f7358SJed Brown if (precision < 0) precision = 6; 2979566063dSJacob Faibussowitsch PetscCall(PetscCommGetNewTag(comm, &tag)); 2989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 3009566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 301552f7358SJed Brown /* VTK only wants the values at cells or vertices */ 3029566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight)); 3039566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd)); 3049566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 305552f7358SJed Brown pStart = PetscMax(PetscMin(cStart, vStart), pStart); 306552f7358SJed Brown pEnd = PetscMin(PetscMax(cEnd, vEnd), pEnd); 3079566063dSJacob Faibussowitsch PetscCall(DMGetStratumSize(dm, "vtk", 1, &numLabelCells)); 3089566063dSJacob Faibussowitsch PetscCall(DMGetStratumSize(dm, "vtk", 2, &numLabelVertices)); 309552f7358SJed Brown hasLabel = numLabelCells > 0 || numLabelVertices > 0 ? PETSC_TRUE : PETSC_FALSE; 310552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 311552f7358SJed Brown /* Reject points not either cells or vertices */ 312552f7358SJed Brown if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue; 313552f7358SJed Brown if (hasLabel) { 314552f7358SJed Brown PetscInt value; 315552f7358SJed Brown 316552f7358SJed Brown if (((p >= cStart) && (p < cEnd) && numLabelCells) || 317552f7358SJed Brown ((p >= vStart) && (p < vEnd) && numLabelVertices)) { 3189566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dm, "vtk", p, &value)); 319552f7358SJed Brown if (value != 1) continue; 320552f7358SJed Brown } 321552f7358SJed Brown } 3229566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &numDof)); 3238865f1eaSKarl Rupp if (numDof) break; 324552f7358SJed Brown } 3251c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, comm)); 326552f7358SJed Brown enforceDof = PetscMax(enforceDof, maxDof); 3279566063dSJacob Faibussowitsch PetscCall(VecGetArray(v, &array)); 328dd400576SPatrick Sanan if (rank == 0) { 32976b700caSToby Isaac PetscVTKReal dval; 33076b700caSToby Isaac PetscScalar val; 331552f7358SJed Brown char formatString[8]; 332552f7358SJed Brown 33363a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(formatString, 8, "%%.%" PetscInt_FMT "e", precision)); 334552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 335552f7358SJed Brown /* Here we lose a way to filter points by keeping them out of the Numbering */ 336552f7358SJed Brown PetscInt dof, off, goff, d; 337552f7358SJed Brown 338552f7358SJed Brown /* Reject points not either cells or vertices */ 339552f7358SJed Brown if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue; 340552f7358SJed Brown if (hasLabel) { 341552f7358SJed Brown PetscInt value; 342552f7358SJed Brown 343552f7358SJed Brown if (((p >= cStart) && (p < cEnd) && numLabelCells) || 344552f7358SJed Brown ((p >= vStart) && (p < vEnd) && numLabelVertices)) { 3459566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dm, "vtk", p, &value)); 346552f7358SJed Brown if (value != 1) continue; 347552f7358SJed Brown } 348552f7358SJed Brown } 3499566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof)); 3509566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off)); 3519566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(globalSection, p, &goff)); 352552f7358SJed Brown if (dof && goff >= 0) { 353552f7358SJed Brown for (d = 0; d < dof; d++) { 354552f7358SJed Brown if (d > 0) { 3559566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " ")); 356552f7358SJed Brown } 35776b700caSToby Isaac val = array[off+d]; 35876b700caSToby Isaac dval = (PetscVTKReal) ((imag ? PetscImaginaryPart(val) : PetscRealPart(val)) * scale); 3599566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, formatString, dval)); 360552f7358SJed Brown } 361552f7358SJed Brown for (d = dof; d < enforceDof; d++) { 3629566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " 0.0")); 363552f7358SJed Brown } 3649566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "\n")); 365552f7358SJed Brown } 366552f7358SJed Brown } 3679852e123SBarry Smith for (proc = 1; proc < size; ++proc) { 368552f7358SJed Brown PetscScalar *remoteValues; 369d892089bSMatthew G. Knepley PetscInt size = 0, d; 370552f7358SJed Brown MPI_Status status; 371552f7358SJed Brown 3729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(&size, 1, MPIU_INT, proc, tag, comm, &status)); 3739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &remoteValues)); 3749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(remoteValues, size, mpiType, proc, tag, comm, &status)); 375552f7358SJed Brown for (p = 0; p < size/maxDof; ++p) { 376552f7358SJed Brown for (d = 0; d < maxDof; ++d) { 377552f7358SJed Brown if (d > 0) { 3789566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " ")); 379552f7358SJed Brown } 38076b700caSToby Isaac val = remoteValues[p*maxDof+d]; 38176b700caSToby Isaac dval = (PetscVTKReal) ((imag ? PetscImaginaryPart(val) : PetscRealPart(val)) * scale); 3829566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, formatString, dval)); 383552f7358SJed Brown } 384552f7358SJed Brown for (d = maxDof; d < enforceDof; ++d) { 3859566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " 0.0")); 386552f7358SJed Brown } 3879566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "\n")); 388552f7358SJed Brown } 3899566063dSJacob Faibussowitsch PetscCall(PetscFree(remoteValues)); 390552f7358SJed Brown } 391552f7358SJed Brown } else { 392552f7358SJed Brown PetscScalar *localValues; 393552f7358SJed Brown PetscInt size, k = 0; 394552f7358SJed Brown 3959566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(section, &size)); 3969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &localValues)); 397552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 398552f7358SJed Brown PetscInt dof, off, goff, d; 399552f7358SJed Brown 400552f7358SJed Brown /* Reject points not either cells or vertices */ 401552f7358SJed Brown if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue; 402552f7358SJed Brown if (hasLabel) { 403552f7358SJed Brown PetscInt value; 404552f7358SJed Brown 405552f7358SJed Brown if (((p >= cStart) && (p < cEnd) && numLabelCells) || 406552f7358SJed Brown ((p >= vStart) && (p < vEnd) && numLabelVertices)) { 4079566063dSJacob Faibussowitsch PetscCall(DMGetLabelValue(dm, "vtk", p, &value)); 408552f7358SJed Brown if (value != 1) continue; 409552f7358SJed Brown } 410552f7358SJed Brown } 4119566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof)); 4129566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off)); 4139566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(globalSection, p, &goff)); 414552f7358SJed Brown if (goff >= 0) { 415552f7358SJed Brown for (d = 0; d < dof; ++d) { 416552f7358SJed Brown localValues[k++] = array[off+d]; 417552f7358SJed Brown } 418552f7358SJed Brown } 419552f7358SJed Brown } 4209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(&k, 1, MPIU_INT, 0, tag, comm)); 4219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(localValues, k, mpiType, 0, tag, comm)); 4229566063dSJacob Faibussowitsch PetscCall(PetscFree(localValues)); 423552f7358SJed Brown } 4249566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v, &array)); 425552f7358SJed Brown PetscFunctionReturn(0); 426552f7358SJed Brown } 427552f7358SJed Brown 42876b700caSToby Isaac static PetscErrorCode DMPlexVTKWriteField_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec field, const char name[], FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale, PetscBool nameComplex, PetscInt imag) 429552f7358SJed Brown { 43082f516ccSBarry Smith MPI_Comm comm; 431552f7358SJed Brown PetscInt numDof = 0, maxDof; 432552f7358SJed Brown PetscInt pStart, pEnd, p; 433552f7358SJed Brown 434552f7358SJed Brown PetscFunctionBegin; 4359566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm,&comm)); 4369566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 437552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 4389566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &numDof)); 4398865f1eaSKarl Rupp if (numDof) break; 440552f7358SJed Brown } 441552f7358SJed Brown numDof = PetscMax(numDof, enforceDof); 4421c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm))); 443552f7358SJed Brown if (!name) name = "Unknown"; 444552f7358SJed Brown if (maxDof == 3) { 44576b700caSToby Isaac if (nameComplex) { 4469566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "VECTORS %s.%s double\n", name, imag ? "Im" : "Re")); 44776b700caSToby Isaac } else { 4489566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "VECTORS %s double\n", name)); 44976b700caSToby Isaac } 45076b700caSToby Isaac } else { 45176b700caSToby Isaac if (nameComplex) { 45263a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "SCALARS %s.%s double %" PetscInt_FMT "\n", name, imag ? "Im" : "Re", maxDof)); 453552f7358SJed Brown } else { 45463a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "SCALARS %s double %" PetscInt_FMT "\n", name, maxDof)); 45576b700caSToby Isaac } 4569566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n")); 457552f7358SJed Brown } 4589566063dSJacob Faibussowitsch PetscCall(DMPlexVTKWriteSection_ASCII(dm, section, globalSection, field, fp, enforceDof, precision, scale, imag)); 459552f7358SJed Brown PetscFunctionReturn(0); 460552f7358SJed Brown } 461552f7358SJed Brown 462552f7358SJed Brown static PetscErrorCode DMPlexVTKWriteAll_ASCII(DM dm, PetscViewer viewer) 463552f7358SJed Brown { 46482f516ccSBarry Smith MPI_Comm comm; 465552f7358SJed Brown PetscViewer_VTK *vtk = (PetscViewer_VTK*) viewer->data; 466552f7358SJed Brown FILE *fp; 467552f7358SJed Brown PetscViewerVTKObjectLink link; 468412e9a14SMatthew G. Knepley PetscInt totVertices, totCells = 0, loops_per_scalar, l; 46976b700caSToby Isaac PetscBool hasPoint = PETSC_FALSE, hasCell = PETSC_FALSE, writePartition = PETSC_FALSE, localized, writeComplex; 47052b05773SJed Brown const char *dmname; 471552f7358SJed Brown 472552f7358SJed Brown PetscFunctionBegin; 47376b700caSToby Isaac #if defined(PETSC_USE_COMPLEX) 47476b700caSToby Isaac loops_per_scalar = 2; 47576b700caSToby Isaac writeComplex = PETSC_TRUE; 47676b700caSToby Isaac #else 47776b700caSToby Isaac loops_per_scalar = 1; 47876b700caSToby Isaac writeComplex = PETSC_FALSE; 47976b700caSToby Isaac #endif 4809566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalized(dm,&localized)); 4819566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm,&comm)); 4821dca8a05SBarry Smith PetscCheck(!localized,comm,PETSC_ERR_SUP,"VTK output with localized coordinates not yet supported"); 4839566063dSJacob Faibussowitsch PetscCall(PetscFOpen(comm, vtk->filename, "wb", &fp)); 4849566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)dm, &dmname)); 4859566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "# vtk DataFile Version 2.0\n")); 4869566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "%s\n", dmname)); 4879566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "ASCII\n")); 4889566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "DATASET UNSTRUCTURED_GRID\n")); 489552f7358SJed Brown /* Vertices */ 4909120ba30SVaclav Hapla { 4919120ba30SVaclav Hapla PetscSection section, coordSection, globalCoordSection; 4929120ba30SVaclav Hapla Vec coordinates; 4939120ba30SVaclav Hapla PetscReal lengthScale; 4949120ba30SVaclav Hapla DMLabel label; 4959120ba30SVaclav Hapla IS vStratumIS; 4969120ba30SVaclav Hapla PetscLayout vLayout; 4979120ba30SVaclav Hapla 4989566063dSJacob Faibussowitsch PetscCall(DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale)); 4999566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 5009566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(dm, &label)); 5019566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, 0, &vStratumIS)); 5029566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, §ion)); /* This section includes all points */ 5036b57843cSMatthew G. Knepley PetscCall(PetscSectionCreateSubdomainSection(section, vStratumIS, &coordSection)); /* This one includes just vertices */ 5049566063dSJacob Faibussowitsch PetscCall(PetscSectionCreateGlobalSection(coordSection, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalCoordSection)); 5059566063dSJacob Faibussowitsch PetscCall(PetscSectionGetPointLayout(comm, globalCoordSection, &vLayout)); 5069566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(vLayout, &totVertices)); 50763a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "POINTS %" PetscInt_FMT " double\n", totVertices)); 5089566063dSJacob Faibussowitsch PetscCall(DMPlexVTKWriteSection_ASCII(dm, coordSection, globalCoordSection, coordinates, fp, 3, PETSC_DETERMINE, lengthScale, 0)); 5099566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vStratumIS)); 5109566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&vLayout)); 5119566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&coordSection)); 5129566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&globalCoordSection)); 5139120ba30SVaclav Hapla } 514552f7358SJed Brown /* Cells */ 5159566063dSJacob Faibussowitsch PetscCall(DMPlexVTKWriteCells_ASCII(dm, fp, &totCells)); 516552f7358SJed Brown /* Vertex fields */ 517552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 518552f7358SJed Brown if ((link->ft == PETSC_VTK_POINT_FIELD) || (link->ft == PETSC_VTK_POINT_VECTOR_FIELD)) hasPoint = PETSC_TRUE; 519552f7358SJed Brown if ((link->ft == PETSC_VTK_CELL_FIELD) || (link->ft == PETSC_VTK_CELL_VECTOR_FIELD)) hasCell = PETSC_TRUE; 520552f7358SJed Brown } 521552f7358SJed Brown if (hasPoint) { 52263a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "POINT_DATA %" PetscInt_FMT "\n", totVertices)); 523552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 524552f7358SJed Brown Vec X = (Vec) link->vec; 525e630c359SToby Isaac PetscSection section = NULL, globalSection, newSection = NULL; 526e630c359SToby Isaac char namebuf[256]; 527552f7358SJed Brown const char *name; 528552f7358SJed Brown PetscInt enforceDof = PETSC_DETERMINE; 529552f7358SJed Brown 530552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 531552f7358SJed Brown if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) enforceDof = 3; 5329566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName(link->vec, &name)); 5339566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject*) §ion)); 534e630c359SToby Isaac if (!section) { 535e630c359SToby Isaac DM dmX; 536e630c359SToby Isaac 5379566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 538552f7358SJed Brown if (dmX) { 539088580cfSMatthew G Knepley DMLabel subpointMap, subpointMapX; 540839dd189SMatthew G Knepley PetscInt dim, dimX, pStart, pEnd, qStart, qEnd; 541839dd189SMatthew G Knepley 5429566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmX, §ion)); 543839dd189SMatthew G Knepley /* Here is where we check whether dmX is a submesh of dm */ 5449566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 5459566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmX, &dimX)); 5469566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 5479566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dmX, &qStart, &qEnd)); 5489566063dSJacob Faibussowitsch PetscCall(DMPlexGetSubpointMap(dm, &subpointMap)); 5499566063dSJacob Faibussowitsch PetscCall(DMPlexGetSubpointMap(dmX, &subpointMapX)); 550839dd189SMatthew G Knepley if (((dim != dimX) || ((pEnd-pStart) < (qEnd-qStart))) && subpointMap && !subpointMapX) { 551434e42b5SMatthew G. Knepley const PetscInt *ind = NULL; 552e6ccafaeSMatthew G Knepley IS subpointIS; 553434e42b5SMatthew G. Knepley PetscInt n = 0, q; 554839dd189SMatthew G Knepley 5559566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &qStart, &qEnd)); 5569566063dSJacob Faibussowitsch PetscCall(DMPlexGetSubpointIS(dm, &subpointIS)); 557434e42b5SMatthew G. Knepley if (subpointIS) { 5589566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(subpointIS, &n)); 5599566063dSJacob Faibussowitsch PetscCall(ISGetIndices(subpointIS, &ind)); 560434e42b5SMatthew G. Knepley } 5619566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, &newSection)); 5629566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(newSection, pStart, pEnd)); 563839dd189SMatthew G Knepley for (q = qStart; q < qEnd; ++q) { 564839dd189SMatthew G Knepley PetscInt dof, off, p; 565839dd189SMatthew G Knepley 5669566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, q, &dof)); 567839dd189SMatthew G Knepley if (dof) { 5689566063dSJacob Faibussowitsch PetscCall(PetscFindInt(q, n, ind, &p)); 569839dd189SMatthew G Knepley if (p >= pStart) { 5709566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(newSection, p, dof)); 5719566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, q, &off)); 5729566063dSJacob Faibussowitsch PetscCall(PetscSectionSetOffset(newSection, p, off)); 573839dd189SMatthew G Knepley } 574839dd189SMatthew G Knepley } 575839dd189SMatthew G Knepley } 576434e42b5SMatthew G. Knepley if (subpointIS) { 5779566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(subpointIS, &ind)); 578434e42b5SMatthew G. Knepley } 579839dd189SMatthew G Knepley /* No need to setup section */ 580839dd189SMatthew G Knepley section = newSection; 581839dd189SMatthew G Knepley } 582552f7358SJed Brown } 583e630c359SToby Isaac } 58428b400f6SJacob Faibussowitsch PetscCheck(section,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it and could not create one from VecGetDM()", name); 585e630c359SToby Isaac if (link->field >= 0) { 586e630c359SToby Isaac const char *fieldname; 587e630c359SToby Isaac 5889566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(section, link->field, &fieldname)); 5899566063dSJacob Faibussowitsch PetscCall(PetscSectionGetField(section, link->field, §ion)); 590e630c359SToby Isaac if (fieldname) { 5919566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s%s", name, fieldname)); 592e630c359SToby Isaac } else { 59363a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s%" PetscInt_FMT, name, link->field)); 594e630c359SToby Isaac } 5952fb742c9SToby Isaac } else { 5969566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s", name)); 597e630c359SToby Isaac } 5989566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKSanitizeName_Internal(namebuf, sizeof(namebuf))); 5999566063dSJacob Faibussowitsch PetscCall(PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection)); 60076b700caSToby Isaac for (l = 0; l < loops_per_scalar; l++) { 6019566063dSJacob Faibussowitsch PetscCall(DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, namebuf, fp, enforceDof, PETSC_DETERMINE, 1.0, writeComplex, l)); 60276b700caSToby Isaac } 6039566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&globalSection)); 6049566063dSJacob Faibussowitsch if (newSection) PetscCall(PetscSectionDestroy(&newSection)); 605552f7358SJed Brown } 606552f7358SJed Brown } 607552f7358SJed Brown /* Cell Fields */ 6089566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_view_partition", &writePartition, NULL)); 609e8f6f0f6SMatthew G. Knepley if (hasCell || writePartition) { 61063a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "CELL_DATA %" PetscInt_FMT "\n", totCells)); 611552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 612552f7358SJed Brown Vec X = (Vec) link->vec; 613e630c359SToby Isaac PetscSection section = NULL, globalSection; 6142fb742c9SToby Isaac const char *name = ""; 615e630c359SToby Isaac char namebuf[256]; 616552f7358SJed Brown PetscInt enforceDof = PETSC_DETERMINE; 617552f7358SJed Brown 618552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 619552f7358SJed Brown if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) enforceDof = 3; 6209566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName(link->vec, &name)); 6219566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(link->vec, "section", (PetscObject*) §ion)); 622e630c359SToby Isaac if (!section) { 623e630c359SToby Isaac DM dmX; 624e630c359SToby Isaac 6259566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dmX)); 626552f7358SJed Brown if (dmX) { 6279566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmX, §ion)); 628552f7358SJed Brown } 629e630c359SToby Isaac } 63028b400f6SJacob Faibussowitsch PetscCheck(section,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it and could not create one from VecGetDM()", name); 631e630c359SToby Isaac if (link->field >= 0) { 632e630c359SToby Isaac const char *fieldname; 633e630c359SToby Isaac 6349566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(section, link->field, &fieldname)); 6359566063dSJacob Faibussowitsch PetscCall(PetscSectionGetField(section, link->field, §ion)); 636e630c359SToby Isaac if (fieldname) { 6379566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s%s", name, fieldname)); 638e630c359SToby Isaac } else { 63963a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s%" PetscInt_FMT, name, link->field)); 640e630c359SToby Isaac } 6412fb742c9SToby Isaac } else { 6429566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(namebuf, sizeof(namebuf), "%s", name)); 643e630c359SToby Isaac } 6449566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKSanitizeName_Internal(namebuf, sizeof(namebuf))); 6459566063dSJacob Faibussowitsch PetscCall(PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection)); 64676b700caSToby Isaac for (l = 0; l < loops_per_scalar; l++) { 6479566063dSJacob Faibussowitsch PetscCall(DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, namebuf, fp, enforceDof, PETSC_DETERMINE, 1.0, writeComplex, l)); 64876b700caSToby Isaac } 6499566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&globalSection)); 650552f7358SJed Brown } 651e8f6f0f6SMatthew G. Knepley if (writePartition) { 6529566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "SCALARS partition int 1\n")); 6539566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n")); 6549566063dSJacob Faibussowitsch PetscCall(DMPlexVTKWritePartition_ASCII(dm, fp)); 655e8f6f0f6SMatthew G. Knepley } 656552f7358SJed Brown } 657552f7358SJed Brown /* Cleanup */ 6589566063dSJacob Faibussowitsch PetscCall(PetscFClose(comm, fp)); 659552f7358SJed Brown PetscFunctionReturn(0); 660552f7358SJed Brown } 661552f7358SJed Brown 662552f7358SJed Brown /*@C 663552f7358SJed Brown DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer 664552f7358SJed Brown 665552f7358SJed Brown Collective 666552f7358SJed Brown 6674165533cSJose E. Roman Input Parameters: 668552f7358SJed Brown + odm - The DMPlex specifying the mesh, passed as a PetscObject 669552f7358SJed Brown - viewer - viewer of type VTK 670552f7358SJed Brown 671552f7358SJed Brown Level: developer 672552f7358SJed Brown 673552f7358SJed Brown Note: 674552f7358SJed Brown This function is a callback used by the VTK viewer to actually write the file. 675552f7358SJed Brown The reason for this odd model is that the VTK file format does not provide any way to write one field at a time. 676552f7358SJed Brown Instead, metadata for the entire file needs to be available up-front before you can start writing the file. 677552f7358SJed Brown 678db781477SPatrick Sanan .seealso: `PETSCVIEWERVTK` 679552f7358SJed Brown @*/ 680552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer) 681552f7358SJed Brown { 682552f7358SJed Brown DM dm = (DM) odm; 683552f7358SJed Brown PetscBool isvtk; 684552f7358SJed Brown 685552f7358SJed Brown PetscFunctionBegin; 686552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 687552f7358SJed Brown PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 6889566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk)); 68928b400f6SJacob Faibussowitsch PetscCheck(isvtk,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name); 690552f7358SJed Brown switch (viewer->format) { 6918ec8862eSJed Brown case PETSC_VIEWER_ASCII_VTK_DEPRECATED: 6929566063dSJacob Faibussowitsch PetscCall(DMPlexVTKWriteAll_ASCII(dm, viewer)); 693552f7358SJed Brown break; 694552f7358SJed Brown case PETSC_VIEWER_VTK_VTU: 6959566063dSJacob Faibussowitsch PetscCall(DMPlexVTKWriteAll_VTU(dm, viewer)); 696552f7358SJed Brown break; 69798921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]); 698552f7358SJed Brown } 699552f7358SJed Brown PetscFunctionReturn(0); 700552f7358SJed Brown } 701