1*ffeef943SBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 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 } 263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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; 460f2609c8SJed Brown PetscInt dim, mx, my, mz, cdim, bs, maxnnodes, maxbs, i, j, k, r; 470f2609c8SJed Brown PetscInt64 boffset; 480298fd71SBarry Smith PetscInt rloc[6], (*grloc)[6] = NULL; 494061b8bfSJed Brown PetscScalar *array, *array2; 504061b8bfSJed Brown 514061b8bfSJed Brown PetscFunctionBegin; 529566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 531dca8a05SBarry Smith PetscCheck(!PetscDefined(USE_COMPLEX), PETSC_COMM_SELF, PETSC_ERR_SUP, "Complex values not supported"); 549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 569566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, &mx, &my, &mz, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL)); 579566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 589566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &Coords)); 592eaa9ef4SLisandro Dalcin if (Coords) { 602eaa9ef4SLisandro Dalcin PetscInt csize; 619566063dSJacob Faibussowitsch PetscCall(VecGetSize(Coords, &csize)); 6208401ef6SPierre Jolivet PetscCheck(csize % (mx * my * mz) == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Coordinate vector size mismatch"); 632eaa9ef4SLisandro Dalcin cdim = csize / (mx * my * mz); 641dca8a05SBarry Smith PetscCheck(cdim >= dim && cdim <= 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Coordinate vector size mismatch"); 652eaa9ef4SLisandro Dalcin } else { 662eaa9ef4SLisandro Dalcin cdim = dim; 672eaa9ef4SLisandro Dalcin } 684061b8bfSJed Brown 699566063dSJacob Faibussowitsch PetscCall(PetscFOpen(comm, vtk->filename, "wb", &fp)); 709566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "<?xml version=\"1.0\"?>\n")); 714a8c153fSJed Brown PetscCall(PetscFPrintf(comm, fp, "<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"%s\" header_type=\"UInt64\">\n", byte_order)); 7263a3b9bcSJacob 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)); 734061b8bfSJed Brown 749566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc1(size * 6, &grloc)); 754061b8bfSJed Brown rloc[0] = info.xs; 764061b8bfSJed Brown rloc[1] = info.xm; 774061b8bfSJed Brown rloc[2] = info.ys; 784061b8bfSJed Brown rloc[3] = info.ym; 794061b8bfSJed Brown rloc[4] = info.zs; 804061b8bfSJed Brown rloc[5] = info.zm; 81810441c8SPierre Jolivet PetscCallMPI(MPI_Gather(rloc, 6, MPIU_INT, rank == 0 ? grloc[0] : NULL, 6, MPIU_INT, 0, comm)); 824061b8bfSJed Brown 834061b8bfSJed Brown /* Write XML header */ 844061b8bfSJed Brown maxnnodes = 0; /* Used for the temporary array size on rank 0 */ 85ea2d7708SPatrick Sanan maxbs = 0; /* Used for the temporary array size on rank 0 */ 864061b8bfSJed Brown boffset = 0; /* Offset into binary file */ 874061b8bfSJed Brown for (r = 0; r < size; r++) { 884061b8bfSJed Brown PetscInt xs = -1, xm = -1, ys = -1, ym = -1, zs = -1, zm = -1, nnodes = 0; 89dd400576SPatrick Sanan if (rank == 0) { 904061b8bfSJed Brown xs = grloc[r][0]; 914061b8bfSJed Brown xm = grloc[r][1]; 924061b8bfSJed Brown ys = grloc[r][2]; 934061b8bfSJed Brown ym = grloc[r][3]; 944061b8bfSJed Brown zs = grloc[r][4]; 954061b8bfSJed Brown zm = grloc[r][5]; 964061b8bfSJed Brown nnodes = xm * ym * zm; 974061b8bfSJed Brown } 984061b8bfSJed Brown maxnnodes = PetscMax(maxnnodes, nnodes); 9963a3b9bcSJacob 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)); 1009566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <Points>\n")); 1010f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, boffset)); 1024a8c153fSJed Brown boffset += 3 * nnodes * sizeof(PetscScalar) + sizeof(PetscInt64); 1039566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </Points>\n")); 1044061b8bfSJed Brown 1059566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <PointData Scalars=\"ScalarPointData\">\n")); 1064061b8bfSJed Brown for (link = vtk->link; link; link = link->next) { 1074061b8bfSJed Brown Vec X = (Vec)link->vec; 108fb5bd1c2SPatrick Sanan PetscInt bs, f; 109ea2d7708SPatrick Sanan DM daCurr; 110fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 111fb5bd1c2SPatrick Sanan const char *vecname = "Unnamed"; 112ea2d7708SPatrick Sanan 1139566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &daCurr)); 1149566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL)); 115ea2d7708SPatrick Sanan maxbs = PetscMax(maxbs, bs); 116ea2d7708SPatrick Sanan 11748a46eb9SPierre Jolivet if (((PetscObject)X)->name || link != vtk->link) PetscCall(PetscObjectGetName((PetscObject)X, &vecname)); 118fb5bd1c2SPatrick Sanan 119fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 1209566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr, &fieldsnamed)); 121fb5bd1c2SPatrick Sanan if (fieldsnamed) { 122fb5bd1c2SPatrick Sanan for (f = 0; f < bs; f++) { 123fb5bd1c2SPatrick Sanan char buf[256]; 124fb5bd1c2SPatrick Sanan const char *fieldname; 1259566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(daCurr, f, &fieldname)); 126fb5bd1c2SPatrick Sanan if (!fieldname) { 12763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, f)); 128fb5bd1c2SPatrick Sanan fieldname = buf; 129fb5bd1c2SPatrick Sanan } 1300f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 1314a8c153fSJed Brown boffset += nnodes * sizeof(PetscScalar) + sizeof(PetscInt64); 132fb5bd1c2SPatrick Sanan } 133fb5bd1c2SPatrick Sanan } else { 1340f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%" PetscInt_FMT "\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, bs, boffset)); 1354a8c153fSJed Brown boffset += bs * nnodes * sizeof(PetscScalar) + sizeof(PetscInt64); 1364061b8bfSJed Brown } 137fb5bd1c2SPatrick Sanan } 1389566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </PointData>\n")); 1399566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </Piece>\n")); 1404061b8bfSJed Brown } 1419566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </StructuredGrid>\n")); 1429566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <AppendedData encoding=\"raw\">\n")); 1439566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "_")); 1444061b8bfSJed Brown 1454061b8bfSJed Brown /* Now write the arrays. */ 1464061b8bfSJed Brown tag = ((PetscObject)viewer)->tag; 1479566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(maxnnodes * PetscMax(3, maxbs), &array, maxnnodes * PetscMax(3, maxbs), &array2)); 1484061b8bfSJed Brown for (r = 0; r < size; r++) { 1494061b8bfSJed Brown MPI_Status status; 1504061b8bfSJed Brown PetscInt xs = -1, xm = -1, ys = -1, ym = -1, zs = -1, zm = -1, nnodes = 0; 151dd400576SPatrick Sanan if (rank == 0) { 1524061b8bfSJed Brown xs = grloc[r][0]; 1534061b8bfSJed Brown xm = grloc[r][1]; 1544061b8bfSJed Brown ys = grloc[r][2]; 1554061b8bfSJed Brown ym = grloc[r][3]; 1564061b8bfSJed Brown zs = grloc[r][4]; 1574061b8bfSJed Brown zm = grloc[r][5]; 1584061b8bfSJed Brown nnodes = xm * ym * zm; 1594061b8bfSJed Brown } else if (r == rank) { 1604061b8bfSJed Brown nnodes = info.xm * info.ym * info.zm; 1614061b8bfSJed Brown } 1624061b8bfSJed Brown 1632eaa9ef4SLisandro Dalcin /* Write the coordinates */ 1644061b8bfSJed Brown if (Coords) { 1654061b8bfSJed Brown const PetscScalar *coords; 1669566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Coords, &coords)); 167dd400576SPatrick Sanan if (rank == 0) { 1684061b8bfSJed Brown if (r) { 1696622f924SJed Brown PetscMPIInt nn; 1709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array, nnodes * cdim, MPIU_SCALAR, r, tag, comm, &status)); 1719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn)); 17208401ef6SPierre Jolivet PetscCheck(nn == nnodes * cdim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch"); 1731baa6e33SBarry Smith } else PetscCall(PetscArraycpy(array, coords, nnodes * cdim)); 1744061b8bfSJed Brown /* Transpose coordinates to VTK (C-style) ordering */ 1754061b8bfSJed Brown for (k = 0; k < zm; k++) { 1764061b8bfSJed Brown for (j = 0; j < ym; j++) { 1774061b8bfSJed Brown for (i = 0; i < xm; i++) { 1784061b8bfSJed Brown PetscInt Iloc = i + xm * (j + ym * k); 1792eaa9ef4SLisandro Dalcin array2[Iloc * 3 + 0] = array[Iloc * cdim + 0]; 1802eaa9ef4SLisandro Dalcin array2[Iloc * 3 + 1] = cdim > 1 ? array[Iloc * cdim + 1] : 0.0; 1812eaa9ef4SLisandro Dalcin array2[Iloc * 3 + 2] = cdim > 2 ? array[Iloc * cdim + 2] : 0.0; 1824061b8bfSJed Brown } 1834061b8bfSJed Brown } 1844061b8bfSJed Brown } 1854061b8bfSJed Brown } else if (r == rank) { 1869566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void *)coords, nnodes * cdim, MPIU_SCALAR, 0, tag, comm)); 1874061b8bfSJed Brown } 1889566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Coords, &coords)); 1894061b8bfSJed Brown } else { /* Fabricate some coordinates using grid index */ 1904061b8bfSJed Brown for (k = 0; k < zm; k++) { 1914061b8bfSJed Brown for (j = 0; j < ym; j++) { 1924061b8bfSJed Brown for (i = 0; i < xm; i++) { 1934061b8bfSJed Brown PetscInt Iloc = i + xm * (j + ym * k); 1944061b8bfSJed Brown array2[Iloc * 3 + 0] = xs + i; 1954061b8bfSJed Brown array2[Iloc * 3 + 1] = ys + j; 1964061b8bfSJed Brown array2[Iloc * 3 + 2] = zs + k; 1974061b8bfSJed Brown } 1984061b8bfSJed Brown } 1994061b8bfSJed Brown } 2004061b8bfSJed Brown } 2019566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, array2, nnodes * 3, MPIU_SCALAR)); 2024061b8bfSJed Brown 2034061b8bfSJed Brown /* Write each of the objects queued up for this file */ 2044061b8bfSJed Brown for (link = vtk->link; link; link = link->next) { 2054061b8bfSJed Brown Vec X = (Vec)link->vec; 2064061b8bfSJed Brown const PetscScalar *x; 207fb5bd1c2SPatrick Sanan PetscInt bs, f; 208ea2d7708SPatrick Sanan DM daCurr; 209fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 2109566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &daCurr)); 2119566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL)); 2129566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 213dd400576SPatrick Sanan if (rank == 0) { 2144061b8bfSJed Brown if (r) { 2156622f924SJed Brown PetscMPIInt nn; 2169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array, nnodes * bs, MPIU_SCALAR, r, tag, comm, &status)); 2179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn)); 21863a3b9bcSJacob Faibussowitsch PetscCheck(nn == nnodes * bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch receiving from rank %" PetscInt_FMT, r); 2191baa6e33SBarry Smith } else PetscCall(PetscArraycpy(array, x, nnodes * bs)); 220fb5bd1c2SPatrick Sanan 221fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 2229566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr, &fieldsnamed)); 223fb5bd1c2SPatrick Sanan if (fieldsnamed) { 224fb5bd1c2SPatrick Sanan for (f = 0; f < bs; f++) { 225fb5bd1c2SPatrick Sanan /* Extract and transpose the f'th field */ 226fb5bd1c2SPatrick Sanan for (k = 0; k < zm; k++) { 227fb5bd1c2SPatrick Sanan for (j = 0; j < ym; j++) { 228fb5bd1c2SPatrick Sanan for (i = 0; i < xm; i++) { 229fb5bd1c2SPatrick Sanan PetscInt Iloc = i + xm * (j + ym * k); 230fb5bd1c2SPatrick Sanan array2[Iloc] = array[Iloc * bs + f]; 231fb5bd1c2SPatrick Sanan } 232fb5bd1c2SPatrick Sanan } 233fb5bd1c2SPatrick Sanan } 2349566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, array2, nnodes, MPIU_SCALAR)); 235fb5bd1c2SPatrick Sanan } 2361baa6e33SBarry Smith } else PetscCall(PetscViewerVTKFWrite(viewer, fp, array, bs * nnodes, MPIU_SCALAR)); 2371baa6e33SBarry Smith } else if (r == rank) PetscCallMPI(MPI_Send((void *)x, nnodes * bs, MPIU_SCALAR, 0, tag, comm)); 2389566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 2394061b8bfSJed Brown } 2404061b8bfSJed Brown } 2419566063dSJacob Faibussowitsch PetscCall(PetscFree2(array, array2)); 2429566063dSJacob Faibussowitsch PetscCall(PetscFree(grloc)); 2434061b8bfSJed Brown 2449566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "\n </AppendedData>\n")); 2459566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "</VTKFile>\n")); 2469566063dSJacob Faibussowitsch PetscCall(PetscFClose(comm, fp)); 2473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2484061b8bfSJed Brown } 2494061b8bfSJed Brown 250d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAVTKWriteAll_VTR(DM da, PetscViewer viewer) 251d71ae5a4SJacob Faibussowitsch { 25230815ce0SLisandro Dalcin const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian"; 253a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_REAL_SINGLE) 254a13bc4e3SShao-Ching Huang const char precision[] = "Float32"; 255a13bc4e3SShao-Ching Huang #elif defined(PETSC_USE_REAL_DOUBLE) 256a13bc4e3SShao-Ching Huang const char precision[] = "Float64"; 257a13bc4e3SShao-Ching Huang #else 258a13bc4e3SShao-Ching Huang const char precision[] = "UnknownPrecision"; 259a13bc4e3SShao-Ching Huang #endif 260a13bc4e3SShao-Ching Huang MPI_Comm comm; 261a13bc4e3SShao-Ching Huang PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data; 262a13bc4e3SShao-Ching Huang PetscViewerVTKObjectLink link; 263a13bc4e3SShao-Ching Huang FILE *fp; 264a13bc4e3SShao-Ching Huang PetscMPIInt rank, size, tag; 265a13bc4e3SShao-Ching Huang DMDALocalInfo info; 2660f2609c8SJed Brown PetscInt dim, mx, my, mz, maxnnodes, maxbs, i, j, k, r; 2670f2609c8SJed Brown PetscInt64 boffset; 268a13bc4e3SShao-Ching Huang PetscInt rloc[6], (*grloc)[6] = NULL; 269a13bc4e3SShao-Ching Huang PetscScalar *array, *array2; 270a13bc4e3SShao-Ching Huang 271a13bc4e3SShao-Ching Huang PetscFunctionBegin; 2729566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 2731dca8a05SBarry Smith PetscCheck(!PetscDefined(USE_COMPLEX), PETSC_COMM_SELF, PETSC_ERR_SUP, "Complex values not supported"); 2749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2769566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, &mx, &my, &mz, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2779566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 2789566063dSJacob Faibussowitsch PetscCall(PetscFOpen(comm, vtk->filename, "wb", &fp)); 2799566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "<?xml version=\"1.0\"?>\n")); 2804a8c153fSJed Brown PetscCall(PetscFPrintf(comm, fp, "<VTKFile type=\"RectilinearGrid\" version=\"0.1\" byte_order=\"%s\" header_type=\"UInt64\">\n", byte_order)); 28163a3b9bcSJacob 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)); 282a13bc4e3SShao-Ching Huang 2839566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc1(size * 6, &grloc)); 284a13bc4e3SShao-Ching Huang rloc[0] = info.xs; 285a13bc4e3SShao-Ching Huang rloc[1] = info.xm; 286a13bc4e3SShao-Ching Huang rloc[2] = info.ys; 287a13bc4e3SShao-Ching Huang rloc[3] = info.ym; 288a13bc4e3SShao-Ching Huang rloc[4] = info.zs; 289a13bc4e3SShao-Ching Huang rloc[5] = info.zm; 290810441c8SPierre Jolivet PetscCallMPI(MPI_Gather(rloc, 6, MPIU_INT, rank == 0 ? grloc[0] : NULL, 6, MPIU_INT, 0, comm)); 291a13bc4e3SShao-Ching Huang 292a13bc4e3SShao-Ching Huang /* Write XML header */ 293a13bc4e3SShao-Ching Huang maxnnodes = 0; /* Used for the temporary array size on rank 0 */ 294ea2d7708SPatrick Sanan maxbs = 0; /* Used for the temporary array size on rank 0 */ 295a13bc4e3SShao-Ching Huang boffset = 0; /* Offset into binary file */ 296a13bc4e3SShao-Ching Huang for (r = 0; r < size; r++) { 297a13bc4e3SShao-Ching Huang PetscInt xs = -1, xm = -1, ys = -1, ym = -1, zs = -1, zm = -1, nnodes = 0; 298dd400576SPatrick Sanan if (rank == 0) { 299a13bc4e3SShao-Ching Huang xs = grloc[r][0]; 300a13bc4e3SShao-Ching Huang xm = grloc[r][1]; 301a13bc4e3SShao-Ching Huang ys = grloc[r][2]; 302a13bc4e3SShao-Ching Huang ym = grloc[r][3]; 303a13bc4e3SShao-Ching Huang zs = grloc[r][4]; 304a13bc4e3SShao-Ching Huang zm = grloc[r][5]; 305a13bc4e3SShao-Ching Huang nnodes = xm * ym * zm; 306a13bc4e3SShao-Ching Huang } 307a13bc4e3SShao-Ching Huang maxnnodes = PetscMax(maxnnodes, nnodes); 30863a3b9bcSJacob 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)); 3099566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <Coordinates>\n")); 3100f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"Xcoord\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, boffset)); 3114a8c153fSJed Brown boffset += xm * sizeof(PetscScalar) + sizeof(PetscInt64); 3120f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"Ycoord\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, boffset)); 3134a8c153fSJed Brown boffset += ym * sizeof(PetscScalar) + sizeof(PetscInt64); 3140f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"Zcoord\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, boffset)); 3154a8c153fSJed Brown boffset += zm * sizeof(PetscScalar) + sizeof(PetscInt64); 3169566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </Coordinates>\n")); 3179566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <PointData Scalars=\"ScalarPointData\">\n")); 318a13bc4e3SShao-Ching Huang for (link = vtk->link; link; link = link->next) { 319a13bc4e3SShao-Ching Huang Vec X = (Vec)link->vec; 320fb5bd1c2SPatrick Sanan PetscInt bs, f; 321ea2d7708SPatrick Sanan DM daCurr; 322fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 323fb5bd1c2SPatrick Sanan const char *vecname = "Unnamed"; 324ea2d7708SPatrick Sanan 3259566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &daCurr)); 3269566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL)); 327ea2d7708SPatrick Sanan maxbs = PetscMax(maxbs, bs); 32848a46eb9SPierre Jolivet if (((PetscObject)X)->name || link != vtk->link) PetscCall(PetscObjectGetName((PetscObject)X, &vecname)); 329ea2d7708SPatrick Sanan 330fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 3319566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr, &fieldsnamed)); 332fb5bd1c2SPatrick Sanan if (fieldsnamed) { 333fb5bd1c2SPatrick Sanan for (f = 0; f < bs; f++) { 334fb5bd1c2SPatrick Sanan char buf[256]; 335fb5bd1c2SPatrick Sanan const char *fieldname; 3369566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(daCurr, f, &fieldname)); 337fb5bd1c2SPatrick Sanan if (!fieldname) { 33863a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, f)); 339fb5bd1c2SPatrick Sanan fieldname = buf; 340fb5bd1c2SPatrick Sanan } 3410f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset)); 3424a8c153fSJed Brown boffset += nnodes * sizeof(PetscScalar) + sizeof(PetscInt64); 343fb5bd1c2SPatrick Sanan } 344fb5bd1c2SPatrick Sanan } else { 3450f2609c8SJed Brown PetscCall(PetscFPrintf(comm, fp, " <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%" PetscInt_FMT "\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, bs, boffset)); 3464a8c153fSJed Brown boffset += bs * nnodes * sizeof(PetscScalar) + sizeof(PetscInt64); 347a13bc4e3SShao-Ching Huang } 348fb5bd1c2SPatrick Sanan } 3499566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </PointData>\n")); 3509566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </Piece>\n")); 351a13bc4e3SShao-Ching Huang } 3529566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " </RectilinearGrid>\n")); 3539566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, " <AppendedData encoding=\"raw\">\n")); 3549566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "_")); 355a13bc4e3SShao-Ching Huang 356a13bc4e3SShao-Ching Huang /* Now write the arrays. */ 357a13bc4e3SShao-Ching Huang tag = ((PetscObject)viewer)->tag; 3589566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(PetscMax(maxnnodes * maxbs, info.xm + info.ym + info.zm), &array, PetscMax(maxnnodes * maxbs, info.xm + info.ym + info.zm), &array2)); 359a13bc4e3SShao-Ching Huang for (r = 0; r < size; r++) { 360a13bc4e3SShao-Ching Huang MPI_Status status; 361a13bc4e3SShao-Ching Huang PetscInt xs = -1, xm = -1, ys = -1, ym = -1, zs = -1, zm = -1, nnodes = 0; 362dd400576SPatrick Sanan if (rank == 0) { 363a13bc4e3SShao-Ching Huang xs = grloc[r][0]; 364a13bc4e3SShao-Ching Huang xm = grloc[r][1]; 365a13bc4e3SShao-Ching Huang ys = grloc[r][2]; 366a13bc4e3SShao-Ching Huang ym = grloc[r][3]; 367a13bc4e3SShao-Ching Huang zs = grloc[r][4]; 368a13bc4e3SShao-Ching Huang zm = grloc[r][5]; 369a13bc4e3SShao-Ching Huang nnodes = xm * ym * zm; 370a13bc4e3SShao-Ching Huang } else if (r == rank) { 371a13bc4e3SShao-Ching Huang nnodes = info.xm * info.ym * info.zm; 372a13bc4e3SShao-Ching Huang } 373a13bc4e3SShao-Ching Huang { /* Write the coordinates */ 374a13bc4e3SShao-Ching Huang Vec Coords; 375ea2d7708SPatrick Sanan 3769566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &Coords)); 377a13bc4e3SShao-Ching Huang if (Coords) { 378a13bc4e3SShao-Ching Huang const PetscScalar *coords; 3799566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Coords, &coords)); 380dd400576SPatrick Sanan if (rank == 0) { 381a13bc4e3SShao-Ching Huang if (r) { 382a13bc4e3SShao-Ching Huang PetscMPIInt nn; 3839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array, xm + ym + zm, MPIU_SCALAR, r, tag, comm, &status)); 3849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn)); 38508401ef6SPierre Jolivet PetscCheck(nn == xm + ym + zm, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch"); 386a13bc4e3SShao-Ching Huang } else { 387a13bc4e3SShao-Ching Huang /* x coordinates */ 388a13bc4e3SShao-Ching Huang for (j = 0, k = 0, i = 0; i < xm; i++) { 389a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 390a13bc4e3SShao-Ching Huang array[i] = coords[Iloc * dim + 0]; 391a13bc4e3SShao-Ching Huang } 392a13bc4e3SShao-Ching Huang /* y coordinates */ 393a13bc4e3SShao-Ching Huang for (i = 0, k = 0, j = 0; j < ym; j++) { 394a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 395a13bc4e3SShao-Ching Huang array[j + xm] = dim > 1 ? coords[Iloc * dim + 1] : 0; 396a13bc4e3SShao-Ching Huang } 397a13bc4e3SShao-Ching Huang /* z coordinates */ 398a13bc4e3SShao-Ching Huang for (i = 0, j = 0, k = 0; k < zm; k++) { 399a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 400a13bc4e3SShao-Ching Huang array[k + xm + ym] = dim > 2 ? coords[Iloc * dim + 2] : 0; 401a13bc4e3SShao-Ching Huang } 402a13bc4e3SShao-Ching Huang } 403a13bc4e3SShao-Ching Huang } else if (r == rank) { 404a13bc4e3SShao-Ching Huang xm = info.xm; 405a13bc4e3SShao-Ching Huang ym = info.ym; 406a13bc4e3SShao-Ching Huang zm = info.zm; 407a13bc4e3SShao-Ching Huang for (j = 0, k = 0, i = 0; i < xm; i++) { 408a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 409a13bc4e3SShao-Ching Huang array2[i] = coords[Iloc * dim + 0]; 410a13bc4e3SShao-Ching Huang } 411a13bc4e3SShao-Ching Huang for (i = 0, k = 0, j = 0; j < ym; j++) { 412a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 413a13bc4e3SShao-Ching Huang array2[j + xm] = dim > 1 ? coords[Iloc * dim + 1] : 0; 414a13bc4e3SShao-Ching Huang } 415a13bc4e3SShao-Ching Huang for (i = 0, j = 0, k = 0; k < zm; k++) { 416a13bc4e3SShao-Ching Huang PetscInt Iloc = i + xm * (j + ym * k); 417a13bc4e3SShao-Ching Huang array2[k + xm + ym] = dim > 2 ? coords[Iloc * dim + 2] : 0; 418a13bc4e3SShao-Ching Huang } 4199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void *)array2, xm + ym + zm, MPIU_SCALAR, 0, tag, comm)); 420a13bc4e3SShao-Ching Huang } 4219566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Coords, &coords)); 422a13bc4e3SShao-Ching Huang } else { /* Fabricate some coordinates using grid index */ 423ad540459SPierre Jolivet for (i = 0; i < xm; i++) array[i] = xs + i; 424ad540459SPierre Jolivet for (j = 0; j < ym; j++) array[j + xm] = ys + j; 425ad540459SPierre Jolivet for (k = 0; k < zm; k++) array[k + xm + ym] = zs + k; 426a13bc4e3SShao-Ching Huang } 427dd400576SPatrick Sanan if (rank == 0) { 428f4f49eeaSPierre Jolivet PetscCall(PetscViewerVTKFWrite(viewer, fp, &array[0], xm, MPIU_SCALAR)); 429f4f49eeaSPierre Jolivet PetscCall(PetscViewerVTKFWrite(viewer, fp, &array[xm], ym, MPIU_SCALAR)); 430f4f49eeaSPierre Jolivet PetscCall(PetscViewerVTKFWrite(viewer, fp, &array[xm + ym], zm, MPIU_SCALAR)); 431a13bc4e3SShao-Ching Huang } 432a13bc4e3SShao-Ching Huang } 433a13bc4e3SShao-Ching Huang 434a13bc4e3SShao-Ching Huang /* Write each of the objects queued up for this file */ 435a13bc4e3SShao-Ching Huang for (link = vtk->link; link; link = link->next) { 436a13bc4e3SShao-Ching Huang Vec X = (Vec)link->vec; 437a13bc4e3SShao-Ching Huang const PetscScalar *x; 438fb5bd1c2SPatrick Sanan PetscInt bs, f; 439ea2d7708SPatrick Sanan DM daCurr; 440fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 4419566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &daCurr)); 4429566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL)); 443a13bc4e3SShao-Ching Huang 4449566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 445dd400576SPatrick Sanan if (rank == 0) { 446a13bc4e3SShao-Ching Huang if (r) { 447a13bc4e3SShao-Ching Huang PetscMPIInt nn; 4489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array, nnodes * bs, MPIU_SCALAR, r, tag, comm, &status)); 4499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn)); 45063a3b9bcSJacob Faibussowitsch PetscCheck(nn == nnodes * bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch receiving from rank %" PetscInt_FMT, r); 4511baa6e33SBarry Smith } else PetscCall(PetscArraycpy(array, x, nnodes * bs)); 452fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 4539566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr, &fieldsnamed)); 454fb5bd1c2SPatrick Sanan if (fieldsnamed) { 455fb5bd1c2SPatrick Sanan for (f = 0; f < bs; f++) { 456fb5bd1c2SPatrick Sanan /* Extract and transpose the f'th field */ 457fb5bd1c2SPatrick Sanan for (k = 0; k < zm; k++) { 458fb5bd1c2SPatrick Sanan for (j = 0; j < ym; j++) { 459fb5bd1c2SPatrick Sanan for (i = 0; i < xm; i++) { 460fb5bd1c2SPatrick Sanan PetscInt Iloc = i + xm * (j + ym * k); 461fb5bd1c2SPatrick Sanan array2[Iloc] = array[Iloc * bs + f]; 462fb5bd1c2SPatrick Sanan } 463fb5bd1c2SPatrick Sanan } 464fb5bd1c2SPatrick Sanan } 4659566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, array2, nnodes, MPIU_SCALAR)); 466fb5bd1c2SPatrick Sanan } 467fb5bd1c2SPatrick Sanan } 4689566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer, fp, array, nnodes * bs, MPIU_SCALAR)); 469ea2d7708SPatrick Sanan 470a13bc4e3SShao-Ching Huang } else if (r == rank) { 4719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void *)x, nnodes * bs, MPIU_SCALAR, 0, tag, comm)); 472a13bc4e3SShao-Ching Huang } 4739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 474a13bc4e3SShao-Ching Huang } 475a13bc4e3SShao-Ching Huang } 4769566063dSJacob Faibussowitsch PetscCall(PetscFree2(array, array2)); 4779566063dSJacob Faibussowitsch PetscCall(PetscFree(grloc)); 478a13bc4e3SShao-Ching Huang 4799566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "\n </AppendedData>\n")); 4809566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, fp, "</VTKFile>\n")); 4819566063dSJacob Faibussowitsch PetscCall(PetscFClose(comm, fp)); 4823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 483a13bc4e3SShao-Ching Huang } 484a13bc4e3SShao-Ching Huang 485*ffeef943SBarry Smith /*@ 4864061b8bfSJed Brown DMDAVTKWriteAll - Write a file containing all the fields that have been provided to the viewer 4874061b8bfSJed Brown 4884061b8bfSJed Brown Collective 4894061b8bfSJed Brown 4904165533cSJose E. Roman Input Parameters: 491dce8aebaSBarry Smith + odm - `DMDA` specifying the grid layout, passed as a `PetscObject` 492dce8aebaSBarry Smith - viewer - viewer of type `PETSCVIEWERVTK` 4934061b8bfSJed Brown 4944061b8bfSJed Brown Level: developer 4954061b8bfSJed Brown 496fb5bd1c2SPatrick Sanan Notes: 49712b4a537SBarry Smith This function is a callback used by the `PETSCVIEWERVTK` viewer to actually write the file. 4984061b8bfSJed 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. 4994061b8bfSJed Brown Instead, metadata for the entire file needs to be available up-front before you can start writing the file. 5004061b8bfSJed Brown 50112b4a537SBarry Smith If any fields have been named (see e.g. `DMDASetFieldName()`), then individual scalar 502fb5bd1c2SPatrick Sanan fields are written. Otherwise, a single multi-dof (vector) field is written. 503fb5bd1c2SPatrick Sanan 50412b4a537SBarry Smith .seealso: [](sec_struct), `DMDA`, `DM`, `PETSCVIEWERVTK`, `DMDASetFieldName()` 5054061b8bfSJed Brown @*/ 506d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVTKWriteAll(PetscObject odm, PetscViewer viewer) 507d71ae5a4SJacob Faibussowitsch { 5084061b8bfSJed Brown DM dm = (DM)odm; 5094061b8bfSJed Brown PetscBool isvtk; 5104061b8bfSJed Brown 5114061b8bfSJed Brown PetscFunctionBegin; 5124061b8bfSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5134061b8bfSJed Brown PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 5149566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk)); 5157a8be351SBarry Smith PetscCheck(isvtk, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name); 5164061b8bfSJed Brown switch (viewer->format) { 517d71ae5a4SJacob Faibussowitsch case PETSC_VIEWER_VTK_VTS: 518d71ae5a4SJacob Faibussowitsch PetscCall(DMDAVTKWriteAll_VTS(dm, viewer)); 519d71ae5a4SJacob Faibussowitsch break; 520d71ae5a4SJacob Faibussowitsch case PETSC_VIEWER_VTK_VTR: 521d71ae5a4SJacob Faibussowitsch PetscCall(DMDAVTKWriteAll_VTR(dm, viewer)); 522d71ae5a4SJacob Faibussowitsch break; 523d71ae5a4SJacob Faibussowitsch default: 524d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]); 5254061b8bfSJed Brown } 5263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5274061b8bfSJed Brown } 528