147c6ae99SBarry Smith /* 2aa219208SBarry Smith Plots vectors obtained with DMDACreate2d() 347c6ae99SBarry Smith */ 447c6ae99SBarry Smith 5af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 68135c375SStefano Zampini #include <petsc/private/glvisvecimpl.h> 72d95be3dSVaclav Hapla #include <petsc/private/viewerhdf5impl.h> 89804daf3SBarry Smith #include <petscdraw.h> 947c6ae99SBarry Smith 1047c6ae99SBarry Smith /* 1147c6ae99SBarry Smith The data that is passed into the graphics callback 1247c6ae99SBarry Smith */ 1347c6ae99SBarry Smith typedef struct { 14e5ab1681SLisandro Dalcin PetscMPIInt rank; 15e5ab1681SLisandro Dalcin PetscInt m, n, dof, k; 16e5ab1681SLisandro Dalcin PetscReal xmin, xmax, ymin, ymax, min, max; 175edff71fSBarry Smith const PetscScalar *xy, *v; 187fe7c8feSLisandro Dalcin PetscBool showaxis, showgrid; 19109c9344SBarry Smith const char *name0, *name1; 2047c6ae99SBarry Smith } ZoomCtx; 2147c6ae99SBarry Smith 2247c6ae99SBarry Smith /* 2347c6ae99SBarry Smith This does the drawing for one particular field 2447c6ae99SBarry Smith in one particular set of coordinates. It is a callback 2547c6ae99SBarry Smith called from PetscDrawZoom() 2647c6ae99SBarry Smith */ 2766976f2fSJacob Faibussowitsch static PetscErrorCode VecView_MPI_Draw_DA2d_Zoom(PetscDraw draw, void *ctx) 28d71ae5a4SJacob Faibussowitsch { 2947c6ae99SBarry Smith ZoomCtx *zctx = (ZoomCtx *)ctx; 30e5ab1681SLisandro Dalcin PetscInt m, n, i, j, k, dof, id, c1, c2, c3, c4; 31e5ab1681SLisandro Dalcin PetscReal min, max, x1, x2, x3, x4, y_1, y2, y3, y4; 32e5ab1681SLisandro Dalcin const PetscScalar *xy, *v; 3347c6ae99SBarry Smith 3447c6ae99SBarry Smith PetscFunctionBegin; 3547c6ae99SBarry Smith m = zctx->m; 3647c6ae99SBarry Smith n = zctx->n; 37e5ab1681SLisandro Dalcin dof = zctx->dof; 3847c6ae99SBarry Smith k = zctx->k; 3947c6ae99SBarry Smith xy = zctx->xy; 40e5ab1681SLisandro Dalcin v = zctx->v; 4147c6ae99SBarry Smith min = zctx->min; 42f3f0eb19SBarry Smith max = zctx->max; 4347c6ae99SBarry Smith 4447c6ae99SBarry Smith /* PetscDraw the contour plot patch */ 45d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 4647c6ae99SBarry Smith for (j = 0; j < n - 1; j++) { 4747c6ae99SBarry Smith for (i = 0; i < m - 1; i++) { 480e4fe250SBarry Smith id = i + j * m; 490e4fe250SBarry Smith x1 = PetscRealPart(xy[2 * id]); 500e4fe250SBarry Smith y_1 = PetscRealPart(xy[2 * id + 1]); 51e5ab1681SLisandro Dalcin c1 = PetscDrawRealToColor(PetscRealPart(v[k + dof * id]), min, max); 520e4fe250SBarry Smith 530e4fe250SBarry Smith id = i + j * m + 1; 540e4fe250SBarry Smith x2 = PetscRealPart(xy[2 * id]); 55a39a4c7dSLisandro Dalcin y2 = PetscRealPart(xy[2 * id + 1]); 56e5ab1681SLisandro Dalcin c2 = PetscDrawRealToColor(PetscRealPart(v[k + dof * id]), min, max); 570e4fe250SBarry Smith 580e4fe250SBarry Smith id = i + j * m + 1 + m; 59a39a4c7dSLisandro Dalcin x3 = PetscRealPart(xy[2 * id]); 600e4fe250SBarry Smith y3 = PetscRealPart(xy[2 * id + 1]); 61e5ab1681SLisandro Dalcin c3 = PetscDrawRealToColor(PetscRealPart(v[k + dof * id]), min, max); 620e4fe250SBarry Smith 630e4fe250SBarry Smith id = i + j * m + m; 64a39a4c7dSLisandro Dalcin x4 = PetscRealPart(xy[2 * id]); 65a39a4c7dSLisandro Dalcin y4 = PetscRealPart(xy[2 * id + 1]); 66e5ab1681SLisandro Dalcin c4 = PetscDrawRealToColor(PetscRealPart(v[k + dof * id]), min, max); 67f3f0eb19SBarry Smith 689566063dSJacob Faibussowitsch PetscCall(PetscDrawTriangle(draw, x1, y_1, x2, y2, x3, y3, c1, c2, c3)); 699566063dSJacob Faibussowitsch PetscCall(PetscDrawTriangle(draw, x1, y_1, x3, y3, x4, y4, c1, c3, c4)); 7047c6ae99SBarry Smith if (zctx->showgrid) { 719566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, x1, y_1, x2, y2, PETSC_DRAW_BLACK)); 729566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, x2, y2, x3, y3, PETSC_DRAW_BLACK)); 739566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, x3, y3, x4, y4, PETSC_DRAW_BLACK)); 749566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, x4, y4, x1, y_1, PETSC_DRAW_BLACK)); 7547c6ae99SBarry Smith } 7647c6ae99SBarry Smith } 7747c6ae99SBarry Smith } 787fe7c8feSLisandro Dalcin if (zctx->showaxis && !zctx->rank) { 79e5ab1681SLisandro Dalcin if (zctx->name0 || zctx->name1) { 80109c9344SBarry Smith PetscReal xl, yl, xr, yr, x, y; 819566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 82e5ab1681SLisandro Dalcin x = xl + .30 * (xr - xl); 83109c9344SBarry Smith xl = xl + .01 * (xr - xl); 84e5ab1681SLisandro Dalcin y = yr - .30 * (yr - yl); 85109c9344SBarry Smith yl = yl + .01 * (yr - yl); 869566063dSJacob Faibussowitsch if (zctx->name0) PetscCall(PetscDrawString(draw, x, yl, PETSC_DRAW_BLACK, zctx->name0)); 879566063dSJacob Faibussowitsch if (zctx->name1) PetscCall(PetscDrawStringVertical(draw, xl, y, PETSC_DRAW_BLACK, zctx->name1)); 88109c9344SBarry Smith } 890e4fe250SBarry Smith /* 900e4fe250SBarry Smith Ideally we would use the PetscDrawAxis object to manage displaying the coordinate limits 910e4fe250SBarry Smith but that may require some refactoring. 920e4fe250SBarry Smith */ 93e5ab1681SLisandro Dalcin { 947fe7c8feSLisandro Dalcin double xmin = (double)zctx->xmin, ymin = (double)zctx->ymin; 957fe7c8feSLisandro Dalcin double xmax = (double)zctx->xmax, ymax = (double)zctx->ymax; 969371c9d4SSatish Balay char value[16]; 979371c9d4SSatish Balay size_t len; 989371c9d4SSatish Balay PetscReal w; 999566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(value, 16, "%0.2e", xmin)); 1009566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw, xmin, ymin - .05 * (ymax - ymin), PETSC_DRAW_BLACK, value)); 1019566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(value, 16, "%0.2e", xmax)); 1029566063dSJacob Faibussowitsch PetscCall(PetscStrlen(value, &len)); 1039566063dSJacob Faibussowitsch PetscCall(PetscDrawStringGetSize(draw, &w, NULL)); 1049566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw, xmax - len * w, ymin - .05 * (ymax - ymin), PETSC_DRAW_BLACK, value)); 1059566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(value, 16, "%0.2e", ymin)); 1069566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw, xmin - .05 * (xmax - xmin), ymin, PETSC_DRAW_BLACK, value)); 1079566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(value, 16, "%0.2e", ymax)); 1089566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw, xmin - .05 * (xmax - xmin), ymax, PETSC_DRAW_BLACK, value)); 109e5ab1681SLisandro Dalcin } 110e5ab1681SLisandro Dalcin } 111d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11347c6ae99SBarry Smith } 11447c6ae99SBarry Smith 11566976f2fSJacob Faibussowitsch static PetscErrorCode VecView_MPI_Draw_DA2d(Vec xin, PetscViewer viewer) 116d71ae5a4SJacob Faibussowitsch { 1179a42bb27SBarry Smith DM da, dac, dag; 118a74ba6f7SBarry Smith PetscInt N, s, M, w, ncoors = 4; 11947c6ae99SBarry Smith const PetscInt *lx, *ly; 120e5ab1681SLisandro Dalcin PetscReal coors[4]; 12147c6ae99SBarry Smith PetscDraw draw, popup; 12247c6ae99SBarry Smith PetscBool isnull, useports = PETSC_FALSE; 12347c6ae99SBarry Smith MPI_Comm comm; 12447c6ae99SBarry Smith Vec xlocal, xcoor, xcoorl; 125bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by; 126aa219208SBarry Smith DMDAStencilType st; 12747c6ae99SBarry Smith ZoomCtx zctx; 1280298fd71SBarry Smith PetscDrawViewPorts *ports = NULL; 12947c6ae99SBarry Smith PetscViewerFormat format; 13020d0051dSBarry Smith PetscInt *displayfields; 13167dd0837SBarry Smith PetscInt ndisplayfields, i, nbounds; 13267dd0837SBarry Smith const PetscReal *bounds; 13347c6ae99SBarry Smith 13447c6ae99SBarry Smith PetscFunctionBegin; 13547c6ae99SBarry Smith zctx.showgrid = PETSC_FALSE; 1367fe7c8feSLisandro Dalcin zctx.showaxis = PETSC_TRUE; 1378865f1eaSKarl Rupp 1389566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 1399566063dSJacob Faibussowitsch PetscCall(PetscDrawIsNull(draw, &isnull)); 1403ba16761SJacob Faibussowitsch if (isnull) PetscFunctionReturn(PETSC_SUCCESS); 1415b399a63SLisandro Dalcin 1429566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetBounds(viewer, &nbounds, &bounds)); 14347c6ae99SBarry Smith 1449566063dSJacob Faibussowitsch PetscCall(VecGetDM(xin, &da)); 1457a8be351SBarry Smith PetscCheck(da, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA"); 14647c6ae99SBarry Smith 1479566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)xin, &comm)); 1489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &zctx.rank)); 14947c6ae99SBarry Smith 1509566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, &M, &N, NULL, &zctx.m, &zctx.n, NULL, &w, &s, &bx, &by, NULL, &st)); 1519566063dSJacob Faibussowitsch PetscCall(DMDAGetOwnershipRanges(da, &lx, &ly, NULL)); 15247c6ae99SBarry Smith 15347c6ae99SBarry Smith /* 15447c6ae99SBarry Smith Obtain a sequential vector that is going to contain the local values plus ONE layer of 155aa219208SBarry Smith ghosted values to draw the graphics from. We also need its corresponding DMDA (dac) that will 156d5b43468SJose E. Roman update the local values plus ONE layer of ghost values. 15747c6ae99SBarry Smith */ 1589566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)da, "GraphicsGhosted", (PetscObject *)&xlocal)); 15947c6ae99SBarry Smith if (!xlocal) { 160bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_NONE || by != DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) { 16147c6ae99SBarry Smith /* 16247c6ae99SBarry Smith if original da is not of stencil width one, or periodic or not a box stencil then 163aa219208SBarry Smith create a special DMDA to handle one level of ghost points for graphics 16447c6ae99SBarry Smith */ 1659566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, zctx.m, zctx.n, w, 1, lx, ly, &dac)); 1669566063dSJacob Faibussowitsch PetscCall(DMSetUp(dac)); 1676aad120cSJose E. Roman PetscCall(PetscInfo(da, "Creating auxiliary DMDA for managing graphics ghost points\n")); 16847c6ae99SBarry Smith } else { 16947c6ae99SBarry Smith /* otherwise we can use the da we already have */ 17047c6ae99SBarry Smith dac = da; 17147c6ae99SBarry Smith } 17247c6ae99SBarry Smith /* create local vector for holding ghosted values used in graphics */ 1739566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dac, &xlocal)); 17447c6ae99SBarry Smith if (dac != da) { 175*15229ffcSPierre Jolivet /* don't keep any public reference of this DMDA, it is only available through xlocal */ 1769566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dac)); 17747c6ae99SBarry Smith } else { 17847c6ae99SBarry Smith /* remove association between xlocal and da, because below we compose in the opposite 17947c6ae99SBarry Smith direction and if we left this connect we'd get a loop, so the objects could 18047c6ae99SBarry Smith never be destroyed */ 1819566063dSJacob Faibussowitsch PetscCall(PetscObjectRemoveReference((PetscObject)xlocal, "__PETSc_dm")); 18247c6ae99SBarry Smith } 1839566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)da, "GraphicsGhosted", (PetscObject)xlocal)); 1849566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)xlocal)); 18547c6ae99SBarry Smith } else { 186bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_NONE || by != DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) { 1879566063dSJacob Faibussowitsch PetscCall(VecGetDM(xlocal, &dac)); 188f7923d8aSBarry Smith } else { 189f7923d8aSBarry Smith dac = da; 19047c6ae99SBarry Smith } 19147c6ae99SBarry Smith } 19247c6ae99SBarry Smith 19347c6ae99SBarry Smith /* 19447c6ae99SBarry Smith Get local (ghosted) values of vector 19547c6ae99SBarry Smith */ 1969566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dac, xin, INSERT_VALUES, xlocal)); 1979566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dac, xin, INSERT_VALUES, xlocal)); 1989566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xlocal, &zctx.v)); 19947c6ae99SBarry Smith 200832b7cebSLisandro Dalcin /* 201832b7cebSLisandro Dalcin Get coordinates of nodes 202832b7cebSLisandro Dalcin */ 2039566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &xcoor)); 20447c6ae99SBarry Smith if (!xcoor) { 2059566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 0.0)); 2069566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &xcoor)); 20747c6ae99SBarry Smith } 20847c6ae99SBarry Smith 20947c6ae99SBarry Smith /* 21047c6ae99SBarry Smith Determine the min and max coordinates in plot 21147c6ae99SBarry Smith */ 2129566063dSJacob Faibussowitsch PetscCall(VecStrideMin(xcoor, 0, NULL, &zctx.xmin)); 2139566063dSJacob Faibussowitsch PetscCall(VecStrideMax(xcoor, 0, NULL, &zctx.xmax)); 2149566063dSJacob Faibussowitsch PetscCall(VecStrideMin(xcoor, 1, NULL, &zctx.ymin)); 2159566063dSJacob Faibussowitsch PetscCall(VecStrideMax(xcoor, 1, NULL, &zctx.ymax)); 2169566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_contour_axis", &zctx.showaxis, NULL)); 2177fe7c8feSLisandro Dalcin if (zctx.showaxis) { 2189371c9d4SSatish Balay coors[0] = zctx.xmin - .05 * (zctx.xmax - zctx.xmin); 2199371c9d4SSatish Balay coors[1] = zctx.ymin - .05 * (zctx.ymax - zctx.ymin); 2209371c9d4SSatish Balay coors[2] = zctx.xmax + .05 * (zctx.xmax - zctx.xmin); 2219371c9d4SSatish Balay coors[3] = zctx.ymax + .05 * (zctx.ymax - zctx.ymin); 2227fe7c8feSLisandro Dalcin } else { 2239371c9d4SSatish Balay coors[0] = zctx.xmin; 2249371c9d4SSatish Balay coors[1] = zctx.ymin; 2259371c9d4SSatish Balay coors[2] = zctx.xmax; 2269371c9d4SSatish Balay coors[3] = zctx.ymax; 2277fe7c8feSLisandro Dalcin } 2289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-draw_coordinates", coors, &ncoors, NULL)); 2299566063dSJacob Faibussowitsch PetscCall(PetscInfo(da, "Preparing DMDA 2d contour plot coordinates %g %g %g %g\n", (double)coors[0], (double)coors[1], (double)coors[2], (double)coors[3])); 230a74ba6f7SBarry Smith 23147c6ae99SBarry Smith /* 232832b7cebSLisandro Dalcin Get local ghosted version of coordinates 23347c6ae99SBarry Smith */ 2349566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)da, "GraphicsCoordinateGhosted", (PetscObject *)&xcoorl)); 23547c6ae99SBarry Smith if (!xcoorl) { 236aa219208SBarry Smith /* create DMDA to get local version of graphics */ 2379566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, zctx.m, zctx.n, 2, 1, lx, ly, &dag)); 2389566063dSJacob Faibussowitsch PetscCall(DMSetUp(dag)); 2396aad120cSJose E. Roman PetscCall(PetscInfo(dag, "Creating auxiliary DMDA for managing graphics coordinates ghost points\n")); 2409566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dag, &xcoorl)); 2419566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)da, "GraphicsCoordinateGhosted", (PetscObject)xcoorl)); 2429566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dag)); 2439566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)xcoorl)); 2441baa6e33SBarry Smith } else PetscCall(VecGetDM(xcoorl, &dag)); 2459566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dag, xcoor, INSERT_VALUES, xcoorl)); 2469566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dag, xcoor, INSERT_VALUES, xcoorl)); 2479566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xcoorl, &zctx.xy)); 2489566063dSJacob Faibussowitsch PetscCall(DMDAGetCoordinateName(da, 0, &zctx.name0)); 2499566063dSJacob Faibussowitsch PetscCall(DMDAGetCoordinateName(da, 1, &zctx.name1)); 25047c6ae99SBarry Smith 25147c6ae99SBarry Smith /* 25247c6ae99SBarry Smith Get information about size of area each processor must do graphics for 25347c6ae99SBarry Smith */ 2549566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &M, &N, NULL, NULL, NULL, NULL, &zctx.dof, NULL, &bx, &by, NULL, NULL)); 2559566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, NULL, NULL, NULL, &zctx.m, &zctx.n, NULL)); 2569566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_contour_grid", &zctx.showgrid, NULL)); 2574e6118eeSBarry Smith 2589566063dSJacob Faibussowitsch PetscCall(DMDASelectFields(da, &ndisplayfields, &displayfields)); 2599566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 2609566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_ports", &useports, NULL)); 261832b7cebSLisandro Dalcin if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE; 262832b7cebSLisandro Dalcin if (useports) { 2639566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 2649566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 2659566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 2669566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsCreate(draw, ndisplayfields, &ports)); 267e5ab1681SLisandro Dalcin } 26820d0051dSBarry Smith 269832b7cebSLisandro Dalcin /* 270832b7cebSLisandro Dalcin Loop over each field; drawing each in a different window 271832b7cebSLisandro Dalcin */ 27220d0051dSBarry Smith for (i = 0; i < ndisplayfields; i++) { 27320d0051dSBarry Smith zctx.k = displayfields[i]; 2745b399a63SLisandro Dalcin 275e5ab1681SLisandro Dalcin /* determine the min and max value in plot */ 2769566063dSJacob Faibussowitsch PetscCall(VecStrideMin(xin, zctx.k, NULL, &zctx.min)); 2779566063dSJacob Faibussowitsch PetscCall(VecStrideMax(xin, zctx.k, NULL, &zctx.max)); 27867dd0837SBarry Smith if (zctx.k < nbounds) { 279f3f0eb19SBarry Smith zctx.min = bounds[2 * zctx.k]; 280f3f0eb19SBarry Smith zctx.max = bounds[2 * zctx.k + 1]; 28167dd0837SBarry Smith } 28247c6ae99SBarry Smith if (zctx.min == zctx.max) { 28347c6ae99SBarry Smith zctx.min -= 1.e-12; 28447c6ae99SBarry Smith zctx.max += 1.e-12; 28547c6ae99SBarry Smith } 2869566063dSJacob Faibussowitsch PetscCall(PetscInfo(da, "DMDA 2d contour plot min %g max %g\n", (double)zctx.min, (double)zctx.max)); 28747c6ae99SBarry Smith 288832b7cebSLisandro Dalcin if (useports) { 2899566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports, i)); 290832b7cebSLisandro Dalcin } else { 291832b7cebSLisandro Dalcin const char *title; 2929566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, i, &draw)); 2939566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da, zctx.k, &title)); 2949566063dSJacob Faibussowitsch if (title) PetscCall(PetscDrawSetTitle(draw, title)); 295832b7cebSLisandro Dalcin } 296832b7cebSLisandro Dalcin 2979566063dSJacob Faibussowitsch PetscCall(PetscDrawGetPopup(draw, &popup)); 2989566063dSJacob Faibussowitsch PetscCall(PetscDrawScalePopup(popup, zctx.min, zctx.max)); 2999566063dSJacob Faibussowitsch PetscCall(PetscDrawSetCoordinates(draw, coors[0], coors[1], coors[2], coors[3])); 3009566063dSJacob Faibussowitsch PetscCall(PetscDrawZoom(draw, VecView_MPI_Draw_DA2d_Zoom, &zctx)); 3019566063dSJacob Faibussowitsch if (!useports) PetscCall(PetscDrawSave(draw)); 30247c6ae99SBarry Smith } 303832b7cebSLisandro Dalcin if (useports) { 3049566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 3059566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 306832b7cebSLisandro Dalcin } 307832b7cebSLisandro Dalcin 3089566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsDestroy(ports)); 3099566063dSJacob Faibussowitsch PetscCall(PetscFree(displayfields)); 3109566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xcoorl, &zctx.xy)); 3119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xlocal, &zctx.v)); 3123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 31347c6ae99SBarry Smith } 31447c6ae99SBarry Smith 3150da24e51SJuha Jäykkä #if defined(PETSC_HAVE_HDF5) 316d71ae5a4SJacob Faibussowitsch static PetscErrorCode VecGetHDF5ChunkSize(DM_DA *da, Vec xin, PetscInt dimension, PetscInt timestep, hsize_t *chunkDims) 317d71ae5a4SJacob Faibussowitsch { 318d4190030SJed Brown PetscMPIInt comm_size; 319561a82e7SJed Brown hsize_t chunk_size, target_size, dim; 320561a82e7SJed Brown hsize_t vec_size = sizeof(PetscScalar) * da->M * da->N * da->P * da->w; 321561a82e7SJed Brown hsize_t avg_local_vec_size, KiB = 1024, MiB = KiB * KiB, GiB = MiB * KiB, min_size = MiB; 322561a82e7SJed Brown hsize_t max_chunks = 64 * KiB; /* HDF5 internal limitation */ 323561a82e7SJed Brown hsize_t max_chunk_size = 4 * GiB; /* HDF5 internal limitation */ 32484179ffaSJed Brown PetscInt zslices = da->p, yslices = da->n, xslices = da->m; 3250da24e51SJuha Jäykkä 326cb5c4f94SJuha Jäykkä PetscFunctionBegin; 3279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)xin), &comm_size)); 328faa75363SBarry Smith avg_local_vec_size = (hsize_t)PetscCeilInt(vec_size, comm_size); /* we will attempt to use this as the chunk size */ 3290da24e51SJuha Jäykkä 330faa75363SBarry Smith target_size = (hsize_t)PetscMin((PetscInt64)vec_size, PetscMin((PetscInt64)max_chunk_size, PetscMax((PetscInt64)avg_local_vec_size, PetscMax(PetscCeilInt64(vec_size, max_chunks), (PetscInt64)min_size)))); 3317d310018SBarry Smith /* following line uses sizeof(PetscReal) instead of sizeof(PetscScalar) because the last dimension of chunkDims[] captures the 2* when complex numbers are being used */ 3327d310018SBarry Smith chunk_size = (hsize_t)PetscMax(1, chunkDims[0]) * PetscMax(1, chunkDims[1]) * PetscMax(1, chunkDims[2]) * PetscMax(1, chunkDims[3]) * PetscMax(1, chunkDims[4]) * PetscMax(1, chunkDims[5]) * sizeof(PetscReal); 3330da24e51SJuha Jäykkä 334cb5c4f94SJuha Jäykkä /* 335cb5c4f94SJuha Jäykkä if size/rank > max_chunk_size, we need radical measures: even going down to 336cb5c4f94SJuha Jäykkä avg_local_vec_size is not enough, so we simply use chunk size of 4 GiB no matter 337cb5c4f94SJuha Jäykkä what, composed in the most efficient way possible. 338cb5c4f94SJuha Jäykkä N.B. this minimises the number of chunks, which may or may not be the optimal 339cb5c4f94SJuha Jäykkä solution. In a BG, for example, the optimal solution is probably to make # chunks = # 340cb5c4f94SJuha Jäykkä IO nodes involved, but this author has no access to a BG to figure out how to 341cb5c4f94SJuha Jäykkä reliably find the right number. And even then it may or may not be enough. 342cb5c4f94SJuha Jäykkä */ 3430da24e51SJuha Jäykkä if (avg_local_vec_size > max_chunk_size) { 344cb5c4f94SJuha Jäykkä /* check if we can just split local z-axis: is that enough? */ 345faa75363SBarry Smith zslices = PetscCeilInt(vec_size, da->p * max_chunk_size) * zslices; 3460da24e51SJuha Jäykkä if (zslices > da->P) { 347cb5c4f94SJuha Jäykkä /* lattice is too large in xy-directions, splitting z only is not enough */ 3480da24e51SJuha Jäykkä zslices = da->P; 349faa75363SBarry Smith yslices = PetscCeilInt(vec_size, zslices * da->n * max_chunk_size) * yslices; 3500da24e51SJuha Jäykkä if (yslices > da->N) { 351cb5c4f94SJuha Jäykkä /* lattice is too large in x-direction, splitting along z, y is not enough */ 3520da24e51SJuha Jäykkä yslices = da->N; 353faa75363SBarry Smith xslices = PetscCeilInt(vec_size, zslices * yslices * da->m * max_chunk_size) * xslices; 3540da24e51SJuha Jäykkä } 3550da24e51SJuha Jäykkä } 3560da24e51SJuha Jäykkä dim = 0; 357ad540459SPierre Jolivet if (timestep >= 0) ++dim; 358cb5c4f94SJuha Jäykkä /* prefer to split z-axis, even down to planar slices */ 359c73cfb54SMatthew G. Knepley if (dimension == 3) { 360cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t)da->P / zslices; 361cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t)da->N / yslices; 362cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t)da->M / xslices; 3630da24e51SJuha Jäykkä } else { 364cb5c4f94SJuha Jäykkä /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */ 365cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t)da->N / yslices; 366cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t)da->M / xslices; 3670da24e51SJuha Jäykkä } 3680da24e51SJuha Jäykkä chunk_size = (hsize_t)PetscMax(1, chunkDims[0]) * PetscMax(1, chunkDims[1]) * PetscMax(1, chunkDims[2]) * PetscMax(1, chunkDims[3]) * PetscMax(1, chunkDims[4]) * PetscMax(1, chunkDims[5]) * sizeof(double); 3690da24e51SJuha Jäykkä } else { 3700da24e51SJuha Jäykkä if (target_size < chunk_size) { 371cb5c4f94SJuha Jäykkä /* only change the defaults if target_size < chunk_size */ 3720da24e51SJuha Jäykkä dim = 0; 373ad540459SPierre Jolivet if (timestep >= 0) ++dim; 374cb5c4f94SJuha Jäykkä /* prefer to split z-axis, even down to planar slices */ 375c73cfb54SMatthew G. Knepley if (dimension == 3) { 376cb5c4f94SJuha Jäykkä /* try splitting the z-axis to core-size bits, i.e. divide chunk size by # comm_size in z-direction */ 3770da24e51SJuha Jäykkä if (target_size >= chunk_size / da->p) { 378cb5c4f94SJuha Jäykkä /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */ 379faa75363SBarry Smith chunkDims[dim] = (hsize_t)PetscCeilInt(da->P, da->p); 3800da24e51SJuha Jäykkä } else { 381cb5c4f94SJuha Jäykkä /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be 382cb5c4f94SJuha Jäykkä radical and let everyone write all they've got */ 383faa75363SBarry Smith chunkDims[dim++] = (hsize_t)PetscCeilInt(da->P, da->p); 384faa75363SBarry Smith chunkDims[dim++] = (hsize_t)PetscCeilInt(da->N, da->n); 385faa75363SBarry Smith chunkDims[dim++] = (hsize_t)PetscCeilInt(da->M, da->m); 3860da24e51SJuha Jäykkä } 3870da24e51SJuha Jäykkä } else { 388cb5c4f94SJuha Jäykkä /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */ 3890da24e51SJuha Jäykkä if (target_size >= chunk_size / da->n) { 390cb5c4f94SJuha Jäykkä /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */ 391faa75363SBarry Smith chunkDims[dim] = (hsize_t)PetscCeilInt(da->N, da->n); 3920da24e51SJuha Jäykkä } else { 393cb5c4f94SJuha Jäykkä /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be 394cb5c4f94SJuha Jäykkä radical and let everyone write all they've got */ 395faa75363SBarry Smith chunkDims[dim++] = (hsize_t)PetscCeilInt(da->N, da->n); 396faa75363SBarry Smith chunkDims[dim++] = (hsize_t)PetscCeilInt(da->M, da->m); 3970da24e51SJuha Jäykkä } 3980da24e51SJuha Jäykkä } 3990da24e51SJuha Jäykkä chunk_size = (hsize_t)PetscMax(1, chunkDims[0]) * PetscMax(1, chunkDims[1]) * PetscMax(1, chunkDims[2]) * PetscMax(1, chunkDims[3]) * PetscMax(1, chunkDims[4]) * PetscMax(1, chunkDims[5]) * sizeof(double); 4000da24e51SJuha Jäykkä } else { 401cb5c4f94SJuha Jäykkä /* precomputed chunks are fine, we don't need to do anything */ 4020da24e51SJuha Jäykkä } 4030da24e51SJuha Jäykkä } 4043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4050da24e51SJuha Jäykkä } 4060da24e51SJuha Jäykkä #endif 40747c6ae99SBarry Smith 40847c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 40966976f2fSJacob Faibussowitsch static PetscErrorCode VecView_MPI_HDF5_DA(Vec xin, PetscViewer viewer) 410d71ae5a4SJacob Faibussowitsch { 4112d95be3dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 4129b2a5a86SJed Brown DM dm; 4139b2a5a86SJed Brown DM_DA *da; 41447c6ae99SBarry Smith hid_t filespace; /* file dataspace identifier */ 4158e2ae6d7SMichael Kraus hid_t chunkspace; /* chunk dataset property identifier */ 41647c6ae99SBarry Smith hid_t dset_id; /* dataset identifier */ 41747c6ae99SBarry Smith hid_t memspace; /* memory dataspace identifier */ 41847c6ae99SBarry Smith hid_t file_id; 41915214e8eSMatthew G Knepley hid_t group; 4209a0502c6SHåkon Strandenes hid_t memscalartype; /* scalar type for mem (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */ 4219a0502c6SHåkon Strandenes hid_t filescalartype; /* scalar type for file (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */ 422d9a4edebSJed Brown hsize_t dim; 4230da24e51SJuha Jäykkä hsize_t maxDims[6] = {0}, dims[6] = {0}, chunkDims[6] = {0}, count[6] = {0}, offset[6] = {0}; /* we depend on these being sane later on */ 424d7dd068bSVaclav Hapla PetscBool timestepping = PETSC_FALSE, dim2, spoutput; 425d7dd068bSVaclav Hapla PetscInt timestep = PETSC_MIN_INT, dimension; 4265edff71fSBarry Smith const PetscScalar *x; 4278e2ae6d7SMichael Kraus const char *vecname; 42847c6ae99SBarry Smith 42947c6ae99SBarry Smith PetscFunctionBegin; 4303014b61aSVaclav Hapla PetscCall(PetscViewerHDF5OpenGroup(viewer, NULL, &file_id, &group)); 4319566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5IsTimestepping(viewer, ×tepping)); 43248a46eb9SPierre Jolivet if (timestepping) PetscCall(PetscViewerHDF5GetTimestep(viewer, ×tep)); 4339566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetBaseDimension2(viewer, &dim2)); 4349566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetSPOutput(viewer, &spoutput)); 43515214e8eSMatthew G Knepley 4369566063dSJacob Faibussowitsch PetscCall(VecGetDM(xin, &dm)); 4377a8be351SBarry Smith PetscCheck(dm, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA"); 4389b2a5a86SJed Brown da = (DM_DA *)dm->data; 4399566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dimension)); 44047c6ae99SBarry Smith 4418e2ae6d7SMichael Kraus /* Create the dataspace for the dataset. 4428e2ae6d7SMichael Kraus * 4438e2ae6d7SMichael Kraus * dims - holds the current dimensions of the dataset 4448e2ae6d7SMichael Kraus * 4458e2ae6d7SMichael Kraus * maxDims - holds the maximum dimensions of the dataset (unlimited 4468e2ae6d7SMichael Kraus * for the number of time steps with the current dimensions for the 4478e2ae6d7SMichael Kraus * other dimensions; so only additional time steps can be added). 4488e2ae6d7SMichael Kraus * 4498e2ae6d7SMichael Kraus * chunkDims - holds the size of a single time step (required to 4508e2ae6d7SMichael Kraus * permit extending dataset). 4518e2ae6d7SMichael Kraus */ 4528e2ae6d7SMichael Kraus dim = 0; 4538e2ae6d7SMichael Kraus if (timestep >= 0) { 4548e2ae6d7SMichael Kraus dims[dim] = timestep + 1; 4558e2ae6d7SMichael Kraus maxDims[dim] = H5S_UNLIMITED; 4568e2ae6d7SMichael Kraus chunkDims[dim] = 1; 4578e2ae6d7SMichael Kraus ++dim; 4588e2ae6d7SMichael Kraus } 459c73cfb54SMatthew G. Knepley if (dimension == 3) { 4609566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(da->P, dims + dim)); 4618e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4628e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4638e2ae6d7SMichael Kraus ++dim; 4648e2ae6d7SMichael Kraus } 465c73cfb54SMatthew G. Knepley if (dimension > 1) { 4669566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(da->N, dims + dim)); 4678e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4688e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4698e2ae6d7SMichael Kraus ++dim; 4708e2ae6d7SMichael Kraus } 4719566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(da->M, dims + dim)); 4728e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4738e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4748e2ae6d7SMichael Kraus ++dim; 47582ea9b62SBarry Smith if (da->w > 1 || dim2) { 4769566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(da->w, dims + dim)); 4778e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4788e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4798e2ae6d7SMichael Kraus ++dim; 4808e2ae6d7SMichael Kraus } 48147c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 4828e2ae6d7SMichael Kraus dims[dim] = 2; 4838e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4848e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4858e2ae6d7SMichael Kraus ++dim; 48647c6ae99SBarry Smith #endif 4870da24e51SJuha Jäykkä 4889566063dSJacob Faibussowitsch PetscCall(VecGetHDF5ChunkSize(da, xin, dimension, timestep, chunkDims)); 4890da24e51SJuha Jäykkä 490792fecdfSBarry Smith PetscCallHDF5Return(filespace, H5Screate_simple, (dim, dims, maxDims)); 49147c6ae99SBarry Smith 49215214e8eSMatthew G Knepley #if defined(PETSC_USE_REAL_SINGLE) 4939a0502c6SHåkon Strandenes memscalartype = H5T_NATIVE_FLOAT; 4949a0502c6SHåkon Strandenes filescalartype = H5T_NATIVE_FLOAT; 49515214e8eSMatthew G Knepley #elif defined(PETSC_USE_REAL___FLOAT128) 49615214e8eSMatthew G Knepley #error "HDF5 output with 128 bit floats not supported." 497570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16) 498570b7f6dSBarry Smith #error "HDF5 output with 16 bit floats not supported." 49915214e8eSMatthew G Knepley #else 5009a0502c6SHåkon Strandenes memscalartype = H5T_NATIVE_DOUBLE; 5019a0502c6SHåkon Strandenes if (spoutput == PETSC_TRUE) filescalartype = H5T_NATIVE_FLOAT; 5029a0502c6SHåkon Strandenes else filescalartype = H5T_NATIVE_DOUBLE; 50315214e8eSMatthew G Knepley #endif 50415214e8eSMatthew G Knepley 50547c6ae99SBarry Smith /* Create the dataset with default properties and close filespace */ 5069566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)xin, &vecname)); 50715214e8eSMatthew G Knepley if (!H5Lexists(group, vecname, H5P_DEFAULT)) { 5088e2ae6d7SMichael Kraus /* Create chunk */ 509792fecdfSBarry Smith PetscCallHDF5Return(chunkspace, H5Pcreate, (H5P_DATASET_CREATE)); 510792fecdfSBarry Smith PetscCallHDF5(H5Pset_chunk, (chunkspace, dim, chunkDims)); 5118e2ae6d7SMichael Kraus 512792fecdfSBarry Smith PetscCallHDF5Return(dset_id, H5Dcreate2, (group, vecname, filescalartype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT)); 51315214e8eSMatthew G Knepley } else { 514792fecdfSBarry Smith PetscCallHDF5Return(dset_id, H5Dopen2, (group, vecname, H5P_DEFAULT)); 515792fecdfSBarry Smith PetscCallHDF5(H5Dset_extent, (dset_id, dims)); 51615214e8eSMatthew G Knepley } 517792fecdfSBarry Smith PetscCallHDF5(H5Sclose, (filespace)); 51847c6ae99SBarry Smith 51947c6ae99SBarry Smith /* Each process defines a dataset and writes it to the hyperslab in the file */ 5208e2ae6d7SMichael Kraus dim = 0; 5218e2ae6d7SMichael Kraus if (timestep >= 0) { 5228e2ae6d7SMichael Kraus offset[dim] = timestep; 5238e2ae6d7SMichael Kraus ++dim; 5248e2ae6d7SMichael Kraus } 5259566063dSJacob Faibussowitsch if (dimension == 3) PetscCall(PetscHDF5IntCast(da->zs, offset + dim++)); 5269566063dSJacob Faibussowitsch if (dimension > 1) PetscCall(PetscHDF5IntCast(da->ys, offset + dim++)); 5279566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(da->xs / da->w, offset + dim++)); 52882ea9b62SBarry Smith if (da->w > 1 || dim2) offset[dim++] = 0; 52947c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 5308e2ae6d7SMichael Kraus offset[dim++] = 0; 53147c6ae99SBarry Smith #endif 5328e2ae6d7SMichael Kraus dim = 0; 5338e2ae6d7SMichael Kraus if (timestep >= 0) { 5348e2ae6d7SMichael Kraus count[dim] = 1; 5358e2ae6d7SMichael Kraus ++dim; 5368e2ae6d7SMichael Kraus } 5379566063dSJacob Faibussowitsch if (dimension == 3) PetscCall(PetscHDF5IntCast(da->ze - da->zs, count + dim++)); 5389566063dSJacob Faibussowitsch if (dimension > 1) PetscCall(PetscHDF5IntCast(da->ye - da->ys, count + dim++)); 5399566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast((da->xe - da->xs) / da->w, count + dim++)); 5409566063dSJacob Faibussowitsch if (da->w > 1 || dim2) PetscCall(PetscHDF5IntCast(da->w, count + dim++)); 54147c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 5428e2ae6d7SMichael Kraus count[dim++] = 2; 54347c6ae99SBarry Smith #endif 544792fecdfSBarry Smith PetscCallHDF5Return(memspace, H5Screate_simple, (dim, count, NULL)); 545792fecdfSBarry Smith PetscCallHDF5Return(filespace, H5Dget_space, (dset_id)); 546792fecdfSBarry Smith PetscCallHDF5(H5Sselect_hyperslab, (filespace, H5S_SELECT_SET, offset, NULL, count, NULL)); 54747c6ae99SBarry Smith 5489566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xin, &x)); 549792fecdfSBarry Smith PetscCallHDF5(H5Dwrite, (dset_id, memscalartype, memspace, filespace, hdf5->dxpl_id, x)); 550792fecdfSBarry Smith PetscCallHDF5(H5Fflush, (file_id, H5F_SCOPE_GLOBAL)); 5519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xin, &x)); 55247c6ae99SBarry Smith 553e0552f36Ssajid__ali #if defined(PETSC_USE_COMPLEX) 554e0552f36Ssajid__ali { 555e0552f36Ssajid__ali PetscBool tru = PETSC_TRUE; 5569566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)xin, "complex", PETSC_BOOL, &tru)); 557e0552f36Ssajid__ali } 558e0552f36Ssajid__ali #endif 55948a46eb9SPierre Jolivet if (timestepping) PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)xin, "timestepping", PETSC_BOOL, ×tepping)); 560e0552f36Ssajid__ali 56147c6ae99SBarry Smith /* Close/release resources */ 562ad540459SPierre Jolivet if (group != file_id) PetscCallHDF5(H5Gclose, (group)); 563792fecdfSBarry Smith PetscCallHDF5(H5Sclose, (filespace)); 564792fecdfSBarry Smith PetscCallHDF5(H5Sclose, (memspace)); 565792fecdfSBarry Smith PetscCallHDF5(H5Dclose, (dset_id)); 5669566063dSJacob Faibussowitsch PetscCall(PetscInfo(xin, "Wrote Vec object with name %s\n", vecname)); 5673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56847c6ae99SBarry Smith } 56947c6ae99SBarry Smith #endif 57047c6ae99SBarry Smith 57109573ac7SBarry Smith extern PetscErrorCode VecView_MPI_Draw_DA1d(Vec, PetscViewer); 57247c6ae99SBarry Smith 57347c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 574d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAArrayMPIIO(DM da, PetscViewer viewer, Vec xin, PetscBool write) 575d71ae5a4SJacob Faibussowitsch { 57647c6ae99SBarry Smith MPI_File mfdes; 57747c6ae99SBarry Smith PetscMPIInt gsizes[4], lsizes[4], lstarts[4], asiz, dof; 57847c6ae99SBarry Smith MPI_Datatype view; 57947c6ae99SBarry Smith const PetscScalar *array; 58047c6ae99SBarry Smith MPI_Offset off; 58147c6ae99SBarry Smith MPI_Aint ub, ul; 58247c6ae99SBarry Smith PetscInt type, rows, vecrows, tr[2]; 58347c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 5840c169764Sdmay PetscBool skipheader; 58547c6ae99SBarry Smith 58647c6ae99SBarry Smith PetscFunctionBegin; 5879566063dSJacob Faibussowitsch PetscCall(VecGetSize(xin, &vecrows)); 5889566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipheader)); 58947c6ae99SBarry Smith if (!write) { 59047c6ae99SBarry Smith /* Read vector header. */ 5910c169764Sdmay if (!skipheader) { 5929566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, tr, 2, NULL, PETSC_INT)); 59347c6ae99SBarry Smith type = tr[0]; 59447c6ae99SBarry Smith rows = tr[1]; 59508401ef6SPierre Jolivet PetscCheck(type == VEC_FILE_CLASSID, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONG, "Not vector next in file"); 59608401ef6SPierre Jolivet PetscCheck(rows == vecrows, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Vector in file not same size as DMDA vector"); 5970c169764Sdmay } 59847c6ae99SBarry Smith } else { 59947c6ae99SBarry Smith tr[0] = VEC_FILE_CLASSID; 60047c6ae99SBarry Smith tr[1] = vecrows; 60148a46eb9SPierre Jolivet if (!skipheader) PetscCall(PetscViewerBinaryWrite(viewer, tr, 2, PETSC_INT)); 6020c169764Sdmay } 60347c6ae99SBarry Smith 6049566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->w, &dof)); 6054dc2109aSBarry Smith gsizes[0] = dof; 6069566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->M, gsizes + 1)); 6079566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->N, gsizes + 2)); 6089566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->P, gsizes + 3)); 6094dc2109aSBarry Smith lsizes[0] = dof; 6109566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast((dd->xe - dd->xs) / dof, lsizes + 1)); 6119566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->ye - dd->ys, lsizes + 2)); 6129566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->ze - dd->zs, lsizes + 3)); 6134dc2109aSBarry Smith lstarts[0] = 0; 6149566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->xs / dof, lstarts + 1)); 6159566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->ys, lstarts + 2)); 6169566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->zs, lstarts + 3)); 6179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_subarray(da->dim + 1, gsizes, lsizes, lstarts, MPI_ORDER_FORTRAN, MPIU_SCALAR, &view)); 6189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&view)); 61947c6ae99SBarry Smith 6209566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetMPIIODescriptor(viewer, &mfdes)); 6219566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetMPIIOOffset(viewer, &off)); 6229566063dSJacob Faibussowitsch PetscCallMPI(MPI_File_set_view(mfdes, off, MPIU_SCALAR, view, (char *)"native", MPI_INFO_NULL)); 6239566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xin, &array)); 62447c6ae99SBarry Smith asiz = lsizes[1] * (lsizes[2] > 0 ? lsizes[2] : 1) * (lsizes[3] > 0 ? lsizes[3] : 1) * dof; 6251baa6e33SBarry Smith if (write) PetscCall(MPIU_File_write_all(mfdes, (PetscScalar *)array, asiz, MPIU_SCALAR, MPI_STATUS_IGNORE)); 6261baa6e33SBarry Smith else PetscCall(MPIU_File_read_all(mfdes, (PetscScalar *)array, asiz, MPIU_SCALAR, MPI_STATUS_IGNORE)); 6279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_get_extent(view, &ul, &ub)); 6289566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryAddMPIIOOffset(viewer, ub)); 6299566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xin, &array)); 6309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&view)); 6313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 63247c6ae99SBarry Smith } 63347c6ae99SBarry Smith #endif 63447c6ae99SBarry Smith 635d71ae5a4SJacob Faibussowitsch PetscErrorCode VecView_MPI_DA(Vec xin, PetscViewer viewer) 636d71ae5a4SJacob Faibussowitsch { 6379a42bb27SBarry Smith DM da; 63847c6ae99SBarry Smith PetscInt dim; 63947c6ae99SBarry Smith Vec natural; 6408135c375SStefano Zampini PetscBool isdraw, isvtk, isglvis; 64147c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 64247c6ae99SBarry Smith PetscBool ishdf5; 64347c6ae99SBarry Smith #endif 6443f3fd955SJed Brown const char *prefix, *name; 645a261c58fSBarry Smith PetscViewerFormat format; 64647c6ae99SBarry Smith 64747c6ae99SBarry Smith PetscFunctionBegin; 6489566063dSJacob Faibussowitsch PetscCall(VecGetDM(xin, &da)); 6497a8be351SBarry Smith PetscCheck(da, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA"); 6509566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 6519566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk)); 65247c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 6539566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5)); 65447c6ae99SBarry Smith #endif 6559566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis)); 65647c6ae99SBarry Smith if (isdraw) { 6579566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 65847c6ae99SBarry Smith if (dim == 1) { 6599566063dSJacob Faibussowitsch PetscCall(VecView_MPI_Draw_DA1d(xin, viewer)); 66047c6ae99SBarry Smith } else if (dim == 2) { 6619566063dSJacob Faibussowitsch PetscCall(VecView_MPI_Draw_DA2d(xin, viewer)); 66263a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Cannot graphically view vector associated with this dimensional DMDA %" PetscInt_FMT, dim); 66301d7c5c3SPatrick Sanan } else if (isvtk) { /* Duplicate the Vec */ 6644061b8bfSJed Brown Vec Y; 6659566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xin, &Y)); 666b51b94faSJed Brown if (((PetscObject)xin)->name) { 667b51b94faSJed Brown /* If xin was named, copy the name over to Y. The duplicate names are safe because nobody else will ever see Y. */ 6689566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)Y, ((PetscObject)xin)->name)); 669b51b94faSJed Brown } 6709566063dSJacob Faibussowitsch PetscCall(VecCopy(xin, Y)); 671a8f87f1dSPatrick Sanan { 672a8f87f1dSPatrick Sanan PetscObject dmvtk; 673a8f87f1dSPatrick Sanan PetscBool compatible, compatibleSet; 6749566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKGetDM(viewer, &dmvtk)); 675a8f87f1dSPatrick Sanan if (dmvtk) { 676064a246eSJacob Faibussowitsch PetscValidHeaderSpecific((DM)dmvtk, DM_CLASSID, 2); 6779566063dSJacob Faibussowitsch PetscCall(DMGetCompatibility(da, (DM)dmvtk, &compatible, &compatibleSet)); 6787a8be351SBarry Smith PetscCheck(compatibleSet && compatible, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "Cannot confirm compatibility of DMs associated with Vecs viewed in the same VTK file. Check that grids are the same."); 679a8f87f1dSPatrick Sanan } 6809566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKAddField(viewer, (PetscObject)da, DMDAVTKWriteAll, PETSC_DEFAULT, PETSC_VTK_POINT_FIELD, PETSC_FALSE, (PetscObject)Y)); 681a8f87f1dSPatrick Sanan } 68247c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 68347c6ae99SBarry Smith } else if (ishdf5) { 6849566063dSJacob Faibussowitsch PetscCall(VecView_MPI_HDF5_DA(xin, viewer)); 68547c6ae99SBarry Smith #endif 6868135c375SStefano Zampini } else if (isglvis) { 6879566063dSJacob Faibussowitsch PetscCall(VecView_GLVis(xin, viewer)); 68847c6ae99SBarry Smith } else { 68947c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 69047c6ae99SBarry Smith PetscBool isbinary, isMPIIO; 69147c6ae99SBarry Smith 6929566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 69347c6ae99SBarry Smith if (isbinary) { 6949566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &isMPIIO)); 69547c6ae99SBarry Smith if (isMPIIO) { 6969566063dSJacob Faibussowitsch PetscCall(DMDAArrayMPIIO(da, viewer, xin, PETSC_TRUE)); 6973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 69847c6ae99SBarry Smith } 69947c6ae99SBarry Smith } 70047c6ae99SBarry Smith #endif 70147c6ae99SBarry Smith 70247c6ae99SBarry Smith /* call viewer on natural ordering */ 7039566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)xin, &prefix)); 7049566063dSJacob Faibussowitsch PetscCall(DMDACreateNaturalVector(da, &natural)); 7059566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)natural, prefix)); 7069566063dSJacob Faibussowitsch PetscCall(DMDAGlobalToNaturalBegin(da, xin, INSERT_VALUES, natural)); 7079566063dSJacob Faibussowitsch PetscCall(DMDAGlobalToNaturalEnd(da, xin, INSERT_VALUES, natural)); 7089566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)xin, &name)); 7099566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)natural, name)); 710a261c58fSBarry Smith 7119566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 712a261c58fSBarry Smith if (format == PETSC_VIEWER_BINARY_MATLAB) { 713a261c58fSBarry Smith /* temporarily remove viewer format so it won't trigger in the VecView() */ 7149566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_DEFAULT)); 715a261c58fSBarry Smith } 716a261c58fSBarry Smith 71723f88b65SBarry Smith ((PetscObject)natural)->donotPetscObjectPrintClassNamePrefixType = PETSC_TRUE; 7189566063dSJacob Faibussowitsch PetscCall(VecView(natural, viewer)); 71923f88b65SBarry Smith ((PetscObject)natural)->donotPetscObjectPrintClassNamePrefixType = PETSC_FALSE; 720a261c58fSBarry Smith 721a261c58fSBarry Smith if (format == PETSC_VIEWER_BINARY_MATLAB) { 722a261c58fSBarry Smith MPI_Comm comm; 723a261c58fSBarry Smith FILE *info; 724a261c58fSBarry Smith const char *fieldname; 725da88d4d4SJed Brown char fieldbuf[256]; 726a261c58fSBarry Smith PetscInt dim, ni, nj, nk, pi, pj, pk, dof, n; 727a261c58fSBarry Smith 728a261c58fSBarry Smith /* set the viewer format back into the viewer */ 7299566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 7309566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm)); 7319566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetInfoPointer(viewer, &info)); 7329566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, &ni, &nj, &nk, &pi, &pj, &pk, &dof, NULL, NULL, NULL, NULL, NULL)); 7339566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n")); 7349566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, info, "#$$ tmp = PetscBinaryRead(fd); \n")); 73548a46eb9SPierre Jolivet if (dim == 1) PetscCall(PetscFPrintf(comm, info, "#$$ tmp = reshape(tmp,%" PetscInt_FMT ",%" PetscInt_FMT ");\n", dof, ni)); 73648a46eb9SPierre Jolivet if (dim == 2) PetscCall(PetscFPrintf(comm, info, "#$$ tmp = reshape(tmp,%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ");\n", dof, ni, nj)); 73748a46eb9SPierre Jolivet if (dim == 3) PetscCall(PetscFPrintf(comm, info, "#$$ tmp = reshape(tmp,%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ");\n", dof, ni, nj, nk)); 738a261c58fSBarry Smith 739a261c58fSBarry Smith for (n = 0; n < dof; n++) { 7409566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da, n, &fieldname)); 741da88d4d4SJed Brown if (!fieldname || !fieldname[0]) { 74263a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(fieldbuf, sizeof fieldbuf, "field%" PetscInt_FMT, n)); 743da88d4d4SJed Brown fieldname = fieldbuf; 744a261c58fSBarry Smith } 74548a46eb9SPierre Jolivet if (dim == 1) PetscCall(PetscFPrintf(comm, info, "#$$ Set.%s.%s = squeeze(tmp(%" PetscInt_FMT ",:))';\n", name, fieldname, n + 1)); 74648a46eb9SPierre Jolivet if (dim == 2) PetscCall(PetscFPrintf(comm, info, "#$$ Set.%s.%s = squeeze(tmp(%" PetscInt_FMT ",:,:))';\n", name, fieldname, n + 1)); 74748a46eb9SPierre Jolivet if (dim == 3) PetscCall(PetscFPrintf(comm, info, "#$$ Set.%s.%s = permute(squeeze(tmp(%" PetscInt_FMT ",:,:,:)),[2 1 3]);\n", name, fieldname, n + 1)); 748a261c58fSBarry Smith } 7499566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, info, "#$$ clear tmp; \n")); 7509566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n")); 751a261c58fSBarry Smith } 752a261c58fSBarry Smith 7539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&natural)); 75447c6ae99SBarry Smith } 7553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 75647c6ae99SBarry Smith } 75747c6ae99SBarry Smith 75847c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 75966976f2fSJacob Faibussowitsch static PetscErrorCode VecLoad_HDF5_DA(Vec xin, PetscViewer viewer) 760d71ae5a4SJacob Faibussowitsch { 7612d95be3dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 7629a42bb27SBarry Smith DM da; 76315b861d2SVaclav Hapla int dim, rdim; 764ec7ae49cSHåkon Strandenes hsize_t dims[6] = {0}, count[6] = {0}, offset[6] = {0}; 765d7dd068bSVaclav Hapla PetscBool dim2 = PETSC_FALSE, timestepping = PETSC_FALSE; 766d7dd068bSVaclav Hapla PetscInt dimension, timestep = PETSC_MIN_INT, dofInd; 76747c6ae99SBarry Smith PetscScalar *x; 76847c6ae99SBarry Smith const char *vecname; 76947c6ae99SBarry Smith hid_t filespace; /* file dataspace identifier */ 77047c6ae99SBarry Smith hid_t dset_id; /* dataset identifier */ 77147c6ae99SBarry Smith hid_t memspace; /* memory dataspace identifier */ 772bfd264e7SBarry Smith hid_t file_id, group; 7737bcbaff4SBarry Smith hid_t scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */ 7749c7c4993SBarry Smith DM_DA *dd; 77547c6ae99SBarry Smith 77647c6ae99SBarry Smith PetscFunctionBegin; 7777bcbaff4SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 7787bcbaff4SBarry Smith scalartype = H5T_NATIVE_FLOAT; 7797bcbaff4SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128) 7807bcbaff4SBarry Smith #error "HDF5 output with 128 bit floats not supported." 781570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16) 782570b7f6dSBarry Smith #error "HDF5 output with 16 bit floats not supported." 7837bcbaff4SBarry Smith #else 7847bcbaff4SBarry Smith scalartype = H5T_NATIVE_DOUBLE; 7857bcbaff4SBarry Smith #endif 7867bcbaff4SBarry Smith 7873014b61aSVaclav Hapla PetscCall(PetscViewerHDF5OpenGroup(viewer, NULL, &file_id, &group)); 7889566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)xin, &vecname)); 7899566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckTimestepping_Internal(viewer, vecname)); 7909566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5IsTimestepping(viewer, ×tepping)); 79148a46eb9SPierre Jolivet if (timestepping) PetscCall(PetscViewerHDF5GetTimestep(viewer, ×tep)); 7929566063dSJacob Faibussowitsch PetscCall(VecGetDM(xin, &da)); 7939c7c4993SBarry Smith dd = (DM_DA *)da->data; 7949566063dSJacob Faibussowitsch PetscCall(DMGetDimension(da, &dimension)); 79547c6ae99SBarry Smith 796ec7ae49cSHåkon Strandenes /* Open dataset */ 797792fecdfSBarry Smith PetscCallHDF5Return(dset_id, H5Dopen2, (group, vecname, H5P_DEFAULT)); 79847c6ae99SBarry Smith 799ec7ae49cSHåkon Strandenes /* Retrieve the dataspace for the dataset */ 800792fecdfSBarry Smith PetscCallHDF5Return(filespace, H5Dget_space, (dset_id)); 801792fecdfSBarry Smith PetscCallHDF5Return(rdim, H5Sget_simple_extent_dims, (filespace, dims, NULL)); 802ec7ae49cSHåkon Strandenes 803ec7ae49cSHåkon Strandenes /* Expected dimension for holding the dof's */ 80447c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 805ec7ae49cSHåkon Strandenes dofInd = rdim - 2; 806ec7ae49cSHåkon Strandenes #else 807ec7ae49cSHåkon Strandenes dofInd = rdim - 1; 80847c6ae99SBarry Smith #endif 809ec7ae49cSHåkon Strandenes 810ec7ae49cSHåkon Strandenes /* The expected number of dimensions, assuming basedimension2 = false */ 811ec7ae49cSHåkon Strandenes dim = dimension; 812ec7ae49cSHåkon Strandenes if (dd->w > 1) ++dim; 813ec7ae49cSHåkon Strandenes if (timestep >= 0) ++dim; 81447c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 815ec7ae49cSHåkon Strandenes ++dim; 81647c6ae99SBarry Smith #endif 817ec7ae49cSHåkon Strandenes 818ec7ae49cSHåkon Strandenes /* In this case the input dataset have one extra, unexpected dimension. */ 819ec7ae49cSHåkon Strandenes if (rdim == dim + 1) { 820ec7ae49cSHåkon Strandenes /* In this case the block size unity */ 821ec7ae49cSHåkon Strandenes if (dd->w == 1 && dims[dofInd] == 1) dim2 = PETSC_TRUE; 822ec7ae49cSHåkon Strandenes 823ec7ae49cSHåkon Strandenes /* Special error message for the case where dof does not match the input file */ 82463a3b9bcSJacob Faibussowitsch else PetscCheck(dd->w == (PetscInt)dims[dofInd], PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Number of dofs in file is %" PetscInt_FMT ", not %" PetscInt_FMT " as expected", (PetscInt)dims[dofInd], dd->w); 825ec7ae49cSHåkon Strandenes 826ec7ae49cSHåkon Strandenes /* Other cases where rdim != dim cannot be handled currently */ 82763a3b9bcSJacob Faibussowitsch } else PetscCheck(rdim == dim, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Dimension of array in file is %d, not %d as expected with dof = %" PetscInt_FMT, rdim, dim, dd->w); 828ec7ae49cSHåkon Strandenes 829ec7ae49cSHåkon Strandenes /* Set up the hyperslab size */ 830ec7ae49cSHåkon Strandenes dim = 0; 831ec7ae49cSHåkon Strandenes if (timestep >= 0) { 832ec7ae49cSHåkon Strandenes offset[dim] = timestep; 833ec7ae49cSHåkon Strandenes count[dim] = 1; 834ec7ae49cSHåkon Strandenes ++dim; 835ec7ae49cSHåkon Strandenes } 836ec7ae49cSHåkon Strandenes if (dimension == 3) { 8379566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(dd->zs, offset + dim)); 8389566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(dd->ze - dd->zs, count + dim)); 839ec7ae49cSHåkon Strandenes ++dim; 840ec7ae49cSHåkon Strandenes } 841ec7ae49cSHåkon Strandenes if (dimension > 1) { 8429566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(dd->ys, offset + dim)); 8439566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(dd->ye - dd->ys, count + dim)); 844ec7ae49cSHåkon Strandenes ++dim; 845ec7ae49cSHåkon Strandenes } 8469566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(dd->xs / dd->w, offset + dim)); 8479566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast((dd->xe - dd->xs) / dd->w, count + dim)); 848ec7ae49cSHåkon Strandenes ++dim; 849ec7ae49cSHåkon Strandenes if (dd->w > 1 || dim2) { 850ec7ae49cSHåkon Strandenes offset[dim] = 0; 8519566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(dd->w, count + dim)); 852ec7ae49cSHåkon Strandenes ++dim; 853ec7ae49cSHåkon Strandenes } 854ec7ae49cSHåkon Strandenes #if defined(PETSC_USE_COMPLEX) 855ec7ae49cSHåkon Strandenes offset[dim] = 0; 856ec7ae49cSHåkon Strandenes count[dim] = 2; 857ec7ae49cSHåkon Strandenes ++dim; 858ec7ae49cSHåkon Strandenes #endif 859ec7ae49cSHåkon Strandenes 860ec7ae49cSHåkon Strandenes /* Create the memory and filespace */ 861792fecdfSBarry Smith PetscCallHDF5Return(memspace, H5Screate_simple, (dim, count, NULL)); 862792fecdfSBarry Smith PetscCallHDF5(H5Sselect_hyperslab, (filespace, H5S_SELECT_SET, offset, NULL, count, NULL)); 86347c6ae99SBarry Smith 8649566063dSJacob Faibussowitsch PetscCall(VecGetArray(xin, &x)); 865792fecdfSBarry Smith PetscCallHDF5(H5Dread, (dset_id, scalartype, memspace, filespace, hdf5->dxpl_id, x)); 8669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xin, &x)); 86747c6ae99SBarry Smith 86847c6ae99SBarry Smith /* Close/release resources */ 869ad540459SPierre Jolivet if (group != file_id) PetscCallHDF5(H5Gclose, (group)); 870792fecdfSBarry Smith PetscCallHDF5(H5Sclose, (filespace)); 871792fecdfSBarry Smith PetscCallHDF5(H5Sclose, (memspace)); 872792fecdfSBarry Smith PetscCallHDF5(H5Dclose, (dset_id)); 8733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 87447c6ae99SBarry Smith } 87547c6ae99SBarry Smith #endif 87647c6ae99SBarry Smith 87766976f2fSJacob Faibussowitsch static PetscErrorCode VecLoad_Binary_DA(Vec xin, PetscViewer viewer) 878d71ae5a4SJacob Faibussowitsch { 8799a42bb27SBarry Smith DM da; 88047c6ae99SBarry Smith Vec natural; 88147c6ae99SBarry Smith const char *prefix; 88247c6ae99SBarry Smith PetscInt bs; 88347c6ae99SBarry Smith PetscBool flag; 88447c6ae99SBarry Smith DM_DA *dd; 88547c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 88647c6ae99SBarry Smith PetscBool isMPIIO; 88747c6ae99SBarry Smith #endif 88847c6ae99SBarry Smith 88947c6ae99SBarry Smith PetscFunctionBegin; 8909566063dSJacob Faibussowitsch PetscCall(VecGetDM(xin, &da)); 89147c6ae99SBarry Smith dd = (DM_DA *)da->data; 89247c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 8939566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &isMPIIO)); 89447c6ae99SBarry Smith if (isMPIIO) { 8959566063dSJacob Faibussowitsch PetscCall(DMDAArrayMPIIO(da, viewer, xin, PETSC_FALSE)); 8963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 89747c6ae99SBarry Smith } 89847c6ae99SBarry Smith #endif 89947c6ae99SBarry Smith 9009566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)xin, &prefix)); 9019566063dSJacob Faibussowitsch PetscCall(DMDACreateNaturalVector(da, &natural)); 9029566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)natural, ((PetscObject)xin)->name)); 9039566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)natural, prefix)); 9049566063dSJacob Faibussowitsch PetscCall(VecLoad(natural, viewer)); 9059566063dSJacob Faibussowitsch PetscCall(DMDANaturalToGlobalBegin(da, natural, INSERT_VALUES, xin)); 9069566063dSJacob Faibussowitsch PetscCall(DMDANaturalToGlobalEnd(da, natural, INSERT_VALUES, xin)); 9079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&natural)); 9089566063dSJacob Faibussowitsch PetscCall(PetscInfo(xin, "Loading vector from natural ordering into DMDA\n")); 9099566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, ((PetscObject)xin)->prefix, "-vecload_block_size", &bs, &flag)); 91048a46eb9SPierre Jolivet if (flag && bs != dd->w) PetscCall(PetscInfo(xin, "Block size in file %" PetscInt_FMT " not equal to DMDA's dof %" PetscInt_FMT "\n", bs, dd->w)); 9113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 91247c6ae99SBarry Smith } 91347c6ae99SBarry Smith 914d71ae5a4SJacob Faibussowitsch PetscErrorCode VecLoad_Default_DA(Vec xin, PetscViewer viewer) 915d71ae5a4SJacob Faibussowitsch { 9169a42bb27SBarry Smith DM da; 91747c6ae99SBarry Smith PetscBool isbinary; 91847c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 91947c6ae99SBarry Smith PetscBool ishdf5; 92047c6ae99SBarry Smith #endif 92147c6ae99SBarry Smith 92247c6ae99SBarry Smith PetscFunctionBegin; 9239566063dSJacob Faibussowitsch PetscCall(VecGetDM(xin, &da)); 9247a8be351SBarry Smith PetscCheck(da, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA"); 92547c6ae99SBarry Smith 92647c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 9279566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5)); 92847c6ae99SBarry Smith #endif 9299566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 93047c6ae99SBarry Smith 93147c6ae99SBarry Smith if (isbinary) { 9329566063dSJacob Faibussowitsch PetscCall(VecLoad_Binary_DA(xin, viewer)); 93347c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 93447c6ae99SBarry Smith } else if (ishdf5) { 9359566063dSJacob Faibussowitsch PetscCall(VecLoad_HDF5_DA(xin, viewer)); 93647c6ae99SBarry Smith #endif 93798921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Viewer type %s not supported for vector loading", ((PetscObject)viewer)->type_name); 9383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 93947c6ae99SBarry Smith } 940