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; 306497c311SBarry Smith PetscInt m, n, i, j, k, dof, id; 316497c311SBarry Smith int c1, c2, c3, c4; 32e5ab1681SLisandro Dalcin PetscReal min, max, x1, x2, x3, x4, y_1, y2, y3, y4; 33e5ab1681SLisandro Dalcin const PetscScalar *xy, *v; 3447c6ae99SBarry Smith 3547c6ae99SBarry Smith PetscFunctionBegin; 3647c6ae99SBarry Smith m = zctx->m; 3747c6ae99SBarry Smith n = zctx->n; 38e5ab1681SLisandro Dalcin dof = zctx->dof; 3947c6ae99SBarry Smith k = zctx->k; 4047c6ae99SBarry Smith xy = zctx->xy; 41e5ab1681SLisandro Dalcin v = zctx->v; 4247c6ae99SBarry Smith min = zctx->min; 43f3f0eb19SBarry Smith max = zctx->max; 4447c6ae99SBarry Smith 4547c6ae99SBarry Smith /* PetscDraw the contour plot patch */ 46d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 4747c6ae99SBarry Smith for (j = 0; j < n - 1; j++) { 4847c6ae99SBarry Smith for (i = 0; i < m - 1; i++) { 490e4fe250SBarry Smith id = i + j * m; 500e4fe250SBarry Smith x1 = PetscRealPart(xy[2 * id]); 510e4fe250SBarry Smith y_1 = PetscRealPart(xy[2 * id + 1]); 52e5ab1681SLisandro Dalcin c1 = PetscDrawRealToColor(PetscRealPart(v[k + dof * id]), min, max); 530e4fe250SBarry Smith 540e4fe250SBarry Smith id = i + j * m + 1; 550e4fe250SBarry Smith x2 = PetscRealPart(xy[2 * id]); 56a39a4c7dSLisandro Dalcin y2 = PetscRealPart(xy[2 * id + 1]); 57e5ab1681SLisandro Dalcin c2 = PetscDrawRealToColor(PetscRealPart(v[k + dof * id]), min, max); 580e4fe250SBarry Smith 590e4fe250SBarry Smith id = i + j * m + 1 + m; 60a39a4c7dSLisandro Dalcin x3 = PetscRealPart(xy[2 * id]); 610e4fe250SBarry Smith y3 = PetscRealPart(xy[2 * id + 1]); 62e5ab1681SLisandro Dalcin c3 = PetscDrawRealToColor(PetscRealPart(v[k + dof * id]), min, max); 630e4fe250SBarry Smith 640e4fe250SBarry Smith id = i + j * m + m; 65a39a4c7dSLisandro Dalcin x4 = PetscRealPart(xy[2 * id]); 66a39a4c7dSLisandro Dalcin y4 = PetscRealPart(xy[2 * id + 1]); 67e5ab1681SLisandro Dalcin c4 = PetscDrawRealToColor(PetscRealPart(v[k + dof * id]), min, max); 68f3f0eb19SBarry Smith 699566063dSJacob Faibussowitsch PetscCall(PetscDrawTriangle(draw, x1, y_1, x2, y2, x3, y3, c1, c2, c3)); 709566063dSJacob Faibussowitsch PetscCall(PetscDrawTriangle(draw, x1, y_1, x3, y3, x4, y4, c1, c3, c4)); 7147c6ae99SBarry Smith if (zctx->showgrid) { 729566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, x1, y_1, x2, y2, PETSC_DRAW_BLACK)); 739566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, x2, y2, x3, y3, PETSC_DRAW_BLACK)); 749566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, x3, y3, x4, y4, PETSC_DRAW_BLACK)); 759566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, x4, y4, x1, y_1, PETSC_DRAW_BLACK)); 7647c6ae99SBarry Smith } 7747c6ae99SBarry Smith } 7847c6ae99SBarry Smith } 797fe7c8feSLisandro Dalcin if (zctx->showaxis && !zctx->rank) { 80e5ab1681SLisandro Dalcin if (zctx->name0 || zctx->name1) { 81109c9344SBarry Smith PetscReal xl, yl, xr, yr, x, y; 829566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 83e5ab1681SLisandro Dalcin x = xl + .30 * (xr - xl); 84109c9344SBarry Smith xl = xl + .01 * (xr - xl); 85e5ab1681SLisandro Dalcin y = yr - .30 * (yr - yl); 86109c9344SBarry Smith yl = yl + .01 * (yr - yl); 879566063dSJacob Faibussowitsch if (zctx->name0) PetscCall(PetscDrawString(draw, x, yl, PETSC_DRAW_BLACK, zctx->name0)); 889566063dSJacob Faibussowitsch if (zctx->name1) PetscCall(PetscDrawStringVertical(draw, xl, y, PETSC_DRAW_BLACK, zctx->name1)); 89109c9344SBarry Smith } 900e4fe250SBarry Smith /* 910e4fe250SBarry Smith Ideally we would use the PetscDrawAxis object to manage displaying the coordinate limits 920e4fe250SBarry Smith but that may require some refactoring. 930e4fe250SBarry Smith */ 94e5ab1681SLisandro Dalcin { 957fe7c8feSLisandro Dalcin double xmin = (double)zctx->xmin, ymin = (double)zctx->ymin; 967fe7c8feSLisandro Dalcin double xmax = (double)zctx->xmax, ymax = (double)zctx->ymax; 979371c9d4SSatish Balay char value[16]; 989371c9d4SSatish Balay size_t len; 999371c9d4SSatish Balay PetscReal w; 1009566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(value, 16, "%0.2e", xmin)); 1019566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw, xmin, ymin - .05 * (ymax - ymin), PETSC_DRAW_BLACK, value)); 1029566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(value, 16, "%0.2e", xmax)); 1039566063dSJacob Faibussowitsch PetscCall(PetscStrlen(value, &len)); 1049566063dSJacob Faibussowitsch PetscCall(PetscDrawStringGetSize(draw, &w, NULL)); 1059566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw, xmax - len * w, ymin - .05 * (ymax - ymin), PETSC_DRAW_BLACK, value)); 1069566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(value, 16, "%0.2e", ymin)); 1079566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw, xmin - .05 * (xmax - xmin), ymin, PETSC_DRAW_BLACK, value)); 1089566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(value, 16, "%0.2e", ymax)); 1099566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw, xmin - .05 * (xmax - xmin), ymax, PETSC_DRAW_BLACK, value)); 110e5ab1681SLisandro Dalcin } 111e5ab1681SLisandro Dalcin } 112d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11447c6ae99SBarry Smith } 11547c6ae99SBarry Smith 11666976f2fSJacob Faibussowitsch static PetscErrorCode VecView_MPI_Draw_DA2d(Vec xin, PetscViewer viewer) 117d71ae5a4SJacob Faibussowitsch { 1189a42bb27SBarry Smith DM da, dac, dag; 119a74ba6f7SBarry Smith PetscInt N, s, M, w, ncoors = 4; 12047c6ae99SBarry Smith const PetscInt *lx, *ly; 121e5ab1681SLisandro Dalcin PetscReal coors[4]; 12247c6ae99SBarry Smith PetscDraw draw, popup; 12347c6ae99SBarry Smith PetscBool isnull, useports = PETSC_FALSE; 12447c6ae99SBarry Smith MPI_Comm comm; 12547c6ae99SBarry Smith Vec xlocal, xcoor, xcoorl; 126bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by; 127aa219208SBarry Smith DMDAStencilType st; 12847c6ae99SBarry Smith ZoomCtx zctx; 1290298fd71SBarry Smith PetscDrawViewPorts *ports = NULL; 13047c6ae99SBarry Smith PetscViewerFormat format; 13120d0051dSBarry Smith PetscInt *displayfields; 13267dd0837SBarry Smith PetscInt ndisplayfields, i, nbounds; 13367dd0837SBarry Smith const PetscReal *bounds; 13447c6ae99SBarry Smith 13547c6ae99SBarry Smith PetscFunctionBegin; 13647c6ae99SBarry Smith zctx.showgrid = PETSC_FALSE; 1377fe7c8feSLisandro Dalcin zctx.showaxis = PETSC_TRUE; 1388865f1eaSKarl Rupp 1399566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 1409566063dSJacob Faibussowitsch PetscCall(PetscDrawIsNull(draw, &isnull)); 1413ba16761SJacob Faibussowitsch if (isnull) PetscFunctionReturn(PETSC_SUCCESS); 1425b399a63SLisandro Dalcin 1439566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetBounds(viewer, &nbounds, &bounds)); 14447c6ae99SBarry Smith 1459566063dSJacob Faibussowitsch PetscCall(VecGetDM(xin, &da)); 1467a8be351SBarry Smith PetscCheck(da, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA"); 14747c6ae99SBarry Smith 1489566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)xin, &comm)); 1499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &zctx.rank)); 15047c6ae99SBarry Smith 1519566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, &M, &N, NULL, &zctx.m, &zctx.n, NULL, &w, &s, &bx, &by, NULL, &st)); 1529566063dSJacob Faibussowitsch PetscCall(DMDAGetOwnershipRanges(da, &lx, &ly, NULL)); 15347c6ae99SBarry Smith 15447c6ae99SBarry Smith /* 15547c6ae99SBarry Smith Obtain a sequential vector that is going to contain the local values plus ONE layer of 156aa219208SBarry Smith ghosted values to draw the graphics from. We also need its corresponding DMDA (dac) that will 157d5b43468SJose E. Roman update the local values plus ONE layer of ghost values. 15847c6ae99SBarry Smith */ 1599566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)da, "GraphicsGhosted", (PetscObject *)&xlocal)); 16047c6ae99SBarry Smith if (!xlocal) { 161bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_NONE || by != DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) { 16247c6ae99SBarry Smith /* 16347c6ae99SBarry Smith if original da is not of stencil width one, or periodic or not a box stencil then 164aa219208SBarry Smith create a special DMDA to handle one level of ghost points for graphics 16547c6ae99SBarry Smith */ 1669566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, zctx.m, zctx.n, w, 1, lx, ly, &dac)); 1679566063dSJacob Faibussowitsch PetscCall(DMSetUp(dac)); 1686aad120cSJose E. Roman PetscCall(PetscInfo(da, "Creating auxiliary DMDA for managing graphics ghost points\n")); 16947c6ae99SBarry Smith } else { 17047c6ae99SBarry Smith /* otherwise we can use the da we already have */ 17147c6ae99SBarry Smith dac = da; 17247c6ae99SBarry Smith } 17347c6ae99SBarry Smith /* create local vector for holding ghosted values used in graphics */ 1749566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dac, &xlocal)); 17547c6ae99SBarry Smith if (dac != da) { 17615229ffcSPierre Jolivet /* don't keep any public reference of this DMDA, it is only available through xlocal */ 1779566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dac)); 17847c6ae99SBarry Smith } else { 17947c6ae99SBarry Smith /* remove association between xlocal and da, because below we compose in the opposite 18047c6ae99SBarry Smith direction and if we left this connect we'd get a loop, so the objects could 18147c6ae99SBarry Smith never be destroyed */ 1829566063dSJacob Faibussowitsch PetscCall(PetscObjectRemoveReference((PetscObject)xlocal, "__PETSc_dm")); 18347c6ae99SBarry Smith } 1849566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)da, "GraphicsGhosted", (PetscObject)xlocal)); 1859566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)xlocal)); 18647c6ae99SBarry Smith } else { 187bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_NONE || by != DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) { 1889566063dSJacob Faibussowitsch PetscCall(VecGetDM(xlocal, &dac)); 189f7923d8aSBarry Smith } else { 190f7923d8aSBarry Smith dac = da; 19147c6ae99SBarry Smith } 19247c6ae99SBarry Smith } 19347c6ae99SBarry Smith 19447c6ae99SBarry Smith /* 19547c6ae99SBarry Smith Get local (ghosted) values of vector 19647c6ae99SBarry Smith */ 1979566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dac, xin, INSERT_VALUES, xlocal)); 1989566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dac, xin, INSERT_VALUES, xlocal)); 1999566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xlocal, &zctx.v)); 20047c6ae99SBarry Smith 201832b7cebSLisandro Dalcin /* 202832b7cebSLisandro Dalcin Get coordinates of nodes 203832b7cebSLisandro Dalcin */ 2049566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &xcoor)); 20547c6ae99SBarry Smith if (!xcoor) { 2069566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 0.0)); 2079566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &xcoor)); 20847c6ae99SBarry Smith } 20947c6ae99SBarry Smith 21047c6ae99SBarry Smith /* 21147c6ae99SBarry Smith Determine the min and max coordinates in plot 21247c6ae99SBarry Smith */ 2139566063dSJacob Faibussowitsch PetscCall(VecStrideMin(xcoor, 0, NULL, &zctx.xmin)); 2149566063dSJacob Faibussowitsch PetscCall(VecStrideMax(xcoor, 0, NULL, &zctx.xmax)); 2159566063dSJacob Faibussowitsch PetscCall(VecStrideMin(xcoor, 1, NULL, &zctx.ymin)); 2169566063dSJacob Faibussowitsch PetscCall(VecStrideMax(xcoor, 1, NULL, &zctx.ymax)); 2179566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_contour_axis", &zctx.showaxis, NULL)); 2187fe7c8feSLisandro Dalcin if (zctx.showaxis) { 2199371c9d4SSatish Balay coors[0] = zctx.xmin - .05 * (zctx.xmax - zctx.xmin); 2209371c9d4SSatish Balay coors[1] = zctx.ymin - .05 * (zctx.ymax - zctx.ymin); 2219371c9d4SSatish Balay coors[2] = zctx.xmax + .05 * (zctx.xmax - zctx.xmin); 2229371c9d4SSatish Balay coors[3] = zctx.ymax + .05 * (zctx.ymax - zctx.ymin); 2237fe7c8feSLisandro Dalcin } else { 2249371c9d4SSatish Balay coors[0] = zctx.xmin; 2259371c9d4SSatish Balay coors[1] = zctx.ymin; 2269371c9d4SSatish Balay coors[2] = zctx.xmax; 2279371c9d4SSatish Balay coors[3] = zctx.ymax; 2287fe7c8feSLisandro Dalcin } 2299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-draw_coordinates", coors, &ncoors, NULL)); 2309566063dSJacob 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])); 231a74ba6f7SBarry Smith 23247c6ae99SBarry Smith /* 233832b7cebSLisandro Dalcin Get local ghosted version of coordinates 23447c6ae99SBarry Smith */ 2359566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)da, "GraphicsCoordinateGhosted", (PetscObject *)&xcoorl)); 23647c6ae99SBarry Smith if (!xcoorl) { 237aa219208SBarry Smith /* create DMDA to get local version of graphics */ 2389566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, zctx.m, zctx.n, 2, 1, lx, ly, &dag)); 2399566063dSJacob Faibussowitsch PetscCall(DMSetUp(dag)); 2406aad120cSJose E. Roman PetscCall(PetscInfo(dag, "Creating auxiliary DMDA for managing graphics coordinates ghost points\n")); 2419566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dag, &xcoorl)); 2429566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)da, "GraphicsCoordinateGhosted", (PetscObject)xcoorl)); 2439566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dag)); 2449566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)xcoorl)); 2451baa6e33SBarry Smith } else PetscCall(VecGetDM(xcoorl, &dag)); 2469566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dag, xcoor, INSERT_VALUES, xcoorl)); 2479566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dag, xcoor, INSERT_VALUES, xcoorl)); 2489566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xcoorl, &zctx.xy)); 2499566063dSJacob Faibussowitsch PetscCall(DMDAGetCoordinateName(da, 0, &zctx.name0)); 2509566063dSJacob Faibussowitsch PetscCall(DMDAGetCoordinateName(da, 1, &zctx.name1)); 25147c6ae99SBarry Smith 25247c6ae99SBarry Smith /* 25347c6ae99SBarry Smith Get information about size of area each processor must do graphics for 25447c6ae99SBarry Smith */ 2559566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dac, NULL, &M, &N, NULL, NULL, NULL, NULL, &zctx.dof, NULL, &bx, &by, NULL, NULL)); 2569566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dac, NULL, NULL, NULL, &zctx.m, &zctx.n, NULL)); 2579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_contour_grid", &zctx.showgrid, NULL)); 2584e6118eeSBarry Smith 2599566063dSJacob Faibussowitsch PetscCall(DMDASelectFields(da, &ndisplayfields, &displayfields)); 2609566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 2619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_ports", &useports, NULL)); 262832b7cebSLisandro Dalcin if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE; 263832b7cebSLisandro Dalcin if (useports) { 2649566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 2659566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 2669566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 2679566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsCreate(draw, ndisplayfields, &ports)); 268e5ab1681SLisandro Dalcin } 26920d0051dSBarry Smith 270832b7cebSLisandro Dalcin /* 271832b7cebSLisandro Dalcin Loop over each field; drawing each in a different window 272832b7cebSLisandro Dalcin */ 27320d0051dSBarry Smith for (i = 0; i < ndisplayfields; i++) { 27420d0051dSBarry Smith zctx.k = displayfields[i]; 2755b399a63SLisandro Dalcin 276e5ab1681SLisandro Dalcin /* determine the min and max value in plot */ 2779566063dSJacob Faibussowitsch PetscCall(VecStrideMin(xin, zctx.k, NULL, &zctx.min)); 2789566063dSJacob Faibussowitsch PetscCall(VecStrideMax(xin, zctx.k, NULL, &zctx.max)); 27967dd0837SBarry Smith if (zctx.k < nbounds) { 280f3f0eb19SBarry Smith zctx.min = bounds[2 * zctx.k]; 281f3f0eb19SBarry Smith zctx.max = bounds[2 * zctx.k + 1]; 28267dd0837SBarry Smith } 28347c6ae99SBarry Smith if (zctx.min == zctx.max) { 28447c6ae99SBarry Smith zctx.min -= 1.e-12; 28547c6ae99SBarry Smith zctx.max += 1.e-12; 28647c6ae99SBarry Smith } 2879566063dSJacob Faibussowitsch PetscCall(PetscInfo(da, "DMDA 2d contour plot min %g max %g\n", (double)zctx.min, (double)zctx.max)); 28847c6ae99SBarry Smith 289832b7cebSLisandro Dalcin if (useports) { 2909566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports, i)); 291832b7cebSLisandro Dalcin } else { 292832b7cebSLisandro Dalcin const char *title; 2939566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, i, &draw)); 2949566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da, zctx.k, &title)); 2959566063dSJacob Faibussowitsch if (title) PetscCall(PetscDrawSetTitle(draw, title)); 296832b7cebSLisandro Dalcin } 297832b7cebSLisandro Dalcin 2989566063dSJacob Faibussowitsch PetscCall(PetscDrawGetPopup(draw, &popup)); 2999566063dSJacob Faibussowitsch PetscCall(PetscDrawScalePopup(popup, zctx.min, zctx.max)); 3009566063dSJacob Faibussowitsch PetscCall(PetscDrawSetCoordinates(draw, coors[0], coors[1], coors[2], coors[3])); 3019566063dSJacob Faibussowitsch PetscCall(PetscDrawZoom(draw, VecView_MPI_Draw_DA2d_Zoom, &zctx)); 3029566063dSJacob Faibussowitsch if (!useports) PetscCall(PetscDrawSave(draw)); 30347c6ae99SBarry Smith } 304832b7cebSLisandro Dalcin if (useports) { 3059566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 3069566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 307832b7cebSLisandro Dalcin } 308832b7cebSLisandro Dalcin 3099566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsDestroy(ports)); 3109566063dSJacob Faibussowitsch PetscCall(PetscFree(displayfields)); 3119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xcoorl, &zctx.xy)); 3129566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xlocal, &zctx.v)); 3133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 31447c6ae99SBarry Smith } 31547c6ae99SBarry Smith 3160da24e51SJuha Jäykkä #if defined(PETSC_HAVE_HDF5) 317d71ae5a4SJacob Faibussowitsch static PetscErrorCode VecGetHDF5ChunkSize(DM_DA *da, Vec xin, PetscInt dimension, PetscInt timestep, hsize_t *chunkDims) 318d71ae5a4SJacob Faibussowitsch { 319d4190030SJed Brown PetscMPIInt comm_size; 320561a82e7SJed Brown hsize_t chunk_size, target_size, dim; 321561a82e7SJed Brown hsize_t vec_size = sizeof(PetscScalar) * da->M * da->N * da->P * da->w; 322561a82e7SJed Brown hsize_t avg_local_vec_size, KiB = 1024, MiB = KiB * KiB, GiB = MiB * KiB, min_size = MiB; 323561a82e7SJed Brown hsize_t max_chunks = 64 * KiB; /* HDF5 internal limitation */ 324561a82e7SJed Brown hsize_t max_chunk_size = 4 * GiB; /* HDF5 internal limitation */ 32584179ffaSJed Brown PetscInt zslices = da->p, yslices = da->n, xslices = da->m; 3260da24e51SJuha Jäykkä 327cb5c4f94SJuha Jäykkä PetscFunctionBegin; 3289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)xin), &comm_size)); 329faa75363SBarry Smith avg_local_vec_size = (hsize_t)PetscCeilInt(vec_size, comm_size); /* we will attempt to use this as the chunk size */ 3300da24e51SJuha Jäykkä 331faa75363SBarry 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)))); 3327d310018SBarry 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 */ 3337d310018SBarry 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); 3340da24e51SJuha Jäykkä 335cb5c4f94SJuha Jäykkä /* 336cb5c4f94SJuha Jäykkä if size/rank > max_chunk_size, we need radical measures: even going down to 337cb5c4f94SJuha Jäykkä avg_local_vec_size is not enough, so we simply use chunk size of 4 GiB no matter 338cb5c4f94SJuha Jäykkä what, composed in the most efficient way possible. 339cb5c4f94SJuha Jäykkä N.B. this minimises the number of chunks, which may or may not be the optimal 340cb5c4f94SJuha Jäykkä solution. In a BG, for example, the optimal solution is probably to make # chunks = # 341cb5c4f94SJuha Jäykkä IO nodes involved, but this author has no access to a BG to figure out how to 342cb5c4f94SJuha Jäykkä reliably find the right number. And even then it may or may not be enough. 343cb5c4f94SJuha Jäykkä */ 3440da24e51SJuha Jäykkä if (avg_local_vec_size > max_chunk_size) { 345cb5c4f94SJuha Jäykkä /* check if we can just split local z-axis: is that enough? */ 346faa75363SBarry Smith zslices = PetscCeilInt(vec_size, da->p * max_chunk_size) * zslices; 3470da24e51SJuha Jäykkä if (zslices > da->P) { 348cb5c4f94SJuha Jäykkä /* lattice is too large in xy-directions, splitting z only is not enough */ 3490da24e51SJuha Jäykkä zslices = da->P; 350faa75363SBarry Smith yslices = PetscCeilInt(vec_size, zslices * da->n * max_chunk_size) * yslices; 3510da24e51SJuha Jäykkä if (yslices > da->N) { 352cb5c4f94SJuha Jäykkä /* lattice is too large in x-direction, splitting along z, y is not enough */ 3530da24e51SJuha Jäykkä yslices = da->N; 354faa75363SBarry Smith xslices = PetscCeilInt(vec_size, zslices * yslices * da->m * max_chunk_size) * xslices; 3550da24e51SJuha Jäykkä } 3560da24e51SJuha Jäykkä } 3570da24e51SJuha Jäykkä dim = 0; 358ad540459SPierre Jolivet if (timestep >= 0) ++dim; 359cb5c4f94SJuha Jäykkä /* prefer to split z-axis, even down to planar slices */ 360c73cfb54SMatthew G. Knepley if (dimension == 3) { 361cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t)da->P / zslices; 362cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t)da->N / yslices; 363cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t)da->M / xslices; 3640da24e51SJuha Jäykkä } else { 365cb5c4f94SJuha Jäykkä /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */ 366cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t)da->N / yslices; 367cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t)da->M / xslices; 3680da24e51SJuha Jäykkä } 3690da24e51SJuha 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); 3700da24e51SJuha Jäykkä } else { 3710da24e51SJuha Jäykkä if (target_size < chunk_size) { 372cb5c4f94SJuha Jäykkä /* only change the defaults if target_size < chunk_size */ 3730da24e51SJuha Jäykkä dim = 0; 374ad540459SPierre Jolivet if (timestep >= 0) ++dim; 375cb5c4f94SJuha Jäykkä /* prefer to split z-axis, even down to planar slices */ 376c73cfb54SMatthew G. Knepley if (dimension == 3) { 377cb5c4f94SJuha Jäykkä /* try splitting the z-axis to core-size bits, i.e. divide chunk size by # comm_size in z-direction */ 3780da24e51SJuha Jäykkä if (target_size >= chunk_size / da->p) { 379cb5c4f94SJuha Jäykkä /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */ 380faa75363SBarry Smith chunkDims[dim] = (hsize_t)PetscCeilInt(da->P, da->p); 3810da24e51SJuha Jäykkä } else { 382cb5c4f94SJuha Jäykkä /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be 383cb5c4f94SJuha Jäykkä radical and let everyone write all they've got */ 384faa75363SBarry Smith chunkDims[dim++] = (hsize_t)PetscCeilInt(da->P, da->p); 385faa75363SBarry Smith chunkDims[dim++] = (hsize_t)PetscCeilInt(da->N, da->n); 386faa75363SBarry Smith chunkDims[dim++] = (hsize_t)PetscCeilInt(da->M, da->m); 3870da24e51SJuha Jäykkä } 3880da24e51SJuha Jäykkä } else { 389cb5c4f94SJuha Jäykkä /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */ 3900da24e51SJuha Jäykkä if (target_size >= chunk_size / da->n) { 391cb5c4f94SJuha Jäykkä /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */ 392faa75363SBarry Smith chunkDims[dim] = (hsize_t)PetscCeilInt(da->N, da->n); 3930da24e51SJuha Jäykkä } else { 394cb5c4f94SJuha Jäykkä /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be 395cb5c4f94SJuha Jäykkä radical and let everyone write all they've got */ 396faa75363SBarry Smith chunkDims[dim++] = (hsize_t)PetscCeilInt(da->N, da->n); 397faa75363SBarry Smith chunkDims[dim++] = (hsize_t)PetscCeilInt(da->M, da->m); 3980da24e51SJuha Jäykkä } 3990da24e51SJuha Jäykkä } 4000da24e51SJuha 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); 4010da24e51SJuha Jäykkä } else { 402cb5c4f94SJuha Jäykkä /* precomputed chunks are fine, we don't need to do anything */ 4030da24e51SJuha Jäykkä } 4040da24e51SJuha Jäykkä } 4053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4060da24e51SJuha Jäykkä } 4070da24e51SJuha Jäykkä #endif 40847c6ae99SBarry Smith 40947c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 41066976f2fSJacob Faibussowitsch static PetscErrorCode VecView_MPI_HDF5_DA(Vec xin, PetscViewer viewer) 411d71ae5a4SJacob Faibussowitsch { 4122d95be3dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 4139b2a5a86SJed Brown DM dm; 4149b2a5a86SJed Brown DM_DA *da; 41547c6ae99SBarry Smith hid_t filespace; /* file dataspace identifier */ 4168e2ae6d7SMichael Kraus hid_t chunkspace; /* chunk dataset property identifier */ 41747c6ae99SBarry Smith hid_t dset_id; /* dataset identifier */ 41847c6ae99SBarry Smith hid_t memspace; /* memory dataspace identifier */ 41947c6ae99SBarry Smith hid_t file_id; 42015214e8eSMatthew G Knepley hid_t group; 4219a0502c6SHåkon Strandenes hid_t memscalartype; /* scalar type for mem (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */ 4229a0502c6SHåkon Strandenes hid_t filescalartype; /* scalar type for file (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */ 423d9a4edebSJed Brown hsize_t dim; 4240da24e51SJuha 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 */ 425d7dd068bSVaclav Hapla PetscBool timestepping = PETSC_FALSE, dim2, spoutput; 426*1690c2aeSBarry Smith PetscInt timestep = PETSC_INT_MIN, dimension; 4275edff71fSBarry Smith const PetscScalar *x; 4288e2ae6d7SMichael Kraus const char *vecname; 42947c6ae99SBarry Smith 43047c6ae99SBarry Smith PetscFunctionBegin; 4313014b61aSVaclav Hapla PetscCall(PetscViewerHDF5OpenGroup(viewer, NULL, &file_id, &group)); 4329566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5IsTimestepping(viewer, ×tepping)); 43348a46eb9SPierre Jolivet if (timestepping) PetscCall(PetscViewerHDF5GetTimestep(viewer, ×tep)); 4349566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetBaseDimension2(viewer, &dim2)); 4359566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetSPOutput(viewer, &spoutput)); 43615214e8eSMatthew G Knepley 4379566063dSJacob Faibussowitsch PetscCall(VecGetDM(xin, &dm)); 4387a8be351SBarry Smith PetscCheck(dm, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA"); 4399b2a5a86SJed Brown da = (DM_DA *)dm->data; 4409566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dimension)); 44147c6ae99SBarry Smith 4428e2ae6d7SMichael Kraus /* Create the dataspace for the dataset. 4438e2ae6d7SMichael Kraus * 4448e2ae6d7SMichael Kraus * dims - holds the current dimensions of the dataset 4458e2ae6d7SMichael Kraus * 4468e2ae6d7SMichael Kraus * maxDims - holds the maximum dimensions of the dataset (unlimited 4478e2ae6d7SMichael Kraus * for the number of time steps with the current dimensions for the 4488e2ae6d7SMichael Kraus * other dimensions; so only additional time steps can be added). 4498e2ae6d7SMichael Kraus * 4508e2ae6d7SMichael Kraus * chunkDims - holds the size of a single time step (required to 4518e2ae6d7SMichael Kraus * permit extending dataset). 4528e2ae6d7SMichael Kraus */ 4538e2ae6d7SMichael Kraus dim = 0; 4548e2ae6d7SMichael Kraus if (timestep >= 0) { 4558e2ae6d7SMichael Kraus dims[dim] = timestep + 1; 4568e2ae6d7SMichael Kraus maxDims[dim] = H5S_UNLIMITED; 4578e2ae6d7SMichael Kraus chunkDims[dim] = 1; 4588e2ae6d7SMichael Kraus ++dim; 4598e2ae6d7SMichael Kraus } 460c73cfb54SMatthew G. Knepley if (dimension == 3) { 4619566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(da->P, dims + dim)); 4628e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4638e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4648e2ae6d7SMichael Kraus ++dim; 4658e2ae6d7SMichael Kraus } 466c73cfb54SMatthew G. Knepley if (dimension > 1) { 4679566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(da->N, dims + dim)); 4688e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4698e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4708e2ae6d7SMichael Kraus ++dim; 4718e2ae6d7SMichael Kraus } 4729566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(da->M, dims + dim)); 4738e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4748e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4758e2ae6d7SMichael Kraus ++dim; 47682ea9b62SBarry Smith if (da->w > 1 || dim2) { 4779566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(da->w, dims + dim)); 4788e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4798e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4808e2ae6d7SMichael Kraus ++dim; 4818e2ae6d7SMichael Kraus } 48247c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 4838e2ae6d7SMichael Kraus dims[dim] = 2; 4848e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4858e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4868e2ae6d7SMichael Kraus ++dim; 48747c6ae99SBarry Smith #endif 4880da24e51SJuha Jäykkä 4899566063dSJacob Faibussowitsch PetscCall(VecGetHDF5ChunkSize(da, xin, dimension, timestep, chunkDims)); 4900da24e51SJuha Jäykkä 4916497c311SBarry Smith PetscCallHDF5Return(filespace, H5Screate_simple, ((int)dim, dims, maxDims)); 49247c6ae99SBarry Smith 49315214e8eSMatthew G Knepley #if defined(PETSC_USE_REAL_SINGLE) 4949a0502c6SHåkon Strandenes memscalartype = H5T_NATIVE_FLOAT; 4959a0502c6SHåkon Strandenes filescalartype = H5T_NATIVE_FLOAT; 49615214e8eSMatthew G Knepley #elif defined(PETSC_USE_REAL___FLOAT128) 49715214e8eSMatthew G Knepley #error "HDF5 output with 128 bit floats not supported." 498570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16) 499570b7f6dSBarry Smith #error "HDF5 output with 16 bit floats not supported." 50015214e8eSMatthew G Knepley #else 5019a0502c6SHåkon Strandenes memscalartype = H5T_NATIVE_DOUBLE; 5029a0502c6SHåkon Strandenes if (spoutput == PETSC_TRUE) filescalartype = H5T_NATIVE_FLOAT; 5039a0502c6SHåkon Strandenes else filescalartype = H5T_NATIVE_DOUBLE; 50415214e8eSMatthew G Knepley #endif 50515214e8eSMatthew G Knepley 50647c6ae99SBarry Smith /* Create the dataset with default properties and close filespace */ 5079566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)xin, &vecname)); 50815214e8eSMatthew G Knepley if (!H5Lexists(group, vecname, H5P_DEFAULT)) { 5098e2ae6d7SMichael Kraus /* Create chunk */ 510792fecdfSBarry Smith PetscCallHDF5Return(chunkspace, H5Pcreate, (H5P_DATASET_CREATE)); 5116497c311SBarry Smith PetscCallHDF5(H5Pset_chunk, (chunkspace, (int)dim, chunkDims)); 5128e2ae6d7SMichael Kraus 513792fecdfSBarry Smith PetscCallHDF5Return(dset_id, H5Dcreate2, (group, vecname, filescalartype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT)); 51415214e8eSMatthew G Knepley } else { 515792fecdfSBarry Smith PetscCallHDF5Return(dset_id, H5Dopen2, (group, vecname, H5P_DEFAULT)); 516792fecdfSBarry Smith PetscCallHDF5(H5Dset_extent, (dset_id, dims)); 51715214e8eSMatthew G Knepley } 518792fecdfSBarry Smith PetscCallHDF5(H5Sclose, (filespace)); 51947c6ae99SBarry Smith 52047c6ae99SBarry Smith /* Each process defines a dataset and writes it to the hyperslab in the file */ 5218e2ae6d7SMichael Kraus dim = 0; 5228e2ae6d7SMichael Kraus if (timestep >= 0) { 5238e2ae6d7SMichael Kraus offset[dim] = timestep; 5248e2ae6d7SMichael Kraus ++dim; 5258e2ae6d7SMichael Kraus } 5269566063dSJacob Faibussowitsch if (dimension == 3) PetscCall(PetscHDF5IntCast(da->zs, offset + dim++)); 5279566063dSJacob Faibussowitsch if (dimension > 1) PetscCall(PetscHDF5IntCast(da->ys, offset + dim++)); 5289566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(da->xs / da->w, offset + dim++)); 52982ea9b62SBarry Smith if (da->w > 1 || dim2) offset[dim++] = 0; 53047c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 5318e2ae6d7SMichael Kraus offset[dim++] = 0; 53247c6ae99SBarry Smith #endif 5338e2ae6d7SMichael Kraus dim = 0; 5348e2ae6d7SMichael Kraus if (timestep >= 0) { 5358e2ae6d7SMichael Kraus count[dim] = 1; 5368e2ae6d7SMichael Kraus ++dim; 5378e2ae6d7SMichael Kraus } 5389566063dSJacob Faibussowitsch if (dimension == 3) PetscCall(PetscHDF5IntCast(da->ze - da->zs, count + dim++)); 5399566063dSJacob Faibussowitsch if (dimension > 1) PetscCall(PetscHDF5IntCast(da->ye - da->ys, count + dim++)); 5409566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast((da->xe - da->xs) / da->w, count + dim++)); 5419566063dSJacob Faibussowitsch if (da->w > 1 || dim2) PetscCall(PetscHDF5IntCast(da->w, count + dim++)); 54247c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 5438e2ae6d7SMichael Kraus count[dim++] = 2; 54447c6ae99SBarry Smith #endif 5456497c311SBarry Smith PetscCallHDF5Return(memspace, H5Screate_simple, ((int)dim, count, NULL)); 546792fecdfSBarry Smith PetscCallHDF5Return(filespace, H5Dget_space, (dset_id)); 547792fecdfSBarry Smith PetscCallHDF5(H5Sselect_hyperslab, (filespace, H5S_SELECT_SET, offset, NULL, count, NULL)); 54847c6ae99SBarry Smith 5499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xin, &x)); 550792fecdfSBarry Smith PetscCallHDF5(H5Dwrite, (dset_id, memscalartype, memspace, filespace, hdf5->dxpl_id, x)); 551792fecdfSBarry Smith PetscCallHDF5(H5Fflush, (file_id, H5F_SCOPE_GLOBAL)); 5529566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xin, &x)); 55347c6ae99SBarry Smith 554e0552f36Ssajid__ali #if defined(PETSC_USE_COMPLEX) 555e0552f36Ssajid__ali { 556e0552f36Ssajid__ali PetscBool tru = PETSC_TRUE; 5579566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)xin, "complex", PETSC_BOOL, &tru)); 558e0552f36Ssajid__ali } 559e0552f36Ssajid__ali #endif 56048a46eb9SPierre Jolivet if (timestepping) PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)xin, "timestepping", PETSC_BOOL, ×tepping)); 561e0552f36Ssajid__ali 56247c6ae99SBarry Smith /* Close/release resources */ 563ad540459SPierre Jolivet if (group != file_id) PetscCallHDF5(H5Gclose, (group)); 564792fecdfSBarry Smith PetscCallHDF5(H5Sclose, (filespace)); 565792fecdfSBarry Smith PetscCallHDF5(H5Sclose, (memspace)); 566792fecdfSBarry Smith PetscCallHDF5(H5Dclose, (dset_id)); 5679566063dSJacob Faibussowitsch PetscCall(PetscInfo(xin, "Wrote Vec object with name %s\n", vecname)); 5683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56947c6ae99SBarry Smith } 57047c6ae99SBarry Smith #endif 57147c6ae99SBarry Smith 57209573ac7SBarry Smith extern PetscErrorCode VecView_MPI_Draw_DA1d(Vec, PetscViewer); 57347c6ae99SBarry Smith 57447c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 575d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAArrayMPIIO(DM da, PetscViewer viewer, Vec xin, PetscBool write) 576d71ae5a4SJacob Faibussowitsch { 57747c6ae99SBarry Smith MPI_File mfdes; 57847c6ae99SBarry Smith PetscMPIInt gsizes[4], lsizes[4], lstarts[4], asiz, dof; 57947c6ae99SBarry Smith MPI_Datatype view; 58047c6ae99SBarry Smith const PetscScalar *array; 58147c6ae99SBarry Smith MPI_Offset off; 58247c6ae99SBarry Smith MPI_Aint ub, ul; 58347c6ae99SBarry Smith PetscInt type, rows, vecrows, tr[2]; 58447c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 5850c169764Sdmay PetscBool skipheader; 58647c6ae99SBarry Smith 58747c6ae99SBarry Smith PetscFunctionBegin; 5889566063dSJacob Faibussowitsch PetscCall(VecGetSize(xin, &vecrows)); 5899566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipheader)); 59047c6ae99SBarry Smith if (!write) { 59147c6ae99SBarry Smith /* Read vector header. */ 5920c169764Sdmay if (!skipheader) { 5939566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, tr, 2, NULL, PETSC_INT)); 59447c6ae99SBarry Smith type = tr[0]; 59547c6ae99SBarry Smith rows = tr[1]; 59608401ef6SPierre Jolivet PetscCheck(type == VEC_FILE_CLASSID, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONG, "Not vector next in file"); 59708401ef6SPierre Jolivet PetscCheck(rows == vecrows, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Vector in file not same size as DMDA vector"); 5980c169764Sdmay } 59947c6ae99SBarry Smith } else { 60047c6ae99SBarry Smith tr[0] = VEC_FILE_CLASSID; 60147c6ae99SBarry Smith tr[1] = vecrows; 60248a46eb9SPierre Jolivet if (!skipheader) PetscCall(PetscViewerBinaryWrite(viewer, tr, 2, PETSC_INT)); 6030c169764Sdmay } 60447c6ae99SBarry Smith 6059566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->w, &dof)); 6064dc2109aSBarry Smith gsizes[0] = dof; 6079566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->M, gsizes + 1)); 6089566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->N, gsizes + 2)); 6099566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->P, gsizes + 3)); 6104dc2109aSBarry Smith lsizes[0] = dof; 6119566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast((dd->xe - dd->xs) / dof, lsizes + 1)); 6129566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->ye - dd->ys, lsizes + 2)); 6139566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->ze - dd->zs, lsizes + 3)); 6144dc2109aSBarry Smith lstarts[0] = 0; 6159566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->xs / dof, lstarts + 1)); 6169566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->ys, lstarts + 2)); 6179566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(dd->zs, lstarts + 3)); 6186497c311SBarry Smith PetscCallMPI(MPI_Type_create_subarray((PetscMPIInt)(da->dim + 1), gsizes, lsizes, lstarts, MPI_ORDER_FORTRAN, MPIU_SCALAR, &view)); 6199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&view)); 62047c6ae99SBarry Smith 6219566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetMPIIODescriptor(viewer, &mfdes)); 6229566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetMPIIOOffset(viewer, &off)); 6239566063dSJacob Faibussowitsch PetscCallMPI(MPI_File_set_view(mfdes, off, MPIU_SCALAR, view, (char *)"native", MPI_INFO_NULL)); 6249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xin, &array)); 62547c6ae99SBarry Smith asiz = lsizes[1] * (lsizes[2] > 0 ? lsizes[2] : 1) * (lsizes[3] > 0 ? lsizes[3] : 1) * dof; 6261baa6e33SBarry Smith if (write) PetscCall(MPIU_File_write_all(mfdes, (PetscScalar *)array, asiz, MPIU_SCALAR, MPI_STATUS_IGNORE)); 6271baa6e33SBarry Smith else PetscCall(MPIU_File_read_all(mfdes, (PetscScalar *)array, asiz, MPIU_SCALAR, MPI_STATUS_IGNORE)); 6289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_get_extent(view, &ul, &ub)); 6299566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryAddMPIIOOffset(viewer, ub)); 6309566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xin, &array)); 6319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&view)); 6323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 63347c6ae99SBarry Smith } 63447c6ae99SBarry Smith #endif 63547c6ae99SBarry Smith 636d71ae5a4SJacob Faibussowitsch PetscErrorCode VecView_MPI_DA(Vec xin, PetscViewer viewer) 637d71ae5a4SJacob Faibussowitsch { 6389a42bb27SBarry Smith DM da; 63947c6ae99SBarry Smith PetscInt dim; 64047c6ae99SBarry Smith Vec natural; 6418135c375SStefano Zampini PetscBool isdraw, isvtk, isglvis; 64247c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 64347c6ae99SBarry Smith PetscBool ishdf5; 64447c6ae99SBarry Smith #endif 6453f3fd955SJed Brown const char *prefix, *name; 646a261c58fSBarry Smith PetscViewerFormat format; 64747c6ae99SBarry Smith 64847c6ae99SBarry Smith PetscFunctionBegin; 6499566063dSJacob Faibussowitsch PetscCall(VecGetDM(xin, &da)); 6507a8be351SBarry Smith PetscCheck(da, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA"); 6519566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 6529566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk)); 65347c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 6549566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5)); 65547c6ae99SBarry Smith #endif 6569566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis)); 65747c6ae99SBarry Smith if (isdraw) { 6589566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 65947c6ae99SBarry Smith if (dim == 1) { 6609566063dSJacob Faibussowitsch PetscCall(VecView_MPI_Draw_DA1d(xin, viewer)); 66147c6ae99SBarry Smith } else if (dim == 2) { 6629566063dSJacob Faibussowitsch PetscCall(VecView_MPI_Draw_DA2d(xin, viewer)); 66363a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Cannot graphically view vector associated with this dimensional DMDA %" PetscInt_FMT, dim); 66401d7c5c3SPatrick Sanan } else if (isvtk) { /* Duplicate the Vec */ 6654061b8bfSJed Brown Vec Y; 6669566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xin, &Y)); 667b51b94faSJed Brown if (((PetscObject)xin)->name) { 668b51b94faSJed Brown /* If xin was named, copy the name over to Y. The duplicate names are safe because nobody else will ever see Y. */ 6699566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)Y, ((PetscObject)xin)->name)); 670b51b94faSJed Brown } 6719566063dSJacob Faibussowitsch PetscCall(VecCopy(xin, Y)); 672a8f87f1dSPatrick Sanan { 673a8f87f1dSPatrick Sanan PetscObject dmvtk; 674a8f87f1dSPatrick Sanan PetscBool compatible, compatibleSet; 6759566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKGetDM(viewer, &dmvtk)); 676a8f87f1dSPatrick Sanan if (dmvtk) { 677064a246eSJacob Faibussowitsch PetscValidHeaderSpecific((DM)dmvtk, DM_CLASSID, 2); 6789566063dSJacob Faibussowitsch PetscCall(DMGetCompatibility(da, (DM)dmvtk, &compatible, &compatibleSet)); 6797a8be351SBarry 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."); 680a8f87f1dSPatrick Sanan } 6819566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKAddField(viewer, (PetscObject)da, DMDAVTKWriteAll, PETSC_DEFAULT, PETSC_VTK_POINT_FIELD, PETSC_FALSE, (PetscObject)Y)); 682a8f87f1dSPatrick Sanan } 68347c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 68447c6ae99SBarry Smith } else if (ishdf5) { 6859566063dSJacob Faibussowitsch PetscCall(VecView_MPI_HDF5_DA(xin, viewer)); 68647c6ae99SBarry Smith #endif 6878135c375SStefano Zampini } else if (isglvis) { 6889566063dSJacob Faibussowitsch PetscCall(VecView_GLVis(xin, viewer)); 68947c6ae99SBarry Smith } else { 69047c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 69147c6ae99SBarry Smith PetscBool isbinary, isMPIIO; 69247c6ae99SBarry Smith 6939566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 69447c6ae99SBarry Smith if (isbinary) { 6959566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &isMPIIO)); 69647c6ae99SBarry Smith if (isMPIIO) { 6979566063dSJacob Faibussowitsch PetscCall(DMDAArrayMPIIO(da, viewer, xin, PETSC_TRUE)); 6983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 69947c6ae99SBarry Smith } 70047c6ae99SBarry Smith } 70147c6ae99SBarry Smith #endif 70247c6ae99SBarry Smith 70347c6ae99SBarry Smith /* call viewer on natural ordering */ 7049566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)xin, &prefix)); 7059566063dSJacob Faibussowitsch PetscCall(DMDACreateNaturalVector(da, &natural)); 7069566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)natural, prefix)); 7079566063dSJacob Faibussowitsch PetscCall(DMDAGlobalToNaturalBegin(da, xin, INSERT_VALUES, natural)); 7089566063dSJacob Faibussowitsch PetscCall(DMDAGlobalToNaturalEnd(da, xin, INSERT_VALUES, natural)); 7099566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)xin, &name)); 7109566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)natural, name)); 711a261c58fSBarry Smith 7129566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 713a261c58fSBarry Smith if (format == PETSC_VIEWER_BINARY_MATLAB) { 714a261c58fSBarry Smith /* temporarily remove viewer format so it won't trigger in the VecView() */ 7159566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_DEFAULT)); 716a261c58fSBarry Smith } 717a261c58fSBarry Smith 71823f88b65SBarry Smith ((PetscObject)natural)->donotPetscObjectPrintClassNamePrefixType = PETSC_TRUE; 7199566063dSJacob Faibussowitsch PetscCall(VecView(natural, viewer)); 72023f88b65SBarry Smith ((PetscObject)natural)->donotPetscObjectPrintClassNamePrefixType = PETSC_FALSE; 721a261c58fSBarry Smith 722a261c58fSBarry Smith if (format == PETSC_VIEWER_BINARY_MATLAB) { 723a261c58fSBarry Smith MPI_Comm comm; 724a261c58fSBarry Smith FILE *info; 725a261c58fSBarry Smith const char *fieldname; 726da88d4d4SJed Brown char fieldbuf[256]; 727a261c58fSBarry Smith PetscInt dim, ni, nj, nk, pi, pj, pk, dof, n; 728a261c58fSBarry Smith 729a261c58fSBarry Smith /* set the viewer format back into the viewer */ 7309566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 7319566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm)); 7329566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetInfoPointer(viewer, &info)); 7339566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, &ni, &nj, &nk, &pi, &pj, &pk, &dof, NULL, NULL, NULL, NULL, NULL)); 7349566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n")); 7359566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, info, "#$$ tmp = PetscBinaryRead(fd); \n")); 73648a46eb9SPierre Jolivet if (dim == 1) PetscCall(PetscFPrintf(comm, info, "#$$ tmp = reshape(tmp,%" PetscInt_FMT ",%" PetscInt_FMT ");\n", dof, ni)); 73748a46eb9SPierre Jolivet if (dim == 2) PetscCall(PetscFPrintf(comm, info, "#$$ tmp = reshape(tmp,%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ");\n", dof, ni, nj)); 73848a46eb9SPierre Jolivet if (dim == 3) PetscCall(PetscFPrintf(comm, info, "#$$ tmp = reshape(tmp,%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ");\n", dof, ni, nj, nk)); 739a261c58fSBarry Smith 740a261c58fSBarry Smith for (n = 0; n < dof; n++) { 7419566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da, n, &fieldname)); 742da88d4d4SJed Brown if (!fieldname || !fieldname[0]) { 74363a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(fieldbuf, sizeof fieldbuf, "field%" PetscInt_FMT, n)); 744da88d4d4SJed Brown fieldname = fieldbuf; 745a261c58fSBarry Smith } 74648a46eb9SPierre Jolivet if (dim == 1) PetscCall(PetscFPrintf(comm, info, "#$$ Set.%s.%s = squeeze(tmp(%" PetscInt_FMT ",:))';\n", name, fieldname, n + 1)); 74748a46eb9SPierre Jolivet if (dim == 2) PetscCall(PetscFPrintf(comm, info, "#$$ Set.%s.%s = squeeze(tmp(%" PetscInt_FMT ",:,:))';\n", name, fieldname, n + 1)); 74848a46eb9SPierre Jolivet if (dim == 3) PetscCall(PetscFPrintf(comm, info, "#$$ Set.%s.%s = permute(squeeze(tmp(%" PetscInt_FMT ",:,:,:)),[2 1 3]);\n", name, fieldname, n + 1)); 749a261c58fSBarry Smith } 7509566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, info, "#$$ clear tmp; \n")); 7519566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n")); 752a261c58fSBarry Smith } 753a261c58fSBarry Smith 7549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&natural)); 75547c6ae99SBarry Smith } 7563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 75747c6ae99SBarry Smith } 75847c6ae99SBarry Smith 75947c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 76066976f2fSJacob Faibussowitsch static PetscErrorCode VecLoad_HDF5_DA(Vec xin, PetscViewer viewer) 761d71ae5a4SJacob Faibussowitsch { 7622d95be3dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 7639a42bb27SBarry Smith DM da; 76415b861d2SVaclav Hapla int dim, rdim; 765ec7ae49cSHåkon Strandenes hsize_t dims[6] = {0}, count[6] = {0}, offset[6] = {0}; 766d7dd068bSVaclav Hapla PetscBool dim2 = PETSC_FALSE, timestepping = PETSC_FALSE; 767*1690c2aeSBarry Smith PetscInt dimension, timestep = PETSC_INT_MIN, dofInd; 76847c6ae99SBarry Smith PetscScalar *x; 76947c6ae99SBarry Smith const char *vecname; 77047c6ae99SBarry Smith hid_t filespace; /* file dataspace identifier */ 77147c6ae99SBarry Smith hid_t dset_id; /* dataset identifier */ 77247c6ae99SBarry Smith hid_t memspace; /* memory dataspace identifier */ 773bfd264e7SBarry Smith hid_t file_id, group; 7747bcbaff4SBarry Smith hid_t scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */ 7759c7c4993SBarry Smith DM_DA *dd; 77647c6ae99SBarry Smith 77747c6ae99SBarry Smith PetscFunctionBegin; 7787bcbaff4SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 7797bcbaff4SBarry Smith scalartype = H5T_NATIVE_FLOAT; 7807bcbaff4SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128) 7817bcbaff4SBarry Smith #error "HDF5 output with 128 bit floats not supported." 782570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16) 783570b7f6dSBarry Smith #error "HDF5 output with 16 bit floats not supported." 7847bcbaff4SBarry Smith #else 7857bcbaff4SBarry Smith scalartype = H5T_NATIVE_DOUBLE; 7867bcbaff4SBarry Smith #endif 7877bcbaff4SBarry Smith 7883014b61aSVaclav Hapla PetscCall(PetscViewerHDF5OpenGroup(viewer, NULL, &file_id, &group)); 7899566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)xin, &vecname)); 7909566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckTimestepping_Internal(viewer, vecname)); 7919566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5IsTimestepping(viewer, ×tepping)); 79248a46eb9SPierre Jolivet if (timestepping) PetscCall(PetscViewerHDF5GetTimestep(viewer, ×tep)); 7939566063dSJacob Faibussowitsch PetscCall(VecGetDM(xin, &da)); 7949c7c4993SBarry Smith dd = (DM_DA *)da->data; 7959566063dSJacob Faibussowitsch PetscCall(DMGetDimension(da, &dimension)); 79647c6ae99SBarry Smith 797ec7ae49cSHåkon Strandenes /* Open dataset */ 798792fecdfSBarry Smith PetscCallHDF5Return(dset_id, H5Dopen2, (group, vecname, H5P_DEFAULT)); 79947c6ae99SBarry Smith 800ec7ae49cSHåkon Strandenes /* Retrieve the dataspace for the dataset */ 801792fecdfSBarry Smith PetscCallHDF5Return(filespace, H5Dget_space, (dset_id)); 802792fecdfSBarry Smith PetscCallHDF5Return(rdim, H5Sget_simple_extent_dims, (filespace, dims, NULL)); 803ec7ae49cSHåkon Strandenes 804ec7ae49cSHåkon Strandenes /* Expected dimension for holding the dof's */ 80547c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 806ec7ae49cSHåkon Strandenes dofInd = rdim - 2; 807ec7ae49cSHåkon Strandenes #else 808ec7ae49cSHåkon Strandenes dofInd = rdim - 1; 80947c6ae99SBarry Smith #endif 810ec7ae49cSHåkon Strandenes 811ec7ae49cSHåkon Strandenes /* The expected number of dimensions, assuming basedimension2 = false */ 812ec7ae49cSHåkon Strandenes dim = dimension; 813ec7ae49cSHåkon Strandenes if (dd->w > 1) ++dim; 814ec7ae49cSHåkon Strandenes if (timestep >= 0) ++dim; 81547c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 816ec7ae49cSHåkon Strandenes ++dim; 81747c6ae99SBarry Smith #endif 818ec7ae49cSHåkon Strandenes 819ec7ae49cSHåkon Strandenes /* In this case the input dataset have one extra, unexpected dimension. */ 820ec7ae49cSHåkon Strandenes if (rdim == dim + 1) { 821ec7ae49cSHåkon Strandenes /* In this case the block size unity */ 822ec7ae49cSHåkon Strandenes if (dd->w == 1 && dims[dofInd] == 1) dim2 = PETSC_TRUE; 823ec7ae49cSHåkon Strandenes 824ec7ae49cSHåkon Strandenes /* Special error message for the case where dof does not match the input file */ 82563a3b9bcSJacob 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); 826ec7ae49cSHåkon Strandenes 827ec7ae49cSHåkon Strandenes /* Other cases where rdim != dim cannot be handled currently */ 82863a3b9bcSJacob 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); 829ec7ae49cSHåkon Strandenes 830ec7ae49cSHåkon Strandenes /* Set up the hyperslab size */ 831ec7ae49cSHåkon Strandenes dim = 0; 832ec7ae49cSHåkon Strandenes if (timestep >= 0) { 833ec7ae49cSHåkon Strandenes offset[dim] = timestep; 834ec7ae49cSHåkon Strandenes count[dim] = 1; 835ec7ae49cSHåkon Strandenes ++dim; 836ec7ae49cSHåkon Strandenes } 837ec7ae49cSHåkon Strandenes if (dimension == 3) { 8389566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(dd->zs, offset + dim)); 8399566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(dd->ze - dd->zs, count + dim)); 840ec7ae49cSHåkon Strandenes ++dim; 841ec7ae49cSHåkon Strandenes } 842ec7ae49cSHåkon Strandenes if (dimension > 1) { 8439566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(dd->ys, offset + dim)); 8449566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(dd->ye - dd->ys, count + dim)); 845ec7ae49cSHåkon Strandenes ++dim; 846ec7ae49cSHåkon Strandenes } 8479566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(dd->xs / dd->w, offset + dim)); 8489566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast((dd->xe - dd->xs) / dd->w, count + dim)); 849ec7ae49cSHåkon Strandenes ++dim; 850ec7ae49cSHåkon Strandenes if (dd->w > 1 || dim2) { 851ec7ae49cSHåkon Strandenes offset[dim] = 0; 8529566063dSJacob Faibussowitsch PetscCall(PetscHDF5IntCast(dd->w, count + dim)); 853ec7ae49cSHåkon Strandenes ++dim; 854ec7ae49cSHåkon Strandenes } 855ec7ae49cSHåkon Strandenes #if defined(PETSC_USE_COMPLEX) 856ec7ae49cSHåkon Strandenes offset[dim] = 0; 857ec7ae49cSHåkon Strandenes count[dim] = 2; 858ec7ae49cSHåkon Strandenes ++dim; 859ec7ae49cSHåkon Strandenes #endif 860ec7ae49cSHåkon Strandenes 861ec7ae49cSHåkon Strandenes /* Create the memory and filespace */ 862792fecdfSBarry Smith PetscCallHDF5Return(memspace, H5Screate_simple, (dim, count, NULL)); 863792fecdfSBarry Smith PetscCallHDF5(H5Sselect_hyperslab, (filespace, H5S_SELECT_SET, offset, NULL, count, NULL)); 86447c6ae99SBarry Smith 8659566063dSJacob Faibussowitsch PetscCall(VecGetArray(xin, &x)); 866792fecdfSBarry Smith PetscCallHDF5(H5Dread, (dset_id, scalartype, memspace, filespace, hdf5->dxpl_id, x)); 8679566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xin, &x)); 86847c6ae99SBarry Smith 86947c6ae99SBarry Smith /* Close/release resources */ 870ad540459SPierre Jolivet if (group != file_id) PetscCallHDF5(H5Gclose, (group)); 871792fecdfSBarry Smith PetscCallHDF5(H5Sclose, (filespace)); 872792fecdfSBarry Smith PetscCallHDF5(H5Sclose, (memspace)); 873792fecdfSBarry Smith PetscCallHDF5(H5Dclose, (dset_id)); 8743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 87547c6ae99SBarry Smith } 87647c6ae99SBarry Smith #endif 87747c6ae99SBarry Smith 87866976f2fSJacob Faibussowitsch static PetscErrorCode VecLoad_Binary_DA(Vec xin, PetscViewer viewer) 879d71ae5a4SJacob Faibussowitsch { 8809a42bb27SBarry Smith DM da; 88147c6ae99SBarry Smith Vec natural; 88247c6ae99SBarry Smith const char *prefix; 88347c6ae99SBarry Smith PetscInt bs; 88447c6ae99SBarry Smith PetscBool flag; 88547c6ae99SBarry Smith DM_DA *dd; 88647c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 88747c6ae99SBarry Smith PetscBool isMPIIO; 88847c6ae99SBarry Smith #endif 88947c6ae99SBarry Smith 89047c6ae99SBarry Smith PetscFunctionBegin; 8919566063dSJacob Faibussowitsch PetscCall(VecGetDM(xin, &da)); 89247c6ae99SBarry Smith dd = (DM_DA *)da->data; 89347c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 8949566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &isMPIIO)); 89547c6ae99SBarry Smith if (isMPIIO) { 8969566063dSJacob Faibussowitsch PetscCall(DMDAArrayMPIIO(da, viewer, xin, PETSC_FALSE)); 8973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 89847c6ae99SBarry Smith } 89947c6ae99SBarry Smith #endif 90047c6ae99SBarry Smith 9019566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)xin, &prefix)); 9029566063dSJacob Faibussowitsch PetscCall(DMDACreateNaturalVector(da, &natural)); 9039566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)natural, ((PetscObject)xin)->name)); 9049566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)natural, prefix)); 9059566063dSJacob Faibussowitsch PetscCall(VecLoad(natural, viewer)); 9069566063dSJacob Faibussowitsch PetscCall(DMDANaturalToGlobalBegin(da, natural, INSERT_VALUES, xin)); 9079566063dSJacob Faibussowitsch PetscCall(DMDANaturalToGlobalEnd(da, natural, INSERT_VALUES, xin)); 9089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&natural)); 9099566063dSJacob Faibussowitsch PetscCall(PetscInfo(xin, "Loading vector from natural ordering into DMDA\n")); 9109566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, ((PetscObject)xin)->prefix, "-vecload_block_size", &bs, &flag)); 91148a46eb9SPierre 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)); 9123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 91347c6ae99SBarry Smith } 91447c6ae99SBarry Smith 915d71ae5a4SJacob Faibussowitsch PetscErrorCode VecLoad_Default_DA(Vec xin, PetscViewer viewer) 916d71ae5a4SJacob Faibussowitsch { 9179a42bb27SBarry Smith DM da; 91847c6ae99SBarry Smith PetscBool isbinary; 91947c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 92047c6ae99SBarry Smith PetscBool ishdf5; 92147c6ae99SBarry Smith #endif 92247c6ae99SBarry Smith 92347c6ae99SBarry Smith PetscFunctionBegin; 9249566063dSJacob Faibussowitsch PetscCall(VecGetDM(xin, &da)); 9257a8be351SBarry Smith PetscCheck(da, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA"); 92647c6ae99SBarry Smith 92747c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 9289566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5)); 92947c6ae99SBarry Smith #endif 9309566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 93147c6ae99SBarry Smith 93247c6ae99SBarry Smith if (isbinary) { 9339566063dSJacob Faibussowitsch PetscCall(VecLoad_Binary_DA(xin, viewer)); 93447c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 93547c6ae99SBarry Smith } else if (ishdf5) { 9369566063dSJacob Faibussowitsch PetscCall(VecLoad_HDF5_DA(xin, viewer)); 93747c6ae99SBarry Smith #endif 93898921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Viewer type %s not supported for vector loading", ((PetscObject)viewer)->type_name); 9393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 94047c6ae99SBarry Smith } 941