#define PETSC_DLL /* Provides the calling sequences for all the basic PetscDraw routines. */ #include "src/sys/draw/drawimpl.h" /*I "petscdraw.h" I*/ #undef __FUNCT__ #define __FUNCT__ "PetscDrawTriangle" /*@ PetscDrawTriangle - PetscDraws a triangle onto a drawable. Not Collective Input Parameters: + draw - the drawing context . x1,y1,x2,y2,x3,y3 - the coordinates of the vertices - c1,c2,c3 - the colors of the three vertices in the same order as the xi,yi Level: beginner Concepts: drawing^triangle Concepts: graphics^triangle Concepts: triangle @*/ PetscErrorCode PETSC_DLLEXPORT PetscDrawTriangle(PetscDraw draw,PetscReal x1,PetscReal y_1,PetscReal x2,PetscReal y2,PetscReal x3,PetscReal y3, int c1,int c2,int c3) { PetscErrorCode ierr; PetscTruth isnull; PetscFunctionBegin; PetscValidHeaderSpecific(draw,PETSC_DRAW_COOKIE,1); ierr = PetscTypeCompare((PetscObject)draw,PETSC_DRAW_NULL,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); ierr = (*draw->ops->triangle)(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "PetscDrawScalePopup" /*@ PetscDrawScalePopup - PetscDraws a contour scale window. Collective on PetscDraw Input Parameters: + popup - the window (often a window obtained via PetscDrawGetPopup() . min - minimum value being plotted - max - maximum value being plotted Level: intermediate Notes: All processors that share the draw MUST call this routine @*/ PetscErrorCode PETSC_DLLEXPORT PetscDrawScalePopup(PetscDraw popup,PetscReal min,PetscReal max) { PetscReal xl = 0.0,yl = 0.0,xr = 1.0,yr = 1.0,value; PetscErrorCode ierr; int i,c = PETSC_DRAW_BASIC_COLORS,rank; char string[32]; MPI_Comm comm; PetscFunctionBegin; ierr = PetscDrawCheckResizedWindow(popup);CHKERRQ(ierr); ierr = PetscObjectGetComm((PetscObject)popup,&comm);CHKERRQ(ierr); ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); if (rank) PetscFunctionReturn(0); for (i=0; i<10; i++) { ierr = PetscDrawRectangle(popup,xl,yl,xr,yr,c,c,c,c);CHKERRQ(ierr); yl += .1; yr += .1; c = (int)((double)c + (245. - PETSC_DRAW_BASIC_COLORS)/9.); } for (i=0; i<10; i++) { value = min + i*(max-min)/9.0; /* look for a value that should be zero, but is not due to round-off */ if (PetscAbsReal(value) < 1.e-10 && max-min > 1.e-6) value = 0.0; sprintf(string,"%g",(double)value); ierr = PetscDrawString(popup,.2,.02 + i/10.0,PETSC_DRAW_BLACK,string);CHKERRQ(ierr); } ierr = PetscDrawSetTitle(popup,"Contour Scale");CHKERRQ(ierr); ierr = PetscDrawFlush(popup);CHKERRQ(ierr); PetscFunctionReturn(0); } typedef struct { int m,n; PetscReal *x,*y,min,max,*v; PetscTruth showgrid; } ZoomCtx; #undef __FUNCT__ #define __FUNCT__ "PetscDrawTensorContour_Zoom" static PetscErrorCode PetscDrawTensorContour_Zoom(PetscDraw win,void *dctx) { PetscErrorCode ierr; int i; ZoomCtx *ctx = (ZoomCtx*)dctx; PetscFunctionBegin; ierr = PetscDrawTensorContourPatch(win,ctx->m,ctx->n,ctx->x,ctx->y,ctx->max,ctx->min,ctx->v);CHKERRQ(ierr); if (ctx->showgrid) { for (i=0; im; i++) { ierr = PetscDrawLine(win,ctx->x[i],ctx->y[0],ctx->x[i],ctx->y[ctx->n-1],PETSC_DRAW_BLACK);CHKERRQ(ierr); } for (i=0; in; i++) { ierr = PetscDrawLine(win,ctx->x[0],ctx->y[i],ctx->x[ctx->m-1],ctx->y[i],PETSC_DRAW_BLACK);CHKERRQ(ierr); } } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "PetscDrawTensorContour" /*@C PetscDrawTensorContour - PetscDraws a contour plot for a two-dimensional array that is stored as a PETSc vector. Collective on PetscDraw, but PetscDraw must be sequential Input Parameters: + win - the window to draw in . m,n - the global number of mesh points in the x and y directions . xi,yi - the locations of the global mesh points (optional, use PETSC_NULL to indicate uniform spacing on [0,1]) - V - the values Options Database Keys: + -draw_x_shared_colormap - Indicates use of private colormap - -draw_contour_grid - PetscDraws grid contour Level: intermediate Concepts: contour plot Concepts: drawing^contour plot .seealso: PetscDrawTensorContourPatch() @*/ PetscErrorCode PETSC_DLLEXPORT PetscDrawTensorContour(PetscDraw win,int m,int n,const PetscReal xi[],const PetscReal yi[],PetscReal *v) { PetscErrorCode ierr; int N = m*n; PetscTruth isnull; PetscDraw popup; MPI_Comm comm; int xin=1,yin=1,i; PetscMPIInt size; PetscReal h; ZoomCtx ctx; PetscFunctionBegin; ierr = PetscDrawIsNull(win,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); ierr = PetscObjectGetComm((PetscObject)win,&comm);CHKERRQ(ierr); ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"May only be used with single processor PetscDraw"); if (N <= 0) { SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"n %d and m %d must be positive",m,n); } /* create scale window */ ierr = PetscDrawGetPopup(win,&popup);CHKERRQ(ierr); ierr = PetscDrawCheckResizedWindow(win);CHKERRQ(ierr); ctx.v = v; ctx.m = m; ctx.n = n; ctx.max = ctx.min = v[0]; for (i=0; i ctx.v[i]) ctx.min = ctx.v[i]; } if (ctx.max - ctx.min < 1.e-7) {ctx.min -= 5.e-8; ctx.max += 5.e-8;} /* PetscDraw the scale window */ if (popup) {ierr = PetscDrawScalePopup(popup,ctx.min,ctx.max);CHKERRQ(ierr);} ierr = PetscOptionsHasName(PETSC_NULL,"-draw_contour_grid",&ctx.showgrid);CHKERRQ(ierr); /* fill up x and y coordinates */ if (!xi) { xin = 0; ierr = PetscMalloc(ctx.m*sizeof(PetscReal),&ctx.x);CHKERRQ(ierr); h = 1.0/(ctx.m-1); ctx.x[0] = 0.0; for (i=1; i