3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * GROningen Mixture of Alchemy and Childrens' Stories
46 #include "gmx_fatal.h"
52 #define round(a) (int)(a+0.5)
54 static const char mapper[] = "ABCDEFGHIJKLMNOPQRSTUVWXYZabcdefghijklmnopqrstuvwxyz0123456789!@#$%^&*()-_=+{}|;:',<.>/?";
55 #define NMAP (long int)strlen(mapper)
57 #define MAX_XPM_LINELENGTH 4096
59 real **mk_matrix(int nx, int ny, gmx_bool b1D)
77 void done_matrix(int nx, real ***m)
87 void clear_matrix(int nx, int ny, real **m)
96 gmx_bool matelmt_cmp(t_xpmelmt e1, t_xpmelmt e2)
98 return (e1.c1 == e2.c1) && (e1.c2 == e2.c2);
101 t_matelmt searchcmap(int n,t_mapping map[],t_xpmelmt c)
106 if (matelmt_cmp(map[i].code, c))
112 int getcmap(FILE *in,const char *fn,t_mapping **map)
116 char code[STRLEN],desc[STRLEN];
120 if (fgets2(line,STRLEN-1,in) == NULL)
121 gmx_fatal(FARGS,"Not enough lines in colormap file %s"
122 "(just wanted to read number of entries)",fn);
123 sscanf(line,"%d",&n);
125 for(i=0; (i<n); i++) {
126 if (fgets2(line,STRLEN-1,in) == NULL)
127 gmx_fatal(FARGS,"Not enough lines in colormap file %s"
128 "(should be %d, found only %d)",fn,n+1,i);
129 sscanf(line,"%s%s%lf%lf%lf",code,desc,&r,&g,&b);
130 m[i].code.c1=code[0];
132 m[i].desc=strdup(desc);
142 int readcmap(const char *fn,t_mapping **map)
148 n=getcmap(in,fn,map);
154 void printcmap(FILE *out,int n,t_mapping map[])
158 fprintf(out,"%d\n",n);
160 fprintf(out,"%c%c %20s %10g %10g %10g\n",
161 map[i].code.c1?map[i].code.c1:' ',
162 map[i].code.c2?map[i].code.c2:' ',
163 map[i].desc,map[i].rgb.r,map[i].rgb.g,map[i].rgb.b);
166 void writecmap(const char *fn,int n,t_mapping map[])
170 out=gmx_fio_fopen(fn,"w");
171 printcmap(out,n,map);
175 static char *fgetline(char **line,int llmax,int *llalloc,FILE *in)
179 if (llmax > *llalloc)
181 srenew(*line,llmax+1);
184 fg=fgets(*line,llmax,in);
190 void skipstr(char *line)
196 while((line[c] != ' ') && (line[c] != '\0'))
199 while(line[c] != '\0')
207 char *line2string(char **line)
212 while (((*line)[0] != '\"' ) && ( (*line)[0] != '\0' ))
215 if ((*line)[0] != '\"')
220 while (( (*line)[i] != '\"' ) && ( (*line)[i] != '\0' ))
223 if ((*line)[i] != '\"')
232 void parsestring(char *line,const char *label, char *string)
234 if (strstr(line,label)) {
235 if (strstr(line,label) < strchr(line,'\"')) {
242 void read_xpm_entry(FILE *in,t_matrix *mm)
245 char *line_buf=NULL,*line=NULL,*str,buf[256]={0};
246 int i,m,col_len,nch,n_axis_x,n_axis_y,llmax;
250 gmx_bool bGetOnWithIt,bSetLine;
265 while ((NULL != fgetline(&line_buf,llmax,&llalloc,in)) &&
266 (strncmp(line_buf,"static",6) != 0)) {
268 parsestring(line,"title",(mm->title));
269 parsestring(line,"legend",(mm->legend));
270 parsestring(line,"x-label",(mm->label_x));
271 parsestring(line,"y-label",(mm->label_y));
272 parsestring(line,"type",buf);
275 if (!line || strncmp(line,"static",6) != 0)
277 gmx_input("Invalid XPixMap");
280 if (buf[0] && (gmx_strcasecmp(buf,"Discrete")==0))
284 fprintf(debug,"%s %s %s %s\n",
285 mm->title,mm->legend,mm->label_x,mm->label_y);
289 while (!bGetOnWithIt && (NULL != fgetline(&line_buf,llmax,&llalloc,in))) {
291 while (( line[0] != '\"' ) && ( line[0] != '\0' ))
294 if ( line[0] == '\"' ) {
296 sscanf(line,"%d %d %d %d",&(mm->nx),&(mm->ny),&(mm->nmap),&nch);
299 gmx_fatal(FARGS,"Sorry can only read xpm's with at most 2 caracters per pixel\n");
301 if (mm->nx <= 0 || mm->ny <= 0 )
303 gmx_fatal(FARGS,"Dimensions of xpm-file have to be larger than 0\n");
305 llmax = max(STRLEN,mm->nx+10);
310 fprintf(debug,"mm->nx %d mm->ny %d mm->nmap %d nch %d\n",
311 mm->nx,mm->ny,mm->nmap,nch);
316 while ((m < mm->nmap) && (NULL != fgetline(&line_buf,llmax,&llalloc,in))) {
317 line = strchr(line_buf,'\"');
320 /* Read xpm color map entry */
321 map[m].code.c1 = line[0];
325 map[m].code.c2 = line[1];
327 str = strchr(line,'#');
331 while (isxdigit(str[col_len]))
334 sscanf(line,"%*s #%2x%2x%2x",&r,&g,&b);
335 map[m].rgb.r=r/255.0;
336 map[m].rgb.g=g/255.0;
337 map[m].rgb.b=b/255.0;
338 } else if (col_len==12) {
339 sscanf(line,"%*s #%4x%4x%4x",&r,&g,&b);
340 map[m].rgb.r=r/65535.0;
341 map[m].rgb.g=g/65535.0;
342 map[m].rgb.b=b/65535.0;
344 gmx_file("Unsupported or invalid colormap in X PixMap");
346 str = strchr(line,'c');
350 gmx_file("Unsupported or invalid colormap in X PixMap");
351 fprintf(stderr,"Using white for color \"%s",str);
356 line=strchr(line,'\"');
359 map[m].desc = strdup(line);
364 gmx_fatal(FARGS,"Number of read colors map entries (%d) does not match the number in the header (%d)",m,mm->nmap);
367 /* Read axes, if there are any */
376 if (strstr(line,"x-axis")) {
377 line=strstr(line,"x-axis");
379 if (mm->axis_x==NULL)
380 snew(mm->axis_x,mm->nx + 1);
381 while (sscanf(line,"%lf",&u)==1) {
382 if (n_axis_x > mm->nx) {
383 gmx_fatal(FARGS,"Too many x-axis labels in xpm (max %d)",mm->nx);
384 } else if (n_axis_x == mm->nx) {
385 mm->flags |= MAT_SPATIAL_X;
387 mm->axis_x[n_axis_x] = u;
392 else if (strstr(line,"y-axis")) {
393 line=strstr(line,"y-axis");
395 if (mm->axis_y==NULL)
396 snew(mm->axis_y,mm->ny + 1);
397 while (sscanf(line,"%lf",&u)==1) {
398 if (n_axis_y > mm->ny) {
399 gmx_file("Too many y-axis labels in xpm");
400 } else if (n_axis_y == mm->ny) {
401 mm->flags |= MAT_SPATIAL_Y;
403 mm->axis_y[n_axis_y] = u;
408 } while ((line[0] != '\"') && (NULL != fgetline(&line_buf,llmax,&llalloc,in)));
411 snew(mm->matrix,mm->nx);
412 for(i=0; i<mm->nx; i++)
413 snew(mm->matrix[i],mm->ny);
420 if(m%(1+mm->ny/100)==0)
421 fprintf(stderr,"%3d%%\b\b\b\b",(100*(mm->ny-m))/mm->ny);
422 while ((line[0] != '\"') && (line[0] != '\0'))
425 gmx_fatal(FARGS,"Not enough caracters in row %d of the matrix\n",m+1);
428 for(i=0; i<mm->nx; i++) {
434 mm->matrix[i][m]=searchcmap(mm->nmap,mm->map,c);
438 } while ((m>=0) && (NULL != fgetline(&line_buf,llmax,&llalloc,in)));
440 gmx_incons("Not enough rows in the matrix");
445 int read_xpm_matrix(const char *fnm,t_matrix **matrix)
452 in=gmx_fio_fopen(fnm,"r");
455 while (NULL != fgetline(&line,STRLEN,&llalloc,in)) {
456 if (strstr(line,"/* XPM */")) {
457 srenew(*matrix,nmat+1);
458 read_xpm_entry(in,&(*matrix)[nmat]);
465 gmx_file("Invalid XPixMap");
472 real **matrix2real(t_matrix *matrix,real **mat)
483 for(i=0; i<nmap; i++) {
484 if ((map[i].desc==NULL) || (sscanf(map[i].desc,"%lf",&tmp)!=1)) {
485 fprintf(stderr,"Could not convert matrix to reals,\n"
486 "color map entry %d has a non-real description: \"%s\"\n",
495 snew(mat,matrix->nx);
496 for(i=0; i<matrix->nx; i++)
497 snew(mat[i],matrix->ny);
499 for(i=0; i<matrix->nx; i++)
500 for(j=0; j<matrix->ny; j++)
501 mat[i][j]=rmap[matrix->matrix[i][j]];
505 fprintf(stderr,"Converted a %dx%d matrix with %d levels to reals\n",
506 matrix->nx,matrix->ny,nmap);
511 void write_xpm_header(FILE *out,
512 const char *title,const char *legend,
513 const char *label_x,const char *label_y,
516 fprintf(out, "/* XPM */\n");
517 fprintf(out, "/* Generated by %s */\n",Program());
518 fprintf(out, "/* This file can be converted to EPS by the GROMACS program xpm2ps */\n");
519 fprintf(out, "/* title: \"%s\" */\n",title);
520 fprintf(out, "/* legend: \"%s\" */\n",legend);
521 fprintf(out, "/* x-label: \"%s\" */\n",label_x);
522 fprintf(out, "/* y-label: \"%s\" */\n",label_y);
524 fprintf(out,"/* type: \"Discrete\" */\n");
526 fprintf(out,"/* type: \"Continuous\" */\n");
529 static int calc_nmid(int nlevels,real lo,real mid,real hi)
531 /* Take care that we have at least 1 entry in the mid to hi range
533 return min(max(0,((mid-lo)/(hi-lo))*(nlevels-1)),nlevels-1);
536 void write_xpm_map3(FILE *out,int n_x,int n_y,int *nlevels,
537 real lo,real mid,real hi,
538 t_rgb rlo,t_rgb rmid,t_rgb rhi)
541 real r,g,b,clev_lo,clev_hi;
543 if (*nlevels > NMAP*NMAP) {
544 fprintf(stderr,"Warning, too many levels (%d) in matrix, using %d only\n",
545 *nlevels,(int)(NMAP*NMAP));
548 else if (*nlevels < 2) {
549 fprintf(stderr,"Warning, too few levels (%d) in matrix, using 2 instead\n",
553 if (!((mid >= lo) && (mid < hi)))
554 gmx_fatal(FARGS,"Lo: %f, Mid: %f, Hi: %f\n",lo,mid,hi);
556 fprintf(out,"static char *gromacs_xpm[] = {\n");
557 fprintf(out,"\"%d %d %d %d\",\n",
558 n_x,n_y,*nlevels,(*nlevels <= NMAP) ? 1 : 2);
560 nmid = calc_nmid(*nlevels,lo,mid,hi);
562 clev_hi = (*nlevels - 1 - nmid);
563 for(i=0; (i<nmid); i++) {
564 r = rlo.r+(i*(rmid.r-rlo.r)/clev_lo);
565 g = rlo.g+(i*(rmid.g-rlo.g)/clev_lo);
566 b = rlo.b+(i*(rmid.b-rlo.b)/clev_lo);
567 fprintf(out,"\"%c%c c #%02X%02X%02X \" /* \"%.3g\" */,\n",
569 (*nlevels <= NMAP) ? ' ' : mapper[i/NMAP],
570 (unsigned int)round(255*r),
571 (unsigned int)round(255*g),
572 (unsigned int)round(255*b),
573 ((nmid - i)*lo + i*mid)/clev_lo);
575 for(i=0; (i<(*nlevels-nmid)); i++) {
576 r = rmid.r+(i*(rhi.r-rmid.r)/clev_hi);
577 g = rmid.g+(i*(rhi.g-rmid.g)/clev_hi);
578 b = rmid.b+(i*(rhi.b-rmid.b)/clev_hi);
579 fprintf(out,"\"%c%c c #%02X%02X%02X \" /* \"%.3g\" */,\n",
580 mapper[(i+nmid) % NMAP],
581 (*nlevels <= NMAP) ? ' ' : mapper[(i+nmid)/NMAP],
582 (unsigned int)round(255*r),
583 (unsigned int)round(255*g),
584 (unsigned int)round(255*b),
585 ((*nlevels - 1 - nmid - i)*mid + i*hi)/clev_hi);
589 static void pr_simple_cmap(FILE *out,real lo,real hi,int nlevel,t_rgb rlo,
595 for(i=0; (i<nlevel); i++) {
596 fac = (i+1.0)/(nlevel);
597 r = rlo.r+fac*(rhi.r-rlo.r);
598 g = rlo.g+fac*(rhi.g-rlo.g);
599 b = rlo.b+fac*(rhi.b-rlo.b);
600 fprintf(out,"\"%c%c c #%02X%02X%02X \" /* \"%.3g\" */,\n",
601 mapper[(i+i0) % NMAP],
602 (nlevel <= NMAP) ? ' ' : mapper[(i+i0)/NMAP],
603 (unsigned int)round(255*r),
604 (unsigned int)round(255*g),
605 (unsigned int)round(255*b),
610 static void pr_discrete_cmap(FILE *out,int *nlevel,int i0)
613 { 1.0, 1.0, 1.0 }, /* white */
614 { 1.0, 0.0, 0.0 }, /* red */
615 { 1.0, 1.0, 0.0 }, /* yellow */
616 { 0.0, 0.0, 1.0 }, /* blue */
617 { 0.0, 1.0, 0.0 }, /* green */
618 { 1.0, 0.0, 1.0 }, /* purple */
619 { 1.0, 0.4, 0.0 }, /* orange */
620 { 0.0, 1.0, 1.0 }, /* cyan */
621 { 1.0, 0.4, 0.4 }, /* pink */
622 { 1.0, 1.0, 0.0 }, /* yellow */
623 { 0.4, 0.4, 1.0 }, /* lightblue */
624 { 0.4, 1.0, 0.4 }, /* lightgreen */
625 { 1.0, 0.4, 1.0 }, /* lightpurple */
626 { 1.0, 0.7, 0.4 }, /* lightorange */
627 { 0.4, 1.0, 1.0 }, /* lightcyan */
628 { 0.0, 0.0, 0.0 } /* black */
633 *nlevel = min(16,*nlevel);
635 for(i=0; (i<n); i++) {
636 fprintf(out,"\"%c%c c #%02X%02X%02X \" /* \"%3d\" */,\n",
637 mapper[(i+i0) % NMAP],
638 (n <= NMAP) ? ' ' : mapper[(i+i0)/NMAP],
639 (unsigned int)round(255*rgbd[i].r),
640 (unsigned int)round(255*rgbd[i].g),
641 (unsigned int)round(255*rgbd[i].b),
648 void write_xpm_map_split(FILE *out,int n_x,int n_y,
649 int *nlevel_top,real lo_top,real hi_top,
650 t_rgb rlo_top,t_rgb rhi_top,
651 gmx_bool bDiscreteColor,
652 int *nlevel_bot,real lo_bot,real hi_bot,
653 t_rgb rlo_bot,t_rgb rhi_bot)
658 ntot = *nlevel_top + *nlevel_bot;
660 gmx_fatal(FARGS,"Warning, too many levels (%d) in matrix",ntot);
662 fprintf(out,"static char *gromacs_xpm[] = {\n");
663 fprintf(out,"\"%d %d %d %d\",\n",n_x,n_y,ntot,1);
666 pr_discrete_cmap(out,nlevel_bot,0);
668 pr_simple_cmap(out,lo_bot,hi_bot,*nlevel_bot,rlo_bot,rhi_bot,0);
670 pr_simple_cmap(out,lo_top,hi_top,*nlevel_top,rlo_top,rhi_top,*nlevel_bot);
674 void write_xpm_map(FILE *out,int n_x, int n_y,int *nlevels,real lo,real hi,
680 if (*nlevels > NMAP*NMAP) {
681 fprintf(stderr,"Warning, too many levels (%d) in matrix, using %d only\n",
682 *nlevels,(int)(NMAP*NMAP));
685 else if (*nlevels < 2) {
686 fprintf(stderr,"Warning, too few levels (%d) in matrix, using 2 instead\n",*nlevels);
690 fprintf(out,"static char *gromacs_xpm[] = {\n");
691 fprintf(out,"\"%d %d %d %d\",\n",
692 n_x,n_y,*nlevels,(*nlevels <= NMAP) ? 1 : 2);
694 invlevel=1.0/(*nlevels-1);
695 for(i=0; (i<*nlevels); i++) {
697 r=(nlo*rlo.r+i*rhi.r)*invlevel;
698 g=(nlo*rlo.g+i*rhi.g)*invlevel;
699 b=(nlo*rlo.b+i*rhi.b)*invlevel;
700 fprintf(out,"\"%c%c c #%02X%02X%02X \" /* \"%.3g\" */,\n",
701 mapper[i % NMAP],(*nlevels <= NMAP) ? ' ' : mapper[i/NMAP],
702 (unsigned int)round(255*r),
703 (unsigned int)round(255*g),
704 (unsigned int)round(255*b),
705 (nlo*lo+i*hi)*invlevel);
709 void write_xpm_axis(FILE *out, const char *axis, gmx_bool bSpatial, int n,
715 for(i=0;i<(bSpatial ? n+1 : n);i++) {
719 fprintf(out,"/* %s-axis: ",axis);
721 fprintf(out,"%g ",label[i]);
727 void write_xpm_data(FILE *out, int n_x, int n_y, real **matrix,
728 real lo, real hi, int nlevels)
733 invlevel=(nlevels-1)/(hi-lo);
734 for(j=n_y-1; (j>=0); j--) {
736 fprintf(stderr,"%3d%%\b\b\b\b",(100*(n_y-j))/n_y);
738 for(i=0; (i<n_x); i++) {
739 c=gmx_nint((matrix[i][j]-lo)*invlevel);
741 if (c>=nlevels) c=nlevels-1;
743 fprintf(out,"%c",mapper[c]);
745 fprintf(out,"%c%c",mapper[c % NMAP],mapper[c / NMAP]);
748 fprintf(out,"\",\n");
754 void write_xpm_data3(FILE *out,int n_x,int n_y,real **matrix,
755 real lo,real mid,real hi,int nlevels)
758 real invlev_lo,invlev_hi;
760 nmid = calc_nmid(nlevels,lo,mid,hi);
761 invlev_hi=(nlevels-1-nmid)/(hi-mid);
762 invlev_lo=(nmid)/(mid-lo);
764 for(j=n_y-1; (j>=0); j--) {
766 fprintf(stderr,"%3d%%\b\b\b\b",(100*(n_y-j))/n_y);
768 for(i=0; (i<n_x); i++) {
769 if (matrix[i][j] >= mid)
770 c=nmid+gmx_nint((matrix[i][j]-mid)*invlev_hi);
771 else if (matrix[i][j] >= lo)
772 c=gmx_nint((matrix[i][j]-lo)*invlev_lo);
781 fprintf(out,"%c",mapper[c]);
783 fprintf(out,"%c%c",mapper[c % NMAP],mapper[c / NMAP]);
786 fprintf(out,"\",\n");
792 void write_xpm_data_split(FILE *out,int n_x,int n_y,real **matrix,
793 real lo_top,real hi_top,int nlevel_top,
794 real lo_bot,real hi_bot,int nlevel_bot)
797 real invlev_top,invlev_bot;
799 invlev_top=(nlevel_top-1)/(hi_top-lo_top);
800 invlev_bot=(nlevel_bot-1)/(hi_bot-lo_bot);
802 for(j=n_y-1; (j>=0); j--) {
803 if(j % (1+n_y/100)==0)
804 fprintf(stderr,"%3d%%\b\b\b\b",(100*(n_y-j))/n_y);
806 for(i=0; (i<n_x); i++) {
808 c = nlevel_bot+round((matrix[i][j]-lo_top)*invlev_top);
809 if ((c < nlevel_bot) || (c >= nlevel_bot+nlevel_top))
810 gmx_fatal(FARGS,"Range checking i = %d, j = %d, c = %d, bot = %d, top = %d matrix[i,j] = %f",i,j,c,nlevel_bot,nlevel_top,matrix[i][j]);
813 c = round((matrix[i][j]-lo_bot)*invlev_bot);
814 if ((c < 0) || (c >= nlevel_bot+nlevel_bot))
815 gmx_fatal(FARGS,"Range checking i = %d, j = %d, c = %d, bot = %d, top = %d matrix[i,j] = %f",i,j,c,nlevel_bot,nlevel_top,matrix[i][j]);
820 fprintf(out,"%c",mapper[c]);
823 fprintf(out,"\",\n");
829 void write_xpm_m(FILE *out, t_matrix m)
831 /* Writes a t_matrix struct to .xpm file */
837 bOneChar=(m.map[0].code.c2 == 0);
838 write_xpm_header(out,m.title,m.legend,m.label_x,m.label_y,
840 fprintf(out,"static char *gromacs_xpm[] = {\n");
841 fprintf(out,"\"%d %d %d %d\",\n",m.nx,m.ny,m.nmap,bOneChar ? 1 : 2);
842 for(i=0; (i<m.nmap); i++)
843 fprintf(out,"\"%c%c c #%02X%02X%02X \" /* \"%s\" */,\n",
845 bOneChar ? ' ' : m.map[i].code.c2,
846 (unsigned int)round(m.map[i].rgb.r*255),
847 (unsigned int)round(m.map[i].rgb.g*255),
848 (unsigned int)round(m.map[i].rgb.b*255),m.map[i].desc);
849 write_xpm_axis(out,"x",m.flags & MAT_SPATIAL_X,m.nx,m.axis_x);
850 write_xpm_axis(out,"y",m.flags & MAT_SPATIAL_Y,m.ny,m.axis_y);
851 for(j=m.ny-1; (j>=0); j--) {
852 if(j%(1+m.ny/100)==0)
853 fprintf(stderr,"%3d%%\b\b\b\b",(100*(m.ny-j))/m.ny);
856 for(i=0; (i<m.nx); i++)
857 fprintf(out,"%c",m.map[m.matrix[i][j]].code.c1);
859 for(i=0; (i<m.nx); i++) {
860 c=m.map[m.matrix[i][j]].code;
861 fprintf(out,"%c%c",c.c1,c.c2);
864 fprintf(out,"\",\n");
870 void write_xpm3(FILE *out,unsigned int flags,
871 const char *title,const char *legend,
872 const char *label_x,const char *label_y,
873 int n_x,int n_y,real axis_x[],real axis_y[],
874 real *matrix[],real lo,real mid,real hi,
875 t_rgb rlo,t_rgb rmid,t_rgb rhi,int *nlevels)
878 * Writes a colormap varying as rlo -> rmid -> rhi.
882 gmx_fatal(FARGS,"hi (%g) <= lo (%g)",hi,lo);
884 write_xpm_header(out,title,legend,label_x,label_y,FALSE);
885 write_xpm_map3(out,n_x,n_y,nlevels,lo,mid,hi,rlo,rmid,rhi);
886 write_xpm_axis(out,"x",flags & MAT_SPATIAL_X,n_x,axis_x);
887 write_xpm_axis(out,"y",flags & MAT_SPATIAL_Y,n_y,axis_y);
888 write_xpm_data3(out,n_x,n_y,matrix,lo,mid,hi,*nlevels);
891 void write_xpm_split(FILE *out,unsigned int flags,
892 const char *title,const char *legend,
893 const char *label_x,const char *label_y,
894 int n_x,int n_y,real axis_x[],real axis_y[],
896 real lo_top,real hi_top,int *nlevel_top,
897 t_rgb rlo_top,t_rgb rhi_top,
898 real lo_bot,real hi_bot,int *nlevel_bot,
899 gmx_bool bDiscreteColor,
900 t_rgb rlo_bot,t_rgb rhi_bot)
903 * Writes a colormap varying as rlo -> rmid -> rhi.
906 if (hi_top <= lo_top)
907 gmx_fatal(FARGS,"hi_top (%g) <= lo_top (%g)",hi_top,lo_top);
908 if (hi_bot <= lo_bot)
909 gmx_fatal(FARGS,"hi_bot (%g) <= lo_bot (%g)",hi_bot,lo_bot);
910 if (bDiscreteColor && (*nlevel_bot >= 16))
911 gmx_impl("Can not plot more than 16 discrete colors");
913 write_xpm_header(out,title,legend,label_x,label_y,FALSE);
914 write_xpm_map_split(out,n_x,n_y,nlevel_top,lo_top,hi_top,rlo_top,rhi_top,
915 bDiscreteColor,nlevel_bot,lo_bot,hi_bot,rlo_bot,rhi_bot);
916 write_xpm_axis(out,"x",flags & MAT_SPATIAL_X,n_x,axis_x);
917 write_xpm_axis(out,"y",flags & MAT_SPATIAL_Y,n_y,axis_y);
918 write_xpm_data_split(out,n_x,n_y,matrix,lo_top,hi_top,*nlevel_top,
919 lo_bot,hi_bot,*nlevel_bot);
922 void write_xpm(FILE *out,unsigned int flags,
923 const char *title,const char *legend,
924 const char *label_x,const char *label_y,
925 int n_x,int n_y,real axis_x[],real axis_y[],
926 real *matrix[],real lo,real hi,
927 t_rgb rlo,t_rgb rhi,int *nlevels)
931 * legend label for the continuous legend
932 * label_x label for the x-axis
933 * label_y label for the y-axis
934 * n_x, n_y size of the matrix
935 * axis_x[] the x-ticklabels
936 * axis_y[] the y-ticklables
937 * *matrix[] element x,y is matrix[x][y]
938 * lo output lower than lo is set to lo
939 * hi output higher than hi is set to hi
940 * rlo rgb value for level lo
941 * rhi rgb value for level hi
942 * nlevels number of color levels for the output
946 gmx_fatal(FARGS,"hi (%f) <= lo (%f)",hi,lo);
948 write_xpm_header(out,title,legend,label_x,label_y,FALSE);
949 write_xpm_map(out,n_x,n_y,nlevels,lo,hi,rlo,rhi);
950 write_xpm_axis(out,"x",flags & MAT_SPATIAL_X,n_x,axis_x);
951 write_xpm_axis(out,"y",flags & MAT_SPATIAL_Y,n_y,axis_y);
952 write_xpm_data(out,n_x,n_y,matrix,lo,hi,*nlevels);