Actual source code: dtri.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

  2: /*
  3:        Provides the calling sequences for all the basic PetscDraw routines.
  4: */
  5:  #include <petsc/private/drawimpl.h>

  7: /*@
  8:    PetscDrawTriangle - PetscDraws a triangle  onto a drawable.

 10:    Not Collective

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

 17:    Level: beginner

 19:    Concepts: drawing^triangle
 20:    Concepts: graphics^triangle
 21:    Concepts: triangle

 23: .seealso: PetscDrawLine(), PetscDrawRectangle(), PetscDrawEllipse(), PetscDrawMarker(), PetscDrawPoint(), PetscDrawArrow()
 24: @*/
 25: PetscErrorCode  PetscDrawTriangle(PetscDraw draw,PetscReal x1,PetscReal y_1,PetscReal x2,PetscReal y2,PetscReal x3,PetscReal y3,int c1,int c2,int c3)
 26: {

 31:   if (!draw->ops->triangle) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"This draw type %s does not support drawing triangles",((PetscObject)draw)->type_name);
 32:   (*draw->ops->triangle)(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);
 33:   return(0);
 34: }

 36: /*@
 37:    PetscDrawScalePopup - draws a contour scale window.

 39:    Collective on PetscDraw

 41:    Input Parameters:
 42: +  popup - the window (often a window obtained via PetscDrawGetPopup()
 43: .  min - minimum value being plotted
 44: -  max - maximum value being plotted

 46:    Level: intermediate

 48:    Notes:
 49:     All processors that share the draw MUST call this routine

 51: .seealso: PetscDrawGetPopup(), PetscDrawTensorContour()

 53: @*/
 54: PetscErrorCode  PetscDrawScalePopup(PetscDraw popup,PetscReal min,PetscReal max)
 55: {
 56:   PetscBool      isnull;
 57:   PetscReal      xl = 0.0,yl = 0.0,xr = 1.0,yr = 1.0;
 58:   PetscMPIInt    rank;
 60:   int            i;
 61:   char           string[32];

 64:   if (!popup) return(0);
 66:   PetscDrawIsNull(popup,&isnull);
 67:   if (isnull) return(0);
 68:   MPI_Comm_rank(PetscObjectComm((PetscObject)popup),&rank);

 70:   PetscDrawCheckResizedWindow(popup);
 71:   PetscDrawClear(popup);
 72:   PetscDrawSetTitle(popup,"Contour Scale");
 73:   PetscDrawSetCoordinates(popup,xl,yl,xr,yr);
 74:   PetscDrawCollectiveBegin(popup);
 75:   if (!rank) {
 76:     for (i=0; i<10; i++) {
 77:       int c = PetscDrawRealToColor((PetscReal)i/9,0,1);
 78:       PetscDrawRectangle(popup,xl,yl,xr,yr,c,c,c,c);
 79:       yl += 0.1;
 80:     }
 81:     for (i=0; i<10; i++) {
 82:       PetscReal value = min + i*(max-min)/9;
 83:       /* look for a value that should be zero, but is not due to round-off */
 84:       if (PetscAbsReal(value) < 1.e-10 && max-min > 1.e-6) value = 0.0;
 85:       PetscSNPrintf(string,sizeof(string),"%18.16e",(double)value);
 86:       PetscDrawString(popup,0.2,0.02+i/10.0,PETSC_DRAW_BLACK,string);
 87:     }
 88:   }
 89:   PetscDrawCollectiveEnd(popup);
 90:   PetscDrawFlush(popup);
 91:   PetscDrawSave(popup);
 92:   return(0);
 93: }

 95: typedef struct {
 96:   int       m,n;
 97:   PetscReal *x,*y,min,max,*v;
 98:   PetscBool showgrid;
 99: } ZoomCtx;

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

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

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

124:    Collective on PetscDraw, but PetscDraw must be sequential

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

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

137:    Level: intermediate

139:    Concepts: contour plot
140:    Concepts: drawing^contour plot

142: .seealso: PetscDrawTensorContourPatch(), PetscDrawScalePopup()

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

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

164:   ctx.v   = v;
165:   ctx.m   = m;
166:   ctx.n   = n;
167:   ctx.max = ctx.min = v[0];
168:   for (i=0; i<N; i++) {
169:     if (ctx.max < ctx.v[i]) ctx.max = ctx.v[i];
170:     if (ctx.min > ctx.v[i]) ctx.min = ctx.v[i];
171:   }
172:   if (ctx.max - ctx.min < 1.e-7) {ctx.min -= 5.e-8; ctx.max += 5.e-8;}

174:   /* PetscDraw the scale window */
175:   PetscDrawGetPopup(draw,&popup);
176:   PetscDrawScalePopup(popup,ctx.min,ctx.max);

178:   ctx.showgrid = PETSC_FALSE;
179:   PetscOptionsGetBool(((PetscObject)draw)->options,NULL,"-draw_contour_grid",&ctx.showgrid,NULL);

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

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

198:   PetscDrawZoom(draw,PetscDrawTensorContour_Zoom,&ctx);

200:   if (!xin) {PetscFree(ctx.x);}
201:   if (!yin) {PetscFree(ctx.y);}
202:   return(0);
203: }

205: /*@
206:    PetscDrawTensorContourPatch - PetscDraws a rectangular patch of a contour plot
207:    for a two-dimensional array.

209:    Not Collective

211:    Input Parameters:
212: +  draw - the draw context
213: .  m,n - the number of local mesh points in the x and y direction
214: .  x,y - the locations of the local mesh points
215: .  min,max - the minimum and maximum value in the entire contour
216: -  v - the data

218:    Options Database Keys:
219: .  -draw_x_shared_colormap - Activates private colormap

221:    Level: advanced

223:    Note:
224:    This is a lower level support routine, usually the user will call
225:    PetscDrawTensorContour().

227:    Concepts: contour plot

229: .seealso: PetscDrawTensorContour()

231: @*/
232: PetscErrorCode  PetscDrawTensorContourPatch(PetscDraw draw,int m,int n,PetscReal *x,PetscReal *y,PetscReal min,PetscReal max,PetscReal *v)
233: {
235:   int            c1,c2,c3,c4,i,j;
236:   PetscReal      x1,x2,x3,x4,y_1,y2,y3,y4;

240:   /* PetscDraw the contour plot patch */
241:   for (j=0; j<n-1; j++) {
242:     for (i=0; i<m-1; i++) {
243:       x1 = x[i];   y_1 = y[j];  c1 = PetscDrawRealToColor(v[i+j*m],min,max);
244:       x2 = x[i+1]; y2 = y_1;    c2 = PetscDrawRealToColor(v[i+j*m+1],min,max);
245:       x3 = x2;     y3 = y[j+1]; c3 = PetscDrawRealToColor(v[i+j*m+1+m],min,max);
246:       x4 = x1;     y4 = y3;     c4 = PetscDrawRealToColor(v[i+j*m+m],min,max);

248:       PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);
249:       PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4);
250:     }
251:   }
252:   return(0);
253: }