1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> 294fbd55eSBarry Smith /* 394fbd55eSBarry Smith Note that the API for using PETSCVIEWERVTK is totally wrong since its use requires 494fbd55eSBarry Smith including the private vtkvimpl.h file. The code should be refactored. 594fbd55eSBarry Smith */ 65c6c1daeSBarry Smith #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h> 74061b8bfSJed Brown 8fb5bd1c2SPatrick Sanan /* Helper function which determines if any DMDA fields are named. This is used 9fb5bd1c2SPatrick Sanan as a proxy for the user's intention to use DMDA fields as distinct 10fb5bd1c2SPatrick Sanan scalar-valued fields as opposed to a single vector-valued field */ 119371c9d4SSatish Balay static PetscErrorCode DMDAGetFieldsNamed(DM da, PetscBool *fieldsnamed) { 12fb5bd1c2SPatrick Sanan PetscInt f, bs; 13fb5bd1c2SPatrick Sanan 14fb5bd1c2SPatrick Sanan PetscFunctionBegin; 15fb5bd1c2SPatrick Sanan *fieldsnamed = PETSC_FALSE; 169566063dSJacob Faibussowitsch PetscCall(DMDAGetDof(da, &bs)); 17fb5bd1c2SPatrick Sanan for (f = 0; f < bs; ++f) { 18fb5bd1c2SPatrick Sanan const char *fieldname; 199566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da, f, &fieldname)); 20fb5bd1c2SPatrick Sanan if (fieldname) { 21fb5bd1c2SPatrick Sanan *fieldsnamed = PETSC_TRUE; 22fb5bd1c2SPatrick Sanan break; 23fb5bd1c2SPatrick Sanan } 24fb5bd1c2SPatrick Sanan } 25fb5bd1c2SPatrick Sanan PetscFunctionReturn(0); 26fb5bd1c2SPatrick Sanan } 27fb5bd1c2SPatrick Sanan 289371c9d4SSatish Balay static PetscErrorCode DMDAVTKWriteAll_VTS(DM da, PetscViewer viewer) { 2930815ce0SLisandro Dalcin const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian"; 304061b8bfSJed Brown #if defined(PETSC_USE_REAL_SINGLE) 314061b8bfSJed Brown const char precision[] = "Float32"; 324061b8bfSJed Brown #elif defined(PETSC_USE_REAL_DOUBLE) 334061b8bfSJed Brown const char precision[] = "Float64"; 344061b8bfSJed Brown #else 354061b8bfSJed Brown const char precision[] = "UnknownPrecision"; 364061b8bfSJed Brown #endif 37ce94432eSBarry Smith MPI_Comm comm; 382eaa9ef4SLisandro Dalcin Vec Coords; 394061b8bfSJed Brown PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data; 404061b8bfSJed Brown PetscViewerVTKObjectLink link; 414061b8bfSJed Brown FILE *fp; 424061b8bfSJed Brown PetscMPIInt rank, size, tag; 434061b8bfSJed Brown DMDALocalInfo info; 44ea2d7708SPatrick Sanan PetscInt dim, mx, my, mz, cdim, bs, boffset, maxnnodes, maxbs, i, j, k, r; 450298fd71SBarry Smith PetscInt rloc[6], (*grloc)[6] = NULL; 464061b8bfSJed Brown PetscScalar *array, *array2; 474061b8bfSJed Brown 484061b8bfSJed Brown PetscFunctionBegin; 499566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 501dca8a05SBarry Smith PetscCheck(!PetscDefined(USE_COMPLEX), PETSC_COMM_SELF, PETSC_ERR_SUP, "Complex values not supported"); 519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 539566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, &mx, &my, &mz, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL)); 549566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 559566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &Coords)); 562eaa9ef4SLisandro Dalcin if (Coords) { 572eaa9ef4SLisandro Dalcin PetscInt csize; 589566063dSJacob Faibussowitsch PetscCall(VecGetSize(Coords, &csize)); 5908401ef6SPierre Jolivet PetscCheck(csize % (mx * my * mz) == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Coordinate vector size mismatch"); 602eaa9ef4SLisandro Dalcin cdim = csize / (mx * my * mz); 611dca8a05SBarry Smith PetscCheck(cdim >= dim && cdim <= 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Coordinate vector size mismatch"); 622eaa9ef4SLisandro Dalcin } else { 632eaa9ef4SLisandro Dalcin cdim = dim; 642eaa9ef4SLisandro Dalcin } 654061b8bfSJed Brown 669566063dSJacob Faibussowitsch PetscCall(PetscFOpen(comm, vtk->filename, "wb", &fp)); 679566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "<?xml version=\"1.0\"?>\n")); 689566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"%s\">\n", byte_order)); 6963a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <StructuredGrid WholeExtent=\"%d %" PetscInt_FMT " %d %" PetscInt_FMT " %d %" PetscInt_FMT "\">\n", 0, mx - 1, 0, my - 1, 0, mz - 1)); 704061b8bfSJed Brown 719566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc1(size * 6, &grloc)); 724061b8bfSJed Brown rloc[0] = info.xs; 734061b8bfSJed Brown rloc[1] = info.xm; 744061b8bfSJed Brown rloc[2] = info.ys; 754061b8bfSJed Brown rloc[3] = info.ym; 764061b8bfSJed Brown rloc[4] = info.zs; 774061b8bfSJed Brown rloc[5] = info.zm; 789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gather(rloc, 6, MPIU_INT, &grloc[0][0], 6, MPIU_INT, 0, comm)); 794061b8bfSJed Brown 804061b8bfSJed Brown /* Write XML header */ 814061b8bfSJed Brown maxnnodes = 0; /* Used for the temporary array size on rank 0 */ 82ea2d7708SPatrick Sanan maxbs = 0; /* Used for the temporary array size on rank 0 */ 834061b8bfSJed Brown boffset = 0; /* Offset into binary file */ 844061b8bfSJed Brown for (r = 0; r < size; r++) { 854061b8bfSJed Brown PetscInt xs = -1, xm = -1, ys = -1, ym = -1, zs = -1, zm = -1, nnodes = 0; 86dd400576SPatrick Sanan if (rank == 0) { 874061b8bfSJed Brown xs = grloc[r][0]; 884061b8bfSJed Brown xm = grloc[r][1]; 894061b8bfSJed Brown ys = grloc[r][2]; 904061b8bfSJed Brown ym = grloc[r][3]; 914061b8bfSJed Brown zs = grloc[r][4]; 924061b8bfSJed Brown zm = grloc[r][5]; 934061b8bfSJed Brown nnodes = xm * ym * zm; 944061b8bfSJed Brown } 954061b8bfSJed Brown maxnnodes = PetscMax(maxnnodes, nnodes); 9663a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <Piece Extent=\"%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\">\n", xs, xs + xm - 1, ys, ys + ym - 1, zs, zs + zm - 1)); 979566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <Points>\n")); 9863a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, boffset)); 994061b8bfSJed Brown boffset += 3 * nnodes * sizeof(PetscScalar) + sizeof(int); 1009566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </Points>\n")); 1014061b8bfSJed Brown 1029566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <PointData Scalars=\"ScalarPointData\">\n")); 1034061b8bfSJed Brown for (link = vtk->link; link; link = link->next) { 1044061b8bfSJed Brown Vec X = (Vec)link->vec; 105fb5bd1c2SPatrick Sanan PetscInt bs, f; 106ea2d7708SPatrick Sanan DM daCurr; 107fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 108fb5bd1c2SPatrick Sanan const char *vecname = "Unnamed"; 109ea2d7708SPatrick Sanan 1109566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &daCurr)); 1119566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL)); 112ea2d7708SPatrick Sanan maxbs = PetscMax(maxbs, bs); 113ea2d7708SPatrick Sanan 114*48a46eb9SPierre Jolivet if (((PetscObject)X)->name || link != vtk->link) PetscCall(PetscObjectGetName((PetscObject)X, &vecname)); 115fb5bd1c2SPatrick Sanan 116fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 1179566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr, &fieldsnamed)); 118fb5bd1c2SPatrick Sanan if (fieldsnamed) { 119fb5bd1c2SPatrick Sanan for (f = 0; f < bs; f++) { 120fb5bd1c2SPatrick Sanan char buf[256]; 121fb5bd1c2SPatrick Sanan const char *fieldname; 1229566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(daCurr, f, &fieldname)); 123fb5bd1c2SPatrick Sanan if (!fieldname) { 12463a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, f)); 125fb5bd1c2SPatrick Sanan fieldname = buf; 126fb5bd1c2SPatrick Sanan } 12763a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, vecname, fieldname, boffset)); 128fb5bd1c2SPatrick Sanan boffset += nnodes * sizeof(PetscScalar) + sizeof(int); 129fb5bd1c2SPatrick Sanan } 130fb5bd1c2SPatrick Sanan } else { 13163a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%" PetscInt_FMT "\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, vecname, bs, boffset)); 132ea2d7708SPatrick Sanan boffset += bs * nnodes * sizeof(PetscScalar) + sizeof(int); 1334061b8bfSJed Brown } 134fb5bd1c2SPatrick Sanan } 1359566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </PointData>\n")); 1369566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </Piece>\n")); 1374061b8bfSJed Brown } 1389566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </StructuredGrid>\n")); 1399566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <AppendedData encoding=\"raw\">\n")); 1409566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "_")); 1414061b8bfSJed Brown 1424061b8bfSJed Brown /* Now write the arrays. */ 1434061b8bfSJed Brown tag = ((PetscObject)viewer)->tag; 1449566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(maxnnodes * PetscMax(3, maxbs), &array, maxnnodes * PetscMax(3, maxbs), &array2)); 1454061b8bfSJed Brown for (r = 0; r < size; r++) { 1464061b8bfSJed Brown MPI_Status status; 1474061b8bfSJed Brown PetscInt xs = -1, xm = -1, ys = -1, ym = -1, zs = -1, zm = -1, nnodes = 0; 148dd400576SPatrick Sanan if (rank == 0) { 1494061b8bfSJed Brown xs = grloc[r][0]; 1504061b8bfSJed Brown xm = grloc[r][1]; 1514061b8bfSJed Brown ys = grloc[r][2]; 1524061b8bfSJed Brown ym = grloc[r][3]; 1534061b8bfSJed Brown zs = grloc[r][4]; 1544061b8bfSJed Brown zm = grloc[r][5]; 1554061b8bfSJed Brown nnodes = xm * ym * zm; 1564061b8bfSJed Brown } else if (r == rank) { 1574061b8bfSJed Brown nnodes = info.xm * info.ym * info.zm; 1584061b8bfSJed Brown } 1594061b8bfSJed Brown 1602eaa9ef4SLisandro Dalcin /* Write the coordinates */ 1614061b8bfSJed Brown if (Coords) { 1624061b8bfSJed Brown const PetscScalar *coords; 1639566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Coords, &coords)); 164dd400576SPatrick Sanan if (rank == 0) { 1654061b8bfSJed Brown if (r) { 1666622f924SJed Brown PetscMPIInt nn; 1679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array, nnodes * cdim, MPIU_SCALAR, r, tag, comm, &status)); 1689566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn)); 16908401ef6SPierre Jolivet PetscCheck(nn == nnodes * cdim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch"); 1701baa6e33SBarry Smith } else PetscCall(PetscArraycpy(array, coords, nnodes * cdim)); 1714061b8bfSJed Brown /* Transpose coordinates to VTK (C-style) ordering */ 1724061b8bfSJed Brown for (k = 0; k < zm; k++) { 1734061b8bfSJed Brown for (j = 0; j < ym; j++) { 1744061b8bfSJed Brown for (i = 0; i < xm; i++) { 1754061b8bfSJed Brown PetscInt Iloc = i + xm * (j + ym * k); 1762eaa9ef4SLisandro Dalcin array2[Iloc * 3 + 0] = array[Iloc * cdim + 0]; 1772eaa9ef4SLisandro Dalcin array2[Iloc * 3 + 1] = cdim > 1 ? array[Iloc * cdim + 1] : 0.0; 1782eaa9ef4SLisandro Dalcin array2[Iloc * 3 + 2] = cdim > 2 ? array[Iloc * cdim + 2] : 0.0; 1794061b8bfSJed Brown } 1804061b8bfSJed Brown } 1814061b8bfSJed Brown } 1824061b8bfSJed Brown } else if (r == rank) { 1839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void *)coords, nnodes * cdim, MPIU_SCALAR, 0, tag, comm)); 1844061b8bfSJed Brown } 1859566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Coords, &coords)); 1864061b8bfSJed Brown } else { /* Fabricate some coordinates using grid index */ 1874061b8bfSJed Brown for (k = 0; k < zm; k++) { 1884061b8bfSJed Brown for (j = 0; j < ym; j++) { 1894061b8bfSJed Brown for (i = 0; i < xm; i++) { 1904061b8bfSJed Brown PetscInt Iloc = i + xm * (j + ym * k); 1914061b8bfSJed Brown array2[Iloc * 3 + 0] = xs + i; 1924061b8bfSJed Brown array2[Iloc * 3 + 1] = ys + j; 1934061b8bfSJed Brown array2[Iloc * 3 + 2] = zs + k; 1944061b8bfSJed Brown } 1954061b8bfSJed Brown } 1964061b8bfSJed Brown } 1974061b8bfSJed Brown } 1989566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, array2, nnodes * 3, MPIU_SCALAR)); 1994061b8bfSJed Brown 2004061b8bfSJed Brown /* Write each of the objects queued up for this file */ 2014061b8bfSJed Brown for (link = vtk->link; link; link = link->next) { 2024061b8bfSJed Brown Vec X = (Vec)link->vec; 2034061b8bfSJed Brown const PetscScalar *x; 204fb5bd1c2SPatrick Sanan PetscInt bs, f; 205ea2d7708SPatrick Sanan DM daCurr; 206fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 2079566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &daCurr)); 2089566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL)); 2099566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 210dd400576SPatrick Sanan if (rank == 0) { 2114061b8bfSJed Brown if (r) { 2126622f924SJed Brown PetscMPIInt nn; 2139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array, nnodes * bs, MPIU_SCALAR, r, tag, comm, &status)); 2149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn)); 21563a3b9bcSJacob Faibussowitsch PetscCheck(nn == nnodes * bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch receiving from rank %" PetscInt_FMT, r); 2161baa6e33SBarry Smith } else PetscCall(PetscArraycpy(array, x, nnodes * bs)); 217fb5bd1c2SPatrick Sanan 218fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 2199566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr, &fieldsnamed)); 220fb5bd1c2SPatrick Sanan if (fieldsnamed) { 221fb5bd1c2SPatrick Sanan for (f = 0; f < bs; f++) { 222fb5bd1c2SPatrick Sanan /* Extract and transpose the f'th field */ 223fb5bd1c2SPatrick Sanan for (k = 0; k < zm; k++) { 224fb5bd1c2SPatrick Sanan for (j = 0; j < ym; j++) { 225fb5bd1c2SPatrick Sanan for (i = 0; i < xm; i++) { 226fb5bd1c2SPatrick Sanan PetscInt Iloc = i + xm * (j + ym * k); 227fb5bd1c2SPatrick Sanan array2[Iloc] = array[Iloc * bs + f]; 228fb5bd1c2SPatrick Sanan } 229fb5bd1c2SPatrick Sanan } 230fb5bd1c2SPatrick Sanan } 2319566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, array2, nnodes, MPIU_SCALAR)); 232fb5bd1c2SPatrick Sanan } 2331baa6e33SBarry Smith } else PetscCall(PetscViewerVTKFWrite(viewer, fp, array, bs * nnodes, MPIU_SCALAR)); 2341baa6e33SBarry Smith } else if (r == rank) PetscCallMPI(MPI_Send((void *)x, nnodes * bs, MPIU_SCALAR, 0, tag, comm)); 2359566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 2364061b8bfSJed Brown } 2374061b8bfSJed Brown } 2389566063dSJacob Faibussowitsch PetscCall(PetscFree2(array, array2)); 2399566063dSJacob Faibussowitsch PetscCall(PetscFree(grloc)); 2404061b8bfSJed Brown 2419566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "\n </AppendedData>\n")); 2429566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "</VTKFile>\n")); 2439566063dSJacob Faibussowitsch PetscCall(PetscFClose(comm, fp)); 2444061b8bfSJed Brown PetscFunctionReturn(0); 2454061b8bfSJed Brown } 2464061b8bfSJed Brown 2479371c9d4SSatish Balay static PetscErrorCode DMDAVTKWriteAll_VTR(DM da, PetscViewer viewer) { 24830815ce0SLisandro Dalcin const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian"; 249a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_REAL_SINGLE) 250a13bc4e3SShao-Ching Huang const char precision[] = "Float32"; 251a13bc4e3SShao-Ching Huang #elif defined(PETSC_USE_REAL_DOUBLE) 252a13bc4e3SShao-Ching Huang const char precision[] = "Float64"; 253a13bc4e3SShao-Ching Huang #else 254a13bc4e3SShao-Ching Huang const char precision[] = "UnknownPrecision"; 255a13bc4e3SShao-Ching Huang #endif 256a13bc4e3SShao-Ching Huang MPI_Comm comm; 257a13bc4e3SShao-Ching Huang PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data; 258a13bc4e3SShao-Ching Huang PetscViewerVTKObjectLink link; 259a13bc4e3SShao-Ching Huang FILE *fp; 260a13bc4e3SShao-Ching Huang PetscMPIInt rank, size, tag; 261a13bc4e3SShao-Ching Huang DMDALocalInfo info; 262ea2d7708SPatrick Sanan PetscInt dim, mx, my, mz, boffset, maxnnodes, maxbs, i, j, k, r; 263a13bc4e3SShao-Ching Huang PetscInt rloc[6], (*grloc)[6] = NULL; 264a13bc4e3SShao-Ching Huang PetscScalar *array, *array2; 265a13bc4e3SShao-Ching Huang 266a13bc4e3SShao-Ching Huang PetscFunctionBegin; 2679566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 2681dca8a05SBarry Smith PetscCheck(!PetscDefined(USE_COMPLEX), PETSC_COMM_SELF, PETSC_ERR_SUP, "Complex values not supported"); 2699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2719566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, &mx, &my, &mz, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2729566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 2739566063dSJacob Faibussowitsch PetscCall(PetscFOpen(comm, vtk->filename, "wb", &fp)); 2749566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "<?xml version=\"1.0\"?>\n")); 2759566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "<VTKFile type=\"RectilinearGrid\" version=\"0.1\" byte_order=\"%s\">\n", byte_order)); 27663a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <RectilinearGrid WholeExtent=\"%d %" PetscInt_FMT " %d %" PetscInt_FMT " %d %" PetscInt_FMT "\">\n", 0, mx - 1, 0, my - 1, 0, mz - 1)); 277a13bc4e3SShao-Ching Huang 2789566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc1(size * 6, &grloc)); 279a13bc4e3SShao-Ching Huang rloc[0] = info.xs; 280a13bc4e3SShao-Ching Huang rloc[1] = info.xm; 281a13bc4e3SShao-Ching Huang rloc[2] = info.ys; 282a13bc4e3SShao-Ching Huang rloc[3] = info.ym; 283a13bc4e3SShao-Ching Huang rloc[4] = info.zs; 284a13bc4e3SShao-Ching Huang rloc[5] = info.zm; 2859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gather(rloc, 6, MPIU_INT, &grloc[0][0], 6, MPIU_INT, 0, comm)); 286a13bc4e3SShao-Ching Huang 287a13bc4e3SShao-Ching Huang /* Write XML header */ 288a13bc4e3SShao-Ching Huang maxnnodes = 0; /* Used for the temporary array size on rank 0 */ 289ea2d7708SPatrick Sanan maxbs = 0; /* Used for the temporary array size on rank 0 */ 290a13bc4e3SShao-Ching Huang boffset = 0; /* Offset into binary file */ 291a13bc4e3SShao-Ching Huang for (r = 0; r < size; r++) { 292a13bc4e3SShao-Ching Huang PetscInt xs = -1, xm = -1, ys = -1, ym = -1, zs = -1, zm = -1, nnodes = 0; 293dd400576SPatrick Sanan if (rank == 0) { 294a13bc4e3SShao-Ching Huang xs = grloc[r][0]; 295a13bc4e3SShao-Ching Huang xm = grloc[r][1]; 296a13bc4e3SShao-Ching Huang ys = grloc[r][2]; 297a13bc4e3SShao-Ching Huang ym = grloc[r][3]; 298a13bc4e3SShao-Ching Huang zs = grloc[r][4]; 299a13bc4e3SShao-Ching Huang zm = grloc[r][5]; 300a13bc4e3SShao-Ching Huang nnodes = xm * ym * zm; 301a13bc4e3SShao-Ching Huang } 302a13bc4e3SShao-Ching Huang maxnnodes = PetscMax(maxnnodes, nnodes); 30363a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <Piece Extent=\"%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\">\n", xs, xs + xm - 1, ys, ys + ym - 1, zs, zs + zm - 1)); 3049566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <Coordinates>\n")); 30563a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"Xcoord\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, boffset)); 306a13bc4e3SShao-Ching Huang boffset += xm * sizeof(PetscScalar) + sizeof(int); 30763a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"Ycoord\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, boffset)); 308a13bc4e3SShao-Ching Huang boffset += ym * sizeof(PetscScalar) + sizeof(int); 30963a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"Zcoord\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, boffset)); 310a13bc4e3SShao-Ching Huang boffset += zm * sizeof(PetscScalar) + sizeof(int); 3119566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </Coordinates>\n")); 3129566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <PointData Scalars=\"ScalarPointData\">\n")); 313a13bc4e3SShao-Ching Huang for (link = vtk->link; link; link = link->next) { 314a13bc4e3SShao-Ching Huang Vec X = (Vec)link->vec; 315fb5bd1c2SPatrick Sanan PetscInt bs, f; 316ea2d7708SPatrick Sanan DM daCurr; 317fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 318fb5bd1c2SPatrick Sanan const char *vecname = "Unnamed"; 319ea2d7708SPatrick Sanan 3209566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &daCurr)); 3219566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL)); 322ea2d7708SPatrick Sanan maxbs = PetscMax(maxbs, bs); 323*48a46eb9SPierre Jolivet if (((PetscObject)X)->name || link != vtk->link) PetscCall(PetscObjectGetName((PetscObject)X, &vecname)); 324ea2d7708SPatrick Sanan 325fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 3269566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr, &fieldsnamed)); 327fb5bd1c2SPatrick Sanan if (fieldsnamed) { 328fb5bd1c2SPatrick Sanan for (f = 0; f < bs; f++) { 329fb5bd1c2SPatrick Sanan char buf[256]; 330fb5bd1c2SPatrick Sanan const char *fieldname; 3319566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(daCurr, f, &fieldname)); 332fb5bd1c2SPatrick Sanan if (!fieldname) { 33363a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, f)); 334fb5bd1c2SPatrick Sanan fieldname = buf; 335fb5bd1c2SPatrick Sanan } 33663a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, vecname, fieldname, boffset)); 337fb5bd1c2SPatrick Sanan boffset += nnodes * sizeof(PetscScalar) + sizeof(int); 338fb5bd1c2SPatrick Sanan } 339fb5bd1c2SPatrick Sanan } else { 34063a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%" PetscInt_FMT "\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, vecname, bs, boffset)); 341ea2d7708SPatrick Sanan boffset += bs * nnodes * sizeof(PetscScalar) + sizeof(int); 342a13bc4e3SShao-Ching Huang } 343fb5bd1c2SPatrick Sanan } 3449566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </PointData>\n")); 3459566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </Piece>\n")); 346a13bc4e3SShao-Ching Huang } 3479566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </RectilinearGrid>\n")); 3489566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <AppendedData encoding=\"raw\">\n")); 3499566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "_")); 350a13bc4e3SShao-Ching Huang 351a13bc4e3SShao-Ching Huang /* Now write the arrays. */ 352a13bc4e3SShao-Ching Huang tag = ((PetscObject)viewer)->tag; 3539566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(PetscMax(maxnnodes * maxbs, info.xm + info.ym + info.zm), &array, PetscMax(maxnnodes * maxbs, info.xm + info.ym + info.zm), &array2)); 354a13bc4e3SShao-Ching Huang for (r = 0; r < size; r++) { 355a13bc4e3SShao-Ching Huang MPI_Status status; 356a13bc4e3SShao-Ching Huang PetscInt xs = -1, xm = -1, ys = -1, ym = -1, zs = -1, zm = -1, nnodes = 0; 357dd400576SPatrick Sanan if (rank == 0) { 358a13bc4e3SShao-Ching Huang xs = grloc[r][0]; 359a13bc4e3SShao-Ching Huang xm = grloc[r][1]; 360a13bc4e3SShao-Ching Huang ys = grloc[r][2]; 361a13bc4e3SShao-Ching Huang ym = grloc[r][3]; 362a13bc4e3SShao-Ching Huang zs = grloc[r][4]; 363a13bc4e3SShao-Ching Huang zm = grloc[r][5]; 364a13bc4e3SShao-Ching Huang nnodes = xm * ym * zm; 365a13bc4e3SShao-Ching Huang } else if (r == rank) { 366a13bc4e3SShao-Ching Huang nnodes = info.xm * info.ym * info.zm; 367a13bc4e3SShao-Ching Huang } 368a13bc4e3SShao-Ching Huang { /* Write the coordinates */ 369a13bc4e3SShao-Ching Huang Vec Coords; 370ea2d7708SPatrick Sanan 3719566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &Coords)); 372a13bc4e3SShao-Ching Huang if (Coords) { 373a13bc4e3SShao-Ching Huang const PetscScalar *coords; 3749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Coords, &coords)); 375dd400576SPatrick Sanan if (rank == 0) { 376a13bc4e3SShao-Ching Huang if (r) { 377a13bc4e3SShao-Ching Huang PetscMPIInt nn; 3789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array, xm + ym + zm, MPIU_SCALAR, r, tag, comm, &status)); 3799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn)); 38008401ef6SPierre Jolivet PetscCheck(nn == xm + ym + zm, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch"); 381a13bc4e3SShao-Ching Huang } else { 382a13bc4e3SShao-Ching Huang /* x coordinates */ 383a13bc4e3SShao-Ching Huang for (j = 0, k = 0, i = 0; i < xm; i++) { 384a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 385a13bc4e3SShao-Ching Huang array[i] = coords[Iloc * dim + 0]; 386a13bc4e3SShao-Ching Huang } 387a13bc4e3SShao-Ching Huang /* y coordinates */ 388a13bc4e3SShao-Ching Huang for (i = 0, k = 0, j = 0; j < ym; j++) { 389a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 390a13bc4e3SShao-Ching Huang array[j + xm] = dim > 1 ? coords[Iloc * dim + 1] : 0; 391a13bc4e3SShao-Ching Huang } 392a13bc4e3SShao-Ching Huang /* z coordinates */ 393a13bc4e3SShao-Ching Huang for (i = 0, j = 0, k = 0; k < zm; k++) { 394a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 395a13bc4e3SShao-Ching Huang array[k + xm + ym] = dim > 2 ? coords[Iloc * dim + 2] : 0; 396a13bc4e3SShao-Ching Huang } 397a13bc4e3SShao-Ching Huang } 398a13bc4e3SShao-Ching Huang } else if (r == rank) { 399a13bc4e3SShao-Ching Huang xm = info.xm; 400a13bc4e3SShao-Ching Huang ym = info.ym; 401a13bc4e3SShao-Ching Huang zm = info.zm; 402a13bc4e3SShao-Ching Huang for (j = 0, k = 0, i = 0; i < xm; i++) { 403a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 404a13bc4e3SShao-Ching Huang array2[i] = coords[Iloc * dim + 0]; 405a13bc4e3SShao-Ching Huang } 406a13bc4e3SShao-Ching Huang for (i = 0, k = 0, j = 0; j < ym; j++) { 407a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 408a13bc4e3SShao-Ching Huang array2[j + xm] = dim > 1 ? coords[Iloc * dim + 1] : 0; 409a13bc4e3SShao-Ching Huang } 410a13bc4e3SShao-Ching Huang for (i = 0, j = 0, k = 0; k < zm; k++) { 411a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 412a13bc4e3SShao-Ching Huang array2[k + xm + ym] = dim > 2 ? coords[Iloc * dim + 2] : 0; 413a13bc4e3SShao-Ching Huang } 4149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void *)array2, xm + ym + zm, MPIU_SCALAR, 0, tag, comm)); 415a13bc4e3SShao-Ching Huang } 4169566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Coords, &coords)); 417a13bc4e3SShao-Ching Huang } else { /* Fabricate some coordinates using grid index */ 4189371c9d4SSatish Balay for (i = 0; i < xm; i++) { array[i] = xs + i; } 4199371c9d4SSatish Balay for (j = 0; j < ym; j++) { array[j + xm] = ys + j; } 4209371c9d4SSatish Balay for (k = 0; k < zm; k++) { array[k + xm + ym] = zs + k; } 421a13bc4e3SShao-Ching Huang } 422dd400576SPatrick Sanan if (rank == 0) { 4239566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, &(array[0]), xm, MPIU_SCALAR)); 4249566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, &(array[xm]), ym, MPIU_SCALAR)); 4259566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, &(array[xm + ym]), zm, MPIU_SCALAR)); 426a13bc4e3SShao-Ching Huang } 427a13bc4e3SShao-Ching Huang } 428a13bc4e3SShao-Ching Huang 429a13bc4e3SShao-Ching Huang /* Write each of the objects queued up for this file */ 430a13bc4e3SShao-Ching Huang for (link = vtk->link; link; link = link->next) { 431a13bc4e3SShao-Ching Huang Vec X = (Vec)link->vec; 432a13bc4e3SShao-Ching Huang const PetscScalar *x; 433fb5bd1c2SPatrick Sanan PetscInt bs, f; 434ea2d7708SPatrick Sanan DM daCurr; 435fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 4369566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &daCurr)); 4379566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL)); 438a13bc4e3SShao-Ching Huang 4399566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 440dd400576SPatrick Sanan if (rank == 0) { 441a13bc4e3SShao-Ching Huang if (r) { 442a13bc4e3SShao-Ching Huang PetscMPIInt nn; 4439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array, nnodes * bs, MPIU_SCALAR, r, tag, comm, &status)); 4449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn)); 44563a3b9bcSJacob Faibussowitsch PetscCheck(nn == nnodes * bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch receiving from rank %" PetscInt_FMT, r); 4461baa6e33SBarry Smith } else PetscCall(PetscArraycpy(array, x, nnodes * bs)); 447fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 4489566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr, &fieldsnamed)); 449fb5bd1c2SPatrick Sanan if (fieldsnamed) { 450fb5bd1c2SPatrick Sanan for (f = 0; f < bs; f++) { 451fb5bd1c2SPatrick Sanan /* Extract and transpose the f'th field */ 452fb5bd1c2SPatrick Sanan for (k = 0; k < zm; k++) { 453fb5bd1c2SPatrick Sanan for (j = 0; j < ym; j++) { 454fb5bd1c2SPatrick Sanan for (i = 0; i < xm; i++) { 455fb5bd1c2SPatrick Sanan PetscInt Iloc = i + xm * (j + ym * k); 456fb5bd1c2SPatrick Sanan array2[Iloc] = array[Iloc * bs + f]; 457fb5bd1c2SPatrick Sanan } 458fb5bd1c2SPatrick Sanan } 459fb5bd1c2SPatrick Sanan } 4609566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, array2, nnodes, MPIU_SCALAR)); 461fb5bd1c2SPatrick Sanan } 462fb5bd1c2SPatrick Sanan } 4639566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, array, nnodes * bs, MPIU_SCALAR)); 464ea2d7708SPatrick Sanan 465a13bc4e3SShao-Ching Huang } else if (r == rank) { 4669566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void *)x, nnodes * bs, MPIU_SCALAR, 0, tag, comm)); 467a13bc4e3SShao-Ching Huang } 4689566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 469a13bc4e3SShao-Ching Huang } 470a13bc4e3SShao-Ching Huang } 4719566063dSJacob Faibussowitsch PetscCall(PetscFree2(array, array2)); 4729566063dSJacob Faibussowitsch PetscCall(PetscFree(grloc)); 473a13bc4e3SShao-Ching Huang 4749566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "\n </AppendedData>\n")); 4759566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "</VTKFile>\n")); 4769566063dSJacob Faibussowitsch PetscCall(PetscFClose(comm, fp)); 477a13bc4e3SShao-Ching Huang PetscFunctionReturn(0); 478a13bc4e3SShao-Ching Huang } 479a13bc4e3SShao-Ching Huang 4804061b8bfSJed Brown /*@C 4814061b8bfSJed Brown DMDAVTKWriteAll - Write a file containing all the fields that have been provided to the viewer 4824061b8bfSJed Brown 4834061b8bfSJed Brown Collective 4844061b8bfSJed Brown 4854165533cSJose E. Roman Input Parameters: 4864165533cSJose E. Roman + odm - DM specifying the grid layout, passed as a PetscObject 4874165533cSJose E. Roman - viewer - viewer of type VTK 4884061b8bfSJed Brown 4894061b8bfSJed Brown Level: developer 4904061b8bfSJed Brown 491fb5bd1c2SPatrick Sanan Notes: 4924061b8bfSJed Brown This function is a callback used by the VTK viewer to actually write the file. 4934061b8bfSJed 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. 4944061b8bfSJed Brown Instead, metadata for the entire file needs to be available up-front before you can start writing the file. 4954061b8bfSJed Brown 496fb5bd1c2SPatrick Sanan If any fields have been named (see e.g. DMDASetFieldName()), then individual scalar 497fb5bd1c2SPatrick Sanan fields are written. Otherwise, a single multi-dof (vector) field is written. 498fb5bd1c2SPatrick Sanan 499db781477SPatrick Sanan .seealso: `PETSCVIEWERVTK` 5004061b8bfSJed Brown @*/ 5019371c9d4SSatish Balay PetscErrorCode DMDAVTKWriteAll(PetscObject odm, PetscViewer viewer) { 5024061b8bfSJed Brown DM dm = (DM)odm; 5034061b8bfSJed Brown PetscBool isvtk; 5044061b8bfSJed Brown 5054061b8bfSJed Brown PetscFunctionBegin; 5064061b8bfSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5074061b8bfSJed Brown PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 5089566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk)); 5097a8be351SBarry Smith PetscCheck(isvtk, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name); 5104061b8bfSJed Brown switch (viewer->format) { 5119371c9d4SSatish Balay case PETSC_VIEWER_VTK_VTS: PetscCall(DMDAVTKWriteAll_VTS(dm, viewer)); break; 5129371c9d4SSatish Balay case PETSC_VIEWER_VTK_VTR: PetscCall(DMDAVTKWriteAll_VTR(dm, viewer)); break; 51398921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]); 5144061b8bfSJed Brown } 5154061b8bfSJed Brown PetscFunctionReturn(0); 5164061b8bfSJed Brown } 517