/*@ PetscDrawLGSPDraw - Redraws a line graph. Collective on PetscDrawLG Input Parameter: . lg - the line graph context Level: intermediate .seealso: PetscDrawLGDraw(), PetscDrawSPDraw() Developer Notes: This code cheats and uses the fact that the LG and SP structs are the same @*/ PetscErrorCode PetscDrawLGSPDraw(PetscDrawLG lg,PetscDrawSP spin) { PetscDrawLG sp = (PetscDrawLG)spin; PetscReal xmin,xmax,ymin,ymax; PetscErrorCode ierr; PetscBool isnull; PetscMPIInt rank; PetscDraw draw; PetscFunctionBegin; PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1); PetscValidHeaderSpecific(sp,PETSC_DRAWSP_CLASSID,2); ierr = PetscDrawIsNull(lg->win,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)lg),&rank);CHKERRQ(ierr); draw = lg->win; ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr); ierr = PetscDrawClear(draw);CHKERRQ(ierr); xmin = PetscMin(lg->xmin,sp->xmin); ymin = PetscMin(lg->ymin,sp->ymin); xmax = PetscMax(lg->xmax,sp->xmax); ymax = PetscMax(lg->ymax,sp->ymax); ierr = PetscDrawAxisSetLimits(lg->axis,xmin,xmax,ymin,ymax);CHKERRQ(ierr); ierr = PetscDrawAxisDraw(lg->axis);CHKERRQ(ierr); ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); if (!rank) { int i,j,dim,nopts; dim = lg->dim; nopts = lg->nopts; for (i=0; i<dim; i++) { for (j=1; j<nopts; j++) { ierr = PetscDrawLine(draw,lg->x[(j-1)*dim+i],lg->y[(j-1)*dim+i],lg->x[j*dim+i],lg->y[j*dim+i],PETSC_DRAW_BLACK+i);CHKERRQ(ierr); if (lg->use_markers) { ierr = PetscDrawMarker(draw,lg->x[j*dim+i],lg->y[j*dim+i],PETSC_DRAW_RED);CHKERRQ(ierr); } } } dim = sp->dim; nopts = sp->nopts; for (i=0; i<dim; i++) { for (j=0; j<nopts; j++) { ierr = PetscDrawMarker(draw,sp->x[j*dim+i],sp->y[j*dim+i],PETSC_DRAW_RED);CHKERRQ(ierr); } } } ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); ierr = PetscDrawFlush(draw);CHKERRQ(ierr); ierr = PetscDrawPause(draw);CHKERRQ(ierr); PetscFunctionReturn(0); }
PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin,PetscViewer v) { DM da; PetscErrorCode ierr; PetscMPIInt rank,size,tag1,tag2; PetscInt i,n,N,step,istart,isize,j,nbounds; MPI_Status status; PetscReal coors[4],ymin,ymax,min,max,xmin = 0.0,xmax = 0.0,tmp = 0.0,xgtmp = 0.0; const PetscScalar *array,*xg; PetscDraw draw; PetscBool isnull,showmarkers = PETSC_FALSE; MPI_Comm comm; PetscDrawAxis axis; Vec xcoor; DMBoundaryType bx; const PetscReal *bounds; PetscInt *displayfields; PetscInt k,ndisplayfields; PetscBool hold; PetscFunctionBegin; ierr = PetscViewerDrawGetDraw(v,0,&draw);CHKERRQ(ierr); ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); ierr = PetscViewerDrawGetBounds(v,&nbounds,&bounds);CHKERRQ(ierr); ierr = VecGetDM(xin,&da);CHKERRQ(ierr); if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); ierr = PetscOptionsGetBool(NULL,NULL,"-draw_vec_use_markers",&showmarkers,NULL);CHKERRQ(ierr); ierr = DMDAGetInfo(da,0,&N,0,0,0,0,0,&step,0,&bx,0,0,0);CHKERRQ(ierr); ierr = DMDAGetCorners(da,&istart,0,0,&isize,0,0);CHKERRQ(ierr); ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr); ierr = VecGetLocalSize(xin,&n);CHKERRQ(ierr); n = n/step; /* get coordinates of nodes */ ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); if (!xcoor) { ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,0.0,0.0,0.0);CHKERRQ(ierr); ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); } ierr = VecGetArrayRead(xcoor,&xg);CHKERRQ(ierr); ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr); ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); /* Determine the min and max x coordinate in plot */ if (!rank) { xmin = PetscRealPart(xg[0]); } if (rank == size-1) { xmax = PetscRealPart(xg[n-1]); } ierr = MPI_Bcast(&xmin,1,MPIU_REAL,0,comm);CHKERRQ(ierr); ierr = MPI_Bcast(&xmax,1,MPIU_REAL,size-1,comm);CHKERRQ(ierr); ierr = DMDASelectFields(da,&ndisplayfields,&displayfields);CHKERRQ(ierr); #if defined(PETSC_HAVE_SETJMP_H) && defined(PETSC_HAVE_X) if (!setjmp(PetscXIOErrorJumpBuf)) XSetIOErrorHandler((XIOErrorHandler)PetscXIOHandler); else { XSetIOErrorHandler(NULL); ierr = PetscDrawSetType(draw,PETSC_DRAW_NULL);CHKERRQ(ierr); PetscFunctionReturn(0); } #endif for (k=0; k<ndisplayfields; k++) { j = displayfields[k]; ierr = PetscViewerDrawGetDraw(v,k,&draw);CHKERRQ(ierr); ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr); /* Determine the min and max y coordinate in plot */ min = 1.e20; max = -1.e20; for (i=0; i<n; i++) { if (PetscRealPart(array[j+i*step]) < min) min = PetscRealPart(array[j+i*step]); if (PetscRealPart(array[j+i*step]) > max) max = PetscRealPart(array[j+i*step]); } if (min + 1.e-10 > max) { min -= 1.e-5; max += 1.e-5; } if (j < nbounds) { min = PetscMin(min,bounds[2*j]); max = PetscMax(max,bounds[2*j+1]); } ierr = MPI_Reduce(&min,&ymin,1,MPIU_REAL,MPIU_MIN,0,comm);CHKERRQ(ierr); ierr = MPI_Reduce(&max,&ymax,1,MPIU_REAL,MPIU_MAX,0,comm);CHKERRQ(ierr); ierr = PetscViewerDrawGetHold(v,&hold);CHKERRQ(ierr); if (!hold) { ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr); } ierr = PetscViewerDrawGetDrawAxis(v,k,&axis);CHKERRQ(ierr); ierr = PetscLogObjectParent((PetscObject)draw,(PetscObject)axis);CHKERRQ(ierr); if (!rank) { const char *title; ierr = PetscDrawAxisSetLimits(axis,xmin,xmax,ymin,ymax);CHKERRQ(ierr); ierr = PetscDrawAxisDraw(axis);CHKERRQ(ierr); ierr = PetscDrawGetCoordinates(draw,coors,coors+1,coors+2,coors+3);CHKERRQ(ierr); ierr = DMDAGetFieldName(da,j,&title);CHKERRQ(ierr); if (title) {ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr);} } ierr = MPI_Bcast(coors,4,MPIU_REAL,0,comm);CHKERRQ(ierr); if (rank) { ierr = PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr); } /* draw local part of vector */ ierr = PetscObjectGetNewTag((PetscObject)xin,&tag1);CHKERRQ(ierr); ierr = PetscObjectGetNewTag((PetscObject)xin,&tag2);CHKERRQ(ierr); if (rank < size-1) { /*send value to right */ ierr = MPI_Send((void*)&array[j+(n-1)*step],1,MPIU_REAL,rank+1,tag1,comm);CHKERRQ(ierr); ierr = MPI_Send((void*)&xg[n-1],1,MPIU_REAL,rank+1,tag1,comm);CHKERRQ(ierr); } if (!rank && bx == DM_BOUNDARY_PERIODIC && size > 1) { /* first processor sends first value to last */ ierr = MPI_Send((void*)&array[j],1,MPIU_REAL,size-1,tag2,comm);CHKERRQ(ierr); } for (i=1; i<n; i++) { ierr = PetscDrawLine(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+step*(i-1)]),PetscRealPart(xg[i]),PetscRealPart(array[j+step*i]),PETSC_DRAW_RED);CHKERRQ(ierr); if (showmarkers) { ierr = PetscDrawMarker(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+step*(i-1)]),PETSC_DRAW_BLACK);CHKERRQ(ierr); } } if (rank) { /* receive value from left */ ierr = MPI_Recv(&tmp,1,MPIU_REAL,rank-1,tag1,comm,&status);CHKERRQ(ierr); ierr = MPI_Recv(&xgtmp,1,MPIU_REAL,rank-1,tag1,comm,&status);CHKERRQ(ierr); ierr = PetscDrawLine(draw,xgtmp,tmp,PetscRealPart(xg[0]),PetscRealPart(array[j]),PETSC_DRAW_RED);CHKERRQ(ierr); if (showmarkers) { ierr = PetscDrawPoint(draw,xgtmp,tmp,PETSC_DRAW_BLACK);CHKERRQ(ierr); } } if (rank == size-1 && bx == DM_BOUNDARY_PERIODIC && size > 1) { ierr = MPI_Recv(&tmp,1,MPIU_REAL,0,tag2,comm,&status);CHKERRQ(ierr); /* If the mesh is not uniform we do not know the mesh spacing between the last point on the right and the first ghost point */ ierr = PetscDrawLine(draw,PetscRealPart(xg[n-1]),PetscRealPart(array[j+step*(n-1)]),PetscRealPart(xg[n-1]+(xg[n-1]-xg[n-2])),tmp,PETSC_DRAW_RED);CHKERRQ(ierr); if (showmarkers) { ierr = PetscDrawMarker(draw,PetscRealPart(xg[n-2]),PetscRealPart(array[j+step*(n-1)]),PETSC_DRAW_BLACK);CHKERRQ(ierr); } } ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); ierr = PetscDrawPause(draw);CHKERRQ(ierr); } #if defined(PETSC_HAVE_SETJMP_H) && defined(PETSC_HAVE_X) XSetIOErrorHandler(NULL); #endif ierr = PetscFree(displayfields);CHKERRQ(ierr); ierr = VecRestoreArrayRead(xcoor,&xg);CHKERRQ(ierr); ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr); PetscFunctionReturn(0); }
/*@ PetscDrawLGDraw - Redraws a line graph. Collective on PetscDrawLG Input Parameter: . lg - the line graph context Level: intermediate .seealso: PetscDrawSPDraw(), PetscDrawLGSPDraw(), PetscDrawLGReset() @*/ PetscErrorCode PetscDrawLGDraw(PetscDrawLG lg) { PetscReal xmin,xmax,ymin,ymax; PetscErrorCode ierr; PetscMPIInt rank; PetscDraw draw; PetscBool isnull; PetscFunctionBegin; PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,1); ierr = PetscDrawIsNull(lg->win,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)lg),&rank);CHKERRQ(ierr); draw = lg->win; ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr); ierr = PetscDrawClear(draw);CHKERRQ(ierr); xmin = lg->xmin; xmax = lg->xmax; ymin = lg->ymin; ymax = lg->ymax; ierr = PetscDrawAxisSetLimits(lg->axis,xmin,xmax,ymin,ymax);CHKERRQ(ierr); ierr = PetscDrawAxisDraw(lg->axis);CHKERRQ(ierr); ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); if (!rank) { int i,j,dim=lg->dim,nopts=lg->nopts,cl; for (i=0; i<dim; i++) { for (j=1; j<nopts; j++) { cl = lg->colors ? lg->colors[i] : (PETSC_DRAW_BLACK + i); ierr = PetscDrawLine(draw,lg->x[(j-1)*dim+i],lg->y[(j-1)*dim+i],lg->x[j*dim+i],lg->y[j*dim+i],cl);CHKERRQ(ierr); if (lg->use_markers) {ierr = PetscDrawMarker(draw,lg->x[j*dim+i],lg->y[j*dim+i],cl);CHKERRQ(ierr);} } } } if (!rank && lg->legend) { int i,dim=lg->dim,cl; PetscReal xl,yl,xr,yr,tw,th; size_t slen,len=0; ierr = PetscDrawAxisGetLimits(lg->axis,&xl,&xr,&yl,&yr);CHKERRQ(ierr); ierr = PetscDrawStringGetSize(draw,&tw,&th);CHKERRQ(ierr); for (i=0; i<dim; i++) { ierr = PetscStrlen(lg->legend[i],&slen);CHKERRQ(ierr); len = PetscMax(len,slen); } xr = xr - 1.5*tw; xl = xr - (len + 7)*tw; yr = yr - 1.0*th; yl = yr - (dim + 1)*th; ierr = PetscDrawLine(draw,xl,yl,xr,yl,PETSC_DRAW_BLACK);CHKERRQ(ierr); ierr = PetscDrawLine(draw,xr,yl,xr,yr,PETSC_DRAW_BLACK);CHKERRQ(ierr); ierr = PetscDrawLine(draw,xr,yr,xl,yr,PETSC_DRAW_BLACK);CHKERRQ(ierr); ierr = PetscDrawLine(draw,xl,yr,xl,yl,PETSC_DRAW_BLACK);CHKERRQ(ierr); for (i=0; i<dim; i++) { cl = lg->colors ? lg->colors[i] : (PETSC_DRAW_BLACK + i); ierr = PetscDrawLine(draw,xl + 1*tw,yr - (i + 1)*th,xl + 5*tw,yr - (i + 1)*th,cl);CHKERRQ(ierr); ierr = PetscDrawString(draw,xl + 6*tw,yr - (i + 1.5)*th,PETSC_DRAW_BLACK,lg->legend[i]);CHKERRQ(ierr); } } ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); ierr = PetscDrawFlush(draw);CHKERRQ(ierr); ierr = PetscDrawPause(draw);CHKERRQ(ierr); PetscFunctionReturn(0); }