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 */ 11d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAGetFieldsNamed(DM da, PetscBool *fieldsnamed) 12d71ae5a4SJacob Faibussowitsch { 13fb5bd1c2SPatrick Sanan PetscInt f, bs; 14fb5bd1c2SPatrick Sanan 15fb5bd1c2SPatrick Sanan PetscFunctionBegin; 16fb5bd1c2SPatrick Sanan *fieldsnamed = PETSC_FALSE; 179566063dSJacob Faibussowitsch PetscCall(DMDAGetDof(da, &bs)); 18fb5bd1c2SPatrick Sanan for (f = 0; f < bs; ++f) { 19fb5bd1c2SPatrick Sanan const char *fieldname; 209566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da, f, &fieldname)); 21fb5bd1c2SPatrick Sanan if (fieldname) { 22fb5bd1c2SPatrick Sanan *fieldsnamed = PETSC_TRUE; 23fb5bd1c2SPatrick Sanan break; 24fb5bd1c2SPatrick Sanan } 25fb5bd1c2SPatrick Sanan } 26fb5bd1c2SPatrick Sanan PetscFunctionReturn(0); 27fb5bd1c2SPatrick Sanan } 28fb5bd1c2SPatrick Sanan 29d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAVTKWriteAll_VTS(DM da, PetscViewer viewer) 30d71ae5a4SJacob Faibussowitsch { 3130815ce0SLisandro Dalcin const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian"; 324061b8bfSJed Brown #if defined(PETSC_USE_REAL_SINGLE) 334061b8bfSJed Brown const char precision[] = "Float32"; 344061b8bfSJed Brown #elif defined(PETSC_USE_REAL_DOUBLE) 354061b8bfSJed Brown const char precision[] = "Float64"; 364061b8bfSJed Brown #else 374061b8bfSJed Brown const char precision[] = "UnknownPrecision"; 384061b8bfSJed Brown #endif 39ce94432eSBarry Smith MPI_Comm comm; 402eaa9ef4SLisandro Dalcin Vec Coords; 414061b8bfSJed Brown PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data; 424061b8bfSJed Brown PetscViewerVTKObjectLink link; 434061b8bfSJed Brown FILE *fp; 444061b8bfSJed Brown PetscMPIInt rank, size, tag; 454061b8bfSJed Brown DMDALocalInfo info; 46ea2d7708SPatrick Sanan PetscInt dim, mx, my, mz, cdim, bs, boffset, maxnnodes, maxbs, i, j, k, r; 470298fd71SBarry Smith PetscInt rloc[6], (*grloc)[6] = NULL; 484061b8bfSJed Brown PetscScalar *array, *array2; 494061b8bfSJed Brown 504061b8bfSJed Brown PetscFunctionBegin; 519566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 521dca8a05SBarry Smith PetscCheck(!PetscDefined(USE_COMPLEX), PETSC_COMM_SELF, PETSC_ERR_SUP, "Complex values not supported"); 539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 559566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, &mx, &my, &mz, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL)); 569566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 579566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &Coords)); 582eaa9ef4SLisandro Dalcin if (Coords) { 592eaa9ef4SLisandro Dalcin PetscInt csize; 609566063dSJacob Faibussowitsch PetscCall(VecGetSize(Coords, &csize)); 6108401ef6SPierre Jolivet PetscCheck(csize % (mx * my * mz) == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Coordinate vector size mismatch"); 622eaa9ef4SLisandro Dalcin cdim = csize / (mx * my * mz); 631dca8a05SBarry Smith PetscCheck(cdim >= dim && cdim <= 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Coordinate vector size mismatch"); 642eaa9ef4SLisandro Dalcin } else { 652eaa9ef4SLisandro Dalcin cdim = dim; 662eaa9ef4SLisandro Dalcin } 674061b8bfSJed Brown 689566063dSJacob Faibussowitsch PetscCall(PetscFOpen(comm, vtk->filename, "wb", &fp)); 699566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "<?xml version=\"1.0\"?>\n")); 709566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"%s\">\n", byte_order)); 7163a3b9bcSJacob 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)); 724061b8bfSJed Brown 739566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc1(size * 6, &grloc)); 744061b8bfSJed Brown rloc[0] = info.xs; 754061b8bfSJed Brown rloc[1] = info.xm; 764061b8bfSJed Brown rloc[2] = info.ys; 774061b8bfSJed Brown rloc[3] = info.ym; 784061b8bfSJed Brown rloc[4] = info.zs; 794061b8bfSJed Brown rloc[5] = info.zm; 809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gather(rloc, 6, MPIU_INT, &grloc[0][0], 6, MPIU_INT, 0, comm)); 814061b8bfSJed Brown 824061b8bfSJed Brown /* Write XML header */ 834061b8bfSJed Brown maxnnodes = 0; /* Used for the temporary array size on rank 0 */ 84ea2d7708SPatrick Sanan maxbs = 0; /* Used for the temporary array size on rank 0 */ 854061b8bfSJed Brown boffset = 0; /* Offset into binary file */ 864061b8bfSJed Brown for (r = 0; r < size; r++) { 874061b8bfSJed Brown PetscInt xs = -1, xm = -1, ys = -1, ym = -1, zs = -1, zm = -1, nnodes = 0; 88dd400576SPatrick Sanan if (rank == 0) { 894061b8bfSJed Brown xs = grloc[r][0]; 904061b8bfSJed Brown xm = grloc[r][1]; 914061b8bfSJed Brown ys = grloc[r][2]; 924061b8bfSJed Brown ym = grloc[r][3]; 934061b8bfSJed Brown zs = grloc[r][4]; 944061b8bfSJed Brown zm = grloc[r][5]; 954061b8bfSJed Brown nnodes = xm * ym * zm; 964061b8bfSJed Brown } 974061b8bfSJed Brown maxnnodes = PetscMax(maxnnodes, nnodes); 9863a3b9bcSJacob 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)); 999566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <Points>\n")); 10063a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, boffset)); 1014061b8bfSJed Brown boffset += 3 * nnodes * sizeof(PetscScalar) + sizeof(int); 1029566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </Points>\n")); 1034061b8bfSJed Brown 1049566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <PointData Scalars=\"ScalarPointData\">\n")); 1054061b8bfSJed Brown for (link = vtk->link; link; link = link->next) { 1064061b8bfSJed Brown Vec X = (Vec)link->vec; 107fb5bd1c2SPatrick Sanan PetscInt bs, f; 108ea2d7708SPatrick Sanan DM daCurr; 109fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 110fb5bd1c2SPatrick Sanan const char *vecname = "Unnamed"; 111ea2d7708SPatrick Sanan 1129566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &daCurr)); 1139566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL)); 114ea2d7708SPatrick Sanan maxbs = PetscMax(maxbs, bs); 115ea2d7708SPatrick Sanan 11648a46eb9SPierre Jolivet if (((PetscObject)X)->name || link != vtk->link) PetscCall(PetscObjectGetName((PetscObject)X, &vecname)); 117fb5bd1c2SPatrick Sanan 118fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 1199566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr, &fieldsnamed)); 120fb5bd1c2SPatrick Sanan if (fieldsnamed) { 121fb5bd1c2SPatrick Sanan for (f = 0; f < bs; f++) { 122fb5bd1c2SPatrick Sanan char buf[256]; 123fb5bd1c2SPatrick Sanan const char *fieldname; 1249566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(daCurr, f, &fieldname)); 125fb5bd1c2SPatrick Sanan if (!fieldname) { 12663a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, f)); 127fb5bd1c2SPatrick Sanan fieldname = buf; 128fb5bd1c2SPatrick Sanan } 12963a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, vecname, fieldname, boffset)); 130fb5bd1c2SPatrick Sanan boffset += nnodes * sizeof(PetscScalar) + sizeof(int); 131fb5bd1c2SPatrick Sanan } 132fb5bd1c2SPatrick Sanan } else { 13363a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%" PetscInt_FMT "\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, vecname, bs, boffset)); 134ea2d7708SPatrick Sanan boffset += bs * nnodes * sizeof(PetscScalar) + sizeof(int); 1354061b8bfSJed Brown } 136fb5bd1c2SPatrick Sanan } 1379566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </PointData>\n")); 1389566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </Piece>\n")); 1394061b8bfSJed Brown } 1409566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </StructuredGrid>\n")); 1419566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <AppendedData encoding=\"raw\">\n")); 1429566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "_")); 1434061b8bfSJed Brown 1444061b8bfSJed Brown /* Now write the arrays. */ 1454061b8bfSJed Brown tag = ((PetscObject)viewer)->tag; 1469566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(maxnnodes * PetscMax(3, maxbs), &array, maxnnodes * PetscMax(3, maxbs), &array2)); 1474061b8bfSJed Brown for (r = 0; r < size; r++) { 1484061b8bfSJed Brown MPI_Status status; 1494061b8bfSJed Brown PetscInt xs = -1, xm = -1, ys = -1, ym = -1, zs = -1, zm = -1, nnodes = 0; 150dd400576SPatrick Sanan if (rank == 0) { 1514061b8bfSJed Brown xs = grloc[r][0]; 1524061b8bfSJed Brown xm = grloc[r][1]; 1534061b8bfSJed Brown ys = grloc[r][2]; 1544061b8bfSJed Brown ym = grloc[r][3]; 1554061b8bfSJed Brown zs = grloc[r][4]; 1564061b8bfSJed Brown zm = grloc[r][5]; 1574061b8bfSJed Brown nnodes = xm * ym * zm; 1584061b8bfSJed Brown } else if (r == rank) { 1594061b8bfSJed Brown nnodes = info.xm * info.ym * info.zm; 1604061b8bfSJed Brown } 1614061b8bfSJed Brown 1622eaa9ef4SLisandro Dalcin /* Write the coordinates */ 1634061b8bfSJed Brown if (Coords) { 1644061b8bfSJed Brown const PetscScalar *coords; 1659566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Coords, &coords)); 166dd400576SPatrick Sanan if (rank == 0) { 1674061b8bfSJed Brown if (r) { 1686622f924SJed Brown PetscMPIInt nn; 1699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array, nnodes * cdim, MPIU_SCALAR, r, tag, comm, &status)); 1709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn)); 17108401ef6SPierre Jolivet PetscCheck(nn == nnodes * cdim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch"); 1721baa6e33SBarry Smith } else PetscCall(PetscArraycpy(array, coords, nnodes * cdim)); 1734061b8bfSJed Brown /* Transpose coordinates to VTK (C-style) ordering */ 1744061b8bfSJed Brown for (k = 0; k < zm; k++) { 1754061b8bfSJed Brown for (j = 0; j < ym; j++) { 1764061b8bfSJed Brown for (i = 0; i < xm; i++) { 1774061b8bfSJed Brown PetscInt Iloc = i + xm * (j + ym * k); 1782eaa9ef4SLisandro Dalcin array2[Iloc * 3 + 0] = array[Iloc * cdim + 0]; 1792eaa9ef4SLisandro Dalcin array2[Iloc * 3 + 1] = cdim > 1 ? array[Iloc * cdim + 1] : 0.0; 1802eaa9ef4SLisandro Dalcin array2[Iloc * 3 + 2] = cdim > 2 ? array[Iloc * cdim + 2] : 0.0; 1814061b8bfSJed Brown } 1824061b8bfSJed Brown } 1834061b8bfSJed Brown } 1844061b8bfSJed Brown } else if (r == rank) { 1859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void *)coords, nnodes * cdim, MPIU_SCALAR, 0, tag, comm)); 1864061b8bfSJed Brown } 1879566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Coords, &coords)); 1884061b8bfSJed Brown } else { /* Fabricate some coordinates using grid index */ 1894061b8bfSJed Brown for (k = 0; k < zm; k++) { 1904061b8bfSJed Brown for (j = 0; j < ym; j++) { 1914061b8bfSJed Brown for (i = 0; i < xm; i++) { 1924061b8bfSJed Brown PetscInt Iloc = i + xm * (j + ym * k); 1934061b8bfSJed Brown array2[Iloc * 3 + 0] = xs + i; 1944061b8bfSJed Brown array2[Iloc * 3 + 1] = ys + j; 1954061b8bfSJed Brown array2[Iloc * 3 + 2] = zs + k; 1964061b8bfSJed Brown } 1974061b8bfSJed Brown } 1984061b8bfSJed Brown } 1994061b8bfSJed Brown } 2009566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, array2, nnodes * 3, MPIU_SCALAR)); 2014061b8bfSJed Brown 2024061b8bfSJed Brown /* Write each of the objects queued up for this file */ 2034061b8bfSJed Brown for (link = vtk->link; link; link = link->next) { 2044061b8bfSJed Brown Vec X = (Vec)link->vec; 2054061b8bfSJed Brown const PetscScalar *x; 206fb5bd1c2SPatrick Sanan PetscInt bs, f; 207ea2d7708SPatrick Sanan DM daCurr; 208fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 2099566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &daCurr)); 2109566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL)); 2119566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 212dd400576SPatrick Sanan if (rank == 0) { 2134061b8bfSJed Brown if (r) { 2146622f924SJed Brown PetscMPIInt nn; 2159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array, nnodes * bs, MPIU_SCALAR, r, tag, comm, &status)); 2169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn)); 21763a3b9bcSJacob Faibussowitsch PetscCheck(nn == nnodes * bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch receiving from rank %" PetscInt_FMT, r); 2181baa6e33SBarry Smith } else PetscCall(PetscArraycpy(array, x, nnodes * bs)); 219fb5bd1c2SPatrick Sanan 220fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 2219566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr, &fieldsnamed)); 222fb5bd1c2SPatrick Sanan if (fieldsnamed) { 223fb5bd1c2SPatrick Sanan for (f = 0; f < bs; f++) { 224fb5bd1c2SPatrick Sanan /* Extract and transpose the f'th field */ 225fb5bd1c2SPatrick Sanan for (k = 0; k < zm; k++) { 226fb5bd1c2SPatrick Sanan for (j = 0; j < ym; j++) { 227fb5bd1c2SPatrick Sanan for (i = 0; i < xm; i++) { 228fb5bd1c2SPatrick Sanan PetscInt Iloc = i + xm * (j + ym * k); 229fb5bd1c2SPatrick Sanan array2[Iloc] = array[Iloc * bs + f]; 230fb5bd1c2SPatrick Sanan } 231fb5bd1c2SPatrick Sanan } 232fb5bd1c2SPatrick Sanan } 2339566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, array2, nnodes, MPIU_SCALAR)); 234fb5bd1c2SPatrick Sanan } 2351baa6e33SBarry Smith } else PetscCall(PetscViewerVTKFWrite(viewer, fp, array, bs * nnodes, MPIU_SCALAR)); 2361baa6e33SBarry Smith } else if (r == rank) PetscCallMPI(MPI_Send((void *)x, nnodes * bs, MPIU_SCALAR, 0, tag, comm)); 2379566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 2384061b8bfSJed Brown } 2394061b8bfSJed Brown } 2409566063dSJacob Faibussowitsch PetscCall(PetscFree2(array, array2)); 2419566063dSJacob Faibussowitsch PetscCall(PetscFree(grloc)); 2424061b8bfSJed Brown 2439566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "\n </AppendedData>\n")); 2449566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "</VTKFile>\n")); 2459566063dSJacob Faibussowitsch PetscCall(PetscFClose(comm, fp)); 2464061b8bfSJed Brown PetscFunctionReturn(0); 2474061b8bfSJed Brown } 2484061b8bfSJed Brown 249d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAVTKWriteAll_VTR(DM da, PetscViewer viewer) 250d71ae5a4SJacob Faibussowitsch { 25130815ce0SLisandro Dalcin const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian"; 252a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_REAL_SINGLE) 253a13bc4e3SShao-Ching Huang const char precision[] = "Float32"; 254a13bc4e3SShao-Ching Huang #elif defined(PETSC_USE_REAL_DOUBLE) 255a13bc4e3SShao-Ching Huang const char precision[] = "Float64"; 256a13bc4e3SShao-Ching Huang #else 257a13bc4e3SShao-Ching Huang const char precision[] = "UnknownPrecision"; 258a13bc4e3SShao-Ching Huang #endif 259a13bc4e3SShao-Ching Huang MPI_Comm comm; 260a13bc4e3SShao-Ching Huang PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data; 261a13bc4e3SShao-Ching Huang PetscViewerVTKObjectLink link; 262a13bc4e3SShao-Ching Huang FILE *fp; 263a13bc4e3SShao-Ching Huang PetscMPIInt rank, size, tag; 264a13bc4e3SShao-Ching Huang DMDALocalInfo info; 265ea2d7708SPatrick Sanan PetscInt dim, mx, my, mz, boffset, maxnnodes, maxbs, i, j, k, r; 266a13bc4e3SShao-Ching Huang PetscInt rloc[6], (*grloc)[6] = NULL; 267a13bc4e3SShao-Ching Huang PetscScalar *array, *array2; 268a13bc4e3SShao-Ching Huang 269a13bc4e3SShao-Ching Huang PetscFunctionBegin; 2709566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 2711dca8a05SBarry Smith PetscCheck(!PetscDefined(USE_COMPLEX), PETSC_COMM_SELF, PETSC_ERR_SUP, "Complex values not supported"); 2729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2749566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, &mx, &my, &mz, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2759566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 2769566063dSJacob Faibussowitsch PetscCall(PetscFOpen(comm, vtk->filename, "wb", &fp)); 2779566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "<?xml version=\"1.0\"?>\n")); 2789566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "<VTKFile type=\"RectilinearGrid\" version=\"0.1\" byte_order=\"%s\">\n", byte_order)); 27963a3b9bcSJacob 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)); 280a13bc4e3SShao-Ching Huang 2819566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc1(size * 6, &grloc)); 282a13bc4e3SShao-Ching Huang rloc[0] = info.xs; 283a13bc4e3SShao-Ching Huang rloc[1] = info.xm; 284a13bc4e3SShao-Ching Huang rloc[2] = info.ys; 285a13bc4e3SShao-Ching Huang rloc[3] = info.ym; 286a13bc4e3SShao-Ching Huang rloc[4] = info.zs; 287a13bc4e3SShao-Ching Huang rloc[5] = info.zm; 2889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gather(rloc, 6, MPIU_INT, &grloc[0][0], 6, MPIU_INT, 0, comm)); 289a13bc4e3SShao-Ching Huang 290a13bc4e3SShao-Ching Huang /* Write XML header */ 291a13bc4e3SShao-Ching Huang maxnnodes = 0; /* Used for the temporary array size on rank 0 */ 292ea2d7708SPatrick Sanan maxbs = 0; /* Used for the temporary array size on rank 0 */ 293a13bc4e3SShao-Ching Huang boffset = 0; /* Offset into binary file */ 294a13bc4e3SShao-Ching Huang for (r = 0; r < size; r++) { 295a13bc4e3SShao-Ching Huang PetscInt xs = -1, xm = -1, ys = -1, ym = -1, zs = -1, zm = -1, nnodes = 0; 296dd400576SPatrick Sanan if (rank == 0) { 297a13bc4e3SShao-Ching Huang xs = grloc[r][0]; 298a13bc4e3SShao-Ching Huang xm = grloc[r][1]; 299a13bc4e3SShao-Ching Huang ys = grloc[r][2]; 300a13bc4e3SShao-Ching Huang ym = grloc[r][3]; 301a13bc4e3SShao-Ching Huang zs = grloc[r][4]; 302a13bc4e3SShao-Ching Huang zm = grloc[r][5]; 303a13bc4e3SShao-Ching Huang nnodes = xm * ym * zm; 304a13bc4e3SShao-Ching Huang } 305a13bc4e3SShao-Ching Huang maxnnodes = PetscMax(maxnnodes, nnodes); 30663a3b9bcSJacob 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)); 3079566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <Coordinates>\n")); 30863a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"Xcoord\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, boffset)); 309a13bc4e3SShao-Ching Huang boffset += xm * sizeof(PetscScalar) + sizeof(int); 31063a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"Ycoord\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, boffset)); 311a13bc4e3SShao-Ching Huang boffset += ym * sizeof(PetscScalar) + sizeof(int); 31263a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"Zcoord\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, boffset)); 313a13bc4e3SShao-Ching Huang boffset += zm * sizeof(PetscScalar) + sizeof(int); 3149566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </Coordinates>\n")); 3159566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <PointData Scalars=\"ScalarPointData\">\n")); 316a13bc4e3SShao-Ching Huang for (link = vtk->link; link; link = link->next) { 317a13bc4e3SShao-Ching Huang Vec X = (Vec)link->vec; 318fb5bd1c2SPatrick Sanan PetscInt bs, f; 319ea2d7708SPatrick Sanan DM daCurr; 320fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 321fb5bd1c2SPatrick Sanan const char *vecname = "Unnamed"; 322ea2d7708SPatrick Sanan 3239566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &daCurr)); 3249566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL)); 325ea2d7708SPatrick Sanan maxbs = PetscMax(maxbs, bs); 32648a46eb9SPierre Jolivet if (((PetscObject)X)->name || link != vtk->link) PetscCall(PetscObjectGetName((PetscObject)X, &vecname)); 327ea2d7708SPatrick Sanan 328fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 3299566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr, &fieldsnamed)); 330fb5bd1c2SPatrick Sanan if (fieldsnamed) { 331fb5bd1c2SPatrick Sanan for (f = 0; f < bs; f++) { 332fb5bd1c2SPatrick Sanan char buf[256]; 333fb5bd1c2SPatrick Sanan const char *fieldname; 3349566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(daCurr, f, &fieldname)); 335fb5bd1c2SPatrick Sanan if (!fieldname) { 33663a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, f)); 337fb5bd1c2SPatrick Sanan fieldname = buf; 338fb5bd1c2SPatrick Sanan } 33963a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, vecname, fieldname, boffset)); 340fb5bd1c2SPatrick Sanan boffset += nnodes * sizeof(PetscScalar) + sizeof(int); 341fb5bd1c2SPatrick Sanan } 342fb5bd1c2SPatrick Sanan } else { 34363a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%" PetscInt_FMT "\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n", precision, vecname, bs, boffset)); 344ea2d7708SPatrick Sanan boffset += bs * nnodes * sizeof(PetscScalar) + sizeof(int); 345a13bc4e3SShao-Ching Huang } 346fb5bd1c2SPatrick Sanan } 3479566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </PointData>\n")); 3489566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </Piece>\n")); 349a13bc4e3SShao-Ching Huang } 3509566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </RectilinearGrid>\n")); 3519566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <AppendedData encoding=\"raw\">\n")); 3529566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "_")); 353a13bc4e3SShao-Ching Huang 354a13bc4e3SShao-Ching Huang /* Now write the arrays. */ 355a13bc4e3SShao-Ching Huang tag = ((PetscObject)viewer)->tag; 3569566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(PetscMax(maxnnodes * maxbs, info.xm + info.ym + info.zm), &array, PetscMax(maxnnodes * maxbs, info.xm + info.ym + info.zm), &array2)); 357a13bc4e3SShao-Ching Huang for (r = 0; r < size; r++) { 358a13bc4e3SShao-Ching Huang MPI_Status status; 359a13bc4e3SShao-Ching Huang PetscInt xs = -1, xm = -1, ys = -1, ym = -1, zs = -1, zm = -1, nnodes = 0; 360dd400576SPatrick Sanan if (rank == 0) { 361a13bc4e3SShao-Ching Huang xs = grloc[r][0]; 362a13bc4e3SShao-Ching Huang xm = grloc[r][1]; 363a13bc4e3SShao-Ching Huang ys = grloc[r][2]; 364a13bc4e3SShao-Ching Huang ym = grloc[r][3]; 365a13bc4e3SShao-Ching Huang zs = grloc[r][4]; 366a13bc4e3SShao-Ching Huang zm = grloc[r][5]; 367a13bc4e3SShao-Ching Huang nnodes = xm * ym * zm; 368a13bc4e3SShao-Ching Huang } else if (r == rank) { 369a13bc4e3SShao-Ching Huang nnodes = info.xm * info.ym * info.zm; 370a13bc4e3SShao-Ching Huang } 371a13bc4e3SShao-Ching Huang { /* Write the coordinates */ 372a13bc4e3SShao-Ching Huang Vec Coords; 373ea2d7708SPatrick Sanan 3749566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &Coords)); 375a13bc4e3SShao-Ching Huang if (Coords) { 376a13bc4e3SShao-Ching Huang const PetscScalar *coords; 3779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Coords, &coords)); 378dd400576SPatrick Sanan if (rank == 0) { 379a13bc4e3SShao-Ching Huang if (r) { 380a13bc4e3SShao-Ching Huang PetscMPIInt nn; 3819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array, xm + ym + zm, MPIU_SCALAR, r, tag, comm, &status)); 3829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn)); 38308401ef6SPierre Jolivet PetscCheck(nn == xm + ym + zm, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch"); 384a13bc4e3SShao-Ching Huang } else { 385a13bc4e3SShao-Ching Huang /* x coordinates */ 386a13bc4e3SShao-Ching Huang for (j = 0, k = 0, i = 0; i < xm; i++) { 387a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 388a13bc4e3SShao-Ching Huang array[i] = coords[Iloc * dim + 0]; 389a13bc4e3SShao-Ching Huang } 390a13bc4e3SShao-Ching Huang /* y coordinates */ 391a13bc4e3SShao-Ching Huang for (i = 0, k = 0, j = 0; j < ym; j++) { 392a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 393a13bc4e3SShao-Ching Huang array[j + xm] = dim > 1 ? coords[Iloc * dim + 1] : 0; 394a13bc4e3SShao-Ching Huang } 395a13bc4e3SShao-Ching Huang /* z coordinates */ 396a13bc4e3SShao-Ching Huang for (i = 0, j = 0, k = 0; k < zm; k++) { 397a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 398a13bc4e3SShao-Ching Huang array[k + xm + ym] = dim > 2 ? coords[Iloc * dim + 2] : 0; 399a13bc4e3SShao-Ching Huang } 400a13bc4e3SShao-Ching Huang } 401a13bc4e3SShao-Ching Huang } else if (r == rank) { 402a13bc4e3SShao-Ching Huang xm = info.xm; 403a13bc4e3SShao-Ching Huang ym = info.ym; 404a13bc4e3SShao-Ching Huang zm = info.zm; 405a13bc4e3SShao-Ching Huang for (j = 0, k = 0, i = 0; i < xm; i++) { 406a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 407a13bc4e3SShao-Ching Huang array2[i] = coords[Iloc * dim + 0]; 408a13bc4e3SShao-Ching Huang } 409a13bc4e3SShao-Ching Huang for (i = 0, k = 0, j = 0; j < ym; j++) { 410a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 411a13bc4e3SShao-Ching Huang array2[j + xm] = dim > 1 ? coords[Iloc * dim + 1] : 0; 412a13bc4e3SShao-Ching Huang } 413a13bc4e3SShao-Ching Huang for (i = 0, j = 0, k = 0; k < zm; k++) { 414a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 415a13bc4e3SShao-Ching Huang array2[k + xm + ym] = dim > 2 ? coords[Iloc * dim + 2] : 0; 416a13bc4e3SShao-Ching Huang } 4179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void *)array2, xm + ym + zm, MPIU_SCALAR, 0, tag, comm)); 418a13bc4e3SShao-Ching Huang } 4199566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Coords, &coords)); 420a13bc4e3SShao-Ching Huang } else { /* Fabricate some coordinates using grid index */ 421ad540459SPierre Jolivet for (i = 0; i < xm; i++) array[i] = xs + i; 422ad540459SPierre Jolivet for (j = 0; j < ym; j++) array[j + xm] = ys + j; 423ad540459SPierre Jolivet for (k = 0; k < zm; k++) array[k + xm + ym] = zs + k; 424a13bc4e3SShao-Ching Huang } 425dd400576SPatrick Sanan if (rank == 0) { 4269566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, &(array[0]), xm, MPIU_SCALAR)); 4279566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, &(array[xm]), ym, MPIU_SCALAR)); 4289566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, &(array[xm + ym]), zm, MPIU_SCALAR)); 429a13bc4e3SShao-Ching Huang } 430a13bc4e3SShao-Ching Huang } 431a13bc4e3SShao-Ching Huang 432a13bc4e3SShao-Ching Huang /* Write each of the objects queued up for this file */ 433a13bc4e3SShao-Ching Huang for (link = vtk->link; link; link = link->next) { 434a13bc4e3SShao-Ching Huang Vec X = (Vec)link->vec; 435a13bc4e3SShao-Ching Huang const PetscScalar *x; 436fb5bd1c2SPatrick Sanan PetscInt bs, f; 437ea2d7708SPatrick Sanan DM daCurr; 438fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 4399566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &daCurr)); 4409566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL)); 441a13bc4e3SShao-Ching Huang 4429566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 443dd400576SPatrick Sanan if (rank == 0) { 444a13bc4e3SShao-Ching Huang if (r) { 445a13bc4e3SShao-Ching Huang PetscMPIInt nn; 4469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array, nnodes * bs, MPIU_SCALAR, r, tag, comm, &status)); 4479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn)); 44863a3b9bcSJacob Faibussowitsch PetscCheck(nn == nnodes * bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch receiving from rank %" PetscInt_FMT, r); 4491baa6e33SBarry Smith } else PetscCall(PetscArraycpy(array, x, nnodes * bs)); 450fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 4519566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr, &fieldsnamed)); 452fb5bd1c2SPatrick Sanan if (fieldsnamed) { 453fb5bd1c2SPatrick Sanan for (f = 0; f < bs; f++) { 454fb5bd1c2SPatrick Sanan /* Extract and transpose the f'th field */ 455fb5bd1c2SPatrick Sanan for (k = 0; k < zm; k++) { 456fb5bd1c2SPatrick Sanan for (j = 0; j < ym; j++) { 457fb5bd1c2SPatrick Sanan for (i = 0; i < xm; i++) { 458fb5bd1c2SPatrick Sanan PetscInt Iloc = i + xm * (j + ym * k); 459fb5bd1c2SPatrick Sanan array2[Iloc] = array[Iloc * bs + f]; 460fb5bd1c2SPatrick Sanan } 461fb5bd1c2SPatrick Sanan } 462fb5bd1c2SPatrick Sanan } 4639566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, array2, nnodes, MPIU_SCALAR)); 464fb5bd1c2SPatrick Sanan } 465fb5bd1c2SPatrick Sanan } 4669566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, array, nnodes * bs, MPIU_SCALAR)); 467ea2d7708SPatrick Sanan 468a13bc4e3SShao-Ching Huang } else if (r == rank) { 4699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void *)x, nnodes * bs, MPIU_SCALAR, 0, tag, comm)); 470a13bc4e3SShao-Ching Huang } 4719566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 472a13bc4e3SShao-Ching Huang } 473a13bc4e3SShao-Ching Huang } 4749566063dSJacob Faibussowitsch PetscCall(PetscFree2(array, array2)); 4759566063dSJacob Faibussowitsch PetscCall(PetscFree(grloc)); 476a13bc4e3SShao-Ching Huang 4779566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "\n </AppendedData>\n")); 4789566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "</VTKFile>\n")); 4799566063dSJacob Faibussowitsch PetscCall(PetscFClose(comm, fp)); 480a13bc4e3SShao-Ching Huang PetscFunctionReturn(0); 481a13bc4e3SShao-Ching Huang } 482a13bc4e3SShao-Ching Huang 4834061b8bfSJed Brown /*@C 4844061b8bfSJed Brown DMDAVTKWriteAll - Write a file containing all the fields that have been provided to the viewer 4854061b8bfSJed Brown 4864061b8bfSJed Brown Collective 4874061b8bfSJed Brown 4884165533cSJose E. Roman Input Parameters: 489*dce8aebaSBarry Smith + odm - `DMDA` specifying the grid layout, passed as a `PetscObject` 490*dce8aebaSBarry Smith - viewer - viewer of type `PETSCVIEWERVTK` 4914061b8bfSJed Brown 4924061b8bfSJed Brown Level: developer 4934061b8bfSJed Brown 494fb5bd1c2SPatrick Sanan Notes: 4954061b8bfSJed Brown This function is a callback used by the VTK viewer to actually write the file. 4964061b8bfSJed 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. 4974061b8bfSJed Brown Instead, metadata for the entire file needs to be available up-front before you can start writing the file. 4984061b8bfSJed Brown 499fb5bd1c2SPatrick Sanan If any fields have been named (see e.g. DMDASetFieldName()), then individual scalar 500fb5bd1c2SPatrick Sanan fields are written. Otherwise, a single multi-dof (vector) field is written. 501fb5bd1c2SPatrick Sanan 502*dce8aebaSBarry Smith .seealso: `DMDA`, `DM`, `PETSCVIEWERVTK` 5034061b8bfSJed Brown @*/ 504d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVTKWriteAll(PetscObject odm, PetscViewer viewer) 505d71ae5a4SJacob Faibussowitsch { 5064061b8bfSJed Brown DM dm = (DM)odm; 5074061b8bfSJed Brown PetscBool isvtk; 5084061b8bfSJed Brown 5094061b8bfSJed Brown PetscFunctionBegin; 5104061b8bfSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5114061b8bfSJed Brown PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 5129566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk)); 5137a8be351SBarry Smith PetscCheck(isvtk, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name); 5144061b8bfSJed Brown switch (viewer->format) { 515d71ae5a4SJacob Faibussowitsch case PETSC_VIEWER_VTK_VTS: 516d71ae5a4SJacob Faibussowitsch PetscCall(DMDAVTKWriteAll_VTS(dm, viewer)); 517d71ae5a4SJacob Faibussowitsch break; 518d71ae5a4SJacob Faibussowitsch case PETSC_VIEWER_VTK_VTR: 519d71ae5a4SJacob Faibussowitsch PetscCall(DMDAVTKWriteAll_VTR(dm, viewer)); 520d71ae5a4SJacob Faibussowitsch break; 521d71ae5a4SJacob Faibussowitsch default: 522d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]); 5234061b8bfSJed Brown } 5244061b8bfSJed Brown PetscFunctionReturn(0); 5254061b8bfSJed Brown } 526