147c6ae99SBarry Smith 247c6ae99SBarry Smith /* 3aa219208SBarry Smith Plots vectors obtained with DMDACreate2d() 447c6ae99SBarry Smith */ 547c6ae99SBarry Smith 6af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 7af0996ceSBarry Smith #include <petsc/private/vecimpl.h> 89804daf3SBarry Smith #include <petscdraw.h> 9665c2dedSJed Brown #include <petscviewerhdf5.h> 1047c6ae99SBarry Smith 1147c6ae99SBarry Smith /* 1247c6ae99SBarry Smith The data that is passed into the graphics callback 1347c6ae99SBarry Smith */ 1447c6ae99SBarry Smith typedef struct { 15e5ab1681SLisandro Dalcin PetscMPIInt rank; 16e5ab1681SLisandro Dalcin PetscInt m,n,dof,k; 17e5ab1681SLisandro Dalcin PetscReal xmin,xmax,ymin,ymax,min,max; 185edff71fSBarry Smith const PetscScalar *xy,*v; 1947c6ae99SBarry Smith PetscBool showgrid; 20109c9344SBarry Smith const char *name0,*name1; 2147c6ae99SBarry Smith } ZoomCtx; 2247c6ae99SBarry Smith 2347c6ae99SBarry Smith /* 2447c6ae99SBarry Smith This does the drawing for one particular field 2547c6ae99SBarry Smith in one particular set of coordinates. It is a callback 2647c6ae99SBarry Smith called from PetscDrawZoom() 2747c6ae99SBarry Smith */ 2847c6ae99SBarry Smith #undef __FUNCT__ 2947c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_Draw_DA2d_Zoom" 3047c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA2d_Zoom(PetscDraw draw,void *ctx) 3147c6ae99SBarry Smith { 3247c6ae99SBarry Smith ZoomCtx *zctx = (ZoomCtx*)ctx; 3347c6ae99SBarry Smith PetscErrorCode ierr; 34e5ab1681SLisandro Dalcin PetscInt m,n,i,j,k,dof,id,c1,c2,c3,c4; 35e5ab1681SLisandro Dalcin PetscReal min,max,x1,x2,x3,x4,y_1,y2,y3,y4; 36e5ab1681SLisandro Dalcin const PetscScalar *xy,*v; 3747c6ae99SBarry Smith 3847c6ae99SBarry Smith PetscFunctionBegin; 3947c6ae99SBarry Smith m = zctx->m; 4047c6ae99SBarry Smith n = zctx->n; 41e5ab1681SLisandro Dalcin dof = zctx->dof; 4247c6ae99SBarry Smith k = zctx->k; 4347c6ae99SBarry Smith xy = zctx->xy; 44e5ab1681SLisandro Dalcin v = zctx->v; 4547c6ae99SBarry Smith min = zctx->min; 46f3f0eb19SBarry Smith max = zctx->max; 4747c6ae99SBarry Smith 4847c6ae99SBarry Smith /* PetscDraw the contour plot patch */ 49e5ab1681SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 5047c6ae99SBarry Smith for (j=0; j<n-1; j++) { 5147c6ae99SBarry Smith for (i=0; i<m-1; i++) { 520e4fe250SBarry Smith id = i+j*m; 530e4fe250SBarry Smith x1 = PetscRealPart(xy[2*id]); 540e4fe250SBarry Smith y_1 = PetscRealPart(xy[2*id+1]); 55e5ab1681SLisandro Dalcin c1 = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max); 560e4fe250SBarry Smith 570e4fe250SBarry Smith id = i+j*m+1; 580e4fe250SBarry Smith x2 = PetscRealPart(xy[2*id]); 59a39a4c7dSLisandro Dalcin y2 = PetscRealPart(xy[2*id+1]); 60e5ab1681SLisandro Dalcin c2 = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max); 610e4fe250SBarry Smith 620e4fe250SBarry Smith id = i+j*m+1+m; 63a39a4c7dSLisandro Dalcin x3 = PetscRealPart(xy[2*id]); 640e4fe250SBarry Smith y3 = PetscRealPart(xy[2*id+1]); 65e5ab1681SLisandro Dalcin c3 = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max); 660e4fe250SBarry Smith 670e4fe250SBarry Smith id = i+j*m+m; 68a39a4c7dSLisandro Dalcin x4 = PetscRealPart(xy[2*id]); 69a39a4c7dSLisandro Dalcin y4 = PetscRealPart(xy[2*id+1]); 70e5ab1681SLisandro Dalcin c4 = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max); 71f3f0eb19SBarry Smith 7247c6ae99SBarry Smith ierr = PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);CHKERRQ(ierr); 7347c6ae99SBarry Smith ierr = PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4);CHKERRQ(ierr); 7447c6ae99SBarry Smith if (zctx->showgrid) { 7547c6ae99SBarry Smith ierr = PetscDrawLine(draw,x1,y_1,x2,y2,PETSC_DRAW_BLACK);CHKERRQ(ierr); 7647c6ae99SBarry Smith ierr = PetscDrawLine(draw,x2,y2,x3,y3,PETSC_DRAW_BLACK);CHKERRQ(ierr); 7747c6ae99SBarry Smith ierr = PetscDrawLine(draw,x3,y3,x4,y4,PETSC_DRAW_BLACK);CHKERRQ(ierr); 7847c6ae99SBarry Smith ierr = PetscDrawLine(draw,x4,y4,x1,y_1,PETSC_DRAW_BLACK);CHKERRQ(ierr); 7947c6ae99SBarry Smith } 8047c6ae99SBarry Smith } 8147c6ae99SBarry Smith } 82e5ab1681SLisandro Dalcin if (!zctx->rank) { 83e5ab1681SLisandro Dalcin if (zctx->name0 || zctx->name1) { 84109c9344SBarry Smith PetscReal xl,yl,xr,yr,x,y; 85109c9344SBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 86e5ab1681SLisandro Dalcin x = xl + .30*(xr - xl); 87109c9344SBarry Smith xl = xl + .01*(xr - xl); 88e5ab1681SLisandro Dalcin y = yr - .30*(yr - yl); 89109c9344SBarry Smith yl = yl + .01*(yr - yl); 90e5ab1681SLisandro Dalcin if (zctx->name0) {ierr = PetscDrawString(draw,x,yl,PETSC_DRAW_BLACK,zctx->name0);CHKERRQ(ierr);} 91e5ab1681SLisandro Dalcin if (zctx->name1) {ierr = PetscDrawStringVertical(draw,xl,y,PETSC_DRAW_BLACK,zctx->name1);CHKERRQ(ierr);} 92109c9344SBarry Smith } 930e4fe250SBarry Smith /* 940e4fe250SBarry Smith Ideally we would use the PetscDrawAxis object to manage displaying the coordinate limits 950e4fe250SBarry Smith but that may require some refactoring. 960e4fe250SBarry Smith */ 97e5ab1681SLisandro Dalcin { 98e5ab1681SLisandro Dalcin PetscReal xmin = zctx->xmin, ymin = zctx->ymin; 99e5ab1681SLisandro Dalcin PetscReal xmax = zctx->xmax, ymax = zctx->ymax; 100e5ab1681SLisandro Dalcin char value[16]; size_t len; PetscReal w; 101e5ab1681SLisandro Dalcin ierr = PetscSNPrintf(value,16,"%f",xmin);CHKERRQ(ierr); 102e5ab1681SLisandro Dalcin ierr = PetscDrawString(draw,xmin,ymin - .05*(ymax - ymin),PETSC_DRAW_BLACK,value);CHKERRQ(ierr); 103e5ab1681SLisandro Dalcin ierr = PetscSNPrintf(value,16,"%f",xmax);CHKERRQ(ierr); 1040e4fe250SBarry Smith ierr = PetscStrlen(value,&len);CHKERRQ(ierr); 1050298fd71SBarry Smith ierr = PetscDrawStringGetSize(draw,&w,NULL);CHKERRQ(ierr); 106e5ab1681SLisandro Dalcin ierr = PetscDrawString(draw,xmax - len*w,ymin - .05*(ymax - ymin),PETSC_DRAW_BLACK,value);CHKERRQ(ierr); 107e5ab1681SLisandro Dalcin ierr = PetscSNPrintf(value,16,"%f",ymin);CHKERRQ(ierr); 108e5ab1681SLisandro Dalcin ierr = PetscDrawString(draw,xmin - .05*(xmax - xmin),ymin,PETSC_DRAW_BLACK,value);CHKERRQ(ierr); 109e5ab1681SLisandro Dalcin ierr = PetscSNPrintf(value,16,"%f",ymax);CHKERRQ(ierr); 110e5ab1681SLisandro Dalcin ierr = PetscDrawString(draw,xmin - .05*(xmax - xmin),ymax,PETSC_DRAW_BLACK,value);CHKERRQ(ierr); 111e5ab1681SLisandro Dalcin } 112e5ab1681SLisandro Dalcin } 113e5ab1681SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 11447c6ae99SBarry Smith PetscFunctionReturn(0); 11547c6ae99SBarry Smith } 11647c6ae99SBarry Smith 11747c6ae99SBarry Smith #undef __FUNCT__ 11847c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_Draw_DA2d" 11947c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA2d(Vec xin,PetscViewer viewer) 12047c6ae99SBarry Smith { 1219a42bb27SBarry Smith DM da,dac,dag; 12247c6ae99SBarry Smith PetscErrorCode ierr; 123a74ba6f7SBarry Smith PetscInt N,s,M,w,ncoors = 4; 12447c6ae99SBarry Smith const PetscInt *lx,*ly; 125e5ab1681SLisandro Dalcin PetscReal coors[4]; 12647c6ae99SBarry Smith PetscDraw draw,popup; 12747c6ae99SBarry Smith PetscBool isnull,useports = PETSC_FALSE; 12847c6ae99SBarry Smith MPI_Comm comm; 12947c6ae99SBarry Smith Vec xlocal,xcoor,xcoorl; 130bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by; 131aa219208SBarry Smith DMDAStencilType st; 13247c6ae99SBarry Smith ZoomCtx zctx; 1330298fd71SBarry Smith PetscDrawViewPorts *ports = NULL; 13447c6ae99SBarry Smith PetscViewerFormat format; 13520d0051dSBarry Smith PetscInt *displayfields; 13667dd0837SBarry Smith PetscInt ndisplayfields,i,nbounds; 13767dd0837SBarry Smith const PetscReal *bounds; 13847c6ae99SBarry Smith 13947c6ae99SBarry Smith PetscFunctionBegin; 14047c6ae99SBarry Smith zctx.showgrid = PETSC_FALSE; 1418865f1eaSKarl Rupp 14247c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 1435b399a63SLisandro Dalcin ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 1445b399a63SLisandro Dalcin if (isnull) PetscFunctionReturn(0); 1455b399a63SLisandro Dalcin 14603193ff8SBarry Smith ierr = PetscViewerDrawGetBounds(viewer,&nbounds,&bounds);CHKERRQ(ierr); 14747c6ae99SBarry Smith 148c688c046SMatthew G Knepley ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 149ce94432eSBarry Smith if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 15047c6ae99SBarry Smith 15147c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr); 152e5ab1681SLisandro Dalcin ierr = MPI_Comm_rank(comm,&zctx.rank);CHKERRQ(ierr); 15347c6ae99SBarry Smith 1541321219cSEthan Coon ierr = DMDAGetInfo(da,0,&M,&N,0,&zctx.m,&zctx.n,0,&w,&s,&bx,&by,0,&st);CHKERRQ(ierr); 1550298fd71SBarry Smith ierr = DMDAGetOwnershipRanges(da,&lx,&ly,NULL);CHKERRQ(ierr); 15647c6ae99SBarry Smith 15747c6ae99SBarry Smith /* 15847c6ae99SBarry Smith Obtain a sequential vector that is going to contain the local values plus ONE layer of 159aa219208SBarry Smith ghosted values to draw the graphics from. We also need its corresponding DMDA (dac) that will 16047c6ae99SBarry Smith update the local values pluse ONE layer of ghost values. 16147c6ae99SBarry Smith */ 16247c6ae99SBarry Smith ierr = PetscObjectQuery((PetscObject)da,"GraphicsGhosted",(PetscObject*)&xlocal);CHKERRQ(ierr); 16347c6ae99SBarry Smith if (!xlocal) { 164bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_NONE || by != DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) { 16547c6ae99SBarry Smith /* 16647c6ae99SBarry Smith if original da is not of stencil width one, or periodic or not a box stencil then 167aa219208SBarry Smith create a special DMDA to handle one level of ghost points for graphics 16847c6ae99SBarry Smith */ 169bff4a2f0SMatthew G. Knepley ierr = DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,zctx.m,zctx.n,w,1,lx,ly,&dac);CHKERRQ(ierr); 170aa219208SBarry Smith ierr = PetscInfo(da,"Creating auxilary DMDA for managing graphics ghost points\n");CHKERRQ(ierr); 17147c6ae99SBarry Smith } else { 17247c6ae99SBarry Smith /* otherwise we can use the da we already have */ 17347c6ae99SBarry Smith dac = da; 17447c6ae99SBarry Smith } 17547c6ae99SBarry Smith /* create local vector for holding ghosted values used in graphics */ 176564755cdSBarry Smith ierr = DMCreateLocalVector(dac,&xlocal);CHKERRQ(ierr); 17747c6ae99SBarry Smith if (dac != da) { 178aa219208SBarry Smith /* don't keep any public reference of this DMDA, is is only available through xlocal */ 179f7923d8aSBarry Smith ierr = PetscObjectDereference((PetscObject)dac);CHKERRQ(ierr); 18047c6ae99SBarry Smith } else { 18147c6ae99SBarry Smith /* remove association between xlocal and da, because below we compose in the opposite 18247c6ae99SBarry Smith direction and if we left this connect we'd get a loop, so the objects could 18347c6ae99SBarry Smith never be destroyed */ 184c688c046SMatthew G Knepley ierr = PetscObjectRemoveReference((PetscObject)xlocal,"__PETSc_dm");CHKERRQ(ierr); 18547c6ae99SBarry Smith } 18647c6ae99SBarry Smith ierr = PetscObjectCompose((PetscObject)da,"GraphicsGhosted",(PetscObject)xlocal);CHKERRQ(ierr); 18747c6ae99SBarry Smith ierr = PetscObjectDereference((PetscObject)xlocal);CHKERRQ(ierr); 18847c6ae99SBarry Smith } else { 189bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_NONE || by != DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) { 190c688c046SMatthew G Knepley ierr = VecGetDM(xlocal, &dac);CHKERRQ(ierr); 191f7923d8aSBarry Smith } else { 192f7923d8aSBarry Smith dac = da; 19347c6ae99SBarry Smith } 19447c6ae99SBarry Smith } 19547c6ae99SBarry Smith 19647c6ae99SBarry Smith /* 19747c6ae99SBarry Smith Get local (ghosted) values of vector 19847c6ae99SBarry Smith */ 1999a42bb27SBarry Smith ierr = DMGlobalToLocalBegin(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr); 2009a42bb27SBarry Smith ierr = DMGlobalToLocalEnd(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr); 2015edff71fSBarry Smith ierr = VecGetArrayRead(xlocal,&zctx.v);CHKERRQ(ierr); 20247c6ae99SBarry Smith 203*832b7cebSLisandro Dalcin /* 204*832b7cebSLisandro Dalcin Get coordinates of nodes 205*832b7cebSLisandro Dalcin */ 2066636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); 20747c6ae99SBarry Smith if (!xcoor) { 208aa219208SBarry Smith ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,0.0);CHKERRQ(ierr); 2096636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); 21047c6ae99SBarry Smith } 21147c6ae99SBarry Smith 21247c6ae99SBarry Smith /* 21347c6ae99SBarry Smith Determine the min and max coordinates in plot 21447c6ae99SBarry Smith */ 215e5ab1681SLisandro Dalcin ierr = VecStrideMin(xcoor,0,NULL,&zctx.xmin);CHKERRQ(ierr); 216e5ab1681SLisandro Dalcin ierr = VecStrideMax(xcoor,0,NULL,&zctx.xmax);CHKERRQ(ierr); 217e5ab1681SLisandro Dalcin ierr = VecStrideMin(xcoor,1,NULL,&zctx.ymin);CHKERRQ(ierr); 218e5ab1681SLisandro Dalcin ierr = VecStrideMax(xcoor,1,NULL,&zctx.ymax);CHKERRQ(ierr); 219e5ab1681SLisandro Dalcin coors[0] = zctx.xmin - .05*(zctx.xmax- zctx.xmin); coors[2] = zctx.xmax + .05*(zctx.xmax - zctx.xmin); 220e5ab1681SLisandro Dalcin coors[1] = zctx.ymin - .05*(zctx.ymax- zctx.ymin); coors[3] = zctx.ymax + .05*(zctx.ymax - zctx.ymin); 221c5929fdfSBarry Smith ierr = PetscOptionsGetRealArray(NULL,NULL,"-draw_coordinates",coors,&ncoors,NULL);CHKERRQ(ierr); 222e5ab1681SLisandro 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); 223a74ba6f7SBarry Smith 22447c6ae99SBarry Smith /* 225*832b7cebSLisandro Dalcin Get local ghosted version of coordinates 22647c6ae99SBarry Smith */ 22747c6ae99SBarry Smith ierr = PetscObjectQuery((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr); 22847c6ae99SBarry Smith if (!xcoorl) { 229aa219208SBarry Smith /* create DMDA to get local version of graphics */ 230bff4a2f0SMatthew 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); 231aa219208SBarry Smith ierr = PetscInfo(dag,"Creating auxilary DMDA for managing graphics coordinates ghost points\n");CHKERRQ(ierr); 232564755cdSBarry Smith ierr = DMCreateLocalVector(dag,&xcoorl);CHKERRQ(ierr); 23347c6ae99SBarry Smith ierr = PetscObjectCompose((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject)xcoorl);CHKERRQ(ierr); 234f7923d8aSBarry Smith ierr = PetscObjectDereference((PetscObject)dag);CHKERRQ(ierr); 23547c6ae99SBarry Smith ierr = PetscObjectDereference((PetscObject)xcoorl);CHKERRQ(ierr); 23647c6ae99SBarry Smith } else { 237c688c046SMatthew G Knepley ierr = VecGetDM(xcoorl,&dag);CHKERRQ(ierr); 23847c6ae99SBarry Smith } 2399a42bb27SBarry Smith ierr = DMGlobalToLocalBegin(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr); 2409a42bb27SBarry Smith ierr = DMGlobalToLocalEnd(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr); 2415edff71fSBarry Smith ierr = VecGetArrayRead(xcoorl,&zctx.xy);CHKERRQ(ierr); 242*832b7cebSLisandro Dalcin ierr = DMDAGetCoordinateName(da,0,&zctx.name0);CHKERRQ(ierr); 243*832b7cebSLisandro Dalcin ierr = DMDAGetCoordinateName(da,1,&zctx.name1);CHKERRQ(ierr); 24447c6ae99SBarry Smith 24547c6ae99SBarry Smith /* 24647c6ae99SBarry Smith Get information about size of area each processor must do graphics for 24747c6ae99SBarry Smith */ 248e5ab1681SLisandro Dalcin ierr = DMDAGetInfo(dac,NULL,&M,&N,NULL,NULL,NULL,NULL,&zctx.dof,NULL,&bx,&by,NULL,NULL);CHKERRQ(ierr); 249e5ab1681SLisandro Dalcin ierr = DMDAGetGhostCorners(dac,NULL,NULL,NULL,&zctx.m,&zctx.n,NULL);CHKERRQ(ierr); 250c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-draw_contour_grid",&zctx.showgrid,NULL);CHKERRQ(ierr); 2514e6118eeSBarry Smith 2524e6118eeSBarry Smith ierr = DMDASelectFields(da,&ndisplayfields,&displayfields);CHKERRQ(ierr); 25347c6ae99SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 254c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-draw_ports",&useports,NULL);CHKERRQ(ierr); 255*832b7cebSLisandro Dalcin if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE; 256*832b7cebSLisandro Dalcin if (useports) { 257e5ab1681SLisandro Dalcin ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 2585b399a63SLisandro Dalcin ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr); 2595b399a63SLisandro Dalcin ierr = PetscDrawClear(draw);CHKERRQ(ierr); 26020d0051dSBarry Smith ierr = PetscDrawViewPortsCreate(draw,ndisplayfields,&ports);CHKERRQ(ierr); 261e5ab1681SLisandro Dalcin } 26220d0051dSBarry Smith 263*832b7cebSLisandro Dalcin /* 264*832b7cebSLisandro Dalcin Loop over each field; drawing each in a different window 265*832b7cebSLisandro Dalcin */ 26620d0051dSBarry Smith for (i=0; i<ndisplayfields; i++) { 26720d0051dSBarry Smith zctx.k = displayfields[i]; 2685b399a63SLisandro Dalcin 269e5ab1681SLisandro Dalcin /* determine the min and max value in plot */ 2700298fd71SBarry Smith ierr = VecStrideMin(xin,zctx.k,NULL,&zctx.min);CHKERRQ(ierr); 2710298fd71SBarry Smith ierr = VecStrideMax(xin,zctx.k,NULL,&zctx.max);CHKERRQ(ierr); 27267dd0837SBarry Smith if (zctx.k < nbounds) { 273f3f0eb19SBarry Smith zctx.min = bounds[2*zctx.k]; 274f3f0eb19SBarry Smith zctx.max = bounds[2*zctx.k+1]; 27567dd0837SBarry Smith } 27647c6ae99SBarry Smith if (zctx.min == zctx.max) { 27747c6ae99SBarry Smith zctx.min -= 1.e-12; 27847c6ae99SBarry Smith zctx.max += 1.e-12; 27947c6ae99SBarry Smith } 28057622a8eSBarry Smith ierr = PetscInfo2(da,"DMDA 2d contour plot min %g max %g\n",(double)zctx.min,(double)zctx.max);CHKERRQ(ierr); 28147c6ae99SBarry Smith 282*832b7cebSLisandro Dalcin if (useports) { 283*832b7cebSLisandro Dalcin ierr = PetscDrawViewPortsSet(ports,i);CHKERRQ(ierr); 284*832b7cebSLisandro Dalcin } else { 285*832b7cebSLisandro Dalcin const char *title; 286*832b7cebSLisandro Dalcin ierr = PetscViewerDrawGetDraw(viewer,i,&draw);CHKERRQ(ierr); 287*832b7cebSLisandro Dalcin ierr = DMDAGetFieldName(da,zctx.k,&title);CHKERRQ(ierr); 288*832b7cebSLisandro Dalcin if (title) {ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr);} 289*832b7cebSLisandro Dalcin } 290*832b7cebSLisandro Dalcin 29147c6ae99SBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 29247c6ae99SBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,zctx.min,zctx.max);CHKERRQ(ierr);} 293e5ab1681SLisandro Dalcin ierr = PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr); 29447c6ae99SBarry Smith ierr = PetscDrawZoom(draw,VecView_MPI_Draw_DA2d_Zoom,&zctx);CHKERRQ(ierr); 295*832b7cebSLisandro Dalcin if (!useports) {ierr = PetscDrawSave(draw);CHKERRQ(ierr);} 29647c6ae99SBarry Smith } 297*832b7cebSLisandro Dalcin if (useports) { 298*832b7cebSLisandro Dalcin ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 299*832b7cebSLisandro Dalcin ierr = PetscDrawSave(draw);CHKERRQ(ierr); 300*832b7cebSLisandro Dalcin } 301*832b7cebSLisandro Dalcin 3026bf464f9SBarry Smith ierr = PetscDrawViewPortsDestroy(ports);CHKERRQ(ierr); 303e5ab1681SLisandro Dalcin ierr = PetscFree(displayfields);CHKERRQ(ierr); 3045edff71fSBarry Smith ierr = VecRestoreArrayRead(xcoorl,&zctx.xy);CHKERRQ(ierr); 3055edff71fSBarry Smith ierr = VecRestoreArrayRead(xlocal,&zctx.v);CHKERRQ(ierr); 30647c6ae99SBarry Smith PetscFunctionReturn(0); 30747c6ae99SBarry Smith } 30847c6ae99SBarry Smith 3090da24e51SJuha Jäykkä #if defined(PETSC_HAVE_HDF5) 3100da24e51SJuha Jäykkä #undef __FUNCT__ 3110da24e51SJuha Jäykkä #define __FUNCT__ "VecGetHDF5ChunkSize" 312c73cfb54SMatthew G. Knepley static PetscErrorCode VecGetHDF5ChunkSize(DM_DA *da, Vec xin, PetscInt dimension, PetscInt timestep, hsize_t *chunkDims) 3130da24e51SJuha Jäykkä { 314d4190030SJed Brown PetscMPIInt comm_size; 315d4190030SJed Brown PetscErrorCode ierr; 316561a82e7SJed Brown hsize_t chunk_size, target_size, dim; 317561a82e7SJed Brown hsize_t vec_size = sizeof(PetscScalar)*da->M*da->N*da->P*da->w; 318561a82e7SJed Brown hsize_t avg_local_vec_size,KiB = 1024,MiB = KiB*KiB,GiB = MiB*KiB,min_size = MiB; 319561a82e7SJed Brown hsize_t max_chunks = 64*KiB; /* HDF5 internal limitation */ 320561a82e7SJed Brown hsize_t max_chunk_size = 4*GiB; /* HDF5 internal limitation */ 32184179ffaSJed Brown PetscInt zslices=da->p, yslices=da->n, xslices=da->m; 3220da24e51SJuha Jäykkä 323cb5c4f94SJuha Jäykkä PetscFunctionBegin; 324d4190030SJed Brown ierr = MPI_Comm_size(PetscObjectComm((PetscObject)xin), &comm_size);CHKERRQ(ierr); 3250da24e51SJuha 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 */ 3260da24e51SJuha Jäykkä 327794a961bSBarry 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))))); 3287d310018SBarry 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 */ 3297d310018SBarry 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); 3300da24e51SJuha Jäykkä 331cb5c4f94SJuha Jäykkä /* 332cb5c4f94SJuha Jäykkä if size/rank > max_chunk_size, we need radical measures: even going down to 333cb5c4f94SJuha Jäykkä avg_local_vec_size is not enough, so we simply use chunk size of 4 GiB no matter 334cb5c4f94SJuha Jäykkä what, composed in the most efficient way possible. 335cb5c4f94SJuha Jäykkä N.B. this minimises the number of chunks, which may or may not be the optimal 336cb5c4f94SJuha Jäykkä solution. In a BG, for example, the optimal solution is probably to make # chunks = # 337cb5c4f94SJuha Jäykkä IO nodes involved, but this author has no access to a BG to figure out how to 338cb5c4f94SJuha Jäykkä reliably find the right number. And even then it may or may not be enough. 339cb5c4f94SJuha Jäykkä */ 3400da24e51SJuha Jäykkä if (avg_local_vec_size > max_chunk_size) { 341cb5c4f94SJuha Jäykkä /* check if we can just split local z-axis: is that enough? */ 34284179ffaSJed Brown zslices = (PetscInt)ceil(vec_size*1.0/(da->p*max_chunk_size))*zslices; 3430da24e51SJuha Jäykkä if (zslices > da->P) { 344cb5c4f94SJuha Jäykkä /* lattice is too large in xy-directions, splitting z only is not enough */ 3450da24e51SJuha Jäykkä zslices = da->P; 34684179ffaSJed Brown yslices= (PetscInt)ceil(vec_size*1.0/(zslices*da->n*max_chunk_size))*yslices; 3470da24e51SJuha Jäykkä if (yslices > da->N) { 348cb5c4f94SJuha Jäykkä /* lattice is too large in x-direction, splitting along z, y is not enough */ 3490da24e51SJuha Jäykkä yslices = da->N; 35084179ffaSJed Brown xslices= (PetscInt)ceil(vec_size*1.0/(zslices*yslices*da->m*max_chunk_size))*xslices; 3510da24e51SJuha Jäykkä } 3520da24e51SJuha Jäykkä } 3530da24e51SJuha Jäykkä dim = 0; 3540da24e51SJuha Jäykkä if (timestep >= 0) { 3550da24e51SJuha Jäykkä ++dim; 3560da24e51SJuha Jäykkä } 357cb5c4f94SJuha Jäykkä /* prefer to split z-axis, even down to planar slices */ 358c73cfb54SMatthew G. Knepley if (dimension == 3) { 359cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t) da->P/zslices; 360cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t) da->N/yslices; 361cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t) da->M/xslices; 3620da24e51SJuha Jäykkä } else { 363cb5c4f94SJuha Jäykkä /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */ 364cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t) da->N/yslices; 365cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t) da->M/xslices; 3660da24e51SJuha Jäykkä } 3670da24e51SJuha 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); 3680da24e51SJuha Jäykkä } else { 3690da24e51SJuha Jäykkä if (target_size < chunk_size) { 370cb5c4f94SJuha Jäykkä /* only change the defaults if target_size < chunk_size */ 3710da24e51SJuha Jäykkä dim = 0; 3720da24e51SJuha Jäykkä if (timestep >= 0) { 3730da24e51SJuha Jäykkä ++dim; 3740da24e51SJuha Jäykkä } 375cb5c4f94SJuha Jäykkä /* prefer to split z-axis, even down to planar slices */ 376c73cfb54SMatthew G. Knepley if (dimension == 3) { 377cb5c4f94SJuha Jäykkä /* try splitting the z-axis to core-size bits, i.e. divide chunk size by # comm_size in z-direction */ 3780da24e51SJuha Jäykkä if (target_size >= chunk_size/da->p) { 379cb5c4f94SJuha Jäykkä /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */ 3800da24e51SJuha Jäykkä chunkDims[dim] = (hsize_t) ceil(da->P*1.0/da->p); 3810da24e51SJuha Jäykkä } else { 382cb5c4f94SJuha Jäykkä /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be 383cb5c4f94SJuha Jäykkä radical and let everyone write all they've got */ 3840da24e51SJuha Jäykkä chunkDims[dim++] = (hsize_t) ceil(da->P*1.0/da->p); 3850da24e51SJuha Jäykkä chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n); 3860da24e51SJuha Jäykkä chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m); 3870da24e51SJuha Jäykkä } 3880da24e51SJuha Jäykkä } else { 389cb5c4f94SJuha Jäykkä /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */ 3900da24e51SJuha Jäykkä if (target_size >= chunk_size/da->n) { 391cb5c4f94SJuha Jäykkä /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */ 3920da24e51SJuha Jäykkä chunkDims[dim] = (hsize_t) ceil(da->N*1.0/da->n); 3930da24e51SJuha Jäykkä } else { 394cb5c4f94SJuha Jäykkä /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be 395cb5c4f94SJuha Jäykkä radical and let everyone write all they've got */ 3960da24e51SJuha Jäykkä chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n); 3970da24e51SJuha Jäykkä chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m); 3980da24e51SJuha Jäykkä } 3990da24e51SJuha Jäykkä 4000da24e51SJuha Jäykkä } 4010da24e51SJuha 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); 4020da24e51SJuha Jäykkä } else { 403cb5c4f94SJuha Jäykkä /* precomputed chunks are fine, we don't need to do anything */ 4040da24e51SJuha Jäykkä } 4050da24e51SJuha Jäykkä } 4060da24e51SJuha Jäykkä PetscFunctionReturn(0); 4070da24e51SJuha Jäykkä } 4080da24e51SJuha Jäykkä #endif 40947c6ae99SBarry Smith 41047c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 41147c6ae99SBarry Smith #undef __FUNCT__ 41247c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_HDF5_DA" 41347c6ae99SBarry Smith PetscErrorCode VecView_MPI_HDF5_DA(Vec xin,PetscViewer viewer) 41447c6ae99SBarry Smith { 4159b2a5a86SJed Brown DM dm; 4169b2a5a86SJed Brown DM_DA *da; 41747c6ae99SBarry Smith hid_t filespace; /* file dataspace identifier */ 4188e2ae6d7SMichael Kraus hid_t chunkspace; /* chunk dataset property identifier */ 41947c6ae99SBarry Smith hid_t plist_id; /* property list identifier */ 42047c6ae99SBarry Smith hid_t dset_id; /* dataset identifier */ 42147c6ae99SBarry Smith hid_t memspace; /* memory dataspace identifier */ 42247c6ae99SBarry Smith hid_t file_id; 42315214e8eSMatthew G Knepley hid_t group; 4249a0502c6SHåkon Strandenes hid_t memscalartype; /* scalar type for mem (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */ 4259a0502c6SHåkon Strandenes hid_t filescalartype; /* scalar type for file (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */ 426d9a4edebSJed Brown hsize_t dim; 4270da24e51SJuha 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 */ 428c73cfb54SMatthew G. Knepley PetscInt timestep, dimension; 4295edff71fSBarry Smith const PetscScalar *x; 4308e2ae6d7SMichael Kraus const char *vecname; 43115214e8eSMatthew G Knepley PetscErrorCode ierr; 43282ea9b62SBarry Smith PetscBool dim2; 4339a0502c6SHåkon Strandenes PetscBool spoutput; 43447c6ae99SBarry Smith 43547c6ae99SBarry Smith PetscFunctionBegin; 43615214e8eSMatthew G Knepley ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr); 4378e2ae6d7SMichael Kraus ierr = PetscViewerHDF5GetTimestep(viewer, ×tep);CHKERRQ(ierr); 43882ea9b62SBarry Smith ierr = PetscViewerHDF5GetBaseDimension2(viewer,&dim2);CHKERRQ(ierr); 4399a0502c6SHåkon Strandenes ierr = PetscViewerHDF5GetSPOutput(viewer,&spoutput);CHKERRQ(ierr); 44015214e8eSMatthew G Knepley 441c688c046SMatthew G Knepley ierr = VecGetDM(xin,&dm);CHKERRQ(ierr); 442ce94432eSBarry Smith if (!dm) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 4439b2a5a86SJed Brown da = (DM_DA*)dm->data; 444c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dimension);CHKERRQ(ierr); 44547c6ae99SBarry Smith 4468e2ae6d7SMichael Kraus /* Create the dataspace for the dataset. 4478e2ae6d7SMichael Kraus * 4488e2ae6d7SMichael Kraus * dims - holds the current dimensions of the dataset 4498e2ae6d7SMichael Kraus * 4508e2ae6d7SMichael Kraus * maxDims - holds the maximum dimensions of the dataset (unlimited 4518e2ae6d7SMichael Kraus * for the number of time steps with the current dimensions for the 4528e2ae6d7SMichael Kraus * other dimensions; so only additional time steps can be added). 4538e2ae6d7SMichael Kraus * 4548e2ae6d7SMichael Kraus * chunkDims - holds the size of a single time step (required to 4558e2ae6d7SMichael Kraus * permit extending dataset). 4568e2ae6d7SMichael Kraus */ 4578e2ae6d7SMichael Kraus dim = 0; 4588e2ae6d7SMichael Kraus if (timestep >= 0) { 4598e2ae6d7SMichael Kraus dims[dim] = timestep+1; 4608e2ae6d7SMichael Kraus maxDims[dim] = H5S_UNLIMITED; 4618e2ae6d7SMichael Kraus chunkDims[dim] = 1; 4628e2ae6d7SMichael Kraus ++dim; 4638e2ae6d7SMichael Kraus } 464c73cfb54SMatthew G. Knepley if (dimension == 3) { 465acba2ac6SBarry Smith ierr = PetscHDF5IntCast(da->P,dims+dim);CHKERRQ(ierr); 4668e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4678e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4688e2ae6d7SMichael Kraus ++dim; 4698e2ae6d7SMichael Kraus } 470c73cfb54SMatthew G. Knepley if (dimension > 1) { 471acba2ac6SBarry Smith ierr = PetscHDF5IntCast(da->N,dims+dim);CHKERRQ(ierr); 4728e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4738e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4748e2ae6d7SMichael Kraus ++dim; 4758e2ae6d7SMichael Kraus } 476acba2ac6SBarry Smith ierr = PetscHDF5IntCast(da->M,dims+dim);CHKERRQ(ierr); 4778e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4788e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4798e2ae6d7SMichael Kraus ++dim; 48082ea9b62SBarry Smith if (da->w > 1 || dim2) { 481acba2ac6SBarry Smith ierr = PetscHDF5IntCast(da->w,dims+dim);CHKERRQ(ierr); 4828e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4838e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4848e2ae6d7SMichael Kraus ++dim; 4858e2ae6d7SMichael Kraus } 48647c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 4878e2ae6d7SMichael Kraus dims[dim] = 2; 4888e2ae6d7SMichael Kraus maxDims[dim] = dims[dim]; 4898e2ae6d7SMichael Kraus chunkDims[dim] = dims[dim]; 4908e2ae6d7SMichael Kraus ++dim; 49147c6ae99SBarry Smith #endif 4920da24e51SJuha Jäykkä 493c73cfb54SMatthew G. Knepley ierr = VecGetHDF5ChunkSize(da, xin, dimension, timestep, chunkDims); CHKERRQ(ierr); 4940da24e51SJuha Jäykkä 495729ab6d9SBarry Smith PetscStackCallHDF5Return(filespace,H5Screate_simple,(dim, dims, maxDims)); 49647c6ae99SBarry Smith 49715214e8eSMatthew G Knepley #if defined(PETSC_USE_REAL_SINGLE) 4989a0502c6SHåkon Strandenes memscalartype = H5T_NATIVE_FLOAT; 4999a0502c6SHåkon Strandenes filescalartype = H5T_NATIVE_FLOAT; 50015214e8eSMatthew G Knepley #elif defined(PETSC_USE_REAL___FLOAT128) 50115214e8eSMatthew G Knepley #error "HDF5 output with 128 bit floats not supported." 50215214e8eSMatthew G Knepley #else 5039a0502c6SHåkon Strandenes memscalartype = H5T_NATIVE_DOUBLE; 5049a0502c6SHåkon Strandenes if (spoutput == PETSC_TRUE) filescalartype = H5T_NATIVE_FLOAT; 5059a0502c6SHåkon Strandenes else filescalartype = H5T_NATIVE_DOUBLE; 50615214e8eSMatthew G Knepley #endif 50715214e8eSMatthew G Knepley 50847c6ae99SBarry Smith /* Create the dataset with default properties and close filespace */ 50947c6ae99SBarry Smith ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr); 51015214e8eSMatthew G Knepley if (!H5Lexists(group, vecname, H5P_DEFAULT)) { 5118e2ae6d7SMichael Kraus /* Create chunk */ 512729ab6d9SBarry Smith PetscStackCallHDF5Return(chunkspace,H5Pcreate,(H5P_DATASET_CREATE)); 513729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_chunk,(chunkspace, dim, chunkDims)); 5148e2ae6d7SMichael Kraus 51547c6ae99SBarry Smith #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 5169a0502c6SHåkon Strandenes PetscStackCallHDF5Return(dset_id,H5Dcreate2,(group, vecname, filescalartype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT)); 51747c6ae99SBarry Smith #else 5189a0502c6SHåkon Strandenes PetscStackCallHDF5Return(dset_id,H5Dcreate,(group, vecname, filescalartype, filespace, H5P_DEFAULT)); 51947c6ae99SBarry Smith #endif 52015214e8eSMatthew G Knepley } else { 521729ab6d9SBarry Smith PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT)); 522729ab6d9SBarry Smith PetscStackCallHDF5(H5Dset_extent,(dset_id, dims)); 52315214e8eSMatthew G Knepley } 524729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(filespace)); 52547c6ae99SBarry Smith 52647c6ae99SBarry Smith /* Each process defines a dataset and writes it to the hyperslab in the file */ 5278e2ae6d7SMichael Kraus dim = 0; 5288e2ae6d7SMichael Kraus if (timestep >= 0) { 5298e2ae6d7SMichael Kraus offset[dim] = timestep; 5308e2ae6d7SMichael Kraus ++dim; 5318e2ae6d7SMichael Kraus } 532c73cfb54SMatthew G. Knepley if (dimension == 3) {ierr = PetscHDF5IntCast(da->zs,offset + dim++);CHKERRQ(ierr);} 533c73cfb54SMatthew G. Knepley if (dimension > 1) {ierr = PetscHDF5IntCast(da->ys,offset + dim++);CHKERRQ(ierr);} 534acba2ac6SBarry Smith ierr = PetscHDF5IntCast(da->xs/da->w,offset + dim++);CHKERRQ(ierr); 53582ea9b62SBarry Smith if (da->w > 1 || dim2) offset[dim++] = 0; 53647c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 5378e2ae6d7SMichael Kraus offset[dim++] = 0; 53847c6ae99SBarry Smith #endif 5398e2ae6d7SMichael Kraus dim = 0; 5408e2ae6d7SMichael Kraus if (timestep >= 0) { 5418e2ae6d7SMichael Kraus count[dim] = 1; 5428e2ae6d7SMichael Kraus ++dim; 5438e2ae6d7SMichael Kraus } 544c73cfb54SMatthew G. Knepley if (dimension == 3) {ierr = PetscHDF5IntCast(da->ze - da->zs,count + dim++);CHKERRQ(ierr);} 545c73cfb54SMatthew G. Knepley if (dimension > 1) {ierr = PetscHDF5IntCast(da->ye - da->ys,count + dim++);CHKERRQ(ierr);} 546acba2ac6SBarry Smith ierr = PetscHDF5IntCast((da->xe - da->xs)/da->w,count + dim++);CHKERRQ(ierr); 54782ea9b62SBarry Smith if (da->w > 1 || dim2) {ierr = PetscHDF5IntCast(da->w,count + dim++);CHKERRQ(ierr);} 54847c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 5498e2ae6d7SMichael Kraus count[dim++] = 2; 55047c6ae99SBarry Smith #endif 551729ab6d9SBarry Smith PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL)); 552729ab6d9SBarry Smith PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id)); 553729ab6d9SBarry Smith PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL)); 55447c6ae99SBarry Smith 55547c6ae99SBarry Smith /* Create property list for collective dataset write */ 556729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_DATASET_XFER)); 55747c6ae99SBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 558729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_dxpl_mpio,(plist_id, H5FD_MPIO_COLLECTIVE)); 55947c6ae99SBarry Smith #endif 56047c6ae99SBarry Smith /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 56147c6ae99SBarry Smith 5625edff71fSBarry Smith ierr = VecGetArrayRead(xin, &x);CHKERRQ(ierr); 5639a0502c6SHåkon Strandenes PetscStackCallHDF5(H5Dwrite,(dset_id, memscalartype, memspace, filespace, plist_id, x)); 564729ab6d9SBarry Smith PetscStackCallHDF5(H5Fflush,(file_id, H5F_SCOPE_GLOBAL)); 5655edff71fSBarry Smith ierr = VecRestoreArrayRead(xin, &x);CHKERRQ(ierr); 56647c6ae99SBarry Smith 56747c6ae99SBarry Smith /* Close/release resources */ 56815214e8eSMatthew G Knepley if (group != file_id) { 569729ab6d9SBarry Smith PetscStackCallHDF5(H5Gclose,(group)); 57015214e8eSMatthew G Knepley } 571729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 572729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(filespace)); 573729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(memspace)); 574729ab6d9SBarry Smith PetscStackCallHDF5(H5Dclose,(dset_id)); 57547c6ae99SBarry Smith ierr = PetscInfo1(xin,"Wrote Vec object with name %s\n",vecname);CHKERRQ(ierr); 57647c6ae99SBarry Smith PetscFunctionReturn(0); 57747c6ae99SBarry Smith } 57847c6ae99SBarry Smith #endif 57947c6ae99SBarry Smith 58009573ac7SBarry Smith extern PetscErrorCode VecView_MPI_Draw_DA1d(Vec,PetscViewer); 58147c6ae99SBarry Smith 58247c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 58347c6ae99SBarry Smith #undef __FUNCT__ 584aa219208SBarry Smith #define __FUNCT__ "DMDAArrayMPIIO" 585aa219208SBarry Smith static PetscErrorCode DMDAArrayMPIIO(DM da,PetscViewer viewer,Vec xin,PetscBool write) 58647c6ae99SBarry Smith { 58747c6ae99SBarry Smith PetscErrorCode ierr; 58847c6ae99SBarry Smith MPI_File mfdes; 58947c6ae99SBarry Smith PetscMPIInt gsizes[4],lsizes[4],lstarts[4],asiz,dof; 59047c6ae99SBarry Smith MPI_Datatype view; 59147c6ae99SBarry Smith const PetscScalar *array; 59247c6ae99SBarry Smith MPI_Offset off; 59347c6ae99SBarry Smith MPI_Aint ub,ul; 59447c6ae99SBarry Smith PetscInt type,rows,vecrows,tr[2]; 59547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 5960c169764Sdmay PetscBool skipheader; 59747c6ae99SBarry Smith 59847c6ae99SBarry Smith PetscFunctionBegin; 59947c6ae99SBarry Smith ierr = VecGetSize(xin,&vecrows);CHKERRQ(ierr); 6000c169764Sdmay ierr = PetscViewerBinaryGetSkipHeader(viewer,&skipheader);CHKERRQ(ierr); 60147c6ae99SBarry Smith if (!write) { 60247c6ae99SBarry Smith /* Read vector header. */ 6030c169764Sdmay if (!skipheader) { 604060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,tr,2,NULL,PETSC_INT);CHKERRQ(ierr); 60547c6ae99SBarry Smith type = tr[0]; 60647c6ae99SBarry Smith rows = tr[1]; 607ce94432eSBarry Smith if (type != VEC_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Not vector next in file"); 608ce94432eSBarry Smith if (rows != vecrows) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Vector in file not same size as DMDA vector"); 6090c169764Sdmay } 61047c6ae99SBarry Smith } else { 61147c6ae99SBarry Smith tr[0] = VEC_FILE_CLASSID; 61247c6ae99SBarry Smith tr[1] = vecrows; 6130c169764Sdmay if (!skipheader) { 61447c6ae99SBarry Smith ierr = PetscViewerBinaryWrite(viewer,tr,2,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 61547c6ae99SBarry Smith } 6160c169764Sdmay } 61747c6ae99SBarry Smith 6184dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->w,&dof);CHKERRQ(ierr); 6194dc2109aSBarry Smith gsizes[0] = dof; 6204dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->M,gsizes+1);CHKERRQ(ierr); 6214dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->N,gsizes+2);CHKERRQ(ierr); 622334634e2SJed Brown ierr = PetscMPIIntCast(dd->P,gsizes+3);CHKERRQ(ierr); 6234dc2109aSBarry Smith lsizes[0] = dof; 6244dc2109aSBarry Smith ierr = PetscMPIIntCast((dd->xe-dd->xs)/dof,lsizes+1);CHKERRQ(ierr); 6254dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->ye-dd->ys,lsizes+2);CHKERRQ(ierr); 6264dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->ze-dd->zs,lsizes+3);CHKERRQ(ierr); 6274dc2109aSBarry Smith lstarts[0] = 0; 6284dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->xs/dof,lstarts+1);CHKERRQ(ierr); 6294dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->ys,lstarts+2);CHKERRQ(ierr); 6304dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->zs,lstarts+3);CHKERRQ(ierr); 631c73cfb54SMatthew G. Knepley ierr = MPI_Type_create_subarray(da->dim+1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view);CHKERRQ(ierr); 63247c6ae99SBarry Smith ierr = MPI_Type_commit(&view);CHKERRQ(ierr); 63347c6ae99SBarry Smith 63447c6ae99SBarry Smith ierr = PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);CHKERRQ(ierr); 63547c6ae99SBarry Smith ierr = PetscViewerBinaryGetMPIIOOffset(viewer,&off);CHKERRQ(ierr); 63647c6ae99SBarry Smith ierr = MPI_File_set_view(mfdes,off,MPIU_SCALAR,view,(char*)"native",MPI_INFO_NULL);CHKERRQ(ierr); 63747c6ae99SBarry Smith ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr); 63847c6ae99SBarry Smith asiz = lsizes[1]*(lsizes[2] > 0 ? lsizes[2] : 1)*(lsizes[3] > 0 ? lsizes[3] : 1)*dof; 63947c6ae99SBarry Smith if (write) { 64047c6ae99SBarry Smith ierr = MPIU_File_write_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr); 64147c6ae99SBarry Smith } else { 64247c6ae99SBarry Smith ierr = MPIU_File_read_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr); 64347c6ae99SBarry Smith } 64447c6ae99SBarry Smith ierr = MPI_Type_get_extent(view,&ul,&ub);CHKERRQ(ierr); 64547c6ae99SBarry Smith ierr = PetscViewerBinaryAddMPIIOOffset(viewer,ub);CHKERRQ(ierr); 64647c6ae99SBarry Smith ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr); 64747c6ae99SBarry Smith ierr = MPI_Type_free(&view);CHKERRQ(ierr); 64847c6ae99SBarry Smith PetscFunctionReturn(0); 64947c6ae99SBarry Smith } 65047c6ae99SBarry Smith #endif 65147c6ae99SBarry Smith 65247c6ae99SBarry Smith #undef __FUNCT__ 65347c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_DA" 6547087cfbeSBarry Smith PetscErrorCode VecView_MPI_DA(Vec xin,PetscViewer viewer) 65547c6ae99SBarry Smith { 6569a42bb27SBarry Smith DM da; 65747c6ae99SBarry Smith PetscErrorCode ierr; 65847c6ae99SBarry Smith PetscInt dim; 65947c6ae99SBarry Smith Vec natural; 6604061b8bfSJed Brown PetscBool isdraw,isvtk; 66147c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 66247c6ae99SBarry Smith PetscBool ishdf5; 66347c6ae99SBarry Smith #endif 6643f3fd955SJed Brown const char *prefix,*name; 665a261c58fSBarry Smith PetscViewerFormat format; 66647c6ae99SBarry Smith 66747c6ae99SBarry Smith PetscFunctionBegin; 668c688c046SMatthew G Knepley ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 669ce94432eSBarry Smith if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 670251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 671251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr); 67247c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 673251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr); 67447c6ae99SBarry Smith #endif 67547c6ae99SBarry Smith if (isdraw) { 6761321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 67747c6ae99SBarry Smith if (dim == 1) { 67847c6ae99SBarry Smith ierr = VecView_MPI_Draw_DA1d(xin,viewer);CHKERRQ(ierr); 67947c6ae99SBarry Smith } else if (dim == 2) { 68047c6ae99SBarry Smith ierr = VecView_MPI_Draw_DA2d(xin,viewer);CHKERRQ(ierr); 681ce94432eSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot graphically view vector associated with this dimensional DMDA %D",dim); 6824061b8bfSJed Brown } else if (isvtk) { /* Duplicate the Vec and hold a reference to the DM */ 6834061b8bfSJed Brown Vec Y; 6844061b8bfSJed Brown ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr); 6854061b8bfSJed Brown ierr = VecDuplicate(xin,&Y);CHKERRQ(ierr); 686b51b94faSJed Brown if (((PetscObject)xin)->name) { 687b51b94faSJed Brown /* If xin was named, copy the name over to Y. The duplicate names are safe because nobody else will ever see Y. */ 688b51b94faSJed Brown ierr = PetscObjectSetName((PetscObject)Y,((PetscObject)xin)->name);CHKERRQ(ierr); 689b51b94faSJed Brown } 6904061b8bfSJed Brown ierr = VecCopy(xin,Y);CHKERRQ(ierr); 69162b69a3fSMatthew G Knepley ierr = PetscViewerVTKAddField(viewer,(PetscObject)da,DMDAVTKWriteAll,PETSC_VTK_POINT_FIELD,(PetscObject)Y);CHKERRQ(ierr); 69247c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 69347c6ae99SBarry Smith } else if (ishdf5) { 69447c6ae99SBarry Smith ierr = VecView_MPI_HDF5_DA(xin,viewer);CHKERRQ(ierr); 69547c6ae99SBarry Smith #endif 69647c6ae99SBarry Smith } else { 69747c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 69847c6ae99SBarry Smith PetscBool isbinary,isMPIIO; 69947c6ae99SBarry Smith 700251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 70147c6ae99SBarry Smith if (isbinary) { 702bc196f7cSDave May ierr = PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO);CHKERRQ(ierr); 70347c6ae99SBarry Smith if (isMPIIO) { 704aa219208SBarry Smith ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_TRUE);CHKERRQ(ierr); 70547c6ae99SBarry Smith PetscFunctionReturn(0); 70647c6ae99SBarry Smith } 70747c6ae99SBarry Smith } 70847c6ae99SBarry Smith #endif 70947c6ae99SBarry Smith 71047c6ae99SBarry Smith /* call viewer on natural ordering */ 71147c6ae99SBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr); 712aa219208SBarry Smith ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr); 71347c6ae99SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr); 714aa219208SBarry Smith ierr = DMDAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr); 715aa219208SBarry Smith ierr = DMDAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr); 7163f3fd955SJed Brown ierr = PetscObjectGetName((PetscObject)xin,&name);CHKERRQ(ierr); 7173f3fd955SJed Brown ierr = PetscObjectSetName((PetscObject)natural,name);CHKERRQ(ierr); 718a261c58fSBarry Smith 719a261c58fSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 720a261c58fSBarry Smith if (format == PETSC_VIEWER_BINARY_MATLAB) { 721a261c58fSBarry Smith /* temporarily remove viewer format so it won't trigger in the VecView() */ 7226a9046bcSBarry Smith ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr); 723a261c58fSBarry Smith } 724a261c58fSBarry Smith 72547c6ae99SBarry Smith ierr = VecView(natural,viewer);CHKERRQ(ierr); 726a261c58fSBarry Smith 727a261c58fSBarry Smith if (format == PETSC_VIEWER_BINARY_MATLAB) { 728a261c58fSBarry Smith MPI_Comm comm; 729a261c58fSBarry Smith FILE *info; 730a261c58fSBarry Smith const char *fieldname; 731da88d4d4SJed Brown char fieldbuf[256]; 732a261c58fSBarry Smith PetscInt dim,ni,nj,nk,pi,pj,pk,dof,n; 733a261c58fSBarry Smith 734a261c58fSBarry Smith /* set the viewer format back into the viewer */ 7356a9046bcSBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 736a261c58fSBarry Smith ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); 737a261c58fSBarry Smith ierr = PetscViewerBinaryGetInfoPointer(viewer,&info);CHKERRQ(ierr); 738a261c58fSBarry Smith ierr = DMDAGetInfo(da,&dim,&ni,&nj,&nk,&pi,&pj,&pk,&dof,0,0,0,0,0);CHKERRQ(ierr); 739da88d4d4SJed Brown ierr = PetscFPrintf(comm,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n");CHKERRQ(ierr); 740da88d4d4SJed Brown ierr = PetscFPrintf(comm,info,"#$$ tmp = PetscBinaryRead(fd); \n");CHKERRQ(ierr); 741da88d4d4SJed Brown if (dim == 1) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d);\n",dof,ni);CHKERRQ(ierr); } 742da88d4d4SJed Brown if (dim == 2) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d);\n",dof,ni,nj);CHKERRQ(ierr); } 743da88d4d4SJed Brown if (dim == 3) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d,%d);\n",dof,ni,nj,nk);CHKERRQ(ierr); } 744a261c58fSBarry Smith 745a261c58fSBarry Smith for (n=0; n<dof; n++) { 746a261c58fSBarry Smith ierr = DMDAGetFieldName(da,n,&fieldname);CHKERRQ(ierr); 747da88d4d4SJed Brown if (!fieldname || !fieldname[0]) { 748da88d4d4SJed Brown ierr = PetscSNPrintf(fieldbuf,sizeof fieldbuf,"field%D",n);CHKERRQ(ierr); 749da88d4d4SJed Brown fieldname = fieldbuf; 750a261c58fSBarry Smith } 751da88d4d4SJed Brown if (dim == 1) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:))';\n",name,fieldname,n+1);CHKERRQ(ierr); } 752da88d4d4SJed Brown if (dim == 2) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:,:))';\n",name,fieldname,n+1);CHKERRQ(ierr); } 753da88d4d4SJed 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);} 754a261c58fSBarry Smith } 755da88d4d4SJed Brown ierr = PetscFPrintf(comm,info,"#$$ clear tmp; \n");CHKERRQ(ierr); 756da88d4d4SJed Brown ierr = PetscFPrintf(comm,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n");CHKERRQ(ierr); 757a261c58fSBarry Smith } 758a261c58fSBarry Smith 759fcfd50ebSBarry Smith ierr = VecDestroy(&natural);CHKERRQ(ierr); 76047c6ae99SBarry Smith } 76147c6ae99SBarry Smith PetscFunctionReturn(0); 76247c6ae99SBarry Smith } 76347c6ae99SBarry Smith 76447c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 76547c6ae99SBarry Smith #undef __FUNCT__ 76647c6ae99SBarry Smith #define __FUNCT__ "VecLoad_HDF5_DA" 76747c6ae99SBarry Smith PetscErrorCode VecLoad_HDF5_DA(Vec xin, PetscViewer viewer) 76847c6ae99SBarry Smith { 7699a42bb27SBarry Smith DM da; 77047c6ae99SBarry Smith PetscErrorCode ierr; 771ec7ae49cSHåkon Strandenes hsize_t dim,rdim; 772ec7ae49cSHåkon Strandenes hsize_t dims[6]={0},count[6]={0},offset[6]={0}; 773ec7ae49cSHåkon Strandenes PetscInt dimension,timestep,dofInd; 77447c6ae99SBarry Smith PetscScalar *x; 77547c6ae99SBarry Smith const char *vecname; 77647c6ae99SBarry Smith hid_t filespace; /* file dataspace identifier */ 77747c6ae99SBarry Smith hid_t plist_id; /* property list identifier */ 77847c6ae99SBarry Smith hid_t dset_id; /* dataset identifier */ 77947c6ae99SBarry Smith hid_t memspace; /* memory dataspace identifier */ 780bfd264e7SBarry Smith hid_t file_id,group; 7817bcbaff4SBarry Smith hid_t scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */ 7829c7c4993SBarry Smith DM_DA *dd; 783ec7ae49cSHåkon Strandenes PetscBool dim2 = PETSC_FALSE; 78447c6ae99SBarry Smith 78547c6ae99SBarry Smith PetscFunctionBegin; 7867bcbaff4SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 7877bcbaff4SBarry Smith scalartype = H5T_NATIVE_FLOAT; 7887bcbaff4SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128) 7897bcbaff4SBarry Smith #error "HDF5 output with 128 bit floats not supported." 7907bcbaff4SBarry Smith #else 7917bcbaff4SBarry Smith scalartype = H5T_NATIVE_DOUBLE; 7927bcbaff4SBarry Smith #endif 7937bcbaff4SBarry Smith 794bfd264e7SBarry Smith ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr); 795ec7ae49cSHåkon Strandenes ierr = PetscViewerHDF5GetTimestep(viewer, ×tep);CHKERRQ(ierr); 796ec7ae49cSHåkon Strandenes ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr); 797c688c046SMatthew G Knepley ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 7989c7c4993SBarry Smith dd = (DM_DA*)da->data; 799c73cfb54SMatthew G. Knepley ierr = DMGetDimension(da, &dimension);CHKERRQ(ierr); 80047c6ae99SBarry Smith 801ec7ae49cSHåkon Strandenes /* Open dataset */ 80247c6ae99SBarry Smith #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 803729ab6d9SBarry Smith PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT)); 80447c6ae99SBarry Smith #else 805729ab6d9SBarry Smith PetscStackCallHDF5Return(dset_id,H5Dopen,(group, vecname)); 80647c6ae99SBarry Smith #endif 80747c6ae99SBarry Smith 808ec7ae49cSHåkon Strandenes /* Retrieve the dataspace for the dataset */ 809ec7ae49cSHåkon Strandenes PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id)); 810ec7ae49cSHåkon Strandenes PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(filespace, dims, NULL)); 811ec7ae49cSHåkon Strandenes 812ec7ae49cSHåkon Strandenes /* Expected dimension for holding the dof's */ 81347c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 814ec7ae49cSHåkon Strandenes dofInd = rdim-2; 815ec7ae49cSHåkon Strandenes #else 816ec7ae49cSHåkon Strandenes dofInd = rdim-1; 81747c6ae99SBarry Smith #endif 818ec7ae49cSHåkon Strandenes 819ec7ae49cSHåkon Strandenes /* The expected number of dimensions, assuming basedimension2 = false */ 820ec7ae49cSHåkon Strandenes dim = dimension; 821ec7ae49cSHåkon Strandenes if (dd->w > 1) ++dim; 822ec7ae49cSHåkon Strandenes if (timestep >= 0) ++dim; 82347c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 824ec7ae49cSHåkon Strandenes ++dim; 82547c6ae99SBarry Smith #endif 826ec7ae49cSHåkon Strandenes 827ec7ae49cSHåkon Strandenes /* In this case the input dataset have one extra, unexpected dimension. */ 828ec7ae49cSHåkon Strandenes if (rdim == dim+1) { 829ec7ae49cSHåkon Strandenes /* In this case the block size unity */ 830ec7ae49cSHåkon Strandenes if (dd->w == 1 && dims[dofInd] == 1) dim2 = PETSC_TRUE; 831ec7ae49cSHåkon Strandenes 832ec7ae49cSHåkon Strandenes /* Special error message for the case where dof does not match the input file */ 833ec7ae49cSHå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); 834ec7ae49cSHåkon Strandenes 835ec7ae49cSHåkon Strandenes /* Other cases where rdim != dim cannot be handled currently */ 8366c4ed002SBarry 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); 837ec7ae49cSHåkon Strandenes 838ec7ae49cSHåkon Strandenes /* Set up the hyperslab size */ 839ec7ae49cSHåkon Strandenes dim = 0; 840ec7ae49cSHåkon Strandenes if (timestep >= 0) { 841ec7ae49cSHåkon Strandenes offset[dim] = timestep; 842ec7ae49cSHåkon Strandenes count[dim] = 1; 843ec7ae49cSHåkon Strandenes ++dim; 844ec7ae49cSHåkon Strandenes } 845ec7ae49cSHåkon Strandenes if (dimension == 3) { 846ec7ae49cSHåkon Strandenes ierr = PetscHDF5IntCast(dd->zs,offset + dim);CHKERRQ(ierr); 847ec7ae49cSHåkon Strandenes ierr = PetscHDF5IntCast(dd->ze - dd->zs,count + dim);CHKERRQ(ierr); 848ec7ae49cSHåkon Strandenes ++dim; 849ec7ae49cSHåkon Strandenes } 850ec7ae49cSHåkon Strandenes if (dimension > 1) { 851ec7ae49cSHåkon Strandenes ierr = PetscHDF5IntCast(dd->ys,offset + dim);CHKERRQ(ierr); 852ec7ae49cSHåkon Strandenes ierr = PetscHDF5IntCast(dd->ye - dd->ys,count + dim);CHKERRQ(ierr); 853ec7ae49cSHåkon Strandenes ++dim; 854ec7ae49cSHåkon Strandenes } 855ec7ae49cSHåkon Strandenes ierr = PetscHDF5IntCast(dd->xs/dd->w,offset + dim);CHKERRQ(ierr); 856ec7ae49cSHåkon Strandenes ierr = PetscHDF5IntCast((dd->xe - dd->xs)/dd->w,count + dim);CHKERRQ(ierr); 857ec7ae49cSHåkon Strandenes ++dim; 858ec7ae49cSHåkon Strandenes if (dd->w > 1 || dim2) { 859ec7ae49cSHåkon Strandenes offset[dim] = 0; 860ec7ae49cSHåkon Strandenes ierr = PetscHDF5IntCast(dd->w,count + dim);CHKERRQ(ierr); 861ec7ae49cSHåkon Strandenes ++dim; 862ec7ae49cSHåkon Strandenes } 863ec7ae49cSHåkon Strandenes #if defined(PETSC_USE_COMPLEX) 864ec7ae49cSHåkon Strandenes offset[dim] = 0; 865ec7ae49cSHåkon Strandenes count[dim] = 2; 866ec7ae49cSHåkon Strandenes ++dim; 867ec7ae49cSHåkon Strandenes #endif 868ec7ae49cSHåkon Strandenes 869ec7ae49cSHåkon Strandenes /* Create the memory and filespace */ 870729ab6d9SBarry Smith PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL)); 871729ab6d9SBarry Smith PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL)); 87247c6ae99SBarry Smith 87347c6ae99SBarry Smith /* Create property list for collective dataset write */ 874729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_DATASET_XFER)); 87547c6ae99SBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 876729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_dxpl_mpio,(plist_id, H5FD_MPIO_COLLECTIVE)); 87747c6ae99SBarry Smith #endif 878ec7ae49cSHåkon Strandenes /* To read dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 87947c6ae99SBarry Smith 88047c6ae99SBarry Smith ierr = VecGetArray(xin, &x);CHKERRQ(ierr); 881a1dbdba5SBarry Smith PetscStackCallHDF5(H5Dread,(dset_id, scalartype, memspace, filespace, plist_id, x)); 88247c6ae99SBarry Smith ierr = VecRestoreArray(xin, &x);CHKERRQ(ierr); 88347c6ae99SBarry Smith 88447c6ae99SBarry Smith /* Close/release resources */ 885bfd264e7SBarry Smith if (group != file_id) { 886729ab6d9SBarry Smith PetscStackCallHDF5(H5Gclose,(group)); 887bfd264e7SBarry Smith } 888729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 889729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(filespace)); 890729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(memspace)); 891729ab6d9SBarry Smith PetscStackCallHDF5(H5Dclose,(dset_id)); 89247c6ae99SBarry Smith PetscFunctionReturn(0); 89347c6ae99SBarry Smith } 89447c6ae99SBarry Smith #endif 89547c6ae99SBarry Smith 89647c6ae99SBarry Smith #undef __FUNCT__ 89747c6ae99SBarry Smith #define __FUNCT__ "VecLoad_Binary_DA" 89847c6ae99SBarry Smith PetscErrorCode VecLoad_Binary_DA(Vec xin, PetscViewer viewer) 89947c6ae99SBarry Smith { 9009a42bb27SBarry Smith DM da; 90147c6ae99SBarry Smith PetscErrorCode ierr; 90247c6ae99SBarry Smith Vec natural; 90347c6ae99SBarry Smith const char *prefix; 90447c6ae99SBarry Smith PetscInt bs; 90547c6ae99SBarry Smith PetscBool flag; 90647c6ae99SBarry Smith DM_DA *dd; 90747c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 90847c6ae99SBarry Smith PetscBool isMPIIO; 90947c6ae99SBarry Smith #endif 91047c6ae99SBarry Smith 91147c6ae99SBarry Smith PetscFunctionBegin; 912c688c046SMatthew G Knepley ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 91347c6ae99SBarry Smith dd = (DM_DA*)da->data; 91447c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 915bc196f7cSDave May ierr = PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO);CHKERRQ(ierr); 91647c6ae99SBarry Smith if (isMPIIO) { 917aa219208SBarry Smith ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_FALSE);CHKERRQ(ierr); 91847c6ae99SBarry Smith PetscFunctionReturn(0); 91947c6ae99SBarry Smith } 92047c6ae99SBarry Smith #endif 92147c6ae99SBarry Smith 92247c6ae99SBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr); 923aa219208SBarry Smith ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr); 92447c6ae99SBarry Smith ierr = PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name);CHKERRQ(ierr); 92547c6ae99SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr); 9268aea9937SBarry Smith ierr = VecLoad(natural,viewer);CHKERRQ(ierr); 927aa219208SBarry Smith ierr = DMDANaturalToGlobalBegin(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr); 928aa219208SBarry Smith ierr = DMDANaturalToGlobalEnd(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr); 929fcfd50ebSBarry Smith ierr = VecDestroy(&natural);CHKERRQ(ierr); 930aa219208SBarry Smith ierr = PetscInfo(xin,"Loading vector from natural ordering into DMDA\n");CHKERRQ(ierr); 931c5929fdfSBarry Smith ierr = PetscOptionsGetInt(NULL,((PetscObject)xin)->prefix,"-vecload_block_size",&bs,&flag);CHKERRQ(ierr); 93247c6ae99SBarry Smith if (flag && bs != dd->w) { 933aa219208SBarry Smith ierr = PetscInfo2(xin,"Block size in file %D not equal to DMDA's dof %D\n",bs,dd->w);CHKERRQ(ierr); 93447c6ae99SBarry Smith } 93547c6ae99SBarry Smith PetscFunctionReturn(0); 93647c6ae99SBarry Smith } 93747c6ae99SBarry Smith 93847c6ae99SBarry Smith #undef __FUNCT__ 93947c6ae99SBarry Smith #define __FUNCT__ "VecLoad_Default_DA" 9407087cfbeSBarry Smith PetscErrorCode VecLoad_Default_DA(Vec xin, PetscViewer viewer) 94147c6ae99SBarry Smith { 94247c6ae99SBarry Smith PetscErrorCode ierr; 9439a42bb27SBarry Smith DM da; 94447c6ae99SBarry Smith PetscBool isbinary; 94547c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 94647c6ae99SBarry Smith PetscBool ishdf5; 94747c6ae99SBarry Smith #endif 94847c6ae99SBarry Smith 94947c6ae99SBarry Smith PetscFunctionBegin; 950c688c046SMatthew G Knepley ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 951ce94432eSBarry Smith if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 95247c6ae99SBarry Smith 95347c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 954251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr); 95547c6ae99SBarry Smith #endif 956251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 95747c6ae99SBarry Smith 95847c6ae99SBarry Smith if (isbinary) { 95947c6ae99SBarry Smith ierr = VecLoad_Binary_DA(xin,viewer);CHKERRQ(ierr); 96047c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 96147c6ae99SBarry Smith } else if (ishdf5) { 96247c6ae99SBarry Smith ierr = VecLoad_HDF5_DA(xin,viewer);CHKERRQ(ierr); 96347c6ae99SBarry Smith #endif 964d34fcf5fSBarry Smith } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for vector loading", ((PetscObject)viewer)->type_name); 96547c6ae99SBarry Smith PetscFunctionReturn(0); 96647c6ae99SBarry Smith } 967