xref: /petsc/src/dm/impls/da/gr2.c (revision d0609ced746bc51b019815ca91d747429db24893)
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;
3447c6ae99SBarry Smith 
3547c6ae99SBarry Smith   PetscFunctionBegin;
3647c6ae99SBarry Smith   m    = zctx->m;
3747c6ae99SBarry Smith   n    = zctx->n;
38e5ab1681SLisandro Dalcin   dof  = zctx->dof;
3947c6ae99SBarry Smith   k    = zctx->k;
4047c6ae99SBarry Smith   xy   = zctx->xy;
41e5ab1681SLisandro Dalcin   v    = zctx->v;
4247c6ae99SBarry Smith   min  = zctx->min;
43f3f0eb19SBarry Smith   max  = zctx->max;
4447c6ae99SBarry Smith 
4547c6ae99SBarry Smith   /* PetscDraw the contour plot patch */
46*d0609cedSBarry Smith   PetscDrawCollectiveBegin(draw);
4747c6ae99SBarry Smith   for (j=0; j<n-1; j++) {
4847c6ae99SBarry Smith     for (i=0; i<m-1; i++) {
490e4fe250SBarry Smith       id   = i+j*m;
500e4fe250SBarry Smith       x1   = PetscRealPart(xy[2*id]);
510e4fe250SBarry Smith       y_1  = PetscRealPart(xy[2*id+1]);
52e5ab1681SLisandro Dalcin       c1   = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max);
530e4fe250SBarry Smith 
540e4fe250SBarry Smith       id   = i+j*m+1;
550e4fe250SBarry Smith       x2   = PetscRealPart(xy[2*id]);
56a39a4c7dSLisandro Dalcin       y2   = PetscRealPart(xy[2*id+1]);
57e5ab1681SLisandro Dalcin       c2   = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max);
580e4fe250SBarry Smith 
590e4fe250SBarry Smith       id   = i+j*m+1+m;
60a39a4c7dSLisandro Dalcin       x3   = PetscRealPart(xy[2*id]);
610e4fe250SBarry Smith       y3   = PetscRealPart(xy[2*id+1]);
62e5ab1681SLisandro Dalcin       c3   = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max);
630e4fe250SBarry Smith 
640e4fe250SBarry Smith       id   = i+j*m+m;
65a39a4c7dSLisandro Dalcin       x4   = PetscRealPart(xy[2*id]);
66a39a4c7dSLisandro Dalcin       y4   = PetscRealPart(xy[2*id+1]);
67e5ab1681SLisandro Dalcin       c4   = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max);
68f3f0eb19SBarry Smith 
699566063dSJacob Faibussowitsch       PetscCall(PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3));
709566063dSJacob Faibussowitsch       PetscCall(PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4));
7147c6ae99SBarry Smith       if (zctx->showgrid) {
729566063dSJacob Faibussowitsch         PetscCall(PetscDrawLine(draw,x1,y_1,x2,y2,PETSC_DRAW_BLACK));
739566063dSJacob Faibussowitsch         PetscCall(PetscDrawLine(draw,x2,y2,x3,y3,PETSC_DRAW_BLACK));
749566063dSJacob Faibussowitsch         PetscCall(PetscDrawLine(draw,x3,y3,x4,y4,PETSC_DRAW_BLACK));
759566063dSJacob Faibussowitsch         PetscCall(PetscDrawLine(draw,x4,y4,x1,y_1,PETSC_DRAW_BLACK));
7647c6ae99SBarry Smith       }
7747c6ae99SBarry Smith     }
7847c6ae99SBarry Smith   }
797fe7c8feSLisandro Dalcin   if (zctx->showaxis && !zctx->rank) {
80e5ab1681SLisandro Dalcin     if (zctx->name0 || zctx->name1) {
81109c9344SBarry Smith       PetscReal xl,yl,xr,yr,x,y;
829566063dSJacob Faibussowitsch       PetscCall(PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr));
83e5ab1681SLisandro Dalcin       x  = xl + .30*(xr - xl);
84109c9344SBarry Smith       xl = xl + .01*(xr - xl);
85e5ab1681SLisandro Dalcin       y  = yr - .30*(yr - yl);
86109c9344SBarry Smith       yl = yl + .01*(yr - yl);
879566063dSJacob Faibussowitsch       if (zctx->name0) PetscCall(PetscDrawString(draw,x,yl,PETSC_DRAW_BLACK,zctx->name0));
889566063dSJacob Faibussowitsch       if (zctx->name1) PetscCall(PetscDrawStringVertical(draw,xl,y,PETSC_DRAW_BLACK,zctx->name1));
89109c9344SBarry Smith     }
900e4fe250SBarry Smith     /*
910e4fe250SBarry Smith        Ideally we would use the PetscDrawAxis object to manage displaying the coordinate limits
920e4fe250SBarry Smith        but that may require some refactoring.
930e4fe250SBarry Smith     */
94e5ab1681SLisandro Dalcin     {
957fe7c8feSLisandro Dalcin       double xmin = (double)zctx->xmin, ymin = (double)zctx->ymin;
967fe7c8feSLisandro Dalcin       double xmax = (double)zctx->xmax, ymax = (double)zctx->ymax;
97e5ab1681SLisandro Dalcin       char   value[16]; size_t len; PetscReal w;
989566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(value,16,"%0.2e",xmin));
999566063dSJacob Faibussowitsch       PetscCall(PetscDrawString(draw,xmin,ymin - .05*(ymax - ymin),PETSC_DRAW_BLACK,value));
1009566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(value,16,"%0.2e",xmax));
1019566063dSJacob Faibussowitsch       PetscCall(PetscStrlen(value,&len));
1029566063dSJacob Faibussowitsch       PetscCall(PetscDrawStringGetSize(draw,&w,NULL));
1039566063dSJacob Faibussowitsch       PetscCall(PetscDrawString(draw,xmax - len*w,ymin - .05*(ymax - ymin),PETSC_DRAW_BLACK,value));
1049566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(value,16,"%0.2e",ymin));
1059566063dSJacob Faibussowitsch       PetscCall(PetscDrawString(draw,xmin - .05*(xmax - xmin),ymin,PETSC_DRAW_BLACK,value));
1069566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(value,16,"%0.2e",ymax));
1079566063dSJacob Faibussowitsch       PetscCall(PetscDrawString(draw,xmin - .05*(xmax - xmin),ymax,PETSC_DRAW_BLACK,value));
108e5ab1681SLisandro Dalcin     }
109e5ab1681SLisandro Dalcin   }
110*d0609cedSBarry Smith   PetscDrawCollectiveEnd(draw);
11147c6ae99SBarry Smith   PetscFunctionReturn(0);
11247c6ae99SBarry Smith }
11347c6ae99SBarry Smith 
11447c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA2d(Vec xin,PetscViewer viewer)
11547c6ae99SBarry Smith {
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));
1669566063dSJacob Faibussowitsch       PetscCall(PetscInfo(da,"Creating auxilary 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) {
2177fe7c8feSLisandro Dalcin     coors[0] = zctx.xmin - .05*(zctx.xmax - zctx.xmin); coors[1] = zctx.ymin - .05*(zctx.ymax - zctx.ymin);
2187fe7c8feSLisandro Dalcin     coors[2] = zctx.xmax + .05*(zctx.xmax - zctx.xmin); coors[3] = zctx.ymax + .05*(zctx.ymax - zctx.ymin);
2197fe7c8feSLisandro Dalcin   } else {
2207fe7c8feSLisandro Dalcin     coors[0] = zctx.xmin; coors[1] = zctx.ymin; coors[2] = zctx.xmax; coors[3] = zctx.ymax;
2217fe7c8feSLisandro Dalcin   }
2229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetRealArray(NULL,NULL,"-draw_coordinates",coors,&ncoors,NULL));
2239566063dSJacob 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]));
224a74ba6f7SBarry Smith 
22547c6ae99SBarry Smith   /*
226832b7cebSLisandro Dalcin       Get local ghosted version of coordinates
22747c6ae99SBarry Smith   */
2289566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject*)&xcoorl));
22947c6ae99SBarry Smith   if (!xcoorl) {
230aa219208SBarry Smith     /* create DMDA to get local version of graphics */
2319566063dSJacob Faibussowitsch     PetscCall(DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,zctx.m,zctx.n,2,1,lx,ly,&dag));
2329566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dag));
2339566063dSJacob Faibussowitsch     PetscCall(PetscInfo(dag,"Creating auxilary DMDA for managing graphics coordinates ghost points\n"));
2349566063dSJacob Faibussowitsch     PetscCall(DMCreateLocalVector(dag,&xcoorl));
2359566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject)xcoorl));
2369566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)dag));
2379566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)xcoorl));
23847c6ae99SBarry Smith   } else {
2399566063dSJacob Faibussowitsch     PetscCall(VecGetDM(xcoorl,&dag));
24047c6ae99SBarry Smith   }
2419566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dag,xcoor,INSERT_VALUES,xcoorl));
2429566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dag,xcoor,INSERT_VALUES,xcoorl));
2439566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xcoorl,&zctx.xy));
2449566063dSJacob Faibussowitsch   PetscCall(DMDAGetCoordinateName(da,0,&zctx.name0));
2459566063dSJacob Faibussowitsch   PetscCall(DMDAGetCoordinateName(da,1,&zctx.name1));
24647c6ae99SBarry Smith 
24747c6ae99SBarry Smith   /*
24847c6ae99SBarry Smith       Get information about size of area each processor must do graphics for
24947c6ae99SBarry Smith   */
2509566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dac,NULL,&M,&N,NULL,NULL,NULL,NULL,&zctx.dof,NULL,&bx,&by,NULL,NULL));
2519566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dac,NULL,NULL,NULL,&zctx.m,&zctx.n,NULL));
2529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-draw_contour_grid",&zctx.showgrid,NULL));
2534e6118eeSBarry Smith 
2549566063dSJacob Faibussowitsch   PetscCall(DMDASelectFields(da,&ndisplayfields,&displayfields));
2559566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer,&format));
2569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-draw_ports",&useports,NULL));
257832b7cebSLisandro Dalcin   if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE;
258832b7cebSLisandro Dalcin   if (useports) {
2599566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw));
2609566063dSJacob Faibussowitsch     PetscCall(PetscDrawCheckResizedWindow(draw));
2619566063dSJacob Faibussowitsch     PetscCall(PetscDrawClear(draw));
2629566063dSJacob Faibussowitsch     PetscCall(PetscDrawViewPortsCreate(draw,ndisplayfields,&ports));
263e5ab1681SLisandro Dalcin   }
26420d0051dSBarry Smith 
265832b7cebSLisandro Dalcin   /*
266832b7cebSLisandro Dalcin       Loop over each field; drawing each in a different window
267832b7cebSLisandro Dalcin   */
26820d0051dSBarry Smith   for (i=0; i<ndisplayfields; i++) {
26920d0051dSBarry Smith     zctx.k = displayfields[i];
2705b399a63SLisandro Dalcin 
271e5ab1681SLisandro Dalcin     /* determine the min and max value in plot */
2729566063dSJacob Faibussowitsch     PetscCall(VecStrideMin(xin,zctx.k,NULL,&zctx.min));
2739566063dSJacob Faibussowitsch     PetscCall(VecStrideMax(xin,zctx.k,NULL,&zctx.max));
27467dd0837SBarry Smith     if (zctx.k < nbounds) {
275f3f0eb19SBarry Smith       zctx.min = bounds[2*zctx.k];
276f3f0eb19SBarry Smith       zctx.max = bounds[2*zctx.k+1];
27767dd0837SBarry Smith     }
27847c6ae99SBarry Smith     if (zctx.min == zctx.max) {
27947c6ae99SBarry Smith       zctx.min -= 1.e-12;
28047c6ae99SBarry Smith       zctx.max += 1.e-12;
28147c6ae99SBarry Smith     }
2829566063dSJacob Faibussowitsch     PetscCall(PetscInfo(da,"DMDA 2d contour plot min %g max %g\n",(double)zctx.min,(double)zctx.max));
28347c6ae99SBarry Smith 
284832b7cebSLisandro Dalcin     if (useports) {
2859566063dSJacob Faibussowitsch       PetscCall(PetscDrawViewPortsSet(ports,i));
286832b7cebSLisandro Dalcin     } else {
287832b7cebSLisandro Dalcin       const char *title;
2889566063dSJacob Faibussowitsch       PetscCall(PetscViewerDrawGetDraw(viewer,i,&draw));
2899566063dSJacob Faibussowitsch       PetscCall(DMDAGetFieldName(da,zctx.k,&title));
2909566063dSJacob Faibussowitsch       if (title) PetscCall(PetscDrawSetTitle(draw,title));
291832b7cebSLisandro Dalcin     }
292832b7cebSLisandro Dalcin 
2939566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetPopup(draw,&popup));
2949566063dSJacob Faibussowitsch     PetscCall(PetscDrawScalePopup(popup,zctx.min,zctx.max));
2959566063dSJacob Faibussowitsch     PetscCall(PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]));
2969566063dSJacob Faibussowitsch     PetscCall(PetscDrawZoom(draw,VecView_MPI_Draw_DA2d_Zoom,&zctx));
2979566063dSJacob Faibussowitsch     if (!useports) PetscCall(PetscDrawSave(draw));
29847c6ae99SBarry Smith   }
299832b7cebSLisandro Dalcin   if (useports) {
3009566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw));
3019566063dSJacob Faibussowitsch     PetscCall(PetscDrawSave(draw));
302832b7cebSLisandro Dalcin   }
303832b7cebSLisandro Dalcin 
3049566063dSJacob Faibussowitsch   PetscCall(PetscDrawViewPortsDestroy(ports));
3059566063dSJacob Faibussowitsch   PetscCall(PetscFree(displayfields));
3069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xcoorl,&zctx.xy));
3079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xlocal,&zctx.v));
30847c6ae99SBarry Smith   PetscFunctionReturn(0);
30947c6ae99SBarry Smith }
31047c6ae99SBarry Smith 
3110da24e51SJuha Jäykkä #if defined(PETSC_HAVE_HDF5)
312c73cfb54SMatthew G. Knepley static PetscErrorCode VecGetHDF5ChunkSize(DM_DA *da, Vec xin, PetscInt dimension, PetscInt timestep, hsize_t *chunkDims)
3130da24e51SJuha Jäykkä {
314d4190030SJed Brown   PetscMPIInt    comm_size;
315561a82e7SJed Brown   hsize_t        chunk_size, target_size, dim;
316561a82e7SJed Brown   hsize_t        vec_size = sizeof(PetscScalar)*da->M*da->N*da->P*da->w;
317561a82e7SJed Brown   hsize_t        avg_local_vec_size,KiB = 1024,MiB = KiB*KiB,GiB = MiB*KiB,min_size = MiB;
318561a82e7SJed Brown   hsize_t        max_chunks = 64*KiB;                                              /* HDF5 internal limitation */
319561a82e7SJed Brown   hsize_t        max_chunk_size = 4*GiB;                                           /* HDF5 internal limitation */
32084179ffaSJed Brown   PetscInt       zslices=da->p, yslices=da->n, xslices=da->m;
3210da24e51SJuha Jäykkä 
322cb5c4f94SJuha Jäykkä   PetscFunctionBegin;
3239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)xin), &comm_size));
324faa75363SBarry Smith   avg_local_vec_size = (hsize_t) PetscCeilInt(vec_size,comm_size);      /* we will attempt to use this as the chunk size */
3250da24e51SJuha Jäykkä 
326faa75363SBarry 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))));
3277d310018SBarry 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 */
3287d310018SBarry 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);
3290da24e51SJuha Jäykkä 
330cb5c4f94SJuha Jäykkä   /*
331cb5c4f94SJuha Jäykkä    if size/rank > max_chunk_size, we need radical measures: even going down to
332cb5c4f94SJuha Jäykkä    avg_local_vec_size is not enough, so we simply use chunk size of 4 GiB no matter
333cb5c4f94SJuha Jäykkä    what, composed in the most efficient way possible.
334cb5c4f94SJuha Jäykkä    N.B. this minimises the number of chunks, which may or may not be the optimal
335cb5c4f94SJuha Jäykkä    solution. In a BG, for example, the optimal solution is probably to make # chunks = #
336cb5c4f94SJuha Jäykkä    IO nodes involved, but this author has no access to a BG to figure out how to
337cb5c4f94SJuha Jäykkä    reliably find the right number. And even then it may or may not be enough.
338cb5c4f94SJuha Jäykkä    */
3390da24e51SJuha Jäykkä   if (avg_local_vec_size > max_chunk_size) {
340cb5c4f94SJuha Jäykkä     /* check if we can just split local z-axis: is that enough? */
341faa75363SBarry Smith     zslices = PetscCeilInt(vec_size,da->p*max_chunk_size)*zslices;
3420da24e51SJuha Jäykkä     if (zslices > da->P) {
343cb5c4f94SJuha Jäykkä       /* lattice is too large in xy-directions, splitting z only is not enough */
3440da24e51SJuha Jäykkä       zslices = da->P;
345faa75363SBarry Smith       yslices = PetscCeilInt(vec_size,zslices*da->n*max_chunk_size)*yslices;
3460da24e51SJuha Jäykkä       if (yslices > da->N) {
347cb5c4f94SJuha Jäykkä         /* lattice is too large in x-direction, splitting along z, y is not enough */
3480da24e51SJuha Jäykkä         yslices = da->N;
349faa75363SBarry Smith         xslices = PetscCeilInt(vec_size,zslices*yslices*da->m*max_chunk_size)*xslices;
3500da24e51SJuha Jäykkä       }
3510da24e51SJuha Jäykkä     }
3520da24e51SJuha Jäykkä     dim = 0;
3530da24e51SJuha Jäykkä     if (timestep >= 0) {
3540da24e51SJuha Jäykkä       ++dim;
3550da24e51SJuha Jäykkä     }
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;
3710da24e51SJuha Jäykkä       if (timestep >= 0) {
3720da24e51SJuha Jäykkä         ++dim;
3730da24e51SJuha Jäykkä       }
374cb5c4f94SJuha Jäykkä       /* prefer to split z-axis, even down to planar slices */
375c73cfb54SMatthew G. Knepley       if (dimension == 3) {
376cb5c4f94SJuha Jäykkä         /* try splitting the z-axis to core-size bits, i.e. divide chunk size by # comm_size in z-direction */
3770da24e51SJuha Jäykkä         if (target_size >= chunk_size/da->p) {
378cb5c4f94SJuha Jäykkä           /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */
379faa75363SBarry Smith           chunkDims[dim] = (hsize_t) PetscCeilInt(da->P,da->p);
3800da24e51SJuha Jäykkä         } else {
381cb5c4f94SJuha Jäykkä           /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be
382cb5c4f94SJuha Jäykkä            radical and let everyone write all they've got */
383faa75363SBarry Smith           chunkDims[dim++] = (hsize_t) PetscCeilInt(da->P,da->p);
384faa75363SBarry Smith           chunkDims[dim++] = (hsize_t) PetscCeilInt(da->N,da->n);
385faa75363SBarry Smith           chunkDims[dim++] = (hsize_t) PetscCeilInt(da->M,da->m);
3860da24e51SJuha Jäykkä         }
3870da24e51SJuha Jäykkä       } else {
388cb5c4f94SJuha Jäykkä         /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */
3890da24e51SJuha Jäykkä         if (target_size >= chunk_size/da->n) {
390cb5c4f94SJuha Jäykkä           /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */
391faa75363SBarry Smith           chunkDims[dim] = (hsize_t) PetscCeilInt(da->N,da->n);
3920da24e51SJuha Jäykkä         } else {
393cb5c4f94SJuha Jäykkä           /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be
394cb5c4f94SJuha Jäykkä            radical and let everyone write all they've got */
395faa75363SBarry Smith           chunkDims[dim++] = (hsize_t) PetscCeilInt(da->N,da->n);
396faa75363SBarry Smith           chunkDims[dim++] = (hsize_t) PetscCeilInt(da->M,da->m);
3970da24e51SJuha Jäykkä         }
3980da24e51SJuha Jäykkä 
3990da24e51SJuha Jäykkä       }
4000da24e51SJuha Jäykkä       chunk_size = (hsize_t) PetscMax(1,chunkDims[0])*PetscMax(1,chunkDims[1])*PetscMax(1,chunkDims[2])*PetscMax(1,chunkDims[3])*PetscMax(1,chunkDims[4])*PetscMax(1,chunkDims[5])*sizeof(double);
4010da24e51SJuha Jäykkä     } else {
402cb5c4f94SJuha Jäykkä       /* precomputed chunks are fine, we don't need to do anything */
4030da24e51SJuha Jäykkä     }
4040da24e51SJuha Jäykkä   }
4050da24e51SJuha Jäykkä   PetscFunctionReturn(0);
4060da24e51SJuha Jäykkä }
4070da24e51SJuha Jäykkä #endif
40847c6ae99SBarry Smith 
40947c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
41047c6ae99SBarry Smith PetscErrorCode VecView_MPI_HDF5_DA(Vec xin,PetscViewer viewer)
41147c6ae99SBarry Smith {
4122d95be3dSVaclav Hapla   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*) viewer->data;
4139b2a5a86SJed Brown   DM                dm;
4149b2a5a86SJed Brown   DM_DA             *da;
41547c6ae99SBarry Smith   hid_t             filespace;  /* file dataspace identifier */
4168e2ae6d7SMichael Kraus   hid_t             chunkspace; /* chunk dataset property identifier */
41747c6ae99SBarry Smith   hid_t             dset_id;    /* dataset identifier */
41847c6ae99SBarry Smith   hid_t             memspace;   /* memory dataspace identifier */
41947c6ae99SBarry Smith   hid_t             file_id;
42015214e8eSMatthew G Knepley   hid_t             group;
4219a0502c6SHåkon Strandenes   hid_t             memscalartype; /* scalar type for mem (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
4229a0502c6SHåkon Strandenes   hid_t             filescalartype; /* scalar type for file (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
423d9a4edebSJed Brown   hsize_t           dim;
4240da24e51SJuha Jäykkä   hsize_t           maxDims[6]={0}, dims[6]={0}, chunkDims[6]={0}, count[6]={0}, offset[6]={0}; /* we depend on these being sane later on  */
425d7dd068bSVaclav Hapla   PetscBool         timestepping=PETSC_FALSE, dim2, spoutput;
426d7dd068bSVaclav Hapla   PetscInt          timestep=PETSC_MIN_INT, dimension;
4275edff71fSBarry Smith   const PetscScalar *x;
4288e2ae6d7SMichael Kraus   const char        *vecname;
42947c6ae99SBarry Smith 
43047c6ae99SBarry Smith   PetscFunctionBegin;
4319566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5OpenGroup(viewer, &file_id, &group));
4329566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5IsTimestepping(viewer, &timestepping));
433d7dd068bSVaclav Hapla   if (timestepping) {
4349566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5GetTimestep(viewer, &timestep));
435d7dd068bSVaclav Hapla   }
4369566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetBaseDimension2(viewer,&dim2));
4379566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetSPOutput(viewer,&spoutput));
43815214e8eSMatthew G Knepley 
4399566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin,&dm));
4407a8be351SBarry Smith   PetscCheck(dm,PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
4419b2a5a86SJed Brown   da = (DM_DA*)dm->data;
4429566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dimension));
44347c6ae99SBarry Smith 
4448e2ae6d7SMichael Kraus   /* Create the dataspace for the dataset.
4458e2ae6d7SMichael Kraus    *
4468e2ae6d7SMichael Kraus    * dims - holds the current dimensions of the dataset
4478e2ae6d7SMichael Kraus    *
4488e2ae6d7SMichael Kraus    * maxDims - holds the maximum dimensions of the dataset (unlimited
4498e2ae6d7SMichael Kraus    * for the number of time steps with the current dimensions for the
4508e2ae6d7SMichael Kraus    * other dimensions; so only additional time steps can be added).
4518e2ae6d7SMichael Kraus    *
4528e2ae6d7SMichael Kraus    * chunkDims - holds the size of a single time step (required to
4538e2ae6d7SMichael Kraus    * permit extending dataset).
4548e2ae6d7SMichael Kraus    */
4558e2ae6d7SMichael Kraus   dim = 0;
4568e2ae6d7SMichael Kraus   if (timestep >= 0) {
4578e2ae6d7SMichael Kraus     dims[dim]      = timestep+1;
4588e2ae6d7SMichael Kraus     maxDims[dim]   = H5S_UNLIMITED;
4598e2ae6d7SMichael Kraus     chunkDims[dim] = 1;
4608e2ae6d7SMichael Kraus     ++dim;
4618e2ae6d7SMichael Kraus   }
462c73cfb54SMatthew G. Knepley   if (dimension == 3) {
4639566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(da->P,dims+dim));
4648e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4658e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4668e2ae6d7SMichael Kraus     ++dim;
4678e2ae6d7SMichael Kraus   }
468c73cfb54SMatthew G. Knepley   if (dimension > 1) {
4699566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(da->N,dims+dim));
4708e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4718e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4728e2ae6d7SMichael Kraus     ++dim;
4738e2ae6d7SMichael Kraus   }
4749566063dSJacob Faibussowitsch   PetscCall(PetscHDF5IntCast(da->M,dims+dim));
4758e2ae6d7SMichael Kraus   maxDims[dim]   = dims[dim];
4768e2ae6d7SMichael Kraus   chunkDims[dim] = dims[dim];
4778e2ae6d7SMichael Kraus   ++dim;
47882ea9b62SBarry Smith   if (da->w > 1 || dim2) {
4799566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(da->w,dims+dim));
4808e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4818e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4828e2ae6d7SMichael Kraus     ++dim;
4838e2ae6d7SMichael Kraus   }
48447c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
4858e2ae6d7SMichael Kraus   dims[dim]      = 2;
4868e2ae6d7SMichael Kraus   maxDims[dim]   = dims[dim];
4878e2ae6d7SMichael Kraus   chunkDims[dim] = dims[dim];
4888e2ae6d7SMichael Kraus   ++dim;
48947c6ae99SBarry Smith #endif
4900da24e51SJuha Jäykkä 
4919566063dSJacob Faibussowitsch   PetscCall(VecGetHDF5ChunkSize(da, xin, dimension, timestep, chunkDims));
4920da24e51SJuha Jäykkä 
493729ab6d9SBarry Smith   PetscStackCallHDF5Return(filespace,H5Screate_simple,(dim, dims, maxDims));
49447c6ae99SBarry Smith 
49515214e8eSMatthew G Knepley #if defined(PETSC_USE_REAL_SINGLE)
4969a0502c6SHåkon Strandenes   memscalartype = H5T_NATIVE_FLOAT;
4979a0502c6SHåkon Strandenes   filescalartype = H5T_NATIVE_FLOAT;
49815214e8eSMatthew G Knepley #elif defined(PETSC_USE_REAL___FLOAT128)
49915214e8eSMatthew G Knepley #error "HDF5 output with 128 bit floats not supported."
500570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
501570b7f6dSBarry Smith #error "HDF5 output with 16 bit floats not supported."
50215214e8eSMatthew G Knepley #else
5039a0502c6SHåkon Strandenes   memscalartype = H5T_NATIVE_DOUBLE;
5049a0502c6SHåkon Strandenes   if (spoutput == PETSC_TRUE) filescalartype = H5T_NATIVE_FLOAT;
5059a0502c6SHåkon Strandenes   else filescalartype = H5T_NATIVE_DOUBLE;
50615214e8eSMatthew G Knepley #endif
50715214e8eSMatthew G Knepley 
50847c6ae99SBarry Smith   /* Create the dataset with default properties and close filespace */
5099566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)xin,&vecname));
51015214e8eSMatthew G Knepley   if (!H5Lexists(group, vecname, H5P_DEFAULT)) {
5118e2ae6d7SMichael Kraus     /* Create chunk */
512729ab6d9SBarry Smith     PetscStackCallHDF5Return(chunkspace,H5Pcreate,(H5P_DATASET_CREATE));
513729ab6d9SBarry Smith     PetscStackCallHDF5(H5Pset_chunk,(chunkspace, dim, chunkDims));
5148e2ae6d7SMichael Kraus 
5159a0502c6SHåkon Strandenes     PetscStackCallHDF5Return(dset_id,H5Dcreate2,(group, vecname, filescalartype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT));
51615214e8eSMatthew G Knepley   } else {
517729ab6d9SBarry Smith     PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT));
518729ab6d9SBarry Smith     PetscStackCallHDF5(H5Dset_extent,(dset_id, dims));
51915214e8eSMatthew G Knepley   }
520729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
52147c6ae99SBarry Smith 
52247c6ae99SBarry Smith   /* Each process defines a dataset and writes it to the hyperslab in the file */
5238e2ae6d7SMichael Kraus   dim = 0;
5248e2ae6d7SMichael Kraus   if (timestep >= 0) {
5258e2ae6d7SMichael Kraus     offset[dim] = timestep;
5268e2ae6d7SMichael Kraus     ++dim;
5278e2ae6d7SMichael Kraus   }
5289566063dSJacob Faibussowitsch   if (dimension == 3) PetscCall(PetscHDF5IntCast(da->zs,offset + dim++));
5299566063dSJacob Faibussowitsch   if (dimension > 1)  PetscCall(PetscHDF5IntCast(da->ys,offset + dim++));
5309566063dSJacob Faibussowitsch   PetscCall(PetscHDF5IntCast(da->xs/da->w,offset + dim++));
53182ea9b62SBarry Smith   if (da->w > 1 || dim2) offset[dim++] = 0;
53247c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5338e2ae6d7SMichael Kraus   offset[dim++] = 0;
53447c6ae99SBarry Smith #endif
5358e2ae6d7SMichael Kraus   dim = 0;
5368e2ae6d7SMichael Kraus   if (timestep >= 0) {
5378e2ae6d7SMichael Kraus     count[dim] = 1;
5388e2ae6d7SMichael Kraus     ++dim;
5398e2ae6d7SMichael Kraus   }
5409566063dSJacob Faibussowitsch   if (dimension == 3) PetscCall(PetscHDF5IntCast(da->ze - da->zs,count + dim++));
5419566063dSJacob Faibussowitsch   if (dimension > 1)  PetscCall(PetscHDF5IntCast(da->ye - da->ys,count + dim++));
5429566063dSJacob Faibussowitsch   PetscCall(PetscHDF5IntCast((da->xe - da->xs)/da->w,count + dim++));
5439566063dSJacob Faibussowitsch   if (da->w > 1 || dim2) PetscCall(PetscHDF5IntCast(da->w,count + dim++));
54447c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5458e2ae6d7SMichael Kraus   count[dim++] = 2;
54647c6ae99SBarry Smith #endif
547729ab6d9SBarry Smith   PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL));
548729ab6d9SBarry Smith   PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id));
549729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL));
55047c6ae99SBarry Smith 
5519566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xin, &x));
5522d95be3dSVaclav Hapla   PetscStackCallHDF5(H5Dwrite,(dset_id, memscalartype, memspace, filespace, hdf5->dxpl_id, x));
553729ab6d9SBarry Smith   PetscStackCallHDF5(H5Fflush,(file_id, H5F_SCOPE_GLOBAL));
5549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xin, &x));
55547c6ae99SBarry Smith 
556e0552f36Ssajid__ali   #if defined(PETSC_USE_COMPLEX)
557e0552f36Ssajid__ali   {
558e0552f36Ssajid__ali     PetscBool tru = PETSC_TRUE;
5599566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer,(PetscObject)xin,"complex",PETSC_BOOL,&tru));
560e0552f36Ssajid__ali   }
561e0552f36Ssajid__ali   #endif
562c7a34aceSVaclav Hapla   if (timestepping) {
5639566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer,(PetscObject)xin,"timestepping",PETSC_BOOL,&timestepping));
564c7a34aceSVaclav Hapla   }
565e0552f36Ssajid__ali 
56647c6ae99SBarry Smith   /* Close/release resources */
56715214e8eSMatthew G Knepley   if (group != file_id) {
568729ab6d9SBarry Smith     PetscStackCallHDF5(H5Gclose,(group));
56915214e8eSMatthew G Knepley   }
570729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
571729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(memspace));
572729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dset_id));
5739566063dSJacob Faibussowitsch   PetscCall(PetscInfo(xin,"Wrote Vec object with name %s\n",vecname));
57447c6ae99SBarry Smith   PetscFunctionReturn(0);
57547c6ae99SBarry Smith }
57647c6ae99SBarry Smith #endif
57747c6ae99SBarry Smith 
57809573ac7SBarry Smith extern PetscErrorCode VecView_MPI_Draw_DA1d(Vec,PetscViewer);
57947c6ae99SBarry Smith 
58047c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
581aa219208SBarry Smith static PetscErrorCode DMDAArrayMPIIO(DM da,PetscViewer viewer,Vec xin,PetscBool write)
58247c6ae99SBarry Smith {
58347c6ae99SBarry Smith   MPI_File          mfdes;
58447c6ae99SBarry Smith   PetscMPIInt       gsizes[4],lsizes[4],lstarts[4],asiz,dof;
58547c6ae99SBarry Smith   MPI_Datatype      view;
58647c6ae99SBarry Smith   const PetscScalar *array;
58747c6ae99SBarry Smith   MPI_Offset        off;
58847c6ae99SBarry Smith   MPI_Aint          ub,ul;
58947c6ae99SBarry Smith   PetscInt          type,rows,vecrows,tr[2];
59047c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
5910c169764Sdmay   PetscBool         skipheader;
59247c6ae99SBarry Smith 
59347c6ae99SBarry Smith   PetscFunctionBegin;
5949566063dSJacob Faibussowitsch   PetscCall(VecGetSize(xin,&vecrows));
5959566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetSkipHeader(viewer,&skipheader));
59647c6ae99SBarry Smith   if (!write) {
59747c6ae99SBarry Smith     /* Read vector header. */
5980c169764Sdmay     if (!skipheader) {
5999566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryRead(viewer,tr,2,NULL,PETSC_INT));
60047c6ae99SBarry Smith       type = tr[0];
60147c6ae99SBarry Smith       rows = tr[1];
60208401ef6SPierre Jolivet       PetscCheck(type == VEC_FILE_CLASSID,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Not vector next in file");
60308401ef6SPierre Jolivet       PetscCheck(rows == vecrows,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Vector in file not same size as DMDA vector");
6040c169764Sdmay     }
60547c6ae99SBarry Smith   } else {
60647c6ae99SBarry Smith     tr[0] = VEC_FILE_CLASSID;
60747c6ae99SBarry Smith     tr[1] = vecrows;
6080c169764Sdmay     if (!skipheader) {
6099566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWrite(viewer,tr,2,PETSC_INT));
61047c6ae99SBarry Smith     }
6110c169764Sdmay   }
61247c6ae99SBarry Smith 
6139566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->w,&dof));
6144dc2109aSBarry Smith   gsizes[0]  = dof;
6159566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->M,gsizes+1));
6169566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->N,gsizes+2));
6179566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->P,gsizes+3));
6184dc2109aSBarry Smith   lsizes[0]  = dof;
6199566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast((dd->xe-dd->xs)/dof,lsizes+1));
6209566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->ye-dd->ys,lsizes+2));
6219566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->ze-dd->zs,lsizes+3));
6224dc2109aSBarry Smith   lstarts[0] = 0;
6239566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->xs/dof,lstarts+1));
6249566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->ys,lstarts+2));
6259566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(dd->zs,lstarts+3));
6269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_create_subarray(da->dim+1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view));
6279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&view));
62847c6ae99SBarry Smith 
6299566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes));
6309566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetMPIIOOffset(viewer,&off));
6319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_File_set_view(mfdes,off,MPIU_SCALAR,view,(char*)"native",MPI_INFO_NULL));
6329566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xin,&array));
63347c6ae99SBarry Smith   asiz = lsizes[1]*(lsizes[2] > 0 ? lsizes[2] : 1)*(lsizes[3] > 0 ? lsizes[3] : 1)*dof;
63447c6ae99SBarry Smith   if (write) {
6359566063dSJacob Faibussowitsch     PetscCall(MPIU_File_write_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE));
63647c6ae99SBarry Smith   } else {
6379566063dSJacob Faibussowitsch     PetscCall(MPIU_File_read_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE));
63847c6ae99SBarry Smith   }
6399566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_get_extent(view,&ul,&ub));
6409566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryAddMPIIOOffset(viewer,ub));
6419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xin,&array));
6429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&view));
64347c6ae99SBarry Smith   PetscFunctionReturn(0);
64447c6ae99SBarry Smith }
64547c6ae99SBarry Smith #endif
64647c6ae99SBarry Smith 
6477087cfbeSBarry Smith PetscErrorCode  VecView_MPI_DA(Vec xin,PetscViewer viewer)
64847c6ae99SBarry Smith {
6499a42bb27SBarry Smith   DM                da;
65047c6ae99SBarry Smith   PetscInt          dim;
65147c6ae99SBarry Smith   Vec               natural;
6528135c375SStefano Zampini   PetscBool         isdraw,isvtk,isglvis;
65347c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
65447c6ae99SBarry Smith   PetscBool         ishdf5;
65547c6ae99SBarry Smith #endif
6563f3fd955SJed Brown   const char        *prefix,*name;
657a261c58fSBarry Smith   PetscViewerFormat format;
65847c6ae99SBarry Smith 
65947c6ae99SBarry Smith   PetscFunctionBegin;
6609566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin,&da));
6617a8be351SBarry Smith   PetscCheck(da,PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
6629566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw));
6639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk));
66447c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
6659566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5));
66647c6ae99SBarry Smith #endif
6679566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis));
66847c6ae99SBarry Smith   if (isdraw) {
6699566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
67047c6ae99SBarry Smith     if (dim == 1) {
6719566063dSJacob Faibussowitsch       PetscCall(VecView_MPI_Draw_DA1d(xin,viewer));
67247c6ae99SBarry Smith     } else if (dim == 2) {
6739566063dSJacob Faibussowitsch       PetscCall(VecView_MPI_Draw_DA2d(xin,viewer));
67498921bdaSJacob Faibussowitsch     } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot graphically view vector associated with this dimensional DMDA %D",dim);
67501d7c5c3SPatrick Sanan   } else if (isvtk) {           /* Duplicate the Vec */
6764061b8bfSJed Brown     Vec Y;
6779566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xin,&Y));
678b51b94faSJed Brown     if (((PetscObject)xin)->name) {
679b51b94faSJed Brown       /* If xin was named, copy the name over to Y. The duplicate names are safe because nobody else will ever see Y. */
6809566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)Y,((PetscObject)xin)->name));
681b51b94faSJed Brown     }
6829566063dSJacob Faibussowitsch     PetscCall(VecCopy(xin,Y));
683a8f87f1dSPatrick Sanan     {
684a8f87f1dSPatrick Sanan       PetscObject dmvtk;
685a8f87f1dSPatrick Sanan       PetscBool   compatible,compatibleSet;
6869566063dSJacob Faibussowitsch       PetscCall(PetscViewerVTKGetDM(viewer,&dmvtk));
687a8f87f1dSPatrick Sanan       if (dmvtk) {
688064a246eSJacob Faibussowitsch         PetscValidHeaderSpecific((DM)dmvtk,DM_CLASSID,2);
6899566063dSJacob Faibussowitsch         PetscCall(DMGetCompatibility(da,(DM)dmvtk,&compatible,&compatibleSet));
6907a8be351SBarry 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.");
691a8f87f1dSPatrick Sanan       }
6929566063dSJacob Faibussowitsch       PetscCall(PetscViewerVTKAddField(viewer,(PetscObject)da,DMDAVTKWriteAll,PETSC_DEFAULT,PETSC_VTK_POINT_FIELD,PETSC_FALSE,(PetscObject)Y));
693a8f87f1dSPatrick Sanan     }
69447c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
69547c6ae99SBarry Smith   } else if (ishdf5) {
6969566063dSJacob Faibussowitsch     PetscCall(VecView_MPI_HDF5_DA(xin,viewer));
69747c6ae99SBarry Smith #endif
6988135c375SStefano Zampini   } else if (isglvis) {
6999566063dSJacob Faibussowitsch     PetscCall(VecView_GLVis(xin,viewer));
70047c6ae99SBarry Smith   } else {
70147c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
70247c6ae99SBarry Smith     PetscBool isbinary,isMPIIO;
70347c6ae99SBarry Smith 
7049566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary));
70547c6ae99SBarry Smith     if (isbinary) {
7069566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO));
70747c6ae99SBarry Smith       if (isMPIIO) {
7089566063dSJacob Faibussowitsch         PetscCall(DMDAArrayMPIIO(da,viewer,xin,PETSC_TRUE));
70947c6ae99SBarry Smith         PetscFunctionReturn(0);
71047c6ae99SBarry Smith       }
71147c6ae99SBarry Smith     }
71247c6ae99SBarry Smith #endif
71347c6ae99SBarry Smith 
71447c6ae99SBarry Smith     /* call viewer on natural ordering */
7159566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix));
7169566063dSJacob Faibussowitsch     PetscCall(DMDACreateNaturalVector(da,&natural));
7179566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)natural,prefix));
7189566063dSJacob Faibussowitsch     PetscCall(DMDAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural));
7199566063dSJacob Faibussowitsch     PetscCall(DMDAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural));
7209566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)xin,&name));
7219566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)natural,name));
722a261c58fSBarry Smith 
7239566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer,&format));
724a261c58fSBarry Smith     if (format == PETSC_VIEWER_BINARY_MATLAB) {
725a261c58fSBarry Smith       /* temporarily remove viewer format so it won't trigger in the VecView() */
7269566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT));
727a261c58fSBarry Smith     }
728a261c58fSBarry Smith 
72923f88b65SBarry Smith     ((PetscObject)natural)->donotPetscObjectPrintClassNamePrefixType = PETSC_TRUE;
7309566063dSJacob Faibussowitsch     PetscCall(VecView(natural,viewer));
73123f88b65SBarry Smith     ((PetscObject)natural)->donotPetscObjectPrintClassNamePrefixType = PETSC_FALSE;
732a261c58fSBarry Smith 
733a261c58fSBarry Smith     if (format == PETSC_VIEWER_BINARY_MATLAB) {
734a261c58fSBarry Smith       MPI_Comm    comm;
735a261c58fSBarry Smith       FILE        *info;
736a261c58fSBarry Smith       const char  *fieldname;
737da88d4d4SJed Brown       char        fieldbuf[256];
738a261c58fSBarry Smith       PetscInt    dim,ni,nj,nk,pi,pj,pk,dof,n;
739a261c58fSBarry Smith 
740a261c58fSBarry Smith       /* set the viewer format back into the viewer */
7419566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(viewer));
7429566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetComm((PetscObject)viewer,&comm));
7439566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryGetInfoPointer(viewer,&info));
7449566063dSJacob Faibussowitsch       PetscCall(DMDAGetInfo(da,&dim,&ni,&nj,&nk,&pi,&pj,&pk,&dof,NULL,NULL,NULL,NULL,NULL));
7459566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n"));
7469566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm,info,"#$$ tmp = PetscBinaryRead(fd); \n"));
7479566063dSJacob Faibussowitsch       if (dim == 1) { PetscCall(PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d);\n",dof,ni)); }
7489566063dSJacob Faibussowitsch       if (dim == 2) { PetscCall(PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d);\n",dof,ni,nj)); }
7499566063dSJacob Faibussowitsch       if (dim == 3) { PetscCall(PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d,%d);\n",dof,ni,nj,nk)); }
750a261c58fSBarry Smith 
751a261c58fSBarry Smith       for (n=0; n<dof; n++) {
7529566063dSJacob Faibussowitsch         PetscCall(DMDAGetFieldName(da,n,&fieldname));
753da88d4d4SJed Brown         if (!fieldname || !fieldname[0]) {
7549566063dSJacob Faibussowitsch           PetscCall(PetscSNPrintf(fieldbuf,sizeof fieldbuf,"field%D",n));
755da88d4d4SJed Brown           fieldname = fieldbuf;
756a261c58fSBarry Smith         }
7579566063dSJacob Faibussowitsch         if (dim == 1) { PetscCall(PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:))';\n",name,fieldname,n+1)); }
7589566063dSJacob Faibussowitsch         if (dim == 2) { PetscCall(PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:,:))';\n",name,fieldname,n+1)); }
7599566063dSJacob Faibussowitsch         if (dim == 3) { PetscCall(PetscFPrintf(comm,info,"#$$ Set.%s.%s = permute(squeeze(tmp(%d,:,:,:)),[2 1 3]);\n",name,fieldname,n+1));}
760a261c58fSBarry Smith       }
7619566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm,info,"#$$ clear tmp; \n"));
7629566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n"));
763a261c58fSBarry Smith     }
764a261c58fSBarry Smith 
7659566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&natural));
76647c6ae99SBarry Smith   }
76747c6ae99SBarry Smith   PetscFunctionReturn(0);
76847c6ae99SBarry Smith }
76947c6ae99SBarry Smith 
77047c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
77147c6ae99SBarry Smith PetscErrorCode VecLoad_HDF5_DA(Vec xin, PetscViewer viewer)
77247c6ae99SBarry Smith {
7732d95be3dSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7749a42bb27SBarry Smith   DM             da;
77515b861d2SVaclav Hapla   int            dim,rdim;
776ec7ae49cSHåkon Strandenes   hsize_t        dims[6]={0},count[6]={0},offset[6]={0};
777d7dd068bSVaclav Hapla   PetscBool      dim2=PETSC_FALSE,timestepping=PETSC_FALSE;
778d7dd068bSVaclav Hapla   PetscInt       dimension,timestep=PETSC_MIN_INT,dofInd;
77947c6ae99SBarry Smith   PetscScalar    *x;
78047c6ae99SBarry Smith   const char     *vecname;
78147c6ae99SBarry Smith   hid_t          filespace; /* file dataspace identifier */
78247c6ae99SBarry Smith   hid_t          dset_id;   /* dataset identifier */
78347c6ae99SBarry Smith   hid_t          memspace;  /* memory dataspace identifier */
784bfd264e7SBarry Smith   hid_t          file_id,group;
7857bcbaff4SBarry Smith   hid_t          scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
7869c7c4993SBarry Smith   DM_DA          *dd;
78747c6ae99SBarry Smith 
78847c6ae99SBarry Smith   PetscFunctionBegin;
7897bcbaff4SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
7907bcbaff4SBarry Smith   scalartype = H5T_NATIVE_FLOAT;
7917bcbaff4SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
7927bcbaff4SBarry Smith #error "HDF5 output with 128 bit floats not supported."
793570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
794570b7f6dSBarry Smith #error "HDF5 output with 16 bit floats not supported."
7957bcbaff4SBarry Smith #else
7967bcbaff4SBarry Smith   scalartype = H5T_NATIVE_DOUBLE;
7977bcbaff4SBarry Smith #endif
7987bcbaff4SBarry Smith 
7999566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5OpenGroup(viewer, &file_id, &group));
8009566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)xin,&vecname));
8019566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckTimestepping_Internal(viewer, vecname));
8029566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5IsTimestepping(viewer, &timestepping));
803d7dd068bSVaclav Hapla   if (timestepping) {
8049566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5GetTimestep(viewer, &timestep));
805d7dd068bSVaclav Hapla   }
8069566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin,&da));
8079c7c4993SBarry Smith   dd   = (DM_DA*)da->data;
8089566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(da, &dimension));
80947c6ae99SBarry Smith 
810ec7ae49cSHåkon Strandenes   /* Open dataset */
811729ab6d9SBarry Smith   PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT));
81247c6ae99SBarry Smith 
813ec7ae49cSHåkon Strandenes   /* Retrieve the dataspace for the dataset */
814ec7ae49cSHåkon Strandenes   PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id));
815ec7ae49cSHåkon Strandenes   PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(filespace, dims, NULL));
816ec7ae49cSHåkon Strandenes 
817ec7ae49cSHåkon Strandenes   /* Expected dimension for holding the dof's */
81847c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
819ec7ae49cSHåkon Strandenes   dofInd = rdim-2;
820ec7ae49cSHåkon Strandenes #else
821ec7ae49cSHåkon Strandenes   dofInd = rdim-1;
82247c6ae99SBarry Smith #endif
823ec7ae49cSHåkon Strandenes 
824ec7ae49cSHåkon Strandenes   /* The expected number of dimensions, assuming basedimension2 = false */
825ec7ae49cSHåkon Strandenes   dim = dimension;
826ec7ae49cSHåkon Strandenes   if (dd->w > 1) ++dim;
827ec7ae49cSHåkon Strandenes   if (timestep >= 0) ++dim;
82847c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
829ec7ae49cSHåkon Strandenes   ++dim;
83047c6ae99SBarry Smith #endif
831ec7ae49cSHåkon Strandenes 
832ec7ae49cSHåkon Strandenes   /* In this case the input dataset have one extra, unexpected dimension. */
833ec7ae49cSHåkon Strandenes   if (rdim == dim+1) {
834ec7ae49cSHåkon Strandenes     /* In this case the block size unity */
835ec7ae49cSHåkon Strandenes     if (dd->w == 1 && dims[dofInd] == 1) dim2 = PETSC_TRUE;
836ec7ae49cSHåkon Strandenes 
837ec7ae49cSHåkon Strandenes     /* Special error message for the case where dof does not match the input file */
83808401ef6SPierre 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);
839ec7ae49cSHåkon Strandenes 
840ec7ae49cSHåkon Strandenes   /* Other cases where rdim != dim cannot be handled currently */
84108401ef6SPierre 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);
842ec7ae49cSHåkon Strandenes 
843ec7ae49cSHåkon Strandenes   /* Set up the hyperslab size */
844ec7ae49cSHåkon Strandenes   dim = 0;
845ec7ae49cSHåkon Strandenes   if (timestep >= 0) {
846ec7ae49cSHåkon Strandenes     offset[dim] = timestep;
847ec7ae49cSHåkon Strandenes     count[dim] = 1;
848ec7ae49cSHåkon Strandenes     ++dim;
849ec7ae49cSHåkon Strandenes   }
850ec7ae49cSHåkon Strandenes   if (dimension == 3) {
8519566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(dd->zs,offset + dim));
8529566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(dd->ze - dd->zs,count + dim));
853ec7ae49cSHåkon Strandenes     ++dim;
854ec7ae49cSHåkon Strandenes   }
855ec7ae49cSHåkon Strandenes   if (dimension > 1) {
8569566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(dd->ys,offset + dim));
8579566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(dd->ye - dd->ys,count + dim));
858ec7ae49cSHåkon Strandenes     ++dim;
859ec7ae49cSHåkon Strandenes   }
8609566063dSJacob Faibussowitsch   PetscCall(PetscHDF5IntCast(dd->xs/dd->w,offset + dim));
8619566063dSJacob Faibussowitsch   PetscCall(PetscHDF5IntCast((dd->xe - dd->xs)/dd->w,count + dim));
862ec7ae49cSHåkon Strandenes   ++dim;
863ec7ae49cSHåkon Strandenes   if (dd->w > 1 || dim2) {
864ec7ae49cSHåkon Strandenes     offset[dim] = 0;
8659566063dSJacob Faibussowitsch     PetscCall(PetscHDF5IntCast(dd->w,count + dim));
866ec7ae49cSHåkon Strandenes     ++dim;
867ec7ae49cSHåkon Strandenes   }
868ec7ae49cSHåkon Strandenes #if defined(PETSC_USE_COMPLEX)
869ec7ae49cSHåkon Strandenes   offset[dim] = 0;
870ec7ae49cSHåkon Strandenes   count[dim] = 2;
871ec7ae49cSHåkon Strandenes   ++dim;
872ec7ae49cSHåkon Strandenes #endif
873ec7ae49cSHåkon Strandenes 
874ec7ae49cSHåkon Strandenes   /* Create the memory and filespace */
875729ab6d9SBarry Smith   PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL));
876729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL));
87747c6ae99SBarry Smith 
8789566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xin, &x));
8792d95be3dSVaclav Hapla   PetscStackCallHDF5(H5Dread,(dset_id, scalartype, memspace, filespace, hdf5->dxpl_id, x));
8809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xin, &x));
88147c6ae99SBarry Smith 
88247c6ae99SBarry Smith   /* Close/release resources */
883bfd264e7SBarry Smith   if (group != file_id) {
884729ab6d9SBarry Smith     PetscStackCallHDF5(H5Gclose,(group));
885bfd264e7SBarry Smith   }
886729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
887729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(memspace));
888729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dset_id));
88947c6ae99SBarry Smith   PetscFunctionReturn(0);
89047c6ae99SBarry Smith }
89147c6ae99SBarry Smith #endif
89247c6ae99SBarry Smith 
89347c6ae99SBarry Smith PetscErrorCode VecLoad_Binary_DA(Vec xin, PetscViewer viewer)
89447c6ae99SBarry Smith {
8959a42bb27SBarry Smith   DM             da;
89647c6ae99SBarry Smith   Vec            natural;
89747c6ae99SBarry Smith   const char     *prefix;
89847c6ae99SBarry Smith   PetscInt       bs;
89947c6ae99SBarry Smith   PetscBool      flag;
90047c6ae99SBarry Smith   DM_DA          *dd;
90147c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
90247c6ae99SBarry Smith   PetscBool isMPIIO;
90347c6ae99SBarry Smith #endif
90447c6ae99SBarry Smith 
90547c6ae99SBarry Smith   PetscFunctionBegin;
9069566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin,&da));
90747c6ae99SBarry Smith   dd   = (DM_DA*)da->data;
90847c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
9099566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO));
91047c6ae99SBarry Smith   if (isMPIIO) {
9119566063dSJacob Faibussowitsch     PetscCall(DMDAArrayMPIIO(da,viewer,xin,PETSC_FALSE));
91247c6ae99SBarry Smith     PetscFunctionReturn(0);
91347c6ae99SBarry Smith   }
91447c6ae99SBarry Smith #endif
91547c6ae99SBarry Smith 
9169566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix));
9179566063dSJacob Faibussowitsch   PetscCall(DMDACreateNaturalVector(da,&natural));
9189566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name));
9199566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)natural,prefix));
9209566063dSJacob Faibussowitsch   PetscCall(VecLoad(natural,viewer));
9219566063dSJacob Faibussowitsch   PetscCall(DMDANaturalToGlobalBegin(da,natural,INSERT_VALUES,xin));
9229566063dSJacob Faibussowitsch   PetscCall(DMDANaturalToGlobalEnd(da,natural,INSERT_VALUES,xin));
9239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&natural));
9249566063dSJacob Faibussowitsch   PetscCall(PetscInfo(xin,"Loading vector from natural ordering into DMDA\n"));
9259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,((PetscObject)xin)->prefix,"-vecload_block_size",&bs,&flag));
92647c6ae99SBarry Smith   if (flag && bs != dd->w) {
9279566063dSJacob Faibussowitsch     PetscCall(PetscInfo(xin,"Block size in file %D not equal to DMDA's dof %D\n",bs,dd->w));
92847c6ae99SBarry Smith   }
92947c6ae99SBarry Smith   PetscFunctionReturn(0);
93047c6ae99SBarry Smith }
93147c6ae99SBarry Smith 
9327087cfbeSBarry Smith PetscErrorCode  VecLoad_Default_DA(Vec xin, PetscViewer viewer)
93347c6ae99SBarry Smith {
9349a42bb27SBarry Smith   DM             da;
93547c6ae99SBarry Smith   PetscBool      isbinary;
93647c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
93747c6ae99SBarry Smith   PetscBool ishdf5;
93847c6ae99SBarry Smith #endif
93947c6ae99SBarry Smith 
94047c6ae99SBarry Smith   PetscFunctionBegin;
9419566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin,&da));
9427a8be351SBarry Smith   PetscCheck(da,PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
94347c6ae99SBarry Smith 
94447c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
9459566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5));
94647c6ae99SBarry Smith #endif
9479566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary));
94847c6ae99SBarry Smith 
94947c6ae99SBarry Smith   if (isbinary) {
9509566063dSJacob Faibussowitsch     PetscCall(VecLoad_Binary_DA(xin,viewer));
95147c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
95247c6ae99SBarry Smith   } else if (ishdf5) {
9539566063dSJacob Faibussowitsch     PetscCall(VecLoad_HDF5_DA(xin,viewer));
95447c6ae99SBarry Smith #endif
95598921bdaSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for vector loading", ((PetscObject)viewer)->type_name);
95647c6ae99SBarry Smith   PetscFunctionReturn(0);
95747c6ae99SBarry Smith }
958