xref: /petsc/src/dm/impls/da/gr2.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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 */
2847c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA2d_Zoom(PetscDraw draw,void *ctx)
2947c6ae99SBarry Smith {
3047c6ae99SBarry Smith   ZoomCtx           *zctx = (ZoomCtx*)ctx;
31e5ab1681SLisandro Dalcin   PetscInt          m,n,i,j,k,dof,id,c1,c2,c3,c4;
32e5ab1681SLisandro Dalcin   PetscReal         min,max,x1,x2,x3,x4,y_1,y2,y3,y4;
33e5ab1681SLisandro Dalcin   const PetscScalar *xy,*v;
345f80ce2aSJacob Faibussowitsch   PetscErrorCode    ierr;
3547c6ae99SBarry Smith 
3647c6ae99SBarry Smith   PetscFunctionBegin;
3747c6ae99SBarry Smith   m    = zctx->m;
3847c6ae99SBarry Smith   n    = zctx->n;
39e5ab1681SLisandro Dalcin   dof  = zctx->dof;
4047c6ae99SBarry Smith   k    = zctx->k;
4147c6ae99SBarry Smith   xy   = zctx->xy;
42e5ab1681SLisandro Dalcin   v    = zctx->v;
4347c6ae99SBarry Smith   min  = zctx->min;
44f3f0eb19SBarry Smith   max  = zctx->max;
4547c6ae99SBarry Smith 
4647c6ae99SBarry Smith   /* PetscDraw the contour plot patch */
479566063dSJacob Faibussowitsch   ierr = PetscDrawCollectiveBegin(draw);PetscCall(ierr);
4847c6ae99SBarry Smith   for (j=0; j<n-1; j++) {
4947c6ae99SBarry Smith     for (i=0; i<m-1; i++) {
500e4fe250SBarry Smith       id   = i+j*m;
510e4fe250SBarry Smith       x1   = PetscRealPart(xy[2*id]);
520e4fe250SBarry Smith       y_1  = PetscRealPart(xy[2*id+1]);
53e5ab1681SLisandro Dalcin       c1   = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max);
540e4fe250SBarry Smith 
550e4fe250SBarry Smith       id   = i+j*m+1;
560e4fe250SBarry Smith       x2   = PetscRealPart(xy[2*id]);
57a39a4c7dSLisandro Dalcin       y2   = PetscRealPart(xy[2*id+1]);
58e5ab1681SLisandro Dalcin       c2   = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max);
590e4fe250SBarry Smith 
600e4fe250SBarry Smith       id   = i+j*m+1+m;
61a39a4c7dSLisandro Dalcin       x3   = PetscRealPart(xy[2*id]);
620e4fe250SBarry Smith       y3   = PetscRealPart(xy[2*id+1]);
63e5ab1681SLisandro Dalcin       c3   = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max);
640e4fe250SBarry Smith 
650e4fe250SBarry Smith       id   = i+j*m+m;
66a39a4c7dSLisandro Dalcin       x4   = PetscRealPart(xy[2*id]);
67a39a4c7dSLisandro Dalcin       y4   = PetscRealPart(xy[2*id+1]);
68e5ab1681SLisandro Dalcin       c4   = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max);
69f3f0eb19SBarry Smith 
709566063dSJacob Faibussowitsch       PetscCall(PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3));
719566063dSJacob Faibussowitsch       PetscCall(PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4));
7247c6ae99SBarry Smith       if (zctx->showgrid) {
739566063dSJacob Faibussowitsch         PetscCall(PetscDrawLine(draw,x1,y_1,x2,y2,PETSC_DRAW_BLACK));
749566063dSJacob Faibussowitsch         PetscCall(PetscDrawLine(draw,x2,y2,x3,y3,PETSC_DRAW_BLACK));
759566063dSJacob Faibussowitsch         PetscCall(PetscDrawLine(draw,x3,y3,x4,y4,PETSC_DRAW_BLACK));
769566063dSJacob Faibussowitsch         PetscCall(PetscDrawLine(draw,x4,y4,x1,y_1,PETSC_DRAW_BLACK));
7747c6ae99SBarry Smith       }
7847c6ae99SBarry Smith     }
7947c6ae99SBarry Smith   }
807fe7c8feSLisandro Dalcin   if (zctx->showaxis && !zctx->rank) {
81e5ab1681SLisandro Dalcin     if (zctx->name0 || zctx->name1) {
82109c9344SBarry Smith       PetscReal xl,yl,xr,yr,x,y;
839566063dSJacob Faibussowitsch       PetscCall(PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr));
84e5ab1681SLisandro Dalcin       x  = xl + .30*(xr - xl);
85109c9344SBarry Smith       xl = xl + .01*(xr - xl);
86e5ab1681SLisandro Dalcin       y  = yr - .30*(yr - yl);
87109c9344SBarry Smith       yl = yl + .01*(yr - yl);
889566063dSJacob Faibussowitsch       if (zctx->name0) PetscCall(PetscDrawString(draw,x,yl,PETSC_DRAW_BLACK,zctx->name0));
899566063dSJacob Faibussowitsch       if (zctx->name1) PetscCall(PetscDrawStringVertical(draw,xl,y,PETSC_DRAW_BLACK,zctx->name1));
90109c9344SBarry Smith     }
910e4fe250SBarry Smith     /*
920e4fe250SBarry Smith        Ideally we would use the PetscDrawAxis object to manage displaying the coordinate limits
930e4fe250SBarry Smith        but that may require some refactoring.
940e4fe250SBarry Smith     */
95e5ab1681SLisandro Dalcin     {
967fe7c8feSLisandro Dalcin       double xmin = (double)zctx->xmin, ymin = (double)zctx->ymin;
977fe7c8feSLisandro Dalcin       double xmax = (double)zctx->xmax, ymax = (double)zctx->ymax;
98e5ab1681SLisandro Dalcin       char   value[16]; size_t len; 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   }
1119566063dSJacob Faibussowitsch   ierr = PetscDrawCollectiveEnd(draw);PetscCall(ierr);
11247c6ae99SBarry Smith   PetscFunctionReturn(0);
11347c6ae99SBarry Smith }
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA2d(Vec xin,PetscViewer viewer)
11647c6ae99SBarry Smith {
1179a42bb27SBarry Smith   DM                 da,dac,dag;
118a74ba6f7SBarry Smith   PetscInt           N,s,M,w,ncoors = 4;
11947c6ae99SBarry Smith   const PetscInt     *lx,*ly;
120e5ab1681SLisandro Dalcin   PetscReal          coors[4];
12147c6ae99SBarry Smith   PetscDraw          draw,popup;
12247c6ae99SBarry Smith   PetscBool          isnull,useports = PETSC_FALSE;
12347c6ae99SBarry Smith   MPI_Comm           comm;
12447c6ae99SBarry Smith   Vec                xlocal,xcoor,xcoorl;
125bff4a2f0SMatthew G. Knepley   DMBoundaryType     bx,by;
126aa219208SBarry Smith   DMDAStencilType    st;
12747c6ae99SBarry Smith   ZoomCtx            zctx;
1280298fd71SBarry Smith   PetscDrawViewPorts *ports = NULL;
12947c6ae99SBarry Smith   PetscViewerFormat  format;
13020d0051dSBarry Smith   PetscInt           *displayfields;
13167dd0837SBarry Smith   PetscInt           ndisplayfields,i,nbounds;
13267dd0837SBarry Smith   const PetscReal    *bounds;
13347c6ae99SBarry Smith 
13447c6ae99SBarry Smith   PetscFunctionBegin;
13547c6ae99SBarry Smith   zctx.showgrid = PETSC_FALSE;
1367fe7c8feSLisandro Dalcin   zctx.showaxis = PETSC_TRUE;
1378865f1eaSKarl Rupp 
1389566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw));
1399566063dSJacob Faibussowitsch   PetscCall(PetscDrawIsNull(draw,&isnull));
1405b399a63SLisandro Dalcin   if (isnull) PetscFunctionReturn(0);
1415b399a63SLisandro Dalcin 
1429566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetBounds(viewer,&nbounds,&bounds));
14347c6ae99SBarry Smith 
1449566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin,&da));
1457a8be351SBarry Smith   PetscCheck(da,PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
14647c6ae99SBarry Smith 
1479566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)xin,&comm));
1489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&zctx.rank));
14947c6ae99SBarry Smith 
1509566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da,NULL,&M,&N,NULL,&zctx.m,&zctx.n,NULL,&w,&s,&bx,&by,NULL,&st));
1519566063dSJacob Faibussowitsch   PetscCall(DMDAGetOwnershipRanges(da,&lx,&ly,NULL));
15247c6ae99SBarry Smith 
15347c6ae99SBarry Smith   /*
15447c6ae99SBarry Smith      Obtain a sequential vector that is going to contain the local values plus ONE layer of
155aa219208SBarry Smith      ghosted values to draw the graphics from. We also need its corresponding DMDA (dac) that will
15647c6ae99SBarry Smith      update the local values pluse ONE layer of ghost values.
15747c6ae99SBarry Smith   */
1589566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)da,"GraphicsGhosted",(PetscObject*)&xlocal));
15947c6ae99SBarry Smith   if (!xlocal) {
160bff4a2f0SMatthew G. Knepley     if (bx !=  DM_BOUNDARY_NONE || by !=  DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) {
16147c6ae99SBarry Smith       /*
16247c6ae99SBarry Smith          if original da is not of stencil width one, or periodic or not a box stencil then
163aa219208SBarry Smith          create a special DMDA to handle one level of ghost points for graphics
16447c6ae99SBarry Smith       */
1659566063dSJacob Faibussowitsch       PetscCall(DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,zctx.m,zctx.n,w,1,lx,ly,&dac));
1669566063dSJacob Faibussowitsch       PetscCall(DMSetUp(dac));
1679566063dSJacob Faibussowitsch       PetscCall(PetscInfo(da,"Creating auxilary DMDA for managing graphics ghost points\n"));
16847c6ae99SBarry Smith     } else {
16947c6ae99SBarry Smith       /* otherwise we can use the da we already have */
17047c6ae99SBarry Smith       dac = da;
17147c6ae99SBarry Smith     }
17247c6ae99SBarry Smith     /* create local vector for holding ghosted values used in graphics */
1739566063dSJacob Faibussowitsch     PetscCall(DMCreateLocalVector(dac,&xlocal));
17447c6ae99SBarry Smith     if (dac != da) {
175aa219208SBarry Smith       /* don't keep any public reference of this DMDA, is is only available through xlocal */
1769566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)dac));
17747c6ae99SBarry Smith     } else {
17847c6ae99SBarry Smith       /* remove association between xlocal and da, because below we compose in the opposite
17947c6ae99SBarry Smith          direction and if we left this connect we'd get a loop, so the objects could
18047c6ae99SBarry Smith          never be destroyed */
1819566063dSJacob Faibussowitsch       PetscCall(PetscObjectRemoveReference((PetscObject)xlocal,"__PETSc_dm"));
18247c6ae99SBarry Smith     }
1839566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)da,"GraphicsGhosted",(PetscObject)xlocal));
1849566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)xlocal));
18547c6ae99SBarry Smith   } else {
186bff4a2f0SMatthew G. Knepley     if (bx !=  DM_BOUNDARY_NONE || by !=  DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) {
1879566063dSJacob Faibussowitsch       PetscCall(VecGetDM(xlocal, &dac));
188f7923d8aSBarry Smith     } else {
189f7923d8aSBarry Smith       dac = da;
19047c6ae99SBarry Smith     }
19147c6ae99SBarry Smith   }
19247c6ae99SBarry Smith 
19347c6ae99SBarry Smith   /*
19447c6ae99SBarry Smith       Get local (ghosted) values of vector
19547c6ae99SBarry Smith   */
1969566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dac,xin,INSERT_VALUES,xlocal));
1979566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dac,xin,INSERT_VALUES,xlocal));
1989566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xlocal,&zctx.v));
19947c6ae99SBarry Smith 
200832b7cebSLisandro Dalcin   /*
201832b7cebSLisandro Dalcin       Get coordinates of nodes
202832b7cebSLisandro Dalcin   */
2039566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da,&xcoor));
20447c6ae99SBarry Smith   if (!xcoor) {
2059566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,0.0));
2069566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(da,&xcoor));
20747c6ae99SBarry Smith   }
20847c6ae99SBarry Smith 
20947c6ae99SBarry Smith   /*
21047c6ae99SBarry Smith       Determine the min and max coordinates in plot
21147c6ae99SBarry Smith   */
2129566063dSJacob Faibussowitsch   PetscCall(VecStrideMin(xcoor,0,NULL,&zctx.xmin));
2139566063dSJacob Faibussowitsch   PetscCall(VecStrideMax(xcoor,0,NULL,&zctx.xmax));
2149566063dSJacob Faibussowitsch   PetscCall(VecStrideMin(xcoor,1,NULL,&zctx.ymin));
2159566063dSJacob Faibussowitsch   PetscCall(VecStrideMax(xcoor,1,NULL,&zctx.ymax));
2169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-draw_contour_axis",&zctx.showaxis,NULL));
2177fe7c8feSLisandro Dalcin   if (zctx.showaxis) {
2187fe7c8feSLisandro Dalcin     coors[0] = zctx.xmin - .05*(zctx.xmax - zctx.xmin); coors[1] = zctx.ymin - .05*(zctx.ymax - zctx.ymin);
2197fe7c8feSLisandro Dalcin     coors[2] = zctx.xmax + .05*(zctx.xmax - zctx.xmin); coors[3] = zctx.ymax + .05*(zctx.ymax - zctx.ymin);
2207fe7c8feSLisandro Dalcin   } else {
2217fe7c8feSLisandro Dalcin     coors[0] = zctx.xmin; coors[1] = zctx.ymin; coors[2] = zctx.xmax; coors[3] = zctx.ymax;
2227fe7c8feSLisandro Dalcin   }
2239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetRealArray(NULL,NULL,"-draw_coordinates",coors,&ncoors,NULL));
2249566063dSJacob 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]));
225a74ba6f7SBarry Smith 
22647c6ae99SBarry Smith   /*
227832b7cebSLisandro Dalcin       Get local ghosted version of coordinates
22847c6ae99SBarry Smith   */
2299566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject*)&xcoorl));
23047c6ae99SBarry Smith   if (!xcoorl) {
231aa219208SBarry Smith     /* create DMDA to get local version of graphics */
2329566063dSJacob Faibussowitsch     PetscCall(DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,zctx.m,zctx.n,2,1,lx,ly,&dag));
2339566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dag));
2349566063dSJacob Faibussowitsch     PetscCall(PetscInfo(dag,"Creating auxilary DMDA for managing graphics coordinates ghost points\n"));
2359566063dSJacob Faibussowitsch     PetscCall(DMCreateLocalVector(dag,&xcoorl));
2369566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject)xcoorl));
2379566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)dag));
2389566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)xcoorl));
23947c6ae99SBarry Smith   } else {
2409566063dSJacob Faibussowitsch     PetscCall(VecGetDM(xcoorl,&dag));
24147c6ae99SBarry Smith   }
2429566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dag,xcoor,INSERT_VALUES,xcoorl));
2439566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dag,xcoor,INSERT_VALUES,xcoorl));
2449566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xcoorl,&zctx.xy));
2459566063dSJacob Faibussowitsch   PetscCall(DMDAGetCoordinateName(da,0,&zctx.name0));
2469566063dSJacob Faibussowitsch   PetscCall(DMDAGetCoordinateName(da,1,&zctx.name1));
24747c6ae99SBarry Smith 
24847c6ae99SBarry Smith   /*
24947c6ae99SBarry Smith       Get information about size of area each processor must do graphics for
25047c6ae99SBarry Smith   */
2519566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac,NULL,&M,&N,NULL,NULL,NULL,NULL,&zctx.dof,NULL,&bx,&by,NULL,NULL));
2529566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dac,NULL,NULL,NULL,&zctx.m,&zctx.n,NULL));
2539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-draw_contour_grid",&zctx.showgrid,NULL));
2544e6118eeSBarry Smith 
2559566063dSJacob Faibussowitsch   PetscCall(DMDASelectFields(da,&ndisplayfields,&displayfields));
2569566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer,&format));
2579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-draw_ports",&useports,NULL));
258832b7cebSLisandro Dalcin   if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE;
259832b7cebSLisandro Dalcin   if (useports) {
2609566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw));
2619566063dSJacob Faibussowitsch     PetscCall(PetscDrawCheckResizedWindow(draw));
2629566063dSJacob Faibussowitsch     PetscCall(PetscDrawClear(draw));
2639566063dSJacob Faibussowitsch     PetscCall(PetscDrawViewPortsCreate(draw,ndisplayfields,&ports));
264e5ab1681SLisandro Dalcin   }
26520d0051dSBarry Smith 
266832b7cebSLisandro Dalcin   /*
267832b7cebSLisandro Dalcin       Loop over each field; drawing each in a different window
268832b7cebSLisandro Dalcin   */
26920d0051dSBarry Smith   for (i=0; i<ndisplayfields; i++) {
27020d0051dSBarry Smith     zctx.k = displayfields[i];
2715b399a63SLisandro Dalcin 
272e5ab1681SLisandro Dalcin     /* determine the min and max value in plot */
2739566063dSJacob Faibussowitsch     PetscCall(VecStrideMin(xin,zctx.k,NULL,&zctx.min));
2749566063dSJacob Faibussowitsch     PetscCall(VecStrideMax(xin,zctx.k,NULL,&zctx.max));
27567dd0837SBarry Smith     if (zctx.k < nbounds) {
276f3f0eb19SBarry Smith       zctx.min = bounds[2*zctx.k];
277f3f0eb19SBarry Smith       zctx.max = bounds[2*zctx.k+1];
27867dd0837SBarry Smith     }
27947c6ae99SBarry Smith     if (zctx.min == zctx.max) {
28047c6ae99SBarry Smith       zctx.min -= 1.e-12;
28147c6ae99SBarry Smith       zctx.max += 1.e-12;
28247c6ae99SBarry Smith     }
2839566063dSJacob Faibussowitsch     PetscCall(PetscInfo(da,"DMDA 2d contour plot min %g max %g\n",(double)zctx.min,(double)zctx.max));
28447c6ae99SBarry Smith 
285832b7cebSLisandro Dalcin     if (useports) {
2869566063dSJacob Faibussowitsch       PetscCall(PetscDrawViewPortsSet(ports,i));
287832b7cebSLisandro Dalcin     } else {
288832b7cebSLisandro Dalcin       const char *title;
2899566063dSJacob Faibussowitsch       PetscCall(PetscViewerDrawGetDraw(viewer,i,&draw));
2909566063dSJacob Faibussowitsch       PetscCall(DMDAGetFieldName(da,zctx.k,&title));
2919566063dSJacob Faibussowitsch       if (title) PetscCall(PetscDrawSetTitle(draw,title));
292832b7cebSLisandro Dalcin     }
293832b7cebSLisandro Dalcin 
2949566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetPopup(draw,&popup));
2959566063dSJacob Faibussowitsch     PetscCall(PetscDrawScalePopup(popup,zctx.min,zctx.max));
2969566063dSJacob Faibussowitsch     PetscCall(PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]));
2979566063dSJacob Faibussowitsch     PetscCall(PetscDrawZoom(draw,VecView_MPI_Draw_DA2d_Zoom,&zctx));
2989566063dSJacob Faibussowitsch     if (!useports) PetscCall(PetscDrawSave(draw));
29947c6ae99SBarry Smith   }
300832b7cebSLisandro Dalcin   if (useports) {
3019566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw));
3029566063dSJacob Faibussowitsch     PetscCall(PetscDrawSave(draw));
303832b7cebSLisandro Dalcin   }
304832b7cebSLisandro Dalcin 
3059566063dSJacob Faibussowitsch   PetscCall(PetscDrawViewPortsDestroy(ports));
3069566063dSJacob Faibussowitsch   PetscCall(PetscFree(displayfields));
3079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xcoorl,&zctx.xy));
3089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xlocal,&zctx.v));
30947c6ae99SBarry Smith   PetscFunctionReturn(0);
31047c6ae99SBarry Smith }
31147c6ae99SBarry Smith 
3120da24e51SJuha Jäykkä #if defined(PETSC_HAVE_HDF5)
313c73cfb54SMatthew G. Knepley static PetscErrorCode VecGetHDF5ChunkSize(DM_DA *da, Vec xin, PetscInt dimension, PetscInt timestep, hsize_t *chunkDims)
3140da24e51SJuha Jäykkä {
315d4190030SJed Brown   PetscMPIInt    comm_size;
316561a82e7SJed Brown   hsize_t        chunk_size, target_size, dim;
317561a82e7SJed Brown   hsize_t        vec_size = sizeof(PetscScalar)*da->M*da->N*da->P*da->w;
318561a82e7SJed Brown   hsize_t        avg_local_vec_size,KiB = 1024,MiB = KiB*KiB,GiB = MiB*KiB,min_size = MiB;
319561a82e7SJed Brown   hsize_t        max_chunks = 64*KiB;                                              /* HDF5 internal limitation */
320561a82e7SJed Brown   hsize_t        max_chunk_size = 4*GiB;                                           /* HDF5 internal limitation */
32184179ffaSJed Brown   PetscInt       zslices=da->p, yslices=da->n, xslices=da->m;
3220da24e51SJuha Jäykkä 
323cb5c4f94SJuha Jäykkä   PetscFunctionBegin;
3249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)xin), &comm_size));
325faa75363SBarry Smith   avg_local_vec_size = (hsize_t) PetscCeilInt(vec_size,comm_size);      /* we will attempt to use this as the chunk size */
3260da24e51SJuha Jäykkä 
327faa75363SBarry 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))));
3287d310018SBarry 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 */
3297d310018SBarry 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);
3300da24e51SJuha Jäykkä 
331cb5c4f94SJuha Jäykkä   /*
332cb5c4f94SJuha Jäykkä    if size/rank > max_chunk_size, we need radical measures: even going down to
333cb5c4f94SJuha Jäykkä    avg_local_vec_size is not enough, so we simply use chunk size of 4 GiB no matter
334cb5c4f94SJuha Jäykkä    what, composed in the most efficient way possible.
335cb5c4f94SJuha Jäykkä    N.B. this minimises the number of chunks, which may or may not be the optimal
336cb5c4f94SJuha Jäykkä    solution. In a BG, for example, the optimal solution is probably to make # chunks = #
337cb5c4f94SJuha Jäykkä    IO nodes involved, but this author has no access to a BG to figure out how to
338cb5c4f94SJuha Jäykkä    reliably find the right number. And even then it may or may not be enough.
339cb5c4f94SJuha Jäykkä    */
3400da24e51SJuha Jäykkä   if (avg_local_vec_size > max_chunk_size) {
341cb5c4f94SJuha Jäykkä     /* check if we can just split local z-axis: is that enough? */
342faa75363SBarry Smith     zslices = PetscCeilInt(vec_size,da->p*max_chunk_size)*zslices;
3430da24e51SJuha Jäykkä     if (zslices > da->P) {
344cb5c4f94SJuha Jäykkä       /* lattice is too large in xy-directions, splitting z only is not enough */
3450da24e51SJuha Jäykkä       zslices = da->P;
346faa75363SBarry Smith       yslices = PetscCeilInt(vec_size,zslices*da->n*max_chunk_size)*yslices;
3470da24e51SJuha Jäykkä       if (yslices > da->N) {
348cb5c4f94SJuha Jäykkä         /* lattice is too large in x-direction, splitting along z, y is not enough */
3490da24e51SJuha Jäykkä         yslices = da->N;
350faa75363SBarry Smith         xslices = PetscCeilInt(vec_size,zslices*yslices*da->m*max_chunk_size)*xslices;
3510da24e51SJuha Jäykkä       }
3520da24e51SJuha Jäykkä     }
3530da24e51SJuha Jäykkä     dim = 0;
3540da24e51SJuha Jäykkä     if (timestep >= 0) {
3550da24e51SJuha Jäykkä       ++dim;
3560da24e51SJuha Jäykkä     }
357cb5c4f94SJuha Jäykkä     /* prefer to split z-axis, even down to planar slices */
358c73cfb54SMatthew G. Knepley     if (dimension == 3) {
359cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->P/zslices;
360cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->N/yslices;
361cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->M/xslices;
3620da24e51SJuha Jäykkä     } else {
363cb5c4f94SJuha Jäykkä       /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */
364cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->N/yslices;
365cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->M/xslices;
3660da24e51SJuha Jäykkä     }
3670da24e51SJuha 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);
3680da24e51SJuha Jäykkä   } else {
3690da24e51SJuha Jäykkä     if (target_size < chunk_size) {
370cb5c4f94SJuha Jäykkä       /* only change the defaults if target_size < chunk_size */
3710da24e51SJuha Jäykkä       dim = 0;
3720da24e51SJuha Jäykkä       if (timestep >= 0) {
3730da24e51SJuha Jäykkä         ++dim;
3740da24e51SJuha Jäykkä       }
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ä       }
4010da24e51SJuha 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);
4020da24e51SJuha Jäykkä     } else {
403cb5c4f94SJuha Jäykkä       /* precomputed chunks are fine, we don't need to do anything */
4040da24e51SJuha Jäykkä     }
4050da24e51SJuha Jäykkä   }
4060da24e51SJuha Jäykkä   PetscFunctionReturn(0);
4070da24e51SJuha Jäykkä }
4080da24e51SJuha Jäykkä #endif
40947c6ae99SBarry Smith 
41047c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
41147c6ae99SBarry Smith PetscErrorCode VecView_MPI_HDF5_DA(Vec xin,PetscViewer viewer)
41247c6ae99SBarry Smith {
4132d95be3dSVaclav Hapla   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*) viewer->data;
4149b2a5a86SJed Brown   DM                dm;
4159b2a5a86SJed Brown   DM_DA             *da;
41647c6ae99SBarry Smith   hid_t             filespace;  /* file dataspace identifier */
4178e2ae6d7SMichael Kraus   hid_t             chunkspace; /* chunk dataset property identifier */
41847c6ae99SBarry Smith   hid_t             dset_id;    /* dataset identifier */
41947c6ae99SBarry Smith   hid_t             memspace;   /* memory dataspace identifier */
42047c6ae99SBarry Smith   hid_t             file_id;
42115214e8eSMatthew G Knepley   hid_t             group;
4229a0502c6SHåkon Strandenes   hid_t             memscalartype; /* scalar type for mem (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
4239a0502c6SHåkon Strandenes   hid_t             filescalartype; /* scalar type for file (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
424d9a4edebSJed Brown   hsize_t           dim;
4250da24e51SJuha 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  */
426d7dd068bSVaclav Hapla   PetscBool         timestepping=PETSC_FALSE, dim2, spoutput;
427d7dd068bSVaclav Hapla   PetscInt          timestep=PETSC_MIN_INT, dimension;
4285edff71fSBarry Smith   const PetscScalar *x;
4298e2ae6d7SMichael Kraus   const char        *vecname;
43047c6ae99SBarry Smith 
43147c6ae99SBarry Smith   PetscFunctionBegin;
4329566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5OpenGroup(viewer, &file_id, &group));
4339566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5IsTimestepping(viewer, &timestepping));
434d7dd068bSVaclav Hapla   if (timestepping) {
4359566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5GetTimestep(viewer, &timestep));
436d7dd068bSVaclav Hapla   }
4379566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetBaseDimension2(viewer,&dim2));
4389566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetSPOutput(viewer,&spoutput));
43915214e8eSMatthew G Knepley 
4409566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin,&dm));
4417a8be351SBarry Smith   PetscCheck(dm,PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
4429b2a5a86SJed Brown   da = (DM_DA*)dm->data;
4439566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dimension));
44447c6ae99SBarry Smith 
4458e2ae6d7SMichael Kraus   /* Create the dataspace for the dataset.
4468e2ae6d7SMichael Kraus    *
4478e2ae6d7SMichael Kraus    * dims - holds the current dimensions of the dataset
4488e2ae6d7SMichael Kraus    *
4498e2ae6d7SMichael Kraus    * maxDims - holds the maximum dimensions of the dataset (unlimited
4508e2ae6d7SMichael Kraus    * for the number of time steps with the current dimensions for the
4518e2ae6d7SMichael Kraus    * other dimensions; so only additional time steps can be added).
4528e2ae6d7SMichael Kraus    *
4538e2ae6d7SMichael Kraus    * chunkDims - holds the size of a single time step (required to
4548e2ae6d7SMichael Kraus    * permit extending dataset).
4558e2ae6d7SMichael Kraus    */
4568e2ae6d7SMichael Kraus   dim = 0;
4578e2ae6d7SMichael Kraus   if (timestep >= 0) {
4588e2ae6d7SMichael Kraus     dims[dim]      = timestep+1;
4598e2ae6d7SMichael Kraus     maxDims[dim]   = H5S_UNLIMITED;
4608e2ae6d7SMichael Kraus     chunkDims[dim] = 1;
4618e2ae6d7SMichael Kraus     ++dim;
4628e2ae6d7SMichael Kraus   }
463c73cfb54SMatthew G. Knepley   if (dimension == 3) {
4649566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(da->P,dims+dim));
4658e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4668e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4678e2ae6d7SMichael Kraus     ++dim;
4688e2ae6d7SMichael Kraus   }
469c73cfb54SMatthew G. Knepley   if (dimension > 1) {
4709566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(da->N,dims+dim));
4718e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4728e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4738e2ae6d7SMichael Kraus     ++dim;
4748e2ae6d7SMichael Kraus   }
4759566063dSJacob Faibussowitsch   PetscCall(PetscHDF5IntCast(da->M,dims+dim));
4768e2ae6d7SMichael Kraus   maxDims[dim]   = dims[dim];
4778e2ae6d7SMichael Kraus   chunkDims[dim] = dims[dim];
4788e2ae6d7SMichael Kraus   ++dim;
47982ea9b62SBarry Smith   if (da->w > 1 || dim2) {
4809566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(da->w,dims+dim));
4818e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4828e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4838e2ae6d7SMichael Kraus     ++dim;
4848e2ae6d7SMichael Kraus   }
48547c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
4868e2ae6d7SMichael Kraus   dims[dim]      = 2;
4878e2ae6d7SMichael Kraus   maxDims[dim]   = dims[dim];
4888e2ae6d7SMichael Kraus   chunkDims[dim] = dims[dim];
4898e2ae6d7SMichael Kraus   ++dim;
49047c6ae99SBarry Smith #endif
4910da24e51SJuha Jäykkä 
4929566063dSJacob Faibussowitsch   PetscCall(VecGetHDF5ChunkSize(da, xin, dimension, timestep, chunkDims));
4930da24e51SJuha Jäykkä 
494729ab6d9SBarry Smith   PetscStackCallHDF5Return(filespace,H5Screate_simple,(dim, dims, maxDims));
49547c6ae99SBarry Smith 
49615214e8eSMatthew G Knepley #if defined(PETSC_USE_REAL_SINGLE)
4979a0502c6SHåkon Strandenes   memscalartype = H5T_NATIVE_FLOAT;
4989a0502c6SHåkon Strandenes   filescalartype = H5T_NATIVE_FLOAT;
49915214e8eSMatthew G Knepley #elif defined(PETSC_USE_REAL___FLOAT128)
50015214e8eSMatthew G Knepley #error "HDF5 output with 128 bit floats not supported."
501570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
502570b7f6dSBarry Smith #error "HDF5 output with 16 bit floats not supported."
50315214e8eSMatthew G Knepley #else
5049a0502c6SHåkon Strandenes   memscalartype = H5T_NATIVE_DOUBLE;
5059a0502c6SHåkon Strandenes   if (spoutput == PETSC_TRUE) filescalartype = H5T_NATIVE_FLOAT;
5069a0502c6SHåkon Strandenes   else filescalartype = H5T_NATIVE_DOUBLE;
50715214e8eSMatthew G Knepley #endif
50815214e8eSMatthew G Knepley 
50947c6ae99SBarry Smith   /* Create the dataset with default properties and close filespace */
5109566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)xin,&vecname));
51115214e8eSMatthew G Knepley   if (!H5Lexists(group, vecname, H5P_DEFAULT)) {
5128e2ae6d7SMichael Kraus     /* Create chunk */
513729ab6d9SBarry Smith     PetscStackCallHDF5Return(chunkspace,H5Pcreate,(H5P_DATASET_CREATE));
514729ab6d9SBarry Smith     PetscStackCallHDF5(H5Pset_chunk,(chunkspace, dim, chunkDims));
5158e2ae6d7SMichael Kraus 
5169a0502c6SHåkon Strandenes     PetscStackCallHDF5Return(dset_id,H5Dcreate2,(group, vecname, filescalartype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT));
51715214e8eSMatthew G Knepley   } else {
518729ab6d9SBarry Smith     PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT));
519729ab6d9SBarry Smith     PetscStackCallHDF5(H5Dset_extent,(dset_id, dims));
52015214e8eSMatthew G Knepley   }
521729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
52247c6ae99SBarry Smith 
52347c6ae99SBarry Smith   /* Each process defines a dataset and writes it to the hyperslab in the file */
5248e2ae6d7SMichael Kraus   dim = 0;
5258e2ae6d7SMichael Kraus   if (timestep >= 0) {
5268e2ae6d7SMichael Kraus     offset[dim] = timestep;
5278e2ae6d7SMichael Kraus     ++dim;
5288e2ae6d7SMichael Kraus   }
5299566063dSJacob Faibussowitsch   if (dimension == 3) PetscCall(PetscHDF5IntCast(da->zs,offset + dim++));
5309566063dSJacob Faibussowitsch   if (dimension > 1)  PetscCall(PetscHDF5IntCast(da->ys,offset + dim++));
5319566063dSJacob Faibussowitsch   PetscCall(PetscHDF5IntCast(da->xs/da->w,offset + dim++));
53282ea9b62SBarry Smith   if (da->w > 1 || dim2) offset[dim++] = 0;
53347c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5348e2ae6d7SMichael Kraus   offset[dim++] = 0;
53547c6ae99SBarry Smith #endif
5368e2ae6d7SMichael Kraus   dim = 0;
5378e2ae6d7SMichael Kraus   if (timestep >= 0) {
5388e2ae6d7SMichael Kraus     count[dim] = 1;
5398e2ae6d7SMichael Kraus     ++dim;
5408e2ae6d7SMichael Kraus   }
5419566063dSJacob Faibussowitsch   if (dimension == 3) PetscCall(PetscHDF5IntCast(da->ze - da->zs,count + dim++));
5429566063dSJacob Faibussowitsch   if (dimension > 1)  PetscCall(PetscHDF5IntCast(da->ye - da->ys,count + dim++));
5439566063dSJacob Faibussowitsch   PetscCall(PetscHDF5IntCast((da->xe - da->xs)/da->w,count + dim++));
5449566063dSJacob Faibussowitsch   if (da->w > 1 || dim2) PetscCall(PetscHDF5IntCast(da->w,count + dim++));
54547c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5468e2ae6d7SMichael Kraus   count[dim++] = 2;
54747c6ae99SBarry Smith #endif
548729ab6d9SBarry Smith   PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL));
549729ab6d9SBarry Smith   PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id));
550729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL));
55147c6ae99SBarry Smith 
5529566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xin, &x));
5532d95be3dSVaclav Hapla   PetscStackCallHDF5(H5Dwrite,(dset_id, memscalartype, memspace, filespace, hdf5->dxpl_id, x));
554729ab6d9SBarry Smith   PetscStackCallHDF5(H5Fflush,(file_id, H5F_SCOPE_GLOBAL));
5559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xin, &x));
55647c6ae99SBarry Smith 
557e0552f36Ssajid__ali   #if defined(PETSC_USE_COMPLEX)
558e0552f36Ssajid__ali   {
559e0552f36Ssajid__ali     PetscBool tru = PETSC_TRUE;
5609566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer,(PetscObject)xin,"complex",PETSC_BOOL,&tru));
561e0552f36Ssajid__ali   }
562e0552f36Ssajid__ali   #endif
563c7a34aceSVaclav Hapla   if (timestepping) {
5649566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer,(PetscObject)xin,"timestepping",PETSC_BOOL,&timestepping));
565c7a34aceSVaclav Hapla   }
566e0552f36Ssajid__ali 
56747c6ae99SBarry Smith   /* Close/release resources */
56815214e8eSMatthew G Knepley   if (group != file_id) {
569729ab6d9SBarry Smith     PetscStackCallHDF5(H5Gclose,(group));
57015214e8eSMatthew G Knepley   }
571729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
572729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(memspace));
573729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dset_id));
5749566063dSJacob Faibussowitsch   PetscCall(PetscInfo(xin,"Wrote Vec object with name %s\n",vecname));
57547c6ae99SBarry Smith   PetscFunctionReturn(0);
57647c6ae99SBarry Smith }
57747c6ae99SBarry Smith #endif
57847c6ae99SBarry Smith 
57909573ac7SBarry Smith extern PetscErrorCode VecView_MPI_Draw_DA1d(Vec,PetscViewer);
58047c6ae99SBarry Smith 
58147c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
582aa219208SBarry Smith static PetscErrorCode DMDAArrayMPIIO(DM da,PetscViewer viewer,Vec xin,PetscBool write)
58347c6ae99SBarry Smith {
58447c6ae99SBarry Smith   MPI_File          mfdes;
58547c6ae99SBarry Smith   PetscMPIInt       gsizes[4],lsizes[4],lstarts[4],asiz,dof;
58647c6ae99SBarry Smith   MPI_Datatype      view;
58747c6ae99SBarry Smith   const PetscScalar *array;
58847c6ae99SBarry Smith   MPI_Offset        off;
58947c6ae99SBarry Smith   MPI_Aint          ub,ul;
59047c6ae99SBarry Smith   PetscInt          type,rows,vecrows,tr[2];
59147c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
5920c169764Sdmay   PetscBool         skipheader;
59347c6ae99SBarry Smith 
59447c6ae99SBarry Smith   PetscFunctionBegin;
5959566063dSJacob Faibussowitsch   PetscCall(VecGetSize(xin,&vecrows));
5969566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetSkipHeader(viewer,&skipheader));
59747c6ae99SBarry Smith   if (!write) {
59847c6ae99SBarry Smith     /* Read vector header. */
5990c169764Sdmay     if (!skipheader) {
6009566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryRead(viewer,tr,2,NULL,PETSC_INT));
60147c6ae99SBarry Smith       type = tr[0];
60247c6ae99SBarry Smith       rows = tr[1];
603*08401ef6SPierre Jolivet       PetscCheck(type == VEC_FILE_CLASSID,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Not vector next in file");
604*08401ef6SPierre Jolivet       PetscCheck(rows == vecrows,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Vector in file not same size as DMDA vector");
6050c169764Sdmay     }
60647c6ae99SBarry Smith   } else {
60747c6ae99SBarry Smith     tr[0] = VEC_FILE_CLASSID;
60847c6ae99SBarry Smith     tr[1] = vecrows;
6090c169764Sdmay     if (!skipheader) {
6109566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWrite(viewer,tr,2,PETSC_INT));
61147c6ae99SBarry Smith     }
6120c169764Sdmay   }
61347c6ae99SBarry Smith 
6149566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->w,&dof));
6154dc2109aSBarry Smith   gsizes[0]  = dof;
6169566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->M,gsizes+1));
6179566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->N,gsizes+2));
6189566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->P,gsizes+3));
6194dc2109aSBarry Smith   lsizes[0]  = dof;
6209566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast((dd->xe-dd->xs)/dof,lsizes+1));
6219566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->ye-dd->ys,lsizes+2));
6229566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->ze-dd->zs,lsizes+3));
6234dc2109aSBarry Smith   lstarts[0] = 0;
6249566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->xs/dof,lstarts+1));
6259566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->ys,lstarts+2));
6269566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->zs,lstarts+3));
6279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_create_subarray(da->dim+1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view));
6289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&view));
62947c6ae99SBarry Smith 
6309566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes));
6319566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetMPIIOOffset(viewer,&off));
6329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_File_set_view(mfdes,off,MPIU_SCALAR,view,(char*)"native",MPI_INFO_NULL));
6339566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xin,&array));
63447c6ae99SBarry Smith   asiz = lsizes[1]*(lsizes[2] > 0 ? lsizes[2] : 1)*(lsizes[3] > 0 ? lsizes[3] : 1)*dof;
63547c6ae99SBarry Smith   if (write) {
6369566063dSJacob Faibussowitsch     PetscCall(MPIU_File_write_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE));
63747c6ae99SBarry Smith   } else {
6389566063dSJacob Faibussowitsch     PetscCall(MPIU_File_read_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE));
63947c6ae99SBarry Smith   }
6409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_get_extent(view,&ul,&ub));
6419566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryAddMPIIOOffset(viewer,ub));
6429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xin,&array));
6439566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&view));
64447c6ae99SBarry Smith   PetscFunctionReturn(0);
64547c6ae99SBarry Smith }
64647c6ae99SBarry Smith #endif
64747c6ae99SBarry Smith 
6487087cfbeSBarry Smith PetscErrorCode  VecView_MPI_DA(Vec xin,PetscViewer viewer)
64947c6ae99SBarry Smith {
6509a42bb27SBarry Smith   DM                da;
65147c6ae99SBarry Smith   PetscInt          dim;
65247c6ae99SBarry Smith   Vec               natural;
6538135c375SStefano Zampini   PetscBool         isdraw,isvtk,isglvis;
65447c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
65547c6ae99SBarry Smith   PetscBool         ishdf5;
65647c6ae99SBarry Smith #endif
6573f3fd955SJed Brown   const char        *prefix,*name;
658a261c58fSBarry Smith   PetscViewerFormat format;
65947c6ae99SBarry Smith 
66047c6ae99SBarry Smith   PetscFunctionBegin;
6619566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin,&da));
6627a8be351SBarry Smith   PetscCheck(da,PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
6639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw));
6649566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk));
66547c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
6669566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5));
66747c6ae99SBarry Smith #endif
6689566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis));
66947c6ae99SBarry Smith   if (isdraw) {
6709566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
67147c6ae99SBarry Smith     if (dim == 1) {
6729566063dSJacob Faibussowitsch       PetscCall(VecView_MPI_Draw_DA1d(xin,viewer));
67347c6ae99SBarry Smith     } else if (dim == 2) {
6749566063dSJacob Faibussowitsch       PetscCall(VecView_MPI_Draw_DA2d(xin,viewer));
67598921bdaSJacob Faibussowitsch     } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot graphically view vector associated with this dimensional DMDA %D",dim);
67601d7c5c3SPatrick Sanan   } else if (isvtk) {           /* Duplicate the Vec */
6774061b8bfSJed Brown     Vec Y;
6789566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xin,&Y));
679b51b94faSJed Brown     if (((PetscObject)xin)->name) {
680b51b94faSJed Brown       /* If xin was named, copy the name over to Y. The duplicate names are safe because nobody else will ever see Y. */
6819566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)Y,((PetscObject)xin)->name));
682b51b94faSJed Brown     }
6839566063dSJacob Faibussowitsch     PetscCall(VecCopy(xin,Y));
684a8f87f1dSPatrick Sanan     {
685a8f87f1dSPatrick Sanan       PetscObject dmvtk;
686a8f87f1dSPatrick Sanan       PetscBool   compatible,compatibleSet;
6879566063dSJacob Faibussowitsch       PetscCall(PetscViewerVTKGetDM(viewer,&dmvtk));
688a8f87f1dSPatrick Sanan       if (dmvtk) {
689064a246eSJacob Faibussowitsch         PetscValidHeaderSpecific((DM)dmvtk,DM_CLASSID,2);
6909566063dSJacob Faibussowitsch         PetscCall(DMGetCompatibility(da,(DM)dmvtk,&compatible,&compatibleSet));
6917a8be351SBarry 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.");
692a8f87f1dSPatrick Sanan       }
6939566063dSJacob Faibussowitsch       PetscCall(PetscViewerVTKAddField(viewer,(PetscObject)da,DMDAVTKWriteAll,PETSC_DEFAULT,PETSC_VTK_POINT_FIELD,PETSC_FALSE,(PetscObject)Y));
694a8f87f1dSPatrick Sanan     }
69547c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
69647c6ae99SBarry Smith   } else if (ishdf5) {
6979566063dSJacob Faibussowitsch     PetscCall(VecView_MPI_HDF5_DA(xin,viewer));
69847c6ae99SBarry Smith #endif
6998135c375SStefano Zampini   } else if (isglvis) {
7009566063dSJacob Faibussowitsch     PetscCall(VecView_GLVis(xin,viewer));
70147c6ae99SBarry Smith   } else {
70247c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
70347c6ae99SBarry Smith     PetscBool isbinary,isMPIIO;
70447c6ae99SBarry Smith 
7059566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary));
70647c6ae99SBarry Smith     if (isbinary) {
7079566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO));
70847c6ae99SBarry Smith       if (isMPIIO) {
7099566063dSJacob Faibussowitsch         PetscCall(DMDAArrayMPIIO(da,viewer,xin,PETSC_TRUE));
71047c6ae99SBarry Smith         PetscFunctionReturn(0);
71147c6ae99SBarry Smith       }
71247c6ae99SBarry Smith     }
71347c6ae99SBarry Smith #endif
71447c6ae99SBarry Smith 
71547c6ae99SBarry Smith     /* call viewer on natural ordering */
7169566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix));
7179566063dSJacob Faibussowitsch     PetscCall(DMDACreateNaturalVector(da,&natural));
7189566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)natural,prefix));
7199566063dSJacob Faibussowitsch     PetscCall(DMDAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural));
7209566063dSJacob Faibussowitsch     PetscCall(DMDAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural));
7219566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)xin,&name));
7229566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)natural,name));
723a261c58fSBarry Smith 
7249566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer,&format));
725a261c58fSBarry Smith     if (format == PETSC_VIEWER_BINARY_MATLAB) {
726a261c58fSBarry Smith       /* temporarily remove viewer format so it won't trigger in the VecView() */
7279566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT));
728a261c58fSBarry Smith     }
729a261c58fSBarry Smith 
73023f88b65SBarry Smith     ((PetscObject)natural)->donotPetscObjectPrintClassNamePrefixType = PETSC_TRUE;
7319566063dSJacob Faibussowitsch     PetscCall(VecView(natural,viewer));
73223f88b65SBarry Smith     ((PetscObject)natural)->donotPetscObjectPrintClassNamePrefixType = PETSC_FALSE;
733a261c58fSBarry Smith 
734a261c58fSBarry Smith     if (format == PETSC_VIEWER_BINARY_MATLAB) {
735a261c58fSBarry Smith       MPI_Comm    comm;
736a261c58fSBarry Smith       FILE        *info;
737a261c58fSBarry Smith       const char  *fieldname;
738da88d4d4SJed Brown       char        fieldbuf[256];
739a261c58fSBarry Smith       PetscInt    dim,ni,nj,nk,pi,pj,pk,dof,n;
740a261c58fSBarry Smith 
741a261c58fSBarry Smith       /* set the viewer format back into the viewer */
7429566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(viewer));
7439566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetComm((PetscObject)viewer,&comm));
7449566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryGetInfoPointer(viewer,&info));
7459566063dSJacob Faibussowitsch       PetscCall(DMDAGetInfo(da,&dim,&ni,&nj,&nk,&pi,&pj,&pk,&dof,NULL,NULL,NULL,NULL,NULL));
7469566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
7479566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm,info,"#$$ tmp = PetscBinaryRead(fd); \n"));
7489566063dSJacob Faibussowitsch       if (dim == 1) { PetscCall(PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d);\n",dof,ni)); }
7499566063dSJacob Faibussowitsch       if (dim == 2) { PetscCall(PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d);\n",dof,ni,nj)); }
7509566063dSJacob Faibussowitsch       if (dim == 3) { PetscCall(PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d,%d);\n",dof,ni,nj,nk)); }
751a261c58fSBarry Smith 
752a261c58fSBarry Smith       for (n=0; n<dof; n++) {
7539566063dSJacob Faibussowitsch         PetscCall(DMDAGetFieldName(da,n,&fieldname));
754da88d4d4SJed Brown         if (!fieldname || !fieldname[0]) {
7559566063dSJacob Faibussowitsch           PetscCall(PetscSNPrintf(fieldbuf,sizeof fieldbuf,"field%D",n));
756da88d4d4SJed Brown           fieldname = fieldbuf;
757a261c58fSBarry Smith         }
7589566063dSJacob Faibussowitsch         if (dim == 1) { PetscCall(PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:))';\n",name,fieldname,n+1)); }
7599566063dSJacob Faibussowitsch         if (dim == 2) { PetscCall(PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:,:))';\n",name,fieldname,n+1)); }
7609566063dSJacob Faibussowitsch         if (dim == 3) { PetscCall(PetscFPrintf(comm,info,"#$$ Set.%s.%s = permute(squeeze(tmp(%d,:,:,:)),[2 1 3]);\n",name,fieldname,n+1));}
761a261c58fSBarry Smith       }
7629566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm,info,"#$$ clear tmp; \n"));
7639566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
764a261c58fSBarry Smith     }
765a261c58fSBarry Smith 
7669566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&natural));
76747c6ae99SBarry Smith   }
76847c6ae99SBarry Smith   PetscFunctionReturn(0);
76947c6ae99SBarry Smith }
77047c6ae99SBarry Smith 
77147c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
77247c6ae99SBarry Smith PetscErrorCode VecLoad_HDF5_DA(Vec xin, PetscViewer viewer)
77347c6ae99SBarry Smith {
7742d95be3dSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7759a42bb27SBarry Smith   DM             da;
77615b861d2SVaclav Hapla   int            dim,rdim;
777ec7ae49cSHåkon Strandenes   hsize_t        dims[6]={0},count[6]={0},offset[6]={0};
778d7dd068bSVaclav Hapla   PetscBool      dim2=PETSC_FALSE,timestepping=PETSC_FALSE;
779d7dd068bSVaclav Hapla   PetscInt       dimension,timestep=PETSC_MIN_INT,dofInd;
78047c6ae99SBarry Smith   PetscScalar    *x;
78147c6ae99SBarry Smith   const char     *vecname;
78247c6ae99SBarry Smith   hid_t          filespace; /* file dataspace identifier */
78347c6ae99SBarry Smith   hid_t          dset_id;   /* dataset identifier */
78447c6ae99SBarry Smith   hid_t          memspace;  /* memory dataspace identifier */
785bfd264e7SBarry Smith   hid_t          file_id,group;
7867bcbaff4SBarry Smith   hid_t          scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
7879c7c4993SBarry Smith   DM_DA          *dd;
78847c6ae99SBarry Smith 
78947c6ae99SBarry Smith   PetscFunctionBegin;
7907bcbaff4SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
7917bcbaff4SBarry Smith   scalartype = H5T_NATIVE_FLOAT;
7927bcbaff4SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
7937bcbaff4SBarry Smith #error "HDF5 output with 128 bit floats not supported."
794570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
795570b7f6dSBarry Smith #error "HDF5 output with 16 bit floats not supported."
7967bcbaff4SBarry Smith #else
7977bcbaff4SBarry Smith   scalartype = H5T_NATIVE_DOUBLE;
7987bcbaff4SBarry Smith #endif
7997bcbaff4SBarry Smith 
8009566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5OpenGroup(viewer, &file_id, &group));
8019566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)xin,&vecname));
8029566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckTimestepping_Internal(viewer, vecname));
8039566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5IsTimestepping(viewer, &timestepping));
804d7dd068bSVaclav Hapla   if (timestepping) {
8059566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5GetTimestep(viewer, &timestep));
806d7dd068bSVaclav Hapla   }
8079566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin,&da));
8089c7c4993SBarry Smith   dd   = (DM_DA*)da->data;
8099566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(da, &dimension));
81047c6ae99SBarry Smith 
811ec7ae49cSHåkon Strandenes   /* Open dataset */
812729ab6d9SBarry Smith   PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT));
81347c6ae99SBarry Smith 
814ec7ae49cSHåkon Strandenes   /* Retrieve the dataspace for the dataset */
815ec7ae49cSHåkon Strandenes   PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id));
816ec7ae49cSHåkon Strandenes   PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(filespace, dims, NULL));
817ec7ae49cSHåkon Strandenes 
818ec7ae49cSHåkon Strandenes   /* Expected dimension for holding the dof's */
81947c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
820ec7ae49cSHåkon Strandenes   dofInd = rdim-2;
821ec7ae49cSHåkon Strandenes #else
822ec7ae49cSHåkon Strandenes   dofInd = rdim-1;
82347c6ae99SBarry Smith #endif
824ec7ae49cSHåkon Strandenes 
825ec7ae49cSHåkon Strandenes   /* The expected number of dimensions, assuming basedimension2 = false */
826ec7ae49cSHåkon Strandenes   dim = dimension;
827ec7ae49cSHåkon Strandenes   if (dd->w > 1) ++dim;
828ec7ae49cSHåkon Strandenes   if (timestep >= 0) ++dim;
82947c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
830ec7ae49cSHåkon Strandenes   ++dim;
83147c6ae99SBarry Smith #endif
832ec7ae49cSHåkon Strandenes 
833ec7ae49cSHåkon Strandenes   /* In this case the input dataset have one extra, unexpected dimension. */
834ec7ae49cSHåkon Strandenes   if (rdim == dim+1) {
835ec7ae49cSHåkon Strandenes     /* In this case the block size unity */
836ec7ae49cSHåkon Strandenes     if (dd->w == 1 && dims[dofInd] == 1) dim2 = PETSC_TRUE;
837ec7ae49cSHåkon Strandenes 
838ec7ae49cSHåkon Strandenes     /* Special error message for the case where dof does not match the input file */
839*08401ef6SPierre Jolivet     else PetscCheck(dd->w == (PetscInt) dims[dofInd],PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Number of dofs in file is %D, not %D as expected",(PetscInt)dims[dofInd],dd->w);
840ec7ae49cSHåkon Strandenes 
841ec7ae49cSHåkon Strandenes   /* Other cases where rdim != dim cannot be handled currently */
842*08401ef6SPierre Jolivet   } else PetscCheck(rdim == dim,PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Dimension of array in file is %d, not %d as expected with dof = %D",rdim,dim,dd->w);
843ec7ae49cSHåkon Strandenes 
844ec7ae49cSHåkon Strandenes   /* Set up the hyperslab size */
845ec7ae49cSHåkon Strandenes   dim = 0;
846ec7ae49cSHåkon Strandenes   if (timestep >= 0) {
847ec7ae49cSHåkon Strandenes     offset[dim] = timestep;
848ec7ae49cSHåkon Strandenes     count[dim] = 1;
849ec7ae49cSHåkon Strandenes     ++dim;
850ec7ae49cSHåkon Strandenes   }
851ec7ae49cSHåkon Strandenes   if (dimension == 3) {
8529566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(dd->zs,offset + dim));
8539566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(dd->ze - dd->zs,count + dim));
854ec7ae49cSHåkon Strandenes     ++dim;
855ec7ae49cSHåkon Strandenes   }
856ec7ae49cSHåkon Strandenes   if (dimension > 1) {
8579566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(dd->ys,offset + dim));
8589566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(dd->ye - dd->ys,count + dim));
859ec7ae49cSHåkon Strandenes     ++dim;
860ec7ae49cSHåkon Strandenes   }
8619566063dSJacob Faibussowitsch   PetscCall(PetscHDF5IntCast(dd->xs/dd->w,offset + dim));
8629566063dSJacob Faibussowitsch   PetscCall(PetscHDF5IntCast((dd->xe - dd->xs)/dd->w,count + dim));
863ec7ae49cSHåkon Strandenes   ++dim;
864ec7ae49cSHåkon Strandenes   if (dd->w > 1 || dim2) {
865ec7ae49cSHåkon Strandenes     offset[dim] = 0;
8669566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(dd->w,count + dim));
867ec7ae49cSHåkon Strandenes     ++dim;
868ec7ae49cSHåkon Strandenes   }
869ec7ae49cSHåkon Strandenes #if defined(PETSC_USE_COMPLEX)
870ec7ae49cSHåkon Strandenes   offset[dim] = 0;
871ec7ae49cSHåkon Strandenes   count[dim] = 2;
872ec7ae49cSHåkon Strandenes   ++dim;
873ec7ae49cSHåkon Strandenes #endif
874ec7ae49cSHåkon Strandenes 
875ec7ae49cSHåkon Strandenes   /* Create the memory and filespace */
876729ab6d9SBarry Smith   PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL));
877729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL));
87847c6ae99SBarry Smith 
8799566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xin, &x));
8802d95be3dSVaclav Hapla   PetscStackCallHDF5(H5Dread,(dset_id, scalartype, memspace, filespace, hdf5->dxpl_id, x));
8819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xin, &x));
88247c6ae99SBarry Smith 
88347c6ae99SBarry Smith   /* Close/release resources */
884bfd264e7SBarry Smith   if (group != file_id) {
885729ab6d9SBarry Smith     PetscStackCallHDF5(H5Gclose,(group));
886bfd264e7SBarry Smith   }
887729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
888729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(memspace));
889729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dset_id));
89047c6ae99SBarry Smith   PetscFunctionReturn(0);
89147c6ae99SBarry Smith }
89247c6ae99SBarry Smith #endif
89347c6ae99SBarry Smith 
89447c6ae99SBarry Smith PetscErrorCode VecLoad_Binary_DA(Vec xin, PetscViewer viewer)
89547c6ae99SBarry Smith {
8969a42bb27SBarry Smith   DM             da;
89747c6ae99SBarry Smith   Vec            natural;
89847c6ae99SBarry Smith   const char     *prefix;
89947c6ae99SBarry Smith   PetscInt       bs;
90047c6ae99SBarry Smith   PetscBool      flag;
90147c6ae99SBarry Smith   DM_DA          *dd;
90247c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
90347c6ae99SBarry Smith   PetscBool isMPIIO;
90447c6ae99SBarry Smith #endif
90547c6ae99SBarry Smith 
90647c6ae99SBarry Smith   PetscFunctionBegin;
9079566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin,&da));
90847c6ae99SBarry Smith   dd   = (DM_DA*)da->data;
90947c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
9109566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO));
91147c6ae99SBarry Smith   if (isMPIIO) {
9129566063dSJacob Faibussowitsch     PetscCall(DMDAArrayMPIIO(da,viewer,xin,PETSC_FALSE));
91347c6ae99SBarry Smith     PetscFunctionReturn(0);
91447c6ae99SBarry Smith   }
91547c6ae99SBarry Smith #endif
91647c6ae99SBarry Smith 
9179566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix));
9189566063dSJacob Faibussowitsch   PetscCall(DMDACreateNaturalVector(da,&natural));
9199566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name));
9209566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)natural,prefix));
9219566063dSJacob Faibussowitsch   PetscCall(VecLoad(natural,viewer));
9229566063dSJacob Faibussowitsch   PetscCall(DMDANaturalToGlobalBegin(da,natural,INSERT_VALUES,xin));
9239566063dSJacob Faibussowitsch   PetscCall(DMDANaturalToGlobalEnd(da,natural,INSERT_VALUES,xin));
9249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&natural));
9259566063dSJacob Faibussowitsch   PetscCall(PetscInfo(xin,"Loading vector from natural ordering into DMDA\n"));
9269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,((PetscObject)xin)->prefix,"-vecload_block_size",&bs,&flag));
92747c6ae99SBarry Smith   if (flag && bs != dd->w) {
9289566063dSJacob Faibussowitsch     PetscCall(PetscInfo(xin,"Block size in file %D not equal to DMDA's dof %D\n",bs,dd->w));
92947c6ae99SBarry Smith   }
93047c6ae99SBarry Smith   PetscFunctionReturn(0);
93147c6ae99SBarry Smith }
93247c6ae99SBarry Smith 
9337087cfbeSBarry Smith PetscErrorCode  VecLoad_Default_DA(Vec xin, PetscViewer viewer)
93447c6ae99SBarry Smith {
9359a42bb27SBarry Smith   DM             da;
93647c6ae99SBarry Smith   PetscBool      isbinary;
93747c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
93847c6ae99SBarry Smith   PetscBool ishdf5;
93947c6ae99SBarry Smith #endif
94047c6ae99SBarry Smith 
94147c6ae99SBarry Smith   PetscFunctionBegin;
9429566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin,&da));
9437a8be351SBarry Smith   PetscCheck(da,PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
94447c6ae99SBarry Smith 
94547c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
9469566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5));
94747c6ae99SBarry Smith #endif
9489566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary));
94947c6ae99SBarry Smith 
95047c6ae99SBarry Smith   if (isbinary) {
9519566063dSJacob Faibussowitsch     PetscCall(VecLoad_Binary_DA(xin,viewer));
95247c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
95347c6ae99SBarry Smith   } else if (ishdf5) {
9549566063dSJacob Faibussowitsch     PetscCall(VecLoad_HDF5_DA(xin,viewer));
95547c6ae99SBarry Smith #endif
95698921bdaSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for vector loading", ((PetscObject)viewer)->type_name);
95747c6ae99SBarry Smith   PetscFunctionReturn(0);
95847c6ae99SBarry Smith }
959