147c6ae99SBarry Smith 247c6ae99SBarry Smith /* 3aa219208SBarry Smith Plots vectors obtained with DMDACreate2d() 447c6ae99SBarry Smith */ 547c6ae99SBarry Smith 64035e84dSBarry Smith #include <petsc-private/dmdaimpl.h> /*I "petscdmda.h" I*/ 7b45d2f2cSJed Brown #include <petsc-private/vecimpl.h> 89804daf3SBarry Smith #include <petscdraw.h> 9665c2dedSJed Brown #include <petscviewerhdf5.h> 1047c6ae99SBarry Smith 1147c6ae99SBarry Smith /* 1247c6ae99SBarry Smith The data that is passed into the graphics callback 1347c6ae99SBarry Smith */ 1447c6ae99SBarry Smith typedef struct { 1547c6ae99SBarry Smith PetscInt m,n,step,k; 1647c6ae99SBarry Smith PetscReal min,max,scale; 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; 340e4fe250SBarry Smith PetscReal s,min,max,x1,x2,x3,x4,y_1,y2,y3,y4,xmin = PETSC_MAX_REAL,xmax = PETSC_MIN_REAL,ymin = PETSC_MAX_REAL,ymax = PETSC_MIN_REAL; 350e4fe250SBarry Smith PetscReal xminf,xmaxf,yminf,ymaxf,w; 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 s = zctx->scale; 4847c6ae99SBarry Smith min = zctx->min; 49f3f0eb19SBarry Smith max = zctx->max; 5047c6ae99SBarry Smith 5147c6ae99SBarry Smith /* PetscDraw the contour plot patch */ 5247c6ae99SBarry Smith for (j=0; j<n-1; j++) { 5347c6ae99SBarry Smith for (i=0; i<m-1; i++) { 540e4fe250SBarry Smith id = i+j*m; 550e4fe250SBarry Smith x1 = PetscRealPart(xy[2*id]); 560e4fe250SBarry Smith y_1 = PetscRealPart(xy[2*id+1]); 570e4fe250SBarry Smith c1 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscClipInterval(PetscRealPart(v[k+step*id]),min,max)-min)); 580e4fe250SBarry Smith xmin = PetscMin(xmin,x1); 590e4fe250SBarry Smith ymin = PetscMin(ymin,y_1); 600e4fe250SBarry Smith xmax = PetscMax(xmax,x1); 610e4fe250SBarry Smith ymax = PetscMax(ymax,y_1); 620e4fe250SBarry Smith 630e4fe250SBarry Smith id = i+j*m+1; 640e4fe250SBarry Smith x2 = PetscRealPart(xy[2*id]); 65a39a4c7dSLisandro Dalcin y2 = PetscRealPart(xy[2*id+1]); 660e4fe250SBarry Smith c2 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscClipInterval(PetscRealPart(v[k+step*id]),min,max)-min)); 670e4fe250SBarry Smith xmin = PetscMin(xmin,x2); 68a39a4c7dSLisandro Dalcin ymin = PetscMin(ymin,y2); 690e4fe250SBarry Smith xmax = PetscMax(xmax,x2); 70a39a4c7dSLisandro Dalcin ymax = PetscMax(ymax,y2); 710e4fe250SBarry Smith 720e4fe250SBarry Smith id = i+j*m+1+m; 73a39a4c7dSLisandro Dalcin x3 = PetscRealPart(xy[2*id]); 740e4fe250SBarry Smith y3 = PetscRealPart(xy[2*id+1]); 750e4fe250SBarry Smith c3 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscClipInterval(PetscRealPart(v[k+step*id]),min,max)-min)); 76a39a4c7dSLisandro Dalcin xmin = PetscMin(xmin,x3); 770e4fe250SBarry Smith ymin = PetscMin(ymin,y3); 78a39a4c7dSLisandro Dalcin xmax = PetscMax(xmax,x3); 790e4fe250SBarry Smith ymax = PetscMax(ymax,y3); 800e4fe250SBarry Smith 810e4fe250SBarry Smith id = i+j*m+m; 82a39a4c7dSLisandro Dalcin x4 = PetscRealPart(xy[2*id]); 83a39a4c7dSLisandro Dalcin y4 = PetscRealPart(xy[2*id+1]); 840e4fe250SBarry Smith c4 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscClipInterval(PetscRealPart(v[k+step*id]),min,max)-min)); 85a39a4c7dSLisandro Dalcin xmin = PetscMin(xmin,x4); 86a39a4c7dSLisandro Dalcin ymin = PetscMin(ymin,y4); 87a39a4c7dSLisandro Dalcin xmax = PetscMax(xmax,x4); 88a39a4c7dSLisandro Dalcin ymax = PetscMax(ymax,y4); 89f3f0eb19SBarry Smith 9047c6ae99SBarry Smith ierr = PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);CHKERRQ(ierr); 9147c6ae99SBarry Smith ierr = PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4);CHKERRQ(ierr); 9247c6ae99SBarry Smith if (zctx->showgrid) { 9347c6ae99SBarry Smith ierr = PetscDrawLine(draw,x1,y_1,x2,y2,PETSC_DRAW_BLACK);CHKERRQ(ierr); 9447c6ae99SBarry Smith ierr = PetscDrawLine(draw,x2,y2,x3,y3,PETSC_DRAW_BLACK);CHKERRQ(ierr); 9547c6ae99SBarry Smith ierr = PetscDrawLine(draw,x3,y3,x4,y4,PETSC_DRAW_BLACK);CHKERRQ(ierr); 9647c6ae99SBarry Smith ierr = PetscDrawLine(draw,x4,y4,x1,y_1,PETSC_DRAW_BLACK);CHKERRQ(ierr); 9747c6ae99SBarry Smith } 9847c6ae99SBarry Smith } 9947c6ae99SBarry Smith } 100109c9344SBarry Smith if (zctx->name0) { 101109c9344SBarry Smith PetscReal xl,yl,xr,yr,x,y; 102109c9344SBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 103109c9344SBarry Smith x = xl + .3*(xr - xl); 104109c9344SBarry Smith xl = xl + .01*(xr - xl); 105109c9344SBarry Smith y = yr - .3*(yr - yl); 106109c9344SBarry Smith yl = yl + .01*(yr - yl); 107109c9344SBarry Smith ierr = PetscDrawString(draw,x,yl,PETSC_DRAW_BLACK,zctx->name0);CHKERRQ(ierr); 108109c9344SBarry Smith ierr = PetscDrawStringVertical(draw,xl,y,PETSC_DRAW_BLACK,zctx->name1);CHKERRQ(ierr); 109109c9344SBarry Smith } 1100e4fe250SBarry Smith /* 1110e4fe250SBarry Smith Ideally we would use the PetscDrawAxis object to manage displaying the coordinate limits 1120e4fe250SBarry Smith but that may require some refactoring. 1130e4fe250SBarry Smith */ 114953e6c96SLisandro Dalcin ierr = MPI_Allreduce(&xmin,&xminf,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)draw));CHKERRQ(ierr); 115ce94432eSBarry Smith ierr = MPI_Allreduce(&xmax,&xmaxf,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)draw));CHKERRQ(ierr); 116953e6c96SLisandro Dalcin ierr = MPI_Allreduce(&ymin,&yminf,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)draw));CHKERRQ(ierr); 117ce94432eSBarry Smith ierr = MPI_Allreduce(&ymax,&ymaxf,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)draw));CHKERRQ(ierr); 1180e4fe250SBarry Smith ierr = PetscSNPrintf(value,16,"%f",xminf);CHKERRQ(ierr); 1190e4fe250SBarry Smith ierr = PetscDrawString(draw,xminf,yminf - .05*(ymaxf - yminf),PETSC_DRAW_BLACK,value);CHKERRQ(ierr); 1200e4fe250SBarry Smith ierr = PetscSNPrintf(value,16,"%f",xmaxf);CHKERRQ(ierr); 1210e4fe250SBarry Smith ierr = PetscStrlen(value,&len);CHKERRQ(ierr); 1220298fd71SBarry Smith ierr = PetscDrawStringGetSize(draw,&w,NULL);CHKERRQ(ierr); 1230e4fe250SBarry Smith ierr = PetscDrawString(draw,xmaxf - len*w,yminf - .05*(ymaxf - yminf),PETSC_DRAW_BLACK,value);CHKERRQ(ierr); 1240e4fe250SBarry Smith ierr = PetscSNPrintf(value,16,"%f",yminf);CHKERRQ(ierr); 1250e4fe250SBarry Smith ierr = PetscDrawString(draw,xminf - .05*(xmaxf - xminf),yminf,PETSC_DRAW_BLACK,value);CHKERRQ(ierr); 1260e4fe250SBarry Smith ierr = PetscSNPrintf(value,16,"%f",ymaxf);CHKERRQ(ierr); 1270e4fe250SBarry Smith ierr = PetscDrawString(draw,xminf - .05*(xmaxf - xminf),ymaxf,PETSC_DRAW_BLACK,value);CHKERRQ(ierr); 12847c6ae99SBarry Smith PetscFunctionReturn(0); 12947c6ae99SBarry Smith } 13047c6ae99SBarry Smith 13147c6ae99SBarry Smith #undef __FUNCT__ 13247c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_Draw_DA2d" 13347c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA2d(Vec xin,PetscViewer viewer) 13447c6ae99SBarry Smith { 1359a42bb27SBarry Smith DM da,dac,dag; 13647c6ae99SBarry Smith PetscErrorCode ierr; 13747c6ae99SBarry Smith PetscMPIInt rank; 138a74ba6f7SBarry Smith PetscInt N,s,M,w,ncoors = 4; 13947c6ae99SBarry Smith const PetscInt *lx,*ly; 14047c6ae99SBarry Smith PetscReal coors[4],ymin,ymax,xmin,xmax; 14147c6ae99SBarry Smith PetscDraw draw,popup; 14247c6ae99SBarry Smith PetscBool isnull,useports = PETSC_FALSE; 14347c6ae99SBarry Smith MPI_Comm comm; 14447c6ae99SBarry Smith Vec xlocal,xcoor,xcoorl; 145bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by; 146aa219208SBarry Smith DMDAStencilType st; 14747c6ae99SBarry Smith ZoomCtx zctx; 1480298fd71SBarry Smith PetscDrawViewPorts *ports = NULL; 14947c6ae99SBarry Smith PetscViewerFormat format; 15020d0051dSBarry Smith PetscInt *displayfields; 15167dd0837SBarry Smith PetscInt ndisplayfields,i,nbounds; 15267dd0837SBarry Smith const PetscReal *bounds; 15347c6ae99SBarry Smith 15447c6ae99SBarry Smith PetscFunctionBegin; 15547c6ae99SBarry Smith zctx.showgrid = PETSC_FALSE; 1568865f1eaSKarl Rupp 15747c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 15847c6ae99SBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 15903193ff8SBarry Smith ierr = PetscViewerDrawGetBounds(viewer,&nbounds,&bounds);CHKERRQ(ierr); 16047c6ae99SBarry Smith 161c688c046SMatthew G Knepley ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 162ce94432eSBarry Smith if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 16347c6ae99SBarry Smith 16447c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr); 16547c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 16647c6ae99SBarry Smith 1671321219cSEthan Coon ierr = DMDAGetInfo(da,0,&M,&N,0,&zctx.m,&zctx.n,0,&w,&s,&bx,&by,0,&st);CHKERRQ(ierr); 1680298fd71SBarry Smith ierr = DMDAGetOwnershipRanges(da,&lx,&ly,NULL);CHKERRQ(ierr); 16947c6ae99SBarry Smith 17047c6ae99SBarry Smith /* 17147c6ae99SBarry Smith Obtain a sequential vector that is going to contain the local values plus ONE layer of 172aa219208SBarry Smith ghosted values to draw the graphics from. We also need its corresponding DMDA (dac) that will 17347c6ae99SBarry Smith update the local values pluse ONE layer of ghost values. 17447c6ae99SBarry Smith */ 17547c6ae99SBarry Smith ierr = PetscObjectQuery((PetscObject)da,"GraphicsGhosted",(PetscObject*)&xlocal);CHKERRQ(ierr); 17647c6ae99SBarry Smith if (!xlocal) { 177bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_NONE || by != DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) { 17847c6ae99SBarry Smith /* 17947c6ae99SBarry Smith if original da is not of stencil width one, or periodic or not a box stencil then 180aa219208SBarry Smith create a special DMDA to handle one level of ghost points for graphics 18147c6ae99SBarry Smith */ 182bff4a2f0SMatthew 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); 183aa219208SBarry Smith ierr = PetscInfo(da,"Creating auxilary DMDA for managing graphics ghost points\n");CHKERRQ(ierr); 18447c6ae99SBarry Smith } else { 18547c6ae99SBarry Smith /* otherwise we can use the da we already have */ 18647c6ae99SBarry Smith dac = da; 18747c6ae99SBarry Smith } 18847c6ae99SBarry Smith /* create local vector for holding ghosted values used in graphics */ 189564755cdSBarry Smith ierr = DMCreateLocalVector(dac,&xlocal);CHKERRQ(ierr); 19047c6ae99SBarry Smith if (dac != da) { 191aa219208SBarry Smith /* don't keep any public reference of this DMDA, is is only available through xlocal */ 192f7923d8aSBarry Smith ierr = PetscObjectDereference((PetscObject)dac);CHKERRQ(ierr); 19347c6ae99SBarry Smith } else { 19447c6ae99SBarry Smith /* remove association between xlocal and da, because below we compose in the opposite 19547c6ae99SBarry Smith direction and if we left this connect we'd get a loop, so the objects could 19647c6ae99SBarry Smith never be destroyed */ 197c688c046SMatthew G Knepley ierr = PetscObjectRemoveReference((PetscObject)xlocal,"__PETSc_dm");CHKERRQ(ierr); 19847c6ae99SBarry Smith } 19947c6ae99SBarry Smith ierr = PetscObjectCompose((PetscObject)da,"GraphicsGhosted",(PetscObject)xlocal);CHKERRQ(ierr); 20047c6ae99SBarry Smith ierr = PetscObjectDereference((PetscObject)xlocal);CHKERRQ(ierr); 20147c6ae99SBarry Smith } else { 202bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_NONE || by != DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) { 203c688c046SMatthew G Knepley ierr = VecGetDM(xlocal, &dac);CHKERRQ(ierr); 204f7923d8aSBarry Smith } else { 205f7923d8aSBarry Smith dac = da; 20647c6ae99SBarry Smith } 20747c6ae99SBarry Smith } 20847c6ae99SBarry Smith 20947c6ae99SBarry Smith /* 21047c6ae99SBarry Smith Get local (ghosted) values of vector 21147c6ae99SBarry Smith */ 2129a42bb27SBarry Smith ierr = DMGlobalToLocalBegin(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr); 2139a42bb27SBarry Smith ierr = DMGlobalToLocalEnd(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr); 2145edff71fSBarry Smith ierr = VecGetArrayRead(xlocal,&zctx.v);CHKERRQ(ierr); 21547c6ae99SBarry Smith 21647c6ae99SBarry Smith /* get coordinates of nodes */ 2176636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); 21847c6ae99SBarry Smith if (!xcoor) { 219aa219208SBarry Smith ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,0.0);CHKERRQ(ierr); 2206636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); 22147c6ae99SBarry Smith } 22247c6ae99SBarry Smith 22347c6ae99SBarry Smith /* 22447c6ae99SBarry Smith Determine the min and max coordinates in plot 22547c6ae99SBarry Smith */ 2260298fd71SBarry Smith ierr = VecStrideMin(xcoor,0,NULL,&xmin);CHKERRQ(ierr); 2270298fd71SBarry Smith ierr = VecStrideMax(xcoor,0,NULL,&xmax);CHKERRQ(ierr); 2280298fd71SBarry Smith ierr = VecStrideMin(xcoor,1,NULL,&ymin);CHKERRQ(ierr); 2290298fd71SBarry Smith ierr = VecStrideMax(xcoor,1,NULL,&ymax);CHKERRQ(ierr); 23047c6ae99SBarry Smith coors[0] = xmin - .05*(xmax- xmin); coors[2] = xmax + .05*(xmax - xmin); 23147c6ae99SBarry Smith coors[1] = ymin - .05*(ymax- ymin); coors[3] = ymax + .05*(ymax - ymin); 23257622a8eSBarry 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); 23347c6ae99SBarry Smith 234a74ba6f7SBarry Smith ierr = PetscOptionsGetRealArray(NULL,"-draw_coordinates",coors,&ncoors,NULL);CHKERRQ(ierr); 235a74ba6f7SBarry Smith 23647c6ae99SBarry Smith /* 23747c6ae99SBarry Smith get local ghosted version of coordinates 23847c6ae99SBarry Smith */ 23947c6ae99SBarry Smith ierr = PetscObjectQuery((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr); 24047c6ae99SBarry Smith if (!xcoorl) { 241aa219208SBarry Smith /* create DMDA to get local version of graphics */ 242bff4a2f0SMatthew 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); 243aa219208SBarry Smith ierr = PetscInfo(dag,"Creating auxilary DMDA for managing graphics coordinates ghost points\n");CHKERRQ(ierr); 244564755cdSBarry Smith ierr = DMCreateLocalVector(dag,&xcoorl);CHKERRQ(ierr); 24547c6ae99SBarry Smith ierr = PetscObjectCompose((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject)xcoorl);CHKERRQ(ierr); 246f7923d8aSBarry Smith ierr = PetscObjectDereference((PetscObject)dag);CHKERRQ(ierr); 24747c6ae99SBarry Smith ierr = PetscObjectDereference((PetscObject)xcoorl);CHKERRQ(ierr); 24847c6ae99SBarry Smith } else { 249c688c046SMatthew G Knepley ierr = VecGetDM(xcoorl,&dag);CHKERRQ(ierr); 25047c6ae99SBarry Smith } 2519a42bb27SBarry Smith ierr = DMGlobalToLocalBegin(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr); 2529a42bb27SBarry Smith ierr = DMGlobalToLocalEnd(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr); 2535edff71fSBarry Smith ierr = VecGetArrayRead(xcoorl,&zctx.xy);CHKERRQ(ierr); 25447c6ae99SBarry Smith 25547c6ae99SBarry Smith /* 25647c6ae99SBarry Smith Get information about size of area each processor must do graphics for 25747c6ae99SBarry Smith */ 2581321219cSEthan Coon ierr = DMDAGetInfo(dac,0,&M,&N,0,0,0,0,&zctx.step,0,&bx,&by,0,0);CHKERRQ(ierr); 259f7923d8aSBarry Smith ierr = DMDAGetGhostCorners(dac,0,0,0,&zctx.m,&zctx.n,0);CHKERRQ(ierr); 26047c6ae99SBarry Smith 2610298fd71SBarry Smith ierr = PetscOptionsGetBool(NULL,"-draw_contour_grid",&zctx.showgrid,NULL);CHKERRQ(ierr); 2624e6118eeSBarry Smith 2634e6118eeSBarry Smith ierr = DMDASelectFields(da,&ndisplayfields,&displayfields);CHKERRQ(ierr); 26447c6ae99SBarry Smith 26547c6ae99SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 2660298fd71SBarry Smith ierr = PetscOptionsGetBool(NULL,"-draw_ports",&useports,NULL);CHKERRQ(ierr); 26747c6ae99SBarry Smith if (useports || format == PETSC_VIEWER_DRAW_PORTS) { 26847c6ae99SBarry Smith ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr); 26920d0051dSBarry Smith ierr = PetscDrawViewPortsCreate(draw,ndisplayfields,&ports);CHKERRQ(ierr); 270109c9344SBarry Smith zctx.name0 = 0; 271109c9344SBarry Smith zctx.name1 = 0; 272109c9344SBarry Smith } else { 273109c9344SBarry Smith ierr = DMDAGetCoordinateName(da,0,&zctx.name0);CHKERRQ(ierr); 274109c9344SBarry Smith ierr = DMDAGetCoordinateName(da,1,&zctx.name1);CHKERRQ(ierr); 27547c6ae99SBarry Smith } 27620d0051dSBarry Smith 27747c6ae99SBarry Smith /* 27847c6ae99SBarry Smith Loop over each field; drawing each in a different window 27947c6ae99SBarry Smith */ 28020d0051dSBarry Smith for (i=0; i<ndisplayfields; i++) { 28120d0051dSBarry Smith zctx.k = displayfields[i]; 28247c6ae99SBarry Smith if (useports) { 28320d0051dSBarry Smith ierr = PetscDrawViewPortsSet(ports,i);CHKERRQ(ierr); 2849332fd86SBarry Smith } else { 2859332fd86SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,i,&draw);CHKERRQ(ierr); 28647c6ae99SBarry Smith } 28747c6ae99SBarry Smith 28847c6ae99SBarry Smith /* 28947c6ae99SBarry Smith Determine the min and max color in plot 29047c6ae99SBarry Smith */ 2910298fd71SBarry Smith ierr = VecStrideMin(xin,zctx.k,NULL,&zctx.min);CHKERRQ(ierr); 2920298fd71SBarry Smith ierr = VecStrideMax(xin,zctx.k,NULL,&zctx.max);CHKERRQ(ierr); 29367dd0837SBarry Smith if (zctx.k < nbounds) { 294f3f0eb19SBarry Smith zctx.min = bounds[2*zctx.k]; 295f3f0eb19SBarry Smith zctx.max = bounds[2*zctx.k+1]; 29667dd0837SBarry Smith } 29747c6ae99SBarry Smith if (zctx.min == zctx.max) { 29847c6ae99SBarry Smith zctx.min -= 1.e-12; 29947c6ae99SBarry Smith zctx.max += 1.e-12; 30047c6ae99SBarry Smith } 30147c6ae99SBarry Smith 30247c6ae99SBarry Smith if (!rank) { 30347c6ae99SBarry Smith const char *title; 30447c6ae99SBarry Smith 305aa219208SBarry Smith ierr = DMDAGetFieldName(da,zctx.k,&title);CHKERRQ(ierr); 30647c6ae99SBarry Smith if (title) { 30747c6ae99SBarry Smith ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr); 30847c6ae99SBarry Smith } 30947c6ae99SBarry Smith } 31047c6ae99SBarry Smith ierr = PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr); 31157622a8eSBarry Smith ierr = PetscInfo2(da,"DMDA 2d contour plot min %g max %g\n",(double)zctx.min,(double)zctx.max);CHKERRQ(ierr); 31247c6ae99SBarry Smith 31347c6ae99SBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 31447c6ae99SBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,zctx.min,zctx.max);CHKERRQ(ierr);} 31547c6ae99SBarry Smith 31647c6ae99SBarry Smith zctx.scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/(zctx.max - zctx.min); 31747c6ae99SBarry Smith 31847c6ae99SBarry Smith ierr = PetscDrawZoom(draw,VecView_MPI_Draw_DA2d_Zoom,&zctx);CHKERRQ(ierr); 31947c6ae99SBarry Smith } 32020d0051dSBarry Smith ierr = PetscFree(displayfields);CHKERRQ(ierr); 3216bf464f9SBarry Smith ierr = PetscDrawViewPortsDestroy(ports);CHKERRQ(ierr); 32247c6ae99SBarry Smith 3235edff71fSBarry Smith ierr = VecRestoreArrayRead(xcoorl,&zctx.xy);CHKERRQ(ierr); 3245edff71fSBarry Smith ierr = VecRestoreArrayRead(xlocal,&zctx.v);CHKERRQ(ierr); 32547c6ae99SBarry Smith PetscFunctionReturn(0); 32647c6ae99SBarry Smith } 32747c6ae99SBarry Smith 3280da24e51SJuha Jäykkä #if defined(PETSC_HAVE_HDF5) 3290da24e51SJuha Jäykkä #undef __FUNCT__ 3300da24e51SJuha Jäykkä #define __FUNCT__ "VecGetHDF5ChunkSize" 331c73cfb54SMatthew G. Knepley static PetscErrorCode VecGetHDF5ChunkSize(DM_DA *da, Vec xin, PetscInt dimension, PetscInt timestep, hsize_t *chunkDims) 3320da24e51SJuha Jäykkä { 333d4190030SJed Brown PetscMPIInt comm_size; 334d4190030SJed Brown PetscErrorCode ierr; 335561a82e7SJed Brown hsize_t chunk_size, target_size, dim; 336561a82e7SJed Brown hsize_t vec_size = sizeof(PetscScalar)*da->M*da->N*da->P*da->w; 337561a82e7SJed Brown hsize_t avg_local_vec_size,KiB = 1024,MiB = KiB*KiB,GiB = MiB*KiB,min_size = MiB; 338561a82e7SJed Brown hsize_t max_chunks = 64*KiB; /* HDF5 internal limitation */ 339561a82e7SJed Brown hsize_t max_chunk_size = 4*GiB; /* HDF5 internal limitation */ 34084179ffaSJed Brown PetscInt zslices=da->p, yslices=da->n, xslices=da->m; 3410da24e51SJuha Jäykkä 342cb5c4f94SJuha Jäykkä PetscFunctionBegin; 343d4190030SJed Brown ierr = MPI_Comm_size(PetscObjectComm((PetscObject)xin), &comm_size);CHKERRQ(ierr); 3440da24e51SJuha 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 */ 3450da24e51SJuha Jäykkä 346794a961bSBarry 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))))); 3477d310018SBarry 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 */ 3487d310018SBarry 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); 3490da24e51SJuha Jäykkä 350cb5c4f94SJuha Jäykkä /* 351cb5c4f94SJuha Jäykkä if size/rank > max_chunk_size, we need radical measures: even going down to 352cb5c4f94SJuha Jäykkä avg_local_vec_size is not enough, so we simply use chunk size of 4 GiB no matter 353cb5c4f94SJuha Jäykkä what, composed in the most efficient way possible. 354cb5c4f94SJuha Jäykkä N.B. this minimises the number of chunks, which may or may not be the optimal 355cb5c4f94SJuha Jäykkä solution. In a BG, for example, the optimal solution is probably to make # chunks = # 356cb5c4f94SJuha Jäykkä IO nodes involved, but this author has no access to a BG to figure out how to 357cb5c4f94SJuha Jäykkä reliably find the right number. And even then it may or may not be enough. 358cb5c4f94SJuha Jäykkä */ 3590da24e51SJuha Jäykkä if (avg_local_vec_size > max_chunk_size) { 360cb5c4f94SJuha Jäykkä /* check if we can just split local z-axis: is that enough? */ 36184179ffaSJed Brown zslices = (PetscInt)ceil(vec_size*1.0/(da->p*max_chunk_size))*zslices; 3620da24e51SJuha Jäykkä if (zslices > da->P) { 363cb5c4f94SJuha Jäykkä /* lattice is too large in xy-directions, splitting z only is not enough */ 3640da24e51SJuha Jäykkä zslices = da->P; 36584179ffaSJed Brown yslices= (PetscInt)ceil(vec_size*1.0/(zslices*da->n*max_chunk_size))*yslices; 3660da24e51SJuha Jäykkä if (yslices > da->N) { 367cb5c4f94SJuha Jäykkä /* lattice is too large in x-direction, splitting along z, y is not enough */ 3680da24e51SJuha Jäykkä yslices = da->N; 36984179ffaSJed Brown xslices= (PetscInt)ceil(vec_size*1.0/(zslices*yslices*da->m*max_chunk_size))*xslices; 3700da24e51SJuha Jäykkä } 3710da24e51SJuha Jäykkä } 3720da24e51SJuha Jäykkä dim = 0; 3730da24e51SJuha Jäykkä if (timestep >= 0) { 3740da24e51SJuha Jäykkä ++dim; 3750da24e51SJuha Jäykkä } 376cb5c4f94SJuha Jäykkä /* prefer to split z-axis, even down to planar slices */ 377c73cfb54SMatthew G. Knepley if (dimension == 3) { 378cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t) da->P/zslices; 379cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t) da->N/yslices; 380cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t) da->M/xslices; 3810da24e51SJuha Jäykkä } else { 382cb5c4f94SJuha Jäykkä /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */ 383cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t) da->N/yslices; 384cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t) da->M/xslices; 3850da24e51SJuha Jäykkä } 3860da24e51SJuha 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); 3870da24e51SJuha Jäykkä } else { 3880da24e51SJuha Jäykkä if (target_size < chunk_size) { 389cb5c4f94SJuha Jäykkä /* only change the defaults if target_size < chunk_size */ 3900da24e51SJuha Jäykkä dim = 0; 3910da24e51SJuha Jäykkä if (timestep >= 0) { 3920da24e51SJuha Jäykkä ++dim; 3930da24e51SJuha Jäykkä } 394cb5c4f94SJuha Jäykkä /* prefer to split z-axis, even down to planar slices */ 395c73cfb54SMatthew G. Knepley if (dimension == 3) { 396cb5c4f94SJuha Jäykkä /* try splitting the z-axis to core-size bits, i.e. divide chunk size by # comm_size in z-direction */ 3970da24e51SJuha Jäykkä if (target_size >= chunk_size/da->p) { 398cb5c4f94SJuha Jäykkä /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */ 3990da24e51SJuha Jäykkä chunkDims[dim] = (hsize_t) ceil(da->P*1.0/da->p); 4000da24e51SJuha Jäykkä } else { 401cb5c4f94SJuha Jäykkä /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be 402cb5c4f94SJuha Jäykkä radical and let everyone write all they've got */ 4030da24e51SJuha Jäykkä chunkDims[dim++] = (hsize_t) ceil(da->P*1.0/da->p); 4040da24e51SJuha Jäykkä chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n); 4050da24e51SJuha Jäykkä chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m); 4060da24e51SJuha Jäykkä } 4070da24e51SJuha Jäykkä } else { 408cb5c4f94SJuha Jäykkä /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */ 4090da24e51SJuha Jäykkä if (target_size >= chunk_size/da->n) { 410cb5c4f94SJuha Jäykkä /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */ 4110da24e51SJuha Jäykkä chunkDims[dim] = (hsize_t) ceil(da->N*1.0/da->n); 4120da24e51SJuha Jäykkä } else { 413cb5c4f94SJuha Jäykkä /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be 414cb5c4f94SJuha Jäykkä radical and let everyone write all they've got */ 4150da24e51SJuha Jäykkä chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n); 4160da24e51SJuha Jäykkä chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m); 4170da24e51SJuha Jäykkä } 4180da24e51SJuha Jäykkä 4190da24e51SJuha Jäykkä } 4200da24e51SJuha 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); 4210da24e51SJuha Jäykkä } else { 422cb5c4f94SJuha Jäykkä /* precomputed chunks are fine, we don't need to do anything */ 4230da24e51SJuha Jäykkä } 4240da24e51SJuha Jäykkä } 4250da24e51SJuha Jäykkä PetscFunctionReturn(0); 4260da24e51SJuha Jäykkä } 4270da24e51SJuha Jäykkä #endif 42847c6ae99SBarry Smith 42947c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 43047c6ae99SBarry Smith #undef __FUNCT__ 43147c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_HDF5_DA" 43247c6ae99SBarry Smith PetscErrorCode VecView_MPI_HDF5_DA(Vec xin,PetscViewer viewer) 43347c6ae99SBarry Smith { 4349b2a5a86SJed Brown DM dm; 4359b2a5a86SJed Brown DM_DA *da; 43647c6ae99SBarry Smith hid_t filespace; /* file dataspace identifier */ 4378e2ae6d7SMichael Kraus hid_t chunkspace; /* chunk dataset property identifier */ 43847c6ae99SBarry Smith hid_t plist_id; /* property list identifier */ 43947c6ae99SBarry Smith hid_t dset_id; /* dataset identifier */ 44047c6ae99SBarry Smith hid_t memspace; /* memory dataspace identifier */ 44147c6ae99SBarry Smith hid_t file_id; 44215214e8eSMatthew G Knepley hid_t group; 4438e2ae6d7SMichael Kraus hid_t scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */ 444d9a4edebSJed Brown hsize_t dim; 4450da24e51SJuha 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 */ 446c73cfb54SMatthew G. Knepley PetscInt timestep, dimension; 4475edff71fSBarry Smith const PetscScalar *x; 4488e2ae6d7SMichael Kraus const char *vecname; 44915214e8eSMatthew G Knepley PetscErrorCode ierr; 45047c6ae99SBarry Smith 45147c6ae99SBarry Smith PetscFunctionBegin; 45215214e8eSMatthew G Knepley ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr); 4538e2ae6d7SMichael Kraus ierr = PetscViewerHDF5GetTimestep(viewer, ×tep);CHKERRQ(ierr); 45415214e8eSMatthew G Knepley 455c688c046SMatthew G Knepley ierr = VecGetDM(xin,&dm);CHKERRQ(ierr); 456ce94432eSBarry Smith if (!dm) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 4579b2a5a86SJed Brown da = (DM_DA*)dm->data; 458c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dimension);CHKERRQ(ierr); 45947c6ae99SBarry Smith 4608e2ae6d7SMichael Kraus /* Create the dataspace for the dataset. 4618e2ae6d7SMichael Kraus * 4628e2ae6d7SMichael Kraus * dims - holds the current dimensions of the dataset 4638e2ae6d7SMichael Kraus * 4648e2ae6d7SMichael Kraus * maxDims - holds the maximum dimensions of the dataset (unlimited 4658e2ae6d7SMichael Kraus * for the number of time steps with the current dimensions for the 4668e2ae6d7SMichael Kraus * other dimensions; so only additional time steps can be added). 4678e2ae6d7SMichael Kraus * 4688e2ae6d7SMichael Kraus * chunkDims - holds the size of a single time step (required to 4698e2ae6d7SMichael Kraus * permit extending dataset). 4708e2ae6d7SMichael Kraus */ 4718e2ae6d7SMichael Kraus dim = 0; 4728e2ae6d7SMichael Kraus if (timestep >= 0) { 4738e2ae6d7SMichael Kraus dims[dim] = timestep+1; 4748e2ae6d7SMichael Kraus maxDims[dim] = H5S_UNLIMITED; 4758e2ae6d7SMichael Kraus chunkDims[dim] = 1; 4768e2ae6d7SMichael Kraus ++dim; 4778e2ae6d7SMichael Kraus } 478c73cfb54SMatthew G. Knepley if (dimension == 3) { 479acba2ac6SBarry Smith ierr = PetscHDF5IntCast(da->P,dims+dim);CHKERRQ(ierr); 4808e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4818e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4828e2ae6d7SMichael Kraus ++dim; 4838e2ae6d7SMichael Kraus } 484c73cfb54SMatthew G. Knepley if (dimension > 1) { 485acba2ac6SBarry Smith ierr = PetscHDF5IntCast(da->N,dims+dim);CHKERRQ(ierr); 4868e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4878e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4888e2ae6d7SMichael Kraus ++dim; 4898e2ae6d7SMichael Kraus } 490acba2ac6SBarry Smith ierr = PetscHDF5IntCast(da->M,dims+dim);CHKERRQ(ierr); 4918e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4928e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4938e2ae6d7SMichael Kraus ++dim; 4948e2ae6d7SMichael Kraus if (da->w > 1) { 495acba2ac6SBarry Smith ierr = PetscHDF5IntCast(da->w,dims+dim);CHKERRQ(ierr); 4968e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4978e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4988e2ae6d7SMichael Kraus ++dim; 4998e2ae6d7SMichael Kraus } 50047c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 5018e2ae6d7SMichael Kraus dims[dim] = 2; 5028e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 5038e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 5048e2ae6d7SMichael Kraus ++dim; 50547c6ae99SBarry Smith #endif 5060da24e51SJuha Jäykkä 507c73cfb54SMatthew G. Knepley ierr = VecGetHDF5ChunkSize(da, xin, dimension, timestep, chunkDims); CHKERRQ(ierr); 5080da24e51SJuha Jäykkä 509729ab6d9SBarry Smith PetscStackCallHDF5Return(filespace,H5Screate_simple,(dim, dims, maxDims)); 51047c6ae99SBarry Smith 51115214e8eSMatthew G Knepley #if defined(PETSC_USE_REAL_SINGLE) 51215214e8eSMatthew G Knepley scalartype = H5T_NATIVE_FLOAT; 51315214e8eSMatthew G Knepley #elif defined(PETSC_USE_REAL___FLOAT128) 51415214e8eSMatthew G Knepley #error "HDF5 output with 128 bit floats not supported." 51515214e8eSMatthew G Knepley #else 51615214e8eSMatthew G Knepley scalartype = H5T_NATIVE_DOUBLE; 51715214e8eSMatthew G Knepley #endif 51815214e8eSMatthew G Knepley 51947c6ae99SBarry Smith /* Create the dataset with default properties and close filespace */ 52047c6ae99SBarry Smith ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr); 52115214e8eSMatthew G Knepley if (!H5Lexists(group, vecname, H5P_DEFAULT)) { 5228e2ae6d7SMichael Kraus /* Create chunk */ 523729ab6d9SBarry Smith PetscStackCallHDF5Return(chunkspace,H5Pcreate,(H5P_DATASET_CREATE)); 524729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_chunk,(chunkspace, dim, chunkDims)); 5258e2ae6d7SMichael Kraus 52647c6ae99SBarry Smith #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 527729ab6d9SBarry Smith PetscStackCallHDF5Return(dset_id,H5Dcreate2,(group, vecname, scalartype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT)); 52847c6ae99SBarry Smith #else 529729ab6d9SBarry Smith PetscStackCallHDF5Return(dset_id,H5Dcreate,(group, vecname, scalartype, filespace, H5P_DEFAULT)); 53047c6ae99SBarry Smith #endif 53115214e8eSMatthew G Knepley } else { 532729ab6d9SBarry Smith PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT)); 533729ab6d9SBarry Smith PetscStackCallHDF5(H5Dset_extent,(dset_id, dims)); 53415214e8eSMatthew G Knepley } 535729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(filespace)); 53647c6ae99SBarry Smith 53747c6ae99SBarry Smith /* Each process defines a dataset and writes it to the hyperslab in the file */ 5388e2ae6d7SMichael Kraus dim = 0; 5398e2ae6d7SMichael Kraus if (timestep >= 0) { 5408e2ae6d7SMichael Kraus offset[dim] = timestep; 5418e2ae6d7SMichael Kraus ++dim; 5428e2ae6d7SMichael Kraus } 543c73cfb54SMatthew G. Knepley if (dimension == 3) {ierr = PetscHDF5IntCast(da->zs,offset + dim++);CHKERRQ(ierr);} 544c73cfb54SMatthew G. Knepley if (dimension > 1) {ierr = PetscHDF5IntCast(da->ys,offset + dim++);CHKERRQ(ierr);} 545acba2ac6SBarry Smith ierr = PetscHDF5IntCast(da->xs/da->w,offset + dim++);CHKERRQ(ierr); 5468e2ae6d7SMichael Kraus if (da->w > 1) offset[dim++] = 0; 54747c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 5488e2ae6d7SMichael Kraus offset[dim++] = 0; 54947c6ae99SBarry Smith #endif 5508e2ae6d7SMichael Kraus dim = 0; 5518e2ae6d7SMichael Kraus if (timestep >= 0) { 5528e2ae6d7SMichael Kraus count[dim] = 1; 5538e2ae6d7SMichael Kraus ++dim; 5548e2ae6d7SMichael Kraus } 555c73cfb54SMatthew G. Knepley if (dimension == 3) {ierr = PetscHDF5IntCast(da->ze - da->zs,count + dim++);CHKERRQ(ierr);} 556c73cfb54SMatthew G. Knepley if (dimension > 1) {ierr = PetscHDF5IntCast(da->ye - da->ys,count + dim++);CHKERRQ(ierr);} 557acba2ac6SBarry Smith ierr = PetscHDF5IntCast((da->xe - da->xs)/da->w,count + dim++);CHKERRQ(ierr); 558acba2ac6SBarry Smith if (da->w > 1) {ierr = PetscHDF5IntCast(da->w,count + dim++);CHKERRQ(ierr);} 55947c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 5608e2ae6d7SMichael Kraus count[dim++] = 2; 56147c6ae99SBarry Smith #endif 562729ab6d9SBarry Smith PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL)); 563729ab6d9SBarry Smith PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id)); 564729ab6d9SBarry Smith PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL)); 56547c6ae99SBarry Smith 56647c6ae99SBarry Smith /* Create property list for collective dataset write */ 567729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_DATASET_XFER)); 56847c6ae99SBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 569729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_dxpl_mpio,(plist_id, H5FD_MPIO_COLLECTIVE)); 57047c6ae99SBarry Smith #endif 57147c6ae99SBarry Smith /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 57247c6ae99SBarry Smith 5735edff71fSBarry Smith ierr = VecGetArrayRead(xin, &x);CHKERRQ(ierr); 574729ab6d9SBarry Smith PetscStackCallHDF5(H5Dwrite,(dset_id, scalartype, memspace, filespace, plist_id, x)); 575729ab6d9SBarry Smith PetscStackCallHDF5(H5Fflush,(file_id, H5F_SCOPE_GLOBAL)); 5765edff71fSBarry Smith ierr = VecRestoreArrayRead(xin, &x);CHKERRQ(ierr); 57747c6ae99SBarry Smith 57847c6ae99SBarry Smith /* Close/release resources */ 57915214e8eSMatthew G Knepley if (group != file_id) { 580729ab6d9SBarry Smith PetscStackCallHDF5(H5Gclose,(group)); 58115214e8eSMatthew G Knepley } 582729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 583729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(filespace)); 584729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(memspace)); 585729ab6d9SBarry Smith PetscStackCallHDF5(H5Dclose,(dset_id)); 58647c6ae99SBarry Smith ierr = PetscInfo1(xin,"Wrote Vec object with name %s\n",vecname);CHKERRQ(ierr); 58747c6ae99SBarry Smith PetscFunctionReturn(0); 58847c6ae99SBarry Smith } 58947c6ae99SBarry Smith #endif 59047c6ae99SBarry Smith 59109573ac7SBarry Smith extern PetscErrorCode VecView_MPI_Draw_DA1d(Vec,PetscViewer); 59247c6ae99SBarry Smith 59347c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 59447c6ae99SBarry Smith #undef __FUNCT__ 595aa219208SBarry Smith #define __FUNCT__ "DMDAArrayMPIIO" 596aa219208SBarry Smith static PetscErrorCode DMDAArrayMPIIO(DM da,PetscViewer viewer,Vec xin,PetscBool write) 59747c6ae99SBarry Smith { 59847c6ae99SBarry Smith PetscErrorCode ierr; 59947c6ae99SBarry Smith MPI_File mfdes; 60047c6ae99SBarry Smith PetscMPIInt gsizes[4],lsizes[4],lstarts[4],asiz,dof; 60147c6ae99SBarry Smith MPI_Datatype view; 60247c6ae99SBarry Smith const PetscScalar *array; 60347c6ae99SBarry Smith MPI_Offset off; 60447c6ae99SBarry Smith MPI_Aint ub,ul; 60547c6ae99SBarry Smith PetscInt type,rows,vecrows,tr[2]; 60647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 6070c169764Sdmay PetscBool skipheader; 60847c6ae99SBarry Smith 60947c6ae99SBarry Smith PetscFunctionBegin; 61047c6ae99SBarry Smith ierr = VecGetSize(xin,&vecrows);CHKERRQ(ierr); 6110c169764Sdmay ierr = PetscViewerBinaryGetSkipHeader(viewer,&skipheader);CHKERRQ(ierr); 61247c6ae99SBarry Smith if (!write) { 61347c6ae99SBarry Smith /* Read vector header. */ 6140c169764Sdmay if (!skipheader) { 615*060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,tr,2,NULL,PETSC_INT);CHKERRQ(ierr); 61647c6ae99SBarry Smith type = tr[0]; 61747c6ae99SBarry Smith rows = tr[1]; 618ce94432eSBarry Smith if (type != VEC_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Not vector next in file"); 619ce94432eSBarry Smith if (rows != vecrows) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Vector in file not same size as DMDA vector"); 6200c169764Sdmay } 62147c6ae99SBarry Smith } else { 62247c6ae99SBarry Smith tr[0] = VEC_FILE_CLASSID; 62347c6ae99SBarry Smith tr[1] = vecrows; 6240c169764Sdmay if (!skipheader) { 62547c6ae99SBarry Smith ierr = PetscViewerBinaryWrite(viewer,tr,2,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 62647c6ae99SBarry Smith } 6270c169764Sdmay } 62847c6ae99SBarry Smith 6294dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->w,&dof);CHKERRQ(ierr); 6304dc2109aSBarry Smith gsizes[0] = dof; 6314dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->M,gsizes+1);CHKERRQ(ierr); 6324dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->N,gsizes+2);CHKERRQ(ierr); 633334634e2SJed Brown ierr = PetscMPIIntCast(dd->P,gsizes+3);CHKERRQ(ierr); 6344dc2109aSBarry Smith lsizes[0] = dof; 6354dc2109aSBarry Smith ierr = PetscMPIIntCast((dd->xe-dd->xs)/dof,lsizes+1);CHKERRQ(ierr); 6364dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->ye-dd->ys,lsizes+2);CHKERRQ(ierr); 6374dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->ze-dd->zs,lsizes+3);CHKERRQ(ierr); 6384dc2109aSBarry Smith lstarts[0] = 0; 6394dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->xs/dof,lstarts+1);CHKERRQ(ierr); 6404dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->ys,lstarts+2);CHKERRQ(ierr); 6414dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->zs,lstarts+3);CHKERRQ(ierr); 642c73cfb54SMatthew G. Knepley ierr = MPI_Type_create_subarray(da->dim+1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view);CHKERRQ(ierr); 64347c6ae99SBarry Smith ierr = MPI_Type_commit(&view);CHKERRQ(ierr); 64447c6ae99SBarry Smith 64547c6ae99SBarry Smith ierr = PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);CHKERRQ(ierr); 64647c6ae99SBarry Smith ierr = PetscViewerBinaryGetMPIIOOffset(viewer,&off);CHKERRQ(ierr); 64747c6ae99SBarry Smith ierr = MPI_File_set_view(mfdes,off,MPIU_SCALAR,view,(char*)"native",MPI_INFO_NULL);CHKERRQ(ierr); 64847c6ae99SBarry Smith ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr); 64947c6ae99SBarry Smith asiz = lsizes[1]*(lsizes[2] > 0 ? lsizes[2] : 1)*(lsizes[3] > 0 ? lsizes[3] : 1)*dof; 65047c6ae99SBarry Smith if (write) { 65147c6ae99SBarry Smith ierr = MPIU_File_write_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr); 65247c6ae99SBarry Smith } else { 65347c6ae99SBarry Smith ierr = MPIU_File_read_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr); 65447c6ae99SBarry Smith } 65547c6ae99SBarry Smith ierr = MPI_Type_get_extent(view,&ul,&ub);CHKERRQ(ierr); 65647c6ae99SBarry Smith ierr = PetscViewerBinaryAddMPIIOOffset(viewer,ub);CHKERRQ(ierr); 65747c6ae99SBarry Smith ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr); 65847c6ae99SBarry Smith ierr = MPI_Type_free(&view);CHKERRQ(ierr); 65947c6ae99SBarry Smith PetscFunctionReturn(0); 66047c6ae99SBarry Smith } 66147c6ae99SBarry Smith #endif 66247c6ae99SBarry Smith 66347c6ae99SBarry Smith #undef __FUNCT__ 66447c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_DA" 6657087cfbeSBarry Smith PetscErrorCode VecView_MPI_DA(Vec xin,PetscViewer viewer) 66647c6ae99SBarry Smith { 6679a42bb27SBarry Smith DM da; 66847c6ae99SBarry Smith PetscErrorCode ierr; 66947c6ae99SBarry Smith PetscInt dim; 67047c6ae99SBarry Smith Vec natural; 6714061b8bfSJed Brown PetscBool isdraw,isvtk; 67247c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 67347c6ae99SBarry Smith PetscBool ishdf5; 67447c6ae99SBarry Smith #endif 6753f3fd955SJed Brown const char *prefix,*name; 676a261c58fSBarry Smith PetscViewerFormat format; 67747c6ae99SBarry Smith 67847c6ae99SBarry Smith PetscFunctionBegin; 679c688c046SMatthew G Knepley ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 680ce94432eSBarry Smith if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 681251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 682251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr); 68347c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 684251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr); 68547c6ae99SBarry Smith #endif 68647c6ae99SBarry Smith if (isdraw) { 6871321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 68847c6ae99SBarry Smith if (dim == 1) { 68947c6ae99SBarry Smith ierr = VecView_MPI_Draw_DA1d(xin,viewer);CHKERRQ(ierr); 69047c6ae99SBarry Smith } else if (dim == 2) { 69147c6ae99SBarry Smith ierr = VecView_MPI_Draw_DA2d(xin,viewer);CHKERRQ(ierr); 692ce94432eSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot graphically view vector associated with this dimensional DMDA %D",dim); 6934061b8bfSJed Brown } else if (isvtk) { /* Duplicate the Vec and hold a reference to the DM */ 6944061b8bfSJed Brown Vec Y; 6954061b8bfSJed Brown ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr); 6964061b8bfSJed Brown ierr = VecDuplicate(xin,&Y);CHKERRQ(ierr); 697b51b94faSJed Brown if (((PetscObject)xin)->name) { 698b51b94faSJed Brown /* If xin was named, copy the name over to Y. The duplicate names are safe because nobody else will ever see Y. */ 699b51b94faSJed Brown ierr = PetscObjectSetName((PetscObject)Y,((PetscObject)xin)->name);CHKERRQ(ierr); 700b51b94faSJed Brown } 7014061b8bfSJed Brown ierr = VecCopy(xin,Y);CHKERRQ(ierr); 70262b69a3fSMatthew G Knepley ierr = PetscViewerVTKAddField(viewer,(PetscObject)da,DMDAVTKWriteAll,PETSC_VTK_POINT_FIELD,(PetscObject)Y);CHKERRQ(ierr); 70347c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 70447c6ae99SBarry Smith } else if (ishdf5) { 70547c6ae99SBarry Smith ierr = VecView_MPI_HDF5_DA(xin,viewer);CHKERRQ(ierr); 70647c6ae99SBarry Smith #endif 70747c6ae99SBarry Smith } else { 70847c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 70947c6ae99SBarry Smith PetscBool isbinary,isMPIIO; 71047c6ae99SBarry Smith 711251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 71247c6ae99SBarry Smith if (isbinary) { 713bc196f7cSDave May ierr = PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO);CHKERRQ(ierr); 71447c6ae99SBarry Smith if (isMPIIO) { 715aa219208SBarry Smith ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_TRUE);CHKERRQ(ierr); 71647c6ae99SBarry Smith PetscFunctionReturn(0); 71747c6ae99SBarry Smith } 71847c6ae99SBarry Smith } 71947c6ae99SBarry Smith #endif 72047c6ae99SBarry Smith 72147c6ae99SBarry Smith /* call viewer on natural ordering */ 72247c6ae99SBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr); 723aa219208SBarry Smith ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr); 72447c6ae99SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr); 725aa219208SBarry Smith ierr = DMDAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr); 726aa219208SBarry Smith ierr = DMDAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr); 7273f3fd955SJed Brown ierr = PetscObjectGetName((PetscObject)xin,&name);CHKERRQ(ierr); 7283f3fd955SJed Brown ierr = PetscObjectSetName((PetscObject)natural,name);CHKERRQ(ierr); 729a261c58fSBarry Smith 730a261c58fSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 731a261c58fSBarry Smith if (format == PETSC_VIEWER_BINARY_MATLAB) { 732a261c58fSBarry Smith /* temporarily remove viewer format so it won't trigger in the VecView() */ 733a261c58fSBarry Smith ierr = PetscViewerSetFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr); 734a261c58fSBarry Smith } 735a261c58fSBarry Smith 73647c6ae99SBarry Smith ierr = VecView(natural,viewer);CHKERRQ(ierr); 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 */ 746a261c58fSBarry Smith ierr = PetscViewerSetFormat(viewer,format);CHKERRQ(ierr); 747a261c58fSBarry Smith ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); 748a261c58fSBarry Smith ierr = PetscViewerBinaryGetInfoPointer(viewer,&info);CHKERRQ(ierr); 749a261c58fSBarry Smith ierr = DMDAGetInfo(da,&dim,&ni,&nj,&nk,&pi,&pj,&pk,&dof,0,0,0,0,0);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 #undef __FUNCT__ 77747c6ae99SBarry Smith #define __FUNCT__ "VecLoad_HDF5_DA" 77847c6ae99SBarry Smith PetscErrorCode VecLoad_HDF5_DA(Vec xin, PetscViewer viewer) 77947c6ae99SBarry Smith { 7809a42bb27SBarry Smith DM da; 78147c6ae99SBarry Smith PetscErrorCode ierr; 78225578ef6SJed Brown hsize_t dim; 78347c6ae99SBarry Smith hsize_t count[5]; 78447c6ae99SBarry Smith hsize_t offset[5]; 785c73cfb54SMatthew G. Knepley PetscInt cnt = 0, dimension; 78647c6ae99SBarry Smith PetscScalar *x; 78747c6ae99SBarry Smith const char *vecname; 78847c6ae99SBarry Smith hid_t filespace; /* file dataspace identifier */ 78947c6ae99SBarry Smith hid_t plist_id; /* property list identifier */ 79047c6ae99SBarry Smith hid_t dset_id; /* dataset identifier */ 79147c6ae99SBarry Smith hid_t memspace; /* memory dataspace identifier */ 792bfd264e7SBarry Smith hid_t file_id,group; 7937bcbaff4SBarry Smith hid_t scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */ 7949c7c4993SBarry Smith DM_DA *dd; 79547c6ae99SBarry Smith 79647c6ae99SBarry Smith PetscFunctionBegin; 7977bcbaff4SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 7987bcbaff4SBarry Smith scalartype = H5T_NATIVE_FLOAT; 7997bcbaff4SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128) 8007bcbaff4SBarry Smith #error "HDF5 output with 128 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); 806c688c046SMatthew G Knepley ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 8079c7c4993SBarry Smith dd = (DM_DA*)da->data; 808c73cfb54SMatthew G. Knepley ierr = DMGetDimension(da, &dimension);CHKERRQ(ierr); 80947c6ae99SBarry Smith 81047c6ae99SBarry Smith /* Create the dataspace for the dataset */ 811c73cfb54SMatthew G. Knepley ierr = PetscHDF5IntCast(dimension + ((dd->w == 1) ? 0 : 1),&dim);CHKERRQ(ierr); 81247c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 81347c6ae99SBarry Smith dim++; 81447c6ae99SBarry Smith #endif 81547c6ae99SBarry Smith 81647c6ae99SBarry Smith /* Create the dataset with default properties and close filespace */ 81747c6ae99SBarry Smith ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr); 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 823729ab6d9SBarry Smith PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id)); 82447c6ae99SBarry Smith 82547c6ae99SBarry Smith /* Each process defines a dataset and reads it from the hyperslab in the file */ 82647c6ae99SBarry Smith cnt = 0; 827c73cfb54SMatthew G. Knepley if (dimension == 3) {ierr = PetscHDF5IntCast(dd->zs,offset + cnt++);CHKERRQ(ierr);} 828c73cfb54SMatthew G. Knepley if (dimension > 1) {ierr = PetscHDF5IntCast(dd->ys,offset + cnt++);CHKERRQ(ierr);} 829acba2ac6SBarry Smith ierr = PetscHDF5IntCast(dd->xs/dd->w,offset + cnt++);CHKERRQ(ierr); 83047c6ae99SBarry Smith if (dd->w > 1) offset[cnt++] = 0; 83147c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 83247c6ae99SBarry Smith offset[cnt++] = 0; 83347c6ae99SBarry Smith #endif 83447c6ae99SBarry Smith cnt = 0; 835c73cfb54SMatthew G. Knepley if (dimension == 3) {ierr = PetscHDF5IntCast(dd->ze - dd->zs,count + cnt++);CHKERRQ(ierr);} 836c73cfb54SMatthew G. Knepley if (dimension > 1) {ierr = PetscHDF5IntCast(dd->ye - dd->ys,count + cnt++);CHKERRQ(ierr);} 837acba2ac6SBarry Smith ierr = PetscHDF5IntCast((dd->xe - dd->xs)/dd->w,count + cnt++);CHKERRQ(ierr); 838acba2ac6SBarry Smith if (dd->w > 1) {ierr = PetscHDF5IntCast(dd->w,count + cnt++);CHKERRQ(ierr);} 83947c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 84047c6ae99SBarry Smith count[cnt++] = 2; 84147c6ae99SBarry Smith #endif 842729ab6d9SBarry Smith PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL)); 843729ab6d9SBarry Smith PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL)); 84447c6ae99SBarry Smith 84547c6ae99SBarry Smith /* Create property list for collective dataset write */ 846729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_DATASET_XFER)); 84747c6ae99SBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 848729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_dxpl_mpio,(plist_id, H5FD_MPIO_COLLECTIVE)); 84947c6ae99SBarry Smith #endif 85047c6ae99SBarry Smith /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 85147c6ae99SBarry Smith 85247c6ae99SBarry Smith ierr = VecGetArray(xin, &x);CHKERRQ(ierr); 853a1dbdba5SBarry Smith PetscStackCallHDF5(H5Dread,(dset_id, scalartype, memspace, filespace, plist_id, x)); 85447c6ae99SBarry Smith ierr = VecRestoreArray(xin, &x);CHKERRQ(ierr); 85547c6ae99SBarry Smith 85647c6ae99SBarry Smith /* Close/release resources */ 857bfd264e7SBarry Smith if (group != file_id) { 858729ab6d9SBarry Smith PetscStackCallHDF5(H5Gclose,(group)); 859bfd264e7SBarry Smith } 860729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 861729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(filespace)); 862729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(memspace)); 863729ab6d9SBarry Smith PetscStackCallHDF5(H5Dclose,(dset_id)); 86447c6ae99SBarry Smith PetscFunctionReturn(0); 86547c6ae99SBarry Smith } 86647c6ae99SBarry Smith #endif 86747c6ae99SBarry Smith 86847c6ae99SBarry Smith #undef __FUNCT__ 86947c6ae99SBarry Smith #define __FUNCT__ "VecLoad_Binary_DA" 87047c6ae99SBarry Smith PetscErrorCode VecLoad_Binary_DA(Vec xin, PetscViewer viewer) 87147c6ae99SBarry Smith { 8729a42bb27SBarry Smith DM da; 87347c6ae99SBarry Smith PetscErrorCode ierr; 87447c6ae99SBarry Smith Vec natural; 87547c6ae99SBarry Smith const char *prefix; 87647c6ae99SBarry Smith PetscInt bs; 87747c6ae99SBarry Smith PetscBool flag; 87847c6ae99SBarry Smith DM_DA *dd; 87947c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 88047c6ae99SBarry Smith PetscBool isMPIIO; 88147c6ae99SBarry Smith #endif 88247c6ae99SBarry Smith 88347c6ae99SBarry Smith PetscFunctionBegin; 884c688c046SMatthew G Knepley ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 88547c6ae99SBarry Smith dd = (DM_DA*)da->data; 88647c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 887bc196f7cSDave May ierr = PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO);CHKERRQ(ierr); 88847c6ae99SBarry Smith if (isMPIIO) { 889aa219208SBarry Smith ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_FALSE);CHKERRQ(ierr); 89047c6ae99SBarry Smith PetscFunctionReturn(0); 89147c6ae99SBarry Smith } 89247c6ae99SBarry Smith #endif 89347c6ae99SBarry Smith 89447c6ae99SBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr); 895aa219208SBarry Smith ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr); 89647c6ae99SBarry Smith ierr = PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name);CHKERRQ(ierr); 89747c6ae99SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr); 8988aea9937SBarry Smith ierr = VecLoad(natural,viewer);CHKERRQ(ierr); 899aa219208SBarry Smith ierr = DMDANaturalToGlobalBegin(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr); 900aa219208SBarry Smith ierr = DMDANaturalToGlobalEnd(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr); 901fcfd50ebSBarry Smith ierr = VecDestroy(&natural);CHKERRQ(ierr); 902aa219208SBarry Smith ierr = PetscInfo(xin,"Loading vector from natural ordering into DMDA\n");CHKERRQ(ierr); 90347c6ae99SBarry Smith ierr = PetscOptionsGetInt(((PetscObject)xin)->prefix,"-vecload_block_size",&bs,&flag);CHKERRQ(ierr); 90447c6ae99SBarry Smith if (flag && bs != dd->w) { 905aa219208SBarry Smith ierr = PetscInfo2(xin,"Block size in file %D not equal to DMDA's dof %D\n",bs,dd->w);CHKERRQ(ierr); 90647c6ae99SBarry Smith } 90747c6ae99SBarry Smith PetscFunctionReturn(0); 90847c6ae99SBarry Smith } 90947c6ae99SBarry Smith 91047c6ae99SBarry Smith #undef __FUNCT__ 91147c6ae99SBarry Smith #define __FUNCT__ "VecLoad_Default_DA" 9127087cfbeSBarry Smith PetscErrorCode VecLoad_Default_DA(Vec xin, PetscViewer viewer) 91347c6ae99SBarry Smith { 91447c6ae99SBarry Smith PetscErrorCode ierr; 9159a42bb27SBarry Smith DM da; 91647c6ae99SBarry Smith PetscBool isbinary; 91747c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 91847c6ae99SBarry Smith PetscBool ishdf5; 91947c6ae99SBarry Smith #endif 92047c6ae99SBarry Smith 92147c6ae99SBarry Smith PetscFunctionBegin; 922c688c046SMatthew G Knepley ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 923ce94432eSBarry Smith if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 92447c6ae99SBarry Smith 92547c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 926251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr); 92747c6ae99SBarry Smith #endif 928251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 92947c6ae99SBarry Smith 93047c6ae99SBarry Smith if (isbinary) { 93147c6ae99SBarry Smith ierr = VecLoad_Binary_DA(xin,viewer);CHKERRQ(ierr); 93247c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 93347c6ae99SBarry Smith } else if (ishdf5) { 93447c6ae99SBarry Smith ierr = VecLoad_HDF5_DA(xin,viewer);CHKERRQ(ierr); 93547c6ae99SBarry Smith #endif 936d34fcf5fSBarry Smith } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for vector loading", ((PetscObject)viewer)->type_name); 93747c6ae99SBarry Smith PetscFunctionReturn(0); 93847c6ae99SBarry Smith } 939