1*552f7358SJed Brown #define PETSCDM_DLL 2*552f7358SJed Brown #include <petsc-private/pleximpl.h> /*I "petscdmplex.h" I*/ 3*552f7358SJed Brown #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h> 4*552f7358SJed Brown 5*552f7358SJed Brown #undef __FUNCT__ 6*552f7358SJed Brown #define __FUNCT__ "DMPlexVTKGetCellType" 7*552f7358SJed Brown PetscErrorCode DMPlexVTKGetCellType(DM dm, PetscInt dim, PetscInt corners, PetscInt *cellType) { 8*552f7358SJed Brown PetscFunctionBegin; 9*552f7358SJed Brown *cellType = -1; 10*552f7358SJed Brown switch(dim) { 11*552f7358SJed Brown case 0: 12*552f7358SJed Brown switch(corners) { 13*552f7358SJed Brown case 1: 14*552f7358SJed Brown *cellType = 1; /* VTK_VERTEX */ 15*552f7358SJed Brown break; 16*552f7358SJed Brown default: 17*552f7358SJed Brown break; 18*552f7358SJed Brown } 19*552f7358SJed Brown break; 20*552f7358SJed Brown case 1: 21*552f7358SJed Brown switch(corners) { 22*552f7358SJed Brown case 2: 23*552f7358SJed Brown *cellType = 3; /* VTK_LINE */ 24*552f7358SJed Brown break; 25*552f7358SJed Brown case 3: 26*552f7358SJed Brown *cellType = 21; /* VTK_QUADRATIC_EDGE */ 27*552f7358SJed Brown break; 28*552f7358SJed Brown default: 29*552f7358SJed Brown break; 30*552f7358SJed Brown } 31*552f7358SJed Brown break; 32*552f7358SJed Brown case 2: 33*552f7358SJed Brown switch(corners) { 34*552f7358SJed Brown case 3: 35*552f7358SJed Brown *cellType = 5; /* VTK_TRIANGLE */ 36*552f7358SJed Brown break; 37*552f7358SJed Brown case 4: 38*552f7358SJed Brown *cellType = 9; /* VTK_QUAD */ 39*552f7358SJed Brown break; 40*552f7358SJed Brown case 6: 41*552f7358SJed Brown *cellType = 22; /* VTK_QUADRATIC_TRIANGLE */ 42*552f7358SJed Brown break; 43*552f7358SJed Brown case 9: 44*552f7358SJed Brown *cellType = 23; /* VTK_QUADRATIC_QUAD */ 45*552f7358SJed Brown break; 46*552f7358SJed Brown default: 47*552f7358SJed Brown break; 48*552f7358SJed Brown } 49*552f7358SJed Brown break; 50*552f7358SJed Brown case 3: 51*552f7358SJed Brown switch(corners) { 52*552f7358SJed Brown case 4: 53*552f7358SJed Brown *cellType = 10; /* VTK_TETRA */ 54*552f7358SJed Brown break; 55*552f7358SJed Brown case 8: 56*552f7358SJed Brown *cellType = 12; /* VTK_HEXAHEDRON */ 57*552f7358SJed Brown break; 58*552f7358SJed Brown case 10: 59*552f7358SJed Brown *cellType = 24; /* VTK_QUADRATIC_TETRA */ 60*552f7358SJed Brown break; 61*552f7358SJed Brown case 27: 62*552f7358SJed Brown *cellType = 29; /* VTK_QUADRATIC_HEXAHEDRON */ 63*552f7358SJed Brown break; 64*552f7358SJed Brown default: 65*552f7358SJed Brown break; 66*552f7358SJed Brown } 67*552f7358SJed Brown } 68*552f7358SJed Brown PetscFunctionReturn(0); 69*552f7358SJed Brown } 70*552f7358SJed Brown 71*552f7358SJed Brown #undef __FUNCT__ 72*552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteCells_ASCII" 73*552f7358SJed Brown PetscErrorCode DMPlexVTKWriteCells_ASCII(DM dm, FILE *fp, PetscInt *totalCells) 74*552f7358SJed Brown { 75*552f7358SJed Brown MPI_Comm comm = ((PetscObject) dm)->comm; 76*552f7358SJed Brown IS globalVertexNumbers; 77*552f7358SJed Brown const PetscInt *gvertex; 78*552f7358SJed Brown PetscInt dim; 79*552f7358SJed Brown PetscInt numCorners = 0, totCorners = 0, maxCorners, *corners; 80*552f7358SJed Brown PetscInt numCells = 0, totCells = 0, maxCells, cellHeight; 81*552f7358SJed Brown PetscInt numLabelCells, cMax, cStart, cEnd, c, vStart, vEnd, v; 82*552f7358SJed Brown PetscMPIInt numProcs, rank, proc, tag; 83*552f7358SJed Brown PetscBool hasLabel; 84*552f7358SJed Brown PetscErrorCode ierr; 85*552f7358SJed Brown 86*552f7358SJed Brown PetscFunctionBegin; 87*552f7358SJed Brown ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr); 88*552f7358SJed Brown ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 89*552f7358SJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 90*552f7358SJed Brown ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 91*552f7358SJed Brown ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 92*552f7358SJed Brown ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 93*552f7358SJed Brown ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 94*552f7358SJed Brown ierr = DMPlexGetVTKBounds(dm, &cMax, PETSC_NULL);CHKERRQ(ierr); 95*552f7358SJed Brown if (cMax >= 0) {cEnd = PetscMin(cEnd, cMax);} 96*552f7358SJed Brown ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 97*552f7358SJed Brown hasLabel = numLabelCells > 0 ? PETSC_TRUE : PETSC_FALSE; 98*552f7358SJed Brown for(c = cStart; c < cEnd; ++c) { 99*552f7358SJed Brown PetscInt *closure = PETSC_NULL; 100*552f7358SJed Brown PetscInt closureSize; 101*552f7358SJed Brown 102*552f7358SJed Brown if (hasLabel) { 103*552f7358SJed Brown PetscInt value; 104*552f7358SJed Brown 105*552f7358SJed Brown ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 106*552f7358SJed Brown if (value != 1) continue; 107*552f7358SJed Brown } 108*552f7358SJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 109*552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 110*552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) ++numCorners; 111*552f7358SJed Brown } 112*552f7358SJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 113*552f7358SJed Brown ++numCells; 114*552f7358SJed Brown } 115*552f7358SJed Brown maxCells = numCells; 116*552f7358SJed Brown ierr = MPI_Reduce(&numCells, &totCells, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRQ(ierr); 117*552f7358SJed Brown ierr = MPI_Reduce(&numCells, &maxCells, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr); 118*552f7358SJed Brown ierr = MPI_Reduce(&numCorners, &totCorners, 1, MPIU_INT, MPI_SUM, 0, comm);CHKERRQ(ierr); 119*552f7358SJed Brown ierr = MPI_Reduce(&numCorners, &maxCorners, 1, MPIU_INT, MPI_MAX, 0, comm);CHKERRQ(ierr); 120*552f7358SJed Brown ierr = DMPlexGetVertexNumbering(dm, &globalVertexNumbers);CHKERRQ(ierr); 121*552f7358SJed Brown ierr = ISGetIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr); 122*552f7358SJed Brown ierr = PetscMalloc(maxCells * sizeof(PetscInt), &corners);CHKERRQ(ierr); 123*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "CELLS %d %d\n", totCells, totCorners+totCells);CHKERRQ(ierr); 124*552f7358SJed Brown if (!rank) { 125*552f7358SJed Brown PetscInt *remoteVertices; 126*552f7358SJed Brown 127*552f7358SJed Brown for(c = cStart, numCells = 0; c < cEnd; ++c) { 128*552f7358SJed Brown PetscInt *closure = PETSC_NULL; 129*552f7358SJed Brown PetscInt closureSize, nC = 0; 130*552f7358SJed Brown 131*552f7358SJed Brown if (hasLabel) { 132*552f7358SJed Brown PetscInt value; 133*552f7358SJed Brown 134*552f7358SJed Brown ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 135*552f7358SJed Brown if (value != 1) continue; 136*552f7358SJed Brown } 137*552f7358SJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 138*552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 139*552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 140*552f7358SJed Brown const PetscInt gv = gvertex[closure[v] - vStart]; 141*552f7358SJed Brown closure[nC++] = gv < 0 ? -(gv+1) : gv; 142*552f7358SJed Brown } 143*552f7358SJed Brown } 144*552f7358SJed Brown corners[numCells++] = nC; 145*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "%d ", nC);CHKERRQ(ierr); 146*552f7358SJed Brown for(v = 0; v < nC; ++v) { 147*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, " %d", closure[v]);CHKERRQ(ierr); 148*552f7358SJed Brown } 149*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr); 150*552f7358SJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 151*552f7358SJed Brown } 152*552f7358SJed Brown if (numProcs > 1) {ierr = PetscMalloc((maxCorners+maxCells) * sizeof(PetscInt), &remoteVertices);CHKERRQ(ierr);} 153*552f7358SJed Brown for(proc = 1; proc < numProcs; ++proc) { 154*552f7358SJed Brown MPI_Status status; 155*552f7358SJed Brown 156*552f7358SJed Brown ierr = MPI_Recv(&numCorners, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 157*552f7358SJed Brown ierr = MPI_Recv(remoteVertices, numCorners, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 158*552f7358SJed Brown for(c = 0; c < numCorners;) { 159*552f7358SJed Brown PetscInt nC = remoteVertices[c++]; 160*552f7358SJed Brown 161*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "%d ", nC);CHKERRQ(ierr); 162*552f7358SJed Brown for(v = 0; v < nC; ++v, ++c) { 163*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, " %d", remoteVertices[c]);CHKERRQ(ierr); 164*552f7358SJed Brown } 165*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr); 166*552f7358SJed Brown } 167*552f7358SJed Brown } 168*552f7358SJed Brown if (numProcs > 1) {ierr = PetscFree(remoteVertices);CHKERRQ(ierr);} 169*552f7358SJed Brown } else { 170*552f7358SJed Brown PetscInt *localVertices, numSend = numCells+numCorners, k = 0; 171*552f7358SJed Brown 172*552f7358SJed Brown ierr = PetscMalloc(numSend * sizeof(PetscInt), &localVertices);CHKERRQ(ierr); 173*552f7358SJed Brown for (c = cStart, numCells = 0; c < cEnd; ++c) { 174*552f7358SJed Brown PetscInt *closure = PETSC_NULL; 175*552f7358SJed Brown PetscInt closureSize, nC = 0; 176*552f7358SJed Brown 177*552f7358SJed Brown if (hasLabel) { 178*552f7358SJed Brown PetscInt value; 179*552f7358SJed Brown 180*552f7358SJed Brown ierr = DMPlexGetLabelValue(dm, "vtk", c, &value);CHKERRQ(ierr); 181*552f7358SJed Brown if (value != 1) continue; 182*552f7358SJed Brown } 183*552f7358SJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 184*552f7358SJed Brown for (v = 0; v < closureSize*2; v += 2) { 185*552f7358SJed Brown if ((closure[v] >= vStart) && (closure[v] < vEnd)) { 186*552f7358SJed Brown const PetscInt gv = gvertex[closure[v] - vStart]; 187*552f7358SJed Brown closure[nC++] = gv < 0 ? -(gv+1) : gv; 188*552f7358SJed Brown } 189*552f7358SJed Brown } 190*552f7358SJed Brown corners[numCells++] = nC; 191*552f7358SJed Brown localVertices[k++] = nC; 192*552f7358SJed Brown for(v = 0; v < nC; ++v, ++k) { 193*552f7358SJed Brown localVertices[k] = closure[v]; 194*552f7358SJed Brown } 195*552f7358SJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 196*552f7358SJed Brown } 197*552f7358SJed Brown if (k != numSend) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid number of vertices to send %d should be %d", k, numSend); 198*552f7358SJed Brown ierr = MPI_Send(&numSend, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 199*552f7358SJed Brown ierr = MPI_Send(localVertices, numSend, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 200*552f7358SJed Brown ierr = PetscFree(localVertices);CHKERRQ(ierr); 201*552f7358SJed Brown } 202*552f7358SJed Brown ierr = ISRestoreIndices(globalVertexNumbers, &gvertex);CHKERRQ(ierr); 203*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "CELL_TYPES %d\n", totCells);CHKERRQ(ierr); 204*552f7358SJed Brown if (!rank) { 205*552f7358SJed Brown PetscInt cellType; 206*552f7358SJed Brown 207*552f7358SJed Brown for (c = 0; c < numCells; ++c) { 208*552f7358SJed Brown ierr = DMPlexVTKGetCellType(dm, dim, corners[c], &cellType);CHKERRQ(ierr); 209*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "%d\n", cellType);CHKERRQ(ierr); 210*552f7358SJed Brown } 211*552f7358SJed Brown for(proc = 1; proc < numProcs; ++proc) { 212*552f7358SJed Brown MPI_Status status; 213*552f7358SJed Brown 214*552f7358SJed Brown ierr = MPI_Recv(&numCells, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 215*552f7358SJed Brown ierr = MPI_Recv(corners, numCells, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 216*552f7358SJed Brown for(c = 0; c < numCells; ++c) { 217*552f7358SJed Brown ierr = DMPlexVTKGetCellType(dm, dim, corners[c], &cellType);CHKERRQ(ierr); 218*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "%d\n", cellType);CHKERRQ(ierr); 219*552f7358SJed Brown } 220*552f7358SJed Brown } 221*552f7358SJed Brown } else { 222*552f7358SJed Brown ierr = MPI_Send(&numCells, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 223*552f7358SJed Brown ierr = MPI_Send(corners, numCells, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 224*552f7358SJed Brown } 225*552f7358SJed Brown ierr = PetscFree(corners);CHKERRQ(ierr); 226*552f7358SJed Brown *totalCells = totCells; 227*552f7358SJed Brown PetscFunctionReturn(0); 228*552f7358SJed Brown } 229*552f7358SJed Brown 230*552f7358SJed Brown #undef __FUNCT__ 231*552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteSection_ASCII" 232*552f7358SJed Brown PetscErrorCode DMPlexVTKWriteSection_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec v, FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale) { 233*552f7358SJed Brown MPI_Comm comm = ((PetscObject) dm)->comm; 234*552f7358SJed Brown const MPI_Datatype mpiType = MPIU_SCALAR; 235*552f7358SJed Brown PetscScalar *array; 236*552f7358SJed Brown PetscInt numDof = 0, maxDof; 237*552f7358SJed Brown PetscInt numLabelCells, cellHeight, cMax, cStart, cEnd, numLabelVertices, vMax, vStart, vEnd, pStart, pEnd, p; 238*552f7358SJed Brown PetscMPIInt numProcs, rank, proc, tag; 239*552f7358SJed Brown PetscBool hasLabel; 240*552f7358SJed Brown PetscErrorCode ierr; 241*552f7358SJed Brown 242*552f7358SJed Brown PetscFunctionBegin; 243*552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 244*552f7358SJed Brown PetscValidHeaderSpecific(v,VEC_CLASSID,4); 245*552f7358SJed Brown if (precision < 0) precision = 6; 246*552f7358SJed Brown ierr = PetscCommGetNewTag(comm, &tag);CHKERRQ(ierr); 247*552f7358SJed Brown ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 248*552f7358SJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 249*552f7358SJed Brown ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 250*552f7358SJed Brown /* VTK only wants the values at cells or vertices */ 251*552f7358SJed Brown ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr); 252*552f7358SJed Brown ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 253*552f7358SJed Brown ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 254*552f7358SJed Brown ierr = DMPlexGetVTKBounds(dm, &cMax, &vMax);CHKERRQ(ierr); 255*552f7358SJed Brown if (cMax >= 0) {cEnd = PetscMin(cEnd, cMax);} 256*552f7358SJed Brown if (vMax >= 0) {vEnd = PetscMin(vEnd, vMax);} 257*552f7358SJed Brown pStart = PetscMax(PetscMin(cStart, vStart), pStart); 258*552f7358SJed Brown pEnd = PetscMin(PetscMax(cEnd, vEnd), pEnd); 259*552f7358SJed Brown ierr = DMPlexGetStratumSize(dm, "vtk", 1, &numLabelCells);CHKERRQ(ierr); 260*552f7358SJed Brown ierr = DMPlexGetStratumSize(dm, "vtk", 2, &numLabelVertices);CHKERRQ(ierr); 261*552f7358SJed Brown hasLabel = numLabelCells > 0 || numLabelVertices > 0 ? PETSC_TRUE : PETSC_FALSE; 262*552f7358SJed Brown for(p = pStart; p < pEnd; ++p) { 263*552f7358SJed Brown /* Reject points not either cells or vertices */ 264*552f7358SJed Brown if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue; 265*552f7358SJed Brown if (hasLabel) { 266*552f7358SJed Brown PetscInt value; 267*552f7358SJed Brown 268*552f7358SJed Brown if (((p >= cStart) && (p < cEnd) && numLabelCells) || 269*552f7358SJed Brown ((p >= vStart) && (p < vEnd) && numLabelVertices)) { 270*552f7358SJed Brown ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr); 271*552f7358SJed Brown if (value != 1) continue; 272*552f7358SJed Brown } 273*552f7358SJed Brown } 274*552f7358SJed Brown ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr); 275*552f7358SJed Brown if (numDof) {break;} 276*552f7358SJed Brown } 277*552f7358SJed Brown ierr = MPI_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr); 278*552f7358SJed Brown enforceDof = PetscMax(enforceDof, maxDof); 279*552f7358SJed Brown ierr = VecGetArray(v, &array);CHKERRQ(ierr); 280*552f7358SJed Brown if (!rank) { 281*552f7358SJed Brown char formatString[8]; 282*552f7358SJed Brown 283*552f7358SJed Brown ierr = PetscSNPrintf(formatString, 8, "%%.%de", precision);CHKERRQ(ierr); 284*552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 285*552f7358SJed Brown /* Here we lose a way to filter points by keeping them out of the Numbering */ 286*552f7358SJed Brown PetscInt dof, off, goff, d; 287*552f7358SJed Brown 288*552f7358SJed Brown /* Reject points not either cells or vertices */ 289*552f7358SJed Brown if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue; 290*552f7358SJed Brown if (hasLabel) { 291*552f7358SJed Brown PetscInt value; 292*552f7358SJed Brown 293*552f7358SJed Brown if (((p >= cStart) && (p < cEnd) && numLabelCells) || 294*552f7358SJed Brown ((p >= vStart) && (p < vEnd) && numLabelVertices)) { 295*552f7358SJed Brown ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr); 296*552f7358SJed Brown if (value != 1) continue; 297*552f7358SJed Brown } 298*552f7358SJed Brown } 299*552f7358SJed Brown ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 300*552f7358SJed Brown ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 301*552f7358SJed Brown ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr); 302*552f7358SJed Brown if (dof && goff >= 0) { 303*552f7358SJed Brown for (d = 0; d < dof; d++) { 304*552f7358SJed Brown if (d > 0) { 305*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr); 306*552f7358SJed Brown } 307*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(array[off+d])*scale);CHKERRQ(ierr); 308*552f7358SJed Brown } 309*552f7358SJed Brown for (d = dof; d < enforceDof; d++) { 310*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr); 311*552f7358SJed Brown } 312*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr); 313*552f7358SJed Brown } 314*552f7358SJed Brown } 315*552f7358SJed Brown for (proc = 1; proc < numProcs; ++proc) { 316*552f7358SJed Brown PetscScalar *remoteValues; 317*552f7358SJed Brown PetscInt size, d; 318*552f7358SJed Brown MPI_Status status; 319*552f7358SJed Brown 320*552f7358SJed Brown ierr = MPI_Recv(&size, 1, MPIU_INT, proc, tag, comm, &status);CHKERRQ(ierr); 321*552f7358SJed Brown ierr = PetscMalloc(size * sizeof(PetscScalar), &remoteValues);CHKERRQ(ierr); 322*552f7358SJed Brown ierr = MPI_Recv(remoteValues, size, mpiType, proc, tag, comm, &status);CHKERRQ(ierr); 323*552f7358SJed Brown for (p = 0; p < size/maxDof; ++p) { 324*552f7358SJed Brown for (d = 0; d < maxDof; ++d) { 325*552f7358SJed Brown if (d > 0) { 326*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, " ");CHKERRQ(ierr); 327*552f7358SJed Brown } 328*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, formatString, PetscRealPart(remoteValues[p*maxDof+d])*scale);CHKERRQ(ierr); 329*552f7358SJed Brown } 330*552f7358SJed Brown for (d = maxDof; d < enforceDof; ++d) { 331*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, " 0.0");CHKERRQ(ierr); 332*552f7358SJed Brown } 333*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "\n");CHKERRQ(ierr); 334*552f7358SJed Brown } 335*552f7358SJed Brown ierr = PetscFree(remoteValues);CHKERRQ(ierr); 336*552f7358SJed Brown } 337*552f7358SJed Brown } else { 338*552f7358SJed Brown PetscScalar *localValues; 339*552f7358SJed Brown PetscInt size, k = 0; 340*552f7358SJed Brown 341*552f7358SJed Brown ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr); 342*552f7358SJed Brown ierr = PetscMalloc(size * sizeof(PetscScalar), &localValues);CHKERRQ(ierr); 343*552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 344*552f7358SJed Brown PetscInt dof, off, goff, d; 345*552f7358SJed Brown 346*552f7358SJed Brown /* Reject points not either cells or vertices */ 347*552f7358SJed Brown if (((p < cStart) || (p >= cEnd)) && ((p < vStart) || (p >= vEnd))) continue; 348*552f7358SJed Brown if (hasLabel) { 349*552f7358SJed Brown PetscInt value; 350*552f7358SJed Brown 351*552f7358SJed Brown if (((p >= cStart) && (p < cEnd) && numLabelCells) || 352*552f7358SJed Brown ((p >= vStart) && (p < vEnd) && numLabelVertices)) { 353*552f7358SJed Brown ierr = DMPlexGetLabelValue(dm, "vtk", p, &value);CHKERRQ(ierr); 354*552f7358SJed Brown if (value != 1) continue; 355*552f7358SJed Brown } 356*552f7358SJed Brown } 357*552f7358SJed Brown ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 358*552f7358SJed Brown ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 359*552f7358SJed Brown ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr); 360*552f7358SJed Brown if (goff >= 0) { 361*552f7358SJed Brown for (d = 0; d < dof; ++d) { 362*552f7358SJed Brown localValues[k++] = array[off+d]; 363*552f7358SJed Brown } 364*552f7358SJed Brown } 365*552f7358SJed Brown } 366*552f7358SJed Brown ierr = MPI_Send(&k, 1, MPIU_INT, 0, tag, comm);CHKERRQ(ierr); 367*552f7358SJed Brown ierr = MPI_Send(localValues, k, mpiType, 0, tag, comm);CHKERRQ(ierr); 368*552f7358SJed Brown ierr = PetscFree(localValues);CHKERRQ(ierr); 369*552f7358SJed Brown } 370*552f7358SJed Brown ierr = VecRestoreArray(v, &array);CHKERRQ(ierr); 371*552f7358SJed Brown PetscFunctionReturn(0); 372*552f7358SJed Brown } 373*552f7358SJed Brown 374*552f7358SJed Brown #undef __FUNCT__ 375*552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteField_ASCII" 376*552f7358SJed Brown PetscErrorCode DMPlexVTKWriteField_ASCII(DM dm, PetscSection section, PetscSection globalSection, Vec field, const char name[], FILE *fp, PetscInt enforceDof, PetscInt precision, PetscReal scale) 377*552f7358SJed Brown { 378*552f7358SJed Brown MPI_Comm comm = ((PetscObject) dm)->comm; 379*552f7358SJed Brown PetscInt numDof = 0, maxDof; 380*552f7358SJed Brown PetscInt pStart, pEnd, p; 381*552f7358SJed Brown PetscErrorCode ierr; 382*552f7358SJed Brown 383*552f7358SJed Brown PetscFunctionBegin; 384*552f7358SJed Brown ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 385*552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 386*552f7358SJed Brown ierr = PetscSectionGetDof(section, p, &numDof);CHKERRQ(ierr); 387*552f7358SJed Brown if (numDof) {break;} 388*552f7358SJed Brown } 389*552f7358SJed Brown numDof = PetscMax(numDof, enforceDof); 390*552f7358SJed Brown ierr = MPI_Allreduce(&numDof, &maxDof, 1, MPIU_INT, MPI_MAX, ((PetscObject) dm)->comm);CHKERRQ(ierr); 391*552f7358SJed Brown if (!name) name = "Unknown"; 392*552f7358SJed Brown if (maxDof == 3) { 393*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "VECTORS %s double\n", name);CHKERRQ(ierr); 394*552f7358SJed Brown } else { 395*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "SCALARS %s double %d\n", name, maxDof);CHKERRQ(ierr); 396*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "LOOKUP_TABLE default\n");CHKERRQ(ierr); 397*552f7358SJed Brown } 398*552f7358SJed Brown ierr = DMPlexVTKWriteSection_ASCII(dm, section, globalSection, field, fp, enforceDof, precision, scale);CHKERRQ(ierr); 399*552f7358SJed Brown PetscFunctionReturn(0); 400*552f7358SJed Brown } 401*552f7358SJed Brown 402*552f7358SJed Brown #undef __FUNCT__ 403*552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll_ASCII" 404*552f7358SJed Brown static PetscErrorCode DMPlexVTKWriteAll_ASCII(DM dm, PetscViewer viewer) 405*552f7358SJed Brown { 406*552f7358SJed Brown MPI_Comm comm = ((PetscObject) dm)->comm; 407*552f7358SJed Brown PetscViewer_VTK *vtk = (PetscViewer_VTK *) viewer->data; 408*552f7358SJed Brown FILE *fp; 409*552f7358SJed Brown PetscViewerVTKObjectLink link; 410*552f7358SJed Brown PetscSection coordSection, globalCoordSection; 411*552f7358SJed Brown PetscLayout vLayout; 412*552f7358SJed Brown Vec coordinates; 413*552f7358SJed Brown PetscReal lengthScale; 414*552f7358SJed Brown PetscInt vMax, totVertices, totCells; 415*552f7358SJed Brown PetscBool hasPoint = PETSC_FALSE, hasCell = PETSC_FALSE; 416*552f7358SJed Brown PetscErrorCode ierr; 417*552f7358SJed Brown 418*552f7358SJed Brown PetscFunctionBegin; 419*552f7358SJed Brown ierr = PetscFOpen(comm, vtk->filename, "wb", &fp);CHKERRQ(ierr); 420*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "# vtk DataFile Version 2.0\n");CHKERRQ(ierr); 421*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "Simplicial Mesh Example\n");CHKERRQ(ierr); 422*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "ASCII\n");CHKERRQ(ierr); 423*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "DATASET UNSTRUCTURED_GRID\n");CHKERRQ(ierr); 424*552f7358SJed Brown /* Vertices */ 425*552f7358SJed Brown ierr = DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);CHKERRQ(ierr); 426*552f7358SJed Brown ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 427*552f7358SJed Brown ierr = PetscSectionCreateGlobalSection(coordSection, dm->sf, PETSC_FALSE, &globalCoordSection);CHKERRQ(ierr); 428*552f7358SJed Brown ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 429*552f7358SJed Brown ierr = DMPlexGetVTKBounds(dm, PETSC_NULL, &vMax);CHKERRQ(ierr); 430*552f7358SJed Brown if (vMax >= 0) { 431*552f7358SJed Brown PetscInt pStart, pEnd, p, localSize = 0; 432*552f7358SJed Brown 433*552f7358SJed Brown ierr = PetscSectionGetChart(globalCoordSection, &pStart, &pEnd);CHKERRQ(ierr); 434*552f7358SJed Brown pEnd = PetscMin(pEnd, vMax); 435*552f7358SJed Brown for (p = pStart; p < pEnd; ++p) { 436*552f7358SJed Brown PetscInt dof; 437*552f7358SJed Brown 438*552f7358SJed Brown ierr = PetscSectionGetDof(globalCoordSection, p, &dof);CHKERRQ(ierr); 439*552f7358SJed Brown if (dof > 0) {++localSize;} 440*552f7358SJed Brown } 441*552f7358SJed Brown ierr = PetscLayoutCreate(((PetscObject) dm)->comm, &vLayout);CHKERRQ(ierr); 442*552f7358SJed Brown ierr = PetscLayoutSetLocalSize(vLayout, localSize);CHKERRQ(ierr); 443*552f7358SJed Brown ierr = PetscLayoutSetBlockSize(vLayout, 1);CHKERRQ(ierr); 444*552f7358SJed Brown ierr = PetscLayoutSetUp(vLayout);CHKERRQ(ierr); 445*552f7358SJed Brown } else { 446*552f7358SJed Brown ierr = PetscSectionGetPointLayout(((PetscObject) dm)->comm, globalCoordSection, &vLayout);CHKERRQ(ierr); 447*552f7358SJed Brown } 448*552f7358SJed Brown ierr = PetscLayoutGetSize(vLayout, &totVertices);CHKERRQ(ierr); 449*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "POINTS %d double\n", totVertices);CHKERRQ(ierr); 450*552f7358SJed Brown ierr = DMPlexVTKWriteSection_ASCII(dm, coordSection, globalCoordSection, coordinates, fp, 3, PETSC_DETERMINE, lengthScale);CHKERRQ(ierr); 451*552f7358SJed Brown /* Cells */ 452*552f7358SJed Brown ierr = DMPlexVTKWriteCells_ASCII(dm, fp, &totCells);CHKERRQ(ierr); 453*552f7358SJed Brown /* Vertex fields */ 454*552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 455*552f7358SJed Brown if ((link->ft == PETSC_VTK_POINT_FIELD) || (link->ft == PETSC_VTK_POINT_VECTOR_FIELD)) hasPoint = PETSC_TRUE; 456*552f7358SJed Brown if ((link->ft == PETSC_VTK_CELL_FIELD) || (link->ft == PETSC_VTK_CELL_VECTOR_FIELD)) hasCell = PETSC_TRUE; 457*552f7358SJed Brown } 458*552f7358SJed Brown if (hasPoint) { 459*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "POINT_DATA %d\n", totVertices); 460*552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 461*552f7358SJed Brown Vec X = (Vec) link->vec; 462*552f7358SJed Brown DM dmX; 463*552f7358SJed Brown PetscSection section, globalSection; 464*552f7358SJed Brown const char *name; 465*552f7358SJed Brown PetscInt enforceDof = PETSC_DETERMINE; 466*552f7358SJed Brown 467*552f7358SJed Brown if ((link->ft != PETSC_VTK_POINT_FIELD) && (link->ft != PETSC_VTK_POINT_VECTOR_FIELD)) continue; 468*552f7358SJed Brown if (link->ft == PETSC_VTK_POINT_VECTOR_FIELD) enforceDof = 3; 469*552f7358SJed Brown ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr); 470*552f7358SJed Brown ierr = VecGetDM(X, &dmX);CHKERRQ(ierr); 471*552f7358SJed Brown if (dmX) { 472*552f7358SJed Brown ierr = DMGetDefaultSection(dmX, §ion); 473*552f7358SJed Brown } else { 474*552f7358SJed Brown PetscContainer c; 475*552f7358SJed Brown 476*552f7358SJed Brown ierr = PetscObjectQuery(link->vec, "section", (PetscObject *) &c);CHKERRQ(ierr); 477*552f7358SJed Brown if (!c) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name); 478*552f7358SJed Brown ierr = PetscContainerGetPointer(c, (void **) §ion);CHKERRQ(ierr); 479*552f7358SJed Brown } 480*552f7358SJed Brown if (!section) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name); 481*552f7358SJed Brown ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, &globalSection);CHKERRQ(ierr); 482*552f7358SJed Brown ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr); 483*552f7358SJed Brown ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr); 484*552f7358SJed Brown } 485*552f7358SJed Brown } 486*552f7358SJed Brown /* Cell Fields */ 487*552f7358SJed Brown if (hasCell) { 488*552f7358SJed Brown ierr = PetscFPrintf(comm, fp, "CELL_DATA %d\n", totCells); 489*552f7358SJed Brown for (link = vtk->link; link; link = link->next) { 490*552f7358SJed Brown Vec X = (Vec) link->vec; 491*552f7358SJed Brown DM dmX; 492*552f7358SJed Brown PetscSection section, globalSection; 493*552f7358SJed Brown const char *name; 494*552f7358SJed Brown PetscInt enforceDof = PETSC_DETERMINE; 495*552f7358SJed Brown 496*552f7358SJed Brown if ((link->ft != PETSC_VTK_CELL_FIELD) && (link->ft != PETSC_VTK_CELL_VECTOR_FIELD)) continue; 497*552f7358SJed Brown if (link->ft == PETSC_VTK_CELL_VECTOR_FIELD) enforceDof = 3; 498*552f7358SJed Brown ierr = PetscObjectGetName(link->vec, &name);CHKERRQ(ierr); 499*552f7358SJed Brown ierr = VecGetDM(X, &dmX);CHKERRQ(ierr); 500*552f7358SJed Brown if (dmX) { 501*552f7358SJed Brown ierr = DMGetDefaultSection(dmX, §ion); 502*552f7358SJed Brown } else { 503*552f7358SJed Brown PetscContainer c; 504*552f7358SJed Brown 505*552f7358SJed Brown ierr = PetscObjectQuery(link->vec, "section", (PetscObject *) &c);CHKERRQ(ierr); 506*552f7358SJed Brown if (!c) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name); 507*552f7358SJed Brown ierr = PetscContainerGetPointer(c, (void **) §ion);CHKERRQ(ierr); 508*552f7358SJed Brown } 509*552f7358SJed Brown if (!section) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "Vector %s had no PetscSection composed with it", name); 510*552f7358SJed Brown ierr = PetscSectionCreateGlobalSection(section, dm->sf, PETSC_FALSE, &globalSection);CHKERRQ(ierr); 511*552f7358SJed Brown ierr = DMPlexVTKWriteField_ASCII(dm, section, globalSection, X, name, fp, enforceDof, PETSC_DETERMINE, 1.0);CHKERRQ(ierr); 512*552f7358SJed Brown ierr = PetscSectionDestroy(&globalSection);CHKERRQ(ierr); 513*552f7358SJed Brown } 514*552f7358SJed Brown } 515*552f7358SJed Brown /* Cleanup */ 516*552f7358SJed Brown ierr = PetscSectionDestroy(&globalCoordSection);CHKERRQ(ierr); 517*552f7358SJed Brown ierr = PetscLayoutDestroy(&vLayout);CHKERRQ(ierr); 518*552f7358SJed Brown ierr = PetscFClose(comm, fp);CHKERRQ(ierr); 519*552f7358SJed Brown PetscFunctionReturn(0); 520*552f7358SJed Brown } 521*552f7358SJed Brown 522*552f7358SJed Brown #undef __FUNCT__ 523*552f7358SJed Brown #define __FUNCT__ "DMPlexVTKWriteAll" 524*552f7358SJed Brown /*@C 525*552f7358SJed Brown DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer 526*552f7358SJed Brown 527*552f7358SJed Brown Collective 528*552f7358SJed Brown 529*552f7358SJed Brown Input Arguments: 530*552f7358SJed Brown + odm - The DMPlex specifying the mesh, passed as a PetscObject 531*552f7358SJed Brown - viewer - viewer of type VTK 532*552f7358SJed Brown 533*552f7358SJed Brown Level: developer 534*552f7358SJed Brown 535*552f7358SJed Brown Note: 536*552f7358SJed Brown This function is a callback used by the VTK viewer to actually write the file. 537*552f7358SJed 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. 538*552f7358SJed Brown Instead, metadata for the entire file needs to be available up-front before you can start writing the file. 539*552f7358SJed Brown 540*552f7358SJed Brown .seealso: PETSCVIEWERVTK 541*552f7358SJed Brown @*/ 542*552f7358SJed Brown PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer) 543*552f7358SJed Brown { 544*552f7358SJed Brown DM dm = (DM) odm; 545*552f7358SJed Brown PetscBool isvtk; 546*552f7358SJed Brown PetscErrorCode ierr; 547*552f7358SJed Brown 548*552f7358SJed Brown PetscFunctionBegin; 549*552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 550*552f7358SJed Brown PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 551*552f7358SJed Brown ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr); 552*552f7358SJed Brown if (!isvtk) SETERRQ1(((PetscObject) viewer)->comm, PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name); 553*552f7358SJed Brown switch (viewer->format) { 554*552f7358SJed Brown case PETSC_VIEWER_ASCII_VTK: 555*552f7358SJed Brown ierr = DMPlexVTKWriteAll_ASCII(dm, viewer);CHKERRQ(ierr); 556*552f7358SJed Brown break; 557*552f7358SJed Brown case PETSC_VIEWER_VTK_VTU: 558*552f7358SJed Brown ierr = DMPlexVTKWriteAll_VTU(dm, viewer);CHKERRQ(ierr); 559*552f7358SJed Brown break; 560*552f7358SJed Brown default: SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]); 561*552f7358SJed Brown } 562*552f7358SJed Brown PetscFunctionReturn(0); 563*552f7358SJed Brown } 564