xref: /petsc/src/dm/impls/da/gr2.c (revision 7fe7c8fe9311cf7183c55e2475b3d9594c67015b)
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*/
7af0996ceSBarry Smith #include <petsc/private/vecimpl.h>
89804daf3SBarry Smith #include <petscdraw.h>
9665c2dedSJed Brown #include <petscviewerhdf5.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;
19*7fe7c8feSLisandro 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 #undef __FUNCT__
2947c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_Draw_DA2d_Zoom"
3047c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA2d_Zoom(PetscDraw draw,void *ctx)
3147c6ae99SBarry Smith {
3247c6ae99SBarry Smith   ZoomCtx           *zctx = (ZoomCtx*)ctx;
3347c6ae99SBarry Smith   PetscErrorCode    ierr;
34e5ab1681SLisandro Dalcin   PetscInt          m,n,i,j,k,dof,id,c1,c2,c3,c4;
35e5ab1681SLisandro Dalcin   PetscReal         min,max,x1,x2,x3,x4,y_1,y2,y3,y4;
36e5ab1681SLisandro Dalcin   const PetscScalar *xy,*v;
3747c6ae99SBarry Smith 
3847c6ae99SBarry Smith   PetscFunctionBegin;
3947c6ae99SBarry Smith   m    = zctx->m;
4047c6ae99SBarry Smith   n    = zctx->n;
41e5ab1681SLisandro Dalcin   dof  = zctx->dof;
4247c6ae99SBarry Smith   k    = zctx->k;
4347c6ae99SBarry Smith   xy   = zctx->xy;
44e5ab1681SLisandro Dalcin   v    = zctx->v;
4547c6ae99SBarry Smith   min  = zctx->min;
46f3f0eb19SBarry Smith   max  = zctx->max;
4747c6ae99SBarry Smith 
4847c6ae99SBarry Smith   /* PetscDraw the contour plot patch */
49e5ab1681SLisandro Dalcin   ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
5047c6ae99SBarry Smith   for (j=0; j<n-1; j++) {
5147c6ae99SBarry Smith     for (i=0; i<m-1; i++) {
520e4fe250SBarry Smith       id   = i+j*m;
530e4fe250SBarry Smith       x1   = PetscRealPart(xy[2*id]);
540e4fe250SBarry Smith       y_1  = PetscRealPart(xy[2*id+1]);
55e5ab1681SLisandro Dalcin       c1   = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max);
560e4fe250SBarry Smith 
570e4fe250SBarry Smith       id   = i+j*m+1;
580e4fe250SBarry Smith       x2   = PetscRealPart(xy[2*id]);
59a39a4c7dSLisandro Dalcin       y2   = PetscRealPart(xy[2*id+1]);
60e5ab1681SLisandro Dalcin       c2   = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max);
610e4fe250SBarry Smith 
620e4fe250SBarry Smith       id   = i+j*m+1+m;
63a39a4c7dSLisandro Dalcin       x3   = PetscRealPart(xy[2*id]);
640e4fe250SBarry Smith       y3   = PetscRealPart(xy[2*id+1]);
65e5ab1681SLisandro Dalcin       c3   = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max);
660e4fe250SBarry Smith 
670e4fe250SBarry Smith       id   = i+j*m+m;
68a39a4c7dSLisandro Dalcin       x4   = PetscRealPart(xy[2*id]);
69a39a4c7dSLisandro Dalcin       y4   = PetscRealPart(xy[2*id+1]);
70e5ab1681SLisandro Dalcin       c4   = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max);
71f3f0eb19SBarry Smith 
7247c6ae99SBarry Smith       ierr = PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);CHKERRQ(ierr);
7347c6ae99SBarry Smith       ierr = PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4);CHKERRQ(ierr);
7447c6ae99SBarry Smith       if (zctx->showgrid) {
7547c6ae99SBarry Smith         ierr = PetscDrawLine(draw,x1,y_1,x2,y2,PETSC_DRAW_BLACK);CHKERRQ(ierr);
7647c6ae99SBarry Smith         ierr = PetscDrawLine(draw,x2,y2,x3,y3,PETSC_DRAW_BLACK);CHKERRQ(ierr);
7747c6ae99SBarry Smith         ierr = PetscDrawLine(draw,x3,y3,x4,y4,PETSC_DRAW_BLACK);CHKERRQ(ierr);
7847c6ae99SBarry Smith         ierr = PetscDrawLine(draw,x4,y4,x1,y_1,PETSC_DRAW_BLACK);CHKERRQ(ierr);
7947c6ae99SBarry Smith       }
8047c6ae99SBarry Smith     }
8147c6ae99SBarry Smith   }
82*7fe7c8feSLisandro Dalcin   if (zctx->showaxis && !zctx->rank) {
83e5ab1681SLisandro Dalcin     if (zctx->name0 || zctx->name1) {
84109c9344SBarry Smith       PetscReal xl,yl,xr,yr,x,y;
85109c9344SBarry Smith       ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
86e5ab1681SLisandro Dalcin       x  = xl + .30*(xr - xl);
87109c9344SBarry Smith       xl = xl + .01*(xr - xl);
88e5ab1681SLisandro Dalcin       y  = yr - .30*(yr - yl);
89109c9344SBarry Smith       yl = yl + .01*(yr - yl);
90e5ab1681SLisandro Dalcin       if (zctx->name0) {ierr = PetscDrawString(draw,x,yl,PETSC_DRAW_BLACK,zctx->name0);CHKERRQ(ierr);}
91e5ab1681SLisandro Dalcin       if (zctx->name1) {ierr = PetscDrawStringVertical(draw,xl,y,PETSC_DRAW_BLACK,zctx->name1);CHKERRQ(ierr);}
92109c9344SBarry Smith     }
930e4fe250SBarry Smith     /*
940e4fe250SBarry Smith        Ideally we would use the PetscDrawAxis object to manage displaying the coordinate limits
950e4fe250SBarry Smith        but that may require some refactoring.
960e4fe250SBarry Smith     */
97e5ab1681SLisandro Dalcin     {
98*7fe7c8feSLisandro Dalcin       double xmin = (double)zctx->xmin, ymin = (double)zctx->ymin;
99*7fe7c8feSLisandro Dalcin       double xmax = (double)zctx->xmax, ymax = (double)zctx->ymax;
100e5ab1681SLisandro Dalcin       char   value[16]; size_t len; PetscReal w;
101*7fe7c8feSLisandro Dalcin       ierr = PetscSNPrintf(value,16,"%0.2e",xmin);CHKERRQ(ierr);
102e5ab1681SLisandro Dalcin       ierr = PetscDrawString(draw,xmin,ymin - .05*(ymax - ymin),PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
103*7fe7c8feSLisandro Dalcin       ierr = PetscSNPrintf(value,16,"%0.2e",xmax);CHKERRQ(ierr);
1040e4fe250SBarry Smith       ierr = PetscStrlen(value,&len);CHKERRQ(ierr);
1050298fd71SBarry Smith       ierr = PetscDrawStringGetSize(draw,&w,NULL);CHKERRQ(ierr);
106e5ab1681SLisandro Dalcin       ierr = PetscDrawString(draw,xmax - len*w,ymin - .05*(ymax - ymin),PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
107*7fe7c8feSLisandro Dalcin       ierr = PetscSNPrintf(value,16,"%0.2e",ymin);CHKERRQ(ierr);
108e5ab1681SLisandro Dalcin       ierr = PetscDrawString(draw,xmin - .05*(xmax - xmin),ymin,PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
109*7fe7c8feSLisandro Dalcin       ierr = PetscSNPrintf(value,16,"%0.2e",ymax);CHKERRQ(ierr);
110e5ab1681SLisandro Dalcin       ierr = PetscDrawString(draw,xmin - .05*(xmax - xmin),ymax,PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
111e5ab1681SLisandro Dalcin     }
112e5ab1681SLisandro Dalcin   }
113e5ab1681SLisandro Dalcin   ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
11447c6ae99SBarry Smith   PetscFunctionReturn(0);
11547c6ae99SBarry Smith }
11647c6ae99SBarry Smith 
11747c6ae99SBarry Smith #undef __FUNCT__
11847c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_Draw_DA2d"
11947c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA2d(Vec xin,PetscViewer viewer)
12047c6ae99SBarry Smith {
1219a42bb27SBarry Smith   DM                 da,dac,dag;
12247c6ae99SBarry Smith   PetscErrorCode     ierr;
123a74ba6f7SBarry Smith   PetscInt           N,s,M,w,ncoors = 4;
12447c6ae99SBarry Smith   const PetscInt     *lx,*ly;
125e5ab1681SLisandro Dalcin   PetscReal          coors[4];
12647c6ae99SBarry Smith   PetscDraw          draw,popup;
12747c6ae99SBarry Smith   PetscBool          isnull,useports = PETSC_FALSE;
12847c6ae99SBarry Smith   MPI_Comm           comm;
12947c6ae99SBarry Smith   Vec                xlocal,xcoor,xcoorl;
130bff4a2f0SMatthew G. Knepley   DMBoundaryType     bx,by;
131aa219208SBarry Smith   DMDAStencilType    st;
13247c6ae99SBarry Smith   ZoomCtx            zctx;
1330298fd71SBarry Smith   PetscDrawViewPorts *ports = NULL;
13447c6ae99SBarry Smith   PetscViewerFormat  format;
13520d0051dSBarry Smith   PetscInt           *displayfields;
13667dd0837SBarry Smith   PetscInt           ndisplayfields,i,nbounds;
13767dd0837SBarry Smith   const PetscReal    *bounds;
13847c6ae99SBarry Smith 
13947c6ae99SBarry Smith   PetscFunctionBegin;
14047c6ae99SBarry Smith   zctx.showgrid = PETSC_FALSE;
141*7fe7c8feSLisandro Dalcin   zctx.showaxis = PETSC_TRUE;
1428865f1eaSKarl Rupp 
14347c6ae99SBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1445b399a63SLisandro Dalcin   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
1455b399a63SLisandro Dalcin   if (isnull) PetscFunctionReturn(0);
1465b399a63SLisandro Dalcin 
14703193ff8SBarry Smith   ierr = PetscViewerDrawGetBounds(viewer,&nbounds,&bounds);CHKERRQ(ierr);
14847c6ae99SBarry Smith 
149c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
150ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
15147c6ae99SBarry Smith 
15247c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr);
153e5ab1681SLisandro Dalcin   ierr = MPI_Comm_rank(comm,&zctx.rank);CHKERRQ(ierr);
15447c6ae99SBarry Smith 
1551321219cSEthan Coon   ierr = DMDAGetInfo(da,0,&M,&N,0,&zctx.m,&zctx.n,0,&w,&s,&bx,&by,0,&st);CHKERRQ(ierr);
1560298fd71SBarry Smith   ierr = DMDAGetOwnershipRanges(da,&lx,&ly,NULL);CHKERRQ(ierr);
15747c6ae99SBarry Smith 
15847c6ae99SBarry Smith   /*
15947c6ae99SBarry Smith      Obtain a sequential vector that is going to contain the local values plus ONE layer of
160aa219208SBarry Smith      ghosted values to draw the graphics from. We also need its corresponding DMDA (dac) that will
16147c6ae99SBarry Smith      update the local values pluse ONE layer of ghost values.
16247c6ae99SBarry Smith   */
16347c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)da,"GraphicsGhosted",(PetscObject*)&xlocal);CHKERRQ(ierr);
16447c6ae99SBarry Smith   if (!xlocal) {
165bff4a2f0SMatthew G. Knepley     if (bx !=  DM_BOUNDARY_NONE || by !=  DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) {
16647c6ae99SBarry Smith       /*
16747c6ae99SBarry Smith          if original da is not of stencil width one, or periodic or not a box stencil then
168aa219208SBarry Smith          create a special DMDA to handle one level of ghost points for graphics
16947c6ae99SBarry Smith       */
170bff4a2f0SMatthew G. Knepley       ierr = DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,zctx.m,zctx.n,w,1,lx,ly,&dac);CHKERRQ(ierr);
171aa219208SBarry Smith       ierr = PetscInfo(da,"Creating auxilary DMDA for managing graphics ghost points\n");CHKERRQ(ierr);
17247c6ae99SBarry Smith     } else {
17347c6ae99SBarry Smith       /* otherwise we can use the da we already have */
17447c6ae99SBarry Smith       dac = da;
17547c6ae99SBarry Smith     }
17647c6ae99SBarry Smith     /* create local vector for holding ghosted values used in graphics */
177564755cdSBarry Smith     ierr = DMCreateLocalVector(dac,&xlocal);CHKERRQ(ierr);
17847c6ae99SBarry Smith     if (dac != da) {
179aa219208SBarry Smith       /* don't keep any public reference of this DMDA, is is only available through xlocal */
180f7923d8aSBarry Smith       ierr = PetscObjectDereference((PetscObject)dac);CHKERRQ(ierr);
18147c6ae99SBarry Smith     } else {
18247c6ae99SBarry Smith       /* remove association between xlocal and da, because below we compose in the opposite
18347c6ae99SBarry Smith          direction and if we left this connect we'd get a loop, so the objects could
18447c6ae99SBarry Smith          never be destroyed */
185c688c046SMatthew G Knepley       ierr = PetscObjectRemoveReference((PetscObject)xlocal,"__PETSc_dm");CHKERRQ(ierr);
18647c6ae99SBarry Smith     }
18747c6ae99SBarry Smith     ierr = PetscObjectCompose((PetscObject)da,"GraphicsGhosted",(PetscObject)xlocal);CHKERRQ(ierr);
18847c6ae99SBarry Smith     ierr = PetscObjectDereference((PetscObject)xlocal);CHKERRQ(ierr);
18947c6ae99SBarry Smith   } else {
190bff4a2f0SMatthew G. Knepley     if (bx !=  DM_BOUNDARY_NONE || by !=  DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) {
191c688c046SMatthew G Knepley       ierr = VecGetDM(xlocal, &dac);CHKERRQ(ierr);
192f7923d8aSBarry Smith     } else {
193f7923d8aSBarry Smith       dac = da;
19447c6ae99SBarry Smith     }
19547c6ae99SBarry Smith   }
19647c6ae99SBarry Smith 
19747c6ae99SBarry Smith   /*
19847c6ae99SBarry Smith       Get local (ghosted) values of vector
19947c6ae99SBarry Smith   */
2009a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr);
2019a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr);
2025edff71fSBarry Smith   ierr = VecGetArrayRead(xlocal,&zctx.v);CHKERRQ(ierr);
20347c6ae99SBarry Smith 
204832b7cebSLisandro Dalcin   /*
205832b7cebSLisandro Dalcin       Get coordinates of nodes
206832b7cebSLisandro Dalcin   */
2076636e97aSMatthew G Knepley   ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
20847c6ae99SBarry Smith   if (!xcoor) {
209aa219208SBarry Smith     ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,0.0);CHKERRQ(ierr);
2106636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
21147c6ae99SBarry Smith   }
21247c6ae99SBarry Smith 
21347c6ae99SBarry Smith   /*
21447c6ae99SBarry Smith       Determine the min and max coordinates in plot
21547c6ae99SBarry Smith   */
216e5ab1681SLisandro Dalcin   ierr = VecStrideMin(xcoor,0,NULL,&zctx.xmin);CHKERRQ(ierr);
217e5ab1681SLisandro Dalcin   ierr = VecStrideMax(xcoor,0,NULL,&zctx.xmax);CHKERRQ(ierr);
218e5ab1681SLisandro Dalcin   ierr = VecStrideMin(xcoor,1,NULL,&zctx.ymin);CHKERRQ(ierr);
219e5ab1681SLisandro Dalcin   ierr = VecStrideMax(xcoor,1,NULL,&zctx.ymax);CHKERRQ(ierr);
220*7fe7c8feSLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-draw_contour_axis",&zctx.showaxis,NULL);CHKERRQ(ierr);
221*7fe7c8feSLisandro Dalcin   if (zctx.showaxis) {
222*7fe7c8feSLisandro Dalcin     coors[0] = zctx.xmin - .05*(zctx.xmax - zctx.xmin); coors[1] = zctx.ymin - .05*(zctx.ymax - zctx.ymin);
223*7fe7c8feSLisandro Dalcin     coors[2] = zctx.xmax + .05*(zctx.xmax - zctx.xmin); coors[3] = zctx.ymax + .05*(zctx.ymax - zctx.ymin);
224*7fe7c8feSLisandro Dalcin   } else {
225*7fe7c8feSLisandro Dalcin     coors[0] = zctx.xmin; coors[1] = zctx.ymin; coors[2] = zctx.xmax; coors[3] = zctx.ymax;
226*7fe7c8feSLisandro Dalcin   }
227c5929fdfSBarry Smith   ierr = PetscOptionsGetRealArray(NULL,NULL,"-draw_coordinates",coors,&ncoors,NULL);CHKERRQ(ierr);
228e5ab1681SLisandro Dalcin   ierr = PetscInfo4(da,"Preparing DMDA 2d contour plot coordinates %g %g %g %g\n",(double)coors[0],(double)coors[1],(double)coors[2],(double)coors[3]);CHKERRQ(ierr);
229a74ba6f7SBarry Smith 
23047c6ae99SBarry Smith   /*
231832b7cebSLisandro Dalcin       Get local ghosted version of coordinates
23247c6ae99SBarry Smith   */
23347c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr);
23447c6ae99SBarry Smith   if (!xcoorl) {
235aa219208SBarry Smith     /* create DMDA to get local version of graphics */
236bff4a2f0SMatthew G. Knepley     ierr = DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,zctx.m,zctx.n,2,1,lx,ly,&dag);CHKERRQ(ierr);
237aa219208SBarry Smith     ierr = PetscInfo(dag,"Creating auxilary DMDA for managing graphics coordinates ghost points\n");CHKERRQ(ierr);
238564755cdSBarry Smith     ierr = DMCreateLocalVector(dag,&xcoorl);CHKERRQ(ierr);
23947c6ae99SBarry Smith     ierr = PetscObjectCompose((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject)xcoorl);CHKERRQ(ierr);
240f7923d8aSBarry Smith     ierr = PetscObjectDereference((PetscObject)dag);CHKERRQ(ierr);
24147c6ae99SBarry Smith     ierr = PetscObjectDereference((PetscObject)xcoorl);CHKERRQ(ierr);
24247c6ae99SBarry Smith   } else {
243c688c046SMatthew G Knepley     ierr = VecGetDM(xcoorl,&dag);CHKERRQ(ierr);
24447c6ae99SBarry Smith   }
2459a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
2469a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
2475edff71fSBarry Smith   ierr = VecGetArrayRead(xcoorl,&zctx.xy);CHKERRQ(ierr);
248832b7cebSLisandro Dalcin   ierr = DMDAGetCoordinateName(da,0,&zctx.name0);CHKERRQ(ierr);
249832b7cebSLisandro Dalcin   ierr = DMDAGetCoordinateName(da,1,&zctx.name1);CHKERRQ(ierr);
25047c6ae99SBarry Smith 
25147c6ae99SBarry Smith   /*
25247c6ae99SBarry Smith       Get information about size of area each processor must do graphics for
25347c6ae99SBarry Smith   */
254e5ab1681SLisandro Dalcin   ierr = DMDAGetInfo(dac,NULL,&M,&N,NULL,NULL,NULL,NULL,&zctx.dof,NULL,&bx,&by,NULL,NULL);CHKERRQ(ierr);
255e5ab1681SLisandro Dalcin   ierr = DMDAGetGhostCorners(dac,NULL,NULL,NULL,&zctx.m,&zctx.n,NULL);CHKERRQ(ierr);
256c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-draw_contour_grid",&zctx.showgrid,NULL);CHKERRQ(ierr);
2574e6118eeSBarry Smith 
2584e6118eeSBarry Smith   ierr = DMDASelectFields(da,&ndisplayfields,&displayfields);CHKERRQ(ierr);
25947c6ae99SBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
260c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-draw_ports",&useports,NULL);CHKERRQ(ierr);
261832b7cebSLisandro Dalcin   if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE;
262832b7cebSLisandro Dalcin   if (useports) {
263e5ab1681SLisandro Dalcin     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
2645b399a63SLisandro Dalcin     ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
2655b399a63SLisandro Dalcin     ierr = PetscDrawClear(draw);CHKERRQ(ierr);
26620d0051dSBarry Smith     ierr = PetscDrawViewPortsCreate(draw,ndisplayfields,&ports);CHKERRQ(ierr);
267e5ab1681SLisandro Dalcin   }
26820d0051dSBarry Smith 
269832b7cebSLisandro Dalcin   /*
270832b7cebSLisandro Dalcin       Loop over each field; drawing each in a different window
271832b7cebSLisandro Dalcin   */
27220d0051dSBarry Smith   for (i=0; i<ndisplayfields; i++) {
27320d0051dSBarry Smith     zctx.k = displayfields[i];
2745b399a63SLisandro Dalcin 
275e5ab1681SLisandro Dalcin     /* determine the min and max value in plot */
2760298fd71SBarry Smith     ierr = VecStrideMin(xin,zctx.k,NULL,&zctx.min);CHKERRQ(ierr);
2770298fd71SBarry Smith     ierr = VecStrideMax(xin,zctx.k,NULL,&zctx.max);CHKERRQ(ierr);
27867dd0837SBarry Smith     if (zctx.k < nbounds) {
279f3f0eb19SBarry Smith       zctx.min = bounds[2*zctx.k];
280f3f0eb19SBarry Smith       zctx.max = bounds[2*zctx.k+1];
28167dd0837SBarry Smith     }
28247c6ae99SBarry Smith     if (zctx.min == zctx.max) {
28347c6ae99SBarry Smith       zctx.min -= 1.e-12;
28447c6ae99SBarry Smith       zctx.max += 1.e-12;
28547c6ae99SBarry Smith     }
28657622a8eSBarry Smith     ierr = PetscInfo2(da,"DMDA 2d contour plot min %g max %g\n",(double)zctx.min,(double)zctx.max);CHKERRQ(ierr);
28747c6ae99SBarry Smith 
288832b7cebSLisandro Dalcin     if (useports) {
289832b7cebSLisandro Dalcin       ierr = PetscDrawViewPortsSet(ports,i);CHKERRQ(ierr);
290832b7cebSLisandro Dalcin     } else {
291832b7cebSLisandro Dalcin       const char *title;
292832b7cebSLisandro Dalcin       ierr = PetscViewerDrawGetDraw(viewer,i,&draw);CHKERRQ(ierr);
293832b7cebSLisandro Dalcin       ierr = DMDAGetFieldName(da,zctx.k,&title);CHKERRQ(ierr);
294832b7cebSLisandro Dalcin       if (title) {ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr);}
295832b7cebSLisandro Dalcin     }
296832b7cebSLisandro Dalcin 
29747c6ae99SBarry Smith     ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
29847c6ae99SBarry Smith     if (popup) {ierr = PetscDrawScalePopup(popup,zctx.min,zctx.max);CHKERRQ(ierr);}
299e5ab1681SLisandro Dalcin     ierr = PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr);
30047c6ae99SBarry Smith     ierr = PetscDrawZoom(draw,VecView_MPI_Draw_DA2d_Zoom,&zctx);CHKERRQ(ierr);
301832b7cebSLisandro Dalcin     if (!useports) {ierr = PetscDrawSave(draw);CHKERRQ(ierr);}
30247c6ae99SBarry Smith   }
303832b7cebSLisandro Dalcin   if (useports) {
304832b7cebSLisandro Dalcin     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
305832b7cebSLisandro Dalcin     ierr = PetscDrawSave(draw);CHKERRQ(ierr);
306832b7cebSLisandro Dalcin   }
307832b7cebSLisandro Dalcin 
3086bf464f9SBarry Smith   ierr = PetscDrawViewPortsDestroy(ports);CHKERRQ(ierr);
309e5ab1681SLisandro Dalcin   ierr = PetscFree(displayfields);CHKERRQ(ierr);
3105edff71fSBarry Smith   ierr = VecRestoreArrayRead(xcoorl,&zctx.xy);CHKERRQ(ierr);
3115edff71fSBarry Smith   ierr = VecRestoreArrayRead(xlocal,&zctx.v);CHKERRQ(ierr);
31247c6ae99SBarry Smith   PetscFunctionReturn(0);
31347c6ae99SBarry Smith }
31447c6ae99SBarry Smith 
3150da24e51SJuha Jäykkä #if defined(PETSC_HAVE_HDF5)
3160da24e51SJuha Jäykkä #undef __FUNCT__
3170da24e51SJuha Jäykkä #define __FUNCT__ "VecGetHDF5ChunkSize"
318c73cfb54SMatthew G. Knepley static PetscErrorCode VecGetHDF5ChunkSize(DM_DA *da, Vec xin, PetscInt dimension, PetscInt timestep, hsize_t *chunkDims)
3190da24e51SJuha Jäykkä {
320d4190030SJed Brown   PetscMPIInt    comm_size;
321d4190030SJed Brown   PetscErrorCode ierr;
322561a82e7SJed Brown   hsize_t        chunk_size, target_size, dim;
323561a82e7SJed Brown   hsize_t        vec_size = sizeof(PetscScalar)*da->M*da->N*da->P*da->w;
324561a82e7SJed Brown   hsize_t        avg_local_vec_size,KiB = 1024,MiB = KiB*KiB,GiB = MiB*KiB,min_size = MiB;
325561a82e7SJed Brown   hsize_t        max_chunks = 64*KiB;                                              /* HDF5 internal limitation */
326561a82e7SJed Brown   hsize_t        max_chunk_size = 4*GiB;                                           /* HDF5 internal limitation */
32784179ffaSJed Brown   PetscInt       zslices=da->p, yslices=da->n, xslices=da->m;
3280da24e51SJuha Jäykkä 
329cb5c4f94SJuha Jäykkä   PetscFunctionBegin;
330d4190030SJed Brown   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)xin), &comm_size);CHKERRQ(ierr);
3310da24e51SJuha Jäykkä   avg_local_vec_size = (hsize_t) ceil(vec_size*1.0/comm_size);      /* we will attempt to use this as the chunk size */
3320da24e51SJuha Jäykkä 
333794a961bSBarry Smith   target_size = (hsize_t) ceil(PetscMin(vec_size,PetscMin(max_chunk_size,PetscMax(avg_local_vec_size,PetscMax(ceil(vec_size*1.0/max_chunks),min_size)))));
3347d310018SBarry 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 */
3357d310018SBarry 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);
3360da24e51SJuha Jäykkä 
337cb5c4f94SJuha Jäykkä   /*
338cb5c4f94SJuha Jäykkä    if size/rank > max_chunk_size, we need radical measures: even going down to
339cb5c4f94SJuha Jäykkä    avg_local_vec_size is not enough, so we simply use chunk size of 4 GiB no matter
340cb5c4f94SJuha Jäykkä    what, composed in the most efficient way possible.
341cb5c4f94SJuha Jäykkä    N.B. this minimises the number of chunks, which may or may not be the optimal
342cb5c4f94SJuha Jäykkä    solution. In a BG, for example, the optimal solution is probably to make # chunks = #
343cb5c4f94SJuha Jäykkä    IO nodes involved, but this author has no access to a BG to figure out how to
344cb5c4f94SJuha Jäykkä    reliably find the right number. And even then it may or may not be enough.
345cb5c4f94SJuha Jäykkä    */
3460da24e51SJuha Jäykkä   if (avg_local_vec_size > max_chunk_size) {
347cb5c4f94SJuha Jäykkä     /* check if we can just split local z-axis: is that enough? */
34884179ffaSJed Brown     zslices = (PetscInt)ceil(vec_size*1.0/(da->p*max_chunk_size))*zslices;
3490da24e51SJuha Jäykkä     if (zslices > da->P) {
350cb5c4f94SJuha Jäykkä       /* lattice is too large in xy-directions, splitting z only is not enough */
3510da24e51SJuha Jäykkä       zslices = da->P;
35284179ffaSJed Brown       yslices= (PetscInt)ceil(vec_size*1.0/(zslices*da->n*max_chunk_size))*yslices;
3530da24e51SJuha Jäykkä       if (yslices > da->N) {
354cb5c4f94SJuha Jäykkä 	/* lattice is too large in x-direction, splitting along z, y is not enough */
3550da24e51SJuha Jäykkä 	yslices = da->N;
35684179ffaSJed Brown 	xslices= (PetscInt)ceil(vec_size*1.0/(zslices*yslices*da->m*max_chunk_size))*xslices;
3570da24e51SJuha Jäykkä       }
3580da24e51SJuha Jäykkä     }
3590da24e51SJuha Jäykkä     dim = 0;
3600da24e51SJuha Jäykkä     if (timestep >= 0) {
3610da24e51SJuha Jäykkä       ++dim;
3620da24e51SJuha Jäykkä     }
363cb5c4f94SJuha Jäykkä     /* prefer to split z-axis, even down to planar slices */
364c73cfb54SMatthew G. Knepley     if (dimension == 3) {
365cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->P/zslices;
366cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->N/yslices;
367cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->M/xslices;
3680da24e51SJuha Jäykkä     } else {
369cb5c4f94SJuha Jäykkä       /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */
370cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->N/yslices;
371cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->M/xslices;
3720da24e51SJuha Jäykkä     }
3730da24e51SJuha 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);
3740da24e51SJuha Jäykkä   } else {
3750da24e51SJuha Jäykkä     if (target_size < chunk_size) {
376cb5c4f94SJuha Jäykkä       /* only change the defaults if target_size < chunk_size */
3770da24e51SJuha Jäykkä       dim = 0;
3780da24e51SJuha Jäykkä       if (timestep >= 0) {
3790da24e51SJuha Jäykkä 	++dim;
3800da24e51SJuha Jäykkä       }
381cb5c4f94SJuha Jäykkä       /* prefer to split z-axis, even down to planar slices */
382c73cfb54SMatthew G. Knepley       if (dimension == 3) {
383cb5c4f94SJuha Jäykkä 	/* try splitting the z-axis to core-size bits, i.e. divide chunk size by # comm_size in z-direction */
3840da24e51SJuha Jäykkä 	if (target_size >= chunk_size/da->p) {
385cb5c4f94SJuha Jäykkä 	  /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */
3860da24e51SJuha Jäykkä 	  chunkDims[dim] = (hsize_t) ceil(da->P*1.0/da->p);
3870da24e51SJuha Jäykkä 	} else {
388cb5c4f94SJuha Jäykkä 	  /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be
389cb5c4f94SJuha Jäykkä            radical and let everyone write all they've got */
3900da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->P*1.0/da->p);
3910da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n);
3920da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m);
3930da24e51SJuha Jäykkä 	}
3940da24e51SJuha Jäykkä       } else {
395cb5c4f94SJuha Jäykkä 	/* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */
3960da24e51SJuha Jäykkä 	if (target_size >= chunk_size/da->n) {
397cb5c4f94SJuha Jäykkä 	  /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */
3980da24e51SJuha Jäykkä 	  chunkDims[dim] = (hsize_t) ceil(da->N*1.0/da->n);
3990da24e51SJuha Jäykkä 	} else {
400cb5c4f94SJuha Jäykkä 	  /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be
401cb5c4f94SJuha Jäykkä 	   radical and let everyone write all they've got */
4020da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n);
4030da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m);
4040da24e51SJuha Jäykkä 	}
4050da24e51SJuha Jäykkä 
4060da24e51SJuha Jäykkä       }
4070da24e51SJuha 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);
4080da24e51SJuha Jäykkä     } else {
409cb5c4f94SJuha Jäykkä       /* precomputed chunks are fine, we don't need to do anything */
4100da24e51SJuha Jäykkä     }
4110da24e51SJuha Jäykkä   }
4120da24e51SJuha Jäykkä   PetscFunctionReturn(0);
4130da24e51SJuha Jäykkä }
4140da24e51SJuha Jäykkä #endif
41547c6ae99SBarry Smith 
41647c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
41747c6ae99SBarry Smith #undef __FUNCT__
41847c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_HDF5_DA"
41947c6ae99SBarry Smith PetscErrorCode VecView_MPI_HDF5_DA(Vec xin,PetscViewer viewer)
42047c6ae99SBarry Smith {
4219b2a5a86SJed Brown   DM                dm;
4229b2a5a86SJed Brown   DM_DA             *da;
42347c6ae99SBarry Smith   hid_t             filespace;  /* file dataspace identifier */
4248e2ae6d7SMichael Kraus   hid_t             chunkspace; /* chunk dataset property identifier */
42547c6ae99SBarry Smith   hid_t             plist_id;   /* property list identifier */
42647c6ae99SBarry Smith   hid_t             dset_id;    /* dataset identifier */
42747c6ae99SBarry Smith   hid_t             memspace;   /* memory dataspace identifier */
42847c6ae99SBarry Smith   hid_t             file_id;
42915214e8eSMatthew G Knepley   hid_t             group;
4309a0502c6SHåkon Strandenes   hid_t             memscalartype; /* scalar type for mem (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
4319a0502c6SHåkon Strandenes   hid_t             filescalartype; /* scalar type for file (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
432d9a4edebSJed Brown   hsize_t           dim;
4330da24e51SJuha 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  */
434c73cfb54SMatthew G. Knepley   PetscInt          timestep, dimension;
4355edff71fSBarry Smith   const PetscScalar *x;
4368e2ae6d7SMichael Kraus   const char        *vecname;
43715214e8eSMatthew G Knepley   PetscErrorCode    ierr;
43882ea9b62SBarry Smith   PetscBool         dim2;
4399a0502c6SHåkon Strandenes   PetscBool         spoutput;
44047c6ae99SBarry Smith 
44147c6ae99SBarry Smith   PetscFunctionBegin;
44215214e8eSMatthew G Knepley   ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr);
4438e2ae6d7SMichael Kraus   ierr = PetscViewerHDF5GetTimestep(viewer, &timestep);CHKERRQ(ierr);
44482ea9b62SBarry Smith   ierr = PetscViewerHDF5GetBaseDimension2(viewer,&dim2);CHKERRQ(ierr);
4459a0502c6SHåkon Strandenes   ierr = PetscViewerHDF5GetSPOutput(viewer,&spoutput);CHKERRQ(ierr);
44615214e8eSMatthew G Knepley 
447c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&dm);CHKERRQ(ierr);
448ce94432eSBarry Smith   if (!dm) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
4499b2a5a86SJed Brown   da = (DM_DA*)dm->data;
450c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dimension);CHKERRQ(ierr);
45147c6ae99SBarry Smith 
4528e2ae6d7SMichael Kraus   /* Create the dataspace for the dataset.
4538e2ae6d7SMichael Kraus    *
4548e2ae6d7SMichael Kraus    * dims - holds the current dimensions of the dataset
4558e2ae6d7SMichael Kraus    *
4568e2ae6d7SMichael Kraus    * maxDims - holds the maximum dimensions of the dataset (unlimited
4578e2ae6d7SMichael Kraus    * for the number of time steps with the current dimensions for the
4588e2ae6d7SMichael Kraus    * other dimensions; so only additional time steps can be added).
4598e2ae6d7SMichael Kraus    *
4608e2ae6d7SMichael Kraus    * chunkDims - holds the size of a single time step (required to
4618e2ae6d7SMichael Kraus    * permit extending dataset).
4628e2ae6d7SMichael Kraus    */
4638e2ae6d7SMichael Kraus   dim = 0;
4648e2ae6d7SMichael Kraus   if (timestep >= 0) {
4658e2ae6d7SMichael Kraus     dims[dim]      = timestep+1;
4668e2ae6d7SMichael Kraus     maxDims[dim]   = H5S_UNLIMITED;
4678e2ae6d7SMichael Kraus     chunkDims[dim] = 1;
4688e2ae6d7SMichael Kraus     ++dim;
4698e2ae6d7SMichael Kraus   }
470c73cfb54SMatthew G. Knepley   if (dimension == 3) {
471acba2ac6SBarry Smith     ierr           = PetscHDF5IntCast(da->P,dims+dim);CHKERRQ(ierr);
4728e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4738e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4748e2ae6d7SMichael Kraus     ++dim;
4758e2ae6d7SMichael Kraus   }
476c73cfb54SMatthew G. Knepley   if (dimension > 1) {
477acba2ac6SBarry Smith     ierr           = PetscHDF5IntCast(da->N,dims+dim);CHKERRQ(ierr);
4788e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4798e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4808e2ae6d7SMichael Kraus     ++dim;
4818e2ae6d7SMichael Kraus   }
482acba2ac6SBarry Smith   ierr           = PetscHDF5IntCast(da->M,dims+dim);CHKERRQ(ierr);
4838e2ae6d7SMichael Kraus   maxDims[dim]   = dims[dim];
4848e2ae6d7SMichael Kraus   chunkDims[dim] = dims[dim];
4858e2ae6d7SMichael Kraus   ++dim;
48682ea9b62SBarry Smith   if (da->w > 1 || dim2) {
487acba2ac6SBarry Smith     ierr           = PetscHDF5IntCast(da->w,dims+dim);CHKERRQ(ierr);
4888e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4898e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4908e2ae6d7SMichael Kraus     ++dim;
4918e2ae6d7SMichael Kraus   }
49247c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
4938e2ae6d7SMichael Kraus   dims[dim]      = 2;
4948e2ae6d7SMichael Kraus   maxDims[dim]   = dims[dim];
4958e2ae6d7SMichael Kraus   chunkDims[dim] = dims[dim];
4968e2ae6d7SMichael Kraus   ++dim;
49747c6ae99SBarry Smith #endif
4980da24e51SJuha Jäykkä 
499c73cfb54SMatthew G. Knepley   ierr = VecGetHDF5ChunkSize(da, xin, dimension, timestep, chunkDims); CHKERRQ(ierr);
5000da24e51SJuha Jäykkä 
501729ab6d9SBarry Smith   PetscStackCallHDF5Return(filespace,H5Screate_simple,(dim, dims, maxDims));
50247c6ae99SBarry Smith 
50315214e8eSMatthew G Knepley #if defined(PETSC_USE_REAL_SINGLE)
5049a0502c6SHåkon Strandenes   memscalartype = H5T_NATIVE_FLOAT;
5059a0502c6SHåkon Strandenes   filescalartype = H5T_NATIVE_FLOAT;
50615214e8eSMatthew G Knepley #elif defined(PETSC_USE_REAL___FLOAT128)
50715214e8eSMatthew G Knepley #error "HDF5 output with 128 bit floats not supported."
50815214e8eSMatthew G Knepley #else
5099a0502c6SHåkon Strandenes   memscalartype = H5T_NATIVE_DOUBLE;
5109a0502c6SHåkon Strandenes   if (spoutput == PETSC_TRUE) filescalartype = H5T_NATIVE_FLOAT;
5119a0502c6SHåkon Strandenes   else filescalartype = H5T_NATIVE_DOUBLE;
51215214e8eSMatthew G Knepley #endif
51315214e8eSMatthew G Knepley 
51447c6ae99SBarry Smith   /* Create the dataset with default properties and close filespace */
51547c6ae99SBarry Smith   ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr);
51615214e8eSMatthew G Knepley   if (!H5Lexists(group, vecname, H5P_DEFAULT)) {
5178e2ae6d7SMichael Kraus     /* Create chunk */
518729ab6d9SBarry Smith     PetscStackCallHDF5Return(chunkspace,H5Pcreate,(H5P_DATASET_CREATE));
519729ab6d9SBarry Smith     PetscStackCallHDF5(H5Pset_chunk,(chunkspace, dim, chunkDims));
5208e2ae6d7SMichael Kraus 
52147c6ae99SBarry Smith #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
5229a0502c6SHåkon Strandenes     PetscStackCallHDF5Return(dset_id,H5Dcreate2,(group, vecname, filescalartype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT));
52347c6ae99SBarry Smith #else
5249a0502c6SHåkon Strandenes     PetscStackCallHDF5Return(dset_id,H5Dcreate,(group, vecname, filescalartype, filespace, H5P_DEFAULT));
52547c6ae99SBarry Smith #endif
52615214e8eSMatthew G Knepley   } else {
527729ab6d9SBarry Smith     PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT));
528729ab6d9SBarry Smith     PetscStackCallHDF5(H5Dset_extent,(dset_id, dims));
52915214e8eSMatthew G Knepley   }
530729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
53147c6ae99SBarry Smith 
53247c6ae99SBarry Smith   /* Each process defines a dataset and writes it to the hyperslab in the file */
5338e2ae6d7SMichael Kraus   dim = 0;
5348e2ae6d7SMichael Kraus   if (timestep >= 0) {
5358e2ae6d7SMichael Kraus     offset[dim] = timestep;
5368e2ae6d7SMichael Kraus     ++dim;
5378e2ae6d7SMichael Kraus   }
538c73cfb54SMatthew G. Knepley   if (dimension == 3) {ierr = PetscHDF5IntCast(da->zs,offset + dim++);CHKERRQ(ierr);}
539c73cfb54SMatthew G. Knepley   if (dimension > 1)  {ierr = PetscHDF5IntCast(da->ys,offset + dim++);CHKERRQ(ierr);}
540acba2ac6SBarry Smith   ierr = PetscHDF5IntCast(da->xs/da->w,offset + dim++);CHKERRQ(ierr);
54182ea9b62SBarry Smith   if (da->w > 1 || dim2) offset[dim++] = 0;
54247c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5438e2ae6d7SMichael Kraus   offset[dim++] = 0;
54447c6ae99SBarry Smith #endif
5458e2ae6d7SMichael Kraus   dim = 0;
5468e2ae6d7SMichael Kraus   if (timestep >= 0) {
5478e2ae6d7SMichael Kraus     count[dim] = 1;
5488e2ae6d7SMichael Kraus     ++dim;
5498e2ae6d7SMichael Kraus   }
550c73cfb54SMatthew G. Knepley   if (dimension == 3) {ierr = PetscHDF5IntCast(da->ze - da->zs,count + dim++);CHKERRQ(ierr);}
551c73cfb54SMatthew G. Knepley   if (dimension > 1)  {ierr = PetscHDF5IntCast(da->ye - da->ys,count + dim++);CHKERRQ(ierr);}
552acba2ac6SBarry Smith   ierr = PetscHDF5IntCast((da->xe - da->xs)/da->w,count + dim++);CHKERRQ(ierr);
55382ea9b62SBarry Smith   if (da->w > 1 || dim2) {ierr = PetscHDF5IntCast(da->w,count + dim++);CHKERRQ(ierr);}
55447c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5558e2ae6d7SMichael Kraus   count[dim++] = 2;
55647c6ae99SBarry Smith #endif
557729ab6d9SBarry Smith   PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL));
558729ab6d9SBarry Smith   PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id));
559729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL));
56047c6ae99SBarry Smith 
56147c6ae99SBarry Smith   /* Create property list for collective dataset write */
562729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_DATASET_XFER));
56347c6ae99SBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
564729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_dxpl_mpio,(plist_id, H5FD_MPIO_COLLECTIVE));
56547c6ae99SBarry Smith #endif
56647c6ae99SBarry Smith   /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
56747c6ae99SBarry Smith 
5685edff71fSBarry Smith   ierr   = VecGetArrayRead(xin, &x);CHKERRQ(ierr);
5699a0502c6SHåkon Strandenes   PetscStackCallHDF5(H5Dwrite,(dset_id, memscalartype, memspace, filespace, plist_id, x));
570729ab6d9SBarry Smith   PetscStackCallHDF5(H5Fflush,(file_id, H5F_SCOPE_GLOBAL));
5715edff71fSBarry Smith   ierr   = VecRestoreArrayRead(xin, &x);CHKERRQ(ierr);
57247c6ae99SBarry Smith 
57347c6ae99SBarry Smith   /* Close/release resources */
57415214e8eSMatthew G Knepley   if (group != file_id) {
575729ab6d9SBarry Smith     PetscStackCallHDF5(H5Gclose,(group));
57615214e8eSMatthew G Knepley   }
577729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
578729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
579729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(memspace));
580729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dset_id));
58147c6ae99SBarry Smith   ierr   = PetscInfo1(xin,"Wrote Vec object with name %s\n",vecname);CHKERRQ(ierr);
58247c6ae99SBarry Smith   PetscFunctionReturn(0);
58347c6ae99SBarry Smith }
58447c6ae99SBarry Smith #endif
58547c6ae99SBarry Smith 
58609573ac7SBarry Smith extern PetscErrorCode VecView_MPI_Draw_DA1d(Vec,PetscViewer);
58747c6ae99SBarry Smith 
58847c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
58947c6ae99SBarry Smith #undef __FUNCT__
590aa219208SBarry Smith #define __FUNCT__ "DMDAArrayMPIIO"
591aa219208SBarry Smith static PetscErrorCode DMDAArrayMPIIO(DM da,PetscViewer viewer,Vec xin,PetscBool write)
59247c6ae99SBarry Smith {
59347c6ae99SBarry Smith   PetscErrorCode    ierr;
59447c6ae99SBarry Smith   MPI_File          mfdes;
59547c6ae99SBarry Smith   PetscMPIInt       gsizes[4],lsizes[4],lstarts[4],asiz,dof;
59647c6ae99SBarry Smith   MPI_Datatype      view;
59747c6ae99SBarry Smith   const PetscScalar *array;
59847c6ae99SBarry Smith   MPI_Offset        off;
59947c6ae99SBarry Smith   MPI_Aint          ub,ul;
60047c6ae99SBarry Smith   PetscInt          type,rows,vecrows,tr[2];
60147c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
6020c169764Sdmay   PetscBool         skipheader;
60347c6ae99SBarry Smith 
60447c6ae99SBarry Smith   PetscFunctionBegin;
60547c6ae99SBarry Smith   ierr = VecGetSize(xin,&vecrows);CHKERRQ(ierr);
6060c169764Sdmay   ierr = PetscViewerBinaryGetSkipHeader(viewer,&skipheader);CHKERRQ(ierr);
60747c6ae99SBarry Smith   if (!write) {
60847c6ae99SBarry Smith     /* Read vector header. */
6090c169764Sdmay     if (!skipheader) {
610060da220SMatthew G. Knepley       ierr = PetscViewerBinaryRead(viewer,tr,2,NULL,PETSC_INT);CHKERRQ(ierr);
61147c6ae99SBarry Smith       type = tr[0];
61247c6ae99SBarry Smith       rows = tr[1];
613ce94432eSBarry Smith       if (type != VEC_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Not vector next in file");
614ce94432eSBarry Smith       if (rows != vecrows) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Vector in file not same size as DMDA vector");
6150c169764Sdmay     }
61647c6ae99SBarry Smith   } else {
61747c6ae99SBarry Smith     tr[0] = VEC_FILE_CLASSID;
61847c6ae99SBarry Smith     tr[1] = vecrows;
6190c169764Sdmay     if (!skipheader) {
62047c6ae99SBarry Smith       ierr  = PetscViewerBinaryWrite(viewer,tr,2,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
62147c6ae99SBarry Smith     }
6220c169764Sdmay   }
62347c6ae99SBarry Smith 
6244dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->w,&dof);CHKERRQ(ierr);
6254dc2109aSBarry Smith   gsizes[0]  = dof;
6264dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->M,gsizes+1);CHKERRQ(ierr);
6274dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->N,gsizes+2);CHKERRQ(ierr);
628334634e2SJed Brown   ierr       = PetscMPIIntCast(dd->P,gsizes+3);CHKERRQ(ierr);
6294dc2109aSBarry Smith   lsizes[0]  = dof;
6304dc2109aSBarry Smith   ierr       = PetscMPIIntCast((dd->xe-dd->xs)/dof,lsizes+1);CHKERRQ(ierr);
6314dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->ye-dd->ys,lsizes+2);CHKERRQ(ierr);
6324dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->ze-dd->zs,lsizes+3);CHKERRQ(ierr);
6334dc2109aSBarry Smith   lstarts[0] = 0;
6344dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->xs/dof,lstarts+1);CHKERRQ(ierr);
6354dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->ys,lstarts+2);CHKERRQ(ierr);
6364dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->zs,lstarts+3);CHKERRQ(ierr);
637c73cfb54SMatthew G. Knepley   ierr       = MPI_Type_create_subarray(da->dim+1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view);CHKERRQ(ierr);
63847c6ae99SBarry Smith   ierr       = MPI_Type_commit(&view);CHKERRQ(ierr);
63947c6ae99SBarry Smith 
64047c6ae99SBarry Smith   ierr = PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);CHKERRQ(ierr);
64147c6ae99SBarry Smith   ierr = PetscViewerBinaryGetMPIIOOffset(viewer,&off);CHKERRQ(ierr);
64247c6ae99SBarry Smith   ierr = MPI_File_set_view(mfdes,off,MPIU_SCALAR,view,(char*)"native",MPI_INFO_NULL);CHKERRQ(ierr);
64347c6ae99SBarry Smith   ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr);
64447c6ae99SBarry Smith   asiz = lsizes[1]*(lsizes[2] > 0 ? lsizes[2] : 1)*(lsizes[3] > 0 ? lsizes[3] : 1)*dof;
64547c6ae99SBarry Smith   if (write) {
64647c6ae99SBarry Smith     ierr = MPIU_File_write_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr);
64747c6ae99SBarry Smith   } else {
64847c6ae99SBarry Smith     ierr = MPIU_File_read_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr);
64947c6ae99SBarry Smith   }
65047c6ae99SBarry Smith   ierr = MPI_Type_get_extent(view,&ul,&ub);CHKERRQ(ierr);
65147c6ae99SBarry Smith   ierr = PetscViewerBinaryAddMPIIOOffset(viewer,ub);CHKERRQ(ierr);
65247c6ae99SBarry Smith   ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr);
65347c6ae99SBarry Smith   ierr = MPI_Type_free(&view);CHKERRQ(ierr);
65447c6ae99SBarry Smith   PetscFunctionReturn(0);
65547c6ae99SBarry Smith }
65647c6ae99SBarry Smith #endif
65747c6ae99SBarry Smith 
65847c6ae99SBarry Smith #undef __FUNCT__
65947c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_DA"
6607087cfbeSBarry Smith PetscErrorCode  VecView_MPI_DA(Vec xin,PetscViewer viewer)
66147c6ae99SBarry Smith {
6629a42bb27SBarry Smith   DM                da;
66347c6ae99SBarry Smith   PetscErrorCode    ierr;
66447c6ae99SBarry Smith   PetscInt          dim;
66547c6ae99SBarry Smith   Vec               natural;
6664061b8bfSJed Brown   PetscBool         isdraw,isvtk;
66747c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
66847c6ae99SBarry Smith   PetscBool         ishdf5;
66947c6ae99SBarry Smith #endif
6703f3fd955SJed Brown   const char        *prefix,*name;
671a261c58fSBarry Smith   PetscViewerFormat format;
67247c6ae99SBarry Smith 
67347c6ae99SBarry Smith   PetscFunctionBegin;
674c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
675ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
676251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
677251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr);
67847c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
679251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr);
68047c6ae99SBarry Smith #endif
68147c6ae99SBarry Smith   if (isdraw) {
6821321219cSEthan Coon     ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
68347c6ae99SBarry Smith     if (dim == 1) {
68447c6ae99SBarry Smith       ierr = VecView_MPI_Draw_DA1d(xin,viewer);CHKERRQ(ierr);
68547c6ae99SBarry Smith     } else if (dim == 2) {
68647c6ae99SBarry Smith       ierr = VecView_MPI_Draw_DA2d(xin,viewer);CHKERRQ(ierr);
687ce94432eSBarry Smith     } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot graphically view vector associated with this dimensional DMDA %D",dim);
6884061b8bfSJed Brown   } else if (isvtk) {           /* Duplicate the Vec and hold a reference to the DM */
6894061b8bfSJed Brown     Vec Y;
6904061b8bfSJed Brown     ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr);
6914061b8bfSJed Brown     ierr = VecDuplicate(xin,&Y);CHKERRQ(ierr);
692b51b94faSJed Brown     if (((PetscObject)xin)->name) {
693b51b94faSJed Brown       /* If xin was named, copy the name over to Y. The duplicate names are safe because nobody else will ever see Y. */
694b51b94faSJed Brown       ierr = PetscObjectSetName((PetscObject)Y,((PetscObject)xin)->name);CHKERRQ(ierr);
695b51b94faSJed Brown     }
6964061b8bfSJed Brown     ierr = VecCopy(xin,Y);CHKERRQ(ierr);
69762b69a3fSMatthew G Knepley     ierr = PetscViewerVTKAddField(viewer,(PetscObject)da,DMDAVTKWriteAll,PETSC_VTK_POINT_FIELD,(PetscObject)Y);CHKERRQ(ierr);
69847c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
69947c6ae99SBarry Smith   } else if (ishdf5) {
70047c6ae99SBarry Smith     ierr = VecView_MPI_HDF5_DA(xin,viewer);CHKERRQ(ierr);
70147c6ae99SBarry Smith #endif
70247c6ae99SBarry Smith   } else {
70347c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
70447c6ae99SBarry Smith     PetscBool isbinary,isMPIIO;
70547c6ae99SBarry Smith 
706251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
70747c6ae99SBarry Smith     if (isbinary) {
708bc196f7cSDave May       ierr = PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO);CHKERRQ(ierr);
70947c6ae99SBarry Smith       if (isMPIIO) {
710aa219208SBarry Smith         ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_TRUE);CHKERRQ(ierr);
71147c6ae99SBarry Smith         PetscFunctionReturn(0);
71247c6ae99SBarry Smith       }
71347c6ae99SBarry Smith     }
71447c6ae99SBarry Smith #endif
71547c6ae99SBarry Smith 
71647c6ae99SBarry Smith     /* call viewer on natural ordering */
71747c6ae99SBarry Smith     ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr);
718aa219208SBarry Smith     ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr);
71947c6ae99SBarry Smith     ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr);
720aa219208SBarry Smith     ierr = DMDAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr);
721aa219208SBarry Smith     ierr = DMDAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr);
7223f3fd955SJed Brown     ierr = PetscObjectGetName((PetscObject)xin,&name);CHKERRQ(ierr);
7233f3fd955SJed Brown     ierr = PetscObjectSetName((PetscObject)natural,name);CHKERRQ(ierr);
724a261c58fSBarry Smith 
725a261c58fSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
726a261c58fSBarry Smith     if (format == PETSC_VIEWER_BINARY_MATLAB) {
727a261c58fSBarry Smith       /* temporarily remove viewer format so it won't trigger in the VecView() */
7286a9046bcSBarry Smith       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
729a261c58fSBarry Smith     }
730a261c58fSBarry Smith 
73147c6ae99SBarry Smith     ierr = VecView(natural,viewer);CHKERRQ(ierr);
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 */
7416a9046bcSBarry Smith       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
742a261c58fSBarry Smith       ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
743a261c58fSBarry Smith       ierr = PetscViewerBinaryGetInfoPointer(viewer,&info);CHKERRQ(ierr);
744a261c58fSBarry Smith       ierr = DMDAGetInfo(da,&dim,&ni,&nj,&nk,&pi,&pj,&pk,&dof,0,0,0,0,0);CHKERRQ(ierr);
745da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n");CHKERRQ(ierr);
746da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#$$ tmp = PetscBinaryRead(fd); \n");CHKERRQ(ierr);
747da88d4d4SJed Brown       if (dim == 1) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d);\n",dof,ni);CHKERRQ(ierr); }
748da88d4d4SJed Brown       if (dim == 2) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d);\n",dof,ni,nj);CHKERRQ(ierr); }
749da88d4d4SJed Brown       if (dim == 3) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d,%d);\n",dof,ni,nj,nk);CHKERRQ(ierr); }
750a261c58fSBarry Smith 
751a261c58fSBarry Smith       for (n=0; n<dof; n++) {
752a261c58fSBarry Smith         ierr = DMDAGetFieldName(da,n,&fieldname);CHKERRQ(ierr);
753da88d4d4SJed Brown         if (!fieldname || !fieldname[0]) {
754da88d4d4SJed Brown           ierr = PetscSNPrintf(fieldbuf,sizeof fieldbuf,"field%D",n);CHKERRQ(ierr);
755da88d4d4SJed Brown           fieldname = fieldbuf;
756a261c58fSBarry Smith         }
757da88d4d4SJed Brown         if (dim == 1) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:))';\n",name,fieldname,n+1);CHKERRQ(ierr); }
758da88d4d4SJed Brown         if (dim == 2) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:,:))';\n",name,fieldname,n+1);CHKERRQ(ierr); }
759da88d4d4SJed Brown         if (dim == 3) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = permute(squeeze(tmp(%d,:,:,:)),[2 1 3]);\n",name,fieldname,n+1);CHKERRQ(ierr);}
760a261c58fSBarry Smith       }
761da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#$$ clear tmp; \n");CHKERRQ(ierr);
762da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n");CHKERRQ(ierr);
763a261c58fSBarry Smith     }
764a261c58fSBarry Smith 
765fcfd50ebSBarry Smith     ierr = VecDestroy(&natural);CHKERRQ(ierr);
76647c6ae99SBarry Smith   }
76747c6ae99SBarry Smith   PetscFunctionReturn(0);
76847c6ae99SBarry Smith }
76947c6ae99SBarry Smith 
77047c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
77147c6ae99SBarry Smith #undef __FUNCT__
77247c6ae99SBarry Smith #define __FUNCT__ "VecLoad_HDF5_DA"
77347c6ae99SBarry Smith PetscErrorCode VecLoad_HDF5_DA(Vec xin, PetscViewer viewer)
77447c6ae99SBarry Smith {
7759a42bb27SBarry Smith   DM             da;
77647c6ae99SBarry Smith   PetscErrorCode ierr;
777ec7ae49cSHåkon Strandenes   hsize_t        dim,rdim;
778ec7ae49cSHåkon Strandenes   hsize_t        dims[6]={0},count[6]={0},offset[6]={0};
779ec7ae49cSHåkon Strandenes   PetscInt       dimension,timestep,dofInd;
78047c6ae99SBarry Smith   PetscScalar    *x;
78147c6ae99SBarry Smith   const char     *vecname;
78247c6ae99SBarry Smith   hid_t          filespace; /* file dataspace identifier */
78347c6ae99SBarry Smith   hid_t          plist_id;  /* property list identifier */
78447c6ae99SBarry Smith   hid_t          dset_id;   /* dataset identifier */
78547c6ae99SBarry Smith   hid_t          memspace;  /* memory dataspace identifier */
786bfd264e7SBarry Smith   hid_t          file_id,group;
7877bcbaff4SBarry Smith   hid_t          scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
7889c7c4993SBarry Smith   DM_DA          *dd;
789ec7ae49cSHåkon Strandenes   PetscBool      dim2 = PETSC_FALSE;
79047c6ae99SBarry Smith 
79147c6ae99SBarry Smith   PetscFunctionBegin;
7927bcbaff4SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
7937bcbaff4SBarry Smith   scalartype = H5T_NATIVE_FLOAT;
7947bcbaff4SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
7957bcbaff4SBarry Smith #error "HDF5 output with 128 bit floats not supported."
7967bcbaff4SBarry Smith #else
7977bcbaff4SBarry Smith   scalartype = H5T_NATIVE_DOUBLE;
7987bcbaff4SBarry Smith #endif
7997bcbaff4SBarry Smith 
800bfd264e7SBarry Smith   ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr);
801ec7ae49cSHåkon Strandenes   ierr = PetscViewerHDF5GetTimestep(viewer, &timestep);CHKERRQ(ierr);
802ec7ae49cSHåkon Strandenes   ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr);
803c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
8049c7c4993SBarry Smith   dd   = (DM_DA*)da->data;
805c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(da, &dimension);CHKERRQ(ierr);
80647c6ae99SBarry Smith 
807ec7ae49cSHåkon Strandenes   /* Open dataset */
80847c6ae99SBarry Smith #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
809729ab6d9SBarry Smith   PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT));
81047c6ae99SBarry Smith #else
811729ab6d9SBarry Smith   PetscStackCallHDF5Return(dset_id,H5Dopen,(group, vecname));
81247c6ae99SBarry Smith #endif
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 */
839ec7ae49cSHåkon Strandenes     else if (dd->w != (PetscInt) dims[dofInd]) SETERRQ2(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 */
8426c4ed002SBarry Smith   } else if (rdim != dim) SETERRQ3(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) {
852ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->zs,offset + dim);CHKERRQ(ierr);
853ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->ze - dd->zs,count + dim);CHKERRQ(ierr);
854ec7ae49cSHåkon Strandenes     ++dim;
855ec7ae49cSHåkon Strandenes   }
856ec7ae49cSHåkon Strandenes   if (dimension > 1) {
857ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->ys,offset + dim);CHKERRQ(ierr);
858ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->ye - dd->ys,count + dim);CHKERRQ(ierr);
859ec7ae49cSHåkon Strandenes     ++dim;
860ec7ae49cSHåkon Strandenes   }
861ec7ae49cSHåkon Strandenes   ierr = PetscHDF5IntCast(dd->xs/dd->w,offset + dim);CHKERRQ(ierr);
862ec7ae49cSHåkon Strandenes   ierr = PetscHDF5IntCast((dd->xe - dd->xs)/dd->w,count + dim);CHKERRQ(ierr);
863ec7ae49cSHåkon Strandenes   ++dim;
864ec7ae49cSHåkon Strandenes   if (dd->w > 1 || dim2) {
865ec7ae49cSHåkon Strandenes     offset[dim] = 0;
866ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->w,count + dim);CHKERRQ(ierr);
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 
87947c6ae99SBarry Smith   /* Create property list for collective dataset write */
880729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_DATASET_XFER));
88147c6ae99SBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
882729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_dxpl_mpio,(plist_id, H5FD_MPIO_COLLECTIVE));
88347c6ae99SBarry Smith #endif
884ec7ae49cSHåkon Strandenes   /* To read dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
88547c6ae99SBarry Smith 
88647c6ae99SBarry Smith   ierr   = VecGetArray(xin, &x);CHKERRQ(ierr);
887a1dbdba5SBarry Smith   PetscStackCallHDF5(H5Dread,(dset_id, scalartype, memspace, filespace, plist_id, x));
88847c6ae99SBarry Smith   ierr   = VecRestoreArray(xin, &x);CHKERRQ(ierr);
88947c6ae99SBarry Smith 
89047c6ae99SBarry Smith   /* Close/release resources */
891bfd264e7SBarry Smith   if (group != file_id) {
892729ab6d9SBarry Smith     PetscStackCallHDF5(H5Gclose,(group));
893bfd264e7SBarry Smith   }
894729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
895729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
896729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(memspace));
897729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dset_id));
89847c6ae99SBarry Smith   PetscFunctionReturn(0);
89947c6ae99SBarry Smith }
90047c6ae99SBarry Smith #endif
90147c6ae99SBarry Smith 
90247c6ae99SBarry Smith #undef __FUNCT__
90347c6ae99SBarry Smith #define __FUNCT__ "VecLoad_Binary_DA"
90447c6ae99SBarry Smith PetscErrorCode VecLoad_Binary_DA(Vec xin, PetscViewer viewer)
90547c6ae99SBarry Smith {
9069a42bb27SBarry Smith   DM             da;
90747c6ae99SBarry Smith   PetscErrorCode ierr;
90847c6ae99SBarry Smith   Vec            natural;
90947c6ae99SBarry Smith   const char     *prefix;
91047c6ae99SBarry Smith   PetscInt       bs;
91147c6ae99SBarry Smith   PetscBool      flag;
91247c6ae99SBarry Smith   DM_DA          *dd;
91347c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
91447c6ae99SBarry Smith   PetscBool isMPIIO;
91547c6ae99SBarry Smith #endif
91647c6ae99SBarry Smith 
91747c6ae99SBarry Smith   PetscFunctionBegin;
918c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
91947c6ae99SBarry Smith   dd   = (DM_DA*)da->data;
92047c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
921bc196f7cSDave May   ierr = PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO);CHKERRQ(ierr);
92247c6ae99SBarry Smith   if (isMPIIO) {
923aa219208SBarry Smith     ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_FALSE);CHKERRQ(ierr);
92447c6ae99SBarry Smith     PetscFunctionReturn(0);
92547c6ae99SBarry Smith   }
92647c6ae99SBarry Smith #endif
92747c6ae99SBarry Smith 
92847c6ae99SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr);
929aa219208SBarry Smith   ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr);
93047c6ae99SBarry Smith   ierr = PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name);CHKERRQ(ierr);
93147c6ae99SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr);
9328aea9937SBarry Smith   ierr = VecLoad(natural,viewer);CHKERRQ(ierr);
933aa219208SBarry Smith   ierr = DMDANaturalToGlobalBegin(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr);
934aa219208SBarry Smith   ierr = DMDANaturalToGlobalEnd(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr);
935fcfd50ebSBarry Smith   ierr = VecDestroy(&natural);CHKERRQ(ierr);
936aa219208SBarry Smith   ierr = PetscInfo(xin,"Loading vector from natural ordering into DMDA\n");CHKERRQ(ierr);
937c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(NULL,((PetscObject)xin)->prefix,"-vecload_block_size",&bs,&flag);CHKERRQ(ierr);
93847c6ae99SBarry Smith   if (flag && bs != dd->w) {
939aa219208SBarry Smith     ierr = PetscInfo2(xin,"Block size in file %D not equal to DMDA's dof %D\n",bs,dd->w);CHKERRQ(ierr);
94047c6ae99SBarry Smith   }
94147c6ae99SBarry Smith   PetscFunctionReturn(0);
94247c6ae99SBarry Smith }
94347c6ae99SBarry Smith 
94447c6ae99SBarry Smith #undef __FUNCT__
94547c6ae99SBarry Smith #define __FUNCT__ "VecLoad_Default_DA"
9467087cfbeSBarry Smith PetscErrorCode  VecLoad_Default_DA(Vec xin, PetscViewer viewer)
94747c6ae99SBarry Smith {
94847c6ae99SBarry Smith   PetscErrorCode ierr;
9499a42bb27SBarry Smith   DM             da;
95047c6ae99SBarry Smith   PetscBool      isbinary;
95147c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
95247c6ae99SBarry Smith   PetscBool ishdf5;
95347c6ae99SBarry Smith #endif
95447c6ae99SBarry Smith 
95547c6ae99SBarry Smith   PetscFunctionBegin;
956c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
957ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
95847c6ae99SBarry Smith 
95947c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
960251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr);
96147c6ae99SBarry Smith #endif
962251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
96347c6ae99SBarry Smith 
96447c6ae99SBarry Smith   if (isbinary) {
96547c6ae99SBarry Smith     ierr = VecLoad_Binary_DA(xin,viewer);CHKERRQ(ierr);
96647c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
96747c6ae99SBarry Smith   } else if (ishdf5) {
96847c6ae99SBarry Smith     ierr = VecLoad_HDF5_DA(xin,viewer);CHKERRQ(ierr);
96947c6ae99SBarry Smith #endif
970d34fcf5fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for vector loading", ((PetscObject)viewer)->type_name);
97147c6ae99SBarry Smith   PetscFunctionReturn(0);
97247c6ae99SBarry Smith }
973