147c6ae99SBarry Smith 247c6ae99SBarry Smith /* 3aa219208SBarry Smith Plots vectors obtained with DMDACreate2d() 447c6ae99SBarry Smith */ 547c6ae99SBarry Smith 6af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 78135c375SStefano Zampini #include <petsc/private/glvisvecimpl.h> 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; 197fe7c8feSLisandro Dalcin PetscBool showaxis,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 PetscErrorCode VecView_MPI_Draw_DA2d_Zoom(PetscDraw draw,void *ctx) 2947c6ae99SBarry Smith { 3047c6ae99SBarry Smith ZoomCtx *zctx = (ZoomCtx*)ctx; 3147c6ae99SBarry Smith PetscErrorCode ierr; 32e5ab1681SLisandro Dalcin PetscInt m,n,i,j,k,dof,id,c1,c2,c3,c4; 33e5ab1681SLisandro Dalcin PetscReal min,max,x1,x2,x3,x4,y_1,y2,y3,y4; 34e5ab1681SLisandro Dalcin const PetscScalar *xy,*v; 3547c6ae99SBarry Smith 3647c6ae99SBarry Smith PetscFunctionBegin; 3747c6ae99SBarry Smith m = zctx->m; 3847c6ae99SBarry Smith n = zctx->n; 39e5ab1681SLisandro Dalcin dof = zctx->dof; 4047c6ae99SBarry Smith k = zctx->k; 4147c6ae99SBarry Smith xy = zctx->xy; 42e5ab1681SLisandro Dalcin v = zctx->v; 4347c6ae99SBarry Smith min = zctx->min; 44f3f0eb19SBarry Smith max = zctx->max; 4547c6ae99SBarry Smith 4647c6ae99SBarry Smith /* PetscDraw the contour plot patch */ 47e5ab1681SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 4847c6ae99SBarry Smith for (j=0; j<n-1; j++) { 4947c6ae99SBarry Smith for (i=0; i<m-1; i++) { 500e4fe250SBarry Smith id = i+j*m; 510e4fe250SBarry Smith x1 = PetscRealPart(xy[2*id]); 520e4fe250SBarry Smith y_1 = PetscRealPart(xy[2*id+1]); 53e5ab1681SLisandro Dalcin c1 = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max); 540e4fe250SBarry Smith 550e4fe250SBarry Smith id = i+j*m+1; 560e4fe250SBarry Smith x2 = PetscRealPart(xy[2*id]); 57a39a4c7dSLisandro Dalcin y2 = PetscRealPart(xy[2*id+1]); 58e5ab1681SLisandro Dalcin c2 = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max); 590e4fe250SBarry Smith 600e4fe250SBarry Smith id = i+j*m+1+m; 61a39a4c7dSLisandro Dalcin x3 = PetscRealPart(xy[2*id]); 620e4fe250SBarry Smith y3 = PetscRealPart(xy[2*id+1]); 63e5ab1681SLisandro Dalcin c3 = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max); 640e4fe250SBarry Smith 650e4fe250SBarry Smith id = i+j*m+m; 66a39a4c7dSLisandro Dalcin x4 = PetscRealPart(xy[2*id]); 67a39a4c7dSLisandro Dalcin y4 = PetscRealPart(xy[2*id+1]); 68e5ab1681SLisandro Dalcin c4 = PetscDrawRealToColor(PetscRealPart(v[k+dof*id]),min,max); 69f3f0eb19SBarry Smith 7047c6ae99SBarry Smith ierr = PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);CHKERRQ(ierr); 7147c6ae99SBarry Smith ierr = PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4);CHKERRQ(ierr); 7247c6ae99SBarry Smith if (zctx->showgrid) { 7347c6ae99SBarry Smith ierr = PetscDrawLine(draw,x1,y_1,x2,y2,PETSC_DRAW_BLACK);CHKERRQ(ierr); 7447c6ae99SBarry Smith ierr = PetscDrawLine(draw,x2,y2,x3,y3,PETSC_DRAW_BLACK);CHKERRQ(ierr); 7547c6ae99SBarry Smith ierr = PetscDrawLine(draw,x3,y3,x4,y4,PETSC_DRAW_BLACK);CHKERRQ(ierr); 7647c6ae99SBarry Smith ierr = PetscDrawLine(draw,x4,y4,x1,y_1,PETSC_DRAW_BLACK);CHKERRQ(ierr); 7747c6ae99SBarry Smith } 7847c6ae99SBarry Smith } 7947c6ae99SBarry Smith } 807fe7c8feSLisandro Dalcin if (zctx->showaxis && !zctx->rank) { 81e5ab1681SLisandro Dalcin if (zctx->name0 || zctx->name1) { 82109c9344SBarry Smith PetscReal xl,yl,xr,yr,x,y; 83109c9344SBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 84e5ab1681SLisandro Dalcin x = xl + .30*(xr - xl); 85109c9344SBarry Smith xl = xl + .01*(xr - xl); 86e5ab1681SLisandro Dalcin y = yr - .30*(yr - yl); 87109c9344SBarry Smith yl = yl + .01*(yr - yl); 88e5ab1681SLisandro Dalcin if (zctx->name0) {ierr = PetscDrawString(draw,x,yl,PETSC_DRAW_BLACK,zctx->name0);CHKERRQ(ierr);} 89e5ab1681SLisandro Dalcin if (zctx->name1) {ierr = PetscDrawStringVertical(draw,xl,y,PETSC_DRAW_BLACK,zctx->name1);CHKERRQ(ierr);} 90109c9344SBarry Smith } 910e4fe250SBarry Smith /* 920e4fe250SBarry Smith Ideally we would use the PetscDrawAxis object to manage displaying the coordinate limits 930e4fe250SBarry Smith but that may require some refactoring. 940e4fe250SBarry Smith */ 95e5ab1681SLisandro Dalcin { 967fe7c8feSLisandro Dalcin double xmin = (double)zctx->xmin, ymin = (double)zctx->ymin; 977fe7c8feSLisandro Dalcin double xmax = (double)zctx->xmax, ymax = (double)zctx->ymax; 98e5ab1681SLisandro Dalcin char value[16]; size_t len; PetscReal w; 997fe7c8feSLisandro Dalcin ierr = PetscSNPrintf(value,16,"%0.2e",xmin);CHKERRQ(ierr); 100e5ab1681SLisandro Dalcin ierr = PetscDrawString(draw,xmin,ymin - .05*(ymax - ymin),PETSC_DRAW_BLACK,value);CHKERRQ(ierr); 1017fe7c8feSLisandro Dalcin ierr = PetscSNPrintf(value,16,"%0.2e",xmax);CHKERRQ(ierr); 1020e4fe250SBarry Smith ierr = PetscStrlen(value,&len);CHKERRQ(ierr); 1030298fd71SBarry Smith ierr = PetscDrawStringGetSize(draw,&w,NULL);CHKERRQ(ierr); 104e5ab1681SLisandro Dalcin ierr = PetscDrawString(draw,xmax - len*w,ymin - .05*(ymax - ymin),PETSC_DRAW_BLACK,value);CHKERRQ(ierr); 1057fe7c8feSLisandro Dalcin ierr = PetscSNPrintf(value,16,"%0.2e",ymin);CHKERRQ(ierr); 106e5ab1681SLisandro Dalcin ierr = PetscDrawString(draw,xmin - .05*(xmax - xmin),ymin,PETSC_DRAW_BLACK,value);CHKERRQ(ierr); 1077fe7c8feSLisandro Dalcin ierr = PetscSNPrintf(value,16,"%0.2e",ymax);CHKERRQ(ierr); 108e5ab1681SLisandro Dalcin ierr = PetscDrawString(draw,xmin - .05*(xmax - xmin),ymax,PETSC_DRAW_BLACK,value);CHKERRQ(ierr); 109e5ab1681SLisandro Dalcin } 110e5ab1681SLisandro Dalcin } 111e5ab1681SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 11247c6ae99SBarry Smith PetscFunctionReturn(0); 11347c6ae99SBarry Smith } 11447c6ae99SBarry Smith 11547c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA2d(Vec xin,PetscViewer viewer) 11647c6ae99SBarry Smith { 1179a42bb27SBarry Smith DM da,dac,dag; 11847c6ae99SBarry Smith PetscErrorCode ierr; 119a74ba6f7SBarry Smith PetscInt N,s,M,w,ncoors = 4; 12047c6ae99SBarry Smith const PetscInt *lx,*ly; 121e5ab1681SLisandro Dalcin PetscReal coors[4]; 12247c6ae99SBarry Smith PetscDraw draw,popup; 12347c6ae99SBarry Smith PetscBool isnull,useports = PETSC_FALSE; 12447c6ae99SBarry Smith MPI_Comm comm; 12547c6ae99SBarry Smith Vec xlocal,xcoor,xcoorl; 126bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by; 127aa219208SBarry Smith DMDAStencilType st; 12847c6ae99SBarry Smith ZoomCtx zctx; 1290298fd71SBarry Smith PetscDrawViewPorts *ports = NULL; 13047c6ae99SBarry Smith PetscViewerFormat format; 13120d0051dSBarry Smith PetscInt *displayfields; 13267dd0837SBarry Smith PetscInt ndisplayfields,i,nbounds; 13367dd0837SBarry Smith const PetscReal *bounds; 13447c6ae99SBarry Smith 13547c6ae99SBarry Smith PetscFunctionBegin; 13647c6ae99SBarry Smith zctx.showgrid = PETSC_FALSE; 1377fe7c8feSLisandro Dalcin zctx.showaxis = PETSC_TRUE; 1388865f1eaSKarl Rupp 13947c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 1405b399a63SLisandro Dalcin ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 1415b399a63SLisandro Dalcin if (isnull) PetscFunctionReturn(0); 1425b399a63SLisandro Dalcin 14303193ff8SBarry Smith ierr = PetscViewerDrawGetBounds(viewer,&nbounds,&bounds);CHKERRQ(ierr); 14447c6ae99SBarry Smith 145c688c046SMatthew G Knepley ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 146ce94432eSBarry Smith if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 14747c6ae99SBarry Smith 14847c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr); 149e5ab1681SLisandro Dalcin ierr = MPI_Comm_rank(comm,&zctx.rank);CHKERRQ(ierr); 15047c6ae99SBarry Smith 1511321219cSEthan Coon ierr = DMDAGetInfo(da,0,&M,&N,0,&zctx.m,&zctx.n,0,&w,&s,&bx,&by,0,&st);CHKERRQ(ierr); 1520298fd71SBarry Smith ierr = DMDAGetOwnershipRanges(da,&lx,&ly,NULL);CHKERRQ(ierr); 15347c6ae99SBarry Smith 15447c6ae99SBarry Smith /* 15547c6ae99SBarry Smith Obtain a sequential vector that is going to contain the local values plus ONE layer of 156aa219208SBarry Smith ghosted values to draw the graphics from. We also need its corresponding DMDA (dac) that will 15747c6ae99SBarry Smith update the local values pluse ONE layer of ghost values. 15847c6ae99SBarry Smith */ 15947c6ae99SBarry Smith ierr = PetscObjectQuery((PetscObject)da,"GraphicsGhosted",(PetscObject*)&xlocal);CHKERRQ(ierr); 16047c6ae99SBarry Smith if (!xlocal) { 161bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_NONE || by != DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) { 16247c6ae99SBarry Smith /* 16347c6ae99SBarry Smith if original da is not of stencil width one, or periodic or not a box stencil then 164aa219208SBarry Smith create a special DMDA to handle one level of ghost points for graphics 16547c6ae99SBarry Smith */ 166bff4a2f0SMatthew 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); 167897f7067SBarry Smith ierr = DMSetUp(dac);CHKERRQ(ierr); 168aa219208SBarry Smith ierr = PetscInfo(da,"Creating auxilary DMDA for managing graphics ghost points\n");CHKERRQ(ierr); 16947c6ae99SBarry Smith } else { 17047c6ae99SBarry Smith /* otherwise we can use the da we already have */ 17147c6ae99SBarry Smith dac = da; 17247c6ae99SBarry Smith } 17347c6ae99SBarry Smith /* create local vector for holding ghosted values used in graphics */ 174564755cdSBarry Smith ierr = DMCreateLocalVector(dac,&xlocal);CHKERRQ(ierr); 17547c6ae99SBarry Smith if (dac != da) { 176aa219208SBarry Smith /* don't keep any public reference of this DMDA, is is only available through xlocal */ 177f7923d8aSBarry Smith ierr = PetscObjectDereference((PetscObject)dac);CHKERRQ(ierr); 17847c6ae99SBarry Smith } else { 17947c6ae99SBarry Smith /* remove association between xlocal and da, because below we compose in the opposite 18047c6ae99SBarry Smith direction and if we left this connect we'd get a loop, so the objects could 18147c6ae99SBarry Smith never be destroyed */ 182c688c046SMatthew G Knepley ierr = PetscObjectRemoveReference((PetscObject)xlocal,"__PETSc_dm");CHKERRQ(ierr); 18347c6ae99SBarry Smith } 18447c6ae99SBarry Smith ierr = PetscObjectCompose((PetscObject)da,"GraphicsGhosted",(PetscObject)xlocal);CHKERRQ(ierr); 18547c6ae99SBarry Smith ierr = PetscObjectDereference((PetscObject)xlocal);CHKERRQ(ierr); 18647c6ae99SBarry Smith } else { 187bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_NONE || by != DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) { 188c688c046SMatthew G Knepley ierr = VecGetDM(xlocal, &dac);CHKERRQ(ierr); 189f7923d8aSBarry Smith } else { 190f7923d8aSBarry Smith dac = da; 19147c6ae99SBarry Smith } 19247c6ae99SBarry Smith } 19347c6ae99SBarry Smith 19447c6ae99SBarry Smith /* 19547c6ae99SBarry Smith Get local (ghosted) values of vector 19647c6ae99SBarry Smith */ 1979a42bb27SBarry Smith ierr = DMGlobalToLocalBegin(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr); 1989a42bb27SBarry Smith ierr = DMGlobalToLocalEnd(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr); 1995edff71fSBarry Smith ierr = VecGetArrayRead(xlocal,&zctx.v);CHKERRQ(ierr); 20047c6ae99SBarry Smith 201832b7cebSLisandro Dalcin /* 202832b7cebSLisandro Dalcin Get coordinates of nodes 203832b7cebSLisandro Dalcin */ 2046636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); 20547c6ae99SBarry Smith if (!xcoor) { 206aa219208SBarry Smith ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,0.0);CHKERRQ(ierr); 2076636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); 20847c6ae99SBarry Smith } 20947c6ae99SBarry Smith 21047c6ae99SBarry Smith /* 21147c6ae99SBarry Smith Determine the min and max coordinates in plot 21247c6ae99SBarry Smith */ 213e5ab1681SLisandro Dalcin ierr = VecStrideMin(xcoor,0,NULL,&zctx.xmin);CHKERRQ(ierr); 214e5ab1681SLisandro Dalcin ierr = VecStrideMax(xcoor,0,NULL,&zctx.xmax);CHKERRQ(ierr); 215e5ab1681SLisandro Dalcin ierr = VecStrideMin(xcoor,1,NULL,&zctx.ymin);CHKERRQ(ierr); 216e5ab1681SLisandro Dalcin ierr = VecStrideMax(xcoor,1,NULL,&zctx.ymax);CHKERRQ(ierr); 2177fe7c8feSLisandro Dalcin ierr = PetscOptionsGetBool(NULL,NULL,"-draw_contour_axis",&zctx.showaxis,NULL);CHKERRQ(ierr); 2187fe7c8feSLisandro Dalcin if (zctx.showaxis) { 2197fe7c8feSLisandro Dalcin coors[0] = zctx.xmin - .05*(zctx.xmax - zctx.xmin); coors[1] = zctx.ymin - .05*(zctx.ymax - zctx.ymin); 2207fe7c8feSLisandro Dalcin coors[2] = zctx.xmax + .05*(zctx.xmax - zctx.xmin); coors[3] = zctx.ymax + .05*(zctx.ymax - zctx.ymin); 2217fe7c8feSLisandro Dalcin } else { 2227fe7c8feSLisandro Dalcin coors[0] = zctx.xmin; coors[1] = zctx.ymin; coors[2] = zctx.xmax; coors[3] = zctx.ymax; 2237fe7c8feSLisandro Dalcin } 224c5929fdfSBarry Smith ierr = PetscOptionsGetRealArray(NULL,NULL,"-draw_coordinates",coors,&ncoors,NULL);CHKERRQ(ierr); 225e5ab1681SLisandro 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); 226a74ba6f7SBarry Smith 22747c6ae99SBarry Smith /* 228832b7cebSLisandro Dalcin Get local ghosted version of coordinates 22947c6ae99SBarry Smith */ 23047c6ae99SBarry Smith ierr = PetscObjectQuery((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr); 23147c6ae99SBarry Smith if (!xcoorl) { 232aa219208SBarry Smith /* create DMDA to get local version of graphics */ 233bff4a2f0SMatthew 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); 234897f7067SBarry Smith ierr = DMSetUp(dag);CHKERRQ(ierr); 235aa219208SBarry Smith ierr = PetscInfo(dag,"Creating auxilary DMDA for managing graphics coordinates ghost points\n");CHKERRQ(ierr); 236564755cdSBarry Smith ierr = DMCreateLocalVector(dag,&xcoorl);CHKERRQ(ierr); 23747c6ae99SBarry Smith ierr = PetscObjectCompose((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject)xcoorl);CHKERRQ(ierr); 238f7923d8aSBarry Smith ierr = PetscObjectDereference((PetscObject)dag);CHKERRQ(ierr); 23947c6ae99SBarry Smith ierr = PetscObjectDereference((PetscObject)xcoorl);CHKERRQ(ierr); 24047c6ae99SBarry Smith } else { 241c688c046SMatthew G Knepley ierr = VecGetDM(xcoorl,&dag);CHKERRQ(ierr); 24247c6ae99SBarry Smith } 2439a42bb27SBarry Smith ierr = DMGlobalToLocalBegin(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr); 2449a42bb27SBarry Smith ierr = DMGlobalToLocalEnd(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr); 2455edff71fSBarry Smith ierr = VecGetArrayRead(xcoorl,&zctx.xy);CHKERRQ(ierr); 246832b7cebSLisandro Dalcin ierr = DMDAGetCoordinateName(da,0,&zctx.name0);CHKERRQ(ierr); 247832b7cebSLisandro Dalcin ierr = DMDAGetCoordinateName(da,1,&zctx.name1);CHKERRQ(ierr); 24847c6ae99SBarry Smith 24947c6ae99SBarry Smith /* 25047c6ae99SBarry Smith Get information about size of area each processor must do graphics for 25147c6ae99SBarry Smith */ 252e5ab1681SLisandro Dalcin ierr = DMDAGetInfo(dac,NULL,&M,&N,NULL,NULL,NULL,NULL,&zctx.dof,NULL,&bx,&by,NULL,NULL);CHKERRQ(ierr); 253e5ab1681SLisandro Dalcin ierr = DMDAGetGhostCorners(dac,NULL,NULL,NULL,&zctx.m,&zctx.n,NULL);CHKERRQ(ierr); 254c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-draw_contour_grid",&zctx.showgrid,NULL);CHKERRQ(ierr); 2554e6118eeSBarry Smith 2564e6118eeSBarry Smith ierr = DMDASelectFields(da,&ndisplayfields,&displayfields);CHKERRQ(ierr); 25747c6ae99SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 258c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-draw_ports",&useports,NULL);CHKERRQ(ierr); 259832b7cebSLisandro Dalcin if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE; 260832b7cebSLisandro Dalcin if (useports) { 261e5ab1681SLisandro Dalcin ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 2625b399a63SLisandro Dalcin ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr); 2635b399a63SLisandro Dalcin ierr = PetscDrawClear(draw);CHKERRQ(ierr); 26420d0051dSBarry Smith ierr = PetscDrawViewPortsCreate(draw,ndisplayfields,&ports);CHKERRQ(ierr); 265e5ab1681SLisandro Dalcin } 26620d0051dSBarry Smith 267832b7cebSLisandro Dalcin /* 268832b7cebSLisandro Dalcin Loop over each field; drawing each in a different window 269832b7cebSLisandro Dalcin */ 27020d0051dSBarry Smith for (i=0; i<ndisplayfields; i++) { 27120d0051dSBarry Smith zctx.k = displayfields[i]; 2725b399a63SLisandro Dalcin 273e5ab1681SLisandro Dalcin /* determine the min and max value in plot */ 2740298fd71SBarry Smith ierr = VecStrideMin(xin,zctx.k,NULL,&zctx.min);CHKERRQ(ierr); 2750298fd71SBarry Smith ierr = VecStrideMax(xin,zctx.k,NULL,&zctx.max);CHKERRQ(ierr); 27667dd0837SBarry Smith if (zctx.k < nbounds) { 277f3f0eb19SBarry Smith zctx.min = bounds[2*zctx.k]; 278f3f0eb19SBarry Smith zctx.max = bounds[2*zctx.k+1]; 27967dd0837SBarry Smith } 28047c6ae99SBarry Smith if (zctx.min == zctx.max) { 28147c6ae99SBarry Smith zctx.min -= 1.e-12; 28247c6ae99SBarry Smith zctx.max += 1.e-12; 28347c6ae99SBarry Smith } 28457622a8eSBarry Smith ierr = PetscInfo2(da,"DMDA 2d contour plot min %g max %g\n",(double)zctx.min,(double)zctx.max);CHKERRQ(ierr); 28547c6ae99SBarry Smith 286832b7cebSLisandro Dalcin if (useports) { 287832b7cebSLisandro Dalcin ierr = PetscDrawViewPortsSet(ports,i);CHKERRQ(ierr); 288832b7cebSLisandro Dalcin } else { 289832b7cebSLisandro Dalcin const char *title; 290832b7cebSLisandro Dalcin ierr = PetscViewerDrawGetDraw(viewer,i,&draw);CHKERRQ(ierr); 291832b7cebSLisandro Dalcin ierr = DMDAGetFieldName(da,zctx.k,&title);CHKERRQ(ierr); 292832b7cebSLisandro Dalcin if (title) {ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr);} 293832b7cebSLisandro Dalcin } 294832b7cebSLisandro Dalcin 29547c6ae99SBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 29645f3bb6eSLisandro Dalcin ierr = PetscDrawScalePopup(popup,zctx.min,zctx.max);CHKERRQ(ierr); 297e5ab1681SLisandro Dalcin ierr = PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr); 29847c6ae99SBarry Smith ierr = PetscDrawZoom(draw,VecView_MPI_Draw_DA2d_Zoom,&zctx);CHKERRQ(ierr); 299832b7cebSLisandro Dalcin if (!useports) {ierr = PetscDrawSave(draw);CHKERRQ(ierr);} 30047c6ae99SBarry Smith } 301832b7cebSLisandro Dalcin if (useports) { 302832b7cebSLisandro Dalcin ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 303832b7cebSLisandro Dalcin ierr = PetscDrawSave(draw);CHKERRQ(ierr); 304832b7cebSLisandro Dalcin } 305832b7cebSLisandro Dalcin 3066bf464f9SBarry Smith ierr = PetscDrawViewPortsDestroy(ports);CHKERRQ(ierr); 307e5ab1681SLisandro Dalcin ierr = PetscFree(displayfields);CHKERRQ(ierr); 3085edff71fSBarry Smith ierr = VecRestoreArrayRead(xcoorl,&zctx.xy);CHKERRQ(ierr); 3095edff71fSBarry Smith ierr = VecRestoreArrayRead(xlocal,&zctx.v);CHKERRQ(ierr); 31047c6ae99SBarry Smith PetscFunctionReturn(0); 31147c6ae99SBarry Smith } 31247c6ae99SBarry Smith 3130da24e51SJuha Jäykkä #if defined(PETSC_HAVE_HDF5) 314c73cfb54SMatthew G. Knepley static PetscErrorCode VecGetHDF5ChunkSize(DM_DA *da, Vec xin, PetscInt dimension, PetscInt timestep, hsize_t *chunkDims) 3150da24e51SJuha Jäykkä { 316d4190030SJed Brown PetscMPIInt comm_size; 317d4190030SJed Brown PetscErrorCode ierr; 318561a82e7SJed Brown hsize_t chunk_size, target_size, dim; 319561a82e7SJed Brown hsize_t vec_size = sizeof(PetscScalar)*da->M*da->N*da->P*da->w; 320561a82e7SJed Brown hsize_t avg_local_vec_size,KiB = 1024,MiB = KiB*KiB,GiB = MiB*KiB,min_size = MiB; 321561a82e7SJed Brown hsize_t max_chunks = 64*KiB; /* HDF5 internal limitation */ 322561a82e7SJed Brown hsize_t max_chunk_size = 4*GiB; /* HDF5 internal limitation */ 32384179ffaSJed Brown PetscInt zslices=da->p, yslices=da->n, xslices=da->m; 3240da24e51SJuha Jäykkä 325cb5c4f94SJuha Jäykkä PetscFunctionBegin; 326d4190030SJed Brown ierr = MPI_Comm_size(PetscObjectComm((PetscObject)xin), &comm_size);CHKERRQ(ierr); 3270da24e51SJuha 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 */ 3280da24e51SJuha Jäykkä 329794a961bSBarry 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))))); 3307d310018SBarry 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 */ 3317d310018SBarry 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); 3320da24e51SJuha Jäykkä 333cb5c4f94SJuha Jäykkä /* 334cb5c4f94SJuha Jäykkä if size/rank > max_chunk_size, we need radical measures: even going down to 335cb5c4f94SJuha Jäykkä avg_local_vec_size is not enough, so we simply use chunk size of 4 GiB no matter 336cb5c4f94SJuha Jäykkä what, composed in the most efficient way possible. 337cb5c4f94SJuha Jäykkä N.B. this minimises the number of chunks, which may or may not be the optimal 338cb5c4f94SJuha Jäykkä solution. In a BG, for example, the optimal solution is probably to make # chunks = # 339cb5c4f94SJuha Jäykkä IO nodes involved, but this author has no access to a BG to figure out how to 340cb5c4f94SJuha Jäykkä reliably find the right number. And even then it may or may not be enough. 341cb5c4f94SJuha Jäykkä */ 3420da24e51SJuha Jäykkä if (avg_local_vec_size > max_chunk_size) { 343cb5c4f94SJuha Jäykkä /* check if we can just split local z-axis: is that enough? */ 34484179ffaSJed Brown zslices = (PetscInt)ceil(vec_size*1.0/(da->p*max_chunk_size))*zslices; 3450da24e51SJuha Jäykkä if (zslices > da->P) { 346cb5c4f94SJuha Jäykkä /* lattice is too large in xy-directions, splitting z only is not enough */ 3470da24e51SJuha Jäykkä zslices = da->P; 34884179ffaSJed Brown yslices= (PetscInt)ceil(vec_size*1.0/(zslices*da->n*max_chunk_size))*yslices; 3490da24e51SJuha Jäykkä if (yslices > da->N) { 350cb5c4f94SJuha Jäykkä /* lattice is too large in x-direction, splitting along z, y is not enough */ 3510da24e51SJuha Jäykkä yslices = da->N; 35284179ffaSJed Brown xslices= (PetscInt)ceil(vec_size*1.0/(zslices*yslices*da->m*max_chunk_size))*xslices; 3530da24e51SJuha Jäykkä } 3540da24e51SJuha Jäykkä } 3550da24e51SJuha Jäykkä dim = 0; 3560da24e51SJuha Jäykkä if (timestep >= 0) { 3570da24e51SJuha Jäykkä ++dim; 3580da24e51SJuha Jäykkä } 359cb5c4f94SJuha Jäykkä /* prefer to split z-axis, even down to planar slices */ 360c73cfb54SMatthew G. Knepley if (dimension == 3) { 361cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t) da->P/zslices; 362cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t) da->N/yslices; 363cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t) da->M/xslices; 3640da24e51SJuha Jäykkä } else { 365cb5c4f94SJuha Jäykkä /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */ 366cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t) da->N/yslices; 367cb5c4f94SJuha Jäykkä chunkDims[dim++] = (hsize_t) da->M/xslices; 3680da24e51SJuha Jäykkä } 3690da24e51SJuha 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); 3700da24e51SJuha Jäykkä } else { 3710da24e51SJuha Jäykkä if (target_size < chunk_size) { 372cb5c4f94SJuha Jäykkä /* only change the defaults if target_size < chunk_size */ 3730da24e51SJuha Jäykkä dim = 0; 3740da24e51SJuha Jäykkä if (timestep >= 0) { 3750da24e51SJuha Jäykkä ++dim; 3760da24e51SJuha Jäykkä } 377cb5c4f94SJuha Jäykkä /* prefer to split z-axis, even down to planar slices */ 378c73cfb54SMatthew G. Knepley if (dimension == 3) { 379cb5c4f94SJuha Jäykkä /* try splitting the z-axis to core-size bits, i.e. divide chunk size by # comm_size in z-direction */ 3800da24e51SJuha Jäykkä if (target_size >= chunk_size/da->p) { 381cb5c4f94SJuha Jäykkä /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */ 3820da24e51SJuha Jäykkä chunkDims[dim] = (hsize_t) ceil(da->P*1.0/da->p); 3830da24e51SJuha Jäykkä } else { 384cb5c4f94SJuha Jäykkä /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be 385cb5c4f94SJuha Jäykkä radical and let everyone write all they've got */ 3860da24e51SJuha Jäykkä chunkDims[dim++] = (hsize_t) ceil(da->P*1.0/da->p); 3870da24e51SJuha Jäykkä chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n); 3880da24e51SJuha Jäykkä chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m); 3890da24e51SJuha Jäykkä } 3900da24e51SJuha Jäykkä } else { 391cb5c4f94SJuha Jäykkä /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */ 3920da24e51SJuha Jäykkä if (target_size >= chunk_size/da->n) { 393cb5c4f94SJuha Jäykkä /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */ 3940da24e51SJuha Jäykkä chunkDims[dim] = (hsize_t) ceil(da->N*1.0/da->n); 3950da24e51SJuha Jäykkä } else { 396cb5c4f94SJuha Jäykkä /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be 397cb5c4f94SJuha Jäykkä radical and let everyone write all they've got */ 3980da24e51SJuha Jäykkä chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n); 3990da24e51SJuha Jäykkä chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m); 4000da24e51SJuha Jäykkä } 4010da24e51SJuha Jäykkä 4020da24e51SJuha Jäykkä } 4030da24e51SJuha 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); 4040da24e51SJuha Jäykkä } else { 405cb5c4f94SJuha Jäykkä /* precomputed chunks are fine, we don't need to do anything */ 4060da24e51SJuha Jäykkä } 4070da24e51SJuha Jäykkä } 4080da24e51SJuha Jäykkä PetscFunctionReturn(0); 4090da24e51SJuha Jäykkä } 4100da24e51SJuha Jäykkä #endif 41147c6ae99SBarry Smith 41247c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 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." 502570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16) 503570b7f6dSBarry Smith #error "HDF5 output with 16 bit floats not supported." 50415214e8eSMatthew G Knepley #else 5059a0502c6SHåkon Strandenes memscalartype = H5T_NATIVE_DOUBLE; 5069a0502c6SHåkon Strandenes if (spoutput == PETSC_TRUE) filescalartype = H5T_NATIVE_FLOAT; 5079a0502c6SHåkon Strandenes else filescalartype = H5T_NATIVE_DOUBLE; 50815214e8eSMatthew G Knepley #endif 50915214e8eSMatthew G Knepley 51047c6ae99SBarry Smith /* Create the dataset with default properties and close filespace */ 51147c6ae99SBarry Smith ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr); 51215214e8eSMatthew G Knepley if (!H5Lexists(group, vecname, H5P_DEFAULT)) { 5138e2ae6d7SMichael Kraus /* Create chunk */ 514729ab6d9SBarry Smith PetscStackCallHDF5Return(chunkspace,H5Pcreate,(H5P_DATASET_CREATE)); 515729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_chunk,(chunkspace, dim, chunkDims)); 5168e2ae6d7SMichael Kraus 5179a0502c6SHåkon Strandenes PetscStackCallHDF5Return(dset_id,H5Dcreate2,(group, vecname, filescalartype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT)); 51815214e8eSMatthew G Knepley } else { 519729ab6d9SBarry Smith PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT)); 520729ab6d9SBarry Smith PetscStackCallHDF5(H5Dset_extent,(dset_id, dims)); 52115214e8eSMatthew G Knepley } 522729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(filespace)); 52347c6ae99SBarry Smith 52447c6ae99SBarry Smith /* Each process defines a dataset and writes it to the hyperslab in the file */ 5258e2ae6d7SMichael Kraus dim = 0; 5268e2ae6d7SMichael Kraus if (timestep >= 0) { 5278e2ae6d7SMichael Kraus offset[dim] = timestep; 5288e2ae6d7SMichael Kraus ++dim; 5298e2ae6d7SMichael Kraus } 530c73cfb54SMatthew G. Knepley if (dimension == 3) {ierr = PetscHDF5IntCast(da->zs,offset + dim++);CHKERRQ(ierr);} 531c73cfb54SMatthew G. Knepley if (dimension > 1) {ierr = PetscHDF5IntCast(da->ys,offset + dim++);CHKERRQ(ierr);} 532acba2ac6SBarry Smith ierr = PetscHDF5IntCast(da->xs/da->w,offset + dim++);CHKERRQ(ierr); 53382ea9b62SBarry Smith if (da->w > 1 || dim2) offset[dim++] = 0; 53447c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 5358e2ae6d7SMichael Kraus offset[dim++] = 0; 53647c6ae99SBarry Smith #endif 5378e2ae6d7SMichael Kraus dim = 0; 5388e2ae6d7SMichael Kraus if (timestep >= 0) { 5398e2ae6d7SMichael Kraus count[dim] = 1; 5408e2ae6d7SMichael Kraus ++dim; 5418e2ae6d7SMichael Kraus } 542c73cfb54SMatthew G. Knepley if (dimension == 3) {ierr = PetscHDF5IntCast(da->ze - da->zs,count + dim++);CHKERRQ(ierr);} 543c73cfb54SMatthew G. Knepley if (dimension > 1) {ierr = PetscHDF5IntCast(da->ye - da->ys,count + dim++);CHKERRQ(ierr);} 544acba2ac6SBarry Smith ierr = PetscHDF5IntCast((da->xe - da->xs)/da->w,count + dim++);CHKERRQ(ierr); 54582ea9b62SBarry Smith if (da->w > 1 || dim2) {ierr = PetscHDF5IntCast(da->w,count + dim++);CHKERRQ(ierr);} 54647c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 5478e2ae6d7SMichael Kraus count[dim++] = 2; 54847c6ae99SBarry Smith #endif 549729ab6d9SBarry Smith PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL)); 550729ab6d9SBarry Smith PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id)); 551729ab6d9SBarry Smith PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL)); 55247c6ae99SBarry Smith 55347c6ae99SBarry Smith /* Create property list for collective dataset write */ 554729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_DATASET_XFER)); 55547c6ae99SBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 556729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_dxpl_mpio,(plist_id, H5FD_MPIO_COLLECTIVE)); 55747c6ae99SBarry Smith #endif 55847c6ae99SBarry Smith /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 55947c6ae99SBarry Smith 5605edff71fSBarry Smith ierr = VecGetArrayRead(xin, &x);CHKERRQ(ierr); 5619a0502c6SHåkon Strandenes PetscStackCallHDF5(H5Dwrite,(dset_id, memscalartype, memspace, filespace, plist_id, x)); 562729ab6d9SBarry Smith PetscStackCallHDF5(H5Fflush,(file_id, H5F_SCOPE_GLOBAL)); 5635edff71fSBarry Smith ierr = VecRestoreArrayRead(xin, &x);CHKERRQ(ierr); 56447c6ae99SBarry Smith 56547c6ae99SBarry Smith /* Close/release resources */ 56615214e8eSMatthew G Knepley if (group != file_id) { 567729ab6d9SBarry Smith PetscStackCallHDF5(H5Gclose,(group)); 56815214e8eSMatthew G Knepley } 569729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 570729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(filespace)); 571729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(memspace)); 572729ab6d9SBarry Smith PetscStackCallHDF5(H5Dclose,(dset_id)); 57347c6ae99SBarry Smith ierr = PetscInfo1(xin,"Wrote Vec object with name %s\n",vecname);CHKERRQ(ierr); 57447c6ae99SBarry Smith PetscFunctionReturn(0); 57547c6ae99SBarry Smith } 57647c6ae99SBarry Smith #endif 57747c6ae99SBarry Smith 57809573ac7SBarry Smith extern PetscErrorCode VecView_MPI_Draw_DA1d(Vec,PetscViewer); 57947c6ae99SBarry Smith 58047c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 581aa219208SBarry Smith static PetscErrorCode DMDAArrayMPIIO(DM da,PetscViewer viewer,Vec xin,PetscBool write) 58247c6ae99SBarry Smith { 58347c6ae99SBarry Smith PetscErrorCode ierr; 58447c6ae99SBarry Smith MPI_File mfdes; 58547c6ae99SBarry Smith PetscMPIInt gsizes[4],lsizes[4],lstarts[4],asiz,dof; 58647c6ae99SBarry Smith MPI_Datatype view; 58747c6ae99SBarry Smith const PetscScalar *array; 58847c6ae99SBarry Smith MPI_Offset off; 58947c6ae99SBarry Smith MPI_Aint ub,ul; 59047c6ae99SBarry Smith PetscInt type,rows,vecrows,tr[2]; 59147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 5920c169764Sdmay PetscBool skipheader; 59347c6ae99SBarry Smith 59447c6ae99SBarry Smith PetscFunctionBegin; 59547c6ae99SBarry Smith ierr = VecGetSize(xin,&vecrows);CHKERRQ(ierr); 5960c169764Sdmay ierr = PetscViewerBinaryGetSkipHeader(viewer,&skipheader);CHKERRQ(ierr); 59747c6ae99SBarry Smith if (!write) { 59847c6ae99SBarry Smith /* Read vector header. */ 5990c169764Sdmay if (!skipheader) { 600060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,tr,2,NULL,PETSC_INT);CHKERRQ(ierr); 60147c6ae99SBarry Smith type = tr[0]; 60247c6ae99SBarry Smith rows = tr[1]; 603ce94432eSBarry Smith if (type != VEC_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Not vector next in file"); 604ce94432eSBarry Smith if (rows != vecrows) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Vector in file not same size as DMDA vector"); 6050c169764Sdmay } 60647c6ae99SBarry Smith } else { 60747c6ae99SBarry Smith tr[0] = VEC_FILE_CLASSID; 60847c6ae99SBarry Smith tr[1] = vecrows; 6090c169764Sdmay if (!skipheader) { 61047c6ae99SBarry Smith ierr = PetscViewerBinaryWrite(viewer,tr,2,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 61147c6ae99SBarry Smith } 6120c169764Sdmay } 61347c6ae99SBarry Smith 6144dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->w,&dof);CHKERRQ(ierr); 6154dc2109aSBarry Smith gsizes[0] = dof; 6164dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->M,gsizes+1);CHKERRQ(ierr); 6174dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->N,gsizes+2);CHKERRQ(ierr); 618334634e2SJed Brown ierr = PetscMPIIntCast(dd->P,gsizes+3);CHKERRQ(ierr); 6194dc2109aSBarry Smith lsizes[0] = dof; 6204dc2109aSBarry Smith ierr = PetscMPIIntCast((dd->xe-dd->xs)/dof,lsizes+1);CHKERRQ(ierr); 6214dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->ye-dd->ys,lsizes+2);CHKERRQ(ierr); 6224dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->ze-dd->zs,lsizes+3);CHKERRQ(ierr); 6234dc2109aSBarry Smith lstarts[0] = 0; 6244dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->xs/dof,lstarts+1);CHKERRQ(ierr); 6254dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->ys,lstarts+2);CHKERRQ(ierr); 6264dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->zs,lstarts+3);CHKERRQ(ierr); 627c73cfb54SMatthew G. Knepley ierr = MPI_Type_create_subarray(da->dim+1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view);CHKERRQ(ierr); 62847c6ae99SBarry Smith ierr = MPI_Type_commit(&view);CHKERRQ(ierr); 62947c6ae99SBarry Smith 63047c6ae99SBarry Smith ierr = PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);CHKERRQ(ierr); 63147c6ae99SBarry Smith ierr = PetscViewerBinaryGetMPIIOOffset(viewer,&off);CHKERRQ(ierr); 63247c6ae99SBarry Smith ierr = MPI_File_set_view(mfdes,off,MPIU_SCALAR,view,(char*)"native",MPI_INFO_NULL);CHKERRQ(ierr); 63347c6ae99SBarry Smith ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr); 63447c6ae99SBarry Smith asiz = lsizes[1]*(lsizes[2] > 0 ? lsizes[2] : 1)*(lsizes[3] > 0 ? lsizes[3] : 1)*dof; 63547c6ae99SBarry Smith if (write) { 63647c6ae99SBarry Smith ierr = MPIU_File_write_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr); 63747c6ae99SBarry Smith } else { 63847c6ae99SBarry Smith ierr = MPIU_File_read_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr); 63947c6ae99SBarry Smith } 64047c6ae99SBarry Smith ierr = MPI_Type_get_extent(view,&ul,&ub);CHKERRQ(ierr); 64147c6ae99SBarry Smith ierr = PetscViewerBinaryAddMPIIOOffset(viewer,ub);CHKERRQ(ierr); 64247c6ae99SBarry Smith ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr); 64347c6ae99SBarry Smith ierr = MPI_Type_free(&view);CHKERRQ(ierr); 64447c6ae99SBarry Smith PetscFunctionReturn(0); 64547c6ae99SBarry Smith } 64647c6ae99SBarry Smith #endif 64747c6ae99SBarry Smith 6487087cfbeSBarry Smith PetscErrorCode VecView_MPI_DA(Vec xin,PetscViewer viewer) 64947c6ae99SBarry Smith { 6509a42bb27SBarry Smith DM da; 65147c6ae99SBarry Smith PetscErrorCode ierr; 65247c6ae99SBarry Smith PetscInt dim; 65347c6ae99SBarry Smith Vec natural; 6548135c375SStefano Zampini PetscBool isdraw,isvtk,isglvis; 65547c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 65647c6ae99SBarry Smith PetscBool ishdf5; 65747c6ae99SBarry Smith #endif 6583f3fd955SJed Brown const char *prefix,*name; 659a261c58fSBarry Smith PetscViewerFormat format; 66047c6ae99SBarry Smith 66147c6ae99SBarry Smith PetscFunctionBegin; 662c688c046SMatthew G Knepley ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 663ce94432eSBarry Smith if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 664251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 665251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr); 66647c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 667251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr); 66847c6ae99SBarry Smith #endif 6698135c375SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr); 67047c6ae99SBarry Smith if (isdraw) { 6711321219cSEthan Coon ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 67247c6ae99SBarry Smith if (dim == 1) { 67347c6ae99SBarry Smith ierr = VecView_MPI_Draw_DA1d(xin,viewer);CHKERRQ(ierr); 67447c6ae99SBarry Smith } else if (dim == 2) { 67547c6ae99SBarry Smith ierr = VecView_MPI_Draw_DA2d(xin,viewer);CHKERRQ(ierr); 676ce94432eSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot graphically view vector associated with this dimensional DMDA %D",dim); 67701d7c5c3SPatrick Sanan } else if (isvtk) { /* Duplicate the Vec */ 6784061b8bfSJed Brown Vec Y; 6794061b8bfSJed Brown ierr = VecDuplicate(xin,&Y);CHKERRQ(ierr); 680b51b94faSJed Brown if (((PetscObject)xin)->name) { 681b51b94faSJed Brown /* If xin was named, copy the name over to Y. The duplicate names are safe because nobody else will ever see Y. */ 682b51b94faSJed Brown ierr = PetscObjectSetName((PetscObject)Y,((PetscObject)xin)->name);CHKERRQ(ierr); 683b51b94faSJed Brown } 6844061b8bfSJed Brown ierr = VecCopy(xin,Y);CHKERRQ(ierr); 685a8f87f1dSPatrick Sanan { 686a8f87f1dSPatrick Sanan PetscObject dmvtk; 687a8f87f1dSPatrick Sanan PetscBool compatible,compatibleSet; 688a8f87f1dSPatrick Sanan ierr = PetscViewerVTKGetDM(viewer,&dmvtk);CHKERRQ(ierr); 689a8f87f1dSPatrick Sanan if (dmvtk) { 690a8f87f1dSPatrick Sanan PetscValidHeaderSpecific((DM)dmvtk,DM_CLASSID,-1); 691a8f87f1dSPatrick Sanan ierr = DMGetCompatibility(da,(DM)dmvtk,&compatible,&compatibleSet);CHKERRQ(ierr); 692a8f87f1dSPatrick Sanan if (!compatibleSet || !compatible) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"Cannot confirm compatibility of DMs associated with Vecs viewed in the same VTK file. Check that grids are the same."); 693a8f87f1dSPatrick Sanan } 694a8f87f1dSPatrick Sanan ierr = PetscViewerVTKAddField(viewer,(PetscObject)da,DMDAVTKWriteAll,PETSC_VTK_POINT_FIELD,PETSC_FALSE,(PetscObject)Y);CHKERRQ(ierr); 695a8f87f1dSPatrick Sanan } 69647c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 69747c6ae99SBarry Smith } else if (ishdf5) { 69847c6ae99SBarry Smith ierr = VecView_MPI_HDF5_DA(xin,viewer);CHKERRQ(ierr); 69947c6ae99SBarry Smith #endif 7008135c375SStefano Zampini } else if (isglvis) { 7018135c375SStefano Zampini ierr = VecView_GLVis(xin,viewer);CHKERRQ(ierr); 70247c6ae99SBarry Smith } else { 70347c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 70447c6ae99SBarry Smith PetscBool isbinary,isMPIIO; 70547c6ae99SBarry Smith 706251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 70747c6ae99SBarry Smith if (isbinary) { 708bc196f7cSDave May ierr = PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO);CHKERRQ(ierr); 70947c6ae99SBarry Smith if (isMPIIO) { 710aa219208SBarry Smith ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_TRUE);CHKERRQ(ierr); 71147c6ae99SBarry Smith PetscFunctionReturn(0); 71247c6ae99SBarry Smith } 71347c6ae99SBarry Smith } 71447c6ae99SBarry Smith #endif 71547c6ae99SBarry Smith 71647c6ae99SBarry Smith /* call viewer on natural ordering */ 71747c6ae99SBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr); 718aa219208SBarry Smith ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr); 71947c6ae99SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr); 720aa219208SBarry Smith ierr = DMDAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr); 721aa219208SBarry Smith ierr = DMDAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr); 7223f3fd955SJed Brown ierr = PetscObjectGetName((PetscObject)xin,&name);CHKERRQ(ierr); 7233f3fd955SJed Brown ierr = PetscObjectSetName((PetscObject)natural,name);CHKERRQ(ierr); 724a261c58fSBarry Smith 725a261c58fSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 726a261c58fSBarry Smith if (format == PETSC_VIEWER_BINARY_MATLAB) { 727a261c58fSBarry Smith /* temporarily remove viewer format so it won't trigger in the VecView() */ 7286a9046bcSBarry Smith ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr); 729a261c58fSBarry Smith } 730a261c58fSBarry Smith 73123f88b65SBarry Smith ((PetscObject)natural)->donotPetscObjectPrintClassNamePrefixType = PETSC_TRUE; 73247c6ae99SBarry Smith ierr = VecView(natural,viewer);CHKERRQ(ierr); 73323f88b65SBarry Smith ((PetscObject)natural)->donotPetscObjectPrintClassNamePrefixType = PETSC_FALSE; 734a261c58fSBarry Smith 735a261c58fSBarry Smith if (format == PETSC_VIEWER_BINARY_MATLAB) { 736a261c58fSBarry Smith MPI_Comm comm; 737a261c58fSBarry Smith FILE *info; 738a261c58fSBarry Smith const char *fieldname; 739da88d4d4SJed Brown char fieldbuf[256]; 740a261c58fSBarry Smith PetscInt dim,ni,nj,nk,pi,pj,pk,dof,n; 741a261c58fSBarry Smith 742a261c58fSBarry Smith /* set the viewer format back into the viewer */ 7436a9046bcSBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 744a261c58fSBarry Smith ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); 745a261c58fSBarry Smith ierr = PetscViewerBinaryGetInfoPointer(viewer,&info);CHKERRQ(ierr); 746a261c58fSBarry Smith ierr = DMDAGetInfo(da,&dim,&ni,&nj,&nk,&pi,&pj,&pk,&dof,0,0,0,0,0);CHKERRQ(ierr); 747da88d4d4SJed Brown ierr = PetscFPrintf(comm,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n");CHKERRQ(ierr); 748da88d4d4SJed Brown ierr = PetscFPrintf(comm,info,"#$$ tmp = PetscBinaryRead(fd); \n");CHKERRQ(ierr); 749da88d4d4SJed Brown if (dim == 1) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d);\n",dof,ni);CHKERRQ(ierr); } 750da88d4d4SJed Brown if (dim == 2) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d);\n",dof,ni,nj);CHKERRQ(ierr); } 751da88d4d4SJed Brown if (dim == 3) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d,%d);\n",dof,ni,nj,nk);CHKERRQ(ierr); } 752a261c58fSBarry Smith 753a261c58fSBarry Smith for (n=0; n<dof; n++) { 754a261c58fSBarry Smith ierr = DMDAGetFieldName(da,n,&fieldname);CHKERRQ(ierr); 755da88d4d4SJed Brown if (!fieldname || !fieldname[0]) { 756da88d4d4SJed Brown ierr = PetscSNPrintf(fieldbuf,sizeof fieldbuf,"field%D",n);CHKERRQ(ierr); 757da88d4d4SJed Brown fieldname = fieldbuf; 758a261c58fSBarry Smith } 759da88d4d4SJed Brown if (dim == 1) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:))';\n",name,fieldname,n+1);CHKERRQ(ierr); } 760da88d4d4SJed Brown if (dim == 2) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:,:))';\n",name,fieldname,n+1);CHKERRQ(ierr); } 761da88d4d4SJed 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);} 762a261c58fSBarry Smith } 763da88d4d4SJed Brown ierr = PetscFPrintf(comm,info,"#$$ clear tmp; \n");CHKERRQ(ierr); 764da88d4d4SJed Brown ierr = PetscFPrintf(comm,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n");CHKERRQ(ierr); 765a261c58fSBarry Smith } 766a261c58fSBarry Smith 767fcfd50ebSBarry Smith ierr = VecDestroy(&natural);CHKERRQ(ierr); 76847c6ae99SBarry Smith } 76947c6ae99SBarry Smith PetscFunctionReturn(0); 77047c6ae99SBarry Smith } 77147c6ae99SBarry Smith 77247c6ae99SBarry Smith #if defined(PETSC_HAVE_HDF5) 77347c6ae99SBarry Smith PetscErrorCode VecLoad_HDF5_DA(Vec xin, PetscViewer viewer) 77447c6ae99SBarry Smith { 7759a42bb27SBarry Smith DM da; 77647c6ae99SBarry Smith PetscErrorCode ierr; 777*15b861d2SVaclav Hapla int dim,rdim; 778ec7ae49cSHåkon Strandenes hsize_t dims[6]={0},count[6]={0},offset[6]={0}; 779ec7ae49cSHåkon Strandenes PetscInt dimension,timestep,dofInd; 78047c6ae99SBarry Smith PetscScalar *x; 78147c6ae99SBarry Smith const char *vecname; 78247c6ae99SBarry Smith hid_t filespace; /* file dataspace identifier */ 78347c6ae99SBarry Smith hid_t plist_id; /* property list identifier */ 78447c6ae99SBarry Smith hid_t dset_id; /* dataset identifier */ 78547c6ae99SBarry Smith hid_t memspace; /* memory dataspace identifier */ 786bfd264e7SBarry Smith hid_t file_id,group; 7877bcbaff4SBarry Smith hid_t scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */ 7889c7c4993SBarry Smith DM_DA *dd; 789ec7ae49cSHåkon Strandenes PetscBool dim2 = PETSC_FALSE; 79047c6ae99SBarry Smith 79147c6ae99SBarry Smith PetscFunctionBegin; 7927bcbaff4SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 7937bcbaff4SBarry Smith scalartype = H5T_NATIVE_FLOAT; 7947bcbaff4SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128) 7957bcbaff4SBarry Smith #error "HDF5 output with 128 bit floats not supported." 796570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16) 797570b7f6dSBarry Smith #error "HDF5 output with 16 bit floats not supported." 7987bcbaff4SBarry Smith #else 7997bcbaff4SBarry Smith scalartype = H5T_NATIVE_DOUBLE; 8007bcbaff4SBarry Smith #endif 8017bcbaff4SBarry Smith 802bfd264e7SBarry Smith ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr); 803ec7ae49cSHåkon Strandenes ierr = PetscViewerHDF5GetTimestep(viewer, ×tep);CHKERRQ(ierr); 804ec7ae49cSHåkon Strandenes ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr); 805c688c046SMatthew G Knepley ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 8069c7c4993SBarry Smith dd = (DM_DA*)da->data; 807c73cfb54SMatthew G. Knepley ierr = DMGetDimension(da, &dimension);CHKERRQ(ierr); 80847c6ae99SBarry Smith 809ec7ae49cSHåkon Strandenes /* Open dataset */ 810729ab6d9SBarry Smith PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT)); 81147c6ae99SBarry Smith 812ec7ae49cSHåkon Strandenes /* Retrieve the dataspace for the dataset */ 813ec7ae49cSHåkon Strandenes PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id)); 814ec7ae49cSHåkon Strandenes PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(filespace, dims, NULL)); 815ec7ae49cSHåkon Strandenes 816ec7ae49cSHåkon Strandenes /* Expected dimension for holding the dof's */ 81747c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 818ec7ae49cSHåkon Strandenes dofInd = rdim-2; 819ec7ae49cSHåkon Strandenes #else 820ec7ae49cSHåkon Strandenes dofInd = rdim-1; 82147c6ae99SBarry Smith #endif 822ec7ae49cSHåkon Strandenes 823ec7ae49cSHåkon Strandenes /* The expected number of dimensions, assuming basedimension2 = false */ 824ec7ae49cSHåkon Strandenes dim = dimension; 825ec7ae49cSHåkon Strandenes if (dd->w > 1) ++dim; 826ec7ae49cSHåkon Strandenes if (timestep >= 0) ++dim; 82747c6ae99SBarry Smith #if defined(PETSC_USE_COMPLEX) 828ec7ae49cSHåkon Strandenes ++dim; 82947c6ae99SBarry Smith #endif 830ec7ae49cSHåkon Strandenes 831ec7ae49cSHåkon Strandenes /* In this case the input dataset have one extra, unexpected dimension. */ 832ec7ae49cSHåkon Strandenes if (rdim == dim+1) { 833ec7ae49cSHåkon Strandenes /* In this case the block size unity */ 834ec7ae49cSHåkon Strandenes if (dd->w == 1 && dims[dofInd] == 1) dim2 = PETSC_TRUE; 835ec7ae49cSHåkon Strandenes 836ec7ae49cSHåkon Strandenes /* Special error message for the case where dof does not match the input file */ 837ec7ae49cSHå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); 838ec7ae49cSHåkon Strandenes 839ec7ae49cSHåkon Strandenes /* Other cases where rdim != dim cannot be handled currently */ 8406c4ed002SBarry 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); 841ec7ae49cSHåkon Strandenes 842ec7ae49cSHåkon Strandenes /* Set up the hyperslab size */ 843ec7ae49cSHåkon Strandenes dim = 0; 844ec7ae49cSHåkon Strandenes if (timestep >= 0) { 845ec7ae49cSHåkon Strandenes offset[dim] = timestep; 846ec7ae49cSHåkon Strandenes count[dim] = 1; 847ec7ae49cSHåkon Strandenes ++dim; 848ec7ae49cSHåkon Strandenes } 849ec7ae49cSHåkon Strandenes if (dimension == 3) { 850ec7ae49cSHåkon Strandenes ierr = PetscHDF5IntCast(dd->zs,offset + dim);CHKERRQ(ierr); 851ec7ae49cSHåkon Strandenes ierr = PetscHDF5IntCast(dd->ze - dd->zs,count + dim);CHKERRQ(ierr); 852ec7ae49cSHåkon Strandenes ++dim; 853ec7ae49cSHåkon Strandenes } 854ec7ae49cSHåkon Strandenes if (dimension > 1) { 855ec7ae49cSHåkon Strandenes ierr = PetscHDF5IntCast(dd->ys,offset + dim);CHKERRQ(ierr); 856ec7ae49cSHåkon Strandenes ierr = PetscHDF5IntCast(dd->ye - dd->ys,count + dim);CHKERRQ(ierr); 857ec7ae49cSHåkon Strandenes ++dim; 858ec7ae49cSHåkon Strandenes } 859ec7ae49cSHåkon Strandenes ierr = PetscHDF5IntCast(dd->xs/dd->w,offset + dim);CHKERRQ(ierr); 860ec7ae49cSHåkon Strandenes ierr = PetscHDF5IntCast((dd->xe - dd->xs)/dd->w,count + dim);CHKERRQ(ierr); 861ec7ae49cSHåkon Strandenes ++dim; 862ec7ae49cSHåkon Strandenes if (dd->w > 1 || dim2) { 863ec7ae49cSHåkon Strandenes offset[dim] = 0; 864ec7ae49cSHåkon Strandenes ierr = PetscHDF5IntCast(dd->w,count + dim);CHKERRQ(ierr); 865ec7ae49cSHåkon Strandenes ++dim; 866ec7ae49cSHåkon Strandenes } 867ec7ae49cSHåkon Strandenes #if defined(PETSC_USE_COMPLEX) 868ec7ae49cSHåkon Strandenes offset[dim] = 0; 869ec7ae49cSHåkon Strandenes count[dim] = 2; 870ec7ae49cSHåkon Strandenes ++dim; 871ec7ae49cSHåkon Strandenes #endif 872ec7ae49cSHåkon Strandenes 873ec7ae49cSHåkon Strandenes /* Create the memory and filespace */ 874729ab6d9SBarry Smith PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL)); 875729ab6d9SBarry Smith PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL)); 87647c6ae99SBarry Smith 87747c6ae99SBarry Smith /* Create property list for collective dataset write */ 878729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_DATASET_XFER)); 87947c6ae99SBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 880729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_dxpl_mpio,(plist_id, H5FD_MPIO_COLLECTIVE)); 88147c6ae99SBarry Smith #endif 882ec7ae49cSHåkon Strandenes /* To read dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 88347c6ae99SBarry Smith 88447c6ae99SBarry Smith ierr = VecGetArray(xin, &x);CHKERRQ(ierr); 885a1dbdba5SBarry Smith PetscStackCallHDF5(H5Dread,(dset_id, scalartype, memspace, filespace, plist_id, x)); 88647c6ae99SBarry Smith ierr = VecRestoreArray(xin, &x);CHKERRQ(ierr); 88747c6ae99SBarry Smith 88847c6ae99SBarry Smith /* Close/release resources */ 889bfd264e7SBarry Smith if (group != file_id) { 890729ab6d9SBarry Smith PetscStackCallHDF5(H5Gclose,(group)); 891bfd264e7SBarry Smith } 892729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 893729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(filespace)); 894729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(memspace)); 895729ab6d9SBarry Smith PetscStackCallHDF5(H5Dclose,(dset_id)); 89647c6ae99SBarry Smith PetscFunctionReturn(0); 89747c6ae99SBarry Smith } 89847c6ae99SBarry Smith #endif 89947c6ae99SBarry Smith 90047c6ae99SBarry Smith PetscErrorCode VecLoad_Binary_DA(Vec xin, PetscViewer viewer) 90147c6ae99SBarry Smith { 9029a42bb27SBarry Smith DM da; 90347c6ae99SBarry Smith PetscErrorCode ierr; 90447c6ae99SBarry Smith Vec natural; 90547c6ae99SBarry Smith const char *prefix; 90647c6ae99SBarry Smith PetscInt bs; 90747c6ae99SBarry Smith PetscBool flag; 90847c6ae99SBarry Smith DM_DA *dd; 90947c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 91047c6ae99SBarry Smith PetscBool isMPIIO; 91147c6ae99SBarry Smith #endif 91247c6ae99SBarry Smith 91347c6ae99SBarry Smith PetscFunctionBegin; 914c688c046SMatthew G Knepley ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 91547c6ae99SBarry Smith dd = (DM_DA*)da->data; 91647c6ae99SBarry Smith #if defined(PETSC_HAVE_MPIIO) 917bc196f7cSDave May ierr = PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO);CHKERRQ(ierr); 91847c6ae99SBarry Smith if (isMPIIO) { 919aa219208SBarry Smith ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_FALSE);CHKERRQ(ierr); 92047c6ae99SBarry Smith PetscFunctionReturn(0); 92147c6ae99SBarry Smith } 92247c6ae99SBarry Smith #endif 92347c6ae99SBarry Smith 92447c6ae99SBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr); 925aa219208SBarry Smith ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr); 92647c6ae99SBarry Smith ierr = PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name);CHKERRQ(ierr); 92747c6ae99SBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr); 9288aea9937SBarry Smith ierr = VecLoad(natural,viewer);CHKERRQ(ierr); 929aa219208SBarry Smith ierr = DMDANaturalToGlobalBegin(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr); 930aa219208SBarry Smith ierr = DMDANaturalToGlobalEnd(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr); 931fcfd50ebSBarry Smith ierr = VecDestroy(&natural);CHKERRQ(ierr); 932aa219208SBarry Smith ierr = PetscInfo(xin,"Loading vector from natural ordering into DMDA\n");CHKERRQ(ierr); 933c5929fdfSBarry Smith ierr = PetscOptionsGetInt(NULL,((PetscObject)xin)->prefix,"-vecload_block_size",&bs,&flag);CHKERRQ(ierr); 93447c6ae99SBarry Smith if (flag && bs != dd->w) { 935aa219208SBarry Smith ierr = PetscInfo2(xin,"Block size in file %D not equal to DMDA's dof %D\n",bs,dd->w);CHKERRQ(ierr); 93647c6ae99SBarry Smith } 93747c6ae99SBarry Smith PetscFunctionReturn(0); 93847c6ae99SBarry Smith } 93947c6ae99SBarry Smith 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