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; 54552f7358SJed Brown case 8: 55552f7358SJed Brown *cellType = 12; /* VTK_HEXAHEDRON */ 56552f7358SJed Brown break; 57552f7358SJed Brown case 10: 58552f7358SJed Brown *cellType = 24; /* VTK_QUADRATIC_TETRA */ 59552f7358SJed Brown break; 60552f7358SJed Brown case 27: 61552f7358SJed Brown *cellType = 29; /* VTK_QUADRATIC_HEXAHEDRON */ 62552f7358SJed Brown break; 63552f7358SJed Brown default: 64552f7358SJed Brown break; 65552f7358SJed Brown } 66552f7358SJed Brown } 67552f7358SJed Brown PetscFunctionReturn(0); 68552f7358SJed Brown } 69552f7358SJed Brown 707508e736SMatthew G. Knepley static PetscErrorCode DMPlexVTKWriteCells_ASCII(DM dm, FILE *fp, PetscInt *totalCells) 71552f7358SJed Brown { 7282f516ccSBarry Smith MPI_Comm comm; 7302281ff3SMatthew G. Knepley DMLabel label; 7402281ff3SMatthew G. Knepley IS globalVertexNumbers = NULL; 75552f7358SJed Brown const PetscInt *gvertex; 76552f7358SJed Brown PetscInt dim; 77552f7358SJed Brown PetscInt numCorners = 0, totCorners = 0, maxCorners, *corners; 78552f7358SJed Brown PetscInt numCells = 0, totCells = 0, maxCells, cellHeight; 79452b7979SMatthew G. Knepley PetscInt numLabelCells, maxLabelCells, cMax, cStart, cEnd, c, vStart, vEnd, v; 809852e123SBarry Smith PetscMPIInt size, rank, proc, tag; 81552f7358SJed Brown PetscErrorCode ierr; 82552f7358SJed Brown 83552f7358SJed Brown PetscFunctionBegin; 8482f516ccSBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 85552f7358SJed Brown ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr); 869852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 87552f7358SJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 88c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 89552f7358SJed Brown ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 90552f7358SJed Brown ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 91552f7358SJed Brown ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 920298fd71SBarry Smith ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 938865f1eaSKarl Rupp if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 94c58f1c22SToby Isaac ierr = DMGetLabel(dm, "vtk", &label);CHKERRQ(ierr); 95c58f1c22SToby Isaac ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 96b2566f29SBarry Smith ierr = MPIU_Allreduce(&numLabelCells, &maxLabelCells, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr); 97452b7979SMatthew G. Knepley if (!maxLabelCells) label = NULL; 98552f7358SJed Brown for (c = cStart; c < cEnd; ++c) { 990298fd71SBarry Smith PetscInt *closure = NULL; 10002281ff3SMatthew G. Knepley PetscInt closureSize, value; 101552f7358SJed Brown 10202281ff3SMatthew G. Knepley if (label) { 10302281ff3SMatthew G. Knepley ierr = DMLabelGetValue(label, c, &value);CHKERRQ(ierr); 104552f7358SJed Brown if (value != 1) continue; 105552f7358SJed Brown } 106552f7358SJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 107552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 108552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) ++numCorners; 109552f7358SJed Brown } 110552f7358SJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 111552f7358SJed Brown ++numCells; 112552f7358SJed Brown } 113552f7358SJed Brown maxCells = numCells; 114552f7358SJed Brown ierr = MPI_Reduce(&numCells, &totCells, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRQ(ierr); 115552f7358SJed Brown ierr = MPI_Reduce(&numCells, &maxCells, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr); 116552f7358SJed Brown ierr = MPI_Reduce(&numCorners, &totCorners, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRQ(ierr); 117552f7358SJed Brown ierr = MPI_Reduce(&numCorners, &maxCorners, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr); 118552f7358SJed Brown ierr = DMPlexGetVertexNumbering(dm, &globalVertexNumbers);CHKERRQ(ierr); 119552f7358SJed Brown ierr = ISGetIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr); 120785e854fSJed Brown ierr = PetscMalloc1(maxCells, &corners);CHKERRQ(ierr); 1216b8858a1SBarry Smith ierr = PetscFPrintf(comm, fp, "CELLS %D %D\n", totCells, totCorners+totCells);CHKERRQ(ierr); 122552f7358SJed Brown if (!rank) { 123552f7358SJed Brown PetscInt *remoteVertices; 1240dcf9c70SMatthew G. Knepley int *vertices; 125552f7358SJed Brown 126785e854fSJed Brown ierr = PetscMalloc1(maxCorners, &vertices);CHKERRQ(ierr); 127552f7358SJed Brown for (c = cStart, numCells = 0; c < cEnd; ++c) { 1280298fd71SBarry Smith PetscInt *closure = NULL; 12902281ff3SMatthew G. Knepley PetscInt closureSize, value, nC = 0; 130552f7358SJed Brown 13102281ff3SMatthew G. Knepley if (label) { 13202281ff3SMatthew G. Knepley ierr = DMLabelGetValue(label, c, &value);CHKERRQ(ierr); 133552f7358SJed Brown if (value != 1) continue; 134552f7358SJed Brown } 135552f7358SJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 136552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 137552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 138552f7358SJed Brown const PetscInt gv = gvertex[closure[v] - vStart]; 1390dcf9c70SMatthew G. Knepley vertices[nC++] = gv < 0 ? -(gv+1) : gv; 140552f7358SJed Brown } 141552f7358SJed Brown } 1420dcf9c70SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 143552f7358SJed Brown corners[numCells++] = nC; 1446b8858a1SBarry Smith ierr = PetscFPrintf(comm, fp, "%D ", nC);CHKERRQ(ierr); 1450dcf9c70SMatthew G. Knepley ierr = DMPlexInvertCell(dim, nC, vertices);CHKERRQ(ierr); 146552f7358SJed Brown for (v = 0; v < nC; ++v) { 1476b8858a1SBarry Smith ierr = PetscFPrintf(comm, fp, " %D", vertices[v]);CHKERRQ(ierr); 148552f7358SJed Brown } 149552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr); 150552f7358SJed Brown } 1519852e123SBarry Smith if (size > 1) {ierr = PetscMalloc1(maxCorners+maxCells, &remoteVertices);CHKERRQ(ierr);} 1529852e123SBarry Smith for (proc = 1; proc < size; ++proc) { 153552f7358SJed Brown MPI_Status status; 154552f7358SJed Brown 155552f7358SJed Brown ierr = MPI_Recv(&numCorners, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 156552f7358SJed Brown ierr = MPI_Recv(remoteVertices, numCorners, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 157552f7358SJed Brown for (c = 0; c < numCorners;) { 158552f7358SJed Brown PetscInt nC = remoteVertices[c++]; 159552f7358SJed Brown 160552f7358SJed Brown for (v = 0; v < nC; ++v, ++c) { 1610dcf9c70SMatthew G. Knepley vertices[v] = remoteVertices[c]; 1620dcf9c70SMatthew G. Knepley } 1630dcf9c70SMatthew G. Knepley ierr = DMPlexInvertCell(dim, nC, vertices);CHKERRQ(ierr); 1646b8858a1SBarry Smith ierr = PetscFPrintf(comm, fp, "%D ", nC);CHKERRQ(ierr); 1650dcf9c70SMatthew G. Knepley for (v = 0; v < nC; ++v) { 1666b8858a1SBarry Smith ierr = PetscFPrintf(comm, fp, " %D", vertices[v]);CHKERRQ(ierr); 167552f7358SJed Brown } 168552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr); 169552f7358SJed Brown } 170552f7358SJed Brown } 1719852e123SBarry Smith if (size > 1) {ierr = PetscFree(remoteVertices);CHKERRQ(ierr);} 1720dcf9c70SMatthew G. Knepley ierr = PetscFree(vertices);CHKERRQ(ierr); 173552f7358SJed Brown } else { 174552f7358SJed Brown PetscInt *localVertices, numSend = numCells+numCorners, k = 0; 175552f7358SJed Brown 176785e854fSJed Brown ierr = PetscMalloc1(numSend, &localVertices);CHKERRQ(ierr); 177552f7358SJed Brown for (c = cStart, numCells = 0; c < cEnd; ++c) { 1780298fd71SBarry Smith PetscInt *closure = NULL; 17902281ff3SMatthew G. Knepley PetscInt closureSize, value, nC = 0; 180552f7358SJed Brown 18102281ff3SMatthew G. Knepley if (label) { 18202281ff3SMatthew G. Knepley ierr = DMLabelGetValue(label, c, &value);CHKERRQ(ierr); 183552f7358SJed Brown if (value != 1) continue; 184552f7358SJed Brown } 185552f7358SJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 186552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 187552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 188552f7358SJed Brown const PetscInt gv = gvertex[closure[v] - vStart]; 189552f7358SJed Brown closure[nC++] = gv < 0 ? -(gv+1) : gv; 190552f7358SJed Brown } 191552f7358SJed Brown } 192552f7358SJed Brown corners[numCells++] = nC; 193552f7358SJed Brown localVertices[k++] = nC; 194552f7358SJed Brown for (v = 0; v < nC; ++v, ++k) { 195552f7358SJed Brown localVertices[k] = closure[v]; 196552f7358SJed Brown } 197552f7358SJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 198552f7358SJed Brown } 1996b8858a1SBarry Smith if (k != numSend) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid number of vertices to send %D should be %D", k, numSend); 200552f7358SJed Brown ierr = MPI_Send(&numSend, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 201552f7358SJed Brown ierr = MPI_Send(localVertices, numSend, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 202552f7358SJed Brown ierr = PetscFree(localVertices);CHKERRQ(ierr); 203552f7358SJed Brown } 204552f7358SJed Brown ierr = ISRestoreIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr); 2056b8858a1SBarry Smith ierr = PetscFPrintf(comm, fp, "CELL_TYPES %D\n", totCells);CHKERRQ(ierr); 206552f7358SJed Brown if (!rank) { 207552f7358SJed Brown PetscInt cellType; 208552f7358SJed Brown 209552f7358SJed Brown for (c = 0; c < numCells; ++c) { 210de0a4605SMatthew G. Knepley ierr = DMPlexVTKGetCellType_Internal(dm, dim, corners[c], &cellType);CHKERRQ(ierr); 2116b8858a1SBarry Smith ierr = PetscFPrintf(comm, fp, "%D\n", cellType);CHKERRQ(ierr); 212552f7358SJed Brown } 2139852e123SBarry Smith for (proc = 1; proc < size; ++proc) { 214552f7358SJed Brown MPI_Status status; 215552f7358SJed Brown 216552f7358SJed Brown ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 217552f7358SJed Brown ierr = MPI_Recv(corners, numCells, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 218552f7358SJed Brown for (c = 0; c < numCells; ++c) { 219de0a4605SMatthew G. Knepley ierr = DMPlexVTKGetCellType_Internal(dm, dim, corners[c], &cellType);CHKERRQ(ierr); 2206b8858a1SBarry Smith ierr = PetscFPrintf(comm, fp, "%D\n", cellType);CHKERRQ(ierr); 221552f7358SJed Brown } 222552f7358SJed Brown } 223552f7358SJed Brown } else { 224552f7358SJed Brown ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 225552f7358SJed Brown ierr = MPI_Send(corners, numCells, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 226552f7358SJed Brown } 227552f7358SJed Brown ierr = PetscFree(corners);CHKERRQ(ierr); 228552f7358SJed Brown *totalCells = totCells; 229552f7358SJed Brown PetscFunctionReturn(0); 230552f7358SJed Brown } 231552f7358SJed Brown 2327508e736SMatthew G. Knepley static PetscErrorCode DMPlexVTKWritePartition_ASCII(DM dm, FILE *fp) 233e8f6f0f6SMatthew G. Knepley { 234e8f6f0f6SMatthew G. Knepley MPI_Comm comm; 235e8f6f0f6SMatthew G. Knepley PetscInt numCells = 0, cellHeight; 236e8f6f0f6SMatthew G. Knepley PetscInt numLabelCells, cMax, cStart, cEnd, c; 2379852e123SBarry Smith PetscMPIInt size, rank, proc, tag; 238e8f6f0f6SMatthew G. Knepley PetscBool hasLabel; 239e8f6f0f6SMatthew G. Knepley PetscErrorCode ierr; 240e8f6f0f6SMatthew G. Knepley 241e8f6f0f6SMatthew G. Knepley PetscFunctionBegin; 242e8f6f0f6SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 243e8f6f0f6SMatthew G. Knepley ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr); 2449852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 245e8f6f0f6SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 246e8f6f0f6SMatthew G. Knepley ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 247e8f6f0f6SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 248e8f6f0f6SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 249e8f6f0f6SMatthew G. Knepley if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 250c58f1c22SToby Isaac ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 251e8f6f0f6SMatthew G. Knepley hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 252e8f6f0f6SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 253e8f6f0f6SMatthew G. Knepley if (hasLabel) { 254e8f6f0f6SMatthew G. Knepley PetscInt value; 255e8f6f0f6SMatthew G. Knepley 256c58f1c22SToby Isaac ierr = DMGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 257e8f6f0f6SMatthew G. Knepley if (value != 1) continue; 258e8f6f0f6SMatthew G. Knepley } 259e8f6f0f6SMatthew G. Knepley ++numCells; 260e8f6f0f6SMatthew G. Knepley } 261e8f6f0f6SMatthew G. Knepley if (!rank) { 262e8f6f0f6SMatthew G. Knepley for (c = 0; c < numCells; ++c) {ierr = PetscFPrintf(comm, fp, "%d\n", rank);CHKERRQ(ierr);} 2639852e123SBarry Smith for (proc = 1; proc < size; ++proc) { 264e8f6f0f6SMatthew G. Knepley MPI_Status status; 265e8f6f0f6SMatthew G. Knepley 266e8f6f0f6SMatthew G. Knepley ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 267e8f6f0f6SMatthew G. Knepley for (c = 0; c < numCells; ++c) {ierr = PetscFPrintf(comm, fp, "%d\n", proc);CHKERRQ(ierr);} 268e8f6f0f6SMatthew G. Knepley } 269e8f6f0f6SMatthew G. Knepley } else { 270e8f6f0f6SMatthew G. Knepley ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 271e8f6f0f6SMatthew G. Knepley } 272e8f6f0f6SMatthew G. Knepley PetscFunctionReturn(0); 273e8f6f0f6SMatthew G. Knepley } 274e8f6f0f6SMatthew G. Knepley 2757508e736SMatthew G. Knepley static PetscErrorCode DMPlexVTKWriteSection_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec v, FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale) 2760adebc6cSBarry Smith { 27782f516ccSBarry Smith MPI_Comm comm; 278552f7358SJed Brown const MPI_Datatype mpiType = MPIU_SCALAR; 279552f7358SJed Brown PetscScalar *array; 280552f7358SJed Brown PetscInt numDof = 0, maxDof; 281552f7358SJed Brown PetscInt numLabelCells, cellHeight, cMax, cStart, cEnd, numLabelVertices, vMax, vStart, vEnd, pStart, pEnd, p; 2829852e123SBarry Smith PetscMPIInt size, rank, proc, tag; 283552f7358SJed Brown PetscBool hasLabel; 284552f7358SJed Brown PetscErrorCode ierr; 285552f7358SJed Brown 286552f7358SJed Brown PetscFunctionBegin; 28782f516ccSBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 288552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 289552f7358SJed Brown PetscValidHeaderSpecific(v,VEC_CLASSID,4); 290552f7358SJed Brown if (precision < 0) precision = 6; 291552f7358SJed Brown ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr); 2929852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 293552f7358SJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 294552f7358SJed Brown ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 295552f7358SJed Brown /* VTK only wants the values at cells or vertices */ 296552f7358SJed Brown ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 297552f7358SJed Brown ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 298552f7358SJed Brown ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 2990298fd71SBarry Smith ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, &vMax);CHKERRQ(ierr); 3008865f1eaSKarl Rupp if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 3018865f1eaSKarl Rupp if (vMax >= 0) vEnd = PetscMin(vEnd, vMax); 302552f7358SJed Brown pStart = PetscMax(PetscMin(cStart, vStart), pStart); 303552f7358SJed Brown pEnd = PetscMin(PetscMax(cEnd, vEnd), pEnd); 304c58f1c22SToby Isaac ierr = DMGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 305c58f1c22SToby Isaac ierr = DMGetStratumSize(dm, "vtk", 2, &numLabelVertices);CHKERRQ(ierr); 306552f7358SJed Brown hasLabel = numLabelCells > 0 || numLabelVertices > 0 ? PETSC_TRUE : PETSC_FALSE; 307552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 308552f7358SJed Brown /* Reject points not either cells or vertices */ 309552f7358SJed Brown if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue; 310552f7358SJed Brown if (hasLabel) { 311552f7358SJed Brown PetscInt value; 312552f7358SJed Brown 313552f7358SJed Brown if (((p >= cStart) && (p < cEnd) && numLabelCells) || 314552f7358SJed Brown ((p >= vStart) && (p < vEnd) && numLabelVertices)) { 315c58f1c22SToby Isaac ierr = DMGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr); 316552f7358SJed Brown if (value != 1) continue; 317552f7358SJed Brown } 318552f7358SJed Brown } 319552f7358SJed Brown ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr); 3208865f1eaSKarl Rupp if (numDof) break; 321552f7358SJed Brown } 322b2566f29SBarry Smith ierr = MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr); 323552f7358SJed Brown enforceDof = PetscMax(enforceDof, maxDof); 324552f7358SJed Brown ierr = VecGetArray(v, &array);CHKERRQ(ierr); 325552f7358SJed Brown if (!rank) { 326552f7358SJed Brown char formatString[8]; 327552f7358SJed Brown 328552f7358SJed Brown ierr = PetscSNPrintf(formatString, 8, "%%.%de", precision);CHKERRQ(ierr); 329552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 330552f7358SJed Brown /* Here we lose a way to filter points by keeping them out of the Numbering */ 331552f7358SJed Brown PetscInt dof, off, goff, d; 332552f7358SJed Brown 333552f7358SJed Brown /* Reject points not either cells or vertices */ 334552f7358SJed Brown if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue; 335552f7358SJed Brown if (hasLabel) { 336552f7358SJed Brown PetscInt value; 337552f7358SJed Brown 338552f7358SJed Brown if (((p >= cStart) && (p < cEnd) && numLabelCells) || 339552f7358SJed Brown ((p >= vStart) && (p < vEnd) && numLabelVertices)) { 340c58f1c22SToby Isaac ierr = DMGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr); 341552f7358SJed Brown if (value != 1) continue; 342552f7358SJed Brown } 343552f7358SJed Brown } 344552f7358SJed Brown ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 345552f7358SJed Brown ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 346552f7358SJed Brown ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr); 347552f7358SJed Brown if (dof && goff >= 0) { 348552f7358SJed Brown for (d = 0; d < dof; d++) { 349552f7358SJed Brown if (d > 0) { 350552f7358SJed Brown ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr); 351552f7358SJed Brown } 352552f7358SJed Brown ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(array[off+d])*scale);CHKERRQ(ierr); 353552f7358SJed Brown } 354552f7358SJed Brown for (d = dof; d < enforceDof; d++) { 355552f7358SJed Brown ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr); 356552f7358SJed Brown } 357552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr); 358552f7358SJed Brown } 359552f7358SJed Brown } 3609852e123SBarry Smith for (proc = 1; proc < size; ++proc) { 361552f7358SJed Brown PetscScalar *remoteValues; 362d892089bSMatthew G. Knepley PetscInt size = 0, d; 363552f7358SJed Brown MPI_Status status; 364552f7358SJed Brown 365552f7358SJed Brown ierr = MPI_Recv(&size, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 366785e854fSJed Brown ierr = PetscMalloc1(size, &remoteValues);CHKERRQ(ierr); 367552f7358SJed Brown ierr = MPI_Recv(remoteValues, size, mpiType, proc, tag, comm, &status);CHKERRQ(ierr); 368552f7358SJed Brown for (p = 0; p < size/maxDof; ++p) { 369552f7358SJed Brown for (d = 0; d < maxDof; ++d) { 370552f7358SJed Brown if (d > 0) { 371552f7358SJed Brown ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr); 372552f7358SJed Brown } 373552f7358SJed Brown ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(remoteValues[p*maxDof+d])*scale);CHKERRQ(ierr); 374552f7358SJed Brown } 375552f7358SJed Brown for (d = maxDof; d < enforceDof; ++d) { 376552f7358SJed Brown ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr); 377552f7358SJed Brown } 378552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr); 379552f7358SJed Brown } 380552f7358SJed Brown ierr = PetscFree(remoteValues);CHKERRQ(ierr); 381552f7358SJed Brown } 382552f7358SJed Brown } else { 383552f7358SJed Brown PetscScalar *localValues; 384552f7358SJed Brown PetscInt size, k = 0; 385552f7358SJed Brown 386552f7358SJed Brown ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr); 387785e854fSJed Brown ierr = PetscMalloc1(size, &localValues);CHKERRQ(ierr); 388552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 389552f7358SJed Brown PetscInt dof, off, goff, d; 390552f7358SJed Brown 391552f7358SJed Brown /* Reject points not either cells or vertices */ 392552f7358SJed Brown if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue; 393552f7358SJed Brown if (hasLabel) { 394552f7358SJed Brown PetscInt value; 395552f7358SJed Brown 396552f7358SJed Brown if (((p >= cStart) && (p < cEnd) && numLabelCells) || 397552f7358SJed Brown ((p >= vStart) && (p < vEnd) && numLabelVertices)) { 398c58f1c22SToby Isaac ierr = DMGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr); 399552f7358SJed Brown if (value != 1) continue; 400552f7358SJed Brown } 401552f7358SJed Brown } 402552f7358SJed Brown ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 403552f7358SJed Brown ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 404552f7358SJed Brown ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr); 405552f7358SJed Brown if (goff >= 0) { 406552f7358SJed Brown for (d = 0; d < dof; ++d) { 407552f7358SJed Brown localValues[k++] = array[off+d]; 408552f7358SJed Brown } 409552f7358SJed Brown } 410552f7358SJed Brown } 411552f7358SJed Brown ierr = MPI_Send(&k, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 412552f7358SJed Brown ierr = MPI_Send(localValues, k, mpiType, 0, tag, comm);CHKERRQ(ierr); 413552f7358SJed Brown ierr = PetscFree(localValues);CHKERRQ(ierr); 414552f7358SJed Brown } 415552f7358SJed Brown ierr = VecRestoreArray(v, &array);CHKERRQ(ierr); 416552f7358SJed Brown PetscFunctionReturn(0); 417552f7358SJed Brown } 418552f7358SJed Brown 4197508e736SMatthew G. Knepley static PetscErrorCode DMPlexVTKWriteField_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec field, const char name[], FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale) 420552f7358SJed Brown { 42182f516ccSBarry Smith MPI_Comm comm; 422552f7358SJed Brown PetscInt numDof = 0, maxDof; 423552f7358SJed Brown PetscInt pStart, pEnd, p; 424552f7358SJed Brown PetscErrorCode ierr; 425552f7358SJed Brown 426552f7358SJed Brown PetscFunctionBegin; 42782f516ccSBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 428552f7358SJed Brown ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 429552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 430552f7358SJed Brown ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr); 4318865f1eaSKarl Rupp if (numDof) break; 432552f7358SJed Brown } 433552f7358SJed Brown numDof = PetscMax(numDof, enforceDof); 434b2566f29SBarry Smith ierr = MPIU_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 435552f7358SJed Brown if (!name) name = "Unknown"; 436552f7358SJed Brown if (maxDof == 3) { 437552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "VECTORS %s double\n", name);CHKERRQ(ierr); 438552f7358SJed Brown } else { 4396b8858a1SBarry Smith ierr = PetscFPrintf(comm, fp, "SCALARS %s double %D\n", name, maxDof);CHKERRQ(ierr); 440552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr); 441552f7358SJed Brown } 442552f7358SJed Brown ierr = DMPlexVTKWriteSection_ASCII(dm, section, globalSection, field, fp, enforceDof, precision, scale);CHKERRQ(ierr); 443552f7358SJed Brown PetscFunctionReturn(0); 444552f7358SJed Brown } 445552f7358SJed Brown 446552f7358SJed Brown static PetscErrorCode DMPlexVTKWriteAll_ASCII(DM dm, PetscViewer viewer) 447552f7358SJed Brown { 44882f516ccSBarry Smith MPI_Comm comm; 449552f7358SJed Brown PetscViewer_VTK *vtk = (PetscViewer_VTK*) viewer->data; 450552f7358SJed Brown FILE *fp; 451552f7358SJed Brown PetscViewerVTKObjectLink link; 452552f7358SJed Brown PetscSection coordSection, globalCoordSection; 453552f7358SJed Brown PetscLayout vLayout; 454552f7358SJed Brown Vec coordinates; 455552f7358SJed Brown PetscReal lengthScale; 456b6f3c738SMatthew G. Knepley PetscInt vMax, totVertices, totCells = 0; 457*d72c04d6SStefano Zampini PetscBool hasPoint = PETSC_FALSE, hasCell = PETSC_FALSE, writePartition = PETSC_FALSE, localized; 458552f7358SJed Brown PetscErrorCode ierr; 459552f7358SJed Brown 460552f7358SJed Brown PetscFunctionBegin; 461*d72c04d6SStefano Zampini ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr); 462*d72c04d6SStefano Zampini if (localized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"VTK output with localized coordinates not yet supported"); 46382f516ccSBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 464552f7358SJed Brown ierr = PetscFOpen(comm, vtk->filename, "wb", &fp);CHKERRQ(ierr); 465552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "# vtk DataFile Version 2.0\n");CHKERRQ(ierr); 466552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "Simplicial Mesh Example\n");CHKERRQ(ierr); 467552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "ASCII\n");CHKERRQ(ierr); 468552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "DATASET UNSTRUCTURED_GRID\n");CHKERRQ(ierr); 469552f7358SJed Brown /* Vertices */ 470552f7358SJed Brown ierr = DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);CHKERRQ(ierr); 47169d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 47215b58121SMatthew G. Knepley ierr = PetscSectionCreateGlobalSection(coordSection, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalCoordSection);CHKERRQ(ierr); 473552f7358SJed Brown ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 4740298fd71SBarry Smith ierr = DMPlexGetHybridBounds(dm, NULL, NULL, NULL, &vMax);CHKERRQ(ierr); 475552f7358SJed Brown if (vMax >= 0) { 476552f7358SJed Brown PetscInt pStart, pEnd, p, localSize = 0; 477552f7358SJed Brown 478552f7358SJed Brown ierr = PetscSectionGetChart(globalCoordSection, &pStart, &pEnd);CHKERRQ(ierr); 479552f7358SJed Brown pEnd = PetscMin(pEnd, vMax); 480552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 481552f7358SJed Brown PetscInt dof; 482552f7358SJed Brown 483552f7358SJed Brown ierr = PetscSectionGetDof(globalCoordSection, p, &dof);CHKERRQ(ierr); 4848865f1eaSKarl Rupp if (dof > 0) ++localSize; 485552f7358SJed Brown } 48682f516ccSBarry Smith ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &vLayout);CHKERRQ(ierr); 487552f7358SJed Brown ierr = PetscLayoutSetLocalSize(vLayout, localSize);CHKERRQ(ierr); 488552f7358SJed Brown ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr); 489552f7358SJed Brown ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr); 490552f7358SJed Brown } else { 49182f516ccSBarry Smith ierr = PetscSectionGetPointLayout(PetscObjectComm((PetscObject)dm), globalCoordSection, &vLayout);CHKERRQ(ierr); 492552f7358SJed Brown } 493552f7358SJed Brown ierr = PetscLayoutGetSize(vLayout, &totVertices);CHKERRQ(ierr); 4946b8858a1SBarry Smith ierr = PetscFPrintf(comm, fp, "POINTS %D double\n", totVertices);CHKERRQ(ierr); 495552f7358SJed Brown ierr = DMPlexVTKWriteSection_ASCII(dm, coordSection, globalCoordSection, coordinates, fp, 3, PETSC_DETERMINE, lengthScale);CHKERRQ(ierr); 496552f7358SJed Brown /* Cells */ 497552f7358SJed Brown ierr = DMPlexVTKWriteCells_ASCII(dm, fp, &totCells);CHKERRQ(ierr); 498552f7358SJed Brown /* Vertex fields */ 499552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 500552f7358SJed Brown if ((link->ft == PETSC_VTK_POINT_FIELD) || (link->ft == PETSC_VTK_POINT_VECTOR_FIELD)) hasPoint = PETSC_TRUE; 501552f7358SJed Brown if ((link->ft == PETSC_VTK_CELL_FIELD) || (link->ft == PETSC_VTK_CELL_VECTOR_FIELD)) hasCell = PETSC_TRUE; 502552f7358SJed Brown } 503552f7358SJed Brown if (hasPoint) { 5046b8858a1SBarry Smith ierr = PetscFPrintf(comm, fp, "POINT_DATA %D\n", totVertices);CHKERRQ(ierr); 505552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 506552f7358SJed Brown Vec X = (Vec) link->vec; 507552f7358SJed Brown DM dmX; 5080298fd71SBarry Smith PetscSection section, globalSection, newSection = NULL; 509552f7358SJed Brown const char *name; 510552f7358SJed Brown PetscInt enforceDof = PETSC_DETERMINE; 511552f7358SJed Brown 512552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 513552f7358SJed Brown if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) enforceDof = 3; 514552f7358SJed Brown ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr); 515552f7358SJed Brown ierr = VecGetDM(X, &dmX);CHKERRQ(ierr); 516552f7358SJed Brown if (dmX) { 517088580cfSMatthew G Knepley DMLabel subpointMap, subpointMapX; 518839dd189SMatthew G Knepley PetscInt dim, dimX, pStart, pEnd, qStart, qEnd; 519839dd189SMatthew G Knepley 52022d28d08SBarry Smith ierr = DMGetDefaultSection(dmX, §ion);CHKERRQ(ierr); 521839dd189SMatthew G Knepley /* Here is where we check whether dmX is a submesh of dm */ 522c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 523c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dmX, &dimX);CHKERRQ(ierr); 524839dd189SMatthew G Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 525839dd189SMatthew G Knepley ierr = DMPlexGetChart(dmX, &qStart, &qEnd);CHKERRQ(ierr); 526839dd189SMatthew G Knepley ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr); 527839dd189SMatthew G Knepley ierr = DMPlexGetSubpointMap(dmX, &subpointMapX);CHKERRQ(ierr); 528839dd189SMatthew G Knepley if (((dim != dimX) || ((pEnd-pStart) < (qEnd-qStart))) && subpointMap && !subpointMapX) { 529434e42b5SMatthew G. Knepley const PetscInt *ind = NULL; 530e6ccafaeSMatthew G Knepley IS subpointIS; 531434e42b5SMatthew G. Knepley PetscInt n = 0, q; 532839dd189SMatthew G Knepley 533839dd189SMatthew G Knepley ierr = PetscSectionGetChart(section, &qStart, &qEnd);CHKERRQ(ierr); 534e6ccafaeSMatthew G Knepley ierr = DMPlexCreateSubpointIS(dm, &subpointIS);CHKERRQ(ierr); 535434e42b5SMatthew G. Knepley if (subpointIS) { 536e6ccafaeSMatthew G Knepley ierr = ISGetLocalSize(subpointIS, &n);CHKERRQ(ierr); 537e6ccafaeSMatthew G Knepley ierr = ISGetIndices(subpointIS, &ind);CHKERRQ(ierr); 538434e42b5SMatthew G. Knepley } 539839dd189SMatthew G Knepley ierr = PetscSectionCreate(comm, &newSection);CHKERRQ(ierr); 540839dd189SMatthew G Knepley ierr = PetscSectionSetChart(newSection, pStart, pEnd);CHKERRQ(ierr); 541839dd189SMatthew G Knepley for (q = qStart; q < qEnd; ++q) { 542839dd189SMatthew G Knepley PetscInt dof, off, p; 543839dd189SMatthew G Knepley 544839dd189SMatthew G Knepley ierr = PetscSectionGetDof(section, q, &dof);CHKERRQ(ierr); 545839dd189SMatthew G Knepley if (dof) { 546839dd189SMatthew G Knepley ierr = PetscFindInt(q, n, ind, &p);CHKERRQ(ierr); 547839dd189SMatthew G Knepley if (p >= pStart) { 548839dd189SMatthew G Knepley ierr = PetscSectionSetDof(newSection, p, dof);CHKERRQ(ierr); 549839dd189SMatthew G Knepley ierr = PetscSectionGetOffset(section, q, &off);CHKERRQ(ierr); 550839dd189SMatthew G Knepley ierr = PetscSectionSetOffset(newSection, p, off);CHKERRQ(ierr); 551839dd189SMatthew G Knepley } 552839dd189SMatthew G Knepley } 553839dd189SMatthew G Knepley } 554434e42b5SMatthew G. Knepley if (subpointIS) { 555e6ccafaeSMatthew G Knepley ierr = ISRestoreIndices(subpointIS, &ind);CHKERRQ(ierr); 556e6ccafaeSMatthew G Knepley ierr = ISDestroy(&subpointIS);CHKERRQ(ierr); 557434e42b5SMatthew G. Knepley } 558839dd189SMatthew G Knepley /* No need to setup section */ 559839dd189SMatthew G Knepley section = newSection; 560839dd189SMatthew G Knepley } 561552f7358SJed Brown } else { 562426ae2f1SMatthew G Knepley ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) §ion);CHKERRQ(ierr); 563426ae2f1SMatthew G Knepley if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name); 564552f7358SJed Brown } 56582f516ccSBarry Smith if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name); 56615b58121SMatthew G. Knepley ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection);CHKERRQ(ierr); 567552f7358SJed Brown ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr); 568552f7358SJed Brown ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr); 569839dd189SMatthew G Knepley if (newSection) {ierr = PetscSectionDestroy(&newSection);CHKERRQ(ierr);} 570552f7358SJed Brown } 571552f7358SJed Brown } 572552f7358SJed Brown /* Cell Fields */ 573c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_view_partition", &writePartition, NULL);CHKERRQ(ierr); 574e8f6f0f6SMatthew G. Knepley if (hasCell || writePartition) { 5756b8858a1SBarry Smith ierr = PetscFPrintf(comm, fp, "CELL_DATA %D\n", totCells);CHKERRQ(ierr); 576552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 577552f7358SJed Brown Vec X = (Vec) link->vec; 578552f7358SJed Brown DM dmX; 579552f7358SJed Brown PetscSection section, globalSection; 580552f7358SJed Brown const char *name; 581552f7358SJed Brown PetscInt enforceDof = PETSC_DETERMINE; 582552f7358SJed Brown 583552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 584552f7358SJed Brown if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) enforceDof = 3; 585552f7358SJed Brown ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr); 586552f7358SJed Brown ierr = VecGetDM(X, &dmX);CHKERRQ(ierr); 587552f7358SJed Brown if (dmX) { 58822d28d08SBarry Smith ierr = DMGetDefaultSection(dmX, §ion);CHKERRQ(ierr); 589552f7358SJed Brown } else { 590552f7358SJed Brown PetscContainer c; 591552f7358SJed Brown 592552f7358SJed Brown ierr = PetscObjectQuery(link->vec, "section", (PetscObject*) &c);CHKERRQ(ierr); 59382f516ccSBarry Smith if (!c) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name); 594552f7358SJed Brown ierr = PetscContainerGetPointer(c, (void**) §ion);CHKERRQ(ierr); 595552f7358SJed Brown } 59682f516ccSBarry Smith if (!section) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name); 59715b58121SMatthew G. Knepley ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, PETSC_FALSE, &globalSection);CHKERRQ(ierr); 598552f7358SJed Brown ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr); 599552f7358SJed Brown ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr); 600552f7358SJed Brown } 601e8f6f0f6SMatthew G. Knepley if (writePartition) { 602e8f6f0f6SMatthew G. Knepley ierr = PetscFPrintf(comm, fp, "SCALARS partition int 1\n");CHKERRQ(ierr); 603e8f6f0f6SMatthew G. Knepley ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr); 604e8f6f0f6SMatthew G. Knepley ierr = DMPlexVTKWritePartition_ASCII(dm, fp);CHKERRQ(ierr); 605e8f6f0f6SMatthew G. Knepley } 606552f7358SJed Brown } 607552f7358SJed Brown /* Cleanup */ 608552f7358SJed Brown ierr = PetscSectionDestroy(&globalCoordSection);CHKERRQ(ierr); 609552f7358SJed Brown ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr); 610552f7358SJed Brown ierr = PetscFClose(comm, fp);CHKERRQ(ierr); 611552f7358SJed Brown PetscFunctionReturn(0); 612552f7358SJed Brown } 613552f7358SJed Brown 614552f7358SJed Brown /*@C 615552f7358SJed Brown DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer 616552f7358SJed Brown 617552f7358SJed Brown Collective 618552f7358SJed Brown 619552f7358SJed Brown Input Arguments: 620552f7358SJed Brown + odm - The DMPlex specifying the mesh, passed as a PetscObject 621552f7358SJed Brown - viewer - viewer of type VTK 622552f7358SJed Brown 623552f7358SJed Brown Level: developer 624552f7358SJed Brown 625552f7358SJed Brown Note: 626552f7358SJed Brown This function is a callback used by the VTK viewer to actually write the file. 627552f7358SJed 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. 628552f7358SJed Brown Instead, metadata for the entire file needs to be available up-front before you can start writing the file. 629552f7358SJed Brown 630552f7358SJed Brown .seealso: PETSCVIEWERVTK 631552f7358SJed Brown @*/ 632552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer) 633552f7358SJed Brown { 634552f7358SJed Brown DM dm = (DM) odm; 635552f7358SJed Brown PetscBool isvtk; 636552f7358SJed Brown PetscErrorCode ierr; 637552f7358SJed Brown 638552f7358SJed Brown PetscFunctionBegin; 639552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 640552f7358SJed Brown PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 641552f7358SJed Brown ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr); 64282f516ccSBarry Smith if (!isvtk) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name); 643552f7358SJed Brown switch (viewer->format) { 644552f7358SJed Brown case PETSC_VIEWER_ASCII_VTK: 645552f7358SJed Brown ierr = DMPlexVTKWriteAll_ASCII(dm, viewer);CHKERRQ(ierr); 646552f7358SJed Brown break; 647552f7358SJed Brown case PETSC_VIEWER_VTK_VTU: 648552f7358SJed Brown ierr = DMPlexVTKWriteAll_VTU(dm, viewer);CHKERRQ(ierr); 649552f7358SJed Brown break; 65082f516ccSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]); 651552f7358SJed Brown } 652552f7358SJed Brown PetscFunctionReturn(0); 653552f7358SJed Brown } 654