xref: /petsc/src/dm/impls/da/gr2.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
147c6ae99SBarry Smith 
247c6ae99SBarry Smith /*
3aa219208SBarry Smith    Plots vectors obtained with DMDACreate2d()
447c6ae99SBarry Smith */
547c6ae99SBarry Smith 
6af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I  "petscdmda.h"   I*/
78135c375SStefano Zampini #include <petsc/private/glvisvecimpl.h>
82d95be3dSVaclav Hapla #include <petsc/private/viewerhdf5impl.h>
99804daf3SBarry Smith #include <petscdraw.h>
1047c6ae99SBarry Smith 
1147c6ae99SBarry Smith /*
1247c6ae99SBarry Smith         The data that is passed into the graphics callback
1347c6ae99SBarry Smith */
1447c6ae99SBarry Smith typedef struct {
15e5ab1681SLisandro Dalcin   PetscMPIInt        rank;
16e5ab1681SLisandro Dalcin   PetscInt           m, n, dof, k;
17e5ab1681SLisandro Dalcin   PetscReal          xmin, xmax, ymin, ymax, min, max;
185edff71fSBarry Smith   const PetscScalar *xy, *v;
197fe7c8feSLisandro Dalcin   PetscBool          showaxis, showgrid;
20109c9344SBarry Smith   const char        *name0, *name1;
2147c6ae99SBarry Smith } ZoomCtx;
2247c6ae99SBarry Smith 
2347c6ae99SBarry Smith /*
2447c6ae99SBarry Smith        This does the drawing for one particular field
2547c6ae99SBarry Smith     in one particular set of coordinates. It is a callback
2647c6ae99SBarry Smith     called from PetscDrawZoom()
2747c6ae99SBarry Smith */
289371c9d4SSatish Balay PetscErrorCode VecView_MPI_Draw_DA2d_Zoom(PetscDraw draw, void *ctx) {
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);
11247c6ae99SBarry Smith   PetscFunctionReturn(0);
11347c6ae99SBarry Smith }
11447c6ae99SBarry Smith 
1159371c9d4SSatish Balay PetscErrorCode VecView_MPI_Draw_DA2d(Vec xin, PetscViewer viewer) {
1169a42bb27SBarry Smith   DM                  da, dac, dag;
117a74ba6f7SBarry Smith   PetscInt            N, s, M, w, ncoors = 4;
11847c6ae99SBarry Smith   const PetscInt     *lx, *ly;
119e5ab1681SLisandro Dalcin   PetscReal           coors[4];
12047c6ae99SBarry Smith   PetscDraw           draw, popup;
12147c6ae99SBarry Smith   PetscBool           isnull, useports = PETSC_FALSE;
12247c6ae99SBarry Smith   MPI_Comm            comm;
12347c6ae99SBarry Smith   Vec                 xlocal, xcoor, xcoorl;
124bff4a2f0SMatthew G. Knepley   DMBoundaryType      bx, by;
125aa219208SBarry Smith   DMDAStencilType     st;
12647c6ae99SBarry Smith   ZoomCtx             zctx;
1270298fd71SBarry Smith   PetscDrawViewPorts *ports = NULL;
12847c6ae99SBarry Smith   PetscViewerFormat   format;
12920d0051dSBarry Smith   PetscInt           *displayfields;
13067dd0837SBarry Smith   PetscInt            ndisplayfields, i, nbounds;
13167dd0837SBarry Smith   const PetscReal    *bounds;
13247c6ae99SBarry Smith 
13347c6ae99SBarry Smith   PetscFunctionBegin;
13447c6ae99SBarry Smith   zctx.showgrid = PETSC_FALSE;
1357fe7c8feSLisandro Dalcin   zctx.showaxis = PETSC_TRUE;
1368865f1eaSKarl Rupp 
1379566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
1389566063dSJacob Faibussowitsch   PetscCall(PetscDrawIsNull(draw, &isnull));
1395b399a63SLisandro Dalcin   if (isnull) PetscFunctionReturn(0);
1405b399a63SLisandro Dalcin 
1419566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetBounds(viewer, &nbounds, &bounds));
14247c6ae99SBarry Smith 
1439566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin, &da));
1447a8be351SBarry Smith   PetscCheck(da, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA");
14547c6ae99SBarry Smith 
1469566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)xin, &comm));
1479566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &zctx.rank));
14847c6ae99SBarry Smith 
1499566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, &M, &N, NULL, &zctx.m, &zctx.n, NULL, &w, &s, &bx, &by, NULL, &st));
1509566063dSJacob Faibussowitsch   PetscCall(DMDAGetOwnershipRanges(da, &lx, &ly, NULL));
15147c6ae99SBarry Smith 
15247c6ae99SBarry Smith   /*
15347c6ae99SBarry Smith      Obtain a sequential vector that is going to contain the local values plus ONE layer of
154aa219208SBarry Smith      ghosted values to draw the graphics from. We also need its corresponding DMDA (dac) that will
15547c6ae99SBarry Smith      update the local values pluse ONE layer of ghost values.
15647c6ae99SBarry Smith   */
1579566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)da, "GraphicsGhosted", (PetscObject *)&xlocal));
15847c6ae99SBarry Smith   if (!xlocal) {
159bff4a2f0SMatthew G. Knepley     if (bx != DM_BOUNDARY_NONE || by != DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) {
16047c6ae99SBarry Smith       /*
16147c6ae99SBarry Smith          if original da is not of stencil width one, or periodic or not a box stencil then
162aa219208SBarry Smith          create a special DMDA to handle one level of ghost points for graphics
16347c6ae99SBarry Smith       */
1649566063dSJacob Faibussowitsch       PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, zctx.m, zctx.n, w, 1, lx, ly, &dac));
1659566063dSJacob Faibussowitsch       PetscCall(DMSetUp(dac));
1666aad120cSJose E. Roman       PetscCall(PetscInfo(da, "Creating auxiliary DMDA for managing graphics ghost points\n"));
16747c6ae99SBarry Smith     } else {
16847c6ae99SBarry Smith       /* otherwise we can use the da we already have */
16947c6ae99SBarry Smith       dac = da;
17047c6ae99SBarry Smith     }
17147c6ae99SBarry Smith     /* create local vector for holding ghosted values used in graphics */
1729566063dSJacob Faibussowitsch     PetscCall(DMCreateLocalVector(dac, &xlocal));
17347c6ae99SBarry Smith     if (dac != da) {
174aa219208SBarry Smith       /* don't keep any public reference of this DMDA, is is only available through xlocal */
1759566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)dac));
17647c6ae99SBarry Smith     } else {
17747c6ae99SBarry Smith       /* remove association between xlocal and da, because below we compose in the opposite
17847c6ae99SBarry Smith          direction and if we left this connect we'd get a loop, so the objects could
17947c6ae99SBarry Smith          never be destroyed */
1809566063dSJacob Faibussowitsch       PetscCall(PetscObjectRemoveReference((PetscObject)xlocal, "__PETSc_dm"));
18147c6ae99SBarry Smith     }
1829566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)da, "GraphicsGhosted", (PetscObject)xlocal));
1839566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)xlocal));
18447c6ae99SBarry Smith   } else {
185bff4a2f0SMatthew G. Knepley     if (bx != DM_BOUNDARY_NONE || by != DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) {
1869566063dSJacob Faibussowitsch       PetscCall(VecGetDM(xlocal, &dac));
187f7923d8aSBarry Smith     } else {
188f7923d8aSBarry Smith       dac = da;
18947c6ae99SBarry Smith     }
19047c6ae99SBarry Smith   }
19147c6ae99SBarry Smith 
19247c6ae99SBarry Smith   /*
19347c6ae99SBarry Smith       Get local (ghosted) values of vector
19447c6ae99SBarry Smith   */
1959566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dac, xin, INSERT_VALUES, xlocal));
1969566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dac, xin, INSERT_VALUES, xlocal));
1979566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xlocal, &zctx.v));
19847c6ae99SBarry Smith 
199832b7cebSLisandro Dalcin   /*
200832b7cebSLisandro Dalcin       Get coordinates of nodes
201832b7cebSLisandro Dalcin   */
2029566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &xcoor));
20347c6ae99SBarry Smith   if (!xcoor) {
2049566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 0.0));
2059566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(da, &xcoor));
20647c6ae99SBarry Smith   }
20747c6ae99SBarry Smith 
20847c6ae99SBarry Smith   /*
20947c6ae99SBarry Smith       Determine the min and max coordinates in plot
21047c6ae99SBarry Smith   */
2119566063dSJacob Faibussowitsch   PetscCall(VecStrideMin(xcoor, 0, NULL, &zctx.xmin));
2129566063dSJacob Faibussowitsch   PetscCall(VecStrideMax(xcoor, 0, NULL, &zctx.xmax));
2139566063dSJacob Faibussowitsch   PetscCall(VecStrideMin(xcoor, 1, NULL, &zctx.ymin));
2149566063dSJacob Faibussowitsch   PetscCall(VecStrideMax(xcoor, 1, NULL, &zctx.ymax));
2159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_contour_axis", &zctx.showaxis, NULL));
2167fe7c8feSLisandro Dalcin   if (zctx.showaxis) {
2179371c9d4SSatish Balay     coors[0] = zctx.xmin - .05 * (zctx.xmax - zctx.xmin);
2189371c9d4SSatish Balay     coors[1] = zctx.ymin - .05 * (zctx.ymax - zctx.ymin);
2199371c9d4SSatish Balay     coors[2] = zctx.xmax + .05 * (zctx.xmax - zctx.xmin);
2209371c9d4SSatish Balay     coors[3] = zctx.ymax + .05 * (zctx.ymax - zctx.ymin);
2217fe7c8feSLisandro Dalcin   } else {
2229371c9d4SSatish Balay     coors[0] = zctx.xmin;
2239371c9d4SSatish Balay     coors[1] = zctx.ymin;
2249371c9d4SSatish Balay     coors[2] = zctx.xmax;
2259371c9d4SSatish Balay     coors[3] = zctx.ymax;
2267fe7c8feSLisandro Dalcin   }
2279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-draw_coordinates", coors, &ncoors, NULL));
2289566063dSJacob 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]));
229a74ba6f7SBarry Smith 
23047c6ae99SBarry Smith   /*
231832b7cebSLisandro Dalcin       Get local ghosted version of coordinates
23247c6ae99SBarry Smith   */
2339566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)da, "GraphicsCoordinateGhosted", (PetscObject *)&xcoorl));
23447c6ae99SBarry Smith   if (!xcoorl) {
235aa219208SBarry Smith     /* create DMDA to get local version of graphics */
2369566063dSJacob Faibussowitsch     PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, zctx.m, zctx.n, 2, 1, lx, ly, &dag));
2379566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dag));
2386aad120cSJose E. Roman     PetscCall(PetscInfo(dag, "Creating auxiliary DMDA for managing graphics coordinates ghost points\n"));
2399566063dSJacob Faibussowitsch     PetscCall(DMCreateLocalVector(dag, &xcoorl));
2409566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)da, "GraphicsCoordinateGhosted", (PetscObject)xcoorl));
2419566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)dag));
2429566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)xcoorl));
2431baa6e33SBarry Smith   } else PetscCall(VecGetDM(xcoorl, &dag));
2449566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dag, xcoor, INSERT_VALUES, xcoorl));
2459566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dag, xcoor, INSERT_VALUES, xcoorl));
2469566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xcoorl, &zctx.xy));
2479566063dSJacob Faibussowitsch   PetscCall(DMDAGetCoordinateName(da, 0, &zctx.name0));
2489566063dSJacob Faibussowitsch   PetscCall(DMDAGetCoordinateName(da, 1, &zctx.name1));
24947c6ae99SBarry Smith 
25047c6ae99SBarry Smith   /*
25147c6ae99SBarry Smith       Get information about size of area each processor must do graphics for
25247c6ae99SBarry Smith   */
2539566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac, NULL, &M, &N, NULL, NULL, NULL, NULL, &zctx.dof, NULL, &bx, &by, NULL, NULL));
2549566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dac, NULL, NULL, NULL, &zctx.m, &zctx.n, NULL));
2559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_contour_grid", &zctx.showgrid, NULL));
2564e6118eeSBarry Smith 
2579566063dSJacob Faibussowitsch   PetscCall(DMDASelectFields(da, &ndisplayfields, &displayfields));
2589566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
2599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_ports", &useports, NULL));
260832b7cebSLisandro Dalcin   if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE;
261832b7cebSLisandro Dalcin   if (useports) {
2629566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
2639566063dSJacob Faibussowitsch     PetscCall(PetscDrawCheckResizedWindow(draw));
2649566063dSJacob Faibussowitsch     PetscCall(PetscDrawClear(draw));
2659566063dSJacob Faibussowitsch     PetscCall(PetscDrawViewPortsCreate(draw, ndisplayfields, &ports));
266e5ab1681SLisandro Dalcin   }
26720d0051dSBarry Smith 
268832b7cebSLisandro Dalcin   /*
269832b7cebSLisandro Dalcin       Loop over each field; drawing each in a different window
270832b7cebSLisandro Dalcin   */
27120d0051dSBarry Smith   for (i = 0; i < ndisplayfields; i++) {
27220d0051dSBarry Smith     zctx.k = displayfields[i];
2735b399a63SLisandro Dalcin 
274e5ab1681SLisandro Dalcin     /* determine the min and max value in plot */
2759566063dSJacob Faibussowitsch     PetscCall(VecStrideMin(xin, zctx.k, NULL, &zctx.min));
2769566063dSJacob Faibussowitsch     PetscCall(VecStrideMax(xin, zctx.k, NULL, &zctx.max));
27767dd0837SBarry Smith     if (zctx.k < nbounds) {
278f3f0eb19SBarry Smith       zctx.min = bounds[2 * zctx.k];
279f3f0eb19SBarry Smith       zctx.max = bounds[2 * zctx.k + 1];
28067dd0837SBarry Smith     }
28147c6ae99SBarry Smith     if (zctx.min == zctx.max) {
28247c6ae99SBarry Smith       zctx.min -= 1.e-12;
28347c6ae99SBarry Smith       zctx.max += 1.e-12;
28447c6ae99SBarry Smith     }
2859566063dSJacob Faibussowitsch     PetscCall(PetscInfo(da, "DMDA 2d contour plot min %g max %g\n", (double)zctx.min, (double)zctx.max));
28647c6ae99SBarry Smith 
287832b7cebSLisandro Dalcin     if (useports) {
2889566063dSJacob Faibussowitsch       PetscCall(PetscDrawViewPortsSet(ports, i));
289832b7cebSLisandro Dalcin     } else {
290832b7cebSLisandro Dalcin       const char *title;
2919566063dSJacob Faibussowitsch       PetscCall(PetscViewerDrawGetDraw(viewer, i, &draw));
2929566063dSJacob Faibussowitsch       PetscCall(DMDAGetFieldName(da, zctx.k, &title));
2939566063dSJacob Faibussowitsch       if (title) PetscCall(PetscDrawSetTitle(draw, title));
294832b7cebSLisandro Dalcin     }
295832b7cebSLisandro Dalcin 
2969566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetPopup(draw, &popup));
2979566063dSJacob Faibussowitsch     PetscCall(PetscDrawScalePopup(popup, zctx.min, zctx.max));
2989566063dSJacob Faibussowitsch     PetscCall(PetscDrawSetCoordinates(draw, coors[0], coors[1], coors[2], coors[3]));
2999566063dSJacob Faibussowitsch     PetscCall(PetscDrawZoom(draw, VecView_MPI_Draw_DA2d_Zoom, &zctx));
3009566063dSJacob Faibussowitsch     if (!useports) PetscCall(PetscDrawSave(draw));
30147c6ae99SBarry Smith   }
302832b7cebSLisandro Dalcin   if (useports) {
3039566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
3049566063dSJacob Faibussowitsch     PetscCall(PetscDrawSave(draw));
305832b7cebSLisandro Dalcin   }
306832b7cebSLisandro Dalcin 
3079566063dSJacob Faibussowitsch   PetscCall(PetscDrawViewPortsDestroy(ports));
3089566063dSJacob Faibussowitsch   PetscCall(PetscFree(displayfields));
3099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xcoorl, &zctx.xy));
3109566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xlocal, &zctx.v));
31147c6ae99SBarry Smith   PetscFunctionReturn(0);
31247c6ae99SBarry Smith }
31347c6ae99SBarry Smith 
3140da24e51SJuha Jäykkä #if defined(PETSC_HAVE_HDF5)
3159371c9d4SSatish Balay static PetscErrorCode VecGetHDF5ChunkSize(DM_DA *da, Vec xin, PetscInt dimension, PetscInt timestep, hsize_t *chunkDims) {
316d4190030SJed Brown   PetscMPIInt comm_size;
317561a82e7SJed Brown   hsize_t     chunk_size, target_size, dim;
318561a82e7SJed Brown   hsize_t     vec_size = sizeof(PetscScalar) * da->M * da->N * da->P * da->w;
319561a82e7SJed Brown   hsize_t     avg_local_vec_size, KiB = 1024, MiB = KiB * KiB, GiB = MiB * KiB, min_size = MiB;
320561a82e7SJed Brown   hsize_t     max_chunks     = 64 * KiB; /* HDF5 internal limitation */
321561a82e7SJed Brown   hsize_t     max_chunk_size = 4 * GiB;  /* HDF5 internal limitation */
32284179ffaSJed Brown   PetscInt    zslices = da->p, yslices = da->n, xslices = da->m;
3230da24e51SJuha Jäykkä 
324cb5c4f94SJuha Jäykkä   PetscFunctionBegin;
3259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)xin), &comm_size));
326faa75363SBarry Smith   avg_local_vec_size = (hsize_t)PetscCeilInt(vec_size, comm_size); /* we will attempt to use this as the chunk size */
3270da24e51SJuha Jäykkä 
328faa75363SBarry 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))));
3297d310018SBarry 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 */
3307d310018SBarry 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);
3310da24e51SJuha Jäykkä 
332cb5c4f94SJuha Jäykkä   /*
333cb5c4f94SJuha Jäykkä    if size/rank > max_chunk_size, we need radical measures: even going down to
334cb5c4f94SJuha Jäykkä    avg_local_vec_size is not enough, so we simply use chunk size of 4 GiB no matter
335cb5c4f94SJuha Jäykkä    what, composed in the most efficient way possible.
336cb5c4f94SJuha Jäykkä    N.B. this minimises the number of chunks, which may or may not be the optimal
337cb5c4f94SJuha Jäykkä    solution. In a BG, for example, the optimal solution is probably to make # chunks = #
338cb5c4f94SJuha Jäykkä    IO nodes involved, but this author has no access to a BG to figure out how to
339cb5c4f94SJuha Jäykkä    reliably find the right number. And even then it may or may not be enough.
340cb5c4f94SJuha Jäykkä    */
3410da24e51SJuha Jäykkä   if (avg_local_vec_size > max_chunk_size) {
342cb5c4f94SJuha Jäykkä     /* check if we can just split local z-axis: is that enough? */
343faa75363SBarry Smith     zslices = PetscCeilInt(vec_size, da->p * max_chunk_size) * zslices;
3440da24e51SJuha Jäykkä     if (zslices > da->P) {
345cb5c4f94SJuha Jäykkä       /* lattice is too large in xy-directions, splitting z only is not enough */
3460da24e51SJuha Jäykkä       zslices = da->P;
347faa75363SBarry Smith       yslices = PetscCeilInt(vec_size, zslices * da->n * max_chunk_size) * yslices;
3480da24e51SJuha Jäykkä       if (yslices > da->N) {
349cb5c4f94SJuha Jäykkä         /* lattice is too large in x-direction, splitting along z, y is not enough */
3500da24e51SJuha Jäykkä         yslices = da->N;
351faa75363SBarry Smith         xslices = PetscCeilInt(vec_size, zslices * yslices * da->m * max_chunk_size) * xslices;
3520da24e51SJuha Jäykkä       }
3530da24e51SJuha Jäykkä     }
3540da24e51SJuha Jäykkä     dim = 0;
3559371c9d4SSatish Balay     if (timestep >= 0) { ++dim; }
356cb5c4f94SJuha Jäykkä     /* prefer to split z-axis, even down to planar slices */
357c73cfb54SMatthew G. Knepley     if (dimension == 3) {
358cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t)da->P / zslices;
359cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t)da->N / yslices;
360cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t)da->M / xslices;
3610da24e51SJuha Jäykkä     } else {
362cb5c4f94SJuha Jäykkä       /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */
363cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t)da->N / yslices;
364cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t)da->M / xslices;
3650da24e51SJuha Jäykkä     }
3660da24e51SJuha 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);
3670da24e51SJuha Jäykkä   } else {
3680da24e51SJuha Jäykkä     if (target_size < chunk_size) {
369cb5c4f94SJuha Jäykkä       /* only change the defaults if target_size < chunk_size */
3700da24e51SJuha Jäykkä       dim = 0;
3719371c9d4SSatish Balay       if (timestep >= 0) { ++dim; }
372cb5c4f94SJuha Jäykkä       /* prefer to split z-axis, even down to planar slices */
373c73cfb54SMatthew G. Knepley       if (dimension == 3) {
374cb5c4f94SJuha Jäykkä         /* try splitting the z-axis to core-size bits, i.e. divide chunk size by # comm_size in z-direction */
3750da24e51SJuha Jäykkä         if (target_size >= chunk_size / da->p) {
376cb5c4f94SJuha Jäykkä           /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */
377faa75363SBarry Smith           chunkDims[dim] = (hsize_t)PetscCeilInt(da->P, da->p);
3780da24e51SJuha Jäykkä         } else {
379cb5c4f94SJuha Jäykkä           /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be
380cb5c4f94SJuha Jäykkä            radical and let everyone write all they've got */
381faa75363SBarry Smith           chunkDims[dim++] = (hsize_t)PetscCeilInt(da->P, da->p);
382faa75363SBarry Smith           chunkDims[dim++] = (hsize_t)PetscCeilInt(da->N, da->n);
383faa75363SBarry Smith           chunkDims[dim++] = (hsize_t)PetscCeilInt(da->M, da->m);
3840da24e51SJuha Jäykkä         }
3850da24e51SJuha Jäykkä       } else {
386cb5c4f94SJuha Jäykkä         /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */
3870da24e51SJuha Jäykkä         if (target_size >= chunk_size / da->n) {
388cb5c4f94SJuha Jäykkä           /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */
389faa75363SBarry Smith           chunkDims[dim] = (hsize_t)PetscCeilInt(da->N, da->n);
3900da24e51SJuha Jäykkä         } else {
391cb5c4f94SJuha Jäykkä           /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be
392cb5c4f94SJuha Jäykkä            radical and let everyone write all they've got */
393faa75363SBarry Smith           chunkDims[dim++] = (hsize_t)PetscCeilInt(da->N, da->n);
394faa75363SBarry Smith           chunkDims[dim++] = (hsize_t)PetscCeilInt(da->M, da->m);
3950da24e51SJuha Jäykkä         }
3960da24e51SJuha Jäykkä       }
3970da24e51SJuha 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);
3980da24e51SJuha Jäykkä     } else {
399cb5c4f94SJuha Jäykkä       /* precomputed chunks are fine, we don't need to do anything */
4000da24e51SJuha Jäykkä     }
4010da24e51SJuha Jäykkä   }
4020da24e51SJuha Jäykkä   PetscFunctionReturn(0);
4030da24e51SJuha Jäykkä }
4040da24e51SJuha Jäykkä #endif
40547c6ae99SBarry Smith 
40647c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
4079371c9d4SSatish Balay PetscErrorCode VecView_MPI_HDF5_DA(Vec xin, PetscViewer viewer) {
4082d95be3dSVaclav Hapla   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5 *)viewer->data;
4099b2a5a86SJed Brown   DM                 dm;
4109b2a5a86SJed Brown   DM_DA             *da;
41147c6ae99SBarry Smith   hid_t              filespace;  /* file dataspace identifier */
4128e2ae6d7SMichael Kraus   hid_t              chunkspace; /* chunk dataset property identifier */
41347c6ae99SBarry Smith   hid_t              dset_id;    /* dataset identifier */
41447c6ae99SBarry Smith   hid_t              memspace;   /* memory dataspace identifier */
41547c6ae99SBarry Smith   hid_t              file_id;
41615214e8eSMatthew G Knepley   hid_t              group;
4179a0502c6SHåkon Strandenes   hid_t              memscalartype;  /* scalar type for mem (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
4189a0502c6SHåkon Strandenes   hid_t              filescalartype; /* scalar type for file (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
419d9a4edebSJed Brown   hsize_t            dim;
4200da24e51SJuha 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  */
421d7dd068bSVaclav Hapla   PetscBool          timestepping = PETSC_FALSE, dim2, spoutput;
422d7dd068bSVaclav Hapla   PetscInt           timestep     = PETSC_MIN_INT, dimension;
4235edff71fSBarry Smith   const PetscScalar *x;
4248e2ae6d7SMichael Kraus   const char        *vecname;
42547c6ae99SBarry Smith 
42647c6ae99SBarry Smith   PetscFunctionBegin;
4279566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5OpenGroup(viewer, &file_id, &group));
4289566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5IsTimestepping(viewer, &timestepping));
429*48a46eb9SPierre Jolivet   if (timestepping) PetscCall(PetscViewerHDF5GetTimestep(viewer, &timestep));
4309566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetBaseDimension2(viewer, &dim2));
4319566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetSPOutput(viewer, &spoutput));
43215214e8eSMatthew G Knepley 
4339566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin, &dm));
4347a8be351SBarry Smith   PetscCheck(dm, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA");
4359b2a5a86SJed Brown   da = (DM_DA *)dm->data;
4369566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dimension));
43747c6ae99SBarry Smith 
4388e2ae6d7SMichael Kraus   /* Create the dataspace for the dataset.
4398e2ae6d7SMichael Kraus    *
4408e2ae6d7SMichael Kraus    * dims - holds the current dimensions of the dataset
4418e2ae6d7SMichael Kraus    *
4428e2ae6d7SMichael Kraus    * maxDims - holds the maximum dimensions of the dataset (unlimited
4438e2ae6d7SMichael Kraus    * for the number of time steps with the current dimensions for the
4448e2ae6d7SMichael Kraus    * other dimensions; so only additional time steps can be added).
4458e2ae6d7SMichael Kraus    *
4468e2ae6d7SMichael Kraus    * chunkDims - holds the size of a single time step (required to
4478e2ae6d7SMichael Kraus    * permit extending dataset).
4488e2ae6d7SMichael Kraus    */
4498e2ae6d7SMichael Kraus   dim = 0;
4508e2ae6d7SMichael Kraus   if (timestep >= 0) {
4518e2ae6d7SMichael Kraus     dims[dim]      = timestep + 1;
4528e2ae6d7SMichael Kraus     maxDims[dim]   = H5S_UNLIMITED;
4538e2ae6d7SMichael Kraus     chunkDims[dim] = 1;
4548e2ae6d7SMichael Kraus     ++dim;
4558e2ae6d7SMichael Kraus   }
456c73cfb54SMatthew G. Knepley   if (dimension == 3) {
4579566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(da->P, dims + dim));
4588e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4598e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4608e2ae6d7SMichael Kraus     ++dim;
4618e2ae6d7SMichael Kraus   }
462c73cfb54SMatthew G. Knepley   if (dimension > 1) {
4639566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(da->N, dims + dim));
4648e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4658e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4668e2ae6d7SMichael Kraus     ++dim;
4678e2ae6d7SMichael Kraus   }
4689566063dSJacob Faibussowitsch   PetscCall(PetscHDF5IntCast(da->M, dims + dim));
4698e2ae6d7SMichael Kraus   maxDims[dim]   = dims[dim];
4708e2ae6d7SMichael Kraus   chunkDims[dim] = dims[dim];
4718e2ae6d7SMichael Kraus   ++dim;
47282ea9b62SBarry Smith   if (da->w > 1 || dim2) {
4739566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(da->w, dims + dim));
4748e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4758e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4768e2ae6d7SMichael Kraus     ++dim;
4778e2ae6d7SMichael Kraus   }
47847c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
4798e2ae6d7SMichael Kraus   dims[dim]      = 2;
4808e2ae6d7SMichael Kraus   maxDims[dim]   = dims[dim];
4818e2ae6d7SMichael Kraus   chunkDims[dim] = dims[dim];
4828e2ae6d7SMichael Kraus   ++dim;
48347c6ae99SBarry Smith #endif
4840da24e51SJuha Jäykkä 
4859566063dSJacob Faibussowitsch   PetscCall(VecGetHDF5ChunkSize(da, xin, dimension, timestep, chunkDims));
4860da24e51SJuha Jäykkä 
487792fecdfSBarry Smith   PetscCallHDF5Return(filespace, H5Screate_simple, (dim, dims, maxDims));
48847c6ae99SBarry Smith 
48915214e8eSMatthew G Knepley #if defined(PETSC_USE_REAL_SINGLE)
4909a0502c6SHåkon Strandenes   memscalartype  = H5T_NATIVE_FLOAT;
4919a0502c6SHåkon Strandenes   filescalartype = H5T_NATIVE_FLOAT;
49215214e8eSMatthew G Knepley #elif defined(PETSC_USE_REAL___FLOAT128)
49315214e8eSMatthew G Knepley #error "HDF5 output with 128 bit floats not supported."
494570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
495570b7f6dSBarry Smith #error "HDF5 output with 16 bit floats not supported."
49615214e8eSMatthew G Knepley #else
4979a0502c6SHåkon Strandenes   memscalartype = H5T_NATIVE_DOUBLE;
4989a0502c6SHåkon Strandenes   if (spoutput == PETSC_TRUE) filescalartype = H5T_NATIVE_FLOAT;
4999a0502c6SHåkon Strandenes   else filescalartype = H5T_NATIVE_DOUBLE;
50015214e8eSMatthew G Knepley #endif
50115214e8eSMatthew G Knepley 
50247c6ae99SBarry Smith   /* Create the dataset with default properties and close filespace */
5039566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)xin, &vecname));
50415214e8eSMatthew G Knepley   if (!H5Lexists(group, vecname, H5P_DEFAULT)) {
5058e2ae6d7SMichael Kraus     /* Create chunk */
506792fecdfSBarry Smith     PetscCallHDF5Return(chunkspace, H5Pcreate, (H5P_DATASET_CREATE));
507792fecdfSBarry Smith     PetscCallHDF5(H5Pset_chunk, (chunkspace, dim, chunkDims));
5088e2ae6d7SMichael Kraus 
509792fecdfSBarry Smith     PetscCallHDF5Return(dset_id, H5Dcreate2, (group, vecname, filescalartype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT));
51015214e8eSMatthew G Knepley   } else {
511792fecdfSBarry Smith     PetscCallHDF5Return(dset_id, H5Dopen2, (group, vecname, H5P_DEFAULT));
512792fecdfSBarry Smith     PetscCallHDF5(H5Dset_extent, (dset_id, dims));
51315214e8eSMatthew G Knepley   }
514792fecdfSBarry Smith   PetscCallHDF5(H5Sclose, (filespace));
51547c6ae99SBarry Smith 
51647c6ae99SBarry Smith   /* Each process defines a dataset and writes it to the hyperslab in the file */
5178e2ae6d7SMichael Kraus   dim = 0;
5188e2ae6d7SMichael Kraus   if (timestep >= 0) {
5198e2ae6d7SMichael Kraus     offset[dim] = timestep;
5208e2ae6d7SMichael Kraus     ++dim;
5218e2ae6d7SMichael Kraus   }
5229566063dSJacob Faibussowitsch   if (dimension == 3) PetscCall(PetscHDF5IntCast(da->zs, offset + dim++));
5239566063dSJacob Faibussowitsch   if (dimension > 1) PetscCall(PetscHDF5IntCast(da->ys, offset + dim++));
5249566063dSJacob Faibussowitsch   PetscCall(PetscHDF5IntCast(da->xs / da->w, offset + dim++));
52582ea9b62SBarry Smith   if (da->w > 1 || dim2) offset[dim++] = 0;
52647c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5278e2ae6d7SMichael Kraus   offset[dim++] = 0;
52847c6ae99SBarry Smith #endif
5298e2ae6d7SMichael Kraus   dim = 0;
5308e2ae6d7SMichael Kraus   if (timestep >= 0) {
5318e2ae6d7SMichael Kraus     count[dim] = 1;
5328e2ae6d7SMichael Kraus     ++dim;
5338e2ae6d7SMichael Kraus   }
5349566063dSJacob Faibussowitsch   if (dimension == 3) PetscCall(PetscHDF5IntCast(da->ze - da->zs, count + dim++));
5359566063dSJacob Faibussowitsch   if (dimension > 1) PetscCall(PetscHDF5IntCast(da->ye - da->ys, count + dim++));
5369566063dSJacob Faibussowitsch   PetscCall(PetscHDF5IntCast((da->xe - da->xs) / da->w, count + dim++));
5379566063dSJacob Faibussowitsch   if (da->w > 1 || dim2) PetscCall(PetscHDF5IntCast(da->w, count + dim++));
53847c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5398e2ae6d7SMichael Kraus   count[dim++] = 2;
54047c6ae99SBarry Smith #endif
541792fecdfSBarry Smith   PetscCallHDF5Return(memspace, H5Screate_simple, (dim, count, NULL));
542792fecdfSBarry Smith   PetscCallHDF5Return(filespace, H5Dget_space, (dset_id));
543792fecdfSBarry Smith   PetscCallHDF5(H5Sselect_hyperslab, (filespace, H5S_SELECT_SET, offset, NULL, count, NULL));
54447c6ae99SBarry Smith 
5459566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xin, &x));
546792fecdfSBarry Smith   PetscCallHDF5(H5Dwrite, (dset_id, memscalartype, memspace, filespace, hdf5->dxpl_id, x));
547792fecdfSBarry Smith   PetscCallHDF5(H5Fflush, (file_id, H5F_SCOPE_GLOBAL));
5489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xin, &x));
54947c6ae99SBarry Smith 
550e0552f36Ssajid__ali #if defined(PETSC_USE_COMPLEX)
551e0552f36Ssajid__ali   {
552e0552f36Ssajid__ali     PetscBool tru = PETSC_TRUE;
5539566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)xin, "complex", PETSC_BOOL, &tru));
554e0552f36Ssajid__ali   }
555e0552f36Ssajid__ali #endif
556*48a46eb9SPierre Jolivet   if (timestepping) PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)xin, "timestepping", PETSC_BOOL, &timestepping));
557e0552f36Ssajid__ali 
55847c6ae99SBarry Smith   /* Close/release resources */
5599371c9d4SSatish Balay   if (group != file_id) { PetscCallHDF5(H5Gclose, (group)); }
560792fecdfSBarry Smith   PetscCallHDF5(H5Sclose, (filespace));
561792fecdfSBarry Smith   PetscCallHDF5(H5Sclose, (memspace));
562792fecdfSBarry Smith   PetscCallHDF5(H5Dclose, (dset_id));
5639566063dSJacob Faibussowitsch   PetscCall(PetscInfo(xin, "Wrote Vec object with name %s\n", vecname));
56447c6ae99SBarry Smith   PetscFunctionReturn(0);
56547c6ae99SBarry Smith }
56647c6ae99SBarry Smith #endif
56747c6ae99SBarry Smith 
56809573ac7SBarry Smith extern PetscErrorCode VecView_MPI_Draw_DA1d(Vec, PetscViewer);
56947c6ae99SBarry Smith 
57047c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
5719371c9d4SSatish Balay static PetscErrorCode DMDAArrayMPIIO(DM da, PetscViewer viewer, Vec xin, PetscBool write) {
57247c6ae99SBarry Smith   MPI_File           mfdes;
57347c6ae99SBarry Smith   PetscMPIInt        gsizes[4], lsizes[4], lstarts[4], asiz, dof;
57447c6ae99SBarry Smith   MPI_Datatype       view;
57547c6ae99SBarry Smith   const PetscScalar *array;
57647c6ae99SBarry Smith   MPI_Offset         off;
57747c6ae99SBarry Smith   MPI_Aint           ub, ul;
57847c6ae99SBarry Smith   PetscInt           type, rows, vecrows, tr[2];
57947c6ae99SBarry Smith   DM_DA             *dd = (DM_DA *)da->data;
5800c169764Sdmay   PetscBool          skipheader;
58147c6ae99SBarry Smith 
58247c6ae99SBarry Smith   PetscFunctionBegin;
5839566063dSJacob Faibussowitsch   PetscCall(VecGetSize(xin, &vecrows));
5849566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipheader));
58547c6ae99SBarry Smith   if (!write) {
58647c6ae99SBarry Smith     /* Read vector header. */
5870c169764Sdmay     if (!skipheader) {
5889566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryRead(viewer, tr, 2, NULL, PETSC_INT));
58947c6ae99SBarry Smith       type = tr[0];
59047c6ae99SBarry Smith       rows = tr[1];
59108401ef6SPierre Jolivet       PetscCheck(type == VEC_FILE_CLASSID, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONG, "Not vector next in file");
59208401ef6SPierre Jolivet       PetscCheck(rows == vecrows, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Vector in file not same size as DMDA vector");
5930c169764Sdmay     }
59447c6ae99SBarry Smith   } else {
59547c6ae99SBarry Smith     tr[0] = VEC_FILE_CLASSID;
59647c6ae99SBarry Smith     tr[1] = vecrows;
597*48a46eb9SPierre Jolivet     if (!skipheader) PetscCall(PetscViewerBinaryWrite(viewer, tr, 2, PETSC_INT));
5980c169764Sdmay   }
59947c6ae99SBarry Smith 
6009566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->w, &dof));
6014dc2109aSBarry Smith   gsizes[0] = dof;
6029566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->M, gsizes + 1));
6039566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->N, gsizes + 2));
6049566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->P, gsizes + 3));
6054dc2109aSBarry Smith   lsizes[0] = dof;
6069566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast((dd->xe - dd->xs) / dof, lsizes + 1));
6079566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->ye - dd->ys, lsizes + 2));
6089566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->ze - dd->zs, lsizes + 3));
6094dc2109aSBarry Smith   lstarts[0] = 0;
6109566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->xs / dof, lstarts + 1));
6119566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->ys, lstarts + 2));
6129566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->zs, lstarts + 3));
6139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_create_subarray(da->dim + 1, gsizes, lsizes, lstarts, MPI_ORDER_FORTRAN, MPIU_SCALAR, &view));
6149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&view));
61547c6ae99SBarry Smith 
6169566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetMPIIODescriptor(viewer, &mfdes));
6179566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetMPIIOOffset(viewer, &off));
6189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_File_set_view(mfdes, off, MPIU_SCALAR, view, (char *)"native", MPI_INFO_NULL));
6199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xin, &array));
62047c6ae99SBarry Smith   asiz = lsizes[1] * (lsizes[2] > 0 ? lsizes[2] : 1) * (lsizes[3] > 0 ? lsizes[3] : 1) * dof;
6211baa6e33SBarry Smith   if (write) PetscCall(MPIU_File_write_all(mfdes, (PetscScalar *)array, asiz, MPIU_SCALAR, MPI_STATUS_IGNORE));
6221baa6e33SBarry Smith   else PetscCall(MPIU_File_read_all(mfdes, (PetscScalar *)array, asiz, MPIU_SCALAR, MPI_STATUS_IGNORE));
6239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_get_extent(view, &ul, &ub));
6249566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryAddMPIIOOffset(viewer, ub));
6259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xin, &array));
6269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&view));
62747c6ae99SBarry Smith   PetscFunctionReturn(0);
62847c6ae99SBarry Smith }
62947c6ae99SBarry Smith #endif
63047c6ae99SBarry Smith 
6319371c9d4SSatish Balay PetscErrorCode VecView_MPI_DA(Vec xin, PetscViewer viewer) {
6329a42bb27SBarry Smith   DM        da;
63347c6ae99SBarry Smith   PetscInt  dim;
63447c6ae99SBarry Smith   Vec       natural;
6358135c375SStefano Zampini   PetscBool isdraw, isvtk, isglvis;
63647c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
63747c6ae99SBarry Smith   PetscBool ishdf5;
63847c6ae99SBarry Smith #endif
6393f3fd955SJed Brown   const char       *prefix, *name;
640a261c58fSBarry Smith   PetscViewerFormat format;
64147c6ae99SBarry Smith 
64247c6ae99SBarry Smith   PetscFunctionBegin;
6439566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin, &da));
6447a8be351SBarry Smith   PetscCheck(da, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA");
6459566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
6469566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
64747c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
6489566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
64947c6ae99SBarry Smith #endif
6509566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis));
65147c6ae99SBarry Smith   if (isdraw) {
6529566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
65347c6ae99SBarry Smith     if (dim == 1) {
6549566063dSJacob Faibussowitsch       PetscCall(VecView_MPI_Draw_DA1d(xin, viewer));
65547c6ae99SBarry Smith     } else if (dim == 2) {
6569566063dSJacob Faibussowitsch       PetscCall(VecView_MPI_Draw_DA2d(xin, viewer));
65763a3b9bcSJacob Faibussowitsch     } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Cannot graphically view vector associated with this dimensional DMDA %" PetscInt_FMT, dim);
65801d7c5c3SPatrick Sanan   } else if (isvtk) { /* Duplicate the Vec */
6594061b8bfSJed Brown     Vec Y;
6609566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xin, &Y));
661b51b94faSJed Brown     if (((PetscObject)xin)->name) {
662b51b94faSJed Brown       /* If xin was named, copy the name over to Y. The duplicate names are safe because nobody else will ever see Y. */
6639566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)Y, ((PetscObject)xin)->name));
664b51b94faSJed Brown     }
6659566063dSJacob Faibussowitsch     PetscCall(VecCopy(xin, Y));
666a8f87f1dSPatrick Sanan     {
667a8f87f1dSPatrick Sanan       PetscObject dmvtk;
668a8f87f1dSPatrick Sanan       PetscBool   compatible, compatibleSet;
6699566063dSJacob Faibussowitsch       PetscCall(PetscViewerVTKGetDM(viewer, &dmvtk));
670a8f87f1dSPatrick Sanan       if (dmvtk) {
671064a246eSJacob Faibussowitsch         PetscValidHeaderSpecific((DM)dmvtk, DM_CLASSID, 2);
6729566063dSJacob Faibussowitsch         PetscCall(DMGetCompatibility(da, (DM)dmvtk, &compatible, &compatibleSet));
6737a8be351SBarry 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.");
674a8f87f1dSPatrick Sanan       }
6759566063dSJacob Faibussowitsch       PetscCall(PetscViewerVTKAddField(viewer, (PetscObject)da, DMDAVTKWriteAll, PETSC_DEFAULT, PETSC_VTK_POINT_FIELD, PETSC_FALSE, (PetscObject)Y));
676a8f87f1dSPatrick Sanan     }
67747c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
67847c6ae99SBarry Smith   } else if (ishdf5) {
6799566063dSJacob Faibussowitsch     PetscCall(VecView_MPI_HDF5_DA(xin, viewer));
68047c6ae99SBarry Smith #endif
6818135c375SStefano Zampini   } else if (isglvis) {
6829566063dSJacob Faibussowitsch     PetscCall(VecView_GLVis(xin, viewer));
68347c6ae99SBarry Smith   } else {
68447c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
68547c6ae99SBarry Smith     PetscBool isbinary, isMPIIO;
68647c6ae99SBarry Smith 
6879566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
68847c6ae99SBarry Smith     if (isbinary) {
6899566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &isMPIIO));
69047c6ae99SBarry Smith       if (isMPIIO) {
6919566063dSJacob Faibussowitsch         PetscCall(DMDAArrayMPIIO(da, viewer, xin, PETSC_TRUE));
69247c6ae99SBarry Smith         PetscFunctionReturn(0);
69347c6ae99SBarry Smith       }
69447c6ae99SBarry Smith     }
69547c6ae99SBarry Smith #endif
69647c6ae99SBarry Smith 
69747c6ae99SBarry Smith     /* call viewer on natural ordering */
6989566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)xin, &prefix));
6999566063dSJacob Faibussowitsch     PetscCall(DMDACreateNaturalVector(da, &natural));
7009566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)natural, prefix));
7019566063dSJacob Faibussowitsch     PetscCall(DMDAGlobalToNaturalBegin(da, xin, INSERT_VALUES, natural));
7029566063dSJacob Faibussowitsch     PetscCall(DMDAGlobalToNaturalEnd(da, xin, INSERT_VALUES, natural));
7039566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)xin, &name));
7049566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)natural, name));
705a261c58fSBarry Smith 
7069566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
707a261c58fSBarry Smith     if (format == PETSC_VIEWER_BINARY_MATLAB) {
708a261c58fSBarry Smith       /* temporarily remove viewer format so it won't trigger in the VecView() */
7099566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_DEFAULT));
710a261c58fSBarry Smith     }
711a261c58fSBarry Smith 
71223f88b65SBarry Smith     ((PetscObject)natural)->donotPetscObjectPrintClassNamePrefixType = PETSC_TRUE;
7139566063dSJacob Faibussowitsch     PetscCall(VecView(natural, viewer));
71423f88b65SBarry Smith     ((PetscObject)natural)->donotPetscObjectPrintClassNamePrefixType = PETSC_FALSE;
715a261c58fSBarry Smith 
716a261c58fSBarry Smith     if (format == PETSC_VIEWER_BINARY_MATLAB) {
717a261c58fSBarry Smith       MPI_Comm    comm;
718a261c58fSBarry Smith       FILE       *info;
719a261c58fSBarry Smith       const char *fieldname;
720da88d4d4SJed Brown       char        fieldbuf[256];
721a261c58fSBarry Smith       PetscInt    dim, ni, nj, nk, pi, pj, pk, dof, n;
722a261c58fSBarry Smith 
723a261c58fSBarry Smith       /* set the viewer format back into the viewer */
7249566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(viewer));
7259566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
7269566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryGetInfoPointer(viewer, &info));
7279566063dSJacob Faibussowitsch       PetscCall(DMDAGetInfo(da, &dim, &ni, &nj, &nk, &pi, &pj, &pk, &dof, NULL, NULL, NULL, NULL, NULL));
7289566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, info, "#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
7299566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, info, "#$$ tmp = PetscBinaryRead(fd); \n"));
730*48a46eb9SPierre Jolivet       if (dim == 1) PetscCall(PetscFPrintf(comm, info, "#$$ tmp = reshape(tmp,%" PetscInt_FMT ",%" PetscInt_FMT ");\n", dof, ni));
731*48a46eb9SPierre Jolivet       if (dim == 2) PetscCall(PetscFPrintf(comm, info, "#$$ tmp = reshape(tmp,%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ");\n", dof, ni, nj));
732*48a46eb9SPierre Jolivet       if (dim == 3) PetscCall(PetscFPrintf(comm, info, "#$$ tmp = reshape(tmp,%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ");\n", dof, ni, nj, nk));
733a261c58fSBarry Smith 
734a261c58fSBarry Smith       for (n = 0; n < dof; n++) {
7359566063dSJacob Faibussowitsch         PetscCall(DMDAGetFieldName(da, n, &fieldname));
736da88d4d4SJed Brown         if (!fieldname || !fieldname[0]) {
73763a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(fieldbuf, sizeof fieldbuf, "field%" PetscInt_FMT, n));
738da88d4d4SJed Brown           fieldname = fieldbuf;
739a261c58fSBarry Smith         }
740*48a46eb9SPierre Jolivet         if (dim == 1) PetscCall(PetscFPrintf(comm, info, "#$$ Set.%s.%s = squeeze(tmp(%" PetscInt_FMT ",:))';\n", name, fieldname, n + 1));
741*48a46eb9SPierre Jolivet         if (dim == 2) PetscCall(PetscFPrintf(comm, info, "#$$ Set.%s.%s = squeeze(tmp(%" PetscInt_FMT ",:,:))';\n", name, fieldname, n + 1));
742*48a46eb9SPierre Jolivet         if (dim == 3) PetscCall(PetscFPrintf(comm, info, "#$$ Set.%s.%s = permute(squeeze(tmp(%" PetscInt_FMT ",:,:,:)),[2 1 3]);\n", name, fieldname, n + 1));
743a261c58fSBarry Smith       }
7449566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, info, "#$$ clear tmp; \n"));
7459566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, info, "#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
746a261c58fSBarry Smith     }
747a261c58fSBarry Smith 
7489566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&natural));
74947c6ae99SBarry Smith   }
75047c6ae99SBarry Smith   PetscFunctionReturn(0);
75147c6ae99SBarry Smith }
75247c6ae99SBarry Smith 
75347c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
7549371c9d4SSatish Balay PetscErrorCode VecLoad_HDF5_DA(Vec xin, PetscViewer viewer) {
7552d95be3dSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
7569a42bb27SBarry Smith   DM                da;
75715b861d2SVaclav Hapla   int               dim, rdim;
758ec7ae49cSHåkon Strandenes   hsize_t           dims[6] = {0}, count[6] = {0}, offset[6] = {0};
759d7dd068bSVaclav Hapla   PetscBool         dim2 = PETSC_FALSE, timestepping = PETSC_FALSE;
760d7dd068bSVaclav Hapla   PetscInt          dimension, timestep              = PETSC_MIN_INT, dofInd;
76147c6ae99SBarry Smith   PetscScalar      *x;
76247c6ae99SBarry Smith   const char       *vecname;
76347c6ae99SBarry Smith   hid_t             filespace; /* file dataspace identifier */
76447c6ae99SBarry Smith   hid_t             dset_id;   /* dataset identifier */
76547c6ae99SBarry Smith   hid_t             memspace;  /* memory dataspace identifier */
766bfd264e7SBarry Smith   hid_t             file_id, group;
7677bcbaff4SBarry Smith   hid_t             scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
7689c7c4993SBarry Smith   DM_DA            *dd;
76947c6ae99SBarry Smith 
77047c6ae99SBarry Smith   PetscFunctionBegin;
7717bcbaff4SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
7727bcbaff4SBarry Smith   scalartype = H5T_NATIVE_FLOAT;
7737bcbaff4SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
7747bcbaff4SBarry Smith #error "HDF5 output with 128 bit floats not supported."
775570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
776570b7f6dSBarry Smith #error "HDF5 output with 16 bit floats not supported."
7777bcbaff4SBarry Smith #else
7787bcbaff4SBarry Smith   scalartype = H5T_NATIVE_DOUBLE;
7797bcbaff4SBarry Smith #endif
7807bcbaff4SBarry Smith 
7819566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5OpenGroup(viewer, &file_id, &group));
7829566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)xin, &vecname));
7839566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckTimestepping_Internal(viewer, vecname));
7849566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5IsTimestepping(viewer, &timestepping));
785*48a46eb9SPierre Jolivet   if (timestepping) PetscCall(PetscViewerHDF5GetTimestep(viewer, &timestep));
7869566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin, &da));
7879c7c4993SBarry Smith   dd = (DM_DA *)da->data;
7889566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(da, &dimension));
78947c6ae99SBarry Smith 
790ec7ae49cSHåkon Strandenes   /* Open dataset */
791792fecdfSBarry Smith   PetscCallHDF5Return(dset_id, H5Dopen2, (group, vecname, H5P_DEFAULT));
79247c6ae99SBarry Smith 
793ec7ae49cSHåkon Strandenes   /* Retrieve the dataspace for the dataset */
794792fecdfSBarry Smith   PetscCallHDF5Return(filespace, H5Dget_space, (dset_id));
795792fecdfSBarry Smith   PetscCallHDF5Return(rdim, H5Sget_simple_extent_dims, (filespace, dims, NULL));
796ec7ae49cSHåkon Strandenes 
797ec7ae49cSHåkon Strandenes   /* Expected dimension for holding the dof's */
79847c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
799ec7ae49cSHåkon Strandenes   dofInd = rdim - 2;
800ec7ae49cSHåkon Strandenes #else
801ec7ae49cSHåkon Strandenes   dofInd = rdim - 1;
80247c6ae99SBarry Smith #endif
803ec7ae49cSHåkon Strandenes 
804ec7ae49cSHåkon Strandenes   /* The expected number of dimensions, assuming basedimension2 = false */
805ec7ae49cSHåkon Strandenes   dim = dimension;
806ec7ae49cSHåkon Strandenes   if (dd->w > 1) ++dim;
807ec7ae49cSHåkon Strandenes   if (timestep >= 0) ++dim;
80847c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
809ec7ae49cSHåkon Strandenes   ++dim;
81047c6ae99SBarry Smith #endif
811ec7ae49cSHåkon Strandenes 
812ec7ae49cSHåkon Strandenes   /* In this case the input dataset have one extra, unexpected dimension. */
813ec7ae49cSHåkon Strandenes   if (rdim == dim + 1) {
814ec7ae49cSHåkon Strandenes     /* In this case the block size unity */
815ec7ae49cSHåkon Strandenes     if (dd->w == 1 && dims[dofInd] == 1) dim2 = PETSC_TRUE;
816ec7ae49cSHåkon Strandenes 
817ec7ae49cSHåkon Strandenes     /* Special error message for the case where dof does not match the input file */
81863a3b9bcSJacob 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);
819ec7ae49cSHåkon Strandenes 
820ec7ae49cSHåkon Strandenes     /* Other cases where rdim != dim cannot be handled currently */
82163a3b9bcSJacob 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);
822ec7ae49cSHåkon Strandenes 
823ec7ae49cSHåkon Strandenes   /* Set up the hyperslab size */
824ec7ae49cSHåkon Strandenes   dim = 0;
825ec7ae49cSHåkon Strandenes   if (timestep >= 0) {
826ec7ae49cSHåkon Strandenes     offset[dim] = timestep;
827ec7ae49cSHåkon Strandenes     count[dim]  = 1;
828ec7ae49cSHåkon Strandenes     ++dim;
829ec7ae49cSHåkon Strandenes   }
830ec7ae49cSHåkon Strandenes   if (dimension == 3) {
8319566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(dd->zs, offset + dim));
8329566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(dd->ze - dd->zs, count + dim));
833ec7ae49cSHåkon Strandenes     ++dim;
834ec7ae49cSHåkon Strandenes   }
835ec7ae49cSHåkon Strandenes   if (dimension > 1) {
8369566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(dd->ys, offset + dim));
8379566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(dd->ye - dd->ys, count + dim));
838ec7ae49cSHåkon Strandenes     ++dim;
839ec7ae49cSHåkon Strandenes   }
8409566063dSJacob Faibussowitsch   PetscCall(PetscHDF5IntCast(dd->xs / dd->w, offset + dim));
8419566063dSJacob Faibussowitsch   PetscCall(PetscHDF5IntCast((dd->xe - dd->xs) / dd->w, count + dim));
842ec7ae49cSHåkon Strandenes   ++dim;
843ec7ae49cSHåkon Strandenes   if (dd->w > 1 || dim2) {
844ec7ae49cSHåkon Strandenes     offset[dim] = 0;
8459566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(dd->w, count + dim));
846ec7ae49cSHåkon Strandenes     ++dim;
847ec7ae49cSHåkon Strandenes   }
848ec7ae49cSHåkon Strandenes #if defined(PETSC_USE_COMPLEX)
849ec7ae49cSHåkon Strandenes   offset[dim] = 0;
850ec7ae49cSHåkon Strandenes   count[dim]  = 2;
851ec7ae49cSHåkon Strandenes   ++dim;
852ec7ae49cSHåkon Strandenes #endif
853ec7ae49cSHåkon Strandenes 
854ec7ae49cSHåkon Strandenes   /* Create the memory and filespace */
855792fecdfSBarry Smith   PetscCallHDF5Return(memspace, H5Screate_simple, (dim, count, NULL));
856792fecdfSBarry Smith   PetscCallHDF5(H5Sselect_hyperslab, (filespace, H5S_SELECT_SET, offset, NULL, count, NULL));
85747c6ae99SBarry Smith 
8589566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xin, &x));
859792fecdfSBarry Smith   PetscCallHDF5(H5Dread, (dset_id, scalartype, memspace, filespace, hdf5->dxpl_id, x));
8609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xin, &x));
86147c6ae99SBarry Smith 
86247c6ae99SBarry Smith   /* Close/release resources */
8639371c9d4SSatish Balay   if (group != file_id) { PetscCallHDF5(H5Gclose, (group)); }
864792fecdfSBarry Smith   PetscCallHDF5(H5Sclose, (filespace));
865792fecdfSBarry Smith   PetscCallHDF5(H5Sclose, (memspace));
866792fecdfSBarry Smith   PetscCallHDF5(H5Dclose, (dset_id));
86747c6ae99SBarry Smith   PetscFunctionReturn(0);
86847c6ae99SBarry Smith }
86947c6ae99SBarry Smith #endif
87047c6ae99SBarry Smith 
8719371c9d4SSatish Balay PetscErrorCode VecLoad_Binary_DA(Vec xin, PetscViewer viewer) {
8729a42bb27SBarry Smith   DM          da;
87347c6ae99SBarry Smith   Vec         natural;
87447c6ae99SBarry Smith   const char *prefix;
87547c6ae99SBarry Smith   PetscInt    bs;
87647c6ae99SBarry Smith   PetscBool   flag;
87747c6ae99SBarry Smith   DM_DA      *dd;
87847c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
87947c6ae99SBarry Smith   PetscBool isMPIIO;
88047c6ae99SBarry Smith #endif
88147c6ae99SBarry Smith 
88247c6ae99SBarry Smith   PetscFunctionBegin;
8839566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin, &da));
88447c6ae99SBarry Smith   dd = (DM_DA *)da->data;
88547c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
8869566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &isMPIIO));
88747c6ae99SBarry Smith   if (isMPIIO) {
8889566063dSJacob Faibussowitsch     PetscCall(DMDAArrayMPIIO(da, viewer, xin, PETSC_FALSE));
88947c6ae99SBarry Smith     PetscFunctionReturn(0);
89047c6ae99SBarry Smith   }
89147c6ae99SBarry Smith #endif
89247c6ae99SBarry Smith 
8939566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)xin, &prefix));
8949566063dSJacob Faibussowitsch   PetscCall(DMDACreateNaturalVector(da, &natural));
8959566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)natural, ((PetscObject)xin)->name));
8969566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)natural, prefix));
8979566063dSJacob Faibussowitsch   PetscCall(VecLoad(natural, viewer));
8989566063dSJacob Faibussowitsch   PetscCall(DMDANaturalToGlobalBegin(da, natural, INSERT_VALUES, xin));
8999566063dSJacob Faibussowitsch   PetscCall(DMDANaturalToGlobalEnd(da, natural, INSERT_VALUES, xin));
9009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&natural));
9019566063dSJacob Faibussowitsch   PetscCall(PetscInfo(xin, "Loading vector from natural ordering into DMDA\n"));
9029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, ((PetscObject)xin)->prefix, "-vecload_block_size", &bs, &flag));
903*48a46eb9SPierre 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));
90447c6ae99SBarry Smith   PetscFunctionReturn(0);
90547c6ae99SBarry Smith }
90647c6ae99SBarry Smith 
9079371c9d4SSatish Balay PetscErrorCode VecLoad_Default_DA(Vec xin, PetscViewer viewer) {
9089a42bb27SBarry Smith   DM        da;
90947c6ae99SBarry Smith   PetscBool isbinary;
91047c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
91147c6ae99SBarry Smith   PetscBool ishdf5;
91247c6ae99SBarry Smith #endif
91347c6ae99SBarry Smith 
91447c6ae99SBarry Smith   PetscFunctionBegin;
9159566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin, &da));
9167a8be351SBarry Smith   PetscCheck(da, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA");
91747c6ae99SBarry Smith 
91847c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
9199566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
92047c6ae99SBarry Smith #endif
9219566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
92247c6ae99SBarry Smith 
92347c6ae99SBarry Smith   if (isbinary) {
9249566063dSJacob Faibussowitsch     PetscCall(VecLoad_Binary_DA(xin, viewer));
92547c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
92647c6ae99SBarry Smith   } else if (ishdf5) {
9279566063dSJacob Faibussowitsch     PetscCall(VecLoad_HDF5_DA(xin, viewer));
92847c6ae99SBarry Smith #endif
92998921bdaSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Viewer type %s not supported for vector loading", ((PetscObject)viewer)->type_name);
93047c6ae99SBarry Smith   PetscFunctionReturn(0);
93147c6ae99SBarry Smith }
932