Actual source code: drawimage.c
petsc-3.7.3 2016-08-01
1: #include <../src/sys/classes/draw/impls/image/drawimage.h> /*I "petscdraw.h" I*/
2: #include <petsc/private/drawimpl.h> /*I "petscdraw.h" I*/
4: #if defined(PETSC_USE_DEBUG)
5: #define PetscDrawValidColor(color) \
6: do { if (PetscUnlikely((color)<0||(color)>=256)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Color value %D out of range [0..255]",(PetscInt)(color)); } while (0)
7: #else
8: #define PetscDrawValidColor(color) do {} while (0)
9: #endif
11: #define XTRANS(draw,img,x) ((int)(((img)->w-1)*((draw)->port_xl + ((((x) - (draw)->coor_xl)*((draw)->port_xr - (draw)->port_xl))/((draw)->coor_xr - (draw)->coor_xl)))))
12: #define YTRANS(draw,img,y) (((img)->h-1) - (int)(((img)->h-1)*((draw)->port_yl + ((((y) - (draw)->coor_yl)*((draw)->port_yr - (draw)->port_yl))/((draw)->coor_yr - (draw)->coor_yl)))))
14: #define ITRANS(draw,img,i) ((draw)->coor_xl + (((PetscReal)(i))*((draw)->coor_xr - (draw)->coor_xl)/((img)->w-1) - (draw)->port_xl)/((draw)->port_xr - (draw)->port_xl))
15: #define JTRANS(draw,img,j) ((draw)->coor_yl + (((PetscReal)(j))/((img)->h-1) + (draw)->port_yl - 1)*((draw)->coor_yr - (draw)->coor_yl)/((draw)->port_yl - (draw)->port_yr))
20: static PetscErrorCode PetscDrawSetViewport_Image(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr)
21: {
22: PetscImage img = (PetscImage)draw->data;
24: {
25: int xmax = img->w - 1, ymax = img->h - 1;
26: int xa = (int)(xl*xmax), ya = ymax - (int)(yr*ymax);
27: int xb = (int)(xr*xmax), yb = ymax - (int)(yl*ymax);
28: PetscImageSetClip(img,xa,ya,xb+1-xa,yb+1-ya);
29: }
30: return(0);
31: }
33: /*
36: static PetscErrorCode PetscDrawSetCoordinates_Image(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr)
37: {
39: return(0);
40: }*/
41: #define PetscDrawSetCoordinates_Image NULL
45: static PetscErrorCode PetscDrawCoordinateToPixel_Image(PetscDraw draw,PetscReal x,PetscReal y,int *i,int *j)
46: {
47: PetscImage img = (PetscImage)draw->data;
49: if (i) *i = XTRANS(draw,img,x);
50: if (j) *j = YTRANS(draw,img,y);
51: return(0);
52: }
56: static PetscErrorCode PetscDrawPixelToCoordinate_Image(PetscDraw draw,int i,int j,PetscReal *x,PetscReal *y)
57: {
58: PetscImage img = (PetscImage)draw->data;
60: if (x) *x = ITRANS(draw,img,i);
61: if (y) *y = JTRANS(draw,img,j);
62: return(0);
63: }
65: /*
68: static PetscErrorCode PetscDrawPointSetSize_Image(PetscDraw draw,PetscReal width)
69: {
71: return(0);
72: }*/
73: #define PetscDrawPointSetSize_Image NULL
77: static PetscErrorCode PetscDrawPoint_Image(PetscDraw draw,PetscReal x,PetscReal y,int c)
78: {
79: PetscImage img = (PetscImage)draw->data;
81: PetscDrawValidColor(c);
82: {
83: int j, xx = XTRANS(draw,img,x);
84: int i, yy = YTRANS(draw,img,y);
85: for (i=-1; i<=1; i++)
86: for (j=-1; j<=1; j++)
87: PetscImageDrawPixel(img,xx+j,yy+i,c);
88: }
89: return(0);
90: }
94: static PetscErrorCode PetscDrawPointPixel_Image(PetscDraw draw,int x,int y,int c)
95: {
96: PetscImage img = (PetscImage)draw->data;
98: PetscDrawValidColor(c);
99: {
100: PetscImageDrawPixel(img,x,y,c);
101: }
102: return(0);
103: }
105: /*
108: static PetscErrorCode PetscDrawLineSetWidth_Image(PetscDraw draw,PetscReal width)
109: {
111: return(0);
112: }*/
113: #define PetscDrawLineSetWidth_Image NULL
117: static PetscErrorCode PetscDrawLineGetWidth_Image(PetscDraw draw,PetscReal *width)
118: {
119: PetscImage img = (PetscImage)draw->data;
121: {
122: int lw = 1;
123: *width = lw*(draw->coor_xr - draw->coor_xl)/(img->w*(draw->port_xr - draw->port_xl));
124: }
125: return(0);
126: }
130: static PetscErrorCode PetscDrawLine_Image(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr,int c)
131: {
132: PetscImage img = (PetscImage)draw->data;
134: {
135: int x_1 = XTRANS(draw,img,xl), x_2 = XTRANS(draw,img,xr);
136: int y_1 = YTRANS(draw,img,yl), y_2 = YTRANS(draw,img,yr);
137: PetscImageDrawLine(img,x_1,y_1,x_2,y_2,c);
138: }
139: return(0);
140: }
144: static PetscErrorCode PetscDrawArrow_Image(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr,int c)
145: {
146: PetscImage img = (PetscImage)draw->data;
148: PetscDrawValidColor(c);
149: {
150: int x_1 = XTRANS(draw,img,xl), x_2 = XTRANS(draw,img,xr);
151: int y_1 = YTRANS(draw,img,yl), y_2 = YTRANS(draw,img,yr);
152: if (x_1 == x_2 && y_1 == y_2) return(0);
153: PetscImageDrawLine(img,x_1,y_1,x_2,y_2,c);
154: if (x_1 == x_2 && PetscAbs(y_1 - y_2) > 7) {
155: if (y_2 > y_1) {
156: PetscImageDrawLine(img,x_2,y_2,x_2-3,y_2-3,c);
157: PetscImageDrawLine(img,x_2,y_2,x_2+3,y_2-3,c);
158: } else {
159: PetscImageDrawLine(img,x_2,y_2,x_2-3,y_2+3,c);
160: PetscImageDrawLine(img,x_2,y_2,x_2+3,y_2+3,c);
161: }
162: }
163: if (y_1 == y_2 && PetscAbs(x_1 - x_2) > 7) {
164: if (x_2 > x_1) {
165: PetscImageDrawLine(img,x_2-3,y_2-3,x_2,y_2,c);
166: PetscImageDrawLine(img,x_2-3,y_2+3,x_2,y_2,c);
167: } else {
168: PetscImageDrawLine(img,x_2,y_2,x_2+3,y_2-3,c);
169: PetscImageDrawLine(img,x_2,y_2,x_2+3,y_2+3,c);
170: }
171: }
172: }
173: return(0);
174: }
178: static PetscErrorCode PetscDrawRectangle_Image(PetscDraw draw,PetscReal xl,PetscReal yl,PetscReal xr,PetscReal yr,int c1,int c2,int c3,int c4)
179: {
180: PetscImage img = (PetscImage)draw->data;
182: PetscDrawValidColor(c1);
183: PetscDrawValidColor(c2);
184: PetscDrawValidColor(c3);
185: PetscDrawValidColor(c4);
186: {
187: int x = XTRANS(draw,img,xl), w = XTRANS(draw,img,xr) + 1 - x;
188: int y = YTRANS(draw,img,yr), h = YTRANS(draw,img,yl) + 1 - y;
189: int c = (c1 + c2 + c3 + c4)/4;
190: PetscImageDrawRectangle(img,x,y,w,h,c);
191: }
192: return(0);
193: }
197: static PetscErrorCode PetscDrawEllipse_Image(PetscDraw draw,PetscReal x,PetscReal y,PetscReal a,PetscReal b,int c)
198: {
199: PetscImage img = (PetscImage)draw->data;
201: PetscDrawValidColor(c);
202: a = PetscAbsReal(a);
203: b = PetscAbsReal(b);
204: {
205: int xc = XTRANS(draw,img,x), w = XTRANS(draw,img,x + a/2) + 0 - xc;
206: int yc = YTRANS(draw,img,y), h = YTRANS(draw,img,y - b/2) + 0 - yc;
207: if (PetscAbsReal(a-b) <= 0) w = h = PetscMin(w,h); /* workaround truncation errors */
208: PetscImageDrawEllipse(img,xc,yc,w,h,c);
209: }
210: return(0);
211: }
215: static PetscErrorCode PetscDrawTriangle_Image(PetscDraw draw,PetscReal X_1,PetscReal Y_1,PetscReal X_2,PetscReal Y_2,PetscReal X_3,PetscReal Y_3,int c1,int c2,int c3)
216: {
217: PetscImage img = (PetscImage)draw->data;
219: PetscDrawValidColor(c1);
220: PetscDrawValidColor(c2);
221: PetscDrawValidColor(c3);
222: {
223: int x_1 = XTRANS(draw,img,X_1), x_2 = XTRANS(draw,img,X_2), x_3 = XTRANS(draw,img,X_3);
224: int y_1 = YTRANS(draw,img,Y_1), y_2 = YTRANS(draw,img,Y_2), y_3 = YTRANS(draw,img,Y_3);
225: PetscImageDrawTriangle(img,x_1,y_1,c1,x_2,y_2,c2,x_3,y_3,c3);
226: }
227: return(0);
228: }
230: /*
233: static PetscErrorCode PetscDrawStringSetSize_Image(PetscDraw draw,PetscReal w,PetscReal h)
234: {
236: return(0);
237: }*/
238: #define PetscDrawStringSetSize_Image NULL
242: static PetscErrorCode PetscDrawStringGetSize_Image(PetscDraw draw,PetscReal *w,PetscReal *h)
243: {
244: PetscImage img = (PetscImage)draw->data;
246: {
247: int tw = PetscImageFontWidth;
248: int th = PetscImageFontHeight;
249: if (w) *w = tw*(draw->coor_xr - draw->coor_xl)/(img->w*(draw->port_xr - draw->port_xl));
250: if (h) *h = th*(draw->coor_yr - draw->coor_yl)/(img->h*(draw->port_yr - draw->port_yl));
251: }
252: return(0);
253: }
257: static PetscErrorCode PetscDrawString_Image(PetscDraw draw,PetscReal x,PetscReal y,int c,const char text[])
258: {
259: PetscImage img = (PetscImage)draw->data;
260: PetscToken token;
261: char *subtext;
264: PetscDrawValidColor(c);
265: {
266: int xx = XTRANS(draw,img,x);
267: int yy = YTRANS(draw,img,y);
268: PetscTokenCreate(text,'\n',&token);
269: PetscTokenFind(token,&subtext);
270: while (subtext) {
271: PetscImageDrawText(img,xx,yy,c,subtext);
272: yy += PetscImageFontHeight;
273: PetscTokenFind(token,&subtext);
274: }
275: PetscTokenDestroy(&token);
276: }
277: return(0);
278: }
282: static PetscErrorCode PetscDrawStringVertical_Image(PetscDraw draw,PetscReal x,PetscReal y,int c,const char text[])
283: {
284: PetscImage img = (PetscImage)draw->data;
286: PetscDrawValidColor(c);
287: {
288: char chr[2] = {0, 0};
289: int xx = XTRANS(draw,img,x);
290: int yy = YTRANS(draw,img,y);
291: int offset = PetscImageFontHeight;
292: while ((chr[0] = *text++)) {
293: PetscImageDrawText(img,xx,yy+offset,c,chr);
294: yy += PetscImageFontHeight;
295: }
296: }
297: return(0);
298: }
300: /*
303: static PetscErrorCode PetscDrawStringBoxed_Image(PetscDraw draw,PetscReal sxl,PetscReal syl,int sc,int bc,const char text[],PetscReal *w,PetscReal *h)
304: {
306: if (w) *w = 0;
307: if (h) *h = 0;
308: return(0);
309: */
310: #define PetscDrawStringBoxed_Image NULL
312: /*
315: static PetscErrorCode PetscDrawFlush_Image(PetscDraw draw)
316: {
318: return(0);
319: }*/
320: #define PetscDrawFlush_Image NULL
324: static PetscErrorCode PetscDrawClear_Image(PetscDraw draw)
325: {
326: PetscImage img = (PetscImage)draw->data;
328: {
329: PetscImageClear(img);
330: }
331: return(0);
332: }
334: /*
337: static PetscErrorCode PetscDrawSetDoubleBuffer_Image(PetscDraw draw)
338: {
340: return(0);
341: }*/
342: #define PetscDrawSetDoubleBuffer_Image NULL
346: static PetscErrorCode PetscDrawGetPopup_Image(PetscDraw draw,PetscDraw *popup)
347: {
348: PetscBool flg = PETSC_FALSE;
352: PetscOptionsGetBool(((PetscObject)draw)->options,((PetscObject)draw)->prefix,"-draw_popup",&flg,NULL);
353: if (!flg) {*popup = NULL; return(0);}
354: PetscDrawCreate(PetscObjectComm((PetscObject)draw),NULL,NULL,0,0,220,220,popup);
355: PetscDrawSetType(*popup,PETSC_DRAW_IMAGE);
356: PetscObjectSetOptionsPrefix((PetscObject)*popup,"popup_");
357: PetscObjectAppendOptionsPrefix((PetscObject)*popup,((PetscObject)draw)->prefix);
358: draw->popup = *popup;
359: return(0);
360: }
362: /*
365: static PetscErrorCode PetscDrawSetTitle_Image(PetscDraw draw,const char title[])
366: {
368: return(0);
369: }*/
370: #define PetscDrawSetTitle_Image NULL
372: /*
375: static PetscErrorCode PetscDrawCheckResizedWindow_Image(PetscDraw draw)
376: {
378: return(0);
379: }*/
380: #define PetscDrawCheckResizedWindow_Image NULL
384: static PetscErrorCode PetscDrawResizeWindow_Image(PetscDraw draw,int w,int h)
385: {
386: PetscImage img = (PetscImage)draw->data;
390: if (w == img->w && h == img->h) return(0);
391: PetscFree(img->buffer);
393: img->w = w; img->h = h;
394: PetscCalloc1((size_t)(img->w*img->h),&img->buffer);
395: PetscDrawSetViewport_Image(draw,draw->port_xl,draw->port_yl,draw->port_xr,draw->port_yr);
396: return(0);
397: }
401: static PetscErrorCode PetscDrawDestroy_Image(PetscDraw draw)
402: {
403: PetscImage img = (PetscImage)draw->data;
407: PetscDrawDestroy(&draw->popup);
408: PetscFree(img->buffer);
409: PetscFree(draw->data);
410: return(0);
411: }
413: /*
416: static PetscErrorCode PetscDrawView_Image(PetscDraw draw,PetscViewer viewer)
417: {
419: return(0);
420: }*/
421: #define PetscDrawView_Image NULL
423: /*
426: static PetscErrorCode PetscDrawGetMouseButton_Image(PetscDraw draw,PetscDrawButton *button,PetscReal *x_user,PetscReal *y_user,PetscReal *x_phys,PetscReal *y_phys)
427: {
429: *button = PETSC_BUTTON_NONE;
430: if (x_user) *x_user = 0;
431: if (y_user) *y_user = 0;
432: if (x_phys) *x_phys = 0;
433: if (y_phys) *y_phys = 0;
434: return(0);
435: }*/
436: #define PetscDrawGetMouseButton_Image NULL
438: /*
441: static PetscErrorCode PetscDrawPause_Image(PetscDraw draw)
442: {
444: return(0);
445: }*/
446: #define PetscDrawPause_Image NULL
448: /*
451: static PetscErrorCode PetscDrawBeginPage_Image(PetscDraw draw)
452: {
454: return(0);
455: }*/
456: #define PetscDrawBeginPage_Image NULL
458: /*
461: static PetscErrorCode PetscDrawEndPage_Image(PetscDraw draw)
462: {
464: return(0);
465: }*/
466: #define PetscDrawEndPage_Image NULL
470: static PetscErrorCode PetscDrawGetSingleton_Image(PetscDraw draw,PetscDraw *sdraw)
471: {
472: PetscImage pimg = (PetscImage)draw->data;
473: PetscImage simg;
477: PetscDrawCreate(PETSC_COMM_SELF,NULL,NULL,0,0,draw->w,draw->h,sdraw);
478: PetscDrawSetType(*sdraw,PETSC_DRAW_IMAGE);
479: (*sdraw)->ops->resizewindow = NULL;
480: simg = (PetscImage)(*sdraw)->data;
481: PetscMemcpy(simg->buffer,pimg->buffer,(size_t)(pimg->w*pimg->h));
482: return(0);
483: }
487: static PetscErrorCode PetscDrawRestoreSingleton_Image(PetscDraw draw,PetscDraw *sdraw)
488: {
489: PetscImage pimg = (PetscImage)draw->data;
490: PetscImage simg = (PetscImage)(*sdraw)->data;
494: PetscMemcpy(pimg->buffer,simg->buffer,(size_t)(pimg->w*pimg->h));
495: PetscDrawDestroy(sdraw);
496: return(0);
497: }
499: /*
502: static PetscErrorCode PetscDrawSave_Image(PetscDraw draw)
503: {
505: return(0);
506: }*/
507: #define PetscDrawSave_Image NULL
511: static PetscErrorCode PetscDrawGetImage_Image(PetscDraw draw,unsigned char palette[256][3],unsigned int *w,unsigned int *h,unsigned char *pixels[])
512: {
513: PetscImage img = (PetscImage)draw->data;
514: unsigned char *buffer = NULL;
515: PetscMPIInt rank,size;
519: if (w) *w = (unsigned int)img->w;
520: if (h) *h = (unsigned int)img->h;
521: if (pixels) *pixels = NULL;
522: MPI_Comm_rank(PetscObjectComm((PetscObject)draw),&rank);
523: if (!rank) {
524: PetscMemcpy(palette,img->palette,sizeof(img->palette));
525: PetscMalloc1((size_t)(img->w*img->h),&buffer);
526: if (pixels) *pixels = buffer;
527: }
528: MPI_Comm_size(PetscObjectComm((PetscObject)draw),&size);
529: if (size == 1) {
530: PetscMemcpy(buffer,img->buffer,(size_t)(img->w*img->h));
531: } else {
532: MPI_Reduce(img->buffer,buffer,img->w*img->h,MPI_UNSIGNED_CHAR,MPI_MAX,0,PetscObjectComm((PetscObject)draw));
533: }
534: return(0);
535: }
537: static struct _PetscDrawOps DvOps = {
538: PetscDrawSetDoubleBuffer_Image,
539: PetscDrawFlush_Image,
540: PetscDrawLine_Image,
541: PetscDrawLineSetWidth_Image,
542: PetscDrawLineGetWidth_Image,
543: PetscDrawPoint_Image,
544: PetscDrawPointSetSize_Image,
545: PetscDrawString_Image,
546: PetscDrawStringVertical_Image,
547: PetscDrawStringSetSize_Image,
548: PetscDrawStringGetSize_Image,
549: PetscDrawSetViewport_Image,
550: PetscDrawClear_Image,
551: PetscDrawRectangle_Image,
552: PetscDrawTriangle_Image,
553: PetscDrawEllipse_Image,
554: PetscDrawGetMouseButton_Image,
555: PetscDrawPause_Image,
556: PetscDrawBeginPage_Image,
557: PetscDrawEndPage_Image,
558: PetscDrawGetPopup_Image,
559: PetscDrawSetTitle_Image,
560: PetscDrawCheckResizedWindow_Image,
561: PetscDrawResizeWindow_Image,
562: PetscDrawDestroy_Image,
563: PetscDrawView_Image,
564: PetscDrawGetSingleton_Image,
565: PetscDrawRestoreSingleton_Image,
566: PetscDrawSave_Image,
567: PetscDrawGetImage_Image,
568: PetscDrawSetCoordinates_Image,
569: PetscDrawArrow_Image,
570: PetscDrawCoordinateToPixel_Image,
571: PetscDrawPixelToCoordinate_Image,
572: PetscDrawPointPixel_Image,
573: PetscDrawStringBoxed_Image
574: };
576: static const unsigned char BasicColors[PETSC_DRAW_BASIC_COLORS][3] = {
577: { 255, 255, 255 }, /* white */
578: { 0, 0, 0 }, /* black */
579: { 255, 0, 0 }, /* red */
580: { 0, 255, 0 }, /* green */
581: { 0, 255, 255 }, /* cyan */
582: { 0, 0, 255 }, /* blue */
583: { 255, 0, 255 }, /* magenta */
584: { 127, 255, 212 }, /* aquamarine */
585: { 34, 139, 34 }, /* forestgreen */
586: { 255, 165, 0 }, /* orange */
587: { 238, 130, 238 }, /* violet */
588: { 165, 42, 42 }, /* brown */
589: { 255, 192, 203 }, /* pink */
590: { 255, 127, 80 }, /* coral */
591: { 190, 190, 190 }, /* gray */
592: { 255, 255, 0 }, /* yellow */
593: { 255, 215, 0 }, /* gold */
594: { 255, 182, 193 }, /* lightpink */
595: { 72, 209, 204 }, /* mediumturquoise */
596: { 240, 230, 140 }, /* khaki */
597: { 105, 105, 105 }, /* dimgray */
598: { 54, 205, 50 }, /* yellowgreen */
599: { 135, 206, 235 }, /* skyblue */
600: { 0, 100, 0 }, /* darkgreen */
601: { 0, 0, 128 }, /* navyblue */
602: { 244, 164, 96 }, /* sandybrown */
603: { 95, 158, 160 }, /* cadetblue */
604: { 176, 224, 230 }, /* powderblue */
605: { 255, 20, 147 }, /* deeppink */
606: { 216, 191, 216 }, /* thistle */
607: { 50, 205, 50 }, /* limegreen */
608: { 255, 240, 245 }, /* lavenderblush */
609: { 221, 160, 221 }, /* plum */
610: };
613: /*MC
614: PETSC_DRAW_IMAGE - PETSc graphics device that uses a raster buffer
616: Options Database Keys:
617: . -draw_size w,h - size of image in pixels
619: Level: beginner
621: .seealso: PetscDrawOpenImage(), PetscDrawSetFromOptions()
622: M*/
623: PETSC_EXTERN PetscErrorCode PetscDrawCreate_Image(PetscDraw);
627: PETSC_EXTERN PetscErrorCode PetscDrawCreate_Image(PetscDraw draw)
628: {
629: PetscImage img;
630: int w = draw->w, h = draw->h;
631: PetscInt size[2], nsize = 2;
632: PetscBool set;
636: draw->pause = 0;
637: draw->coor_xl = 0; draw->coor_xr = 1;
638: draw->coor_yl = 0; draw->coor_yr = 1;
639: draw->port_xl = 0; draw->port_xr = 1;
640: draw->port_yl = 0; draw->port_yr = 1;
642: size[0] = w; if (size[0] < 1) size[0] = 300;
643: size[1] = h; if (size[1] < 1) size[1] = size[0];
644: PetscOptionsGetIntArray(((PetscObject)draw)->options,((PetscObject)draw)->prefix,"-draw_size",size,&nsize,&set);
645: if (set && nsize == 1) size[1] = size[0];
646: if (size[0] < 1) size[0] = 300;
647: if (size[1] < 1) size[1] = size[0];
648: draw->w = w = size[0]; draw->x = 0;
649: draw->h = h = size[1]; draw->x = 0;
651: PetscNewLog(draw,&img);
652: PetscMemcpy(draw->ops,&DvOps,sizeof(DvOps));
653: draw->data = (void*)img;
655: img->w = w; img->h = h;
656: PetscCalloc1((size_t)(img->w*img->h),&img->buffer);
657: PetscImageSetClip(img,0,0,img->w,img->h);
658: {
659: int i,k,ncolors = 256-PETSC_DRAW_BASIC_COLORS;
660: unsigned char R[256-PETSC_DRAW_BASIC_COLORS];
661: unsigned char G[256-PETSC_DRAW_BASIC_COLORS];
662: unsigned char B[256-PETSC_DRAW_BASIC_COLORS];
663: PetscDrawUtilitySetCmap(NULL,ncolors,R,G,B);
664: for (k=0; k<PETSC_DRAW_BASIC_COLORS; k++) {
665: img->palette[k][0] = BasicColors[k][0];
666: img->palette[k][1] = BasicColors[k][1];
667: img->palette[k][2] = BasicColors[k][2];
668: }
669: for (i=0; i<ncolors; i++, k++) {
670: img->palette[k][0] = R[i];
671: img->palette[k][1] = G[i];
672: img->palette[k][2] = B[i];
673: }
674: }
676: if (!draw->savefilename ){PetscDrawSetSave(draw,"");}
677: return(0);
678: }
682: /*@C
683: PetscDrawOpenImage - Opens an image for use with the PetscDraw routines.
685: Collective on MPI_Comm
687: Input Parameters:
688: + comm - the communicator that will share image
689: - filename - optional name of the file
690: - w, h - the image width and height in pixels
692: Output Parameters:
693: . draw - the drawing context.
695: Level: beginner
697: .seealso: PetscDrawSetSave(), PetscDrawSetFromOptions(), PetscDrawCreate(), PetscDrawDestroy()
698: @*/
699: PetscErrorCode PetscDrawOpenImage(MPI_Comm comm,const char filename[],int w,int h,PetscDraw *draw)
700: {
704: PetscDrawCreate(comm,NULL,NULL,0,0,w,h,draw);
705: PetscDrawSetType(*draw,PETSC_DRAW_IMAGE);
706: PetscDrawSetSave(*draw,filename);
707: return(0);
708: }