xref: /petsc/src/dm/impls/da/gr2.c (revision 84179ffa140a496dba545a3c72fc4347854e0f4b)
147c6ae99SBarry Smith 
247c6ae99SBarry Smith /*
3aa219208SBarry Smith    Plots vectors obtained with DMDACreate2d()
447c6ae99SBarry Smith */
547c6ae99SBarry Smith 
64035e84dSBarry Smith #include <petsc-private/dmdaimpl.h>      /*I  "petscdmda.h"   I*/
7b45d2f2cSJed Brown #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 {
1547c6ae99SBarry Smith   PetscInt    m,n,step,k;
1647c6ae99SBarry Smith   PetscReal   min,max,scale;
1747c6ae99SBarry Smith   PetscScalar *xy,*v;
1847c6ae99SBarry Smith   PetscBool   showgrid;
19109c9344SBarry Smith   const char  *name0,*name1;
2047c6ae99SBarry Smith } ZoomCtx;
2147c6ae99SBarry Smith 
2247c6ae99SBarry Smith /*
2347c6ae99SBarry Smith        This does the drawing for one particular field
2447c6ae99SBarry Smith     in one particular set of coordinates. It is a callback
2547c6ae99SBarry Smith     called from PetscDrawZoom()
2647c6ae99SBarry Smith */
2747c6ae99SBarry Smith #undef __FUNCT__
2847c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_Draw_DA2d_Zoom"
2947c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA2d_Zoom(PetscDraw draw,void *ctx)
3047c6ae99SBarry Smith {
3147c6ae99SBarry Smith   ZoomCtx        *zctx = (ZoomCtx*)ctx;
3247c6ae99SBarry Smith   PetscErrorCode ierr;
3347c6ae99SBarry Smith   PetscInt       m,n,i,j,k,step,id,c1,c2,c3,c4;
340e4fe250SBarry Smith   PetscReal      s,min,max,x1,x2,x3,x4,y_1,y2,y3,y4,xmin = PETSC_MAX_REAL,xmax = PETSC_MIN_REAL,ymin = PETSC_MAX_REAL,ymax = PETSC_MIN_REAL;
350e4fe250SBarry Smith   PetscReal      xminf,xmaxf,yminf,ymaxf,w;
3647c6ae99SBarry Smith   PetscScalar    *v,*xy;
370e4fe250SBarry Smith   char           value[16];
380e4fe250SBarry Smith   size_t         len;
3947c6ae99SBarry Smith 
4047c6ae99SBarry Smith   PetscFunctionBegin;
4147c6ae99SBarry Smith   m    = zctx->m;
4247c6ae99SBarry Smith   n    = zctx->n;
4347c6ae99SBarry Smith   step = zctx->step;
4447c6ae99SBarry Smith   k    = zctx->k;
4547c6ae99SBarry Smith   v    = zctx->v;
4647c6ae99SBarry Smith   xy   = zctx->xy;
4747c6ae99SBarry Smith   s    = zctx->scale;
4847c6ae99SBarry Smith   min  = zctx->min;
49f3f0eb19SBarry Smith   max  = zctx->max;
5047c6ae99SBarry Smith 
5147c6ae99SBarry Smith   /* PetscDraw the contour plot patch */
5247c6ae99SBarry Smith   for (j=0; j<n-1; j++) {
5347c6ae99SBarry Smith     for (i=0; i<m-1; i++) {
540e4fe250SBarry Smith       id   = i+j*m;
550e4fe250SBarry Smith       x1   = PetscRealPart(xy[2*id]);
560e4fe250SBarry Smith       y_1  = PetscRealPart(xy[2*id+1]);
570e4fe250SBarry Smith       c1   = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscClipInterval(PetscRealPart(v[k+step*id]),min,max)-min));
580e4fe250SBarry Smith       xmin = PetscMin(xmin,x1);
590e4fe250SBarry Smith       ymin = PetscMin(ymin,y_1);
600e4fe250SBarry Smith       xmax = PetscMax(xmax,x1);
610e4fe250SBarry Smith       ymax = PetscMax(ymax,y_1);
620e4fe250SBarry Smith 
630e4fe250SBarry Smith       id   = i+j*m+1;
640e4fe250SBarry Smith       x2   = PetscRealPart(xy[2*id]);
650e4fe250SBarry Smith       y2   = y_1;
660e4fe250SBarry Smith       c2   = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscClipInterval(PetscRealPart(v[k+step*id]),min,max)-min));
670e4fe250SBarry Smith       xmin = PetscMin(xmin,x2);
680e4fe250SBarry Smith       xmax = PetscMax(xmax,x2);
690e4fe250SBarry Smith 
700e4fe250SBarry Smith       id   = i+j*m+1+m;
710e4fe250SBarry Smith       x3   = x2;
720e4fe250SBarry Smith       y3   = PetscRealPart(xy[2*id+1]);
730e4fe250SBarry Smith       c3   = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscClipInterval(PetscRealPart(v[k+step*id]),min,max)-min));
740e4fe250SBarry Smith       ymin = PetscMin(ymin,y3);
750e4fe250SBarry Smith       ymax = PetscMax(ymax,y3);
760e4fe250SBarry Smith 
770e4fe250SBarry Smith       id = i+j*m+m;
780e4fe250SBarry Smith       x4 = x1;
790e4fe250SBarry Smith       y4 = y3;
800e4fe250SBarry Smith       c4 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscClipInterval(PetscRealPart(v[k+step*id]),min,max)-min));
81f3f0eb19SBarry Smith 
8247c6ae99SBarry Smith       ierr = PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);CHKERRQ(ierr);
8347c6ae99SBarry Smith       ierr = PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4);CHKERRQ(ierr);
8447c6ae99SBarry Smith       if (zctx->showgrid) {
8547c6ae99SBarry Smith         ierr = PetscDrawLine(draw,x1,y_1,x2,y2,PETSC_DRAW_BLACK);CHKERRQ(ierr);
8647c6ae99SBarry Smith         ierr = PetscDrawLine(draw,x2,y2,x3,y3,PETSC_DRAW_BLACK);CHKERRQ(ierr);
8747c6ae99SBarry Smith         ierr = PetscDrawLine(draw,x3,y3,x4,y4,PETSC_DRAW_BLACK);CHKERRQ(ierr);
8847c6ae99SBarry Smith         ierr = PetscDrawLine(draw,x4,y4,x1,y_1,PETSC_DRAW_BLACK);CHKERRQ(ierr);
8947c6ae99SBarry Smith       }
9047c6ae99SBarry Smith     }
9147c6ae99SBarry Smith   }
92109c9344SBarry Smith   if (zctx->name0) {
93109c9344SBarry Smith     PetscReal xl,yl,xr,yr,x,y;
94109c9344SBarry Smith     ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
95109c9344SBarry Smith     x    = xl + .3*(xr - xl);
96109c9344SBarry Smith     xl   = xl + .01*(xr - xl);
97109c9344SBarry Smith     y    = yr - .3*(yr - yl);
98109c9344SBarry Smith     yl   = yl + .01*(yr - yl);
99109c9344SBarry Smith     ierr = PetscDrawString(draw,x,yl,PETSC_DRAW_BLACK,zctx->name0);CHKERRQ(ierr);
100109c9344SBarry Smith     ierr = PetscDrawStringVertical(draw,xl,y,PETSC_DRAW_BLACK,zctx->name1);CHKERRQ(ierr);
101109c9344SBarry Smith   }
1020e4fe250SBarry Smith   /*
1030e4fe250SBarry Smith      Ideally we would use the PetscDrawAxis object to manage displaying the coordinate limits
1040e4fe250SBarry Smith      but that may require some refactoring.
1050e4fe250SBarry Smith   */
106ce94432eSBarry Smith   ierr = MPI_Allreduce(&xmin,&xminf,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)draw));CHKERRQ(ierr);
107ce94432eSBarry Smith   ierr = MPI_Allreduce(&xmax,&xmaxf,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)draw));CHKERRQ(ierr);
108ce94432eSBarry Smith   ierr = MPI_Allreduce(&ymin,&yminf,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)draw));CHKERRQ(ierr);
109ce94432eSBarry Smith   ierr = MPI_Allreduce(&ymax,&ymaxf,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)draw));CHKERRQ(ierr);
1100e4fe250SBarry Smith   ierr = PetscSNPrintf(value,16,"%f",xminf);CHKERRQ(ierr);
1110e4fe250SBarry Smith   ierr = PetscDrawString(draw,xminf,yminf - .05*(ymaxf - yminf),PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
1120e4fe250SBarry Smith   ierr = PetscSNPrintf(value,16,"%f",xmaxf);CHKERRQ(ierr);
1130e4fe250SBarry Smith   ierr = PetscStrlen(value,&len);CHKERRQ(ierr);
1140298fd71SBarry Smith   ierr = PetscDrawStringGetSize(draw,&w,NULL);CHKERRQ(ierr);
1150e4fe250SBarry Smith   ierr = PetscDrawString(draw,xmaxf - len*w,yminf - .05*(ymaxf - yminf),PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
1160e4fe250SBarry Smith   ierr = PetscSNPrintf(value,16,"%f",yminf);CHKERRQ(ierr);
1170e4fe250SBarry Smith   ierr = PetscDrawString(draw,xminf - .05*(xmaxf - xminf),yminf,PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
1180e4fe250SBarry Smith   ierr = PetscSNPrintf(value,16,"%f",ymaxf);CHKERRQ(ierr);
1190e4fe250SBarry Smith   ierr = PetscDrawString(draw,xminf - .05*(xmaxf - xminf),ymaxf,PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
12047c6ae99SBarry Smith   PetscFunctionReturn(0);
12147c6ae99SBarry Smith }
12247c6ae99SBarry Smith 
12347c6ae99SBarry Smith #undef __FUNCT__
12447c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_Draw_DA2d"
12547c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA2d(Vec xin,PetscViewer viewer)
12647c6ae99SBarry Smith {
1279a42bb27SBarry Smith   DM                 da,dac,dag;
12847c6ae99SBarry Smith   PetscErrorCode     ierr;
12947c6ae99SBarry Smith   PetscMPIInt        rank;
130a74ba6f7SBarry Smith   PetscInt           N,s,M,w,ncoors = 4;
13147c6ae99SBarry Smith   const PetscInt     *lx,*ly;
13247c6ae99SBarry Smith   PetscReal          coors[4],ymin,ymax,xmin,xmax;
13347c6ae99SBarry Smith   PetscDraw          draw,popup;
13447c6ae99SBarry Smith   PetscBool          isnull,useports = PETSC_FALSE;
13547c6ae99SBarry Smith   MPI_Comm           comm;
13647c6ae99SBarry Smith   Vec                xlocal,xcoor,xcoorl;
1371321219cSEthan Coon   DMDABoundaryType   bx,by;
138aa219208SBarry Smith   DMDAStencilType    st;
13947c6ae99SBarry Smith   ZoomCtx            zctx;
1400298fd71SBarry Smith   PetscDrawViewPorts *ports = NULL;
14147c6ae99SBarry Smith   PetscViewerFormat  format;
14220d0051dSBarry Smith   PetscInt           *displayfields;
14367dd0837SBarry Smith   PetscInt           ndisplayfields,i,nbounds;
14467dd0837SBarry Smith   const PetscReal    *bounds;
14547c6ae99SBarry Smith 
14647c6ae99SBarry Smith   PetscFunctionBegin;
14747c6ae99SBarry Smith   zctx.showgrid = PETSC_FALSE;
1488865f1eaSKarl Rupp 
14947c6ae99SBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
15047c6ae99SBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
15103193ff8SBarry Smith   ierr = PetscViewerDrawGetBounds(viewer,&nbounds,&bounds);CHKERRQ(ierr);
15247c6ae99SBarry Smith 
153c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
154ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
15547c6ae99SBarry Smith 
15647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr);
15747c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
15847c6ae99SBarry Smith 
1591321219cSEthan Coon   ierr = DMDAGetInfo(da,0,&M,&N,0,&zctx.m,&zctx.n,0,&w,&s,&bx,&by,0,&st);CHKERRQ(ierr);
1600298fd71SBarry Smith   ierr = DMDAGetOwnershipRanges(da,&lx,&ly,NULL);CHKERRQ(ierr);
16147c6ae99SBarry Smith 
16247c6ae99SBarry Smith   /*
16347c6ae99SBarry Smith         Obtain a sequential vector that is going to contain the local values plus ONE layer of
164aa219208SBarry Smith      ghosted values to draw the graphics from. We also need its corresponding DMDA (dac) that will
16547c6ae99SBarry Smith      update the local values pluse ONE layer of ghost values.
16647c6ae99SBarry Smith   */
16747c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)da,"GraphicsGhosted",(PetscObject*)&xlocal);CHKERRQ(ierr);
16847c6ae99SBarry Smith   if (!xlocal) {
169f7923d8aSBarry Smith     if (bx !=  DMDA_BOUNDARY_NONE || by !=  DMDA_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) {
17047c6ae99SBarry Smith       /*
17147c6ae99SBarry Smith          if original da is not of stencil width one, or periodic or not a box stencil then
172aa219208SBarry Smith          create a special DMDA to handle one level of ghost points for graphics
17347c6ae99SBarry Smith       */
1741321219cSEthan Coon       ierr = DMDACreate2d(comm,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,zctx.m,zctx.n,w,1,lx,ly,&dac);CHKERRQ(ierr);
175aa219208SBarry Smith       ierr = PetscInfo(da,"Creating auxilary DMDA for managing graphics ghost points\n");CHKERRQ(ierr);
17647c6ae99SBarry Smith     } else {
17747c6ae99SBarry Smith       /* otherwise we can use the da we already have */
17847c6ae99SBarry Smith       dac = da;
17947c6ae99SBarry Smith     }
18047c6ae99SBarry Smith     /* create local vector for holding ghosted values used in graphics */
181564755cdSBarry Smith     ierr = DMCreateLocalVector(dac,&xlocal);CHKERRQ(ierr);
18247c6ae99SBarry Smith     if (dac != da) {
183aa219208SBarry Smith       /* don't keep any public reference of this DMDA, is is only available through xlocal */
184f7923d8aSBarry Smith       ierr = PetscObjectDereference((PetscObject)dac);CHKERRQ(ierr);
18547c6ae99SBarry Smith     } else {
18647c6ae99SBarry Smith       /* remove association between xlocal and da, because below we compose in the opposite
18747c6ae99SBarry Smith          direction and if we left this connect we'd get a loop, so the objects could
18847c6ae99SBarry Smith          never be destroyed */
189c688c046SMatthew G Knepley       ierr = PetscObjectRemoveReference((PetscObject)xlocal,"__PETSc_dm");CHKERRQ(ierr);
19047c6ae99SBarry Smith     }
19147c6ae99SBarry Smith     ierr = PetscObjectCompose((PetscObject)da,"GraphicsGhosted",(PetscObject)xlocal);CHKERRQ(ierr);
19247c6ae99SBarry Smith     ierr = PetscObjectDereference((PetscObject)xlocal);CHKERRQ(ierr);
19347c6ae99SBarry Smith   } else {
194f7923d8aSBarry Smith     if (bx !=  DMDA_BOUNDARY_NONE || by !=  DMDA_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) {
195c688c046SMatthew G Knepley       ierr = VecGetDM(xlocal, &dac);CHKERRQ(ierr);
196f7923d8aSBarry Smith     } else {
197f7923d8aSBarry Smith       dac = da;
19847c6ae99SBarry Smith     }
19947c6ae99SBarry Smith   }
20047c6ae99SBarry Smith 
20147c6ae99SBarry Smith   /*
20247c6ae99SBarry Smith       Get local (ghosted) values of vector
20347c6ae99SBarry Smith   */
2049a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr);
2059a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr);
20647c6ae99SBarry Smith   ierr = VecGetArray(xlocal,&zctx.v);CHKERRQ(ierr);
20747c6ae99SBarry Smith 
20847c6ae99SBarry Smith   /* get coordinates of nodes */
2096636e97aSMatthew G Knepley   ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
21047c6ae99SBarry Smith   if (!xcoor) {
211aa219208SBarry Smith     ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,0.0);CHKERRQ(ierr);
2126636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
21347c6ae99SBarry Smith   }
21447c6ae99SBarry Smith 
21547c6ae99SBarry Smith   /*
21647c6ae99SBarry Smith       Determine the min and max  coordinates in plot
21747c6ae99SBarry Smith   */
2180298fd71SBarry Smith   ierr     = VecStrideMin(xcoor,0,NULL,&xmin);CHKERRQ(ierr);
2190298fd71SBarry Smith   ierr     = VecStrideMax(xcoor,0,NULL,&xmax);CHKERRQ(ierr);
2200298fd71SBarry Smith   ierr     = VecStrideMin(xcoor,1,NULL,&ymin);CHKERRQ(ierr);
2210298fd71SBarry Smith   ierr     = VecStrideMax(xcoor,1,NULL,&ymax);CHKERRQ(ierr);
22247c6ae99SBarry Smith   coors[0] = xmin - .05*(xmax- xmin); coors[2] = xmax + .05*(xmax - xmin);
22347c6ae99SBarry Smith   coors[1] = ymin - .05*(ymax- ymin); coors[3] = ymax + .05*(ymax - ymin);
224aa219208SBarry Smith   ierr     = PetscInfo4(da,"Preparing DMDA 2d contour plot coordinates %G %G %G %G\n",coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr);
22547c6ae99SBarry Smith 
226a74ba6f7SBarry Smith   ierr = PetscOptionsGetRealArray(NULL,"-draw_coordinates",coors,&ncoors,NULL);CHKERRQ(ierr);
227a74ba6f7SBarry Smith 
22847c6ae99SBarry Smith   /*
22947c6ae99SBarry Smith        get local ghosted version of coordinates
23047c6ae99SBarry Smith   */
23147c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr);
23247c6ae99SBarry Smith   if (!xcoorl) {
233aa219208SBarry Smith     /* create DMDA to get local version of graphics */
2341321219cSEthan Coon     ierr = DMDACreate2d(comm,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,zctx.m,zctx.n,2,1,lx,ly,&dag);CHKERRQ(ierr);
235aa219208SBarry Smith     ierr = PetscInfo(dag,"Creating auxilary DMDA for managing graphics coordinates ghost points\n");CHKERRQ(ierr);
236564755cdSBarry Smith     ierr = DMCreateLocalVector(dag,&xcoorl);CHKERRQ(ierr);
23747c6ae99SBarry Smith     ierr = PetscObjectCompose((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject)xcoorl);CHKERRQ(ierr);
238f7923d8aSBarry Smith     ierr = PetscObjectDereference((PetscObject)dag);CHKERRQ(ierr);
23947c6ae99SBarry Smith     ierr = PetscObjectDereference((PetscObject)xcoorl);CHKERRQ(ierr);
24047c6ae99SBarry Smith   } else {
241c688c046SMatthew G Knepley     ierr = VecGetDM(xcoorl,&dag);CHKERRQ(ierr);
24247c6ae99SBarry Smith   }
2439a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
2449a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
24547c6ae99SBarry Smith   ierr = VecGetArray(xcoorl,&zctx.xy);CHKERRQ(ierr);
24647c6ae99SBarry Smith 
24747c6ae99SBarry Smith   /*
24847c6ae99SBarry Smith         Get information about size of area each processor must do graphics for
24947c6ae99SBarry Smith   */
2501321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&M,&N,0,0,0,0,&zctx.step,0,&bx,&by,0,0);CHKERRQ(ierr);
251f7923d8aSBarry Smith   ierr = DMDAGetGhostCorners(dac,0,0,0,&zctx.m,&zctx.n,0);CHKERRQ(ierr);
25247c6ae99SBarry Smith 
2530298fd71SBarry Smith   ierr = PetscOptionsGetBool(NULL,"-draw_contour_grid",&zctx.showgrid,NULL);CHKERRQ(ierr);
2544e6118eeSBarry Smith 
2554e6118eeSBarry Smith   ierr = DMDASelectFields(da,&ndisplayfields,&displayfields);CHKERRQ(ierr);
25647c6ae99SBarry Smith 
25747c6ae99SBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
2580298fd71SBarry Smith   ierr = PetscOptionsGetBool(NULL,"-draw_ports",&useports,NULL);CHKERRQ(ierr);
25947c6ae99SBarry Smith   if (useports || format == PETSC_VIEWER_DRAW_PORTS) {
26047c6ae99SBarry Smith     ierr       = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
26120d0051dSBarry Smith     ierr       = PetscDrawViewPortsCreate(draw,ndisplayfields,&ports);CHKERRQ(ierr);
262109c9344SBarry Smith     zctx.name0 = 0;
263109c9344SBarry Smith     zctx.name1 = 0;
264109c9344SBarry Smith   } else {
265109c9344SBarry Smith     ierr = DMDAGetCoordinateName(da,0,&zctx.name0);CHKERRQ(ierr);
266109c9344SBarry Smith     ierr = DMDAGetCoordinateName(da,1,&zctx.name1);CHKERRQ(ierr);
26747c6ae99SBarry Smith   }
26820d0051dSBarry Smith 
26947c6ae99SBarry Smith   /*
27047c6ae99SBarry Smith      Loop over each field; drawing each in a different window
27147c6ae99SBarry Smith   */
27220d0051dSBarry Smith   for (i=0; i<ndisplayfields; i++) {
27320d0051dSBarry Smith     zctx.k = displayfields[i];
27447c6ae99SBarry Smith     if (useports) {
27520d0051dSBarry Smith       ierr = PetscDrawViewPortsSet(ports,i);CHKERRQ(ierr);
2769332fd86SBarry Smith     } else {
2779332fd86SBarry Smith       ierr = PetscViewerDrawGetDraw(viewer,i,&draw);CHKERRQ(ierr);
27847c6ae99SBarry Smith     }
27947c6ae99SBarry Smith 
28047c6ae99SBarry Smith     /*
28147c6ae99SBarry Smith         Determine the min and max color in plot
28247c6ae99SBarry Smith     */
2830298fd71SBarry Smith     ierr = VecStrideMin(xin,zctx.k,NULL,&zctx.min);CHKERRQ(ierr);
2840298fd71SBarry Smith     ierr = VecStrideMax(xin,zctx.k,NULL,&zctx.max);CHKERRQ(ierr);
28567dd0837SBarry Smith     if (zctx.k < nbounds) {
286f3f0eb19SBarry Smith       zctx.min = bounds[2*zctx.k];
287f3f0eb19SBarry Smith       zctx.max = bounds[2*zctx.k+1];
28867dd0837SBarry Smith     }
28947c6ae99SBarry Smith     if (zctx.min == zctx.max) {
29047c6ae99SBarry Smith       zctx.min -= 1.e-12;
29147c6ae99SBarry Smith       zctx.max += 1.e-12;
29247c6ae99SBarry Smith     }
29347c6ae99SBarry Smith 
29447c6ae99SBarry Smith     if (!rank) {
29547c6ae99SBarry Smith       const char *title;
29647c6ae99SBarry Smith 
297aa219208SBarry Smith       ierr = DMDAGetFieldName(da,zctx.k,&title);CHKERRQ(ierr);
29847c6ae99SBarry Smith       if (title) {
29947c6ae99SBarry Smith         ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr);
30047c6ae99SBarry Smith       }
30147c6ae99SBarry Smith     }
30247c6ae99SBarry Smith     ierr = PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr);
303aa219208SBarry Smith     ierr = PetscInfo2(da,"DMDA 2d contour plot min %G max %G\n",zctx.min,zctx.max);CHKERRQ(ierr);
30447c6ae99SBarry Smith 
30547c6ae99SBarry Smith     ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
30647c6ae99SBarry Smith     if (popup) {ierr = PetscDrawScalePopup(popup,zctx.min,zctx.max);CHKERRQ(ierr);}
30747c6ae99SBarry Smith 
30847c6ae99SBarry Smith     zctx.scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/(zctx.max - zctx.min);
30947c6ae99SBarry Smith 
31047c6ae99SBarry Smith     ierr = PetscDrawZoom(draw,VecView_MPI_Draw_DA2d_Zoom,&zctx);CHKERRQ(ierr);
31147c6ae99SBarry Smith   }
31220d0051dSBarry Smith   ierr = PetscFree(displayfields);CHKERRQ(ierr);
3136bf464f9SBarry Smith   ierr = PetscDrawViewPortsDestroy(ports);CHKERRQ(ierr);
31447c6ae99SBarry Smith 
31547c6ae99SBarry Smith   ierr = VecRestoreArray(xcoorl,&zctx.xy);CHKERRQ(ierr);
31647c6ae99SBarry Smith   ierr = VecRestoreArray(xlocal,&zctx.v);CHKERRQ(ierr);
31747c6ae99SBarry Smith   PetscFunctionReturn(0);
31847c6ae99SBarry Smith }
31947c6ae99SBarry Smith 
3200da24e51SJuha Jäykkä #if defined(PETSC_HAVE_HDF5)
3210da24e51SJuha Jäykkä #undef __FUNCT__
3220da24e51SJuha Jäykkä #define __FUNCT__ "VecGetHDF5ChunkSize"
323d4190030SJed Brown static PetscErrorCode VecGetHDF5ChunkSize(DM_DA *da, Vec xin, PetscInt timestep, hsize_t *chunkDims)
3240da24e51SJuha Jäykkä {
325d4190030SJed Brown   PetscMPIInt comm_size;
326d4190030SJed Brown   PetscErrorCode ierr;
3270da24e51SJuha Jäykkä   hsize_t chunk_size, target_size, dim,
3280da24e51SJuha Jäykkä     vec_size = sizeof(PetscScalar)*da->M*da->N*da->P*da->w,
3290da24e51SJuha Jäykkä     avg_local_vec_size,
3300da24e51SJuha Jäykkä     KiB = 1024,
3310da24e51SJuha Jäykkä     MiB = KiB*KiB,
3320da24e51SJuha Jäykkä     GiB = MiB*KiB,
3330da24e51SJuha Jäykkä     min_size = MiB,
3340da24e51SJuha Jäykkä     max_chunks = 64*KiB,                                              /* HDF5 internal limitation */
335*84179ffaSJed Brown     max_chunk_size = 4*GiB;                                           /* HDF5 internal limitation */
336*84179ffaSJed Brown   PetscInt zslices=da->p, yslices=da->n, xslices=da->m;
3370da24e51SJuha Jäykkä 
338cb5c4f94SJuha Jäykkä   PetscFunctionBegin;
339d4190030SJed Brown   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)xin), &comm_size);CHKERRQ(ierr);
3400da24e51SJuha 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 */
3410da24e51SJuha Jäykkä 
3420da24e51SJuha Jäykkä   target_size = (hsize_t) ceil(PetscMin(vec_size,
3430da24e51SJuha Jäykkä                                         PetscMin(max_chunk_size,
3440da24e51SJuha Jäykkä                                                  PetscMax(avg_local_vec_size,
3450da24e51SJuha Jäykkä                                                           PetscMax(ceil(vec_size*1.0/max_chunks),
3460da24e51SJuha Jäykkä                                                                    min_size)
3470da24e51SJuha Jäykkä                                                           )
3480da24e51SJuha Jäykkä                                                  )
3490da24e51SJuha Jäykkä                                         )
3500da24e51SJuha Jäykkä                                );
3510da24e51SJuha 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(PetscScalar);
3520da24e51SJuha Jäykkä 
353cb5c4f94SJuha Jäykkä   /*
354cb5c4f94SJuha Jäykkä    if size/rank > max_chunk_size, we need radical measures: even going down to
355cb5c4f94SJuha Jäykkä    avg_local_vec_size is not enough, so we simply use chunk size of 4 GiB no matter
356cb5c4f94SJuha Jäykkä    what, composed in the most efficient way possible.
357cb5c4f94SJuha Jäykkä    N.B. this minimises the number of chunks, which may or may not be the optimal
358cb5c4f94SJuha Jäykkä    solution. In a BG, for example, the optimal solution is probably to make # chunks = #
359cb5c4f94SJuha Jäykkä    IO nodes involved, but this author has no access to a BG to figure out how to
360cb5c4f94SJuha Jäykkä    reliably find the right number. And even then it may or may not be enough.
361cb5c4f94SJuha Jäykkä    */
3620da24e51SJuha Jäykkä   if (avg_local_vec_size > max_chunk_size) {
363cb5c4f94SJuha Jäykkä     /* check if we can just split local z-axis: is that enough? */
364*84179ffaSJed Brown     zslices = (PetscInt)ceil(vec_size*1.0/(da->p*max_chunk_size))*zslices;
3650da24e51SJuha Jäykkä     if (zslices > da->P) {
366cb5c4f94SJuha Jäykkä       /* lattice is too large in xy-directions, splitting z only is not enough */
3670da24e51SJuha Jäykkä       zslices = da->P;
368*84179ffaSJed Brown       yslices= (PetscInt)ceil(vec_size*1.0/(zslices*da->n*max_chunk_size))*yslices;
3690da24e51SJuha Jäykkä       if (yslices > da->N) {
370cb5c4f94SJuha Jäykkä 	/* lattice is too large in x-direction, splitting along z, y is not enough */
3710da24e51SJuha Jäykkä 	yslices = da->N;
372*84179ffaSJed Brown 	xslices= (PetscInt)ceil(vec_size*1.0/(zslices*yslices*da->m*max_chunk_size))*xslices;
3730da24e51SJuha Jäykkä       }
3740da24e51SJuha Jäykkä     }
3750da24e51SJuha Jäykkä     dim = 0;
3760da24e51SJuha Jäykkä     if (timestep >= 0) {
3770da24e51SJuha Jäykkä       ++dim;
3780da24e51SJuha Jäykkä     }
379cb5c4f94SJuha Jäykkä     /* prefer to split z-axis, even down to planar slices */
3800da24e51SJuha Jäykkä     if (da->dim == 3) {
381cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->P/zslices;
382cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->N/yslices;
383cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->M/xslices;
3840da24e51SJuha Jäykkä     } else {
385cb5c4f94SJuha Jäykkä       /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */
386cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->N/yslices;
387cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->M/xslices;
3880da24e51SJuha Jäykkä     }
3890da24e51SJuha 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);
3900da24e51SJuha Jäykkä   } else {
3910da24e51SJuha Jäykkä     if (target_size < chunk_size) {
392cb5c4f94SJuha Jäykkä       /* only change the defaults if target_size < chunk_size */
3930da24e51SJuha Jäykkä       dim = 0;
3940da24e51SJuha Jäykkä       if (timestep >= 0) {
3950da24e51SJuha Jäykkä 	++dim;
3960da24e51SJuha Jäykkä       }
397cb5c4f94SJuha Jäykkä       /* prefer to split z-axis, even down to planar slices */
3980da24e51SJuha Jäykkä       if (da->dim == 3) {
399cb5c4f94SJuha Jäykkä 	/* try splitting the z-axis to core-size bits, i.e. divide chunk size by # comm_size in z-direction */
4000da24e51SJuha Jäykkä 	if (target_size >= chunk_size/da->p) {
401cb5c4f94SJuha Jäykkä 	  /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */
4020da24e51SJuha Jäykkä 	  chunkDims[dim] = (hsize_t) ceil(da->P*1.0/da->p);
4030da24e51SJuha Jäykkä 	} else {
404cb5c4f94SJuha Jäykkä 	  /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be
405cb5c4f94SJuha Jäykkä            radical and let everyone write all they've got */
4060da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->P*1.0/da->p);
4070da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n);
4080da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m);
4090da24e51SJuha Jäykkä 	}
4100da24e51SJuha Jäykkä       } else {
411cb5c4f94SJuha Jäykkä 	/* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */
4120da24e51SJuha Jäykkä 	if (target_size >= chunk_size/da->n) {
413cb5c4f94SJuha Jäykkä 	  /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */
4140da24e51SJuha Jäykkä 	  chunkDims[dim] = (hsize_t) ceil(da->N*1.0/da->n);
4150da24e51SJuha Jäykkä 	} else {
416cb5c4f94SJuha Jäykkä 	  /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be
417cb5c4f94SJuha Jäykkä 	   radical and let everyone write all they've got */
4180da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n);
4190da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m);
4200da24e51SJuha Jäykkä 	}
4210da24e51SJuha Jäykkä 
4220da24e51SJuha Jäykkä       }
4230da24e51SJuha 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);
4240da24e51SJuha Jäykkä     } else {
425cb5c4f94SJuha Jäykkä       /* precomputed chunks are fine, we don't need to do anything */
4260da24e51SJuha Jäykkä     }
4270da24e51SJuha Jäykkä   }
4280da24e51SJuha Jäykkä   PetscFunctionReturn(0);
4290da24e51SJuha Jäykkä }
4300da24e51SJuha Jäykkä #endif
43147c6ae99SBarry Smith 
43247c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
43347c6ae99SBarry Smith #undef __FUNCT__
43447c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_HDF5_DA"
43547c6ae99SBarry Smith PetscErrorCode VecView_MPI_HDF5_DA(Vec xin,PetscViewer viewer)
43647c6ae99SBarry Smith {
4379b2a5a86SJed Brown   DM             dm;
4389b2a5a86SJed Brown   DM_DA          *da;
43947c6ae99SBarry Smith   hid_t          filespace;  /* file dataspace identifier */
4408e2ae6d7SMichael Kraus   hid_t          chunkspace; /* chunk dataset property identifier */
44147c6ae99SBarry Smith   hid_t          plist_id;   /* property list identifier */
44247c6ae99SBarry Smith   hid_t          dset_id;    /* dataset identifier */
44347c6ae99SBarry Smith   hid_t          memspace;   /* memory dataspace identifier */
44447c6ae99SBarry Smith   hid_t          file_id;
44515214e8eSMatthew G Knepley   hid_t          group;
4468e2ae6d7SMichael Kraus   hid_t          scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
44747c6ae99SBarry Smith   herr_t         status;
448d9a4edebSJed Brown   hsize_t        dim;
4490da24e51SJuha 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  */
45015214e8eSMatthew G Knepley   PetscInt       timestep;
4518e2ae6d7SMichael Kraus   PetscScalar    *x;
4528e2ae6d7SMichael Kraus   const char     *vecname;
45315214e8eSMatthew G Knepley   PetscErrorCode ierr;
45447c6ae99SBarry Smith 
45547c6ae99SBarry Smith   PetscFunctionBegin;
45615214e8eSMatthew G Knepley   ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr);
4578e2ae6d7SMichael Kraus   ierr = PetscViewerHDF5GetTimestep(viewer, &timestep);CHKERRQ(ierr);
45815214e8eSMatthew G Knepley 
459c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&dm);CHKERRQ(ierr);
460ce94432eSBarry Smith   if (!dm) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
4619b2a5a86SJed Brown   da = (DM_DA*)dm->data;
46247c6ae99SBarry Smith 
4638e2ae6d7SMichael Kraus   /* Create the dataspace for the dataset.
4648e2ae6d7SMichael Kraus    *
4658e2ae6d7SMichael Kraus    * dims - holds the current dimensions of the dataset
4668e2ae6d7SMichael Kraus    *
4678e2ae6d7SMichael Kraus    * maxDims - holds the maximum dimensions of the dataset (unlimited
4688e2ae6d7SMichael Kraus    * for the number of time steps with the current dimensions for the
4698e2ae6d7SMichael Kraus    * other dimensions; so only additional time steps can be added).
4708e2ae6d7SMichael Kraus    *
4718e2ae6d7SMichael Kraus    * chunkDims - holds the size of a single time step (required to
4728e2ae6d7SMichael Kraus    * permit extending dataset).
4738e2ae6d7SMichael Kraus    */
4748e2ae6d7SMichael Kraus   dim = 0;
4758e2ae6d7SMichael Kraus   if (timestep >= 0) {
4768e2ae6d7SMichael Kraus     dims[dim]      = timestep+1;
4778e2ae6d7SMichael Kraus     maxDims[dim]   = H5S_UNLIMITED;
4788e2ae6d7SMichael Kraus     chunkDims[dim] = 1;
4798e2ae6d7SMichael Kraus     ++dim;
4808e2ae6d7SMichael Kraus   }
4818e2ae6d7SMichael Kraus   if (da->dim == 3) {
482acba2ac6SBarry Smith     ierr           = PetscHDF5IntCast(da->P,dims+dim);CHKERRQ(ierr);
4838e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4848e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4858e2ae6d7SMichael Kraus     ++dim;
4868e2ae6d7SMichael Kraus   }
4878e2ae6d7SMichael Kraus   if (da->dim > 1) {
488acba2ac6SBarry Smith     ierr           = PetscHDF5IntCast(da->N,dims+dim);CHKERRQ(ierr);
4898e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4908e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4918e2ae6d7SMichael Kraus     ++dim;
4928e2ae6d7SMichael Kraus   }
493acba2ac6SBarry Smith   ierr           = PetscHDF5IntCast(da->M,dims+dim);CHKERRQ(ierr);
4948e2ae6d7SMichael Kraus   maxDims[dim]   = dims[dim];
4958e2ae6d7SMichael Kraus   chunkDims[dim] = dims[dim];
4968e2ae6d7SMichael Kraus   ++dim;
4978e2ae6d7SMichael Kraus   if (da->w > 1) {
498acba2ac6SBarry Smith     ierr           = PetscHDF5IntCast(da->w,dims+dim);CHKERRQ(ierr);
4998e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
5008e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
5018e2ae6d7SMichael Kraus     ++dim;
5028e2ae6d7SMichael Kraus   }
50347c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5048e2ae6d7SMichael Kraus   dims[dim]      = 2;
5058e2ae6d7SMichael Kraus   maxDims[dim]   = dims[dim];
5068e2ae6d7SMichael Kraus   chunkDims[dim] = dims[dim];
5078e2ae6d7SMichael Kraus   ++dim;
50847c6ae99SBarry Smith #endif
5090da24e51SJuha Jäykkä 
510cb5c4f94SJuha Jäykkä   ierr = VecGetHDF5ChunkSize(da, xin, timestep, chunkDims); CHKERRQ(ierr);
5110da24e51SJuha Jäykkä 
5128e2ae6d7SMichael Kraus   filespace = H5Screate_simple(dim, dims, maxDims);
51347c6ae99SBarry Smith   if (filespace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Screate_simple()");
51447c6ae99SBarry Smith 
51515214e8eSMatthew G Knepley #if defined(PETSC_USE_REAL_SINGLE)
51615214e8eSMatthew G Knepley   scalartype = H5T_NATIVE_FLOAT;
51715214e8eSMatthew G Knepley #elif defined(PETSC_USE_REAL___FLOAT128)
51815214e8eSMatthew G Knepley #error "HDF5 output with 128 bit floats not supported."
51915214e8eSMatthew G Knepley #else
52015214e8eSMatthew G Knepley   scalartype = H5T_NATIVE_DOUBLE;
52115214e8eSMatthew G Knepley #endif
52215214e8eSMatthew G Knepley 
52347c6ae99SBarry Smith   /* Create the dataset with default properties and close filespace */
52447c6ae99SBarry Smith   ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr);
52515214e8eSMatthew G Knepley   if (!H5Lexists(group, vecname, H5P_DEFAULT)) {
5268e2ae6d7SMichael Kraus     /* Create chunk */
5278e2ae6d7SMichael Kraus     chunkspace = H5Pcreate(H5P_DATASET_CREATE);
5288e2ae6d7SMichael Kraus     if (chunkspace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Pcreate()");
5298e2ae6d7SMichael Kraus     status = H5Pset_chunk(chunkspace, dim, chunkDims);CHKERRQ(status);
5308e2ae6d7SMichael Kraus 
53147c6ae99SBarry Smith #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
5328e2ae6d7SMichael Kraus     dset_id = H5Dcreate2(group, vecname, scalartype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT);
53347c6ae99SBarry Smith #else
53415214e8eSMatthew G Knepley     dset_id = H5Dcreate(group, vecname, scalartype, filespace, H5P_DEFAULT);
53547c6ae99SBarry Smith #endif
53647c6ae99SBarry Smith     if (dset_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Dcreate2()");
53715214e8eSMatthew G Knepley   } else {
53815214e8eSMatthew G Knepley     dset_id = H5Dopen2(group, vecname, H5P_DEFAULT);
53915214e8eSMatthew G Knepley     status  = H5Dset_extent(dset_id, dims);CHKERRQ(status);
54015214e8eSMatthew G Knepley   }
54147c6ae99SBarry Smith   status = H5Sclose(filespace);CHKERRQ(status);
54247c6ae99SBarry Smith 
54347c6ae99SBarry Smith   /* Each process defines a dataset and writes it to the hyperslab in the file */
5448e2ae6d7SMichael Kraus   dim = 0;
5458e2ae6d7SMichael Kraus   if (timestep >= 0) {
5468e2ae6d7SMichael Kraus     offset[dim] = timestep;
5478e2ae6d7SMichael Kraus     ++dim;
5488e2ae6d7SMichael Kraus   }
549acba2ac6SBarry Smith   if (da->dim == 3) {ierr = PetscHDF5IntCast(da->zs,offset + dim++);CHKERRQ(ierr);}
550acba2ac6SBarry Smith   if (da->dim > 1)  {ierr = PetscHDF5IntCast(da->ys,offset + dim++);CHKERRQ(ierr);}
551acba2ac6SBarry Smith   ierr = PetscHDF5IntCast(da->xs/da->w,offset + dim++);CHKERRQ(ierr);
5528e2ae6d7SMichael Kraus   if (da->w > 1) offset[dim++] = 0;
55347c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5548e2ae6d7SMichael Kraus   offset[dim++] = 0;
55547c6ae99SBarry Smith #endif
5568e2ae6d7SMichael Kraus   dim = 0;
5578e2ae6d7SMichael Kraus   if (timestep >= 0) {
5588e2ae6d7SMichael Kraus     count[dim] = 1;
5598e2ae6d7SMichael Kraus     ++dim;
5608e2ae6d7SMichael Kraus   }
561acba2ac6SBarry Smith   if (da->dim == 3) {ierr = PetscHDF5IntCast(da->ze - da->zs,count + dim++);CHKERRQ(ierr);}
562acba2ac6SBarry Smith   if (da->dim > 1)  {ierr = PetscHDF5IntCast(da->ye - da->ys,count + dim++);CHKERRQ(ierr);}
563acba2ac6SBarry Smith   ierr = PetscHDF5IntCast((da->xe - da->xs)/da->w,count + dim++);CHKERRQ(ierr);
564acba2ac6SBarry Smith   if (da->w > 1) {ierr = PetscHDF5IntCast(da->w,count + dim++);CHKERRQ(ierr);}
56547c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5668e2ae6d7SMichael Kraus   count[dim++] = 2;
56747c6ae99SBarry Smith #endif
56847c6ae99SBarry Smith   memspace = H5Screate_simple(dim, count, NULL);
56947c6ae99SBarry Smith   if (memspace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Screate_simple()");
57047c6ae99SBarry Smith 
57147c6ae99SBarry Smith   filespace = H5Dget_space(dset_id);
57247c6ae99SBarry Smith   if (filespace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Dget_space()");
57347c6ae99SBarry Smith   status = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);CHKERRQ(status);
57447c6ae99SBarry Smith 
57547c6ae99SBarry Smith   /* Create property list for collective dataset write */
57647c6ae99SBarry Smith   plist_id = H5Pcreate(H5P_DATASET_XFER);
57747c6ae99SBarry Smith   if (plist_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Pcreate()");
57847c6ae99SBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
57947c6ae99SBarry Smith   status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);CHKERRQ(status);
58047c6ae99SBarry Smith #endif
58147c6ae99SBarry Smith   /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
58247c6ae99SBarry Smith 
58347c6ae99SBarry Smith   ierr   = VecGetArray(xin, &x);CHKERRQ(ierr);
58415214e8eSMatthew G Knepley   status = H5Dwrite(dset_id, scalartype, memspace, filespace, plist_id, x);CHKERRQ(status);
58547c6ae99SBarry Smith   status = H5Fflush(file_id, H5F_SCOPE_GLOBAL);CHKERRQ(status);
58647c6ae99SBarry Smith   ierr   = VecRestoreArray(xin, &x);CHKERRQ(ierr);
58747c6ae99SBarry Smith 
58847c6ae99SBarry Smith   /* Close/release resources */
58915214e8eSMatthew G Knepley   if (group != file_id) {
59015214e8eSMatthew G Knepley     status = H5Gclose(group);CHKERRQ(status);
59115214e8eSMatthew G Knepley   }
59247c6ae99SBarry Smith   status = H5Pclose(plist_id);CHKERRQ(status);
59347c6ae99SBarry Smith   status = H5Sclose(filespace);CHKERRQ(status);
59447c6ae99SBarry Smith   status = H5Sclose(memspace);CHKERRQ(status);
59547c6ae99SBarry Smith   status = H5Dclose(dset_id);CHKERRQ(status);
59647c6ae99SBarry Smith   ierr   = PetscInfo1(xin,"Wrote Vec object with name %s\n",vecname);CHKERRQ(ierr);
59747c6ae99SBarry Smith   PetscFunctionReturn(0);
59847c6ae99SBarry Smith }
59947c6ae99SBarry Smith #endif
60047c6ae99SBarry Smith 
60109573ac7SBarry Smith extern PetscErrorCode VecView_MPI_Draw_DA1d(Vec,PetscViewer);
60247c6ae99SBarry Smith 
60347c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
60447c6ae99SBarry Smith #undef __FUNCT__
605aa219208SBarry Smith #define __FUNCT__ "DMDAArrayMPIIO"
606aa219208SBarry Smith static PetscErrorCode DMDAArrayMPIIO(DM da,PetscViewer viewer,Vec xin,PetscBool write)
60747c6ae99SBarry Smith {
60847c6ae99SBarry Smith   PetscErrorCode    ierr;
60947c6ae99SBarry Smith   MPI_File          mfdes;
61047c6ae99SBarry Smith   PetscMPIInt       gsizes[4],lsizes[4],lstarts[4],asiz,dof;
61147c6ae99SBarry Smith   MPI_Datatype      view;
61247c6ae99SBarry Smith   const PetscScalar *array;
61347c6ae99SBarry Smith   MPI_Offset        off;
61447c6ae99SBarry Smith   MPI_Aint          ub,ul;
61547c6ae99SBarry Smith   PetscInt          type,rows,vecrows,tr[2];
61647c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
61747c6ae99SBarry Smith 
61847c6ae99SBarry Smith   PetscFunctionBegin;
61947c6ae99SBarry Smith   ierr = VecGetSize(xin,&vecrows);CHKERRQ(ierr);
62047c6ae99SBarry Smith   if (!write) {
62147c6ae99SBarry Smith     /* Read vector header. */
62247c6ae99SBarry Smith     ierr = PetscViewerBinaryRead(viewer,tr,2,PETSC_INT);CHKERRQ(ierr);
62347c6ae99SBarry Smith     type = tr[0];
62447c6ae99SBarry Smith     rows = tr[1];
625ce94432eSBarry Smith     if (type != VEC_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Not vector next in file");
626ce94432eSBarry Smith     if (rows != vecrows) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Vector in file not same size as DMDA vector");
62747c6ae99SBarry Smith   } else {
62847c6ae99SBarry Smith     tr[0] = VEC_FILE_CLASSID;
62947c6ae99SBarry Smith     tr[1] = vecrows;
63047c6ae99SBarry Smith     ierr  = PetscViewerBinaryWrite(viewer,tr,2,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
63147c6ae99SBarry Smith   }
63247c6ae99SBarry Smith 
6334dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->w,&dof);CHKERRQ(ierr);
6344dc2109aSBarry Smith   gsizes[0]  = dof;
6354dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->M,gsizes+1);CHKERRQ(ierr);
6364dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->N,gsizes+2);CHKERRQ(ierr);
6374dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->P,gsizes+1);CHKERRQ(ierr);
6384dc2109aSBarry Smith   lsizes[0]  = dof;
6394dc2109aSBarry Smith   ierr       = PetscMPIIntCast((dd->xe-dd->xs)/dof,lsizes+1);CHKERRQ(ierr);
6404dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->ye-dd->ys,lsizes+2);CHKERRQ(ierr);
6414dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->ze-dd->zs,lsizes+3);CHKERRQ(ierr);
6424dc2109aSBarry Smith   lstarts[0] = 0;
6434dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->xs/dof,lstarts+1);CHKERRQ(ierr);
6444dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->ys,lstarts+2);CHKERRQ(ierr);
6454dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->zs,lstarts+3);CHKERRQ(ierr);
64647c6ae99SBarry Smith   ierr       = MPI_Type_create_subarray(dd->dim+1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view);CHKERRQ(ierr);
64747c6ae99SBarry Smith   ierr       = MPI_Type_commit(&view);CHKERRQ(ierr);
64847c6ae99SBarry Smith 
64947c6ae99SBarry Smith   ierr = PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);CHKERRQ(ierr);
65047c6ae99SBarry Smith   ierr = PetscViewerBinaryGetMPIIOOffset(viewer,&off);CHKERRQ(ierr);
65147c6ae99SBarry Smith   ierr = MPI_File_set_view(mfdes,off,MPIU_SCALAR,view,(char*)"native",MPI_INFO_NULL);CHKERRQ(ierr);
65247c6ae99SBarry Smith   ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr);
65347c6ae99SBarry Smith   asiz = lsizes[1]*(lsizes[2] > 0 ? lsizes[2] : 1)*(lsizes[3] > 0 ? lsizes[3] : 1)*dof;
65447c6ae99SBarry Smith   if (write) {
65547c6ae99SBarry Smith     ierr = MPIU_File_write_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr);
65647c6ae99SBarry Smith   } else {
65747c6ae99SBarry Smith     ierr = MPIU_File_read_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr);
65847c6ae99SBarry Smith   }
65947c6ae99SBarry Smith   ierr = MPI_Type_get_extent(view,&ul,&ub);CHKERRQ(ierr);
66047c6ae99SBarry Smith   ierr = PetscViewerBinaryAddMPIIOOffset(viewer,ub);CHKERRQ(ierr);
66147c6ae99SBarry Smith   ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr);
66247c6ae99SBarry Smith   ierr = MPI_Type_free(&view);CHKERRQ(ierr);
66347c6ae99SBarry Smith   PetscFunctionReturn(0);
66447c6ae99SBarry Smith }
66547c6ae99SBarry Smith #endif
66647c6ae99SBarry Smith 
66747c6ae99SBarry Smith #undef __FUNCT__
66847c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_DA"
6697087cfbeSBarry Smith PetscErrorCode  VecView_MPI_DA(Vec xin,PetscViewer viewer)
67047c6ae99SBarry Smith {
6719a42bb27SBarry Smith   DM                da;
67247c6ae99SBarry Smith   PetscErrorCode    ierr;
67347c6ae99SBarry Smith   PetscInt          dim;
67447c6ae99SBarry Smith   Vec               natural;
6754061b8bfSJed Brown   PetscBool         isdraw,isvtk;
67647c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
67747c6ae99SBarry Smith   PetscBool         ishdf5;
67847c6ae99SBarry Smith #endif
6793f3fd955SJed Brown   const char        *prefix,*name;
680a261c58fSBarry Smith   PetscViewerFormat format;
68147c6ae99SBarry Smith 
68247c6ae99SBarry Smith   PetscFunctionBegin;
683c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
684ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
685251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
686251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr);
68747c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
688251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr);
68947c6ae99SBarry Smith #endif
69047c6ae99SBarry Smith   if (isdraw) {
6911321219cSEthan Coon     ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
69247c6ae99SBarry Smith     if (dim == 1) {
69347c6ae99SBarry Smith       ierr = VecView_MPI_Draw_DA1d(xin,viewer);CHKERRQ(ierr);
69447c6ae99SBarry Smith     } else if (dim == 2) {
69547c6ae99SBarry Smith       ierr = VecView_MPI_Draw_DA2d(xin,viewer);CHKERRQ(ierr);
696ce94432eSBarry Smith     } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot graphically view vector associated with this dimensional DMDA %D",dim);
6974061b8bfSJed Brown   } else if (isvtk) {           /* Duplicate the Vec and hold a reference to the DM */
6984061b8bfSJed Brown     Vec Y;
6994061b8bfSJed Brown     ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr);
7004061b8bfSJed Brown     ierr = VecDuplicate(xin,&Y);CHKERRQ(ierr);
701b51b94faSJed Brown     if (((PetscObject)xin)->name) {
702b51b94faSJed Brown       /* If xin was named, copy the name over to Y. The duplicate names are safe because nobody else will ever see Y. */
703b51b94faSJed Brown       ierr = PetscObjectSetName((PetscObject)Y,((PetscObject)xin)->name);CHKERRQ(ierr);
704b51b94faSJed Brown     }
7054061b8bfSJed Brown     ierr = VecCopy(xin,Y);CHKERRQ(ierr);
70662b69a3fSMatthew G Knepley     ierr = PetscViewerVTKAddField(viewer,(PetscObject)da,DMDAVTKWriteAll,PETSC_VTK_POINT_FIELD,(PetscObject)Y);CHKERRQ(ierr);
70747c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
70847c6ae99SBarry Smith   } else if (ishdf5) {
70947c6ae99SBarry Smith     ierr = VecView_MPI_HDF5_DA(xin,viewer);CHKERRQ(ierr);
71047c6ae99SBarry Smith #endif
71147c6ae99SBarry Smith   } else {
71247c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
71347c6ae99SBarry Smith     PetscBool isbinary,isMPIIO;
71447c6ae99SBarry Smith 
715251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
71647c6ae99SBarry Smith     if (isbinary) {
71747c6ae99SBarry Smith       ierr = PetscViewerBinaryGetMPIIO(viewer,&isMPIIO);CHKERRQ(ierr);
71847c6ae99SBarry Smith       if (isMPIIO) {
719aa219208SBarry Smith         ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_TRUE);CHKERRQ(ierr);
72047c6ae99SBarry Smith         PetscFunctionReturn(0);
72147c6ae99SBarry Smith       }
72247c6ae99SBarry Smith     }
72347c6ae99SBarry Smith #endif
72447c6ae99SBarry Smith 
72547c6ae99SBarry Smith     /* call viewer on natural ordering */
72647c6ae99SBarry Smith     ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr);
727aa219208SBarry Smith     ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr);
72847c6ae99SBarry Smith     ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr);
729aa219208SBarry Smith     ierr = DMDAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr);
730aa219208SBarry Smith     ierr = DMDAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr);
7313f3fd955SJed Brown     ierr = PetscObjectGetName((PetscObject)xin,&name);CHKERRQ(ierr);
7323f3fd955SJed Brown     ierr = PetscObjectSetName((PetscObject)natural,name);CHKERRQ(ierr);
733a261c58fSBarry Smith 
734a261c58fSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
735a261c58fSBarry Smith     if (format == PETSC_VIEWER_BINARY_MATLAB) {
736a261c58fSBarry Smith       /* temporarily remove viewer format so it won't trigger in the VecView() */
737a261c58fSBarry Smith       ierr = PetscViewerSetFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
738a261c58fSBarry Smith     }
739a261c58fSBarry Smith 
74047c6ae99SBarry Smith     ierr = VecView(natural,viewer);CHKERRQ(ierr);
741a261c58fSBarry Smith 
742a261c58fSBarry Smith     if (format == PETSC_VIEWER_BINARY_MATLAB) {
743a261c58fSBarry Smith       MPI_Comm    comm;
744a261c58fSBarry Smith       FILE        *info;
745a261c58fSBarry Smith       const char  *fieldname;
746da88d4d4SJed Brown       char        fieldbuf[256];
747a261c58fSBarry Smith       PetscInt    dim,ni,nj,nk,pi,pj,pk,dof,n;
748a261c58fSBarry Smith 
749a261c58fSBarry Smith       /* set the viewer format back into the viewer */
750a261c58fSBarry Smith       ierr = PetscViewerSetFormat(viewer,format);CHKERRQ(ierr);
751a261c58fSBarry Smith       ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
752a261c58fSBarry Smith       ierr = PetscViewerBinaryGetInfoPointer(viewer,&info);CHKERRQ(ierr);
753a261c58fSBarry Smith       ierr = DMDAGetInfo(da,&dim,&ni,&nj,&nk,&pi,&pj,&pk,&dof,0,0,0,0,0);CHKERRQ(ierr);
754da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n");CHKERRQ(ierr);
755da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#$$ tmp = PetscBinaryRead(fd); \n");CHKERRQ(ierr);
756da88d4d4SJed Brown       if (dim == 1) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d);\n",dof,ni);CHKERRQ(ierr); }
757da88d4d4SJed Brown       if (dim == 2) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d);\n",dof,ni,nj);CHKERRQ(ierr); }
758da88d4d4SJed Brown       if (dim == 3) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d,%d);\n",dof,ni,nj,nk);CHKERRQ(ierr); }
759a261c58fSBarry Smith 
760a261c58fSBarry Smith       for (n=0; n<dof; n++) {
761a261c58fSBarry Smith         ierr = DMDAGetFieldName(da,n,&fieldname);CHKERRQ(ierr);
762da88d4d4SJed Brown         if (!fieldname || !fieldname[0]) {
763da88d4d4SJed Brown           ierr = PetscSNPrintf(fieldbuf,sizeof fieldbuf,"field%D",n);CHKERRQ(ierr);
764da88d4d4SJed Brown           fieldname = fieldbuf;
765a261c58fSBarry Smith         }
766da88d4d4SJed Brown         if (dim == 1) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:))';\n",name,fieldname,n+1);CHKERRQ(ierr); }
767da88d4d4SJed Brown         if (dim == 2) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:,:))';\n",name,fieldname,n+1);CHKERRQ(ierr); }
768da88d4d4SJed 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);}
769a261c58fSBarry Smith       }
770da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#$$ clear tmp; \n");CHKERRQ(ierr);
771da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n");CHKERRQ(ierr);
772a261c58fSBarry Smith     }
773a261c58fSBarry Smith 
774fcfd50ebSBarry Smith     ierr = VecDestroy(&natural);CHKERRQ(ierr);
77547c6ae99SBarry Smith   }
77647c6ae99SBarry Smith   PetscFunctionReturn(0);
77747c6ae99SBarry Smith }
77847c6ae99SBarry Smith 
77947c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
78047c6ae99SBarry Smith #undef __FUNCT__
78147c6ae99SBarry Smith #define __FUNCT__ "VecLoad_HDF5_DA"
78247c6ae99SBarry Smith PetscErrorCode VecLoad_HDF5_DA(Vec xin, PetscViewer viewer)
78347c6ae99SBarry Smith {
7849a42bb27SBarry Smith   DM             da;
78547c6ae99SBarry Smith   PetscErrorCode ierr;
78625578ef6SJed Brown   hsize_t        dim;
78747c6ae99SBarry Smith   hsize_t        count[5];
78847c6ae99SBarry Smith   hsize_t        offset[5];
78947c6ae99SBarry Smith   PetscInt       cnt = 0;
79047c6ae99SBarry Smith   PetscScalar    *x;
79147c6ae99SBarry Smith   const char     *vecname;
79247c6ae99SBarry Smith   hid_t          filespace; /* file dataspace identifier */
79347c6ae99SBarry Smith   hid_t          plist_id;  /* property list identifier */
79447c6ae99SBarry Smith   hid_t          dset_id;   /* dataset identifier */
79547c6ae99SBarry Smith   hid_t          memspace;  /* memory dataspace identifier */
79647c6ae99SBarry Smith   hid_t          file_id;
79747c6ae99SBarry Smith   herr_t         status;
7989c7c4993SBarry Smith   DM_DA          *dd;
79947c6ae99SBarry Smith 
80047c6ae99SBarry Smith   PetscFunctionBegin;
80147c6ae99SBarry Smith   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
802c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
8039c7c4993SBarry Smith   dd   = (DM_DA*)da->data;
80447c6ae99SBarry Smith 
80547c6ae99SBarry Smith   /* Create the dataspace for the dataset */
806acba2ac6SBarry Smith   ierr = PetscHDF5IntCast(dd->dim + ((dd->w == 1) ? 0 : 1),&dim);CHKERRQ(ierr);
80747c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
80847c6ae99SBarry Smith   dim++;
80947c6ae99SBarry Smith #endif
81047c6ae99SBarry Smith 
81147c6ae99SBarry Smith   /* Create the dataset with default properties and close filespace */
81247c6ae99SBarry Smith   ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr);
81347c6ae99SBarry Smith #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
81447c6ae99SBarry Smith   dset_id = H5Dopen2(file_id, vecname, H5P_DEFAULT);
81547c6ae99SBarry Smith #else
81647c6ae99SBarry Smith   dset_id = H5Dopen(file_id, vecname);
81747c6ae99SBarry Smith #endif
81847c6ae99SBarry Smith   if (dset_id == -1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Dopen2() with Vec named %s",vecname);
81947c6ae99SBarry Smith   filespace = H5Dget_space(dset_id);
82047c6ae99SBarry Smith 
82147c6ae99SBarry Smith   /* Each process defines a dataset and reads it from the hyperslab in the file */
82247c6ae99SBarry Smith   cnt = 0;
823acba2ac6SBarry Smith   if (dd->dim == 3) {ierr = PetscHDF5IntCast(dd->zs,offset + cnt++);CHKERRQ(ierr);}
824acba2ac6SBarry Smith   if (dd->dim > 1)  {ierr = PetscHDF5IntCast(dd->ys,offset + cnt++);CHKERRQ(ierr);}
825acba2ac6SBarry Smith   ierr = PetscHDF5IntCast(dd->xs/dd->w,offset + cnt++);CHKERRQ(ierr);
82647c6ae99SBarry Smith   if (dd->w > 1) offset[cnt++] = 0;
82747c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
82847c6ae99SBarry Smith   offset[cnt++] = 0;
82947c6ae99SBarry Smith #endif
83047c6ae99SBarry Smith   cnt = 0;
831acba2ac6SBarry Smith   if (dd->dim == 3) {ierr = PetscHDF5IntCast(dd->ze - dd->zs,count + cnt++);CHKERRQ(ierr);}
832acba2ac6SBarry Smith   if (dd->dim > 1)  {ierr = PetscHDF5IntCast(dd->ye - dd->ys,count + cnt++);CHKERRQ(ierr);}
833acba2ac6SBarry Smith   ierr = PetscHDF5IntCast((dd->xe - dd->xs)/dd->w,count + cnt++);CHKERRQ(ierr);
834acba2ac6SBarry Smith   if (dd->w > 1) {ierr = PetscHDF5IntCast(dd->w,count + cnt++);CHKERRQ(ierr);}
83547c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
83647c6ae99SBarry Smith   count[cnt++] = 2;
83747c6ae99SBarry Smith #endif
83847c6ae99SBarry Smith   memspace = H5Screate_simple(dim, count, NULL);
83947c6ae99SBarry Smith   if (memspace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Screate_simple()");
84047c6ae99SBarry Smith 
84147c6ae99SBarry Smith   status = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);CHKERRQ(status);
84247c6ae99SBarry Smith 
84347c6ae99SBarry Smith   /* Create property list for collective dataset write */
84447c6ae99SBarry Smith   plist_id = H5Pcreate(H5P_DATASET_XFER);
84547c6ae99SBarry Smith   if (plist_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Pcreate()");
84647c6ae99SBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
84747c6ae99SBarry Smith   status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);CHKERRQ(status);
84847c6ae99SBarry Smith #endif
84947c6ae99SBarry Smith   /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
85047c6ae99SBarry Smith 
85147c6ae99SBarry Smith   ierr   = VecGetArray(xin, &x);CHKERRQ(ierr);
85247c6ae99SBarry Smith   status = H5Dread(dset_id, H5T_NATIVE_DOUBLE, memspace, filespace, plist_id, x);CHKERRQ(status);
85347c6ae99SBarry Smith   ierr   = VecRestoreArray(xin, &x);CHKERRQ(ierr);
85447c6ae99SBarry Smith 
85547c6ae99SBarry Smith   /* Close/release resources */
85647c6ae99SBarry Smith   status = H5Pclose(plist_id);CHKERRQ(status);
85747c6ae99SBarry Smith   status = H5Sclose(filespace);CHKERRQ(status);
85847c6ae99SBarry Smith   status = H5Sclose(memspace);CHKERRQ(status);
85947c6ae99SBarry Smith   status = H5Dclose(dset_id);CHKERRQ(status);
86047c6ae99SBarry Smith   PetscFunctionReturn(0);
86147c6ae99SBarry Smith }
86247c6ae99SBarry Smith #endif
86347c6ae99SBarry Smith 
86447c6ae99SBarry Smith #undef __FUNCT__
86547c6ae99SBarry Smith #define __FUNCT__ "VecLoad_Binary_DA"
86647c6ae99SBarry Smith PetscErrorCode VecLoad_Binary_DA(Vec xin, PetscViewer viewer)
86747c6ae99SBarry Smith {
8689a42bb27SBarry Smith   DM             da;
86947c6ae99SBarry Smith   PetscErrorCode ierr;
87047c6ae99SBarry Smith   Vec            natural;
87147c6ae99SBarry Smith   const char     *prefix;
87247c6ae99SBarry Smith   PetscInt       bs;
87347c6ae99SBarry Smith   PetscBool      flag;
87447c6ae99SBarry Smith   DM_DA          *dd;
87547c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
87647c6ae99SBarry Smith   PetscBool isMPIIO;
87747c6ae99SBarry Smith #endif
87847c6ae99SBarry Smith 
87947c6ae99SBarry Smith   PetscFunctionBegin;
880c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
88147c6ae99SBarry Smith   dd   = (DM_DA*)da->data;
88247c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
88347c6ae99SBarry Smith   ierr = PetscViewerBinaryGetMPIIO(viewer,&isMPIIO);CHKERRQ(ierr);
88447c6ae99SBarry Smith   if (isMPIIO) {
885aa219208SBarry Smith     ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_FALSE);CHKERRQ(ierr);
88647c6ae99SBarry Smith     PetscFunctionReturn(0);
88747c6ae99SBarry Smith   }
88847c6ae99SBarry Smith #endif
88947c6ae99SBarry Smith 
89047c6ae99SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr);
891aa219208SBarry Smith   ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr);
89247c6ae99SBarry Smith   ierr = PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name);CHKERRQ(ierr);
89347c6ae99SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr);
8948aea9937SBarry Smith   ierr = VecLoad(natural,viewer);CHKERRQ(ierr);
895aa219208SBarry Smith   ierr = DMDANaturalToGlobalBegin(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr);
896aa219208SBarry Smith   ierr = DMDANaturalToGlobalEnd(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr);
897fcfd50ebSBarry Smith   ierr = VecDestroy(&natural);CHKERRQ(ierr);
898aa219208SBarry Smith   ierr = PetscInfo(xin,"Loading vector from natural ordering into DMDA\n");CHKERRQ(ierr);
89947c6ae99SBarry Smith   ierr = PetscOptionsGetInt(((PetscObject)xin)->prefix,"-vecload_block_size",&bs,&flag);CHKERRQ(ierr);
90047c6ae99SBarry Smith   if (flag && bs != dd->w) {
901aa219208SBarry Smith     ierr = PetscInfo2(xin,"Block size in file %D not equal to DMDA's dof %D\n",bs,dd->w);CHKERRQ(ierr);
90247c6ae99SBarry Smith   }
90347c6ae99SBarry Smith   PetscFunctionReturn(0);
90447c6ae99SBarry Smith }
90547c6ae99SBarry Smith 
90647c6ae99SBarry Smith #undef __FUNCT__
90747c6ae99SBarry Smith #define __FUNCT__ "VecLoad_Default_DA"
9087087cfbeSBarry Smith PetscErrorCode  VecLoad_Default_DA(Vec xin, PetscViewer viewer)
90947c6ae99SBarry Smith {
91047c6ae99SBarry Smith   PetscErrorCode ierr;
9119a42bb27SBarry Smith   DM             da;
91247c6ae99SBarry Smith   PetscBool      isbinary;
91347c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
91447c6ae99SBarry Smith   PetscBool ishdf5;
91547c6ae99SBarry Smith #endif
91647c6ae99SBarry Smith 
91747c6ae99SBarry Smith   PetscFunctionBegin;
918c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
919ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
92047c6ae99SBarry Smith 
92147c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
922251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr);
92347c6ae99SBarry Smith #endif
924251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
92547c6ae99SBarry Smith 
92647c6ae99SBarry Smith   if (isbinary) {
92747c6ae99SBarry Smith     ierr = VecLoad_Binary_DA(xin,viewer);CHKERRQ(ierr);
92847c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
92947c6ae99SBarry Smith   } else if (ishdf5) {
93047c6ae99SBarry Smith     ierr = VecLoad_HDF5_DA(xin,viewer);CHKERRQ(ierr);
93147c6ae99SBarry Smith #endif
932d34fcf5fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for vector loading", ((PetscObject)viewer)->type_name);
93347c6ae99SBarry Smith   PetscFunctionReturn(0);
93447c6ae99SBarry Smith }
935