xref: /petsc/src/dm/impls/da/gr2.c (revision d7dd068b66a8daa3a37c9e1556b34bd8def54922)
147c6ae99SBarry Smith 
247c6ae99SBarry Smith /*
3aa219208SBarry Smith    Plots vectors obtained with DMDACreate2d()
447c6ae99SBarry Smith */
547c6ae99SBarry Smith 
6af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>      /*I  "petscdmda.h"   I*/
78135c375SStefano Zampini #include <petsc/private/glvisvecimpl.h>
82d95be3dSVaclav Hapla #include <petsc/private/viewerimpl.h>
92d95be3dSVaclav Hapla #include <petsc/private/viewerhdf5impl.h>
109804daf3SBarry Smith #include <petscdraw.h>
1147c6ae99SBarry Smith 
1247c6ae99SBarry Smith /*
1347c6ae99SBarry Smith         The data that is passed into the graphics callback
1447c6ae99SBarry Smith */
1547c6ae99SBarry Smith typedef struct {
16e5ab1681SLisandro Dalcin   PetscMPIInt       rank;
17e5ab1681SLisandro Dalcin   PetscInt          m,n,dof,k;
18e5ab1681SLisandro Dalcin   PetscReal         xmin,xmax,ymin,ymax,min,max;
195edff71fSBarry Smith   const PetscScalar *xy,*v;
207fe7c8feSLisandro Dalcin   PetscBool         showaxis,showgrid;
21109c9344SBarry Smith   const char        *name0,*name1;
2247c6ae99SBarry Smith } ZoomCtx;
2347c6ae99SBarry Smith 
2447c6ae99SBarry Smith /*
2547c6ae99SBarry Smith        This does the drawing for one particular field
2647c6ae99SBarry Smith     in one particular set of coordinates. It is a callback
2747c6ae99SBarry Smith     called from PetscDrawZoom()
2847c6ae99SBarry Smith */
2947c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA2d_Zoom(PetscDraw draw,void *ctx)
3047c6ae99SBarry Smith {
3147c6ae99SBarry Smith   ZoomCtx           *zctx = (ZoomCtx*)ctx;
3247c6ae99SBarry Smith   PetscErrorCode    ierr;
33e5ab1681SLisandro Dalcin   PetscInt          m,n,i,j,k,dof,id,c1,c2,c3,c4;
34e5ab1681SLisandro Dalcin   PetscReal         min,max,x1,x2,x3,x4,y_1,y2,y3,y4;
35e5ab1681SLisandro Dalcin   const PetscScalar *xy,*v;
3647c6ae99SBarry Smith 
3747c6ae99SBarry Smith   PetscFunctionBegin;
3847c6ae99SBarry Smith   m    = zctx->m;
3947c6ae99SBarry Smith   n    = zctx->n;
40e5ab1681SLisandro Dalcin   dof  = zctx->dof;
4147c6ae99SBarry Smith   k    = zctx->k;
4247c6ae99SBarry Smith   xy   = zctx->xy;
43e5ab1681SLisandro Dalcin   v    = zctx->v;
4447c6ae99SBarry Smith   min  = zctx->min;
45f3f0eb19SBarry Smith   max  = zctx->max;
4647c6ae99SBarry Smith 
4747c6ae99SBarry Smith   /* PetscDraw the contour plot patch */
48e5ab1681SLisandro Dalcin   ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
4947c6ae99SBarry Smith   for (j=0; j<n-1; j++) {
5047c6ae99SBarry Smith     for (i=0; i<m-1; i++) {
510e4fe250SBarry Smith       id   = i+j*m;
520e4fe250SBarry Smith       x1   = PetscRealPart(xy[2*id]);
530e4fe250SBarry Smith       y_1  = PetscRealPart(xy[2*id+1]);
54e5ab1681SLisandro Dalcin       c1   = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max);
550e4fe250SBarry Smith 
560e4fe250SBarry Smith       id   = i+j*m+1;
570e4fe250SBarry Smith       x2   = PetscRealPart(xy[2*id]);
58a39a4c7dSLisandro Dalcin       y2   = PetscRealPart(xy[2*id+1]);
59e5ab1681SLisandro Dalcin       c2   = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max);
600e4fe250SBarry Smith 
610e4fe250SBarry Smith       id   = i+j*m+1+m;
62a39a4c7dSLisandro Dalcin       x3   = PetscRealPart(xy[2*id]);
630e4fe250SBarry Smith       y3   = PetscRealPart(xy[2*id+1]);
64e5ab1681SLisandro Dalcin       c3   = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max);
650e4fe250SBarry Smith 
660e4fe250SBarry Smith       id   = i+j*m+m;
67a39a4c7dSLisandro Dalcin       x4   = PetscRealPart(xy[2*id]);
68a39a4c7dSLisandro Dalcin       y4   = PetscRealPart(xy[2*id+1]);
69e5ab1681SLisandro Dalcin       c4   = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max);
70f3f0eb19SBarry Smith 
7147c6ae99SBarry Smith       ierr = PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);CHKERRQ(ierr);
7247c6ae99SBarry Smith       ierr = PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4);CHKERRQ(ierr);
7347c6ae99SBarry Smith       if (zctx->showgrid) {
7447c6ae99SBarry Smith         ierr = PetscDrawLine(draw,x1,y_1,x2,y2,PETSC_DRAW_BLACK);CHKERRQ(ierr);
7547c6ae99SBarry Smith         ierr = PetscDrawLine(draw,x2,y2,x3,y3,PETSC_DRAW_BLACK);CHKERRQ(ierr);
7647c6ae99SBarry Smith         ierr = PetscDrawLine(draw,x3,y3,x4,y4,PETSC_DRAW_BLACK);CHKERRQ(ierr);
7747c6ae99SBarry Smith         ierr = PetscDrawLine(draw,x4,y4,x1,y_1,PETSC_DRAW_BLACK);CHKERRQ(ierr);
7847c6ae99SBarry Smith       }
7947c6ae99SBarry Smith     }
8047c6ae99SBarry Smith   }
817fe7c8feSLisandro Dalcin   if (zctx->showaxis && !zctx->rank) {
82e5ab1681SLisandro Dalcin     if (zctx->name0 || zctx->name1) {
83109c9344SBarry Smith       PetscReal xl,yl,xr,yr,x,y;
84109c9344SBarry Smith       ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
85e5ab1681SLisandro Dalcin       x  = xl + .30*(xr - xl);
86109c9344SBarry Smith       xl = xl + .01*(xr - xl);
87e5ab1681SLisandro Dalcin       y  = yr - .30*(yr - yl);
88109c9344SBarry Smith       yl = yl + .01*(yr - yl);
89e5ab1681SLisandro Dalcin       if (zctx->name0) {ierr = PetscDrawString(draw,x,yl,PETSC_DRAW_BLACK,zctx->name0);CHKERRQ(ierr);}
90e5ab1681SLisandro Dalcin       if (zctx->name1) {ierr = PetscDrawStringVertical(draw,xl,y,PETSC_DRAW_BLACK,zctx->name1);CHKERRQ(ierr);}
91109c9344SBarry Smith     }
920e4fe250SBarry Smith     /*
930e4fe250SBarry Smith        Ideally we would use the PetscDrawAxis object to manage displaying the coordinate limits
940e4fe250SBarry Smith        but that may require some refactoring.
950e4fe250SBarry Smith     */
96e5ab1681SLisandro Dalcin     {
977fe7c8feSLisandro Dalcin       double xmin = (double)zctx->xmin, ymin = (double)zctx->ymin;
987fe7c8feSLisandro Dalcin       double xmax = (double)zctx->xmax, ymax = (double)zctx->ymax;
99e5ab1681SLisandro Dalcin       char   value[16]; size_t len; PetscReal w;
1007fe7c8feSLisandro Dalcin       ierr = PetscSNPrintf(value,16,"%0.2e",xmin);CHKERRQ(ierr);
101e5ab1681SLisandro Dalcin       ierr = PetscDrawString(draw,xmin,ymin - .05*(ymax - ymin),PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
1027fe7c8feSLisandro Dalcin       ierr = PetscSNPrintf(value,16,"%0.2e",xmax);CHKERRQ(ierr);
1030e4fe250SBarry Smith       ierr = PetscStrlen(value,&len);CHKERRQ(ierr);
1040298fd71SBarry Smith       ierr = PetscDrawStringGetSize(draw,&w,NULL);CHKERRQ(ierr);
105e5ab1681SLisandro Dalcin       ierr = PetscDrawString(draw,xmax - len*w,ymin - .05*(ymax - ymin),PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
1067fe7c8feSLisandro Dalcin       ierr = PetscSNPrintf(value,16,"%0.2e",ymin);CHKERRQ(ierr);
107e5ab1681SLisandro Dalcin       ierr = PetscDrawString(draw,xmin - .05*(xmax - xmin),ymin,PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
1087fe7c8feSLisandro Dalcin       ierr = PetscSNPrintf(value,16,"%0.2e",ymax);CHKERRQ(ierr);
109e5ab1681SLisandro Dalcin       ierr = PetscDrawString(draw,xmin - .05*(xmax - xmin),ymax,PETSC_DRAW_BLACK,value);CHKERRQ(ierr);
110e5ab1681SLisandro Dalcin     }
111e5ab1681SLisandro Dalcin   }
112e5ab1681SLisandro Dalcin   ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
11347c6ae99SBarry Smith   PetscFunctionReturn(0);
11447c6ae99SBarry Smith }
11547c6ae99SBarry Smith 
11647c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA2d(Vec xin,PetscViewer viewer)
11747c6ae99SBarry Smith {
1189a42bb27SBarry Smith   DM                 da,dac,dag;
11947c6ae99SBarry Smith   PetscErrorCode     ierr;
120a74ba6f7SBarry Smith   PetscInt           N,s,M,w,ncoors = 4;
12147c6ae99SBarry Smith   const PetscInt     *lx,*ly;
122e5ab1681SLisandro Dalcin   PetscReal          coors[4];
12347c6ae99SBarry Smith   PetscDraw          draw,popup;
12447c6ae99SBarry Smith   PetscBool          isnull,useports = PETSC_FALSE;
12547c6ae99SBarry Smith   MPI_Comm           comm;
12647c6ae99SBarry Smith   Vec                xlocal,xcoor,xcoorl;
127bff4a2f0SMatthew G. Knepley   DMBoundaryType     bx,by;
128aa219208SBarry Smith   DMDAStencilType    st;
12947c6ae99SBarry Smith   ZoomCtx            zctx;
1300298fd71SBarry Smith   PetscDrawViewPorts *ports = NULL;
13147c6ae99SBarry Smith   PetscViewerFormat  format;
13220d0051dSBarry Smith   PetscInt           *displayfields;
13367dd0837SBarry Smith   PetscInt           ndisplayfields,i,nbounds;
13467dd0837SBarry Smith   const PetscReal    *bounds;
13547c6ae99SBarry Smith 
13647c6ae99SBarry Smith   PetscFunctionBegin;
13747c6ae99SBarry Smith   zctx.showgrid = PETSC_FALSE;
1387fe7c8feSLisandro Dalcin   zctx.showaxis = PETSC_TRUE;
1398865f1eaSKarl Rupp 
14047c6ae99SBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1415b399a63SLisandro Dalcin   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
1425b399a63SLisandro Dalcin   if (isnull) PetscFunctionReturn(0);
1435b399a63SLisandro Dalcin 
14403193ff8SBarry Smith   ierr = PetscViewerDrawGetBounds(viewer,&nbounds,&bounds);CHKERRQ(ierr);
14547c6ae99SBarry Smith 
146c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
147ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
14847c6ae99SBarry Smith 
14947c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr);
150ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&zctx.rank);CHKERRMPI(ierr);
15147c6ae99SBarry Smith 
152ea78f98cSLisandro Dalcin   ierr = DMDAGetInfo(da,NULL,&M,&N,NULL,&zctx.m,&zctx.n,NULL,&w,&s,&bx,&by,NULL,&st);CHKERRQ(ierr);
1530298fd71SBarry Smith   ierr = DMDAGetOwnershipRanges(da,&lx,&ly,NULL);CHKERRQ(ierr);
15447c6ae99SBarry Smith 
15547c6ae99SBarry Smith   /*
15647c6ae99SBarry Smith      Obtain a sequential vector that is going to contain the local values plus ONE layer of
157aa219208SBarry Smith      ghosted values to draw the graphics from. We also need its corresponding DMDA (dac) that will
15847c6ae99SBarry Smith      update the local values pluse ONE layer of ghost values.
15947c6ae99SBarry Smith   */
16047c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)da,"GraphicsGhosted",(PetscObject*)&xlocal);CHKERRQ(ierr);
16147c6ae99SBarry Smith   if (!xlocal) {
162bff4a2f0SMatthew G. Knepley     if (bx !=  DM_BOUNDARY_NONE || by !=  DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) {
16347c6ae99SBarry Smith       /*
16447c6ae99SBarry Smith          if original da is not of stencil width one, or periodic or not a box stencil then
165aa219208SBarry Smith          create a special DMDA to handle one level of ghost points for graphics
16647c6ae99SBarry Smith       */
167bff4a2f0SMatthew 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);
168897f7067SBarry Smith       ierr = DMSetUp(dac);CHKERRQ(ierr);
169aa219208SBarry Smith       ierr = PetscInfo(da,"Creating auxilary DMDA for managing graphics ghost points\n");CHKERRQ(ierr);
17047c6ae99SBarry Smith     } else {
17147c6ae99SBarry Smith       /* otherwise we can use the da we already have */
17247c6ae99SBarry Smith       dac = da;
17347c6ae99SBarry Smith     }
17447c6ae99SBarry Smith     /* create local vector for holding ghosted values used in graphics */
175564755cdSBarry Smith     ierr = DMCreateLocalVector(dac,&xlocal);CHKERRQ(ierr);
17647c6ae99SBarry Smith     if (dac != da) {
177aa219208SBarry Smith       /* don't keep any public reference of this DMDA, is is only available through xlocal */
178f7923d8aSBarry Smith       ierr = PetscObjectDereference((PetscObject)dac);CHKERRQ(ierr);
17947c6ae99SBarry Smith     } else {
18047c6ae99SBarry Smith       /* remove association between xlocal and da, because below we compose in the opposite
18147c6ae99SBarry Smith          direction and if we left this connect we'd get a loop, so the objects could
18247c6ae99SBarry Smith          never be destroyed */
183c688c046SMatthew G Knepley       ierr = PetscObjectRemoveReference((PetscObject)xlocal,"__PETSc_dm");CHKERRQ(ierr);
18447c6ae99SBarry Smith     }
18547c6ae99SBarry Smith     ierr = PetscObjectCompose((PetscObject)da,"GraphicsGhosted",(PetscObject)xlocal);CHKERRQ(ierr);
18647c6ae99SBarry Smith     ierr = PetscObjectDereference((PetscObject)xlocal);CHKERRQ(ierr);
18747c6ae99SBarry Smith   } else {
188bff4a2f0SMatthew G. Knepley     if (bx !=  DM_BOUNDARY_NONE || by !=  DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) {
189c688c046SMatthew G Knepley       ierr = VecGetDM(xlocal, &dac);CHKERRQ(ierr);
190f7923d8aSBarry Smith     } else {
191f7923d8aSBarry Smith       dac = da;
19247c6ae99SBarry Smith     }
19347c6ae99SBarry Smith   }
19447c6ae99SBarry Smith 
19547c6ae99SBarry Smith   /*
19647c6ae99SBarry Smith       Get local (ghosted) values of vector
19747c6ae99SBarry Smith   */
1989a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr);
1999a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr);
2005edff71fSBarry Smith   ierr = VecGetArrayRead(xlocal,&zctx.v);CHKERRQ(ierr);
20147c6ae99SBarry Smith 
202832b7cebSLisandro Dalcin   /*
203832b7cebSLisandro Dalcin       Get coordinates of nodes
204832b7cebSLisandro Dalcin   */
2056636e97aSMatthew G Knepley   ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
20647c6ae99SBarry Smith   if (!xcoor) {
207aa219208SBarry Smith     ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,0.0);CHKERRQ(ierr);
2086636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
20947c6ae99SBarry Smith   }
21047c6ae99SBarry Smith 
21147c6ae99SBarry Smith   /*
21247c6ae99SBarry Smith       Determine the min and max coordinates in plot
21347c6ae99SBarry Smith   */
214e5ab1681SLisandro Dalcin   ierr = VecStrideMin(xcoor,0,NULL,&zctx.xmin);CHKERRQ(ierr);
215e5ab1681SLisandro Dalcin   ierr = VecStrideMax(xcoor,0,NULL,&zctx.xmax);CHKERRQ(ierr);
216e5ab1681SLisandro Dalcin   ierr = VecStrideMin(xcoor,1,NULL,&zctx.ymin);CHKERRQ(ierr);
217e5ab1681SLisandro Dalcin   ierr = VecStrideMax(xcoor,1,NULL,&zctx.ymax);CHKERRQ(ierr);
2187fe7c8feSLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-draw_contour_axis",&zctx.showaxis,NULL);CHKERRQ(ierr);
2197fe7c8feSLisandro Dalcin   if (zctx.showaxis) {
2207fe7c8feSLisandro Dalcin     coors[0] = zctx.xmin - .05*(zctx.xmax - zctx.xmin); coors[1] = zctx.ymin - .05*(zctx.ymax - zctx.ymin);
2217fe7c8feSLisandro Dalcin     coors[2] = zctx.xmax + .05*(zctx.xmax - zctx.xmin); coors[3] = zctx.ymax + .05*(zctx.ymax - zctx.ymin);
2227fe7c8feSLisandro Dalcin   } else {
2237fe7c8feSLisandro Dalcin     coors[0] = zctx.xmin; coors[1] = zctx.ymin; coors[2] = zctx.xmax; coors[3] = zctx.ymax;
2247fe7c8feSLisandro Dalcin   }
225c5929fdfSBarry Smith   ierr = PetscOptionsGetRealArray(NULL,NULL,"-draw_coordinates",coors,&ncoors,NULL);CHKERRQ(ierr);
226e5ab1681SLisandro Dalcin   ierr = PetscInfo4(da,"Preparing DMDA 2d contour plot coordinates %g %g %g %g\n",(double)coors[0],(double)coors[1],(double)coors[2],(double)coors[3]);CHKERRQ(ierr);
227a74ba6f7SBarry Smith 
22847c6ae99SBarry Smith   /*
229832b7cebSLisandro Dalcin       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 */
234bff4a2f0SMatthew 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);
235897f7067SBarry Smith     ierr = DMSetUp(dag);CHKERRQ(ierr);
236aa219208SBarry Smith     ierr = PetscInfo(dag,"Creating auxilary DMDA for managing graphics coordinates ghost points\n");CHKERRQ(ierr);
237564755cdSBarry Smith     ierr = DMCreateLocalVector(dag,&xcoorl);CHKERRQ(ierr);
23847c6ae99SBarry Smith     ierr = PetscObjectCompose((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject)xcoorl);CHKERRQ(ierr);
239f7923d8aSBarry Smith     ierr = PetscObjectDereference((PetscObject)dag);CHKERRQ(ierr);
24047c6ae99SBarry Smith     ierr = PetscObjectDereference((PetscObject)xcoorl);CHKERRQ(ierr);
24147c6ae99SBarry Smith   } else {
242c688c046SMatthew G Knepley     ierr = VecGetDM(xcoorl,&dag);CHKERRQ(ierr);
24347c6ae99SBarry Smith   }
2449a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
2459a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
2465edff71fSBarry Smith   ierr = VecGetArrayRead(xcoorl,&zctx.xy);CHKERRQ(ierr);
247832b7cebSLisandro Dalcin   ierr = DMDAGetCoordinateName(da,0,&zctx.name0);CHKERRQ(ierr);
248832b7cebSLisandro Dalcin   ierr = DMDAGetCoordinateName(da,1,&zctx.name1);CHKERRQ(ierr);
24947c6ae99SBarry Smith 
25047c6ae99SBarry Smith   /*
25147c6ae99SBarry Smith       Get information about size of area each processor must do graphics for
25247c6ae99SBarry Smith   */
253e5ab1681SLisandro Dalcin   ierr = DMDAGetInfo(dac,NULL,&M,&N,NULL,NULL,NULL,NULL,&zctx.dof,NULL,&bx,&by,NULL,NULL);CHKERRQ(ierr);
254e5ab1681SLisandro Dalcin   ierr = DMDAGetGhostCorners(dac,NULL,NULL,NULL,&zctx.m,&zctx.n,NULL);CHKERRQ(ierr);
255c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-draw_contour_grid",&zctx.showgrid,NULL);CHKERRQ(ierr);
2564e6118eeSBarry Smith 
2574e6118eeSBarry Smith   ierr = DMDASelectFields(da,&ndisplayfields,&displayfields);CHKERRQ(ierr);
25847c6ae99SBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
259c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-draw_ports",&useports,NULL);CHKERRQ(ierr);
260832b7cebSLisandro Dalcin   if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE;
261832b7cebSLisandro Dalcin   if (useports) {
262e5ab1681SLisandro Dalcin     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
2635b399a63SLisandro Dalcin     ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
2645b399a63SLisandro Dalcin     ierr = PetscDrawClear(draw);CHKERRQ(ierr);
26520d0051dSBarry Smith     ierr = PetscDrawViewPortsCreate(draw,ndisplayfields,&ports);CHKERRQ(ierr);
266e5ab1681SLisandro Dalcin   }
26720d0051dSBarry Smith 
268832b7cebSLisandro Dalcin   /*
269832b7cebSLisandro Dalcin       Loop over each field; drawing each in a different window
270832b7cebSLisandro Dalcin   */
27120d0051dSBarry Smith   for (i=0; i<ndisplayfields; i++) {
27220d0051dSBarry Smith     zctx.k = displayfields[i];
2735b399a63SLisandro Dalcin 
274e5ab1681SLisandro Dalcin     /* determine the min and max value in plot */
2750298fd71SBarry Smith     ierr = VecStrideMin(xin,zctx.k,NULL,&zctx.min);CHKERRQ(ierr);
2760298fd71SBarry Smith     ierr = VecStrideMax(xin,zctx.k,NULL,&zctx.max);CHKERRQ(ierr);
27767dd0837SBarry Smith     if (zctx.k < nbounds) {
278f3f0eb19SBarry Smith       zctx.min = bounds[2*zctx.k];
279f3f0eb19SBarry Smith       zctx.max = bounds[2*zctx.k+1];
28067dd0837SBarry Smith     }
28147c6ae99SBarry Smith     if (zctx.min == zctx.max) {
28247c6ae99SBarry Smith       zctx.min -= 1.e-12;
28347c6ae99SBarry Smith       zctx.max += 1.e-12;
28447c6ae99SBarry Smith     }
28557622a8eSBarry Smith     ierr = PetscInfo2(da,"DMDA 2d contour plot min %g max %g\n",(double)zctx.min,(double)zctx.max);CHKERRQ(ierr);
28647c6ae99SBarry Smith 
287832b7cebSLisandro Dalcin     if (useports) {
288832b7cebSLisandro Dalcin       ierr = PetscDrawViewPortsSet(ports,i);CHKERRQ(ierr);
289832b7cebSLisandro Dalcin     } else {
290832b7cebSLisandro Dalcin       const char *title;
291832b7cebSLisandro Dalcin       ierr = PetscViewerDrawGetDraw(viewer,i,&draw);CHKERRQ(ierr);
292832b7cebSLisandro Dalcin       ierr = DMDAGetFieldName(da,zctx.k,&title);CHKERRQ(ierr);
293832b7cebSLisandro Dalcin       if (title) {ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr);}
294832b7cebSLisandro Dalcin     }
295832b7cebSLisandro Dalcin 
29647c6ae99SBarry Smith     ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
29745f3bb6eSLisandro Dalcin     ierr = PetscDrawScalePopup(popup,zctx.min,zctx.max);CHKERRQ(ierr);
298e5ab1681SLisandro Dalcin     ierr = PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr);
29947c6ae99SBarry Smith     ierr = PetscDrawZoom(draw,VecView_MPI_Draw_DA2d_Zoom,&zctx);CHKERRQ(ierr);
300832b7cebSLisandro Dalcin     if (!useports) {ierr = PetscDrawSave(draw);CHKERRQ(ierr);}
30147c6ae99SBarry Smith   }
302832b7cebSLisandro Dalcin   if (useports) {
303832b7cebSLisandro Dalcin     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
304832b7cebSLisandro Dalcin     ierr = PetscDrawSave(draw);CHKERRQ(ierr);
305832b7cebSLisandro Dalcin   }
306832b7cebSLisandro Dalcin 
3076bf464f9SBarry Smith   ierr = PetscDrawViewPortsDestroy(ports);CHKERRQ(ierr);
308e5ab1681SLisandro Dalcin   ierr = PetscFree(displayfields);CHKERRQ(ierr);
3095edff71fSBarry Smith   ierr = VecRestoreArrayRead(xcoorl,&zctx.xy);CHKERRQ(ierr);
3105edff71fSBarry Smith   ierr = VecRestoreArrayRead(xlocal,&zctx.v);CHKERRQ(ierr);
31147c6ae99SBarry Smith   PetscFunctionReturn(0);
31247c6ae99SBarry Smith }
31347c6ae99SBarry Smith 
3140da24e51SJuha Jäykkä #if defined(PETSC_HAVE_HDF5)
315c73cfb54SMatthew G. Knepley static PetscErrorCode VecGetHDF5ChunkSize(DM_DA *da, Vec xin, PetscInt dimension, PetscInt timestep, hsize_t *chunkDims)
3160da24e51SJuha Jäykkä {
317d4190030SJed Brown   PetscMPIInt    comm_size;
318d4190030SJed Brown   PetscErrorCode ierr;
319561a82e7SJed Brown   hsize_t        chunk_size, target_size, dim;
320561a82e7SJed Brown   hsize_t        vec_size = sizeof(PetscScalar)*da->M*da->N*da->P*da->w;
321561a82e7SJed Brown   hsize_t        avg_local_vec_size,KiB = 1024,MiB = KiB*KiB,GiB = MiB*KiB,min_size = MiB;
322561a82e7SJed Brown   hsize_t        max_chunks = 64*KiB;                                              /* HDF5 internal limitation */
323561a82e7SJed Brown   hsize_t        max_chunk_size = 4*GiB;                                           /* HDF5 internal limitation */
32484179ffaSJed Brown   PetscInt       zslices=da->p, yslices=da->n, xslices=da->m;
3250da24e51SJuha Jäykkä 
326cb5c4f94SJuha Jäykkä   PetscFunctionBegin;
327ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)xin), &comm_size);CHKERRMPI(ierr);
3280da24e51SJuha 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 */
3290da24e51SJuha Jäykkä 
330794a961bSBarry 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)))));
3317d310018SBarry 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 */
3327d310018SBarry 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);
3330da24e51SJuha Jäykkä 
334cb5c4f94SJuha Jäykkä   /*
335cb5c4f94SJuha Jäykkä    if size/rank > max_chunk_size, we need radical measures: even going down to
336cb5c4f94SJuha Jäykkä    avg_local_vec_size is not enough, so we simply use chunk size of 4 GiB no matter
337cb5c4f94SJuha Jäykkä    what, composed in the most efficient way possible.
338cb5c4f94SJuha Jäykkä    N.B. this minimises the number of chunks, which may or may not be the optimal
339cb5c4f94SJuha Jäykkä    solution. In a BG, for example, the optimal solution is probably to make # chunks = #
340cb5c4f94SJuha Jäykkä    IO nodes involved, but this author has no access to a BG to figure out how to
341cb5c4f94SJuha Jäykkä    reliably find the right number. And even then it may or may not be enough.
342cb5c4f94SJuha Jäykkä    */
3430da24e51SJuha Jäykkä   if (avg_local_vec_size > max_chunk_size) {
344cb5c4f94SJuha Jäykkä     /* check if we can just split local z-axis: is that enough? */
34584179ffaSJed Brown     zslices = (PetscInt)ceil(vec_size*1.0/(da->p*max_chunk_size))*zslices;
3460da24e51SJuha Jäykkä     if (zslices > da->P) {
347cb5c4f94SJuha Jäykkä       /* lattice is too large in xy-directions, splitting z only is not enough */
3480da24e51SJuha Jäykkä       zslices = da->P;
34984179ffaSJed Brown       yslices= (PetscInt)ceil(vec_size*1.0/(zslices*da->n*max_chunk_size))*yslices;
3500da24e51SJuha Jäykkä       if (yslices > da->N) {
351cb5c4f94SJuha Jäykkä         /* lattice is too large in x-direction, splitting along z, y is not enough */
3520da24e51SJuha Jäykkä         yslices = da->N;
35384179ffaSJed Brown         xslices= (PetscInt)ceil(vec_size*1.0/(zslices*yslices*da->m*max_chunk_size))*xslices;
3540da24e51SJuha Jäykkä       }
3550da24e51SJuha Jäykkä     }
3560da24e51SJuha Jäykkä     dim = 0;
3570da24e51SJuha Jäykkä     if (timestep >= 0) {
3580da24e51SJuha Jäykkä       ++dim;
3590da24e51SJuha Jäykkä     }
360cb5c4f94SJuha Jäykkä     /* prefer to split z-axis, even down to planar slices */
361c73cfb54SMatthew G. Knepley     if (dimension == 3) {
362cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->P/zslices;
363cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->N/yslices;
364cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->M/xslices;
3650da24e51SJuha Jäykkä     } else {
366cb5c4f94SJuha Jäykkä       /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */
367cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->N/yslices;
368cb5c4f94SJuha Jäykkä       chunkDims[dim++] = (hsize_t) da->M/xslices;
3690da24e51SJuha Jäykkä     }
3700da24e51SJuha 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);
3710da24e51SJuha Jäykkä   } else {
3720da24e51SJuha Jäykkä     if (target_size < chunk_size) {
373cb5c4f94SJuha Jäykkä       /* only change the defaults if target_size < chunk_size */
3740da24e51SJuha Jäykkä       dim = 0;
3750da24e51SJuha Jäykkä       if (timestep >= 0) {
3760da24e51SJuha Jäykkä         ++dim;
3770da24e51SJuha Jäykkä       }
378cb5c4f94SJuha Jäykkä       /* prefer to split z-axis, even down to planar slices */
379c73cfb54SMatthew G. Knepley       if (dimension == 3) {
380cb5c4f94SJuha Jäykkä         /* try splitting the z-axis to core-size bits, i.e. divide chunk size by # comm_size in z-direction */
3810da24e51SJuha Jäykkä         if (target_size >= chunk_size/da->p) {
382cb5c4f94SJuha Jäykkä           /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */
3830da24e51SJuha Jäykkä           chunkDims[dim] = (hsize_t) ceil(da->P*1.0/da->p);
3840da24e51SJuha Jäykkä         } else {
385cb5c4f94SJuha Jäykkä           /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be
386cb5c4f94SJuha Jäykkä            radical and let everyone write all they've got */
3870da24e51SJuha Jäykkä           chunkDims[dim++] = (hsize_t) ceil(da->P*1.0/da->p);
3880da24e51SJuha Jäykkä           chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n);
3890da24e51SJuha Jäykkä           chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m);
3900da24e51SJuha Jäykkä         }
3910da24e51SJuha Jäykkä       } else {
392cb5c4f94SJuha Jäykkä         /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */
3930da24e51SJuha Jäykkä         if (target_size >= chunk_size/da->n) {
394cb5c4f94SJuha Jäykkä           /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */
3950da24e51SJuha Jäykkä           chunkDims[dim] = (hsize_t) ceil(da->N*1.0/da->n);
3960da24e51SJuha Jäykkä         } else {
397cb5c4f94SJuha Jäykkä           /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be
398cb5c4f94SJuha Jäykkä            radical and let everyone write all they've got */
3990da24e51SJuha Jäykkä           chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n);
4000da24e51SJuha Jäykkä           chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m);
4010da24e51SJuha Jäykkä         }
4020da24e51SJuha Jäykkä 
4030da24e51SJuha Jäykkä       }
4040da24e51SJuha 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);
4050da24e51SJuha Jäykkä     } else {
406cb5c4f94SJuha Jäykkä       /* precomputed chunks are fine, we don't need to do anything */
4070da24e51SJuha Jäykkä     }
4080da24e51SJuha Jäykkä   }
4090da24e51SJuha Jäykkä   PetscFunctionReturn(0);
4100da24e51SJuha Jäykkä }
4110da24e51SJuha Jäykkä #endif
41247c6ae99SBarry Smith 
41347c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
41447c6ae99SBarry Smith PetscErrorCode VecView_MPI_HDF5_DA(Vec xin,PetscViewer viewer)
41547c6ae99SBarry Smith {
4162d95be3dSVaclav Hapla   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*) viewer->data;
4179b2a5a86SJed Brown   DM                dm;
4189b2a5a86SJed Brown   DM_DA             *da;
41947c6ae99SBarry Smith   hid_t             filespace;  /* file dataspace identifier */
4208e2ae6d7SMichael Kraus   hid_t             chunkspace; /* chunk dataset property identifier */
42147c6ae99SBarry Smith   hid_t             dset_id;    /* dataset identifier */
42247c6ae99SBarry Smith   hid_t             memspace;   /* memory dataspace identifier */
42347c6ae99SBarry Smith   hid_t             file_id;
42415214e8eSMatthew G Knepley   hid_t             group;
4259a0502c6SHåkon Strandenes   hid_t             memscalartype; /* scalar type for mem (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
4269a0502c6SHåkon Strandenes   hid_t             filescalartype; /* scalar type for file (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
427d9a4edebSJed Brown   hsize_t           dim;
4280da24e51SJuha 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  */
429*d7dd068bSVaclav Hapla   PetscBool         timestepping=PETSC_FALSE, dim2, spoutput;
430*d7dd068bSVaclav Hapla   PetscInt          timestep=PETSC_MIN_INT, dimension;
4315edff71fSBarry Smith   const PetscScalar *x;
4328e2ae6d7SMichael Kraus   const char        *vecname;
43315214e8eSMatthew G Knepley   PetscErrorCode    ierr;
43447c6ae99SBarry Smith 
43547c6ae99SBarry Smith   PetscFunctionBegin;
43615214e8eSMatthew G Knepley   ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr);
437*d7dd068bSVaclav Hapla   ierr = PetscViewerHDF5IsTimestepping(viewer, &timestepping);CHKERRQ(ierr);
438*d7dd068bSVaclav Hapla   if (timestepping) {
4398e2ae6d7SMichael Kraus     ierr = PetscViewerHDF5GetTimestep(viewer, &timestep);CHKERRQ(ierr);
440*d7dd068bSVaclav Hapla   }
44182ea9b62SBarry Smith   ierr = PetscViewerHDF5GetBaseDimension2(viewer,&dim2);CHKERRQ(ierr);
4429a0502c6SHåkon Strandenes   ierr = PetscViewerHDF5GetSPOutput(viewer,&spoutput);CHKERRQ(ierr);
44315214e8eSMatthew G Knepley 
444c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&dm);CHKERRQ(ierr);
445ce94432eSBarry Smith   if (!dm) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
4469b2a5a86SJed Brown   da = (DM_DA*)dm->data;
447c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dimension);CHKERRQ(ierr);
44847c6ae99SBarry Smith 
4498e2ae6d7SMichael Kraus   /* Create the dataspace for the dataset.
4508e2ae6d7SMichael Kraus    *
4518e2ae6d7SMichael Kraus    * dims - holds the current dimensions of the dataset
4528e2ae6d7SMichael Kraus    *
4538e2ae6d7SMichael Kraus    * maxDims - holds the maximum dimensions of the dataset (unlimited
4548e2ae6d7SMichael Kraus    * for the number of time steps with the current dimensions for the
4558e2ae6d7SMichael Kraus    * other dimensions; so only additional time steps can be added).
4568e2ae6d7SMichael Kraus    *
4578e2ae6d7SMichael Kraus    * chunkDims - holds the size of a single time step (required to
4588e2ae6d7SMichael Kraus    * permit extending dataset).
4598e2ae6d7SMichael Kraus    */
4608e2ae6d7SMichael Kraus   dim = 0;
4618e2ae6d7SMichael Kraus   if (timestep >= 0) {
4628e2ae6d7SMichael Kraus     dims[dim]      = timestep+1;
4638e2ae6d7SMichael Kraus     maxDims[dim]   = H5S_UNLIMITED;
4648e2ae6d7SMichael Kraus     chunkDims[dim] = 1;
4658e2ae6d7SMichael Kraus     ++dim;
4668e2ae6d7SMichael Kraus   }
467c73cfb54SMatthew G. Knepley   if (dimension == 3) {
468acba2ac6SBarry Smith     ierr           = PetscHDF5IntCast(da->P,dims+dim);CHKERRQ(ierr);
4698e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4708e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4718e2ae6d7SMichael Kraus     ++dim;
4728e2ae6d7SMichael Kraus   }
473c73cfb54SMatthew G. Knepley   if (dimension > 1) {
474acba2ac6SBarry Smith     ierr           = PetscHDF5IntCast(da->N,dims+dim);CHKERRQ(ierr);
4758e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4768e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4778e2ae6d7SMichael Kraus     ++dim;
4788e2ae6d7SMichael Kraus   }
479acba2ac6SBarry Smith   ierr           = PetscHDF5IntCast(da->M,dims+dim);CHKERRQ(ierr);
4808e2ae6d7SMichael Kraus   maxDims[dim]   = dims[dim];
4818e2ae6d7SMichael Kraus   chunkDims[dim] = dims[dim];
4828e2ae6d7SMichael Kraus   ++dim;
48382ea9b62SBarry Smith   if (da->w > 1 || dim2) {
484acba2ac6SBarry Smith     ierr           = PetscHDF5IntCast(da->w,dims+dim);CHKERRQ(ierr);
4858e2ae6d7SMichael Kraus     maxDims[dim]   = dims[dim];
4868e2ae6d7SMichael Kraus     chunkDims[dim] = dims[dim];
4878e2ae6d7SMichael Kraus     ++dim;
4888e2ae6d7SMichael Kraus   }
48947c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
4908e2ae6d7SMichael Kraus   dims[dim]      = 2;
4918e2ae6d7SMichael Kraus   maxDims[dim]   = dims[dim];
4928e2ae6d7SMichael Kraus   chunkDims[dim] = dims[dim];
4938e2ae6d7SMichael Kraus   ++dim;
49447c6ae99SBarry Smith #endif
4950da24e51SJuha Jäykkä 
496c73cfb54SMatthew G. Knepley   ierr = VecGetHDF5ChunkSize(da, xin, dimension, timestep, chunkDims);CHKERRQ(ierr);
4970da24e51SJuha Jäykkä 
498729ab6d9SBarry Smith   PetscStackCallHDF5Return(filespace,H5Screate_simple,(dim, dims, maxDims));
49947c6ae99SBarry Smith 
50015214e8eSMatthew G Knepley #if defined(PETSC_USE_REAL_SINGLE)
5019a0502c6SHåkon Strandenes   memscalartype = H5T_NATIVE_FLOAT;
5029a0502c6SHåkon Strandenes   filescalartype = H5T_NATIVE_FLOAT;
50315214e8eSMatthew G Knepley #elif defined(PETSC_USE_REAL___FLOAT128)
50415214e8eSMatthew G Knepley #error "HDF5 output with 128 bit floats not supported."
505570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
506570b7f6dSBarry Smith #error "HDF5 output with 16 bit floats not supported."
50715214e8eSMatthew G Knepley #else
5089a0502c6SHåkon Strandenes   memscalartype = H5T_NATIVE_DOUBLE;
5099a0502c6SHåkon Strandenes   if (spoutput == PETSC_TRUE) filescalartype = H5T_NATIVE_FLOAT;
5109a0502c6SHåkon Strandenes   else filescalartype = H5T_NATIVE_DOUBLE;
51115214e8eSMatthew G Knepley #endif
51215214e8eSMatthew G Knepley 
51347c6ae99SBarry Smith   /* Create the dataset with default properties and close filespace */
51447c6ae99SBarry Smith   ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr);
51515214e8eSMatthew G Knepley   if (!H5Lexists(group, vecname, H5P_DEFAULT)) {
5168e2ae6d7SMichael Kraus     /* Create chunk */
517729ab6d9SBarry Smith     PetscStackCallHDF5Return(chunkspace,H5Pcreate,(H5P_DATASET_CREATE));
518729ab6d9SBarry Smith     PetscStackCallHDF5(H5Pset_chunk,(chunkspace, dim, chunkDims));
5198e2ae6d7SMichael Kraus 
5209a0502c6SHåkon Strandenes     PetscStackCallHDF5Return(dset_id,H5Dcreate2,(group, vecname, filescalartype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT));
52115214e8eSMatthew G Knepley   } else {
522729ab6d9SBarry Smith     PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT));
523729ab6d9SBarry Smith     PetscStackCallHDF5(H5Dset_extent,(dset_id, dims));
52415214e8eSMatthew G Knepley   }
525729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
52647c6ae99SBarry Smith 
52747c6ae99SBarry Smith   /* Each process defines a dataset and writes it to the hyperslab in the file */
5288e2ae6d7SMichael Kraus   dim = 0;
5298e2ae6d7SMichael Kraus   if (timestep >= 0) {
5308e2ae6d7SMichael Kraus     offset[dim] = timestep;
5318e2ae6d7SMichael Kraus     ++dim;
5328e2ae6d7SMichael Kraus   }
533c73cfb54SMatthew G. Knepley   if (dimension == 3) {ierr = PetscHDF5IntCast(da->zs,offset + dim++);CHKERRQ(ierr);}
534c73cfb54SMatthew G. Knepley   if (dimension > 1)  {ierr = PetscHDF5IntCast(da->ys,offset + dim++);CHKERRQ(ierr);}
535acba2ac6SBarry Smith   ierr = PetscHDF5IntCast(da->xs/da->w,offset + dim++);CHKERRQ(ierr);
53682ea9b62SBarry Smith   if (da->w > 1 || dim2) offset[dim++] = 0;
53747c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5388e2ae6d7SMichael Kraus   offset[dim++] = 0;
53947c6ae99SBarry Smith #endif
5408e2ae6d7SMichael Kraus   dim = 0;
5418e2ae6d7SMichael Kraus   if (timestep >= 0) {
5428e2ae6d7SMichael Kraus     count[dim] = 1;
5438e2ae6d7SMichael Kraus     ++dim;
5448e2ae6d7SMichael Kraus   }
545c73cfb54SMatthew G. Knepley   if (dimension == 3) {ierr = PetscHDF5IntCast(da->ze - da->zs,count + dim++);CHKERRQ(ierr);}
546c73cfb54SMatthew G. Knepley   if (dimension > 1)  {ierr = PetscHDF5IntCast(da->ye - da->ys,count + dim++);CHKERRQ(ierr);}
547acba2ac6SBarry Smith   ierr = PetscHDF5IntCast((da->xe - da->xs)/da->w,count + dim++);CHKERRQ(ierr);
54882ea9b62SBarry Smith   if (da->w > 1 || dim2) {ierr = PetscHDF5IntCast(da->w,count + dim++);CHKERRQ(ierr);}
54947c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
5508e2ae6d7SMichael Kraus   count[dim++] = 2;
55147c6ae99SBarry Smith #endif
552729ab6d9SBarry Smith   PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL));
553729ab6d9SBarry Smith   PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id));
554729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL));
55547c6ae99SBarry Smith 
5565edff71fSBarry Smith   ierr   = VecGetArrayRead(xin, &x);CHKERRQ(ierr);
5572d95be3dSVaclav Hapla   PetscStackCallHDF5(H5Dwrite,(dset_id, memscalartype, memspace, filespace, hdf5->dxpl_id, x));
558729ab6d9SBarry Smith   PetscStackCallHDF5(H5Fflush,(file_id, H5F_SCOPE_GLOBAL));
5595edff71fSBarry Smith   ierr   = VecRestoreArrayRead(xin, &x);CHKERRQ(ierr);
56047c6ae99SBarry Smith 
561e0552f36Ssajid__ali   #if defined(PETSC_USE_COMPLEX)
562e0552f36Ssajid__ali   {
563e0552f36Ssajid__ali     PetscBool tru = PETSC_TRUE;
564e0552f36Ssajid__ali     ierr = PetscViewerHDF5WriteObjectAttribute(viewer,(PetscObject)xin,"complex",PETSC_BOOL,&tru);CHKERRQ(ierr);
565e0552f36Ssajid__ali   }
566e0552f36Ssajid__ali   #endif
567*d7dd068bSVaclav Hapla   ierr = PetscViewerHDF5WriteObjectAttribute(viewer,(PetscObject)xin,"timestepping",PETSC_BOOL,&timestepping);CHKERRQ(ierr);
568e0552f36Ssajid__ali 
56947c6ae99SBarry Smith   /* Close/release resources */
57015214e8eSMatthew G Knepley   if (group != file_id) {
571729ab6d9SBarry Smith     PetscStackCallHDF5(H5Gclose,(group));
57215214e8eSMatthew G Knepley   }
573729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
574729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(memspace));
575729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dset_id));
57647c6ae99SBarry Smith   ierr   = PetscInfo1(xin,"Wrote Vec object with name %s\n",vecname);CHKERRQ(ierr);
57747c6ae99SBarry Smith   PetscFunctionReturn(0);
57847c6ae99SBarry Smith }
57947c6ae99SBarry Smith #endif
58047c6ae99SBarry Smith 
58109573ac7SBarry Smith extern PetscErrorCode VecView_MPI_Draw_DA1d(Vec,PetscViewer);
58247c6ae99SBarry Smith 
58347c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
584aa219208SBarry Smith static PetscErrorCode DMDAArrayMPIIO(DM da,PetscViewer viewer,Vec xin,PetscBool write)
58547c6ae99SBarry Smith {
58647c6ae99SBarry Smith   PetscErrorCode    ierr;
58747c6ae99SBarry Smith   MPI_File          mfdes;
58847c6ae99SBarry Smith   PetscMPIInt       gsizes[4],lsizes[4],lstarts[4],asiz,dof;
58947c6ae99SBarry Smith   MPI_Datatype      view;
59047c6ae99SBarry Smith   const PetscScalar *array;
59147c6ae99SBarry Smith   MPI_Offset        off;
59247c6ae99SBarry Smith   MPI_Aint          ub,ul;
59347c6ae99SBarry Smith   PetscInt          type,rows,vecrows,tr[2];
59447c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
5950c169764Sdmay   PetscBool         skipheader;
59647c6ae99SBarry Smith 
59747c6ae99SBarry Smith   PetscFunctionBegin;
59847c6ae99SBarry Smith   ierr = VecGetSize(xin,&vecrows);CHKERRQ(ierr);
5990c169764Sdmay   ierr = PetscViewerBinaryGetSkipHeader(viewer,&skipheader);CHKERRQ(ierr);
60047c6ae99SBarry Smith   if (!write) {
60147c6ae99SBarry Smith     /* Read vector header. */
6020c169764Sdmay     if (!skipheader) {
603060da220SMatthew G. Knepley       ierr = PetscViewerBinaryRead(viewer,tr,2,NULL,PETSC_INT);CHKERRQ(ierr);
60447c6ae99SBarry Smith       type = tr[0];
60547c6ae99SBarry Smith       rows = tr[1];
606ce94432eSBarry Smith       if (type != VEC_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Not vector next in file");
607ce94432eSBarry Smith       if (rows != vecrows) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Vector in file not same size as DMDA vector");
6080c169764Sdmay     }
60947c6ae99SBarry Smith   } else {
61047c6ae99SBarry Smith     tr[0] = VEC_FILE_CLASSID;
61147c6ae99SBarry Smith     tr[1] = vecrows;
6120c169764Sdmay     if (!skipheader) {
613f253e43cSLisandro Dalcin       ierr  = PetscViewerBinaryWrite(viewer,tr,2,PETSC_INT);CHKERRQ(ierr);
61447c6ae99SBarry Smith     }
6150c169764Sdmay   }
61647c6ae99SBarry Smith 
6174dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->w,&dof);CHKERRQ(ierr);
6184dc2109aSBarry Smith   gsizes[0]  = dof;
6194dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->M,gsizes+1);CHKERRQ(ierr);
6204dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->N,gsizes+2);CHKERRQ(ierr);
621334634e2SJed Brown   ierr       = PetscMPIIntCast(dd->P,gsizes+3);CHKERRQ(ierr);
6224dc2109aSBarry Smith   lsizes[0]  = dof;
6234dc2109aSBarry Smith   ierr       = PetscMPIIntCast((dd->xe-dd->xs)/dof,lsizes+1);CHKERRQ(ierr);
6244dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->ye-dd->ys,lsizes+2);CHKERRQ(ierr);
6254dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->ze-dd->zs,lsizes+3);CHKERRQ(ierr);
6264dc2109aSBarry Smith   lstarts[0] = 0;
6274dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->xs/dof,lstarts+1);CHKERRQ(ierr);
6284dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->ys,lstarts+2);CHKERRQ(ierr);
6294dc2109aSBarry Smith   ierr       = PetscMPIIntCast(dd->zs,lstarts+3);CHKERRQ(ierr);
63055b25c41SPierre Jolivet   ierr       = MPI_Type_create_subarray(da->dim+1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view);CHKERRMPI(ierr);
63155b25c41SPierre Jolivet   ierr       = MPI_Type_commit(&view);CHKERRMPI(ierr);
63247c6ae99SBarry Smith 
63347c6ae99SBarry Smith   ierr = PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);CHKERRQ(ierr);
63447c6ae99SBarry Smith   ierr = PetscViewerBinaryGetMPIIOOffset(viewer,&off);CHKERRQ(ierr);
635ffc4695bSBarry Smith   ierr = MPI_File_set_view(mfdes,off,MPIU_SCALAR,view,(char*)"native",MPI_INFO_NULL);CHKERRMPI(ierr);
63647c6ae99SBarry Smith   ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr);
63747c6ae99SBarry Smith   asiz = lsizes[1]*(lsizes[2] > 0 ? lsizes[2] : 1)*(lsizes[3] > 0 ? lsizes[3] : 1)*dof;
63847c6ae99SBarry Smith   if (write) {
63947c6ae99SBarry Smith     ierr = MPIU_File_write_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr);
64047c6ae99SBarry Smith   } else {
64147c6ae99SBarry Smith     ierr = MPIU_File_read_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr);
64247c6ae99SBarry Smith   }
643ffc4695bSBarry Smith   ierr = MPI_Type_get_extent(view,&ul,&ub);CHKERRMPI(ierr);
64447c6ae99SBarry Smith   ierr = PetscViewerBinaryAddMPIIOOffset(viewer,ub);CHKERRQ(ierr);
64547c6ae99SBarry Smith   ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr);
646ffc4695bSBarry Smith   ierr = MPI_Type_free(&view);CHKERRMPI(ierr);
64747c6ae99SBarry Smith   PetscFunctionReturn(0);
64847c6ae99SBarry Smith }
64947c6ae99SBarry Smith #endif
65047c6ae99SBarry Smith 
6517087cfbeSBarry Smith PetscErrorCode  VecView_MPI_DA(Vec xin,PetscViewer viewer)
65247c6ae99SBarry Smith {
6539a42bb27SBarry Smith   DM                da;
65447c6ae99SBarry Smith   PetscErrorCode    ierr;
65547c6ae99SBarry Smith   PetscInt          dim;
65647c6ae99SBarry Smith   Vec               natural;
6578135c375SStefano Zampini   PetscBool         isdraw,isvtk,isglvis;
65847c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
65947c6ae99SBarry Smith   PetscBool         ishdf5;
66047c6ae99SBarry Smith #endif
6613f3fd955SJed Brown   const char        *prefix,*name;
662a261c58fSBarry Smith   PetscViewerFormat format;
66347c6ae99SBarry Smith 
66447c6ae99SBarry Smith   PetscFunctionBegin;
665c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
666ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
667251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
668251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr);
66947c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
670251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr);
67147c6ae99SBarry Smith #endif
6728135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr);
67347c6ae99SBarry Smith   if (isdraw) {
674ea78f98cSLisandro Dalcin     ierr = DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
67547c6ae99SBarry Smith     if (dim == 1) {
67647c6ae99SBarry Smith       ierr = VecView_MPI_Draw_DA1d(xin,viewer);CHKERRQ(ierr);
67747c6ae99SBarry Smith     } else if (dim == 2) {
67847c6ae99SBarry Smith       ierr = VecView_MPI_Draw_DA2d(xin,viewer);CHKERRQ(ierr);
679ce94432eSBarry Smith     } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot graphically view vector associated with this dimensional DMDA %D",dim);
68001d7c5c3SPatrick Sanan   } else if (isvtk) {           /* Duplicate the Vec */
6814061b8bfSJed Brown     Vec Y;
6824061b8bfSJed Brown     ierr = VecDuplicate(xin,&Y);CHKERRQ(ierr);
683b51b94faSJed Brown     if (((PetscObject)xin)->name) {
684b51b94faSJed Brown       /* If xin was named, copy the name over to Y. The duplicate names are safe because nobody else will ever see Y. */
685b51b94faSJed Brown       ierr = PetscObjectSetName((PetscObject)Y,((PetscObject)xin)->name);CHKERRQ(ierr);
686b51b94faSJed Brown     }
6874061b8bfSJed Brown     ierr = VecCopy(xin,Y);CHKERRQ(ierr);
688a8f87f1dSPatrick Sanan     {
689a8f87f1dSPatrick Sanan       PetscObject dmvtk;
690a8f87f1dSPatrick Sanan       PetscBool   compatible,compatibleSet;
691a8f87f1dSPatrick Sanan       ierr = PetscViewerVTKGetDM(viewer,&dmvtk);CHKERRQ(ierr);
692a8f87f1dSPatrick Sanan       if (dmvtk) {
693a8f87f1dSPatrick Sanan         PetscValidHeaderSpecific((DM)dmvtk,DM_CLASSID,-1);
694a8f87f1dSPatrick Sanan         ierr = DMGetCompatibility(da,(DM)dmvtk,&compatible,&compatibleSet);CHKERRQ(ierr);
695a8f87f1dSPatrick Sanan         if (!compatibleSet || !compatible) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"Cannot confirm compatibility of DMs associated with Vecs viewed in the same VTK file. Check that grids are the same.");
696a8f87f1dSPatrick Sanan       }
697e630c359SToby Isaac       ierr = PetscViewerVTKAddField(viewer,(PetscObject)da,DMDAVTKWriteAll,PETSC_DEFAULT,PETSC_VTK_POINT_FIELD,PETSC_FALSE,(PetscObject)Y);CHKERRQ(ierr);
698a8f87f1dSPatrick Sanan     }
69947c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
70047c6ae99SBarry Smith   } else if (ishdf5) {
70147c6ae99SBarry Smith     ierr = VecView_MPI_HDF5_DA(xin,viewer);CHKERRQ(ierr);
70247c6ae99SBarry Smith #endif
7038135c375SStefano Zampini   } else if (isglvis) {
7048135c375SStefano Zampini     ierr = VecView_GLVis(xin,viewer);CHKERRQ(ierr);
70547c6ae99SBarry Smith   } else {
70647c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
70747c6ae99SBarry Smith     PetscBool isbinary,isMPIIO;
70847c6ae99SBarry Smith 
709251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
71047c6ae99SBarry Smith     if (isbinary) {
711bc196f7cSDave May       ierr = PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO);CHKERRQ(ierr);
71247c6ae99SBarry Smith       if (isMPIIO) {
713aa219208SBarry Smith         ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_TRUE);CHKERRQ(ierr);
71447c6ae99SBarry Smith         PetscFunctionReturn(0);
71547c6ae99SBarry Smith       }
71647c6ae99SBarry Smith     }
71747c6ae99SBarry Smith #endif
71847c6ae99SBarry Smith 
71947c6ae99SBarry Smith     /* call viewer on natural ordering */
72047c6ae99SBarry Smith     ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr);
721aa219208SBarry Smith     ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr);
72247c6ae99SBarry Smith     ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr);
723aa219208SBarry Smith     ierr = DMDAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr);
724aa219208SBarry Smith     ierr = DMDAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr);
7253f3fd955SJed Brown     ierr = PetscObjectGetName((PetscObject)xin,&name);CHKERRQ(ierr);
7263f3fd955SJed Brown     ierr = PetscObjectSetName((PetscObject)natural,name);CHKERRQ(ierr);
727a261c58fSBarry Smith 
728a261c58fSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
729a261c58fSBarry Smith     if (format == PETSC_VIEWER_BINARY_MATLAB) {
730a261c58fSBarry Smith       /* temporarily remove viewer format so it won't trigger in the VecView() */
7316a9046bcSBarry Smith       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
732a261c58fSBarry Smith     }
733a261c58fSBarry Smith 
73423f88b65SBarry Smith     ((PetscObject)natural)->donotPetscObjectPrintClassNamePrefixType = PETSC_TRUE;
73547c6ae99SBarry Smith     ierr = VecView(natural,viewer);CHKERRQ(ierr);
73623f88b65SBarry Smith     ((PetscObject)natural)->donotPetscObjectPrintClassNamePrefixType = PETSC_FALSE;
737a261c58fSBarry Smith 
738a261c58fSBarry Smith     if (format == PETSC_VIEWER_BINARY_MATLAB) {
739a261c58fSBarry Smith       MPI_Comm    comm;
740a261c58fSBarry Smith       FILE        *info;
741a261c58fSBarry Smith       const char  *fieldname;
742da88d4d4SJed Brown       char        fieldbuf[256];
743a261c58fSBarry Smith       PetscInt    dim,ni,nj,nk,pi,pj,pk,dof,n;
744a261c58fSBarry Smith 
745a261c58fSBarry Smith       /* set the viewer format back into the viewer */
7466a9046bcSBarry Smith       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
747a261c58fSBarry Smith       ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
748a261c58fSBarry Smith       ierr = PetscViewerBinaryGetInfoPointer(viewer,&info);CHKERRQ(ierr);
749ea78f98cSLisandro Dalcin       ierr = DMDAGetInfo(da,&dim,&ni,&nj,&nk,&pi,&pj,&pk,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
750da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n");CHKERRQ(ierr);
751da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#$$ tmp = PetscBinaryRead(fd); \n");CHKERRQ(ierr);
752da88d4d4SJed Brown       if (dim == 1) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d);\n",dof,ni);CHKERRQ(ierr); }
753da88d4d4SJed Brown       if (dim == 2) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d);\n",dof,ni,nj);CHKERRQ(ierr); }
754da88d4d4SJed Brown       if (dim == 3) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d,%d);\n",dof,ni,nj,nk);CHKERRQ(ierr); }
755a261c58fSBarry Smith 
756a261c58fSBarry Smith       for (n=0; n<dof; n++) {
757a261c58fSBarry Smith         ierr = DMDAGetFieldName(da,n,&fieldname);CHKERRQ(ierr);
758da88d4d4SJed Brown         if (!fieldname || !fieldname[0]) {
759da88d4d4SJed Brown           ierr = PetscSNPrintf(fieldbuf,sizeof fieldbuf,"field%D",n);CHKERRQ(ierr);
760da88d4d4SJed Brown           fieldname = fieldbuf;
761a261c58fSBarry Smith         }
762da88d4d4SJed Brown         if (dim == 1) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:))';\n",name,fieldname,n+1);CHKERRQ(ierr); }
763da88d4d4SJed Brown         if (dim == 2) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:,:))';\n",name,fieldname,n+1);CHKERRQ(ierr); }
764da88d4d4SJed 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);}
765a261c58fSBarry Smith       }
766da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#$$ clear tmp; \n");CHKERRQ(ierr);
767da88d4d4SJed Brown       ierr = PetscFPrintf(comm,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n");CHKERRQ(ierr);
768a261c58fSBarry Smith     }
769a261c58fSBarry Smith 
770fcfd50ebSBarry Smith     ierr = VecDestroy(&natural);CHKERRQ(ierr);
77147c6ae99SBarry Smith   }
77247c6ae99SBarry Smith   PetscFunctionReturn(0);
77347c6ae99SBarry Smith }
77447c6ae99SBarry Smith 
77547c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
77647c6ae99SBarry Smith PetscErrorCode VecLoad_HDF5_DA(Vec xin, PetscViewer viewer)
77747c6ae99SBarry Smith {
7782d95be3dSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7799a42bb27SBarry Smith   DM             da;
78047c6ae99SBarry Smith   PetscErrorCode ierr;
78115b861d2SVaclav Hapla   int            dim,rdim;
782ec7ae49cSHåkon Strandenes   hsize_t        dims[6]={0},count[6]={0},offset[6]={0};
783*d7dd068bSVaclav Hapla   PetscBool      dim2=PETSC_FALSE,timestepping=PETSC_FALSE;
784*d7dd068bSVaclav Hapla   PetscInt       dimension,timestep=PETSC_MIN_INT,dofInd;
78547c6ae99SBarry Smith   PetscScalar    *x;
78647c6ae99SBarry Smith   const char     *vecname;
78747c6ae99SBarry Smith   hid_t          filespace; /* file dataspace identifier */
78847c6ae99SBarry Smith   hid_t          dset_id;   /* dataset identifier */
78947c6ae99SBarry Smith   hid_t          memspace;  /* memory dataspace identifier */
790bfd264e7SBarry Smith   hid_t          file_id,group;
7917bcbaff4SBarry Smith   hid_t          scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
7929c7c4993SBarry Smith   DM_DA          *dd;
79347c6ae99SBarry Smith 
79447c6ae99SBarry Smith   PetscFunctionBegin;
7957bcbaff4SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
7967bcbaff4SBarry Smith   scalartype = H5T_NATIVE_FLOAT;
7977bcbaff4SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
7987bcbaff4SBarry Smith #error "HDF5 output with 128 bit floats not supported."
799570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
800570b7f6dSBarry Smith #error "HDF5 output with 16 bit floats not supported."
8017bcbaff4SBarry Smith #else
8027bcbaff4SBarry Smith   scalartype = H5T_NATIVE_DOUBLE;
8037bcbaff4SBarry Smith #endif
8047bcbaff4SBarry Smith 
805bfd264e7SBarry Smith   ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr);
806ec7ae49cSHåkon Strandenes   ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr);
807*d7dd068bSVaclav Hapla   ierr = PetscViewerHDF5CheckTimestepping_Internal(viewer, vecname);CHKERRQ(ierr);
808*d7dd068bSVaclav Hapla   ierr = PetscViewerHDF5IsTimestepping(viewer, &timestepping);CHKERRQ(ierr);
809*d7dd068bSVaclav Hapla   if (timestepping) {
810*d7dd068bSVaclav Hapla     ierr = PetscViewerHDF5GetTimestep(viewer, &timestep);CHKERRQ(ierr);
811*d7dd068bSVaclav Hapla   }
812c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
8139c7c4993SBarry Smith   dd   = (DM_DA*)da->data;
814c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(da, &dimension);CHKERRQ(ierr);
81547c6ae99SBarry Smith 
816ec7ae49cSHåkon Strandenes   /* Open dataset */
817729ab6d9SBarry Smith   PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT));
81847c6ae99SBarry Smith 
819ec7ae49cSHåkon Strandenes   /* Retrieve the dataspace for the dataset */
820ec7ae49cSHåkon Strandenes   PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id));
821ec7ae49cSHåkon Strandenes   PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(filespace, dims, NULL));
822ec7ae49cSHåkon Strandenes 
823ec7ae49cSHåkon Strandenes   /* Expected dimension for holding the dof's */
82447c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
825ec7ae49cSHåkon Strandenes   dofInd = rdim-2;
826ec7ae49cSHåkon Strandenes #else
827ec7ae49cSHåkon Strandenes   dofInd = rdim-1;
82847c6ae99SBarry Smith #endif
829ec7ae49cSHåkon Strandenes 
830ec7ae49cSHåkon Strandenes   /* The expected number of dimensions, assuming basedimension2 = false */
831ec7ae49cSHåkon Strandenes   dim = dimension;
832ec7ae49cSHåkon Strandenes   if (dd->w > 1) ++dim;
833ec7ae49cSHåkon Strandenes   if (timestep >= 0) ++dim;
83447c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX)
835ec7ae49cSHåkon Strandenes   ++dim;
83647c6ae99SBarry Smith #endif
837ec7ae49cSHåkon Strandenes 
838ec7ae49cSHåkon Strandenes   /* In this case the input dataset have one extra, unexpected dimension. */
839ec7ae49cSHåkon Strandenes   if (rdim == dim+1) {
840ec7ae49cSHåkon Strandenes     /* In this case the block size unity */
841ec7ae49cSHåkon Strandenes     if (dd->w == 1 && dims[dofInd] == 1) dim2 = PETSC_TRUE;
842ec7ae49cSHåkon Strandenes 
843ec7ae49cSHåkon Strandenes     /* Special error message for the case where dof does not match the input file */
844ec7ae49cSHå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);
845ec7ae49cSHåkon Strandenes 
846ec7ae49cSHåkon Strandenes   /* Other cases where rdim != dim cannot be handled currently */
8476c4ed002SBarry 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);
848ec7ae49cSHåkon Strandenes 
849ec7ae49cSHåkon Strandenes   /* Set up the hyperslab size */
850ec7ae49cSHåkon Strandenes   dim = 0;
851ec7ae49cSHåkon Strandenes   if (timestep >= 0) {
852ec7ae49cSHåkon Strandenes     offset[dim] = timestep;
853ec7ae49cSHåkon Strandenes     count[dim] = 1;
854ec7ae49cSHåkon Strandenes     ++dim;
855ec7ae49cSHåkon Strandenes   }
856ec7ae49cSHåkon Strandenes   if (dimension == 3) {
857ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->zs,offset + dim);CHKERRQ(ierr);
858ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->ze - dd->zs,count + dim);CHKERRQ(ierr);
859ec7ae49cSHåkon Strandenes     ++dim;
860ec7ae49cSHåkon Strandenes   }
861ec7ae49cSHåkon Strandenes   if (dimension > 1) {
862ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->ys,offset + dim);CHKERRQ(ierr);
863ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->ye - dd->ys,count + dim);CHKERRQ(ierr);
864ec7ae49cSHåkon Strandenes     ++dim;
865ec7ae49cSHåkon Strandenes   }
866ec7ae49cSHåkon Strandenes   ierr = PetscHDF5IntCast(dd->xs/dd->w,offset + dim);CHKERRQ(ierr);
867ec7ae49cSHåkon Strandenes   ierr = PetscHDF5IntCast((dd->xe - dd->xs)/dd->w,count + dim);CHKERRQ(ierr);
868ec7ae49cSHåkon Strandenes   ++dim;
869ec7ae49cSHåkon Strandenes   if (dd->w > 1 || dim2) {
870ec7ae49cSHåkon Strandenes     offset[dim] = 0;
871ec7ae49cSHåkon Strandenes     ierr = PetscHDF5IntCast(dd->w,count + dim);CHKERRQ(ierr);
872ec7ae49cSHåkon Strandenes     ++dim;
873ec7ae49cSHåkon Strandenes   }
874ec7ae49cSHåkon Strandenes #if defined(PETSC_USE_COMPLEX)
875ec7ae49cSHåkon Strandenes   offset[dim] = 0;
876ec7ae49cSHåkon Strandenes   count[dim] = 2;
877ec7ae49cSHåkon Strandenes   ++dim;
878ec7ae49cSHåkon Strandenes #endif
879ec7ae49cSHåkon Strandenes 
880ec7ae49cSHåkon Strandenes   /* Create the memory and filespace */
881729ab6d9SBarry Smith   PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL));
882729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL));
88347c6ae99SBarry Smith 
88447c6ae99SBarry Smith   ierr   = VecGetArray(xin, &x);CHKERRQ(ierr);
8852d95be3dSVaclav Hapla   PetscStackCallHDF5(H5Dread,(dset_id, scalartype, memspace, filespace, hdf5->dxpl_id, x));
88647c6ae99SBarry Smith   ierr   = VecRestoreArray(xin, &x);CHKERRQ(ierr);
88747c6ae99SBarry Smith 
88847c6ae99SBarry Smith   /* Close/release resources */
889bfd264e7SBarry Smith   if (group != file_id) {
890729ab6d9SBarry Smith     PetscStackCallHDF5(H5Gclose,(group));
891bfd264e7SBarry Smith   }
892729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(filespace));
893729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(memspace));
894729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dset_id));
89547c6ae99SBarry Smith   PetscFunctionReturn(0);
89647c6ae99SBarry Smith }
89747c6ae99SBarry Smith #endif
89847c6ae99SBarry Smith 
89947c6ae99SBarry Smith PetscErrorCode VecLoad_Binary_DA(Vec xin, PetscViewer viewer)
90047c6ae99SBarry Smith {
9019a42bb27SBarry Smith   DM             da;
90247c6ae99SBarry Smith   PetscErrorCode ierr;
90347c6ae99SBarry Smith   Vec            natural;
90447c6ae99SBarry Smith   const char     *prefix;
90547c6ae99SBarry Smith   PetscInt       bs;
90647c6ae99SBarry Smith   PetscBool      flag;
90747c6ae99SBarry Smith   DM_DA          *dd;
90847c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
90947c6ae99SBarry Smith   PetscBool isMPIIO;
91047c6ae99SBarry Smith #endif
91147c6ae99SBarry Smith 
91247c6ae99SBarry Smith   PetscFunctionBegin;
913c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
91447c6ae99SBarry Smith   dd   = (DM_DA*)da->data;
91547c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO)
916bc196f7cSDave May   ierr = PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO);CHKERRQ(ierr);
91747c6ae99SBarry Smith   if (isMPIIO) {
918aa219208SBarry Smith     ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_FALSE);CHKERRQ(ierr);
91947c6ae99SBarry Smith     PetscFunctionReturn(0);
92047c6ae99SBarry Smith   }
92147c6ae99SBarry Smith #endif
92247c6ae99SBarry Smith 
92347c6ae99SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr);
924aa219208SBarry Smith   ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr);
92547c6ae99SBarry Smith   ierr = PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name);CHKERRQ(ierr);
92647c6ae99SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr);
9278aea9937SBarry Smith   ierr = VecLoad(natural,viewer);CHKERRQ(ierr);
928aa219208SBarry Smith   ierr = DMDANaturalToGlobalBegin(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr);
929aa219208SBarry Smith   ierr = DMDANaturalToGlobalEnd(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr);
930fcfd50ebSBarry Smith   ierr = VecDestroy(&natural);CHKERRQ(ierr);
931aa219208SBarry Smith   ierr = PetscInfo(xin,"Loading vector from natural ordering into DMDA\n");CHKERRQ(ierr);
932c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(NULL,((PetscObject)xin)->prefix,"-vecload_block_size",&bs,&flag);CHKERRQ(ierr);
93347c6ae99SBarry Smith   if (flag && bs != dd->w) {
934aa219208SBarry Smith     ierr = PetscInfo2(xin,"Block size in file %D not equal to DMDA's dof %D\n",bs,dd->w);CHKERRQ(ierr);
93547c6ae99SBarry Smith   }
93647c6ae99SBarry Smith   PetscFunctionReturn(0);
93747c6ae99SBarry Smith }
93847c6ae99SBarry Smith 
9397087cfbeSBarry Smith PetscErrorCode  VecLoad_Default_DA(Vec xin, PetscViewer viewer)
94047c6ae99SBarry Smith {
94147c6ae99SBarry Smith   PetscErrorCode ierr;
9429a42bb27SBarry Smith   DM             da;
94347c6ae99SBarry Smith   PetscBool      isbinary;
94447c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
94547c6ae99SBarry Smith   PetscBool ishdf5;
94647c6ae99SBarry Smith #endif
94747c6ae99SBarry Smith 
94847c6ae99SBarry Smith   PetscFunctionBegin;
949c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
950ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
95147c6ae99SBarry Smith 
95247c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
953251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr);
95447c6ae99SBarry Smith #endif
955251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
95647c6ae99SBarry Smith 
95747c6ae99SBarry Smith   if (isbinary) {
95847c6ae99SBarry Smith     ierr = VecLoad_Binary_DA(xin,viewer);CHKERRQ(ierr);
95947c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5)
96047c6ae99SBarry Smith   } else if (ishdf5) {
96147c6ae99SBarry Smith     ierr = VecLoad_HDF5_DA(xin,viewer);CHKERRQ(ierr);
96247c6ae99SBarry Smith #endif
963d34fcf5fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for vector loading", ((PetscObject)viewer)->type_name);
96447c6ae99SBarry Smith   PetscFunctionReturn(0);
96547c6ae99SBarry Smith }
966