xref: /petsc/src/dm/impls/da/gr2.c (revision e5ab1681210b093c2ee26b200f61e98eac43a430)
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 {
15*e5ab1681SLisandro Dalcin   PetscMPIInt       rank;
16*e5ab1681SLisandro Dalcin   PetscInt          m,n,dof,k;
17*e5ab1681SLisandro Dalcin   PetscReal         xmin,xmax,ymin,ymax,min,max;
185edff71fSBarry Smith   const PetscScalar *xy,*v;
1947c6ae99SBarry Smith   PetscBool         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;
34*e5ab1681SLisandro Dalcin   PetscInt          m,n,i,j,k,dof,id,c1,c2,c3,c4;
35*e5ab1681SLisandro Dalcin   PetscReal         min,max,x1,x2,x3,x4,y_1,y2,y3,y4;
36*e5ab1681SLisandro Dalcin   const PetscScalar *xy,*v;
3747c6ae99SBarry Smith 
3847c6ae99SBarry Smith   PetscFunctionBegin;
3947c6ae99SBarry Smith   m    = zctx->m;
4047c6ae99SBarry Smith   n    = zctx->n;
41*e5ab1681SLisandro Dalcin   dof  = zctx->dof;
4247c6ae99SBarry Smith   k    = zctx->k;
4347c6ae99SBarry Smith   xy   = zctx->xy;
44*e5ab1681SLisandro Dalcin   v    = zctx->v;
4547c6ae99SBarry Smith   min  = zctx->min;
46f3f0eb19SBarry Smith   max  = zctx->max;
4747c6ae99SBarry Smith 
4847c6ae99SBarry Smith   /* PetscDraw the contour plot patch */
49*e5ab1681SLisandro 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]);
55*e5ab1681SLisandro 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]);
60*e5ab1681SLisandro 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]);
65*e5ab1681SLisandro 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]);
70*e5ab1681SLisandro 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*e5ab1681SLisandro Dalcin   if (!zctx->rank) {
83*e5ab1681SLisandro 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);
86*e5ab1681SLisandro Dalcin       x  = xl + .30*(xr - xl);
87109c9344SBarry Smith       xl = xl + .01*(xr - xl);
88*e5ab1681SLisandro Dalcin       y  = yr - .30*(yr - yl);
89109c9344SBarry Smith       yl = yl + .01*(yr - yl);
90*e5ab1681SLisandro Dalcin       if (zctx->name0) {ierr = PetscDrawString(draw,x,yl,PETSC_DRAW_BLACK,zctx->name0);CHKERRQ(ierr);}
91*e5ab1681SLisandro 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     */
97*e5ab1681SLisandro Dalcin     {
98*e5ab1681SLisandro Dalcin       PetscReal xmin = zctx->xmin, ymin = zctx->ymin;
99*e5ab1681SLisandro Dalcin       PetscReal xmax = zctx->xmax, ymax = zctx->ymax;
100*e5ab1681SLisandro Dalcin       char      value[16]; size_t len; PetscReal w;
101*e5ab1681SLisandro Dalcin       ierr = PetscSNPrintf(value,16,"%f",xmin);CHKERRQ(ierr);
102*e5ab1681SLisandro Dalcin       ierr = PetscDrawString(draw,xmin,ymin - .05*(ymax - ymin),PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
103*e5ab1681SLisandro Dalcin       ierr = PetscSNPrintf(value,16,"%f",xmax);CHKERRQ(ierr);
1040e4fe250SBarry Smith       ierr = PetscStrlen(value,&len);CHKERRQ(ierr);
1050298fd71SBarry Smith       ierr = PetscDrawStringGetSize(draw,&w,NULL);CHKERRQ(ierr);
106*e5ab1681SLisandro Dalcin       ierr = PetscDrawString(draw,xmax - len*w,ymin - .05*(ymax - ymin),PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
107*e5ab1681SLisandro Dalcin       ierr = PetscSNPrintf(value,16,"%f",ymin);CHKERRQ(ierr);
108*e5ab1681SLisandro Dalcin       ierr = PetscDrawString(draw,xmin - .05*(xmax - xmin),ymin,PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
109*e5ab1681SLisandro Dalcin       ierr = PetscSNPrintf(value,16,"%f",ymax);CHKERRQ(ierr);
110*e5ab1681SLisandro Dalcin       ierr = PetscDrawString(draw,xmin - .05*(xmax - xmin),ymax,PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
111*e5ab1681SLisandro Dalcin     }
112*e5ab1681SLisandro Dalcin   }
113*e5ab1681SLisandro 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;
125*e5ab1681SLisandro 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;
1418865f1eaSKarl Rupp 
14247c6ae99SBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1435b399a63SLisandro Dalcin   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
1445b399a63SLisandro Dalcin   if (isnull) PetscFunctionReturn(0);
1455b399a63SLisandro Dalcin 
14603193ff8SBarry Smith   ierr = PetscViewerDrawGetBounds(viewer,&nbounds,&bounds);CHKERRQ(ierr);
14747c6ae99SBarry Smith 
148c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
149ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
15047c6ae99SBarry Smith 
15147c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr);
152*e5ab1681SLisandro Dalcin   ierr = MPI_Comm_rank(comm,&zctx.rank);CHKERRQ(ierr);
15347c6ae99SBarry Smith 
1541321219cSEthan Coon   ierr = DMDAGetInfo(da,0,&M,&N,0,&zctx.m,&zctx.n,0,&w,&s,&bx,&by,0,&st);CHKERRQ(ierr);
1550298fd71SBarry Smith   ierr = DMDAGetOwnershipRanges(da,&lx,&ly,NULL);CHKERRQ(ierr);
15647c6ae99SBarry Smith 
15747c6ae99SBarry Smith   /*
15847c6ae99SBarry Smith         Obtain a sequential vector that is going to contain the local values plus ONE layer of
159aa219208SBarry Smith      ghosted values to draw the graphics from. We also need its corresponding DMDA (dac) that will
16047c6ae99SBarry Smith      update the local values pluse ONE layer of ghost values.
16147c6ae99SBarry Smith   */
16247c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)da,"GraphicsGhosted",(PetscObject*)&xlocal);CHKERRQ(ierr);
16347c6ae99SBarry Smith   if (!xlocal) {
164bff4a2f0SMatthew G. Knepley     if (bx !=  DM_BOUNDARY_NONE || by !=  DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) {
16547c6ae99SBarry Smith       /*
16647c6ae99SBarry Smith          if original da is not of stencil width one, or periodic or not a box stencil then
167aa219208SBarry Smith          create a special DMDA to handle one level of ghost points for graphics
16847c6ae99SBarry Smith       */
169bff4a2f0SMatthew 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);
170aa219208SBarry Smith       ierr = PetscInfo(da,"Creating auxilary DMDA for managing graphics ghost points\n");CHKERRQ(ierr);
17147c6ae99SBarry Smith     } else {
17247c6ae99SBarry Smith       /* otherwise we can use the da we already have */
17347c6ae99SBarry Smith       dac = da;
17447c6ae99SBarry Smith     }
17547c6ae99SBarry Smith     /* create local vector for holding ghosted values used in graphics */
176564755cdSBarry Smith     ierr = DMCreateLocalVector(dac,&xlocal);CHKERRQ(ierr);
17747c6ae99SBarry Smith     if (dac != da) {
178aa219208SBarry Smith       /* don't keep any public reference of this DMDA, is is only available through xlocal */
179f7923d8aSBarry Smith       ierr = PetscObjectDereference((PetscObject)dac);CHKERRQ(ierr);
18047c6ae99SBarry Smith     } else {
18147c6ae99SBarry Smith       /* remove association between xlocal and da, because below we compose in the opposite
18247c6ae99SBarry Smith          direction and if we left this connect we'd get a loop, so the objects could
18347c6ae99SBarry Smith          never be destroyed */
184c688c046SMatthew G Knepley       ierr = PetscObjectRemoveReference((PetscObject)xlocal,"__PETSc_dm");CHKERRQ(ierr);
18547c6ae99SBarry Smith     }
18647c6ae99SBarry Smith     ierr = PetscObjectCompose((PetscObject)da,"GraphicsGhosted",(PetscObject)xlocal);CHKERRQ(ierr);
18747c6ae99SBarry Smith     ierr = PetscObjectDereference((PetscObject)xlocal);CHKERRQ(ierr);
18847c6ae99SBarry Smith   } else {
189bff4a2f0SMatthew G. Knepley     if (bx !=  DM_BOUNDARY_NONE || by !=  DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) {
190c688c046SMatthew G Knepley       ierr = VecGetDM(xlocal, &dac);CHKERRQ(ierr);
191f7923d8aSBarry Smith     } else {
192f7923d8aSBarry Smith       dac = da;
19347c6ae99SBarry Smith     }
19447c6ae99SBarry Smith   }
19547c6ae99SBarry Smith 
19647c6ae99SBarry Smith   /*
19747c6ae99SBarry Smith       Get local (ghosted) values of vector
19847c6ae99SBarry Smith   */
1999a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr);
2009a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr);
2015edff71fSBarry Smith   ierr = VecGetArrayRead(xlocal,&zctx.v);CHKERRQ(ierr);
20247c6ae99SBarry Smith 
20347c6ae99SBarry Smith   /* get coordinates of nodes */
2046636e97aSMatthew G Knepley   ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
20547c6ae99SBarry Smith   if (!xcoor) {
206aa219208SBarry Smith     ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,0.0);CHKERRQ(ierr);
2076636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
20847c6ae99SBarry Smith   }
20947c6ae99SBarry Smith 
21047c6ae99SBarry Smith   /*
21147c6ae99SBarry Smith       Determine the min and max coordinates in plot
21247c6ae99SBarry Smith   */
213*e5ab1681SLisandro Dalcin   ierr = VecStrideMin(xcoor,0,NULL,&zctx.xmin);CHKERRQ(ierr);
214*e5ab1681SLisandro Dalcin   ierr = VecStrideMax(xcoor,0,NULL,&zctx.xmax);CHKERRQ(ierr);
215*e5ab1681SLisandro Dalcin   ierr = VecStrideMin(xcoor,1,NULL,&zctx.ymin);CHKERRQ(ierr);
216*e5ab1681SLisandro Dalcin   ierr = VecStrideMax(xcoor,1,NULL,&zctx.ymax);CHKERRQ(ierr);
217*e5ab1681SLisandro Dalcin   coors[0] = zctx.xmin - .05*(zctx.xmax- zctx.xmin); coors[2] = zctx.xmax + .05*(zctx.xmax - zctx.xmin);
218*e5ab1681SLisandro Dalcin   coors[1] = zctx.ymin - .05*(zctx.ymax- zctx.ymin); coors[3] = zctx.ymax + .05*(zctx.ymax - zctx.ymin);
219c5929fdfSBarry Smith   ierr = PetscOptionsGetRealArray(NULL,NULL,"-draw_coordinates",coors,&ncoors,NULL);CHKERRQ(ierr);
220*e5ab1681SLisandro 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);
221a74ba6f7SBarry Smith 
22247c6ae99SBarry Smith   /*
22347c6ae99SBarry Smith        get local ghosted version of coordinates
22447c6ae99SBarry Smith   */
22547c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr);
22647c6ae99SBarry Smith   if (!xcoorl) {
227aa219208SBarry Smith     /* create DMDA to get local version of graphics */
228bff4a2f0SMatthew 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);
229aa219208SBarry Smith     ierr = PetscInfo(dag,"Creating auxilary DMDA for managing graphics coordinates ghost points\n");CHKERRQ(ierr);
230564755cdSBarry Smith     ierr = DMCreateLocalVector(dag,&xcoorl);CHKERRQ(ierr);
23147c6ae99SBarry Smith     ierr = PetscObjectCompose((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject)xcoorl);CHKERRQ(ierr);
232f7923d8aSBarry Smith     ierr = PetscObjectDereference((PetscObject)dag);CHKERRQ(ierr);
23347c6ae99SBarry Smith     ierr = PetscObjectDereference((PetscObject)xcoorl);CHKERRQ(ierr);
23447c6ae99SBarry Smith   } else {
235c688c046SMatthew G Knepley     ierr = VecGetDM(xcoorl,&dag);CHKERRQ(ierr);
23647c6ae99SBarry Smith   }
2379a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
2389a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
2395edff71fSBarry Smith   ierr = VecGetArrayRead(xcoorl,&zctx.xy);CHKERRQ(ierr);
24047c6ae99SBarry Smith 
24147c6ae99SBarry Smith   /*
24247c6ae99SBarry Smith         Get information about size of area each processor must do graphics for
24347c6ae99SBarry Smith   */
244*e5ab1681SLisandro Dalcin   ierr = DMDAGetInfo(dac,NULL,&M,&N,NULL,NULL,NULL,NULL,&zctx.dof,NULL,&bx,&by,NULL,NULL);CHKERRQ(ierr);
245*e5ab1681SLisandro Dalcin   ierr = DMDAGetGhostCorners(dac,NULL,NULL,NULL,&zctx.m,&zctx.n,NULL);CHKERRQ(ierr);
246c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-draw_contour_grid",&zctx.showgrid,NULL);CHKERRQ(ierr);
2474e6118eeSBarry Smith 
2484e6118eeSBarry Smith   ierr = DMDASelectFields(da,&ndisplayfields,&displayfields);CHKERRQ(ierr);
24947c6ae99SBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
250c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-draw_ports",&useports,NULL);CHKERRQ(ierr);
25147c6ae99SBarry Smith   if (useports || format == PETSC_VIEWER_DRAW_PORTS) {
252*e5ab1681SLisandro Dalcin     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
2535b399a63SLisandro Dalcin     ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
2545b399a63SLisandro Dalcin     ierr = PetscDrawClear(draw);CHKERRQ(ierr);
25520d0051dSBarry Smith     ierr = PetscDrawViewPortsCreate(draw,ndisplayfields,&ports);CHKERRQ(ierr);
256*e5ab1681SLisandro Dalcin   }
257109c9344SBarry Smith   ierr = DMDAGetCoordinateName(da,0,&zctx.name0);CHKERRQ(ierr);
258109c9344SBarry Smith   ierr = DMDAGetCoordinateName(da,1,&zctx.name1);CHKERRQ(ierr);
25920d0051dSBarry Smith 
260*e5ab1681SLisandro Dalcin   /* loop over each field; drawing each in a different window */
26120d0051dSBarry Smith   for (i=0; i<ndisplayfields; i++) {
26220d0051dSBarry Smith     zctx.k = displayfields[i];
263*e5ab1681SLisandro Dalcin     if (useports || format == PETSC_VIEWER_DRAW_PORTS) {
26420d0051dSBarry Smith       ierr = PetscDrawViewPortsSet(ports,i);CHKERRQ(ierr);
2659332fd86SBarry Smith     } else {
266*e5ab1681SLisandro Dalcin       const char *title;
2679332fd86SBarry Smith       ierr = PetscViewerDrawGetDraw(viewer,i,&draw);CHKERRQ(ierr);
2685b399a63SLisandro Dalcin       ierr = DMDAGetFieldName(da,zctx.k,&title);CHKERRQ(ierr);
2695b399a63SLisandro Dalcin       if (title) {ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr);}
270*e5ab1681SLisandro Dalcin     }
2715b399a63SLisandro Dalcin 
272*e5ab1681SLisandro Dalcin     /* determine the min and max value in plot */
2730298fd71SBarry Smith     ierr = VecStrideMin(xin,zctx.k,NULL,&zctx.min);CHKERRQ(ierr);
2740298fd71SBarry Smith     ierr = VecStrideMax(xin,zctx.k,NULL,&zctx.max);CHKERRQ(ierr);
27567dd0837SBarry Smith     if (zctx.k < nbounds) {
276f3f0eb19SBarry Smith       zctx.min = bounds[2*zctx.k];
277f3f0eb19SBarry Smith       zctx.max = bounds[2*zctx.k+1];
27867dd0837SBarry Smith     }
27947c6ae99SBarry Smith     if (zctx.min == zctx.max) {
28047c6ae99SBarry Smith       zctx.min -= 1.e-12;
28147c6ae99SBarry Smith       zctx.max += 1.e-12;
28247c6ae99SBarry Smith     }
28357622a8eSBarry Smith     ierr = PetscInfo2(da,"DMDA 2d contour plot min %g max %g\n",(double)zctx.min,(double)zctx.max);CHKERRQ(ierr);
28447c6ae99SBarry Smith 
28547c6ae99SBarry Smith     ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
28647c6ae99SBarry Smith     if (popup) {ierr = PetscDrawScalePopup(popup,zctx.min,zctx.max);CHKERRQ(ierr);}
287*e5ab1681SLisandro Dalcin     ierr = PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr);
28847c6ae99SBarry Smith     ierr = PetscDrawZoom(draw,VecView_MPI_Draw_DA2d_Zoom,&zctx);CHKERRQ(ierr);
28947c6ae99SBarry Smith   }
2906bf464f9SBarry Smith   ierr = PetscDrawViewPortsDestroy(ports);CHKERRQ(ierr);
291*e5ab1681SLisandro Dalcin   ierr = PetscFree(displayfields);CHKERRQ(ierr);
29247c6ae99SBarry Smith 
2935edff71fSBarry Smith   ierr = VecRestoreArrayRead(xcoorl,&zctx.xy);CHKERRQ(ierr);
2945edff71fSBarry Smith   ierr = VecRestoreArrayRead(xlocal,&zctx.v);CHKERRQ(ierr);
29547c6ae99SBarry Smith   PetscFunctionReturn(0);
29647c6ae99SBarry Smith }
29747c6ae99SBarry Smith 
2980da24e51SJuha Jäykkä #if defined(PETSC_HAVE_HDF5)
2990da24e51SJuha Jäykkä #undef __FUNCT__
3000da24e51SJuha Jäykkä #define __FUNCT__ "VecGetHDF5ChunkSize"
301c73cfb54SMatthew G. Knepley static PetscErrorCode VecGetHDF5ChunkSize(DM_DA *da, Vec xin, PetscInt dimension, PetscInt timestep, hsize_t *chunkDims)
3020da24e51SJuha Jäykkä {
303d4190030SJed Brown   PetscMPIInt    comm_size;
304d4190030SJed Brown   PetscErrorCode ierr;
305561a82e7SJed Brown   hsize_t        chunk_size, target_size, dim;
306561a82e7SJed Brown   hsize_t        vec_size = sizeof(PetscScalar)*da->M*da->N*da->P*da->w;
307561a82e7SJed Brown   hsize_t        avg_local_vec_size,KiB = 1024,MiB = KiB*KiB,GiB = MiB*KiB,min_size = MiB;
308561a82e7SJed Brown   hsize_t        max_chunks = 64*KiB;                                              /* HDF5 internal limitation */
309561a82e7SJed Brown   hsize_t        max_chunk_size = 4*GiB;                                           /* HDF5 internal limitation */
31084179ffaSJed Brown   PetscInt       zslices=da->p, yslices=da->n, xslices=da->m;
3110da24e51SJuha Jäykkä 
312cb5c4f94SJuha Jäykkä   PetscFunctionBegin;
313d4190030SJed Brown   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)xin), &comm_size);CHKERRQ(ierr);
3140da24e51SJuha 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 */
3150da24e51SJuha Jäykkä 
316794a961bSBarry 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)))));
3177d310018SBarry 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 */
3187d310018SBarry 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);
3190da24e51SJuha Jäykkä 
320cb5c4f94SJuha Jäykkä   /*
321cb5c4f94SJuha Jäykkä    if size/rank > max_chunk_size, we need radical measures: even going down to
322cb5c4f94SJuha Jäykkä    avg_local_vec_size is not enough, so we simply use chunk size of 4 GiB no matter
323cb5c4f94SJuha Jäykkä    what, composed in the most efficient way possible.
324cb5c4f94SJuha Jäykkä    N.B. this minimises the number of chunks, which may or may not be the optimal
325cb5c4f94SJuha Jäykkä    solution. In a BG, for example, the optimal solution is probably to make # chunks = #
326cb5c4f94SJuha Jäykkä    IO nodes involved, but this author has no access to a BG to figure out how to
327cb5c4f94SJuha Jäykkä    reliably find the right number. And even then it may or may not be enough.
328cb5c4f94SJuha Jäykkä    */
3290da24e51SJuha Jäykkä   if (avg_local_vec_size > max_chunk_size) {
330cb5c4f94SJuha Jäykkä     /* check if we can just split local z-axis: is that enough? */
33184179ffaSJed Brown     zslices = (PetscInt)ceil(vec_size*1.0/(da->p*max_chunk_size))*zslices;
3320da24e51SJuha Jäykkä     if (zslices > da->P) {
333cb5c4f94SJuha Jäykkä       /* lattice is too large in xy-directions, splitting z only is not enough */
3340da24e51SJuha Jäykkä       zslices = da->P;
33584179ffaSJed Brown       yslices= (PetscInt)ceil(vec_size*1.0/(zslices*da->n*max_chunk_size))*yslices;
3360da24e51SJuha Jäykkä       if (yslices > da->N) {
337cb5c4f94SJuha Jäykkä 	/* lattice is too large in x-direction, splitting along z, y is not enough */
3380da24e51SJuha Jäykkä 	yslices = da->N;
33984179ffaSJed Brown 	xslices= (PetscInt)ceil(vec_size*1.0/(zslices*yslices*da->m*max_chunk_size))*xslices;
3400da24e51SJuha Jäykkä       }
3410da24e51SJuha Jäykkä     }
3420da24e51SJuha Jäykkä     dim = 0;
3430da24e51SJuha Jäykkä     if (timestep >= 0) {
3440da24e51SJuha Jäykkä       ++dim;
3450da24e51SJuha Jäykkä     }
346cb5c4f94SJuha Jäykkä     /* prefer to split z-axis, even down to planar slices */
347c73cfb54SMatthew G. Knepley     if (dimension == 3) {
348cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->P/zslices;
349cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->N/yslices;
350cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->M/xslices;
3510da24e51SJuha Jäykkä     } else {
352cb5c4f94SJuha Jäykkä       /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */
353cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->N/yslices;
354cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->M/xslices;
3550da24e51SJuha Jäykkä     }
3560da24e51SJuha 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);
3570da24e51SJuha Jäykkä   } else {
3580da24e51SJuha Jäykkä     if (target_size < chunk_size) {
359cb5c4f94SJuha Jäykkä       /* only change the defaults if target_size < chunk_size */
3600da24e51SJuha Jäykkä       dim = 0;
3610da24e51SJuha Jäykkä       if (timestep >= 0) {
3620da24e51SJuha Jäykkä 	++dim;
3630da24e51SJuha Jäykkä       }
364cb5c4f94SJuha Jäykkä       /* prefer to split z-axis, even down to planar slices */
365c73cfb54SMatthew G. Knepley       if (dimension == 3) {
366cb5c4f94SJuha Jäykkä 	/* try splitting the z-axis to core-size bits, i.e. divide chunk size by # comm_size in z-direction */
3670da24e51SJuha Jäykkä 	if (target_size >= chunk_size/da->p) {
368cb5c4f94SJuha Jäykkä 	  /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */
3690da24e51SJuha Jäykkä 	  chunkDims[dim] = (hsize_t) ceil(da->P*1.0/da->p);
3700da24e51SJuha Jäykkä 	} else {
371cb5c4f94SJuha Jäykkä 	  /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be
372cb5c4f94SJuha Jäykkä            radical and let everyone write all they've got */
3730da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->P*1.0/da->p);
3740da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n);
3750da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m);
3760da24e51SJuha Jäykkä 	}
3770da24e51SJuha Jäykkä       } else {
378cb5c4f94SJuha Jäykkä 	/* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */
3790da24e51SJuha Jäykkä 	if (target_size >= chunk_size/da->n) {
380cb5c4f94SJuha Jäykkä 	  /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */
3810da24e51SJuha Jäykkä 	  chunkDims[dim] = (hsize_t) ceil(da->N*1.0/da->n);
3820da24e51SJuha Jäykkä 	} else {
383cb5c4f94SJuha Jäykkä 	  /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be
384cb5c4f94SJuha Jäykkä 	   radical and let everyone write all they've got */
3850da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n);
3860da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m);
3870da24e51SJuha Jäykkä 	}
3880da24e51SJuha Jäykkä 
3890da24e51SJuha Jäykkä       }
3900da24e51SJuha 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);
3910da24e51SJuha Jäykkä     } else {
392cb5c4f94SJuha Jäykkä       /* precomputed chunks are fine, we don't need to do anything */
3930da24e51SJuha Jäykkä     }
3940da24e51SJuha Jäykkä   }
3950da24e51SJuha Jäykkä   PetscFunctionReturn(0);
3960da24e51SJuha Jäykkä }
3970da24e51SJuha Jäykkä #endif
39847c6ae99SBarry Smith 
39947c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
40047c6ae99SBarry Smith #undef __FUNCT__
40147c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_HDF5_DA"
40247c6ae99SBarry Smith PetscErrorCode VecView_MPI_HDF5_DA(Vec xin,PetscViewer viewer)
40347c6ae99SBarry Smith {
4049b2a5a86SJed Brown   DM                dm;
4059b2a5a86SJed Brown   DM_DA             *da;
40647c6ae99SBarry Smith   hid_t             filespace;  /* file dataspace identifier */
4078e2ae6d7SMichael Kraus   hid_t             chunkspace; /* chunk dataset property identifier */
40847c6ae99SBarry Smith   hid_t             plist_id;   /* property list identifier */
40947c6ae99SBarry Smith   hid_t             dset_id;    /* dataset identifier */
41047c6ae99SBarry Smith   hid_t             memspace;   /* memory dataspace identifier */
41147c6ae99SBarry Smith   hid_t             file_id;
41215214e8eSMatthew G Knepley   hid_t             group;
4139a0502c6SHåkon Strandenes   hid_t             memscalartype; /* scalar type for mem (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
4149a0502c6SHåkon Strandenes   hid_t             filescalartype; /* scalar type for file (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
415d9a4edebSJed Brown   hsize_t           dim;
4160da24e51SJuha 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  */
417c73cfb54SMatthew G. Knepley   PetscInt          timestep, dimension;
4185edff71fSBarry Smith   const PetscScalar *x;
4198e2ae6d7SMichael Kraus   const char        *vecname;
42015214e8eSMatthew G Knepley   PetscErrorCode    ierr;
42182ea9b62SBarry Smith   PetscBool         dim2;
4229a0502c6SHåkon Strandenes   PetscBool         spoutput;
42347c6ae99SBarry Smith 
42447c6ae99SBarry Smith   PetscFunctionBegin;
42515214e8eSMatthew G Knepley   ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr);
4268e2ae6d7SMichael Kraus   ierr = PetscViewerHDF5GetTimestep(viewer, &timestep);CHKERRQ(ierr);
42782ea9b62SBarry Smith   ierr = PetscViewerHDF5GetBaseDimension2(viewer,&dim2);CHKERRQ(ierr);
4289a0502c6SHåkon Strandenes   ierr = PetscViewerHDF5GetSPOutput(viewer,&spoutput);CHKERRQ(ierr);
42915214e8eSMatthew G Knepley 
430c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&dm);CHKERRQ(ierr);
431ce94432eSBarry Smith   if (!dm) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
4329b2a5a86SJed Brown   da = (DM_DA*)dm->data;
433c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dimension);CHKERRQ(ierr);
43447c6ae99SBarry Smith 
4358e2ae6d7SMichael Kraus   /* Create the dataspace for the dataset.
4368e2ae6d7SMichael Kraus    *
4378e2ae6d7SMichael Kraus    * dims - holds the current dimensions of the dataset
4388e2ae6d7SMichael Kraus    *
4398e2ae6d7SMichael Kraus    * maxDims - holds the maximum dimensions of the dataset (unlimited
4408e2ae6d7SMichael Kraus    * for the number of time steps with the current dimensions for the
4418e2ae6d7SMichael Kraus    * other dimensions; so only additional time steps can be added).
4428e2ae6d7SMichael Kraus    *
4438e2ae6d7SMichael Kraus    * chunkDims - holds the size of a single time step (required to
4448e2ae6d7SMichael Kraus    * permit extending dataset).
4458e2ae6d7SMichael Kraus    */
4468e2ae6d7SMichael Kraus   dim = 0;
4478e2ae6d7SMichael Kraus   if (timestep >= 0) {
4488e2ae6d7SMichael Kraus     dims[dim]      = timestep+1;
4498e2ae6d7SMichael Kraus     maxDims[dim]   = H5S_UNLIMITED;
4508e2ae6d7SMichael Kraus     chunkDims[dim] = 1;
4518e2ae6d7SMichael Kraus     ++dim;
4528e2ae6d7SMichael Kraus   }
453c73cfb54SMatthew G. Knepley   if (dimension == 3) {
454acba2ac6SBarry Smith     ierr           = PetscHDF5IntCast(da->P,dims+dim);CHKERRQ(ierr);
4558e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4568e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4578e2ae6d7SMichael Kraus     ++dim;
4588e2ae6d7SMichael Kraus   }
459c73cfb54SMatthew G. Knepley   if (dimension > 1) {
460acba2ac6SBarry Smith     ierr           = PetscHDF5IntCast(da->N,dims+dim);CHKERRQ(ierr);
4618e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4628e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4638e2ae6d7SMichael Kraus     ++dim;
4648e2ae6d7SMichael Kraus   }
465acba2ac6SBarry Smith   ierr           = PetscHDF5IntCast(da->M,dims+dim);CHKERRQ(ierr);
4668e2ae6d7SMichael Kraus   maxDims[dim]   = dims[dim];
4678e2ae6d7SMichael Kraus   chunkDims[dim] = dims[dim];
4688e2ae6d7SMichael Kraus   ++dim;
46982ea9b62SBarry Smith   if (da->w > 1 || dim2) {
470acba2ac6SBarry Smith     ierr           = PetscHDF5IntCast(da->w,dims+dim);CHKERRQ(ierr);
4718e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4728e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4738e2ae6d7SMichael Kraus     ++dim;
4748e2ae6d7SMichael Kraus   }
47547c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
4768e2ae6d7SMichael Kraus   dims[dim]      = 2;
4778e2ae6d7SMichael Kraus   maxDims[dim]   = dims[dim];
4788e2ae6d7SMichael Kraus   chunkDims[dim] = dims[dim];
4798e2ae6d7SMichael Kraus   ++dim;
48047c6ae99SBarry Smith #endif
4810da24e51SJuha Jäykkä 
482c73cfb54SMatthew G. Knepley   ierr = VecGetHDF5ChunkSize(da, xin, dimension, timestep, chunkDims); CHKERRQ(ierr);
4830da24e51SJuha Jäykkä 
484729ab6d9SBarry Smith   PetscStackCallHDF5Return(filespace,H5Screate_simple,(dim, dims, maxDims));
48547c6ae99SBarry Smith 
48615214e8eSMatthew G Knepley #if defined(PETSC_USE_REAL_SINGLE)
4879a0502c6SHåkon Strandenes   memscalartype = H5T_NATIVE_FLOAT;
4889a0502c6SHåkon Strandenes   filescalartype = H5T_NATIVE_FLOAT;
48915214e8eSMatthew G Knepley #elif defined(PETSC_USE_REAL___FLOAT128)
49015214e8eSMatthew G Knepley #error "HDF5 output with 128 bit floats not supported."
49115214e8eSMatthew G Knepley #else
4929a0502c6SHåkon Strandenes   memscalartype = H5T_NATIVE_DOUBLE;
4939a0502c6SHåkon Strandenes   if (spoutput == PETSC_TRUE) filescalartype = H5T_NATIVE_FLOAT;
4949a0502c6SHåkon Strandenes   else filescalartype = H5T_NATIVE_DOUBLE;
49515214e8eSMatthew G Knepley #endif
49615214e8eSMatthew G Knepley 
49747c6ae99SBarry Smith   /* Create the dataset with default properties and close filespace */
49847c6ae99SBarry Smith   ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr);
49915214e8eSMatthew G Knepley   if (!H5Lexists(group, vecname, H5P_DEFAULT)) {
5008e2ae6d7SMichael Kraus     /* Create chunk */
501729ab6d9SBarry Smith     PetscStackCallHDF5Return(chunkspace,H5Pcreate,(H5P_DATASET_CREATE));
502729ab6d9SBarry Smith     PetscStackCallHDF5(H5Pset_chunk,(chunkspace, dim, chunkDims));
5038e2ae6d7SMichael Kraus 
50447c6ae99SBarry Smith #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
5059a0502c6SHåkon Strandenes     PetscStackCallHDF5Return(dset_id,H5Dcreate2,(group, vecname, filescalartype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT));
50647c6ae99SBarry Smith #else
5079a0502c6SHåkon Strandenes     PetscStackCallHDF5Return(dset_id,H5Dcreate,(group, vecname, filescalartype, filespace, H5P_DEFAULT));
50847c6ae99SBarry Smith #endif
50915214e8eSMatthew G Knepley   } else {
510729ab6d9SBarry Smith     PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT));
511729ab6d9SBarry Smith     PetscStackCallHDF5(H5Dset_extent,(dset_id, dims));
51215214e8eSMatthew G Knepley   }
513729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
51447c6ae99SBarry Smith 
51547c6ae99SBarry Smith   /* Each process defines a dataset and writes it to the hyperslab in the file */
5168e2ae6d7SMichael Kraus   dim = 0;
5178e2ae6d7SMichael Kraus   if (timestep >= 0) {
5188e2ae6d7SMichael Kraus     offset[dim] = timestep;
5198e2ae6d7SMichael Kraus     ++dim;
5208e2ae6d7SMichael Kraus   }
521c73cfb54SMatthew G. Knepley   if (dimension == 3) {ierr = PetscHDF5IntCast(da->zs,offset + dim++);CHKERRQ(ierr);}
522c73cfb54SMatthew G. Knepley   if (dimension > 1)  {ierr = PetscHDF5IntCast(da->ys,offset + dim++);CHKERRQ(ierr);}
523acba2ac6SBarry Smith   ierr = PetscHDF5IntCast(da->xs/da->w,offset + dim++);CHKERRQ(ierr);
52482ea9b62SBarry Smith   if (da->w > 1 || dim2) offset[dim++] = 0;
52547c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5268e2ae6d7SMichael Kraus   offset[dim++] = 0;
52747c6ae99SBarry Smith #endif
5288e2ae6d7SMichael Kraus   dim = 0;
5298e2ae6d7SMichael Kraus   if (timestep >= 0) {
5308e2ae6d7SMichael Kraus     count[dim] = 1;
5318e2ae6d7SMichael Kraus     ++dim;
5328e2ae6d7SMichael Kraus   }
533c73cfb54SMatthew G. Knepley   if (dimension == 3) {ierr = PetscHDF5IntCast(da->ze - da->zs,count + dim++);CHKERRQ(ierr);}
534c73cfb54SMatthew G. Knepley   if (dimension > 1)  {ierr = PetscHDF5IntCast(da->ye - da->ys,count + dim++);CHKERRQ(ierr);}
535acba2ac6SBarry Smith   ierr = PetscHDF5IntCast((da->xe - da->xs)/da->w,count + dim++);CHKERRQ(ierr);
53682ea9b62SBarry Smith   if (da->w > 1 || dim2) {ierr = PetscHDF5IntCast(da->w,count + dim++);CHKERRQ(ierr);}
53747c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5388e2ae6d7SMichael Kraus   count[dim++] = 2;
53947c6ae99SBarry Smith #endif
540729ab6d9SBarry Smith   PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL));
541729ab6d9SBarry Smith   PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id));
542729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL));
54347c6ae99SBarry Smith 
54447c6ae99SBarry Smith   /* Create property list for collective dataset write */
545729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_DATASET_XFER));
54647c6ae99SBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
547729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_dxpl_mpio,(plist_id, H5FD_MPIO_COLLECTIVE));
54847c6ae99SBarry Smith #endif
54947c6ae99SBarry Smith   /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
55047c6ae99SBarry Smith 
5515edff71fSBarry Smith   ierr   = VecGetArrayRead(xin, &x);CHKERRQ(ierr);
5529a0502c6SHåkon Strandenes   PetscStackCallHDF5(H5Dwrite,(dset_id, memscalartype, memspace, filespace, plist_id, x));
553729ab6d9SBarry Smith   PetscStackCallHDF5(H5Fflush,(file_id, H5F_SCOPE_GLOBAL));
5545edff71fSBarry Smith   ierr   = VecRestoreArrayRead(xin, &x);CHKERRQ(ierr);
55547c6ae99SBarry Smith 
55647c6ae99SBarry Smith   /* Close/release resources */
55715214e8eSMatthew G Knepley   if (group != file_id) {
558729ab6d9SBarry Smith     PetscStackCallHDF5(H5Gclose,(group));
55915214e8eSMatthew G Knepley   }
560729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
561729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
562729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(memspace));
563729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dset_id));
56447c6ae99SBarry Smith   ierr   = PetscInfo1(xin,"Wrote Vec object with name %s\n",vecname);CHKERRQ(ierr);
56547c6ae99SBarry Smith   PetscFunctionReturn(0);
56647c6ae99SBarry Smith }
56747c6ae99SBarry Smith #endif
56847c6ae99SBarry Smith 
56909573ac7SBarry Smith extern PetscErrorCode VecView_MPI_Draw_DA1d(Vec,PetscViewer);
57047c6ae99SBarry Smith 
57147c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
57247c6ae99SBarry Smith #undef __FUNCT__
573aa219208SBarry Smith #define __FUNCT__ "DMDAArrayMPIIO"
574aa219208SBarry Smith static PetscErrorCode DMDAArrayMPIIO(DM da,PetscViewer viewer,Vec xin,PetscBool write)
57547c6ae99SBarry Smith {
57647c6ae99SBarry Smith   PetscErrorCode    ierr;
57747c6ae99SBarry Smith   MPI_File          mfdes;
57847c6ae99SBarry Smith   PetscMPIInt       gsizes[4],lsizes[4],lstarts[4],asiz,dof;
57947c6ae99SBarry Smith   MPI_Datatype      view;
58047c6ae99SBarry Smith   const PetscScalar *array;
58147c6ae99SBarry Smith   MPI_Offset        off;
58247c6ae99SBarry Smith   MPI_Aint          ub,ul;
58347c6ae99SBarry Smith   PetscInt          type,rows,vecrows,tr[2];
58447c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
5850c169764Sdmay   PetscBool         skipheader;
58647c6ae99SBarry Smith 
58747c6ae99SBarry Smith   PetscFunctionBegin;
58847c6ae99SBarry Smith   ierr = VecGetSize(xin,&vecrows);CHKERRQ(ierr);
5890c169764Sdmay   ierr = PetscViewerBinaryGetSkipHeader(viewer,&skipheader);CHKERRQ(ierr);
59047c6ae99SBarry Smith   if (!write) {
59147c6ae99SBarry Smith     /* Read vector header. */
5920c169764Sdmay     if (!skipheader) {
593060da220SMatthew G. Knepley       ierr = PetscViewerBinaryRead(viewer,tr,2,NULL,PETSC_INT);CHKERRQ(ierr);
59447c6ae99SBarry Smith       type = tr[0];
59547c6ae99SBarry Smith       rows = tr[1];
596ce94432eSBarry Smith       if (type != VEC_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Not vector next in file");
597ce94432eSBarry Smith       if (rows != vecrows) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Vector in file not same size as DMDA vector");
5980c169764Sdmay     }
59947c6ae99SBarry Smith   } else {
60047c6ae99SBarry Smith     tr[0] = VEC_FILE_CLASSID;
60147c6ae99SBarry Smith     tr[1] = vecrows;
6020c169764Sdmay     if (!skipheader) {
60347c6ae99SBarry Smith       ierr  = PetscViewerBinaryWrite(viewer,tr,2,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
60447c6ae99SBarry Smith     }
6050c169764Sdmay   }
60647c6ae99SBarry Smith 
6074dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->w,&dof);CHKERRQ(ierr);
6084dc2109aSBarry Smith   gsizes[0]  = dof;
6094dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->M,gsizes+1);CHKERRQ(ierr);
6104dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->N,gsizes+2);CHKERRQ(ierr);
611334634e2SJed Brown   ierr       = PetscMPIIntCast(dd->P,gsizes+3);CHKERRQ(ierr);
6124dc2109aSBarry Smith   lsizes[0]  = dof;
6134dc2109aSBarry Smith   ierr       = PetscMPIIntCast((dd->xe-dd->xs)/dof,lsizes+1);CHKERRQ(ierr);
6144dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->ye-dd->ys,lsizes+2);CHKERRQ(ierr);
6154dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->ze-dd->zs,lsizes+3);CHKERRQ(ierr);
6164dc2109aSBarry Smith   lstarts[0] = 0;
6174dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->xs/dof,lstarts+1);CHKERRQ(ierr);
6184dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->ys,lstarts+2);CHKERRQ(ierr);
6194dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->zs,lstarts+3);CHKERRQ(ierr);
620c73cfb54SMatthew G. Knepley   ierr       = MPI_Type_create_subarray(da->dim+1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view);CHKERRQ(ierr);
62147c6ae99SBarry Smith   ierr       = MPI_Type_commit(&view);CHKERRQ(ierr);
62247c6ae99SBarry Smith 
62347c6ae99SBarry Smith   ierr = PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);CHKERRQ(ierr);
62447c6ae99SBarry Smith   ierr = PetscViewerBinaryGetMPIIOOffset(viewer,&off);CHKERRQ(ierr);
62547c6ae99SBarry Smith   ierr = MPI_File_set_view(mfdes,off,MPIU_SCALAR,view,(char*)"native",MPI_INFO_NULL);CHKERRQ(ierr);
62647c6ae99SBarry Smith   ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr);
62747c6ae99SBarry Smith   asiz = lsizes[1]*(lsizes[2] > 0 ? lsizes[2] : 1)*(lsizes[3] > 0 ? lsizes[3] : 1)*dof;
62847c6ae99SBarry Smith   if (write) {
62947c6ae99SBarry Smith     ierr = MPIU_File_write_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr);
63047c6ae99SBarry Smith   } else {
63147c6ae99SBarry Smith     ierr = MPIU_File_read_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr);
63247c6ae99SBarry Smith   }
63347c6ae99SBarry Smith   ierr = MPI_Type_get_extent(view,&ul,&ub);CHKERRQ(ierr);
63447c6ae99SBarry Smith   ierr = PetscViewerBinaryAddMPIIOOffset(viewer,ub);CHKERRQ(ierr);
63547c6ae99SBarry Smith   ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr);
63647c6ae99SBarry Smith   ierr = MPI_Type_free(&view);CHKERRQ(ierr);
63747c6ae99SBarry Smith   PetscFunctionReturn(0);
63847c6ae99SBarry Smith }
63947c6ae99SBarry Smith #endif
64047c6ae99SBarry Smith 
64147c6ae99SBarry Smith #undef __FUNCT__
64247c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_DA"
6437087cfbeSBarry Smith PetscErrorCode  VecView_MPI_DA(Vec xin,PetscViewer viewer)
64447c6ae99SBarry Smith {
6459a42bb27SBarry Smith   DM                da;
64647c6ae99SBarry Smith   PetscErrorCode    ierr;
64747c6ae99SBarry Smith   PetscInt          dim;
64847c6ae99SBarry Smith   Vec               natural;
6494061b8bfSJed Brown   PetscBool         isdraw,isvtk;
65047c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
65147c6ae99SBarry Smith   PetscBool         ishdf5;
65247c6ae99SBarry Smith #endif
6533f3fd955SJed Brown   const char        *prefix,*name;
654a261c58fSBarry Smith   PetscViewerFormat format;
65547c6ae99SBarry Smith 
65647c6ae99SBarry Smith   PetscFunctionBegin;
657c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
658ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
659251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
660251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr);
66147c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
662251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr);
66347c6ae99SBarry Smith #endif
66447c6ae99SBarry Smith   if (isdraw) {
6651321219cSEthan Coon     ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
66647c6ae99SBarry Smith     if (dim == 1) {
66747c6ae99SBarry Smith       ierr = VecView_MPI_Draw_DA1d(xin,viewer);CHKERRQ(ierr);
66847c6ae99SBarry Smith     } else if (dim == 2) {
66947c6ae99SBarry Smith       ierr = VecView_MPI_Draw_DA2d(xin,viewer);CHKERRQ(ierr);
670ce94432eSBarry Smith     } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot graphically view vector associated with this dimensional DMDA %D",dim);
6714061b8bfSJed Brown   } else if (isvtk) {           /* Duplicate the Vec and hold a reference to the DM */
6724061b8bfSJed Brown     Vec Y;
6734061b8bfSJed Brown     ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr);
6744061b8bfSJed Brown     ierr = VecDuplicate(xin,&Y);CHKERRQ(ierr);
675b51b94faSJed Brown     if (((PetscObject)xin)->name) {
676b51b94faSJed Brown       /* If xin was named, copy the name over to Y. The duplicate names are safe because nobody else will ever see Y. */
677b51b94faSJed Brown       ierr = PetscObjectSetName((PetscObject)Y,((PetscObject)xin)->name);CHKERRQ(ierr);
678b51b94faSJed Brown     }
6794061b8bfSJed Brown     ierr = VecCopy(xin,Y);CHKERRQ(ierr);
68062b69a3fSMatthew G Knepley     ierr = PetscViewerVTKAddField(viewer,(PetscObject)da,DMDAVTKWriteAll,PETSC_VTK_POINT_FIELD,(PetscObject)Y);CHKERRQ(ierr);
68147c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
68247c6ae99SBarry Smith   } else if (ishdf5) {
68347c6ae99SBarry Smith     ierr = VecView_MPI_HDF5_DA(xin,viewer);CHKERRQ(ierr);
68447c6ae99SBarry Smith #endif
68547c6ae99SBarry Smith   } else {
68647c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
68747c6ae99SBarry Smith     PetscBool isbinary,isMPIIO;
68847c6ae99SBarry Smith 
689251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
69047c6ae99SBarry Smith     if (isbinary) {
691bc196f7cSDave May       ierr = PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO);CHKERRQ(ierr);
69247c6ae99SBarry Smith       if (isMPIIO) {
693aa219208SBarry Smith         ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_TRUE);CHKERRQ(ierr);
69447c6ae99SBarry Smith         PetscFunctionReturn(0);
69547c6ae99SBarry Smith       }
69647c6ae99SBarry Smith     }
69747c6ae99SBarry Smith #endif
69847c6ae99SBarry Smith 
69947c6ae99SBarry Smith     /* call viewer on natural ordering */
70047c6ae99SBarry Smith     ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr);
701aa219208SBarry Smith     ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr);
70247c6ae99SBarry Smith     ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr);
703aa219208SBarry Smith     ierr = DMDAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr);
704aa219208SBarry Smith     ierr = DMDAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr);
7053f3fd955SJed Brown     ierr = PetscObjectGetName((PetscObject)xin,&name);CHKERRQ(ierr);
7063f3fd955SJed Brown     ierr = PetscObjectSetName((PetscObject)natural,name);CHKERRQ(ierr);
707a261c58fSBarry Smith 
708a261c58fSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
709a261c58fSBarry Smith     if (format == PETSC_VIEWER_BINARY_MATLAB) {
710a261c58fSBarry Smith       /* temporarily remove viewer format so it won't trigger in the VecView() */
7116a9046bcSBarry Smith       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
712a261c58fSBarry Smith     }
713a261c58fSBarry Smith 
71447c6ae99SBarry Smith     ierr = VecView(natural,viewer);CHKERRQ(ierr);
715a261c58fSBarry Smith 
716a261c58fSBarry Smith     if (format == PETSC_VIEWER_BINARY_MATLAB) {
717a261c58fSBarry Smith       MPI_Comm    comm;
718a261c58fSBarry Smith       FILE        *info;
719a261c58fSBarry Smith       const char  *fieldname;
720da88d4d4SJed Brown       char        fieldbuf[256];
721a261c58fSBarry Smith       PetscInt    dim,ni,nj,nk,pi,pj,pk,dof,n;
722a261c58fSBarry Smith 
723a261c58fSBarry Smith       /* set the viewer format back into the viewer */
7246a9046bcSBarry Smith       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
725a261c58fSBarry Smith       ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
726a261c58fSBarry Smith       ierr = PetscViewerBinaryGetInfoPointer(viewer,&info);CHKERRQ(ierr);
727a261c58fSBarry Smith       ierr = DMDAGetInfo(da,&dim,&ni,&nj,&nk,&pi,&pj,&pk,&dof,0,0,0,0,0);CHKERRQ(ierr);
728da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n");CHKERRQ(ierr);
729da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#$$ tmp = PetscBinaryRead(fd); \n");CHKERRQ(ierr);
730da88d4d4SJed Brown       if (dim == 1) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d);\n",dof,ni);CHKERRQ(ierr); }
731da88d4d4SJed Brown       if (dim == 2) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d);\n",dof,ni,nj);CHKERRQ(ierr); }
732da88d4d4SJed Brown       if (dim == 3) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d,%d);\n",dof,ni,nj,nk);CHKERRQ(ierr); }
733a261c58fSBarry Smith 
734a261c58fSBarry Smith       for (n=0; n<dof; n++) {
735a261c58fSBarry Smith         ierr = DMDAGetFieldName(da,n,&fieldname);CHKERRQ(ierr);
736da88d4d4SJed Brown         if (!fieldname || !fieldname[0]) {
737da88d4d4SJed Brown           ierr = PetscSNPrintf(fieldbuf,sizeof fieldbuf,"field%D",n);CHKERRQ(ierr);
738da88d4d4SJed Brown           fieldname = fieldbuf;
739a261c58fSBarry Smith         }
740da88d4d4SJed Brown         if (dim == 1) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:))';\n",name,fieldname,n+1);CHKERRQ(ierr); }
741da88d4d4SJed Brown         if (dim == 2) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:,:))';\n",name,fieldname,n+1);CHKERRQ(ierr); }
742da88d4d4SJed 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);}
743a261c58fSBarry Smith       }
744da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#$$ clear tmp; \n");CHKERRQ(ierr);
745da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n");CHKERRQ(ierr);
746a261c58fSBarry Smith     }
747a261c58fSBarry Smith 
748fcfd50ebSBarry Smith     ierr = VecDestroy(&natural);CHKERRQ(ierr);
74947c6ae99SBarry Smith   }
75047c6ae99SBarry Smith   PetscFunctionReturn(0);
75147c6ae99SBarry Smith }
75247c6ae99SBarry Smith 
75347c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
75447c6ae99SBarry Smith #undef __FUNCT__
75547c6ae99SBarry Smith #define __FUNCT__ "VecLoad_HDF5_DA"
75647c6ae99SBarry Smith PetscErrorCode VecLoad_HDF5_DA(Vec xin, PetscViewer viewer)
75747c6ae99SBarry Smith {
7589a42bb27SBarry Smith   DM             da;
75947c6ae99SBarry Smith   PetscErrorCode ierr;
760ec7ae49cSHåkon Strandenes   hsize_t        dim,rdim;
761ec7ae49cSHåkon Strandenes   hsize_t        dims[6]={0},count[6]={0},offset[6]={0};
762ec7ae49cSHåkon Strandenes   PetscInt       dimension,timestep,dofInd;
76347c6ae99SBarry Smith   PetscScalar    *x;
76447c6ae99SBarry Smith   const char     *vecname;
76547c6ae99SBarry Smith   hid_t          filespace; /* file dataspace identifier */
76647c6ae99SBarry Smith   hid_t          plist_id;  /* property list identifier */
76747c6ae99SBarry Smith   hid_t          dset_id;   /* dataset identifier */
76847c6ae99SBarry Smith   hid_t          memspace;  /* memory dataspace identifier */
769bfd264e7SBarry Smith   hid_t          file_id,group;
7707bcbaff4SBarry Smith   hid_t          scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
7719c7c4993SBarry Smith   DM_DA          *dd;
772ec7ae49cSHåkon Strandenes   PetscBool      dim2 = PETSC_FALSE;
77347c6ae99SBarry Smith 
77447c6ae99SBarry Smith   PetscFunctionBegin;
7757bcbaff4SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
7767bcbaff4SBarry Smith   scalartype = H5T_NATIVE_FLOAT;
7777bcbaff4SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
7787bcbaff4SBarry Smith #error "HDF5 output with 128 bit floats not supported."
7797bcbaff4SBarry Smith #else
7807bcbaff4SBarry Smith   scalartype = H5T_NATIVE_DOUBLE;
7817bcbaff4SBarry Smith #endif
7827bcbaff4SBarry Smith 
783bfd264e7SBarry Smith   ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr);
784ec7ae49cSHåkon Strandenes   ierr = PetscViewerHDF5GetTimestep(viewer, &timestep);CHKERRQ(ierr);
785ec7ae49cSHåkon Strandenes   ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr);
786c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
7879c7c4993SBarry Smith   dd   = (DM_DA*)da->data;
788c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(da, &dimension);CHKERRQ(ierr);
78947c6ae99SBarry Smith 
790ec7ae49cSHåkon Strandenes   /* Open dataset */
79147c6ae99SBarry Smith #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
792729ab6d9SBarry Smith   PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT));
79347c6ae99SBarry Smith #else
794729ab6d9SBarry Smith   PetscStackCallHDF5Return(dset_id,H5Dopen,(group, vecname));
79547c6ae99SBarry Smith #endif
79647c6ae99SBarry Smith 
797ec7ae49cSHåkon Strandenes   /* Retrieve the dataspace for the dataset */
798ec7ae49cSHåkon Strandenes   PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id));
799ec7ae49cSHåkon Strandenes   PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(filespace, dims, NULL));
800ec7ae49cSHåkon Strandenes 
801ec7ae49cSHåkon Strandenes   /* Expected dimension for holding the dof's */
80247c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
803ec7ae49cSHåkon Strandenes   dofInd = rdim-2;
804ec7ae49cSHåkon Strandenes #else
805ec7ae49cSHåkon Strandenes   dofInd = rdim-1;
80647c6ae99SBarry Smith #endif
807ec7ae49cSHåkon Strandenes 
808ec7ae49cSHåkon Strandenes   /* The expected number of dimensions, assuming basedimension2 = false */
809ec7ae49cSHåkon Strandenes   dim = dimension;
810ec7ae49cSHåkon Strandenes   if (dd->w > 1) ++dim;
811ec7ae49cSHåkon Strandenes   if (timestep >= 0) ++dim;
81247c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
813ec7ae49cSHåkon Strandenes   ++dim;
81447c6ae99SBarry Smith #endif
815ec7ae49cSHåkon Strandenes 
816ec7ae49cSHåkon Strandenes   /* In this case the input dataset have one extra, unexpected dimension. */
817ec7ae49cSHåkon Strandenes   if (rdim == dim+1) {
818ec7ae49cSHåkon Strandenes     /* In this case the block size unity */
819ec7ae49cSHåkon Strandenes     if (dd->w == 1 && dims[dofInd] == 1) dim2 = PETSC_TRUE;
820ec7ae49cSHåkon Strandenes 
821ec7ae49cSHåkon Strandenes     /* Special error message for the case where dof does not match the input file */
822ec7ae49cSHå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);
823ec7ae49cSHåkon Strandenes 
824ec7ae49cSHåkon Strandenes   /* Other cases where rdim != dim cannot be handled currently */
8256c4ed002SBarry 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);
826ec7ae49cSHåkon Strandenes 
827ec7ae49cSHåkon Strandenes   /* Set up the hyperslab size */
828ec7ae49cSHåkon Strandenes   dim = 0;
829ec7ae49cSHåkon Strandenes   if (timestep >= 0) {
830ec7ae49cSHåkon Strandenes     offset[dim] = timestep;
831ec7ae49cSHåkon Strandenes     count[dim] = 1;
832ec7ae49cSHåkon Strandenes     ++dim;
833ec7ae49cSHåkon Strandenes   }
834ec7ae49cSHåkon Strandenes   if (dimension == 3) {
835ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->zs,offset + dim);CHKERRQ(ierr);
836ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->ze - dd->zs,count + dim);CHKERRQ(ierr);
837ec7ae49cSHåkon Strandenes     ++dim;
838ec7ae49cSHåkon Strandenes   }
839ec7ae49cSHåkon Strandenes   if (dimension > 1) {
840ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->ys,offset + dim);CHKERRQ(ierr);
841ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->ye - dd->ys,count + dim);CHKERRQ(ierr);
842ec7ae49cSHåkon Strandenes     ++dim;
843ec7ae49cSHåkon Strandenes   }
844ec7ae49cSHåkon Strandenes   ierr = PetscHDF5IntCast(dd->xs/dd->w,offset + dim);CHKERRQ(ierr);
845ec7ae49cSHåkon Strandenes   ierr = PetscHDF5IntCast((dd->xe - dd->xs)/dd->w,count + dim);CHKERRQ(ierr);
846ec7ae49cSHåkon Strandenes   ++dim;
847ec7ae49cSHåkon Strandenes   if (dd->w > 1 || dim2) {
848ec7ae49cSHåkon Strandenes     offset[dim] = 0;
849ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->w,count + dim);CHKERRQ(ierr);
850ec7ae49cSHåkon Strandenes     ++dim;
851ec7ae49cSHåkon Strandenes   }
852ec7ae49cSHåkon Strandenes #if defined(PETSC_USE_COMPLEX)
853ec7ae49cSHåkon Strandenes   offset[dim] = 0;
854ec7ae49cSHåkon Strandenes   count[dim] = 2;
855ec7ae49cSHåkon Strandenes   ++dim;
856ec7ae49cSHåkon Strandenes #endif
857ec7ae49cSHåkon Strandenes 
858ec7ae49cSHåkon Strandenes   /* Create the memory and filespace */
859729ab6d9SBarry Smith   PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL));
860729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL));
86147c6ae99SBarry Smith 
86247c6ae99SBarry Smith   /* Create property list for collective dataset write */
863729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_DATASET_XFER));
86447c6ae99SBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
865729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_dxpl_mpio,(plist_id, H5FD_MPIO_COLLECTIVE));
86647c6ae99SBarry Smith #endif
867ec7ae49cSHåkon Strandenes   /* To read dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
86847c6ae99SBarry Smith 
86947c6ae99SBarry Smith   ierr   = VecGetArray(xin, &x);CHKERRQ(ierr);
870a1dbdba5SBarry Smith   PetscStackCallHDF5(H5Dread,(dset_id, scalartype, memspace, filespace, plist_id, x));
87147c6ae99SBarry Smith   ierr   = VecRestoreArray(xin, &x);CHKERRQ(ierr);
87247c6ae99SBarry Smith 
87347c6ae99SBarry Smith   /* Close/release resources */
874bfd264e7SBarry Smith   if (group != file_id) {
875729ab6d9SBarry Smith     PetscStackCallHDF5(H5Gclose,(group));
876bfd264e7SBarry Smith   }
877729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
878729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
879729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(memspace));
880729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dset_id));
88147c6ae99SBarry Smith   PetscFunctionReturn(0);
88247c6ae99SBarry Smith }
88347c6ae99SBarry Smith #endif
88447c6ae99SBarry Smith 
88547c6ae99SBarry Smith #undef __FUNCT__
88647c6ae99SBarry Smith #define __FUNCT__ "VecLoad_Binary_DA"
88747c6ae99SBarry Smith PetscErrorCode VecLoad_Binary_DA(Vec xin, PetscViewer viewer)
88847c6ae99SBarry Smith {
8899a42bb27SBarry Smith   DM             da;
89047c6ae99SBarry Smith   PetscErrorCode ierr;
89147c6ae99SBarry Smith   Vec            natural;
89247c6ae99SBarry Smith   const char     *prefix;
89347c6ae99SBarry Smith   PetscInt       bs;
89447c6ae99SBarry Smith   PetscBool      flag;
89547c6ae99SBarry Smith   DM_DA          *dd;
89647c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
89747c6ae99SBarry Smith   PetscBool isMPIIO;
89847c6ae99SBarry Smith #endif
89947c6ae99SBarry Smith 
90047c6ae99SBarry Smith   PetscFunctionBegin;
901c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
90247c6ae99SBarry Smith   dd   = (DM_DA*)da->data;
90347c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
904bc196f7cSDave May   ierr = PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO);CHKERRQ(ierr);
90547c6ae99SBarry Smith   if (isMPIIO) {
906aa219208SBarry Smith     ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_FALSE);CHKERRQ(ierr);
90747c6ae99SBarry Smith     PetscFunctionReturn(0);
90847c6ae99SBarry Smith   }
90947c6ae99SBarry Smith #endif
91047c6ae99SBarry Smith 
91147c6ae99SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr);
912aa219208SBarry Smith   ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr);
91347c6ae99SBarry Smith   ierr = PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name);CHKERRQ(ierr);
91447c6ae99SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr);
9158aea9937SBarry Smith   ierr = VecLoad(natural,viewer);CHKERRQ(ierr);
916aa219208SBarry Smith   ierr = DMDANaturalToGlobalBegin(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr);
917aa219208SBarry Smith   ierr = DMDANaturalToGlobalEnd(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr);
918fcfd50ebSBarry Smith   ierr = VecDestroy(&natural);CHKERRQ(ierr);
919aa219208SBarry Smith   ierr = PetscInfo(xin,"Loading vector from natural ordering into DMDA\n");CHKERRQ(ierr);
920c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(NULL,((PetscObject)xin)->prefix,"-vecload_block_size",&bs,&flag);CHKERRQ(ierr);
92147c6ae99SBarry Smith   if (flag && bs != dd->w) {
922aa219208SBarry Smith     ierr = PetscInfo2(xin,"Block size in file %D not equal to DMDA's dof %D\n",bs,dd->w);CHKERRQ(ierr);
92347c6ae99SBarry Smith   }
92447c6ae99SBarry Smith   PetscFunctionReturn(0);
92547c6ae99SBarry Smith }
92647c6ae99SBarry Smith 
92747c6ae99SBarry Smith #undef __FUNCT__
92847c6ae99SBarry Smith #define __FUNCT__ "VecLoad_Default_DA"
9297087cfbeSBarry Smith PetscErrorCode  VecLoad_Default_DA(Vec xin, PetscViewer viewer)
93047c6ae99SBarry Smith {
93147c6ae99SBarry Smith   PetscErrorCode ierr;
9329a42bb27SBarry Smith   DM             da;
93347c6ae99SBarry Smith   PetscBool      isbinary;
93447c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
93547c6ae99SBarry Smith   PetscBool ishdf5;
93647c6ae99SBarry Smith #endif
93747c6ae99SBarry Smith 
93847c6ae99SBarry Smith   PetscFunctionBegin;
939c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
940ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
94147c6ae99SBarry Smith 
94247c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
943251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr);
94447c6ae99SBarry Smith #endif
945251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
94647c6ae99SBarry Smith 
94747c6ae99SBarry Smith   if (isbinary) {
94847c6ae99SBarry Smith     ierr = VecLoad_Binary_DA(xin,viewer);CHKERRQ(ierr);
94947c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
95047c6ae99SBarry Smith   } else if (ishdf5) {
95147c6ae99SBarry Smith     ierr = VecLoad_HDF5_DA(xin,viewer);CHKERRQ(ierr);
95247c6ae99SBarry Smith #endif
953d34fcf5fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for vector loading", ((PetscObject)viewer)->type_name);
95447c6ae99SBarry Smith   PetscFunctionReturn(0);
95547c6ae99SBarry Smith }
956