Actual source code: dtri.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  2: /*
  3:        Provides the calling sequences for all the basic PetscDraw routines.
  4: */
  5: #include <petsc/private/drawimpl.h>  /*I "petscdraw.h" I*/

  9: /*@
 10:    PetscDrawTriangle - PetscDraws a triangle  onto a drawable.

 12:    Not Collective

 14:    Input Parameters:
 15: +  draw - the drawing context
 16: .  x1,y1,x2,y2,x3,y3 - the coordinates of the vertices
 17: -  c1,c2,c3 - the colors of the three vertices in the same order as the xi,yi

 19:    Level: beginner

 21:    Concepts: drawing^triangle
 22:    Concepts: graphics^triangle
 23:    Concepts: triangle

 25: @*/
 26: PetscErrorCode  PetscDrawTriangle(PetscDraw draw,PetscReal x1,PetscReal y_1,PetscReal x2,PetscReal y2,PetscReal x3,PetscReal y3,int c1,int c2,int c3)
 27: {
 29:   PetscBool      isnull;

 33:   PetscObjectTypeCompare((PetscObject)draw,PETSC_DRAW_NULL,&isnull);
 34:   if (isnull) return(0);
 35:   if (!draw->ops->triangle) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No triangles with this PetscDrawType");
 36:   (*draw->ops->triangle)(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);
 37:   return(0);
 38: }


 43: /*@
 44:        PetscDrawScalePopup - PetscDraws a contour scale window.

 46:      Collective on PetscDraw

 48:   Input Parameters:
 49: +    popup - the window (often a window obtained via PetscDrawGetPopup()
 50: .    min - minimum value being plotted
 51: -    max - maximum value being plotted

 53:   Level: intermediate

 55:   Notes:
 56:      All processors that share the draw MUST call this routine

 58: @*/
 59: PetscErrorCode  PetscDrawScalePopup(PetscDraw popup,PetscReal min,PetscReal max)
 60: {
 61:   PetscReal      xl = 0.0,yl = 0.0,xr = 2.0,yr = 1.0,value;
 63:   int            i,c = PETSC_DRAW_BASIC_COLORS,rank;
 64:   char           string[32];
 65:   MPI_Comm       comm;

 68:   PetscDrawCheckResizedWindow(popup);
 69:   PetscDrawSynchronizedClear(popup);
 70:   PetscObjectGetComm((PetscObject)popup,&comm);
 71:   MPI_Comm_rank(comm,&rank);
 72:   if (rank) return(0);

 74:   for (i=0; i<10; i++) {
 75:     PetscDrawRectangle(popup,xl,yl,xr,yr,c,c,c,c);

 77:     yl += .1; yr += .1;
 78:     c = (int)((double)c + (245. - PETSC_DRAW_BASIC_COLORS)/9.);
 79:   }
 80:   for (i=0; i<10; i++) {
 81:     value = min + i*(max-min)/9.0;
 82:     /* look for a value that should be zero, but is not due to round-off */
 83:     if (PetscAbsReal(value) < 1.e-10 && max-min > 1.e-6) value = 0.0;
 84:     sprintf(string,"%18.16e",(double)value);
 85:     PetscDrawString(popup,.2,.02 + i/10.0,PETSC_DRAW_BLACK,string);
 86:   }
 87:   PetscDrawSetTitle(popup,"Contour Scale");
 88:   PetscDrawFlush(popup);
 89:   return(0);
 90: }

 92: typedef struct {
 93:   int       m,n;
 94:   PetscReal *x,*y,min,max,*v;
 95:   PetscBool showgrid;
 96: } ZoomCtx;

100: static PetscErrorCode PetscDrawTensorContour_Zoom(PetscDraw win,void *dctx)
101: {
103:   int            i;
104:   ZoomCtx        *ctx = (ZoomCtx*)dctx;

107:   PetscDrawTensorContourPatch(win,ctx->m,ctx->n,ctx->x,ctx->y,ctx->max,ctx->min,ctx->v);
108:   if (ctx->showgrid) {
109:     for (i=0; i<ctx->m; i++) {
110:       PetscDrawLine(win,ctx->x[i],ctx->y[0],ctx->x[i],ctx->y[ctx->n-1],PETSC_DRAW_BLACK);
111:     }
112:     for (i=0; i<ctx->n; i++) {
113:       PetscDrawLine(win,ctx->x[0],ctx->y[i],ctx->x[ctx->m-1],ctx->y[i],PETSC_DRAW_BLACK);
114:     }
115:   }
116:   return(0);
117: }

121: /*@C
122:    PetscDrawTensorContour - PetscDraws a contour plot for a two-dimensional array
123:    that is stored as a PETSc vector.

125:    Collective on PetscDraw, but PetscDraw must be sequential

127:    Input Parameters:
128: +   win   - the window to draw in
129: .   m,n   - the global number of mesh points in the x and y directions
130: .   xi,yi - the locations of the global mesh points (optional, use NULL
131:             to indicate uniform spacing on [0,1])
132: -   V     - the values

134:    Options Database Keys:
135: +  -draw_x_shared_colormap - Indicates use of private colormap
136: -  -draw_contour_grid - PetscDraws grid contour

138:    Level: intermediate

140:    Concepts: contour plot
141:    Concepts: drawing^contour plot

143: .seealso: PetscDrawTensorContourPatch()

145: @*/
146: PetscErrorCode  PetscDrawTensorContour(PetscDraw win,int m,int n,const PetscReal xi[],const PetscReal yi[],PetscReal *v)
147: {
149:   int            N = m*n;
150:   PetscBool      isnull;
151:   PetscDraw      popup;
152:   MPI_Comm       comm;
153:   int            xin=1,yin=1,i;
154:   PetscMPIInt    size;
155:   PetscReal      h;
156:   ZoomCtx        ctx;

159:   PetscDrawIsNull(win,&isnull); if (isnull) return(0);
160:   PetscObjectGetComm((PetscObject)win,&comm);
161:   MPI_Comm_size(comm,&size);
162:   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"May only be used with single processor PetscDraw");
163:   if (N <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"n %d and m %d must be positive",m,n);

165:   /* create scale window */
166:   PetscDrawGetPopup(win,&popup);
167:   PetscDrawCheckResizedWindow(win);

169:   ctx.v   = v;
170:   ctx.m   = m;
171:   ctx.n   = n;
172:   ctx.max = ctx.min = v[0];
173:   for (i=0; i<N; i++) {
174:     if (ctx.max < ctx.v[i]) ctx.max = ctx.v[i];
175:     if (ctx.min > ctx.v[i]) ctx.min = ctx.v[i];
176:   }
177:   if (ctx.max - ctx.min < 1.e-7) {ctx.min -= 5.e-8; ctx.max += 5.e-8;}

179:   /* PetscDraw the scale window */
180:   if (popup) {PetscDrawScalePopup(popup,ctx.min,ctx.max);}

182:   ctx.showgrid = PETSC_FALSE;
183:   PetscOptionsGetBool(NULL,"-draw_contour_grid",&ctx.showgrid,NULL);

185:   /* fill up x and y coordinates */
186:   if (!xi) {
187:     xin      = 0;
188:     PetscMalloc1(ctx.m,&ctx.x);
189:     h        = 1.0/(ctx.m-1);
190:     ctx.x[0] = 0.0;
191:     for (i=1; i<ctx.m; i++) ctx.x[i] = ctx.x[i-1] + h;
192:   } else ctx.x = (PetscReal*)xi;

194:   if (!yi) {
195:     yin      = 0;
196:     PetscMalloc1(ctx.n,&ctx.y);
197:     h        = 1.0/(ctx.n-1);
198:     ctx.y[0] = 0.0;
199:     for (i=1; i<ctx.n; i++) ctx.y[i] = ctx.y[i-1] + h;
200:   } else ctx.y = (PetscReal*)yi;

202:   PetscDrawZoom(win,(PetscErrorCode (*)(PetscDraw,void*))PetscDrawTensorContour_Zoom,&ctx);

204:   if (!xin) {PetscFree(ctx.x);}
205:   if (!yin) {PetscFree(ctx.y);}
206:   return(0);
207: }

211: /*@
212:    PetscDrawTensorContourPatch - PetscDraws a rectangular patch of a contour plot
213:    for a two-dimensional array.

215:    Not Collective

217:    Input Parameters:
218: +  win - the window to draw in
219: .  m,n - the number of local mesh points in the x and y direction
220: .  x,y - the locations of the local mesh points
221: .  max,min - the maximum and minimum value in the entire contour
222: -  v - the data

224:    Options Database Keys:
225: .  -draw_x_shared_colormap - Activates private colormap

227:    Level: advanced

229:    Note:
230:    This is a lower level support routine, usually the user will call
231:    PetscDrawTensorContour().

233:    Concepts: contour plot

235: .seealso: PetscDrawTensorContour()

237: @*/
238: PetscErrorCode  PetscDrawTensorContourPatch(PetscDraw draw,int m,int n,PetscReal *x,PetscReal *y,PetscReal max,PetscReal min,PetscReal *v)
239: {
241:   int            c1,c2,c3,c4,i,j;
242:   PetscReal      x1,x2,x3,x4,y_1,y2,y3,y4,scale;

245:   scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/(max - min);

247:   /* PetscDraw the contour plot patch */
248:   for (j=0; j<n-1; j++) {
249:     for (i=0; i<m-1; i++) {
250:       x1 = x[i];  y_1 = y[j]; c1 = (int)(PETSC_DRAW_BASIC_COLORS + scale*(v[i+j*m] - min));
251:       x2 = x[i+1];y2 = y_1;   c2 = (int)(PETSC_DRAW_BASIC_COLORS + scale*(v[i+j*m+1]-min));
252:       x3 = x2;    y3 = y[j+1];c3 = (int)(PETSC_DRAW_BASIC_COLORS + scale*(v[i+j*m+1+m]-min));
253:       x4 = x1;    y4 = y3;    c4 = (int)(PETSC_DRAW_BASIC_COLORS + scale*(v[i+j*m+m]-min));

255:       PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);
256:       PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4);
257:     }
258:   }
259:   return(0);
260: }