xref: /petsc/src/dm/impls/da/gr2.c (revision b05fc00097cdf6db9066065a4de0654ad03e32cd)
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 {
1547c6ae99SBarry Smith   PetscInt          m,n,step,k;
16*b05fc000SLisandro Dalcin   PetscReal         min,max;
175edff71fSBarry Smith   const 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;
34*b05fc000SLisandro Dalcin   PetscReal         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;
365edff71fSBarry Smith   const 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   min  = zctx->min;
48f3f0eb19SBarry Smith   max  = zctx->max;
4947c6ae99SBarry Smith 
5047c6ae99SBarry Smith   /* PetscDraw the contour plot patch */
5147c6ae99SBarry Smith   for (j=0; j<n-1; j++) {
5247c6ae99SBarry Smith     for (i=0; i<m-1; i++) {
530e4fe250SBarry Smith       id   = i+j*m;
540e4fe250SBarry Smith       x1   = PetscRealPart(xy[2*id]);
550e4fe250SBarry Smith       y_1  = PetscRealPart(xy[2*id+1]);
56*b05fc000SLisandro Dalcin       c1   = PetscDrawRealToColor(PetscRealPart(v[k+step*id]),min,max);
570e4fe250SBarry Smith       xmin = PetscMin(xmin,x1);
580e4fe250SBarry Smith       ymin = PetscMin(ymin,y_1);
590e4fe250SBarry Smith       xmax = PetscMax(xmax,x1);
600e4fe250SBarry Smith       ymax = PetscMax(ymax,y_1);
610e4fe250SBarry Smith 
620e4fe250SBarry Smith       id   = i+j*m+1;
630e4fe250SBarry Smith       x2   = PetscRealPart(xy[2*id]);
64a39a4c7dSLisandro Dalcin       y2   = PetscRealPart(xy[2*id+1]);
65*b05fc000SLisandro Dalcin       c2   = PetscDrawRealToColor(PetscRealPart(v[k+step*id]),min,max);
660e4fe250SBarry Smith       xmin = PetscMin(xmin,x2);
67a39a4c7dSLisandro Dalcin       ymin = PetscMin(ymin,y2);
680e4fe250SBarry Smith       xmax = PetscMax(xmax,x2);
69a39a4c7dSLisandro Dalcin       ymax = PetscMax(ymax,y2);
700e4fe250SBarry Smith 
710e4fe250SBarry Smith       id   = i+j*m+1+m;
72a39a4c7dSLisandro Dalcin       x3   = PetscRealPart(xy[2*id]);
730e4fe250SBarry Smith       y3   = PetscRealPart(xy[2*id+1]);
74*b05fc000SLisandro Dalcin       c3   = PetscDrawRealToColor(PetscRealPart(v[k+step*id]),min,max);
75a39a4c7dSLisandro Dalcin       xmin = PetscMin(xmin,x3);
760e4fe250SBarry Smith       ymin = PetscMin(ymin,y3);
77a39a4c7dSLisandro Dalcin       xmax = PetscMax(xmax,x3);
780e4fe250SBarry Smith       ymax = PetscMax(ymax,y3);
790e4fe250SBarry Smith 
800e4fe250SBarry Smith       id = i+j*m+m;
81a39a4c7dSLisandro Dalcin       x4 = PetscRealPart(xy[2*id]);
82a39a4c7dSLisandro Dalcin       y4 = PetscRealPart(xy[2*id+1]);
83*b05fc000SLisandro Dalcin       c4 = PetscDrawRealToColor(PetscRealPart(v[k+step*id]),min,max);
84a39a4c7dSLisandro Dalcin       xmin = PetscMin(xmin,x4);
85a39a4c7dSLisandro Dalcin       ymin = PetscMin(ymin,y4);
86a39a4c7dSLisandro Dalcin       xmax = PetscMax(xmax,x4);
87a39a4c7dSLisandro Dalcin       ymax = PetscMax(ymax,y4);
88f3f0eb19SBarry Smith 
8947c6ae99SBarry Smith       ierr = PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);CHKERRQ(ierr);
9047c6ae99SBarry Smith       ierr = PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4);CHKERRQ(ierr);
9147c6ae99SBarry Smith       if (zctx->showgrid) {
9247c6ae99SBarry Smith         ierr = PetscDrawLine(draw,x1,y_1,x2,y2,PETSC_DRAW_BLACK);CHKERRQ(ierr);
9347c6ae99SBarry Smith         ierr = PetscDrawLine(draw,x2,y2,x3,y3,PETSC_DRAW_BLACK);CHKERRQ(ierr);
9447c6ae99SBarry Smith         ierr = PetscDrawLine(draw,x3,y3,x4,y4,PETSC_DRAW_BLACK);CHKERRQ(ierr);
9547c6ae99SBarry Smith         ierr = PetscDrawLine(draw,x4,y4,x1,y_1,PETSC_DRAW_BLACK);CHKERRQ(ierr);
9647c6ae99SBarry Smith       }
9747c6ae99SBarry Smith     }
9847c6ae99SBarry Smith   }
99109c9344SBarry Smith   if (zctx->name0) {
100109c9344SBarry Smith     PetscReal xl,yl,xr,yr,x,y;
101109c9344SBarry Smith     ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
102109c9344SBarry Smith     x    = xl + .3*(xr - xl);
103109c9344SBarry Smith     xl   = xl + .01*(xr - xl);
104109c9344SBarry Smith     y    = yr - .3*(yr - yl);
105109c9344SBarry Smith     yl   = yl + .01*(yr - yl);
106109c9344SBarry Smith     ierr = PetscDrawString(draw,x,yl,PETSC_DRAW_BLACK,zctx->name0);CHKERRQ(ierr);
107109c9344SBarry Smith     ierr = PetscDrawStringVertical(draw,xl,y,PETSC_DRAW_BLACK,zctx->name1);CHKERRQ(ierr);
108109c9344SBarry Smith   }
1090e4fe250SBarry Smith   /*
1100e4fe250SBarry Smith      Ideally we would use the PetscDrawAxis object to manage displaying the coordinate limits
1110e4fe250SBarry Smith      but that may require some refactoring.
1120e4fe250SBarry Smith   */
113b2566f29SBarry Smith   ierr = MPIU_Allreduce(&xmin,&xminf,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)draw));CHKERRQ(ierr);
114b2566f29SBarry Smith   ierr = MPIU_Allreduce(&xmax,&xmaxf,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)draw));CHKERRQ(ierr);
115b2566f29SBarry Smith   ierr = MPIU_Allreduce(&ymin,&yminf,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)draw));CHKERRQ(ierr);
116b2566f29SBarry Smith   ierr = MPIU_Allreduce(&ymax,&ymaxf,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)draw));CHKERRQ(ierr);
1170e4fe250SBarry Smith   ierr = PetscSNPrintf(value,16,"%f",xminf);CHKERRQ(ierr);
1180e4fe250SBarry Smith   ierr = PetscDrawString(draw,xminf,yminf - .05*(ymaxf - yminf),PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
1190e4fe250SBarry Smith   ierr = PetscSNPrintf(value,16,"%f",xmaxf);CHKERRQ(ierr);
1200e4fe250SBarry Smith   ierr = PetscStrlen(value,&len);CHKERRQ(ierr);
1210298fd71SBarry Smith   ierr = PetscDrawStringGetSize(draw,&w,NULL);CHKERRQ(ierr);
1220e4fe250SBarry Smith   ierr = PetscDrawString(draw,xmaxf - len*w,yminf - .05*(ymaxf - yminf),PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
1230e4fe250SBarry Smith   ierr = PetscSNPrintf(value,16,"%f",yminf);CHKERRQ(ierr);
1240e4fe250SBarry Smith   ierr = PetscDrawString(draw,xminf - .05*(xmaxf - xminf),yminf,PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
1250e4fe250SBarry Smith   ierr = PetscSNPrintf(value,16,"%f",ymaxf);CHKERRQ(ierr);
1260e4fe250SBarry Smith   ierr = PetscDrawString(draw,xminf - .05*(xmaxf - xminf),ymaxf,PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
12747c6ae99SBarry Smith   PetscFunctionReturn(0);
12847c6ae99SBarry Smith }
12947c6ae99SBarry Smith 
13047c6ae99SBarry Smith #undef __FUNCT__
13147c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_Draw_DA2d"
13247c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA2d(Vec xin,PetscViewer viewer)
13347c6ae99SBarry Smith {
1349a42bb27SBarry Smith   DM                 da,dac,dag;
13547c6ae99SBarry Smith   PetscErrorCode     ierr;
13647c6ae99SBarry Smith   PetscMPIInt        rank;
137a74ba6f7SBarry Smith   PetscInt           N,s,M,w,ncoors = 4;
13847c6ae99SBarry Smith   const PetscInt     *lx,*ly;
13947c6ae99SBarry Smith   PetscReal          coors[4],ymin,ymax,xmin,xmax;
14047c6ae99SBarry Smith   PetscDraw          draw,popup;
14147c6ae99SBarry Smith   PetscBool          isnull,useports = PETSC_FALSE;
14247c6ae99SBarry Smith   MPI_Comm           comm;
14347c6ae99SBarry Smith   Vec                xlocal,xcoor,xcoorl;
144bff4a2f0SMatthew G. Knepley   DMBoundaryType     bx,by;
145aa219208SBarry Smith   DMDAStencilType    st;
14647c6ae99SBarry Smith   ZoomCtx            zctx;
1470298fd71SBarry Smith   PetscDrawViewPorts *ports = NULL;
14847c6ae99SBarry Smith   PetscViewerFormat  format;
14920d0051dSBarry Smith   PetscInt           *displayfields;
15067dd0837SBarry Smith   PetscInt           ndisplayfields,i,nbounds;
15167dd0837SBarry Smith   const PetscReal    *bounds;
15247c6ae99SBarry Smith 
15347c6ae99SBarry Smith   PetscFunctionBegin;
15447c6ae99SBarry Smith   zctx.showgrid = PETSC_FALSE;
1558865f1eaSKarl Rupp 
15647c6ae99SBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
15747c6ae99SBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
15803193ff8SBarry Smith   ierr = PetscViewerDrawGetBounds(viewer,&nbounds,&bounds);CHKERRQ(ierr);
15947c6ae99SBarry Smith 
160c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
161ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
16247c6ae99SBarry Smith 
16347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr);
16447c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
16547c6ae99SBarry Smith 
1661321219cSEthan Coon   ierr = DMDAGetInfo(da,0,&M,&N,0,&zctx.m,&zctx.n,0,&w,&s,&bx,&by,0,&st);CHKERRQ(ierr);
1670298fd71SBarry Smith   ierr = DMDAGetOwnershipRanges(da,&lx,&ly,NULL);CHKERRQ(ierr);
16847c6ae99SBarry Smith 
16947c6ae99SBarry Smith   /*
17047c6ae99SBarry Smith         Obtain a sequential vector that is going to contain the local values plus ONE layer of
171aa219208SBarry Smith      ghosted values to draw the graphics from. We also need its corresponding DMDA (dac) that will
17247c6ae99SBarry Smith      update the local values pluse ONE layer of ghost values.
17347c6ae99SBarry Smith   */
17447c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)da,"GraphicsGhosted",(PetscObject*)&xlocal);CHKERRQ(ierr);
17547c6ae99SBarry Smith   if (!xlocal) {
176bff4a2f0SMatthew G. Knepley     if (bx !=  DM_BOUNDARY_NONE || by !=  DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) {
17747c6ae99SBarry Smith       /*
17847c6ae99SBarry Smith          if original da is not of stencil width one, or periodic or not a box stencil then
179aa219208SBarry Smith          create a special DMDA to handle one level of ghost points for graphics
18047c6ae99SBarry Smith       */
181bff4a2f0SMatthew 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);
182aa219208SBarry Smith       ierr = PetscInfo(da,"Creating auxilary DMDA for managing graphics ghost points\n");CHKERRQ(ierr);
18347c6ae99SBarry Smith     } else {
18447c6ae99SBarry Smith       /* otherwise we can use the da we already have */
18547c6ae99SBarry Smith       dac = da;
18647c6ae99SBarry Smith     }
18747c6ae99SBarry Smith     /* create local vector for holding ghosted values used in graphics */
188564755cdSBarry Smith     ierr = DMCreateLocalVector(dac,&xlocal);CHKERRQ(ierr);
18947c6ae99SBarry Smith     if (dac != da) {
190aa219208SBarry Smith       /* don't keep any public reference of this DMDA, is is only available through xlocal */
191f7923d8aSBarry Smith       ierr = PetscObjectDereference((PetscObject)dac);CHKERRQ(ierr);
19247c6ae99SBarry Smith     } else {
19347c6ae99SBarry Smith       /* remove association between xlocal and da, because below we compose in the opposite
19447c6ae99SBarry Smith          direction and if we left this connect we'd get a loop, so the objects could
19547c6ae99SBarry Smith          never be destroyed */
196c688c046SMatthew G Knepley       ierr = PetscObjectRemoveReference((PetscObject)xlocal,"__PETSc_dm");CHKERRQ(ierr);
19747c6ae99SBarry Smith     }
19847c6ae99SBarry Smith     ierr = PetscObjectCompose((PetscObject)da,"GraphicsGhosted",(PetscObject)xlocal);CHKERRQ(ierr);
19947c6ae99SBarry Smith     ierr = PetscObjectDereference((PetscObject)xlocal);CHKERRQ(ierr);
20047c6ae99SBarry Smith   } else {
201bff4a2f0SMatthew G. Knepley     if (bx !=  DM_BOUNDARY_NONE || by !=  DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) {
202c688c046SMatthew G Knepley       ierr = VecGetDM(xlocal, &dac);CHKERRQ(ierr);
203f7923d8aSBarry Smith     } else {
204f7923d8aSBarry Smith       dac = da;
20547c6ae99SBarry Smith     }
20647c6ae99SBarry Smith   }
20747c6ae99SBarry Smith 
20847c6ae99SBarry Smith   /*
20947c6ae99SBarry Smith       Get local (ghosted) values of vector
21047c6ae99SBarry Smith   */
2119a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr);
2129a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr);
2135edff71fSBarry Smith   ierr = VecGetArrayRead(xlocal,&zctx.v);CHKERRQ(ierr);
21447c6ae99SBarry Smith 
21547c6ae99SBarry Smith   /* get coordinates of nodes */
2166636e97aSMatthew G Knepley   ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
21747c6ae99SBarry Smith   if (!xcoor) {
218aa219208SBarry Smith     ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,0.0);CHKERRQ(ierr);
2196636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
22047c6ae99SBarry Smith   }
22147c6ae99SBarry Smith 
22247c6ae99SBarry Smith   /*
22347c6ae99SBarry Smith       Determine the min and max  coordinates in plot
22447c6ae99SBarry Smith   */
2250298fd71SBarry Smith   ierr     = VecStrideMin(xcoor,0,NULL,&xmin);CHKERRQ(ierr);
2260298fd71SBarry Smith   ierr     = VecStrideMax(xcoor,0,NULL,&xmax);CHKERRQ(ierr);
2270298fd71SBarry Smith   ierr     = VecStrideMin(xcoor,1,NULL,&ymin);CHKERRQ(ierr);
2280298fd71SBarry Smith   ierr     = VecStrideMax(xcoor,1,NULL,&ymax);CHKERRQ(ierr);
22947c6ae99SBarry Smith   coors[0] = xmin - .05*(xmax- xmin); coors[2] = xmax + .05*(xmax - xmin);
23047c6ae99SBarry Smith   coors[1] = ymin - .05*(ymax- ymin); coors[3] = ymax + .05*(ymax - ymin);
23157622a8eSBarry Smith   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);
23247c6ae99SBarry Smith 
233c5929fdfSBarry Smith   ierr = PetscOptionsGetRealArray(NULL,NULL,"-draw_coordinates",coors,&ncoors,NULL);CHKERRQ(ierr);
234a74ba6f7SBarry Smith 
23547c6ae99SBarry Smith   /*
23647c6ae99SBarry Smith        get local ghosted version of coordinates
23747c6ae99SBarry Smith   */
23847c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr);
23947c6ae99SBarry Smith   if (!xcoorl) {
240aa219208SBarry Smith     /* create DMDA to get local version of graphics */
241bff4a2f0SMatthew 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);
242aa219208SBarry Smith     ierr = PetscInfo(dag,"Creating auxilary DMDA for managing graphics coordinates ghost points\n");CHKERRQ(ierr);
243564755cdSBarry Smith     ierr = DMCreateLocalVector(dag,&xcoorl);CHKERRQ(ierr);
24447c6ae99SBarry Smith     ierr = PetscObjectCompose((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject)xcoorl);CHKERRQ(ierr);
245f7923d8aSBarry Smith     ierr = PetscObjectDereference((PetscObject)dag);CHKERRQ(ierr);
24647c6ae99SBarry Smith     ierr = PetscObjectDereference((PetscObject)xcoorl);CHKERRQ(ierr);
24747c6ae99SBarry Smith   } else {
248c688c046SMatthew G Knepley     ierr = VecGetDM(xcoorl,&dag);CHKERRQ(ierr);
24947c6ae99SBarry Smith   }
2509a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
2519a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
2525edff71fSBarry Smith   ierr = VecGetArrayRead(xcoorl,&zctx.xy);CHKERRQ(ierr);
25347c6ae99SBarry Smith 
25447c6ae99SBarry Smith   /*
25547c6ae99SBarry Smith         Get information about size of area each processor must do graphics for
25647c6ae99SBarry Smith   */
2571321219cSEthan Coon   ierr = DMDAGetInfo(dac,0,&M,&N,0,0,0,0,&zctx.step,0,&bx,&by,0,0);CHKERRQ(ierr);
258f7923d8aSBarry Smith   ierr = DMDAGetGhostCorners(dac,0,0,0,&zctx.m,&zctx.n,0);CHKERRQ(ierr);
25947c6ae99SBarry Smith 
260c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-draw_contour_grid",&zctx.showgrid,NULL);CHKERRQ(ierr);
2614e6118eeSBarry Smith 
2624e6118eeSBarry Smith   ierr = DMDASelectFields(da,&ndisplayfields,&displayfields);CHKERRQ(ierr);
26347c6ae99SBarry Smith 
26447c6ae99SBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
265c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-draw_ports",&useports,NULL);CHKERRQ(ierr);
26647c6ae99SBarry Smith   if (useports || format == PETSC_VIEWER_DRAW_PORTS) {
26747c6ae99SBarry Smith     ierr       = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
26820d0051dSBarry Smith     ierr       = PetscDrawViewPortsCreate(draw,ndisplayfields,&ports);CHKERRQ(ierr);
269109c9344SBarry Smith     zctx.name0 = 0;
270109c9344SBarry Smith     zctx.name1 = 0;
271109c9344SBarry Smith   } else {
272109c9344SBarry Smith     ierr = DMDAGetCoordinateName(da,0,&zctx.name0);CHKERRQ(ierr);
273109c9344SBarry Smith     ierr = DMDAGetCoordinateName(da,1,&zctx.name1);CHKERRQ(ierr);
27447c6ae99SBarry Smith   }
27520d0051dSBarry Smith 
27647c6ae99SBarry Smith   /*
27747c6ae99SBarry Smith      Loop over each field; drawing each in a different window
27847c6ae99SBarry Smith   */
27920d0051dSBarry Smith   for (i=0; i<ndisplayfields; i++) {
28020d0051dSBarry Smith     zctx.k = displayfields[i];
28147c6ae99SBarry Smith     if (useports) {
28220d0051dSBarry Smith       ierr = PetscDrawViewPortsSet(ports,i);CHKERRQ(ierr);
2839332fd86SBarry Smith     } else {
2849332fd86SBarry Smith       ierr = PetscViewerDrawGetDraw(viewer,i,&draw);CHKERRQ(ierr);
28547c6ae99SBarry Smith     }
28647c6ae99SBarry Smith 
28747c6ae99SBarry Smith     /*
28847c6ae99SBarry Smith         Determine the min and max color in plot
28947c6ae99SBarry Smith     */
2900298fd71SBarry Smith     ierr = VecStrideMin(xin,zctx.k,NULL,&zctx.min);CHKERRQ(ierr);
2910298fd71SBarry Smith     ierr = VecStrideMax(xin,zctx.k,NULL,&zctx.max);CHKERRQ(ierr);
29267dd0837SBarry Smith     if (zctx.k < nbounds) {
293f3f0eb19SBarry Smith       zctx.min = bounds[2*zctx.k];
294f3f0eb19SBarry Smith       zctx.max = bounds[2*zctx.k+1];
29567dd0837SBarry Smith     }
29647c6ae99SBarry Smith     if (zctx.min == zctx.max) {
29747c6ae99SBarry Smith       zctx.min -= 1.e-12;
29847c6ae99SBarry Smith       zctx.max += 1.e-12;
29947c6ae99SBarry Smith     }
30047c6ae99SBarry Smith 
30147c6ae99SBarry Smith     if (!rank) {
30247c6ae99SBarry Smith       const char *title;
30347c6ae99SBarry Smith 
304aa219208SBarry Smith       ierr = DMDAGetFieldName(da,zctx.k,&title);CHKERRQ(ierr);
30547c6ae99SBarry Smith       if (title) {
30647c6ae99SBarry Smith         ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr);
30747c6ae99SBarry Smith       }
30847c6ae99SBarry Smith     }
30947c6ae99SBarry Smith     ierr = PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr);
31057622a8eSBarry Smith     ierr = PetscInfo2(da,"DMDA 2d contour plot min %g max %g\n",(double)zctx.min,(double)zctx.max);CHKERRQ(ierr);
31147c6ae99SBarry Smith 
31247c6ae99SBarry Smith     ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
31347c6ae99SBarry Smith     if (popup) {ierr = PetscDrawScalePopup(popup,zctx.min,zctx.max);CHKERRQ(ierr);}
31447c6ae99SBarry Smith 
31547c6ae99SBarry Smith     ierr = PetscDrawZoom(draw,VecView_MPI_Draw_DA2d_Zoom,&zctx);CHKERRQ(ierr);
31647c6ae99SBarry Smith   }
31720d0051dSBarry Smith   ierr = PetscFree(displayfields);CHKERRQ(ierr);
3186bf464f9SBarry Smith   ierr = PetscDrawViewPortsDestroy(ports);CHKERRQ(ierr);
31947c6ae99SBarry Smith 
3205edff71fSBarry Smith   ierr = VecRestoreArrayRead(xcoorl,&zctx.xy);CHKERRQ(ierr);
3215edff71fSBarry Smith   ierr = VecRestoreArrayRead(xlocal,&zctx.v);CHKERRQ(ierr);
32247c6ae99SBarry Smith   PetscFunctionReturn(0);
32347c6ae99SBarry Smith }
32447c6ae99SBarry Smith 
3250da24e51SJuha Jäykkä #if defined(PETSC_HAVE_HDF5)
3260da24e51SJuha Jäykkä #undef __FUNCT__
3270da24e51SJuha Jäykkä #define __FUNCT__ "VecGetHDF5ChunkSize"
328c73cfb54SMatthew G. Knepley static PetscErrorCode VecGetHDF5ChunkSize(DM_DA *da, Vec xin, PetscInt dimension, PetscInt timestep, hsize_t *chunkDims)
3290da24e51SJuha Jäykkä {
330d4190030SJed Brown   PetscMPIInt    comm_size;
331d4190030SJed Brown   PetscErrorCode ierr;
332561a82e7SJed Brown   hsize_t        chunk_size, target_size, dim;
333561a82e7SJed Brown   hsize_t        vec_size = sizeof(PetscScalar)*da->M*da->N*da->P*da->w;
334561a82e7SJed Brown   hsize_t        avg_local_vec_size,KiB = 1024,MiB = KiB*KiB,GiB = MiB*KiB,min_size = MiB;
335561a82e7SJed Brown   hsize_t        max_chunks = 64*KiB;                                              /* HDF5 internal limitation */
336561a82e7SJed Brown   hsize_t        max_chunk_size = 4*GiB;                                           /* HDF5 internal limitation */
33784179ffaSJed Brown   PetscInt       zslices=da->p, yslices=da->n, xslices=da->m;
3380da24e51SJuha Jäykkä 
339cb5c4f94SJuha Jäykkä   PetscFunctionBegin;
340d4190030SJed Brown   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)xin), &comm_size);CHKERRQ(ierr);
3410da24e51SJuha 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 */
3420da24e51SJuha Jäykkä 
343794a961bSBarry 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)))));
3447d310018SBarry 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 */
3457d310018SBarry 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);
3460da24e51SJuha Jäykkä 
347cb5c4f94SJuha Jäykkä   /*
348cb5c4f94SJuha Jäykkä    if size/rank > max_chunk_size, we need radical measures: even going down to
349cb5c4f94SJuha Jäykkä    avg_local_vec_size is not enough, so we simply use chunk size of 4 GiB no matter
350cb5c4f94SJuha Jäykkä    what, composed in the most efficient way possible.
351cb5c4f94SJuha Jäykkä    N.B. this minimises the number of chunks, which may or may not be the optimal
352cb5c4f94SJuha Jäykkä    solution. In a BG, for example, the optimal solution is probably to make # chunks = #
353cb5c4f94SJuha Jäykkä    IO nodes involved, but this author has no access to a BG to figure out how to
354cb5c4f94SJuha Jäykkä    reliably find the right number. And even then it may or may not be enough.
355cb5c4f94SJuha Jäykkä    */
3560da24e51SJuha Jäykkä   if (avg_local_vec_size > max_chunk_size) {
357cb5c4f94SJuha Jäykkä     /* check if we can just split local z-axis: is that enough? */
35884179ffaSJed Brown     zslices = (PetscInt)ceil(vec_size*1.0/(da->p*max_chunk_size))*zslices;
3590da24e51SJuha Jäykkä     if (zslices > da->P) {
360cb5c4f94SJuha Jäykkä       /* lattice is too large in xy-directions, splitting z only is not enough */
3610da24e51SJuha Jäykkä       zslices = da->P;
36284179ffaSJed Brown       yslices= (PetscInt)ceil(vec_size*1.0/(zslices*da->n*max_chunk_size))*yslices;
3630da24e51SJuha Jäykkä       if (yslices > da->N) {
364cb5c4f94SJuha Jäykkä 	/* lattice is too large in x-direction, splitting along z, y is not enough */
3650da24e51SJuha Jäykkä 	yslices = da->N;
36684179ffaSJed Brown 	xslices= (PetscInt)ceil(vec_size*1.0/(zslices*yslices*da->m*max_chunk_size))*xslices;
3670da24e51SJuha Jäykkä       }
3680da24e51SJuha Jäykkä     }
3690da24e51SJuha Jäykkä     dim = 0;
3700da24e51SJuha Jäykkä     if (timestep >= 0) {
3710da24e51SJuha Jäykkä       ++dim;
3720da24e51SJuha Jäykkä     }
373cb5c4f94SJuha Jäykkä     /* prefer to split z-axis, even down to planar slices */
374c73cfb54SMatthew G. Knepley     if (dimension == 3) {
375cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->P/zslices;
376cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->N/yslices;
377cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->M/xslices;
3780da24e51SJuha Jäykkä     } else {
379cb5c4f94SJuha Jäykkä       /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */
380cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->N/yslices;
381cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->M/xslices;
3820da24e51SJuha Jäykkä     }
3830da24e51SJuha 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);
3840da24e51SJuha Jäykkä   } else {
3850da24e51SJuha Jäykkä     if (target_size < chunk_size) {
386cb5c4f94SJuha Jäykkä       /* only change the defaults if target_size < chunk_size */
3870da24e51SJuha Jäykkä       dim = 0;
3880da24e51SJuha Jäykkä       if (timestep >= 0) {
3890da24e51SJuha Jäykkä 	++dim;
3900da24e51SJuha Jäykkä       }
391cb5c4f94SJuha Jäykkä       /* prefer to split z-axis, even down to planar slices */
392c73cfb54SMatthew G. Knepley       if (dimension == 3) {
393cb5c4f94SJuha Jäykkä 	/* try splitting the z-axis to core-size bits, i.e. divide chunk size by # comm_size in z-direction */
3940da24e51SJuha Jäykkä 	if (target_size >= chunk_size/da->p) {
395cb5c4f94SJuha Jäykkä 	  /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */
3960da24e51SJuha Jäykkä 	  chunkDims[dim] = (hsize_t) ceil(da->P*1.0/da->p);
3970da24e51SJuha Jäykkä 	} else {
398cb5c4f94SJuha Jäykkä 	  /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be
399cb5c4f94SJuha Jäykkä            radical and let everyone write all they've got */
4000da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->P*1.0/da->p);
4010da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n);
4020da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m);
4030da24e51SJuha Jäykkä 	}
4040da24e51SJuha Jäykkä       } else {
405cb5c4f94SJuha Jäykkä 	/* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */
4060da24e51SJuha Jäykkä 	if (target_size >= chunk_size/da->n) {
407cb5c4f94SJuha Jäykkä 	  /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */
4080da24e51SJuha Jäykkä 	  chunkDims[dim] = (hsize_t) ceil(da->N*1.0/da->n);
4090da24e51SJuha Jäykkä 	} else {
410cb5c4f94SJuha Jäykkä 	  /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be
411cb5c4f94SJuha Jäykkä 	   radical and let everyone write all they've got */
4120da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n);
4130da24e51SJuha Jäykkä 	  chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m);
4140da24e51SJuha Jäykkä 	}
4150da24e51SJuha Jäykkä 
4160da24e51SJuha Jäykkä       }
4170da24e51SJuha 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);
4180da24e51SJuha Jäykkä     } else {
419cb5c4f94SJuha Jäykkä       /* precomputed chunks are fine, we don't need to do anything */
4200da24e51SJuha Jäykkä     }
4210da24e51SJuha Jäykkä   }
4220da24e51SJuha Jäykkä   PetscFunctionReturn(0);
4230da24e51SJuha Jäykkä }
4240da24e51SJuha Jäykkä #endif
42547c6ae99SBarry Smith 
42647c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
42747c6ae99SBarry Smith #undef __FUNCT__
42847c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_HDF5_DA"
42947c6ae99SBarry Smith PetscErrorCode VecView_MPI_HDF5_DA(Vec xin,PetscViewer viewer)
43047c6ae99SBarry Smith {
4319b2a5a86SJed Brown   DM                dm;
4329b2a5a86SJed Brown   DM_DA             *da;
43347c6ae99SBarry Smith   hid_t             filespace;  /* file dataspace identifier */
4348e2ae6d7SMichael Kraus   hid_t             chunkspace; /* chunk dataset property identifier */
43547c6ae99SBarry Smith   hid_t             plist_id;   /* property list identifier */
43647c6ae99SBarry Smith   hid_t             dset_id;    /* dataset identifier */
43747c6ae99SBarry Smith   hid_t             memspace;   /* memory dataspace identifier */
43847c6ae99SBarry Smith   hid_t             file_id;
43915214e8eSMatthew G Knepley   hid_t             group;
4409a0502c6SHåkon Strandenes   hid_t             memscalartype; /* scalar type for mem (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
4419a0502c6SHåkon Strandenes   hid_t             filescalartype; /* scalar type for file (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
442d9a4edebSJed Brown   hsize_t           dim;
4430da24e51SJuha 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  */
444c73cfb54SMatthew G. Knepley   PetscInt          timestep, dimension;
4455edff71fSBarry Smith   const PetscScalar *x;
4468e2ae6d7SMichael Kraus   const char        *vecname;
44715214e8eSMatthew G Knepley   PetscErrorCode    ierr;
44882ea9b62SBarry Smith   PetscBool         dim2;
4499a0502c6SHåkon Strandenes   PetscBool         spoutput;
45047c6ae99SBarry Smith 
45147c6ae99SBarry Smith   PetscFunctionBegin;
45215214e8eSMatthew G Knepley   ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr);
4538e2ae6d7SMichael Kraus   ierr = PetscViewerHDF5GetTimestep(viewer, &timestep);CHKERRQ(ierr);
45482ea9b62SBarry Smith   ierr = PetscViewerHDF5GetBaseDimension2(viewer,&dim2);CHKERRQ(ierr);
4559a0502c6SHåkon Strandenes   ierr = PetscViewerHDF5GetSPOutput(viewer,&spoutput);CHKERRQ(ierr);
45615214e8eSMatthew G Knepley 
457c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&dm);CHKERRQ(ierr);
458ce94432eSBarry Smith   if (!dm) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
4599b2a5a86SJed Brown   da = (DM_DA*)dm->data;
460c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dimension);CHKERRQ(ierr);
46147c6ae99SBarry Smith 
4628e2ae6d7SMichael Kraus   /* Create the dataspace for the dataset.
4638e2ae6d7SMichael Kraus    *
4648e2ae6d7SMichael Kraus    * dims - holds the current dimensions of the dataset
4658e2ae6d7SMichael Kraus    *
4668e2ae6d7SMichael Kraus    * maxDims - holds the maximum dimensions of the dataset (unlimited
4678e2ae6d7SMichael Kraus    * for the number of time steps with the current dimensions for the
4688e2ae6d7SMichael Kraus    * other dimensions; so only additional time steps can be added).
4698e2ae6d7SMichael Kraus    *
4708e2ae6d7SMichael Kraus    * chunkDims - holds the size of a single time step (required to
4718e2ae6d7SMichael Kraus    * permit extending dataset).
4728e2ae6d7SMichael Kraus    */
4738e2ae6d7SMichael Kraus   dim = 0;
4748e2ae6d7SMichael Kraus   if (timestep >= 0) {
4758e2ae6d7SMichael Kraus     dims[dim]      = timestep+1;
4768e2ae6d7SMichael Kraus     maxDims[dim]   = H5S_UNLIMITED;
4778e2ae6d7SMichael Kraus     chunkDims[dim] = 1;
4788e2ae6d7SMichael Kraus     ++dim;
4798e2ae6d7SMichael Kraus   }
480c73cfb54SMatthew G. Knepley   if (dimension == 3) {
481acba2ac6SBarry Smith     ierr           = PetscHDF5IntCast(da->P,dims+dim);CHKERRQ(ierr);
4828e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4838e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4848e2ae6d7SMichael Kraus     ++dim;
4858e2ae6d7SMichael Kraus   }
486c73cfb54SMatthew G. Knepley   if (dimension > 1) {
487acba2ac6SBarry Smith     ierr           = PetscHDF5IntCast(da->N,dims+dim);CHKERRQ(ierr);
4888e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4898e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4908e2ae6d7SMichael Kraus     ++dim;
4918e2ae6d7SMichael Kraus   }
492acba2ac6SBarry Smith   ierr           = PetscHDF5IntCast(da->M,dims+dim);CHKERRQ(ierr);
4938e2ae6d7SMichael Kraus   maxDims[dim]   = dims[dim];
4948e2ae6d7SMichael Kraus   chunkDims[dim] = dims[dim];
4958e2ae6d7SMichael Kraus   ++dim;
49682ea9b62SBarry Smith   if (da->w > 1 || dim2) {
497acba2ac6SBarry Smith     ierr           = PetscHDF5IntCast(da->w,dims+dim);CHKERRQ(ierr);
4988e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4998e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
5008e2ae6d7SMichael Kraus     ++dim;
5018e2ae6d7SMichael Kraus   }
50247c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5038e2ae6d7SMichael Kraus   dims[dim]      = 2;
5048e2ae6d7SMichael Kraus   maxDims[dim]   = dims[dim];
5058e2ae6d7SMichael Kraus   chunkDims[dim] = dims[dim];
5068e2ae6d7SMichael Kraus   ++dim;
50747c6ae99SBarry Smith #endif
5080da24e51SJuha Jäykkä 
509c73cfb54SMatthew G. Knepley   ierr = VecGetHDF5ChunkSize(da, xin, dimension, timestep, chunkDims); CHKERRQ(ierr);
5100da24e51SJuha Jäykkä 
511729ab6d9SBarry Smith   PetscStackCallHDF5Return(filespace,H5Screate_simple,(dim, dims, maxDims));
51247c6ae99SBarry Smith 
51315214e8eSMatthew G Knepley #if defined(PETSC_USE_REAL_SINGLE)
5149a0502c6SHåkon Strandenes   memscalartype = H5T_NATIVE_FLOAT;
5159a0502c6SHåkon Strandenes   filescalartype = H5T_NATIVE_FLOAT;
51615214e8eSMatthew G Knepley #elif defined(PETSC_USE_REAL___FLOAT128)
51715214e8eSMatthew G Knepley #error "HDF5 output with 128 bit floats not supported."
51815214e8eSMatthew G Knepley #else
5199a0502c6SHåkon Strandenes   memscalartype = H5T_NATIVE_DOUBLE;
5209a0502c6SHåkon Strandenes   if (spoutput == PETSC_TRUE) filescalartype = H5T_NATIVE_FLOAT;
5219a0502c6SHåkon Strandenes   else filescalartype = H5T_NATIVE_DOUBLE;
52215214e8eSMatthew G Knepley #endif
52315214e8eSMatthew G Knepley 
52447c6ae99SBarry Smith   /* Create the dataset with default properties and close filespace */
52547c6ae99SBarry Smith   ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr);
52615214e8eSMatthew G Knepley   if (!H5Lexists(group, vecname, H5P_DEFAULT)) {
5278e2ae6d7SMichael Kraus     /* Create chunk */
528729ab6d9SBarry Smith     PetscStackCallHDF5Return(chunkspace,H5Pcreate,(H5P_DATASET_CREATE));
529729ab6d9SBarry Smith     PetscStackCallHDF5(H5Pset_chunk,(chunkspace, dim, chunkDims));
5308e2ae6d7SMichael Kraus 
53147c6ae99SBarry Smith #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
5329a0502c6SHåkon Strandenes     PetscStackCallHDF5Return(dset_id,H5Dcreate2,(group, vecname, filescalartype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT));
53347c6ae99SBarry Smith #else
5349a0502c6SHåkon Strandenes     PetscStackCallHDF5Return(dset_id,H5Dcreate,(group, vecname, filescalartype, filespace, H5P_DEFAULT));
53547c6ae99SBarry Smith #endif
53615214e8eSMatthew G Knepley   } else {
537729ab6d9SBarry Smith     PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT));
538729ab6d9SBarry Smith     PetscStackCallHDF5(H5Dset_extent,(dset_id, dims));
53915214e8eSMatthew G Knepley   }
540729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
54147c6ae99SBarry Smith 
54247c6ae99SBarry Smith   /* Each process defines a dataset and writes it to the hyperslab in the file */
5438e2ae6d7SMichael Kraus   dim = 0;
5448e2ae6d7SMichael Kraus   if (timestep >= 0) {
5458e2ae6d7SMichael Kraus     offset[dim] = timestep;
5468e2ae6d7SMichael Kraus     ++dim;
5478e2ae6d7SMichael Kraus   }
548c73cfb54SMatthew G. Knepley   if (dimension == 3) {ierr = PetscHDF5IntCast(da->zs,offset + dim++);CHKERRQ(ierr);}
549c73cfb54SMatthew G. Knepley   if (dimension > 1)  {ierr = PetscHDF5IntCast(da->ys,offset + dim++);CHKERRQ(ierr);}
550acba2ac6SBarry Smith   ierr = PetscHDF5IntCast(da->xs/da->w,offset + dim++);CHKERRQ(ierr);
55182ea9b62SBarry Smith   if (da->w > 1 || dim2) offset[dim++] = 0;
55247c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5538e2ae6d7SMichael Kraus   offset[dim++] = 0;
55447c6ae99SBarry Smith #endif
5558e2ae6d7SMichael Kraus   dim = 0;
5568e2ae6d7SMichael Kraus   if (timestep >= 0) {
5578e2ae6d7SMichael Kraus     count[dim] = 1;
5588e2ae6d7SMichael Kraus     ++dim;
5598e2ae6d7SMichael Kraus   }
560c73cfb54SMatthew G. Knepley   if (dimension == 3) {ierr = PetscHDF5IntCast(da->ze - da->zs,count + dim++);CHKERRQ(ierr);}
561c73cfb54SMatthew G. Knepley   if (dimension > 1)  {ierr = PetscHDF5IntCast(da->ye - da->ys,count + dim++);CHKERRQ(ierr);}
562acba2ac6SBarry Smith   ierr = PetscHDF5IntCast((da->xe - da->xs)/da->w,count + dim++);CHKERRQ(ierr);
56382ea9b62SBarry Smith   if (da->w > 1 || dim2) {ierr = PetscHDF5IntCast(da->w,count + dim++);CHKERRQ(ierr);}
56447c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5658e2ae6d7SMichael Kraus   count[dim++] = 2;
56647c6ae99SBarry Smith #endif
567729ab6d9SBarry Smith   PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL));
568729ab6d9SBarry Smith   PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id));
569729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL));
57047c6ae99SBarry Smith 
57147c6ae99SBarry Smith   /* Create property list for collective dataset write */
572729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_DATASET_XFER));
57347c6ae99SBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
574729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_dxpl_mpio,(plist_id, H5FD_MPIO_COLLECTIVE));
57547c6ae99SBarry Smith #endif
57647c6ae99SBarry Smith   /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
57747c6ae99SBarry Smith 
5785edff71fSBarry Smith   ierr   = VecGetArrayRead(xin, &x);CHKERRQ(ierr);
5799a0502c6SHåkon Strandenes   PetscStackCallHDF5(H5Dwrite,(dset_id, memscalartype, memspace, filespace, plist_id, x));
580729ab6d9SBarry Smith   PetscStackCallHDF5(H5Fflush,(file_id, H5F_SCOPE_GLOBAL));
5815edff71fSBarry Smith   ierr   = VecRestoreArrayRead(xin, &x);CHKERRQ(ierr);
58247c6ae99SBarry Smith 
58347c6ae99SBarry Smith   /* Close/release resources */
58415214e8eSMatthew G Knepley   if (group != file_id) {
585729ab6d9SBarry Smith     PetscStackCallHDF5(H5Gclose,(group));
58615214e8eSMatthew G Knepley   }
587729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
588729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
589729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(memspace));
590729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dset_id));
59147c6ae99SBarry Smith   ierr   = PetscInfo1(xin,"Wrote Vec object with name %s\n",vecname);CHKERRQ(ierr);
59247c6ae99SBarry Smith   PetscFunctionReturn(0);
59347c6ae99SBarry Smith }
59447c6ae99SBarry Smith #endif
59547c6ae99SBarry Smith 
59609573ac7SBarry Smith extern PetscErrorCode VecView_MPI_Draw_DA1d(Vec,PetscViewer);
59747c6ae99SBarry Smith 
59847c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
59947c6ae99SBarry Smith #undef __FUNCT__
600aa219208SBarry Smith #define __FUNCT__ "DMDAArrayMPIIO"
601aa219208SBarry Smith static PetscErrorCode DMDAArrayMPIIO(DM da,PetscViewer viewer,Vec xin,PetscBool write)
60247c6ae99SBarry Smith {
60347c6ae99SBarry Smith   PetscErrorCode    ierr;
60447c6ae99SBarry Smith   MPI_File          mfdes;
60547c6ae99SBarry Smith   PetscMPIInt       gsizes[4],lsizes[4],lstarts[4],asiz,dof;
60647c6ae99SBarry Smith   MPI_Datatype      view;
60747c6ae99SBarry Smith   const PetscScalar *array;
60847c6ae99SBarry Smith   MPI_Offset        off;
60947c6ae99SBarry Smith   MPI_Aint          ub,ul;
61047c6ae99SBarry Smith   PetscInt          type,rows,vecrows,tr[2];
61147c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
6120c169764Sdmay   PetscBool         skipheader;
61347c6ae99SBarry Smith 
61447c6ae99SBarry Smith   PetscFunctionBegin;
61547c6ae99SBarry Smith   ierr = VecGetSize(xin,&vecrows);CHKERRQ(ierr);
6160c169764Sdmay   ierr = PetscViewerBinaryGetSkipHeader(viewer,&skipheader);CHKERRQ(ierr);
61747c6ae99SBarry Smith   if (!write) {
61847c6ae99SBarry Smith     /* Read vector header. */
6190c169764Sdmay     if (!skipheader) {
620060da220SMatthew G. Knepley       ierr = PetscViewerBinaryRead(viewer,tr,2,NULL,PETSC_INT);CHKERRQ(ierr);
62147c6ae99SBarry Smith       type = tr[0];
62247c6ae99SBarry Smith       rows = tr[1];
623ce94432eSBarry Smith       if (type != VEC_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Not vector next in file");
624ce94432eSBarry Smith       if (rows != vecrows) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Vector in file not same size as DMDA vector");
6250c169764Sdmay     }
62647c6ae99SBarry Smith   } else {
62747c6ae99SBarry Smith     tr[0] = VEC_FILE_CLASSID;
62847c6ae99SBarry Smith     tr[1] = vecrows;
6290c169764Sdmay     if (!skipheader) {
63047c6ae99SBarry Smith       ierr  = PetscViewerBinaryWrite(viewer,tr,2,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
63147c6ae99SBarry Smith     }
6320c169764Sdmay   }
63347c6ae99SBarry Smith 
6344dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->w,&dof);CHKERRQ(ierr);
6354dc2109aSBarry Smith   gsizes[0]  = dof;
6364dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->M,gsizes+1);CHKERRQ(ierr);
6374dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->N,gsizes+2);CHKERRQ(ierr);
638334634e2SJed Brown   ierr       = PetscMPIIntCast(dd->P,gsizes+3);CHKERRQ(ierr);
6394dc2109aSBarry Smith   lsizes[0]  = dof;
6404dc2109aSBarry Smith   ierr       = PetscMPIIntCast((dd->xe-dd->xs)/dof,lsizes+1);CHKERRQ(ierr);
6414dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->ye-dd->ys,lsizes+2);CHKERRQ(ierr);
6424dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->ze-dd->zs,lsizes+3);CHKERRQ(ierr);
6434dc2109aSBarry Smith   lstarts[0] = 0;
6444dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->xs/dof,lstarts+1);CHKERRQ(ierr);
6454dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->ys,lstarts+2);CHKERRQ(ierr);
6464dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->zs,lstarts+3);CHKERRQ(ierr);
647c73cfb54SMatthew G. Knepley   ierr       = MPI_Type_create_subarray(da->dim+1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view);CHKERRQ(ierr);
64847c6ae99SBarry Smith   ierr       = MPI_Type_commit(&view);CHKERRQ(ierr);
64947c6ae99SBarry Smith 
65047c6ae99SBarry Smith   ierr = PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);CHKERRQ(ierr);
65147c6ae99SBarry Smith   ierr = PetscViewerBinaryGetMPIIOOffset(viewer,&off);CHKERRQ(ierr);
65247c6ae99SBarry Smith   ierr = MPI_File_set_view(mfdes,off,MPIU_SCALAR,view,(char*)"native",MPI_INFO_NULL);CHKERRQ(ierr);
65347c6ae99SBarry Smith   ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr);
65447c6ae99SBarry Smith   asiz = lsizes[1]*(lsizes[2] > 0 ? lsizes[2] : 1)*(lsizes[3] > 0 ? lsizes[3] : 1)*dof;
65547c6ae99SBarry Smith   if (write) {
65647c6ae99SBarry Smith     ierr = MPIU_File_write_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr);
65747c6ae99SBarry Smith   } else {
65847c6ae99SBarry Smith     ierr = MPIU_File_read_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr);
65947c6ae99SBarry Smith   }
66047c6ae99SBarry Smith   ierr = MPI_Type_get_extent(view,&ul,&ub);CHKERRQ(ierr);
66147c6ae99SBarry Smith   ierr = PetscViewerBinaryAddMPIIOOffset(viewer,ub);CHKERRQ(ierr);
66247c6ae99SBarry Smith   ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr);
66347c6ae99SBarry Smith   ierr = MPI_Type_free(&view);CHKERRQ(ierr);
66447c6ae99SBarry Smith   PetscFunctionReturn(0);
66547c6ae99SBarry Smith }
66647c6ae99SBarry Smith #endif
66747c6ae99SBarry Smith 
66847c6ae99SBarry Smith #undef __FUNCT__
66947c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_DA"
6707087cfbeSBarry Smith PetscErrorCode  VecView_MPI_DA(Vec xin,PetscViewer viewer)
67147c6ae99SBarry Smith {
6729a42bb27SBarry Smith   DM                da;
67347c6ae99SBarry Smith   PetscErrorCode    ierr;
67447c6ae99SBarry Smith   PetscInt          dim;
67547c6ae99SBarry Smith   Vec               natural;
6764061b8bfSJed Brown   PetscBool         isdraw,isvtk;
67747c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
67847c6ae99SBarry Smith   PetscBool         ishdf5;
67947c6ae99SBarry Smith #endif
6803f3fd955SJed Brown   const char        *prefix,*name;
681a261c58fSBarry Smith   PetscViewerFormat format;
68247c6ae99SBarry Smith 
68347c6ae99SBarry Smith   PetscFunctionBegin;
684c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
685ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
686251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
687251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr);
68847c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
689251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr);
69047c6ae99SBarry Smith #endif
69147c6ae99SBarry Smith   if (isdraw) {
6921321219cSEthan Coon     ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
69347c6ae99SBarry Smith     if (dim == 1) {
69447c6ae99SBarry Smith       ierr = VecView_MPI_Draw_DA1d(xin,viewer);CHKERRQ(ierr);
69547c6ae99SBarry Smith     } else if (dim == 2) {
69647c6ae99SBarry Smith       ierr = VecView_MPI_Draw_DA2d(xin,viewer);CHKERRQ(ierr);
697ce94432eSBarry Smith     } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot graphically view vector associated with this dimensional DMDA %D",dim);
6984061b8bfSJed Brown   } else if (isvtk) {           /* Duplicate the Vec and hold a reference to the DM */
6994061b8bfSJed Brown     Vec Y;
7004061b8bfSJed Brown     ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr);
7014061b8bfSJed Brown     ierr = VecDuplicate(xin,&Y);CHKERRQ(ierr);
702b51b94faSJed Brown     if (((PetscObject)xin)->name) {
703b51b94faSJed Brown       /* If xin was named, copy the name over to Y. The duplicate names are safe because nobody else will ever see Y. */
704b51b94faSJed Brown       ierr = PetscObjectSetName((PetscObject)Y,((PetscObject)xin)->name);CHKERRQ(ierr);
705b51b94faSJed Brown     }
7064061b8bfSJed Brown     ierr = VecCopy(xin,Y);CHKERRQ(ierr);
70762b69a3fSMatthew G Knepley     ierr = PetscViewerVTKAddField(viewer,(PetscObject)da,DMDAVTKWriteAll,PETSC_VTK_POINT_FIELD,(PetscObject)Y);CHKERRQ(ierr);
70847c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
70947c6ae99SBarry Smith   } else if (ishdf5) {
71047c6ae99SBarry Smith     ierr = VecView_MPI_HDF5_DA(xin,viewer);CHKERRQ(ierr);
71147c6ae99SBarry Smith #endif
71247c6ae99SBarry Smith   } else {
71347c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
71447c6ae99SBarry Smith     PetscBool isbinary,isMPIIO;
71547c6ae99SBarry Smith 
716251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
71747c6ae99SBarry Smith     if (isbinary) {
718bc196f7cSDave May       ierr = PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO);CHKERRQ(ierr);
71947c6ae99SBarry Smith       if (isMPIIO) {
720aa219208SBarry Smith         ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_TRUE);CHKERRQ(ierr);
72147c6ae99SBarry Smith         PetscFunctionReturn(0);
72247c6ae99SBarry Smith       }
72347c6ae99SBarry Smith     }
72447c6ae99SBarry Smith #endif
72547c6ae99SBarry Smith 
72647c6ae99SBarry Smith     /* call viewer on natural ordering */
72747c6ae99SBarry Smith     ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr);
728aa219208SBarry Smith     ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr);
72947c6ae99SBarry Smith     ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr);
730aa219208SBarry Smith     ierr = DMDAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr);
731aa219208SBarry Smith     ierr = DMDAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr);
7323f3fd955SJed Brown     ierr = PetscObjectGetName((PetscObject)xin,&name);CHKERRQ(ierr);
7333f3fd955SJed Brown     ierr = PetscObjectSetName((PetscObject)natural,name);CHKERRQ(ierr);
734a261c58fSBarry Smith 
735a261c58fSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
736a261c58fSBarry Smith     if (format == PETSC_VIEWER_BINARY_MATLAB) {
737a261c58fSBarry Smith       /* temporarily remove viewer format so it won't trigger in the VecView() */
7386a9046bcSBarry Smith       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
739a261c58fSBarry Smith     }
740a261c58fSBarry Smith 
74147c6ae99SBarry Smith     ierr = VecView(natural,viewer);CHKERRQ(ierr);
742a261c58fSBarry Smith 
743a261c58fSBarry Smith     if (format == PETSC_VIEWER_BINARY_MATLAB) {
744a261c58fSBarry Smith       MPI_Comm    comm;
745a261c58fSBarry Smith       FILE        *info;
746a261c58fSBarry Smith       const char  *fieldname;
747da88d4d4SJed Brown       char        fieldbuf[256];
748a261c58fSBarry Smith       PetscInt    dim,ni,nj,nk,pi,pj,pk,dof,n;
749a261c58fSBarry Smith 
750a261c58fSBarry Smith       /* set the viewer format back into the viewer */
7516a9046bcSBarry Smith       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
752a261c58fSBarry Smith       ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
753a261c58fSBarry Smith       ierr = PetscViewerBinaryGetInfoPointer(viewer,&info);CHKERRQ(ierr);
754a261c58fSBarry Smith       ierr = DMDAGetInfo(da,&dim,&ni,&nj,&nk,&pi,&pj,&pk,&dof,0,0,0,0,0);CHKERRQ(ierr);
755da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n");CHKERRQ(ierr);
756da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#$$ tmp = PetscBinaryRead(fd); \n");CHKERRQ(ierr);
757da88d4d4SJed Brown       if (dim == 1) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d);\n",dof,ni);CHKERRQ(ierr); }
758da88d4d4SJed Brown       if (dim == 2) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d);\n",dof,ni,nj);CHKERRQ(ierr); }
759da88d4d4SJed Brown       if (dim == 3) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d,%d);\n",dof,ni,nj,nk);CHKERRQ(ierr); }
760a261c58fSBarry Smith 
761a261c58fSBarry Smith       for (n=0; n<dof; n++) {
762a261c58fSBarry Smith         ierr = DMDAGetFieldName(da,n,&fieldname);CHKERRQ(ierr);
763da88d4d4SJed Brown         if (!fieldname || !fieldname[0]) {
764da88d4d4SJed Brown           ierr = PetscSNPrintf(fieldbuf,sizeof fieldbuf,"field%D",n);CHKERRQ(ierr);
765da88d4d4SJed Brown           fieldname = fieldbuf;
766a261c58fSBarry Smith         }
767da88d4d4SJed Brown         if (dim == 1) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:))';\n",name,fieldname,n+1);CHKERRQ(ierr); }
768da88d4d4SJed Brown         if (dim == 2) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:,:))';\n",name,fieldname,n+1);CHKERRQ(ierr); }
769da88d4d4SJed 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);}
770a261c58fSBarry Smith       }
771da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#$$ clear tmp; \n");CHKERRQ(ierr);
772da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n");CHKERRQ(ierr);
773a261c58fSBarry Smith     }
774a261c58fSBarry Smith 
775fcfd50ebSBarry Smith     ierr = VecDestroy(&natural);CHKERRQ(ierr);
77647c6ae99SBarry Smith   }
77747c6ae99SBarry Smith   PetscFunctionReturn(0);
77847c6ae99SBarry Smith }
77947c6ae99SBarry Smith 
78047c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
78147c6ae99SBarry Smith #undef __FUNCT__
78247c6ae99SBarry Smith #define __FUNCT__ "VecLoad_HDF5_DA"
78347c6ae99SBarry Smith PetscErrorCode VecLoad_HDF5_DA(Vec xin, PetscViewer viewer)
78447c6ae99SBarry Smith {
7859a42bb27SBarry Smith   DM             da;
78647c6ae99SBarry Smith   PetscErrorCode ierr;
787ec7ae49cSHåkon Strandenes   hsize_t        dim,rdim;
788ec7ae49cSHåkon Strandenes   hsize_t        dims[6]={0},count[6]={0},offset[6]={0};
789ec7ae49cSHåkon Strandenes   PetscInt       dimension,timestep,dofInd;
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 */
796bfd264e7SBarry Smith   hid_t          file_id,group;
7977bcbaff4SBarry Smith   hid_t          scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
7989c7c4993SBarry Smith   DM_DA          *dd;
799ec7ae49cSHåkon Strandenes   PetscBool      dim2 = PETSC_FALSE;
80047c6ae99SBarry Smith 
80147c6ae99SBarry Smith   PetscFunctionBegin;
8027bcbaff4SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
8037bcbaff4SBarry Smith   scalartype = H5T_NATIVE_FLOAT;
8047bcbaff4SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
8057bcbaff4SBarry Smith #error "HDF5 output with 128 bit floats not supported."
8067bcbaff4SBarry Smith #else
8077bcbaff4SBarry Smith   scalartype = H5T_NATIVE_DOUBLE;
8087bcbaff4SBarry Smith #endif
8097bcbaff4SBarry Smith 
810bfd264e7SBarry Smith   ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr);
811ec7ae49cSHåkon Strandenes   ierr = PetscViewerHDF5GetTimestep(viewer, &timestep);CHKERRQ(ierr);
812ec7ae49cSHåkon Strandenes   ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr);
813c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
8149c7c4993SBarry Smith   dd   = (DM_DA*)da->data;
815c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(da, &dimension);CHKERRQ(ierr);
81647c6ae99SBarry Smith 
817ec7ae49cSHåkon Strandenes   /* Open dataset */
81847c6ae99SBarry Smith #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
819729ab6d9SBarry Smith   PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT));
82047c6ae99SBarry Smith #else
821729ab6d9SBarry Smith   PetscStackCallHDF5Return(dset_id,H5Dopen,(group, vecname));
82247c6ae99SBarry Smith #endif
82347c6ae99SBarry Smith 
824ec7ae49cSHåkon Strandenes   /* Retrieve the dataspace for the dataset */
825ec7ae49cSHåkon Strandenes   PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id));
826ec7ae49cSHåkon Strandenes   PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(filespace, dims, NULL));
827ec7ae49cSHåkon Strandenes 
828ec7ae49cSHåkon Strandenes   /* Expected dimension for holding the dof's */
82947c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
830ec7ae49cSHåkon Strandenes   dofInd = rdim-2;
831ec7ae49cSHåkon Strandenes #else
832ec7ae49cSHåkon Strandenes   dofInd = rdim-1;
83347c6ae99SBarry Smith #endif
834ec7ae49cSHåkon Strandenes 
835ec7ae49cSHåkon Strandenes   /* The expected number of dimensions, assuming basedimension2 = false */
836ec7ae49cSHåkon Strandenes   dim = dimension;
837ec7ae49cSHåkon Strandenes   if (dd->w > 1) ++dim;
838ec7ae49cSHåkon Strandenes   if (timestep >= 0) ++dim;
83947c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
840ec7ae49cSHåkon Strandenes   ++dim;
84147c6ae99SBarry Smith #endif
842ec7ae49cSHåkon Strandenes 
843ec7ae49cSHåkon Strandenes   /* In this case the input dataset have one extra, unexpected dimension. */
844ec7ae49cSHåkon Strandenes   if (rdim == dim+1) {
845ec7ae49cSHåkon Strandenes     /* In this case the block size unity */
846ec7ae49cSHåkon Strandenes     if (dd->w == 1 && dims[dofInd] == 1) dim2 = PETSC_TRUE;
847ec7ae49cSHåkon Strandenes 
848ec7ae49cSHåkon Strandenes     /* Special error message for the case where dof does not match the input file */
849ec7ae49cSHå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);
850ec7ae49cSHåkon Strandenes 
851ec7ae49cSHåkon Strandenes   /* Other cases where rdim != dim cannot be handled currently */
8526c4ed002SBarry 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);
853ec7ae49cSHåkon Strandenes 
854ec7ae49cSHåkon Strandenes   /* Set up the hyperslab size */
855ec7ae49cSHåkon Strandenes   dim = 0;
856ec7ae49cSHåkon Strandenes   if (timestep >= 0) {
857ec7ae49cSHåkon Strandenes     offset[dim] = timestep;
858ec7ae49cSHåkon Strandenes     count[dim] = 1;
859ec7ae49cSHåkon Strandenes     ++dim;
860ec7ae49cSHåkon Strandenes   }
861ec7ae49cSHåkon Strandenes   if (dimension == 3) {
862ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->zs,offset + dim);CHKERRQ(ierr);
863ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->ze - dd->zs,count + dim);CHKERRQ(ierr);
864ec7ae49cSHåkon Strandenes     ++dim;
865ec7ae49cSHåkon Strandenes   }
866ec7ae49cSHåkon Strandenes   if (dimension > 1) {
867ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->ys,offset + dim);CHKERRQ(ierr);
868ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->ye - dd->ys,count + dim);CHKERRQ(ierr);
869ec7ae49cSHåkon Strandenes     ++dim;
870ec7ae49cSHåkon Strandenes   }
871ec7ae49cSHåkon Strandenes   ierr = PetscHDF5IntCast(dd->xs/dd->w,offset + dim);CHKERRQ(ierr);
872ec7ae49cSHåkon Strandenes   ierr = PetscHDF5IntCast((dd->xe - dd->xs)/dd->w,count + dim);CHKERRQ(ierr);
873ec7ae49cSHåkon Strandenes   ++dim;
874ec7ae49cSHåkon Strandenes   if (dd->w > 1 || dim2) {
875ec7ae49cSHåkon Strandenes     offset[dim] = 0;
876ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->w,count + dim);CHKERRQ(ierr);
877ec7ae49cSHåkon Strandenes     ++dim;
878ec7ae49cSHåkon Strandenes   }
879ec7ae49cSHåkon Strandenes #if defined(PETSC_USE_COMPLEX)
880ec7ae49cSHåkon Strandenes   offset[dim] = 0;
881ec7ae49cSHåkon Strandenes   count[dim] = 2;
882ec7ae49cSHåkon Strandenes   ++dim;
883ec7ae49cSHåkon Strandenes #endif
884ec7ae49cSHåkon Strandenes 
885ec7ae49cSHåkon Strandenes   /* Create the memory and filespace */
886729ab6d9SBarry Smith   PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL));
887729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL));
88847c6ae99SBarry Smith 
88947c6ae99SBarry Smith   /* Create property list for collective dataset write */
890729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_DATASET_XFER));
89147c6ae99SBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
892729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_dxpl_mpio,(plist_id, H5FD_MPIO_COLLECTIVE));
89347c6ae99SBarry Smith #endif
894ec7ae49cSHåkon Strandenes   /* To read dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
89547c6ae99SBarry Smith 
89647c6ae99SBarry Smith   ierr   = VecGetArray(xin, &x);CHKERRQ(ierr);
897a1dbdba5SBarry Smith   PetscStackCallHDF5(H5Dread,(dset_id, scalartype, memspace, filespace, plist_id, x));
89847c6ae99SBarry Smith   ierr   = VecRestoreArray(xin, &x);CHKERRQ(ierr);
89947c6ae99SBarry Smith 
90047c6ae99SBarry Smith   /* Close/release resources */
901bfd264e7SBarry Smith   if (group != file_id) {
902729ab6d9SBarry Smith     PetscStackCallHDF5(H5Gclose,(group));
903bfd264e7SBarry Smith   }
904729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
905729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
906729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(memspace));
907729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dset_id));
90847c6ae99SBarry Smith   PetscFunctionReturn(0);
90947c6ae99SBarry Smith }
91047c6ae99SBarry Smith #endif
91147c6ae99SBarry Smith 
91247c6ae99SBarry Smith #undef __FUNCT__
91347c6ae99SBarry Smith #define __FUNCT__ "VecLoad_Binary_DA"
91447c6ae99SBarry Smith PetscErrorCode VecLoad_Binary_DA(Vec xin, PetscViewer viewer)
91547c6ae99SBarry Smith {
9169a42bb27SBarry Smith   DM             da;
91747c6ae99SBarry Smith   PetscErrorCode ierr;
91847c6ae99SBarry Smith   Vec            natural;
91947c6ae99SBarry Smith   const char     *prefix;
92047c6ae99SBarry Smith   PetscInt       bs;
92147c6ae99SBarry Smith   PetscBool      flag;
92247c6ae99SBarry Smith   DM_DA          *dd;
92347c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
92447c6ae99SBarry Smith   PetscBool isMPIIO;
92547c6ae99SBarry Smith #endif
92647c6ae99SBarry Smith 
92747c6ae99SBarry Smith   PetscFunctionBegin;
928c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
92947c6ae99SBarry Smith   dd   = (DM_DA*)da->data;
93047c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
931bc196f7cSDave May   ierr = PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO);CHKERRQ(ierr);
93247c6ae99SBarry Smith   if (isMPIIO) {
933aa219208SBarry Smith     ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_FALSE);CHKERRQ(ierr);
93447c6ae99SBarry Smith     PetscFunctionReturn(0);
93547c6ae99SBarry Smith   }
93647c6ae99SBarry Smith #endif
93747c6ae99SBarry Smith 
93847c6ae99SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr);
939aa219208SBarry Smith   ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr);
94047c6ae99SBarry Smith   ierr = PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name);CHKERRQ(ierr);
94147c6ae99SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr);
9428aea9937SBarry Smith   ierr = VecLoad(natural,viewer);CHKERRQ(ierr);
943aa219208SBarry Smith   ierr = DMDANaturalToGlobalBegin(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr);
944aa219208SBarry Smith   ierr = DMDANaturalToGlobalEnd(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr);
945fcfd50ebSBarry Smith   ierr = VecDestroy(&natural);CHKERRQ(ierr);
946aa219208SBarry Smith   ierr = PetscInfo(xin,"Loading vector from natural ordering into DMDA\n");CHKERRQ(ierr);
947c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(NULL,((PetscObject)xin)->prefix,"-vecload_block_size",&bs,&flag);CHKERRQ(ierr);
94847c6ae99SBarry Smith   if (flag && bs != dd->w) {
949aa219208SBarry Smith     ierr = PetscInfo2(xin,"Block size in file %D not equal to DMDA's dof %D\n",bs,dd->w);CHKERRQ(ierr);
95047c6ae99SBarry Smith   }
95147c6ae99SBarry Smith   PetscFunctionReturn(0);
95247c6ae99SBarry Smith }
95347c6ae99SBarry Smith 
95447c6ae99SBarry Smith #undef __FUNCT__
95547c6ae99SBarry Smith #define __FUNCT__ "VecLoad_Default_DA"
9567087cfbeSBarry Smith PetscErrorCode  VecLoad_Default_DA(Vec xin, PetscViewer viewer)
95747c6ae99SBarry Smith {
95847c6ae99SBarry Smith   PetscErrorCode ierr;
9599a42bb27SBarry Smith   DM             da;
96047c6ae99SBarry Smith   PetscBool      isbinary;
96147c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
96247c6ae99SBarry Smith   PetscBool ishdf5;
96347c6ae99SBarry Smith #endif
96447c6ae99SBarry Smith 
96547c6ae99SBarry Smith   PetscFunctionBegin;
966c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
967ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
96847c6ae99SBarry Smith 
96947c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
970251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr);
97147c6ae99SBarry Smith #endif
972251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
97347c6ae99SBarry Smith 
97447c6ae99SBarry Smith   if (isbinary) {
97547c6ae99SBarry Smith     ierr = VecLoad_Binary_DA(xin,viewer);CHKERRQ(ierr);
97647c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
97747c6ae99SBarry Smith   } else if (ishdf5) {
97847c6ae99SBarry Smith     ierr = VecLoad_HDF5_DA(xin,viewer);CHKERRQ(ierr);
97947c6ae99SBarry Smith #endif
980d34fcf5fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for vector loading", ((PetscObject)viewer)->type_name);
98147c6ae99SBarry Smith   PetscFunctionReturn(0);
98247c6ae99SBarry Smith }
983