1552f7358SJed Brown #define PETSCDM_DLL 234541f0dSBarry Smith #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3552f7358SJed Brown #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h> 4552f7358SJed Brown 5552f7358SJed Brown #undef __FUNCT__ 6552f7358SJed Brown #define __FUNCT__ "DMPlexVTKGetCellType" 70adebc6cSBarry Smith PetscErrorCode DMPlexVTKGetCellType(DM dm, PetscInt dim, PetscInt corners, PetscInt *cellType) 80adebc6cSBarry Smith { 9552f7358SJed Brown PetscFunctionBegin; 10552f7358SJed Brown *cellType = -1; 11552f7358SJed Brown switch (dim) { 12552f7358SJed Brown case 0: 13552f7358SJed Brown switch (corners) { 14552f7358SJed Brown case 1: 15552f7358SJed Brown *cellType = 1; /* VTK_VERTEX */ 16552f7358SJed Brown break; 17552f7358SJed Brown default: 18552f7358SJed Brown break; 19552f7358SJed Brown } 20552f7358SJed Brown break; 21552f7358SJed Brown case 1: 22552f7358SJed Brown switch (corners) { 23552f7358SJed Brown case 2: 24552f7358SJed Brown *cellType = 3; /* VTK_LINE */ 25552f7358SJed Brown break; 26552f7358SJed Brown case 3: 27552f7358SJed Brown *cellType = 21; /* VTK_QUADRATIC_EDGE */ 28552f7358SJed Brown break; 29552f7358SJed Brown default: 30552f7358SJed Brown break; 31552f7358SJed Brown } 32552f7358SJed Brown break; 33552f7358SJed Brown case 2: 34552f7358SJed Brown switch (corners) { 35552f7358SJed Brown case 3: 36552f7358SJed Brown *cellType = 5; /* VTK_TRIANGLE */ 37552f7358SJed Brown break; 38552f7358SJed Brown case 4: 39552f7358SJed Brown *cellType = 9; /* VTK_QUAD */ 40552f7358SJed Brown break; 41552f7358SJed Brown case 6: 42552f7358SJed Brown *cellType = 22; /* VTK_QUADRATIC_TRIANGLE */ 43552f7358SJed Brown break; 44552f7358SJed Brown case 9: 45552f7358SJed Brown *cellType = 23; /* VTK_QUADRATIC_QUAD */ 46552f7358SJed Brown break; 47552f7358SJed Brown default: 48552f7358SJed Brown break; 49552f7358SJed Brown } 50552f7358SJed Brown break; 51552f7358SJed Brown case 3: 52552f7358SJed Brown switch (corners) { 53552f7358SJed Brown case 4: 54552f7358SJed Brown *cellType = 10; /* VTK_TETRA */ 55552f7358SJed Brown break; 56552f7358SJed Brown case 8: 57552f7358SJed Brown *cellType = 12; /* VTK_HEXAHEDRON */ 58552f7358SJed Brown break; 59552f7358SJed Brown case 10: 60552f7358SJed Brown *cellType = 24; /* VTK_QUADRATIC_TETRA */ 61552f7358SJed Brown break; 62552f7358SJed Brown case 27: 63552f7358SJed Brown *cellType = 29; /* VTK_QUADRATIC_HEXAHEDRON */ 64552f7358SJed Brown break; 65552f7358SJed Brown default: 66552f7358SJed Brown break; 67552f7358SJed Brown } 68552f7358SJed Brown } 69552f7358SJed Brown PetscFunctionReturn(0); 70552f7358SJed Brown } 71552f7358SJed Brown 72552f7358SJed Brown #undef __FUNCT__ 73552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteCells_ASCII" 74552f7358SJed Brown PetscErrorCode DMPlexVTKWriteCells_ASCII(DM dm, FILE *fp, PetscInt *totalCells) 75552f7358SJed Brown { 7682f516ccSBarry Smith MPI_Comm comm; 77552f7358SJed Brown IS globalVertexNumbers; 78552f7358SJed Brown const PetscInt *gvertex; 79552f7358SJed Brown PetscInt dim; 80552f7358SJed Brown PetscInt numCorners = 0, totCorners = 0, maxCorners, *corners; 81552f7358SJed Brown PetscInt numCells = 0, totCells = 0, maxCells, cellHeight; 82552f7358SJed Brown PetscInt numLabelCells, cMax, cStart, cEnd, c, vStart, vEnd, v; 83552f7358SJed Brown PetscMPIInt numProcs, rank, proc, tag; 84552f7358SJed Brown PetscBool hasLabel; 85552f7358SJed Brown PetscErrorCode ierr; 86552f7358SJed Brown 87552f7358SJed Brown PetscFunctionBegin; 8882f516ccSBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 89552f7358SJed Brown ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr); 90552f7358SJed Brown ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 91552f7358SJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 92552f7358SJed Brown ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 93552f7358SJed Brown ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 94552f7358SJed Brown ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 95552f7358SJed Brown ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 960298fd71SBarry Smith ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 978865f1eaSKarl Rupp if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 98552f7358SJed Brown ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 998865f1eaSKarl Rupp 100552f7358SJed Brown hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 101552f7358SJed Brown for (c = cStart; c < cEnd; ++c) { 1020298fd71SBarry Smith PetscInt *closure = NULL; 103552f7358SJed Brown PetscInt closureSize; 104552f7358SJed Brown 105552f7358SJed Brown if (hasLabel) { 106552f7358SJed Brown PetscInt value; 107552f7358SJed Brown 108552f7358SJed Brown ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 109552f7358SJed Brown if (value != 1) continue; 110552f7358SJed Brown } 111552f7358SJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 112552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 113552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) ++numCorners; 114552f7358SJed Brown } 115552f7358SJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 116552f7358SJed Brown ++numCells; 117552f7358SJed Brown } 118552f7358SJed Brown maxCells = numCells; 119552f7358SJed Brown ierr = MPI_Reduce(&numCells, &totCells, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRQ(ierr); 120552f7358SJed Brown ierr = MPI_Reduce(&numCells, &maxCells, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr); 121552f7358SJed Brown ierr = MPI_Reduce(&numCorners, &totCorners, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRQ(ierr); 122552f7358SJed Brown ierr = MPI_Reduce(&numCorners, &maxCorners, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr); 123552f7358SJed Brown ierr = DMPlexGetVertexNumbering(dm, &globalVertexNumbers);CHKERRQ(ierr); 124552f7358SJed Brown ierr = ISGetIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr); 125552f7358SJed Brown ierr = PetscMalloc(maxCells * sizeof(PetscInt), &corners);CHKERRQ(ierr); 126552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "CELLS %d %d\n", totCells, totCorners+totCells);CHKERRQ(ierr); 127552f7358SJed Brown if (!rank) { 128552f7358SJed Brown PetscInt *remoteVertices; 129552f7358SJed Brown 130552f7358SJed Brown for (c = cStart, numCells = 0; c < cEnd; ++c) { 1310298fd71SBarry Smith PetscInt *closure = NULL; 132baf8153dSMatthew G. Knepley PetscInt closureSize, nC = 0, tmp; 133552f7358SJed Brown 134552f7358SJed Brown if (hasLabel) { 135552f7358SJed Brown PetscInt value; 136552f7358SJed Brown 137552f7358SJed Brown ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 138552f7358SJed Brown if (value != 1) continue; 139552f7358SJed Brown } 140552f7358SJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 141552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 142552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 143552f7358SJed Brown const PetscInt gv = gvertex[closure[v] - vStart]; 144552f7358SJed Brown closure[nC++] = gv < 0 ? -(gv+1) : gv; 145552f7358SJed Brown } 146552f7358SJed Brown } 147552f7358SJed Brown corners[numCells++] = nC; 148552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "%d ", nC);CHKERRQ(ierr); 149baf8153dSMatthew G. Knepley tmp = closure[0]; 150baf8153dSMatthew G. Knepley closure[0] = closure[1]; 151baf8153dSMatthew G. Knepley closure[1] = tmp; 152552f7358SJed Brown for (v = 0; v < nC; ++v) { 153552f7358SJed Brown ierr = PetscFPrintf(comm, fp, " %d", closure[v]);CHKERRQ(ierr); 154552f7358SJed Brown } 155552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr); 156552f7358SJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 157552f7358SJed Brown } 158552f7358SJed Brown if (numProcs > 1) {ierr = PetscMalloc((maxCorners+maxCells) * sizeof(PetscInt), &remoteVertices);CHKERRQ(ierr);} 159552f7358SJed Brown for (proc = 1; proc < numProcs; ++proc) { 160552f7358SJed Brown MPI_Status status; 161552f7358SJed Brown 162552f7358SJed Brown ierr = MPI_Recv(&numCorners, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 163552f7358SJed Brown ierr = MPI_Recv(remoteVertices, numCorners, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 164552f7358SJed Brown for (c = 0; c < numCorners;) { 165552f7358SJed Brown PetscInt nC = remoteVertices[c++]; 166552f7358SJed Brown 167552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "%d ", nC);CHKERRQ(ierr); 168552f7358SJed Brown for (v = 0; v < nC; ++v, ++c) { 169552f7358SJed Brown ierr = PetscFPrintf(comm, fp, " %d", remoteVertices[c]);CHKERRQ(ierr); 170552f7358SJed Brown } 171552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr); 172552f7358SJed Brown } 173552f7358SJed Brown } 174552f7358SJed Brown if (numProcs > 1) {ierr = PetscFree(remoteVertices);CHKERRQ(ierr);} 175552f7358SJed Brown } else { 176552f7358SJed Brown PetscInt *localVertices, numSend = numCells+numCorners, k = 0; 177552f7358SJed Brown 178552f7358SJed Brown ierr = PetscMalloc(numSend * sizeof(PetscInt), &localVertices);CHKERRQ(ierr); 179552f7358SJed Brown for (c = cStart, numCells = 0; c < cEnd; ++c) { 1800298fd71SBarry Smith PetscInt *closure = NULL; 181552f7358SJed Brown PetscInt closureSize, nC = 0; 182552f7358SJed Brown 183552f7358SJed Brown if (hasLabel) { 184552f7358SJed Brown PetscInt value; 185552f7358SJed Brown 186552f7358SJed Brown ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 187552f7358SJed Brown if (value != 1) continue; 188552f7358SJed Brown } 189552f7358SJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 190552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 191552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 192552f7358SJed Brown const PetscInt gv = gvertex[closure[v] - vStart]; 193552f7358SJed Brown closure[nC++] = gv < 0 ? -(gv+1) : gv; 194552f7358SJed Brown } 195552f7358SJed Brown } 196552f7358SJed Brown corners[numCells++] = nC; 197552f7358SJed Brown localVertices[k++] = nC; 198552f7358SJed Brown for (v = 0; v < nC; ++v, ++k) { 199552f7358SJed Brown localVertices[k] = closure[v]; 200552f7358SJed Brown } 201552f7358SJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 202552f7358SJed Brown } 203552f7358SJed Brown if (k != numSend) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid number of vertices to send %d should be %d", k, numSend); 204552f7358SJed Brown ierr = MPI_Send(&numSend, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 205552f7358SJed Brown ierr = MPI_Send(localVertices, numSend, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 206552f7358SJed Brown ierr = PetscFree(localVertices);CHKERRQ(ierr); 207552f7358SJed Brown } 208552f7358SJed Brown ierr = ISRestoreIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr); 209552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "CELL_TYPES %d\n", totCells);CHKERRQ(ierr); 210552f7358SJed Brown if (!rank) { 211552f7358SJed Brown PetscInt cellType; 212552f7358SJed Brown 213552f7358SJed Brown for (c = 0; c < numCells; ++c) { 214552f7358SJed Brown ierr = DMPlexVTKGetCellType(dm, dim, corners[c], &cellType);CHKERRQ(ierr); 215552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "%d\n", cellType);CHKERRQ(ierr); 216552f7358SJed Brown } 217552f7358SJed Brown for (proc = 1; proc < numProcs; ++proc) { 218552f7358SJed Brown MPI_Status status; 219552f7358SJed Brown 220552f7358SJed Brown ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 221552f7358SJed Brown ierr = MPI_Recv(corners, numCells, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 222552f7358SJed Brown for (c = 0; c < numCells; ++c) { 223552f7358SJed Brown ierr = DMPlexVTKGetCellType(dm, dim, corners[c], &cellType);CHKERRQ(ierr); 224552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "%d\n", cellType);CHKERRQ(ierr); 225552f7358SJed Brown } 226552f7358SJed Brown } 227552f7358SJed Brown } else { 228552f7358SJed Brown ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 229552f7358SJed Brown ierr = MPI_Send(corners, numCells, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 230552f7358SJed Brown } 231552f7358SJed Brown ierr = PetscFree(corners);CHKERRQ(ierr); 232552f7358SJed Brown *totalCells = totCells; 233552f7358SJed Brown PetscFunctionReturn(0); 234552f7358SJed Brown } 235552f7358SJed Brown 236552f7358SJed Brown #undef __FUNCT__ 237*e8f6f0f6SMatthew G. Knepley #define __FUNCT__ "DMPlexVTKWritePartition_ASCII" 238*e8f6f0f6SMatthew G. Knepley PetscErrorCode DMPlexVTKWritePartition_ASCII(DM dm, FILE *fp) 239*e8f6f0f6SMatthew G. Knepley { 240*e8f6f0f6SMatthew G. Knepley MPI_Comm comm; 241*e8f6f0f6SMatthew G. Knepley PetscInt numCells = 0, cellHeight; 242*e8f6f0f6SMatthew G. Knepley PetscInt numLabelCells, cMax, cStart, cEnd, c; 243*e8f6f0f6SMatthew G. Knepley PetscMPIInt numProcs, rank, proc, tag; 244*e8f6f0f6SMatthew G. Knepley PetscBool hasLabel; 245*e8f6f0f6SMatthew G. Knepley PetscErrorCode ierr; 246*e8f6f0f6SMatthew G. Knepley 247*e8f6f0f6SMatthew G. Knepley PetscFunctionBegin; 248*e8f6f0f6SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 249*e8f6f0f6SMatthew G. Knepley ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr); 250*e8f6f0f6SMatthew G. Knepley ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 251*e8f6f0f6SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 252*e8f6f0f6SMatthew G. Knepley ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 253*e8f6f0f6SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 254*e8f6f0f6SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 255*e8f6f0f6SMatthew G. Knepley if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 256*e8f6f0f6SMatthew G. Knepley ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 257*e8f6f0f6SMatthew G. Knepley hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 258*e8f6f0f6SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 259*e8f6f0f6SMatthew G. Knepley if (hasLabel) { 260*e8f6f0f6SMatthew G. Knepley PetscInt value; 261*e8f6f0f6SMatthew G. Knepley 262*e8f6f0f6SMatthew G. Knepley ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 263*e8f6f0f6SMatthew G. Knepley if (value != 1) continue; 264*e8f6f0f6SMatthew G. Knepley } 265*e8f6f0f6SMatthew G. Knepley ++numCells; 266*e8f6f0f6SMatthew G. Knepley } 267*e8f6f0f6SMatthew G. Knepley if (!rank) { 268*e8f6f0f6SMatthew G. Knepley for (c = 0; c < numCells; ++c) {ierr = PetscFPrintf(comm, fp, "%d\n", rank);CHKERRQ(ierr);} 269*e8f6f0f6SMatthew G. Knepley for (proc = 1; proc < numProcs; ++proc) { 270*e8f6f0f6SMatthew G. Knepley MPI_Status status; 271*e8f6f0f6SMatthew G. Knepley 272*e8f6f0f6SMatthew G. Knepley ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 273*e8f6f0f6SMatthew G. Knepley for (c = 0; c < numCells; ++c) {ierr = PetscFPrintf(comm, fp, "%d\n", proc);CHKERRQ(ierr);} 274*e8f6f0f6SMatthew G. Knepley } 275*e8f6f0f6SMatthew G. Knepley } else { 276*e8f6f0f6SMatthew G. Knepley ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 277*e8f6f0f6SMatthew G. Knepley } 278*e8f6f0f6SMatthew G. Knepley PetscFunctionReturn(0); 279*e8f6f0f6SMatthew G. Knepley } 280*e8f6f0f6SMatthew G. Knepley 281*e8f6f0f6SMatthew G. Knepley #undef __FUNCT__ 282552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteSection_ASCII" 2830adebc6cSBarry Smith PetscErrorCode DMPlexVTKWriteSection_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec v, FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale) 2840adebc6cSBarry Smith { 28582f516ccSBarry Smith MPI_Comm comm; 286552f7358SJed Brown const MPI_Datatype mpiType = MPIU_SCALAR; 287552f7358SJed Brown PetscScalar *array; 288552f7358SJed Brown PetscInt numDof = 0, maxDof; 289552f7358SJed Brown PetscInt numLabelCells, cellHeight, cMax, cStart, cEnd, numLabelVertices, vMax, vStart, vEnd, pStart, pEnd, p; 290552f7358SJed Brown PetscMPIInt numProcs, rank, proc, tag; 291552f7358SJed Brown PetscBool hasLabel; 292552f7358SJed Brown PetscErrorCode ierr; 293552f7358SJed Brown 294552f7358SJed Brown PetscFunctionBegin; 29582f516ccSBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 296552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 297552f7358SJed Brown PetscValidHeaderSpecific(v,VEC_CLASSID,4); 298552f7358SJed Brown if (precision < 0) precision = 6; 299552f7358SJed Brown ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr); 300552f7358SJed Brown ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 301552f7358SJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 302552f7358SJed Brown ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 303552f7358SJed Brown /* VTK only wants the values at cells or vertices */ 304552f7358SJed Brown ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 305552f7358SJed Brown ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 306552f7358SJed Brown ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 3070298fd71SBarry Smith ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, &vMax);CHKERRQ(ierr); 3088865f1eaSKarl Rupp if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 3098865f1eaSKarl Rupp if (vMax >= 0) vEnd = PetscMin(vEnd, vMax); 310552f7358SJed Brown pStart = PetscMax(PetscMin(cStart, vStart), pStart); 311552f7358SJed Brown pEnd = PetscMin(PetscMax(cEnd, vEnd), pEnd); 312552f7358SJed Brown ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 313552f7358SJed Brown ierr = DMPlexGetStratumSize(dm, "vtk", 2, &numLabelVertices);CHKERRQ(ierr); 314552f7358SJed Brown hasLabel = numLabelCells > 0 || numLabelVertices > 0 ? PETSC_TRUE : PETSC_FALSE; 315552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 316552f7358SJed Brown /* Reject points not either cells or vertices */ 317552f7358SJed Brown if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue; 318552f7358SJed Brown if (hasLabel) { 319552f7358SJed Brown PetscInt value; 320552f7358SJed Brown 321552f7358SJed Brown if (((p >= cStart) && (p < cEnd) && numLabelCells) || 322552f7358SJed Brown ((p >= vStart) && (p < vEnd) && numLabelVertices)) { 323552f7358SJed Brown ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr); 324552f7358SJed Brown if (value != 1) continue; 325552f7358SJed Brown } 326552f7358SJed Brown } 327552f7358SJed Brown ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr); 3288865f1eaSKarl Rupp if (numDof) break; 329552f7358SJed Brown } 330552f7358SJed Brown ierr = MPI_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr); 331552f7358SJed Brown enforceDof = PetscMax(enforceDof, maxDof); 332552f7358SJed Brown ierr = VecGetArray(v, &array);CHKERRQ(ierr); 333552f7358SJed Brown if (!rank) { 334552f7358SJed Brown char formatString[8]; 335552f7358SJed Brown 336552f7358SJed Brown ierr = PetscSNPrintf(formatString, 8, "%%.%de", precision);CHKERRQ(ierr); 337552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 338552f7358SJed Brown /* Here we lose a way to filter points by keeping them out of the Numbering */ 339552f7358SJed Brown PetscInt dof, off, goff, d; 340552f7358SJed Brown 341552f7358SJed Brown /* Reject points not either cells or vertices */ 342552f7358SJed Brown if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue; 343552f7358SJed Brown if (hasLabel) { 344552f7358SJed Brown PetscInt value; 345552f7358SJed Brown 346552f7358SJed Brown if (((p >= cStart) && (p < cEnd) && numLabelCells) || 347552f7358SJed Brown ((p >= vStart) && (p < vEnd) && numLabelVertices)) { 348552f7358SJed Brown ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr); 349552f7358SJed Brown if (value != 1) continue; 350552f7358SJed Brown } 351552f7358SJed Brown } 352552f7358SJed Brown ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 353552f7358SJed Brown ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 354552f7358SJed Brown ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr); 355552f7358SJed Brown if (dof && goff >= 0) { 356552f7358SJed Brown for (d = 0; d < dof; d++) { 357552f7358SJed Brown if (d > 0) { 358552f7358SJed Brown ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr); 359552f7358SJed Brown } 360552f7358SJed Brown ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(array[off+d])*scale);CHKERRQ(ierr); 361552f7358SJed Brown } 362552f7358SJed Brown for (d = dof; d < enforceDof; d++) { 363552f7358SJed Brown ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr); 364552f7358SJed Brown } 365552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr); 366552f7358SJed Brown } 367552f7358SJed Brown } 368552f7358SJed Brown for (proc = 1; proc < numProcs; ++proc) { 369552f7358SJed Brown PetscScalar *remoteValues; 370552f7358SJed Brown PetscInt size, d; 371552f7358SJed Brown MPI_Status status; 372552f7358SJed Brown 373552f7358SJed Brown ierr = MPI_Recv(&size, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 374552f7358SJed Brown ierr = PetscMalloc(size * sizeof(PetscScalar), &remoteValues);CHKERRQ(ierr); 375552f7358SJed Brown ierr = MPI_Recv(remoteValues, size, mpiType, proc, tag, comm, &status);CHKERRQ(ierr); 376552f7358SJed Brown for (p = 0; p < size/maxDof; ++p) { 377552f7358SJed Brown for (d = 0; d < maxDof; ++d) { 378552f7358SJed Brown if (d > 0) { 379552f7358SJed Brown ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr); 380552f7358SJed Brown } 381552f7358SJed Brown ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(remoteValues[p*maxDof+d])*scale);CHKERRQ(ierr); 382552f7358SJed Brown } 383552f7358SJed Brown for (d = maxDof; d < enforceDof; ++d) { 384552f7358SJed Brown ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr); 385552f7358SJed Brown } 386552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr); 387552f7358SJed Brown } 388552f7358SJed Brown ierr = PetscFree(remoteValues);CHKERRQ(ierr); 389552f7358SJed Brown } 390552f7358SJed Brown } else { 391552f7358SJed Brown PetscScalar *localValues; 392552f7358SJed Brown PetscInt size, k = 0; 393552f7358SJed Brown 394552f7358SJed Brown ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr); 395552f7358SJed Brown ierr = PetscMalloc(size * sizeof(PetscScalar), &localValues);CHKERRQ(ierr); 396552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 397552f7358SJed Brown PetscInt dof, off, goff, d; 398552f7358SJed Brown 399552f7358SJed Brown /* Reject points not either cells or vertices */ 400552f7358SJed Brown if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue; 401552f7358SJed Brown if (hasLabel) { 402552f7358SJed Brown PetscInt value; 403552f7358SJed Brown 404552f7358SJed Brown if (((p >= cStart) && (p < cEnd) && numLabelCells) || 405552f7358SJed Brown ((p >= vStart) && (p < vEnd) && numLabelVertices)) { 406552f7358SJed Brown ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr); 407552f7358SJed Brown if (value != 1) continue; 408552f7358SJed Brown } 409552f7358SJed Brown } 410552f7358SJed Brown ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 411552f7358SJed Brown ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 412552f7358SJed Brown ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr); 413552f7358SJed Brown if (goff >= 0) { 414552f7358SJed Brown for (d = 0; d < dof; ++d) { 415552f7358SJed Brown localValues[k++] = array[off+d]; 416552f7358SJed Brown } 417552f7358SJed Brown } 418552f7358SJed Brown } 419552f7358SJed Brown ierr = MPI_Send(&k, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 420552f7358SJed Brown ierr = MPI_Send(localValues, k, mpiType, 0, tag, comm);CHKERRQ(ierr); 421552f7358SJed Brown ierr = PetscFree(localValues);CHKERRQ(ierr); 422552f7358SJed Brown } 423552f7358SJed Brown ierr = VecRestoreArray(v, &array);CHKERRQ(ierr); 424552f7358SJed Brown PetscFunctionReturn(0); 425552f7358SJed Brown } 426552f7358SJed Brown 427552f7358SJed Brown #undef __FUNCT__ 428552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteField_ASCII" 429552f7358SJed Brown PetscErrorCode DMPlexVTKWriteField_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec field, const char name[], FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale) 430552f7358SJed Brown { 43182f516ccSBarry Smith MPI_Comm comm; 432552f7358SJed Brown PetscInt numDof = 0, maxDof; 433552f7358SJed Brown PetscInt pStart, pEnd, p; 434552f7358SJed Brown PetscErrorCode ierr; 435552f7358SJed Brown 436552f7358SJed Brown PetscFunctionBegin; 43782f516ccSBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 438552f7358SJed Brown ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 439552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 440552f7358SJed Brown ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr); 4418865f1eaSKarl Rupp if (numDof) break; 442552f7358SJed Brown } 443552f7358SJed Brown numDof = PetscMax(numDof, enforceDof); 44482f516ccSBarry Smith ierr = MPI_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 445552f7358SJed Brown if (!name) name = "Unknown"; 446552f7358SJed Brown if (maxDof == 3) { 447552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "VECTORS %s double\n", name);CHKERRQ(ierr); 448552f7358SJed Brown } else { 449552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "SCALARS %s double %d\n", name, maxDof);CHKERRQ(ierr); 450552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr); 451552f7358SJed Brown } 452552f7358SJed Brown ierr = DMPlexVTKWriteSection_ASCII(dm, section, globalSection, field, fp, enforceDof, precision, scale);CHKERRQ(ierr); 453552f7358SJed Brown PetscFunctionReturn(0); 454552f7358SJed Brown } 455552f7358SJed Brown 456552f7358SJed Brown #undef __FUNCT__ 457552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll_ASCII" 458552f7358SJed Brown static PetscErrorCode DMPlexVTKWriteAll_ASCII(DM dm, PetscViewer viewer) 459552f7358SJed Brown { 46082f516ccSBarry Smith MPI_Comm comm; 461552f7358SJed Brown PetscViewer_VTK *vtk = (PetscViewer_VTK*) viewer->data; 462552f7358SJed Brown FILE *fp; 463552f7358SJed Brown PetscViewerVTKObjectLink link; 464552f7358SJed Brown PetscSection coordSection, globalCoordSection; 465552f7358SJed Brown PetscLayout vLayout; 466552f7358SJed Brown Vec coordinates; 467552f7358SJed Brown PetscReal lengthScale; 468552f7358SJed Brown PetscInt vMax, totVertices, totCells; 469*e8f6f0f6SMatthew G. Knepley PetscBool hasPoint = PETSC_FALSE, hasCell = PETSC_FALSE, writePartition = PETSC_FALSE; 470552f7358SJed Brown PetscErrorCode ierr; 471552f7358SJed Brown 472552f7358SJed Brown PetscFunctionBegin; 47382f516ccSBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 474552f7358SJed Brown ierr = PetscFOpen(comm, vtk->filename, "wb", &fp);CHKERRQ(ierr); 475552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "# vtk DataFile Version 2.0\n");CHKERRQ(ierr); 476552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "Simplicial Mesh Example\n");CHKERRQ(ierr); 477552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "ASCII\n");CHKERRQ(ierr); 478552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "DATASET UNSTRUCTURED_GRID\n");CHKERRQ(ierr); 479552f7358SJed Brown /* Vertices */ 480552f7358SJed Brown ierr = DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);CHKERRQ(ierr); 481552f7358SJed Brown ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 482552f7358SJed Brown ierr = PetscSectionCreateGlobalSection(coordSection, dm->sf, PETSC_FALSE, &globalCoordSection);CHKERRQ(ierr); 483552f7358SJed Brown ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 4840298fd71SBarry Smith ierr = DMPlexGetHybridBounds(dm, NULL, NULL, NULL, &vMax);CHKERRQ(ierr); 485552f7358SJed Brown if (vMax >= 0) { 486552f7358SJed Brown PetscInt pStart, pEnd, p, localSize = 0; 487552f7358SJed Brown 488552f7358SJed Brown ierr = PetscSectionGetChart(globalCoordSection, &pStart, &pEnd);CHKERRQ(ierr); 489552f7358SJed Brown pEnd = PetscMin(pEnd, vMax); 490552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 491552f7358SJed Brown PetscInt dof; 492552f7358SJed Brown 493552f7358SJed Brown ierr = PetscSectionGetDof(globalCoordSection, p, &dof);CHKERRQ(ierr); 4948865f1eaSKarl Rupp if (dof > 0) ++localSize; 495552f7358SJed Brown } 49682f516ccSBarry Smith ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &vLayout);CHKERRQ(ierr); 497552f7358SJed Brown ierr = PetscLayoutSetLocalSize(vLayout, localSize);CHKERRQ(ierr); 498552f7358SJed Brown ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr); 499552f7358SJed Brown ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr); 500552f7358SJed Brown } else { 50182f516ccSBarry Smith ierr = PetscSectionGetPointLayout(PetscObjectComm((PetscObject)dm), globalCoordSection, &vLayout);CHKERRQ(ierr); 502552f7358SJed Brown } 503552f7358SJed Brown ierr = PetscLayoutGetSize(vLayout, &totVertices);CHKERRQ(ierr); 504552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "POINTS %d double\n", totVertices);CHKERRQ(ierr); 505552f7358SJed Brown ierr = DMPlexVTKWriteSection_ASCII(dm, coordSection, globalCoordSection, coordinates, fp, 3, PETSC_DETERMINE, lengthScale);CHKERRQ(ierr); 506552f7358SJed Brown /* Cells */ 507552f7358SJed Brown ierr = DMPlexVTKWriteCells_ASCII(dm, fp, &totCells);CHKERRQ(ierr); 508552f7358SJed Brown /* Vertex fields */ 509552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 510552f7358SJed Brown if ((link->ft == PETSC_VTK_POINT_FIELD) || (link->ft == PETSC_VTK_POINT_VECTOR_FIELD)) hasPoint = PETSC_TRUE; 511552f7358SJed Brown if ((link->ft == PETSC_VTK_CELL_FIELD) || (link->ft == PETSC_VTK_CELL_VECTOR_FIELD)) hasCell = PETSC_TRUE; 512552f7358SJed Brown } 513552f7358SJed Brown if (hasPoint) { 51422d28d08SBarry Smith ierr = PetscFPrintf(comm, fp, "POINT_DATA %d\n", totVertices);CHKERRQ(ierr); 515552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 516552f7358SJed Brown Vec X = (Vec) link->vec; 517552f7358SJed Brown DM dmX; 5180298fd71SBarry Smith PetscSection section, globalSection, newSection = NULL; 519552f7358SJed Brown const char *name; 520552f7358SJed Brown PetscInt enforceDof = PETSC_DETERMINE; 521552f7358SJed Brown 522552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 523552f7358SJed Brown if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) enforceDof = 3; 524552f7358SJed Brown ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr); 525552f7358SJed Brown ierr = VecGetDM(X, &dmX);CHKERRQ(ierr); 526552f7358SJed Brown if (dmX) { 527088580cfSMatthew G Knepley DMLabel subpointMap, subpointMapX; 528839dd189SMatthew G Knepley PetscInt dim, dimX, pStart, pEnd, qStart, qEnd; 529839dd189SMatthew G Knepley 53022d28d08SBarry Smith ierr = DMGetDefaultSection(dmX, §ion);CHKERRQ(ierr); 531839dd189SMatthew G Knepley /* Here is where we check whether dmX is a submesh of dm */ 532839dd189SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 533839dd189SMatthew G Knepley ierr = DMPlexGetDimension(dmX, &dimX);CHKERRQ(ierr); 534839dd189SMatthew G Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 535839dd189SMatthew G Knepley ierr = DMPlexGetChart(dmX, &qStart, &qEnd);CHKERRQ(ierr); 536839dd189SMatthew G Knepley ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr); 537839dd189SMatthew G Knepley ierr = DMPlexGetSubpointMap(dmX, &subpointMapX);CHKERRQ(ierr); 538839dd189SMatthew G Knepley if (((dim != dimX) || ((pEnd-pStart) < (qEnd-qStart))) && subpointMap && !subpointMapX) { 539839dd189SMatthew G Knepley const PetscInt *ind; 540e6ccafaeSMatthew G Knepley IS subpointIS; 541839dd189SMatthew G Knepley PetscInt n, q; 542839dd189SMatthew G Knepley 543839dd189SMatthew G Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "Making translation PetscSection\n");CHKERRQ(ierr); 544839dd189SMatthew G Knepley ierr = PetscSectionGetChart(section, &qStart, &qEnd);CHKERRQ(ierr); 545e6ccafaeSMatthew G Knepley ierr = DMPlexCreateSubpointIS(dm, &subpointIS);CHKERRQ(ierr); 546e6ccafaeSMatthew G Knepley ierr = ISGetLocalSize(subpointIS, &n);CHKERRQ(ierr); 547e6ccafaeSMatthew G Knepley ierr = ISGetIndices(subpointIS, &ind);CHKERRQ(ierr); 548839dd189SMatthew G Knepley ierr = PetscSectionCreate(comm, &newSection);CHKERRQ(ierr); 549839dd189SMatthew G Knepley ierr = PetscSectionSetChart(newSection, pStart, pEnd);CHKERRQ(ierr); 550839dd189SMatthew G Knepley for (q = qStart; q < qEnd; ++q) { 551839dd189SMatthew G Knepley PetscInt dof, off, p; 552839dd189SMatthew G Knepley 553839dd189SMatthew G Knepley ierr = PetscSectionGetDof(section, q, &dof);CHKERRQ(ierr); 554839dd189SMatthew G Knepley if (dof) { 555839dd189SMatthew G Knepley ierr = PetscFindInt(q, n, ind, &p);CHKERRQ(ierr); 556839dd189SMatthew G Knepley if (p >= pStart) { 557839dd189SMatthew G Knepley ierr = PetscSectionSetDof(newSection, p, dof);CHKERRQ(ierr); 558839dd189SMatthew G Knepley ierr = PetscSectionGetOffset(section, q, &off);CHKERRQ(ierr); 559839dd189SMatthew G Knepley ierr = PetscSectionSetOffset(newSection, p, off);CHKERRQ(ierr); 560839dd189SMatthew G Knepley } 561839dd189SMatthew G Knepley } 562839dd189SMatthew G Knepley } 563e6ccafaeSMatthew G Knepley ierr = ISRestoreIndices(subpointIS, &ind);CHKERRQ(ierr); 564e6ccafaeSMatthew G Knepley ierr = ISDestroy(&subpointIS);CHKERRQ(ierr); 565839dd189SMatthew G Knepley /* No need to setup section */ 566839dd189SMatthew G Knepley section = newSection; 567839dd189SMatthew G Knepley } 568552f7358SJed Brown } else { 569426ae2f1SMatthew G Knepley ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) §ion);CHKERRQ(ierr); 570426ae2f1SMatthew G Knepley if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name); 571552f7358SJed Brown } 57282f516ccSBarry Smith if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name); 573552f7358SJed Brown ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, &globalSection);CHKERRQ(ierr); 574552f7358SJed Brown ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr); 575552f7358SJed Brown ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr); 576839dd189SMatthew G Knepley if (newSection) {ierr = PetscSectionDestroy(&newSection);CHKERRQ(ierr);} 577552f7358SJed Brown } 578552f7358SJed Brown } 579552f7358SJed Brown /* Cell Fields */ 580*e8f6f0f6SMatthew G. Knepley ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_view_partition", &writePartition, NULL);CHKERRQ(ierr); 581*e8f6f0f6SMatthew G. Knepley if (hasCell || writePartition) { 58222d28d08SBarry Smith ierr = PetscFPrintf(comm, fp, "CELL_DATA %d\n", totCells);CHKERRQ(ierr); 583552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 584552f7358SJed Brown Vec X = (Vec) link->vec; 585552f7358SJed Brown DM dmX; 586552f7358SJed Brown PetscSection section, globalSection; 587552f7358SJed Brown const char *name; 588552f7358SJed Brown PetscInt enforceDof = PETSC_DETERMINE; 589552f7358SJed Brown 590552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 591552f7358SJed Brown if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) enforceDof = 3; 592552f7358SJed Brown ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr); 593552f7358SJed Brown ierr = VecGetDM(X, &dmX);CHKERRQ(ierr); 594552f7358SJed Brown if (dmX) { 59522d28d08SBarry Smith ierr = DMGetDefaultSection(dmX, §ion);CHKERRQ(ierr); 596552f7358SJed Brown } else { 597552f7358SJed Brown PetscContainer c; 598552f7358SJed Brown 599552f7358SJed Brown ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &c);CHKERRQ(ierr); 60082f516ccSBarry Smith if (!c) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name); 601552f7358SJed Brown ierr = PetscContainerGetPointer(c, (void**) §ion);CHKERRQ(ierr); 602552f7358SJed Brown } 60382f516ccSBarry Smith if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name); 604552f7358SJed Brown ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, &globalSection);CHKERRQ(ierr); 605552f7358SJed Brown ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr); 606552f7358SJed Brown ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr); 607552f7358SJed Brown } 608*e8f6f0f6SMatthew G. Knepley if (writePartition) { 609*e8f6f0f6SMatthew G. Knepley ierr = PetscFPrintf(comm, fp, "SCALARS partition int 1\n");CHKERRQ(ierr); 610*e8f6f0f6SMatthew G. Knepley ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr); 611*e8f6f0f6SMatthew G. Knepley ierr = DMPlexVTKWritePartition_ASCII(dm, fp);CHKERRQ(ierr); 612*e8f6f0f6SMatthew G. Knepley } 613552f7358SJed Brown } 614552f7358SJed Brown /* Cleanup */ 615552f7358SJed Brown ierr = PetscSectionDestroy(&globalCoordSection);CHKERRQ(ierr); 616552f7358SJed Brown ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr); 617552f7358SJed Brown ierr = PetscFClose(comm, fp);CHKERRQ(ierr); 618552f7358SJed Brown PetscFunctionReturn(0); 619552f7358SJed Brown } 620552f7358SJed Brown 621552f7358SJed Brown #undef __FUNCT__ 622552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll" 623552f7358SJed Brown /*@C 624552f7358SJed Brown DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer 625552f7358SJed Brown 626552f7358SJed Brown Collective 627552f7358SJed Brown 628552f7358SJed Brown Input Arguments: 629552f7358SJed Brown + odm - The DMPlex specifying the mesh, passed as a PetscObject 630552f7358SJed Brown - viewer - viewer of type VTK 631552f7358SJed Brown 632552f7358SJed Brown Level: developer 633552f7358SJed Brown 634552f7358SJed Brown Note: 635552f7358SJed Brown This function is a callback used by the VTK viewer to actually write the file. 636552f7358SJed 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. 637552f7358SJed Brown Instead, metadata for the entire file needs to be available up-front before you can start writing the file. 638552f7358SJed Brown 639552f7358SJed Brown .seealso: PETSCVIEWERVTK 640552f7358SJed Brown @*/ 641552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer) 642552f7358SJed Brown { 643552f7358SJed Brown DM dm = (DM) odm; 644552f7358SJed Brown PetscBool isvtk; 645552f7358SJed Brown PetscErrorCode ierr; 646552f7358SJed Brown 647552f7358SJed Brown PetscFunctionBegin; 648552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 649552f7358SJed Brown PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 650552f7358SJed Brown ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr); 65182f516ccSBarry Smith if (!isvtk) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name); 652552f7358SJed Brown switch (viewer->format) { 653552f7358SJed Brown case PETSC_VIEWER_ASCII_VTK: 654552f7358SJed Brown ierr = DMPlexVTKWriteAll_ASCII(dm, viewer);CHKERRQ(ierr); 655552f7358SJed Brown break; 656552f7358SJed Brown case PETSC_VIEWER_VTK_VTU: 657552f7358SJed Brown ierr = DMPlexVTKWriteAll_VTU(dm, viewer);CHKERRQ(ierr); 658552f7358SJed Brown break; 65982f516ccSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]); 660552f7358SJed Brown } 661552f7358SJed Brown PetscFunctionReturn(0); 662552f7358SJed Brown } 663