Actual source code: dtri.c
petsc-3.7.7 2017-09-25
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: {
32: if (!draw->ops->triangle) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"This draw type %s does not support drawing triangles",((PetscObject)draw)->type_name);
33: (*draw->ops->triangle)(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);
34: return(0);
35: }
39: /*@
40: PetscDrawScalePopup - PetscDraws a contour scale window.
42: Collective on PetscDraw
44: Input Parameters:
45: + popup - the window (often a window obtained via PetscDrawGetPopup()
46: . min - minimum value being plotted
47: - max - maximum value being plotted
49: Level: intermediate
51: Notes: All processors that share the draw MUST call this routine
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;
103: static PetscErrorCode PetscDrawTensorContour_Zoom(PetscDraw win,void *dctx)
104: {
106: int i;
107: ZoomCtx *ctx = (ZoomCtx*)dctx;
110: PetscDrawTensorContourPatch(win,ctx->m,ctx->n,ctx->x,ctx->y,ctx->min,ctx->max,ctx->v);
111: if (ctx->showgrid) {
112: for (i=0; i<ctx->m; i++) {
113: PetscDrawLine(win,ctx->x[i],ctx->y[0],ctx->x[i],ctx->y[ctx->n-1],PETSC_DRAW_BLACK);
114: }
115: for (i=0; i<ctx->n; i++) {
116: PetscDrawLine(win,ctx->x[0],ctx->y[i],ctx->x[ctx->m-1],ctx->y[i],PETSC_DRAW_BLACK);
117: }
118: }
119: return(0);
120: }
124: /*@C
125: PetscDrawTensorContour - PetscDraws a contour plot for a two-dimensional array
126: that is stored as a PETSc vector.
128: Collective on PetscDraw, but PetscDraw must be sequential
130: Input Parameters:
131: + draw - the draw context
132: . m,n - the global number of mesh points in the x and y directions
133: . xi,yi - the locations of the global mesh points (optional, use NULL
134: to indicate uniform spacing on [0,1])
135: - V - the values
137: Options Database Keys:
138: + -draw_x_shared_colormap - Indicates use of private colormap
139: - -draw_contour_grid - PetscDraws grid contour
141: Level: intermediate
143: Concepts: contour plot
144: Concepts: drawing^contour plot
146: .seealso: PetscDrawTensorContourPatch()
148: @*/
149: PetscErrorCode PetscDrawTensorContour(PetscDraw draw,int m,int n,const PetscReal xi[],const PetscReal yi[],PetscReal *v)
150: {
152: int N = m*n;
153: PetscBool isnull;
154: PetscDraw popup;
155: int xin=1,yin=1,i;
156: PetscMPIInt size;
157: PetscReal h;
158: ZoomCtx ctx;
162: PetscDrawIsNull(draw,&isnull);
163: if (isnull) return(0);
164: MPI_Comm_size(PetscObjectComm((PetscObject)draw),&size);
165: if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"May only be used with single processor PetscDraw");
166: if (N <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"n %d and m %d must be positive",m,n);
168: ctx.v = v;
169: ctx.m = m;
170: ctx.n = n;
171: ctx.max = ctx.min = v[0];
172: for (i=0; i<N; i++) {
173: if (ctx.max < ctx.v[i]) ctx.max = ctx.v[i];
174: if (ctx.min > ctx.v[i]) ctx.min = ctx.v[i];
175: }
176: if (ctx.max - ctx.min < 1.e-7) {ctx.min -= 5.e-8; ctx.max += 5.e-8;}
178: /* PetscDraw the scale window */
179: PetscDrawGetPopup(draw,&popup);
180: PetscDrawScalePopup(popup,ctx.min,ctx.max);
182: ctx.showgrid = PETSC_FALSE;
183: PetscOptionsGetBool(((PetscObject)draw)->options,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(draw,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: + draw - the draw context
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: . min,max - the minimum and maximum 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 min,PetscReal max,PetscReal *v)
239: {
241: int c1,c2,c3,c4,i,j;
242: PetscReal x1,x2,x3,x4,y_1,y2,y3,y4;
246: /* PetscDraw the contour plot patch */
247: for (j=0; j<n-1; j++) {
248: for (i=0; i<m-1; i++) {
249: x1 = x[i]; y_1 = y[j]; c1 = PetscDrawRealToColor(v[i+j*m],min,max);
250: x2 = x[i+1]; y2 = y_1; c2 = PetscDrawRealToColor(v[i+j*m+1],min,max);
251: x3 = x2; y3 = y[j+1]; c3 = PetscDrawRealToColor(v[i+j*m+1+m],min,max);
252: x4 = x1; y4 = y3; c4 = PetscDrawRealToColor(v[i+j*m+m],min,max);
254: PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);
255: PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4);
256: }
257: }
258: return(0);
259: }