2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * Copyright (c) 2013, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
51 #include "gmx_fatal.h"
58 #include "gromacs/utility/exceptions.h"
59 #include "gromacs/utility/programinfo.h"
61 #define round(a) (int)(a+0.5)
63 static const char mapper[] = "ABCDEFGHIJKLMNOPQRSTUVWXYZabcdefghijklmnopqrstuvwxyz0123456789!@#$%^&*()-_=+{}|;:',<.>/?";
64 #define NMAP (long int)strlen(mapper)
66 #define MAX_XPM_LINELENGTH 4096
68 real **mk_matrix(int nx, int ny, gmx_bool b1D)
79 for (i = 0; (i < nx); i++)
94 void done_matrix(int nx, real ***m)
98 for (i = 0; (i < nx); i++)
106 void clear_matrix(int nx, int ny, real **m)
110 for (x = 0; x < nx; x++)
112 for (y = 0; y < ny; y++)
119 gmx_bool matelmt_cmp(t_xpmelmt e1, t_xpmelmt e2)
121 return (e1.c1 == e2.c1) && (e1.c2 == e2.c2);
124 t_matelmt searchcmap(int n, t_mapping map[], t_xpmelmt c)
128 for (i = 0; (i < n); i++)
130 if (matelmt_cmp(map[i].code, c))
139 int getcmap(FILE *in, const char *fn, t_mapping **map)
143 char code[STRLEN], desc[STRLEN];
147 if (fgets2(line, STRLEN-1, in) == NULL)
149 gmx_fatal(FARGS, "Not enough lines in colormap file %s"
150 "(just wanted to read number of entries)", fn);
152 sscanf(line, "%d", &n);
154 for (i = 0; (i < n); i++)
156 if (fgets2(line, STRLEN-1, in) == NULL)
158 gmx_fatal(FARGS, "Not enough lines in colormap file %s"
159 "(should be %d, found only %d)", fn, n+1, i);
161 sscanf(line, "%s%s%lf%lf%lf", code, desc, &r, &g, &b);
162 m[i].code.c1 = code[0];
164 m[i].desc = strdup(desc);
174 int readcmap(const char *fn, t_mapping **map)
180 n = getcmap(in, fn, map);
186 void printcmap(FILE *out, int n, t_mapping map[])
190 fprintf(out, "%d\n", n);
191 for (i = 0; (i < n); i++)
193 fprintf(out, "%c%c %20s %10g %10g %10g\n",
194 map[i].code.c1 ? map[i].code.c1 : ' ',
195 map[i].code.c2 ? map[i].code.c2 : ' ',
196 map[i].desc, map[i].rgb.r, map[i].rgb.g, map[i].rgb.b);
200 void writecmap(const char *fn, int n, t_mapping map[])
204 out = gmx_fio_fopen(fn, "w");
205 printcmap(out, n, map);
209 static char *fgetline(char **line, int llmax, int *llalloc, FILE *in)
213 if (llmax > *llalloc)
215 srenew(*line, llmax+1);
218 fg = fgets(*line, llmax, in);
224 void skipstr(char *line)
230 while ((line[c] != ' ') && (line[c] != '\0'))
235 while (line[c] != '\0')
243 char *line2string(char **line)
249 while (((*line)[0] != '\"' ) && ( (*line)[0] != '\0' ))
254 if ((*line)[0] != '\"')
261 while (( (*line)[i] != '\"' ) && ( (*line)[i] != '\0' ))
266 if ((*line)[i] != '\"')
279 void parsestring(char *line, const char *label, char *string)
281 if (strstr(line, label))
283 if (strstr(line, label) < strchr(line, '\"'))
286 strcpy(string, line);
291 void read_xpm_entry(FILE *in, t_matrix *mm)
294 char *line_buf = NULL, *line = NULL, *str, buf[256] = {0};
295 int i, m, col_len, nch = 0, n_axis_x, n_axis_y, llmax;
297 unsigned int r, g, b;
299 gmx_bool bGetOnWithIt, bSetLine;
310 mm->bDiscrete = FALSE;
314 while ((NULL != fgetline(&line_buf, llmax, &llalloc, in)) &&
315 (strncmp(line_buf, "static", 6) != 0))
318 parsestring(line, "title", (mm->title));
319 parsestring(line, "legend", (mm->legend));
320 parsestring(line, "x-label", (mm->label_x));
321 parsestring(line, "y-label", (mm->label_y));
322 parsestring(line, "type", buf);
325 if (!line || strncmp(line, "static", 6) != 0)
327 gmx_input("Invalid XPixMap");
330 if (buf[0] && (gmx_strcasecmp(buf, "Discrete") == 0))
332 mm->bDiscrete = TRUE;
337 fprintf(debug, "%s %s %s %s\n",
338 mm->title, mm->legend, mm->label_x, mm->label_y);
342 bGetOnWithIt = FALSE;
343 while (!bGetOnWithIt && (NULL != fgetline(&line_buf, llmax, &llalloc, in)))
346 while (( line[0] != '\"' ) && ( line[0] != '\0' ))
354 sscanf(line, "%d %d %d %d", &(mm->nx), &(mm->ny), &(mm->nmap), &nch);
357 gmx_fatal(FARGS, "Sorry can only read xpm's with at most 2 caracters per pixel\n");
359 if (mm->nx <= 0 || mm->ny <= 0)
361 gmx_fatal(FARGS, "Dimensions of xpm-file have to be larger than 0\n");
363 llmax = std::max(STRLEN, mm->nx+10);
369 fprintf(debug, "mm->nx %d mm->ny %d mm->nmap %d nch %d\n",
370 mm->nx, mm->ny, mm->nmap, nch);
374 gmx_fatal(FARGS, "Number of characters per pixel not found in xpm\n");
380 while ((m < mm->nmap) && (NULL != fgetline(&line_buf, llmax, &llalloc, in)))
382 line = strchr(line_buf, '\"');
386 /* Read xpm color map entry */
387 map[m].code.c1 = line[0];
394 map[m].code.c2 = line[1];
397 str = strchr(line, '#');
402 while (isxdigit(str[col_len]))
408 sscanf(line, "%*s #%2x%2x%2x", &r, &g, &b);
409 map[m].rgb.r = r/255.0;
410 map[m].rgb.g = g/255.0;
411 map[m].rgb.b = b/255.0;
413 else if (col_len == 12)
415 sscanf(line, "%*s #%4x%4x%4x", &r, &g, &b);
416 map[m].rgb.r = r/65535.0;
417 map[m].rgb.g = g/65535.0;
418 map[m].rgb.b = b/65535.0;
422 gmx_file("Unsupported or invalid colormap in X PixMap");
427 str = strchr(line, 'c');
434 gmx_file("Unsupported or invalid colormap in X PixMap");
436 fprintf(stderr, "Using white for color \"%s", str);
441 line = strchr(line, '\"');
444 map[m].desc = strdup(line);
450 gmx_fatal(FARGS, "Number of read colors map entries (%d) does not match the number in the header (%d)", m, mm->nmap);
454 /* Read axes, if there are any */
465 if (strstr(line, "x-axis"))
467 line = strstr(line, "x-axis");
469 if (mm->axis_x == NULL)
471 snew(mm->axis_x, mm->nx + 1);
473 while (sscanf(line, "%lf", &u) == 1)
475 if (n_axis_x > mm->nx)
477 gmx_fatal(FARGS, "Too many x-axis labels in xpm (max %d)", mm->nx);
479 else if (n_axis_x == mm->nx)
481 mm->flags |= MAT_SPATIAL_X;
483 mm->axis_x[n_axis_x] = u;
488 else if (strstr(line, "y-axis"))
490 line = strstr(line, "y-axis");
492 if (mm->axis_y == NULL)
494 snew(mm->axis_y, mm->ny + 1);
496 while (sscanf(line, "%lf", &u) == 1)
498 if (n_axis_y > mm->ny)
500 gmx_file("Too many y-axis labels in xpm");
502 else if (n_axis_y == mm->ny)
504 mm->flags |= MAT_SPATIAL_Y;
506 mm->axis_y[n_axis_y] = u;
512 while ((line[0] != '\"') && (NULL != fgetline(&line_buf, llmax, &llalloc, in)));
515 snew(mm->matrix, mm->nx);
516 for (i = 0; i < mm->nx; i++)
518 snew(mm->matrix[i], mm->ny);
529 if (m%(1+mm->ny/100) == 0)
531 fprintf(stderr, "%3d%%\b\b\b\b", (100*(mm->ny-m))/mm->ny);
533 while ((line[0] != '\"') && (line[0] != '\0'))
539 gmx_fatal(FARGS, "Not enough caracters in row %d of the matrix\n", m+1);
544 for (i = 0; i < mm->nx; i++)
553 c.c2 = line[nch*i+1];
555 mm->matrix[i][m] = searchcmap(mm->nmap, mm->map, c);
560 while ((m >= 0) && (NULL != fgetline(&line_buf, llmax, &llalloc, in)));
563 gmx_incons("Not enough rows in the matrix");
569 int read_xpm_matrix(const char *fnm, t_matrix **matrix)
576 in = gmx_fio_fopen(fnm, "r");
579 while (NULL != fgetline(&line, STRLEN, &llalloc, in))
581 if (strstr(line, "/* XPM */"))
583 srenew(*matrix, nmat+1);
584 read_xpm_entry(in, &(*matrix)[nmat]);
592 gmx_file("Invalid XPixMap");
600 real **matrix2real(t_matrix *matrix, real **mat)
611 for (i = 0; i < nmap; i++)
613 if ((map[i].desc == NULL) || (sscanf(map[i].desc, "%lf", &tmp) != 1))
615 fprintf(stderr, "Could not convert matrix to reals,\n"
616 "color map entry %d has a non-real description: \"%s\"\n",
626 snew(mat, matrix->nx);
627 for (i = 0; i < matrix->nx; i++)
629 snew(mat[i], matrix->ny);
632 for (i = 0; i < matrix->nx; i++)
634 for (j = 0; j < matrix->ny; j++)
636 mat[i][j] = rmap[matrix->matrix[i][j]];
642 fprintf(stderr, "Converted a %dx%d matrix with %d levels to reals\n",
643 matrix->nx, matrix->ny, nmap);
648 void write_xpm_header(FILE *out,
649 const char *title, const char *legend,
650 const char *label_x, const char *label_y,
653 fprintf(out, "/* XPM */\n");
656 gmx::BinaryInformationSettings settings;
657 settings.generatedByHeader(true);
658 settings.linePrefix("/* ");
659 settings.lineSuffix(" */");
660 gmx::printBinaryInformation(out, gmx::ProgramInfo::getInstance(),
663 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
664 fprintf(out, "/* This file can be converted to EPS by the GROMACS program xpm2ps */\n");
665 fprintf(out, "/* title: \"%s\" */\n", title);
666 fprintf(out, "/* legend: \"%s\" */\n", legend);
667 fprintf(out, "/* x-label: \"%s\" */\n", label_x);
668 fprintf(out, "/* y-label: \"%s\" */\n", label_y);
671 fprintf(out, "/* type: \"Discrete\" */\n");
675 fprintf(out, "/* type: \"Continuous\" */\n");
679 static int calc_nmid(int nlevels, real lo, real mid, real hi)
681 /* Take care that we have at least 1 entry in the mid to hi range
683 return std::min(std::max(0, static_cast<int>(((mid-lo)/(hi-lo))*(nlevels-1))),
687 void write_xpm_map3(FILE *out, int n_x, int n_y, int *nlevels,
688 real lo, real mid, real hi,
689 t_rgb rlo, t_rgb rmid, t_rgb rhi)
692 real r, g, b, clev_lo, clev_hi;
694 if (*nlevels > NMAP*NMAP)
696 fprintf(stderr, "Warning, too many levels (%d) in matrix, using %d only\n",
697 *nlevels, (int)(NMAP*NMAP));
698 *nlevels = NMAP*NMAP;
700 else if (*nlevels < 2)
702 fprintf(stderr, "Warning, too few levels (%d) in matrix, using 2 instead\n",
706 if (!((mid >= lo) && (mid < hi)))
708 gmx_fatal(FARGS, "Lo: %f, Mid: %f, Hi: %f\n", lo, mid, hi);
711 fprintf(out, "static char *gromacs_xpm[] = {\n");
712 fprintf(out, "\"%d %d %d %d\",\n",
713 n_x, n_y, *nlevels, (*nlevels <= NMAP) ? 1 : 2);
715 nmid = calc_nmid(*nlevels, lo, mid, hi);
717 clev_hi = (*nlevels - 1 - nmid);
718 for (i = 0; (i < nmid); i++)
720 r = rlo.r+(i*(rmid.r-rlo.r)/clev_lo);
721 g = rlo.g+(i*(rmid.g-rlo.g)/clev_lo);
722 b = rlo.b+(i*(rmid.b-rlo.b)/clev_lo);
723 fprintf(out, "\"%c%c c #%02X%02X%02X \" /* \"%.3g\" */,\n",
725 (*nlevels <= NMAP) ? ' ' : mapper[i/NMAP],
726 (unsigned int)round(255*r),
727 (unsigned int)round(255*g),
728 (unsigned int)round(255*b),
729 ((nmid - i)*lo + i*mid)/clev_lo);
731 for (i = 0; (i < (*nlevels-nmid)); i++)
733 r = rmid.r+(i*(rhi.r-rmid.r)/clev_hi);
734 g = rmid.g+(i*(rhi.g-rmid.g)/clev_hi);
735 b = rmid.b+(i*(rhi.b-rmid.b)/clev_hi);
736 fprintf(out, "\"%c%c c #%02X%02X%02X \" /* \"%.3g\" */,\n",
737 mapper[(i+nmid) % NMAP],
738 (*nlevels <= NMAP) ? ' ' : mapper[(i+nmid)/NMAP],
739 (unsigned int)round(255*r),
740 (unsigned int)round(255*g),
741 (unsigned int)round(255*b),
742 ((*nlevels - 1 - nmid - i)*mid + i*hi)/clev_hi);
746 static void pr_simple_cmap(FILE *out, real lo, real hi, int nlevel, t_rgb rlo,
752 for (i = 0; (i < nlevel); i++)
754 fac = (i+1.0)/(nlevel);
755 r = rlo.r+fac*(rhi.r-rlo.r);
756 g = rlo.g+fac*(rhi.g-rlo.g);
757 b = rlo.b+fac*(rhi.b-rlo.b);
758 fprintf(out, "\"%c%c c #%02X%02X%02X \" /* \"%.3g\" */,\n",
759 mapper[(i+i0) % NMAP],
760 (nlevel <= NMAP) ? ' ' : mapper[(i+i0)/NMAP],
761 (unsigned int)round(255*r),
762 (unsigned int)round(255*g),
763 (unsigned int)round(255*b),
768 static void pr_discrete_cmap(FILE *out, int *nlevel, int i0)
771 { 1.0, 1.0, 1.0 }, /* white */
772 { 1.0, 0.0, 0.0 }, /* red */
773 { 1.0, 1.0, 0.0 }, /* yellow */
774 { 0.0, 0.0, 1.0 }, /* blue */
775 { 0.0, 1.0, 0.0 }, /* green */
776 { 1.0, 0.0, 1.0 }, /* purple */
777 { 1.0, 0.4, 0.0 }, /* orange */
778 { 0.0, 1.0, 1.0 }, /* cyan */
779 { 1.0, 0.4, 0.4 }, /* pink */
780 { 1.0, 1.0, 0.0 }, /* yellow */
781 { 0.4, 0.4, 1.0 }, /* lightblue */
782 { 0.4, 1.0, 0.4 }, /* lightgreen */
783 { 1.0, 0.4, 1.0 }, /* lightpurple */
784 { 1.0, 0.7, 0.4 }, /* lightorange */
785 { 0.4, 1.0, 1.0 }, /* lightcyan */
786 { 0.0, 0.0, 0.0 } /* black */
791 *nlevel = std::min(16, *nlevel);
793 for (i = 0; (i < n); i++)
795 fprintf(out, "\"%c%c c #%02X%02X%02X \" /* \"%3d\" */,\n",
796 mapper[(i+i0) % NMAP],
797 (n <= NMAP) ? ' ' : mapper[(i+i0)/NMAP],
798 (unsigned int)round(255*rgbd[i].r),
799 (unsigned int)round(255*rgbd[i].g),
800 (unsigned int)round(255*rgbd[i].b),
807 void write_xpm_map_split(FILE *out, int n_x, int n_y,
808 int *nlevel_top, real lo_top, real hi_top,
809 t_rgb rlo_top, t_rgb rhi_top,
810 gmx_bool bDiscreteColor,
811 int *nlevel_bot, real lo_bot, real hi_bot,
812 t_rgb rlo_bot, t_rgb rhi_bot)
816 ntot = *nlevel_top + *nlevel_bot;
819 gmx_fatal(FARGS, "Warning, too many levels (%d) in matrix", ntot);
822 fprintf(out, "static char *gromacs_xpm[] = {\n");
823 fprintf(out, "\"%d %d %d %d\",\n", n_x, n_y, ntot, 1);
827 pr_discrete_cmap(out, nlevel_bot, 0);
831 pr_simple_cmap(out, lo_bot, hi_bot, *nlevel_bot, rlo_bot, rhi_bot, 0);
834 pr_simple_cmap(out, lo_top, hi_top, *nlevel_top, rlo_top, rhi_top, *nlevel_bot);
838 void write_xpm_map(FILE *out, int n_x, int n_y, int *nlevels, real lo, real hi,
839 t_rgb rlo, t_rgb rhi)
842 real invlevel, r, g, b;
844 if (*nlevels > NMAP*NMAP)
846 fprintf(stderr, "Warning, too many levels (%d) in matrix, using %d only\n",
847 *nlevels, (int)(NMAP*NMAP));
848 *nlevels = NMAP*NMAP;
850 else if (*nlevels < 2)
852 fprintf(stderr, "Warning, too few levels (%d) in matrix, using 2 instead\n", *nlevels);
856 fprintf(out, "static char *gromacs_xpm[] = {\n");
857 fprintf(out, "\"%d %d %d %d\",\n",
858 n_x, n_y, *nlevels, (*nlevels <= NMAP) ? 1 : 2);
860 invlevel = 1.0/(*nlevels-1);
861 for (i = 0; (i < *nlevels); i++)
864 r = (nlo*rlo.r+i*rhi.r)*invlevel;
865 g = (nlo*rlo.g+i*rhi.g)*invlevel;
866 b = (nlo*rlo.b+i*rhi.b)*invlevel;
867 fprintf(out, "\"%c%c c #%02X%02X%02X \" /* \"%.3g\" */,\n",
868 mapper[i % NMAP], (*nlevels <= NMAP) ? ' ' : mapper[i/NMAP],
869 (unsigned int)round(255*r),
870 (unsigned int)round(255*g),
871 (unsigned int)round(255*b),
872 (nlo*lo+i*hi)*invlevel);
876 void write_xpm_axis(FILE *out, const char *axis, gmx_bool bSpatial, int n,
883 for (i = 0; i < (bSpatial ? n+1 : n); i++)
889 fprintf(out, "*/\n");
891 fprintf(out, "/* %s-axis: ", axis);
893 fprintf(out, "%g ", label[i]);
895 fprintf(out, "*/\n");
899 void write_xpm_data(FILE *out, int n_x, int n_y, real **matrix,
900 real lo, real hi, int nlevels)
905 invlevel = (nlevels-1)/(hi-lo);
906 for (j = n_y-1; (j >= 0); j--)
908 if (j%(1+n_y/100) == 0)
910 fprintf(stderr, "%3d%%\b\b\b\b", (100*(n_y-j))/n_y);
913 for (i = 0; (i < n_x); i++)
915 c = gmx_nint((matrix[i][j]-lo)*invlevel);
926 fprintf(out, "%c", mapper[c]);
930 fprintf(out, "%c%c", mapper[c % NMAP], mapper[c / NMAP]);
935 fprintf(out, "\",\n");
939 fprintf(out, "\"\n");
944 void write_xpm_data3(FILE *out, int n_x, int n_y, real **matrix,
945 real lo, real mid, real hi, int nlevels)
947 int i, j, c = 0, nmid;
948 real invlev_lo, invlev_hi;
950 nmid = calc_nmid(nlevels, lo, mid, hi);
951 invlev_hi = (nlevels-1-nmid)/(hi-mid);
952 invlev_lo = (nmid)/(mid-lo);
954 for (j = n_y-1; (j >= 0); j--)
956 if (j%(1+n_y/100) == 0)
958 fprintf(stderr, "%3d%%\b\b\b\b", (100*(n_y-j))/n_y);
961 for (i = 0; (i < n_x); i++)
963 if (matrix[i][j] >= mid)
965 c = nmid+gmx_nint((matrix[i][j]-mid)*invlev_hi);
967 else if (matrix[i][j] >= lo)
969 c = gmx_nint((matrix[i][j]-lo)*invlev_lo);
986 fprintf(out, "%c", mapper[c]);
990 fprintf(out, "%c%c", mapper[c % NMAP], mapper[c / NMAP]);
995 fprintf(out, "\",\n");
999 fprintf(out, "\"\n");
1004 void write_xpm_data_split(FILE *out, int n_x, int n_y, real **matrix,
1005 real lo_top, real hi_top, int nlevel_top,
1006 real lo_bot, real hi_bot, int nlevel_bot)
1009 real invlev_top, invlev_bot;
1011 invlev_top = (nlevel_top-1)/(hi_top-lo_top);
1012 invlev_bot = (nlevel_bot-1)/(hi_bot-lo_bot);
1014 for (j = n_y-1; (j >= 0); j--)
1016 if (j % (1+n_y/100) == 0)
1018 fprintf(stderr, "%3d%%\b\b\b\b", (100*(n_y-j))/n_y);
1021 for (i = 0; (i < n_x); i++)
1025 c = nlevel_bot+round((matrix[i][j]-lo_top)*invlev_top);
1026 if ((c < nlevel_bot) || (c >= nlevel_bot+nlevel_top))
1028 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]);
1033 c = round((matrix[i][j]-lo_bot)*invlev_bot);
1034 if ((c < 0) || (c >= nlevel_bot+nlevel_bot))
1036 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]);
1044 fprintf(out, "%c", mapper[c]);
1048 fprintf(out, "\",\n");
1052 fprintf(out, "\"\n");
1057 void write_xpm_m(FILE *out, t_matrix m)
1059 /* Writes a t_matrix struct to .xpm file */
1065 bOneChar = (m.map[0].code.c2 == 0);
1066 write_xpm_header(out, m.title, m.legend, m.label_x, m.label_y,
1068 fprintf(out, "static char *gromacs_xpm[] = {\n");
1069 fprintf(out, "\"%d %d %d %d\",\n", m.nx, m.ny, m.nmap, bOneChar ? 1 : 2);
1070 for (i = 0; (i < m.nmap); i++)
1072 fprintf(out, "\"%c%c c #%02X%02X%02X \" /* \"%s\" */,\n",
1074 bOneChar ? ' ' : m.map[i].code.c2,
1075 (unsigned int)round(m.map[i].rgb.r*255),
1076 (unsigned int)round(m.map[i].rgb.g*255),
1077 (unsigned int)round(m.map[i].rgb.b*255), m.map[i].desc);
1079 write_xpm_axis(out, "x", m.flags & MAT_SPATIAL_X, m.nx, m.axis_x);
1080 write_xpm_axis(out, "y", m.flags & MAT_SPATIAL_Y, m.ny, m.axis_y);
1081 for (j = m.ny-1; (j >= 0); j--)
1083 if (j%(1+m.ny/100) == 0)
1085 fprintf(stderr, "%3d%%\b\b\b\b", (100*(m.ny-j))/m.ny);
1090 for (i = 0; (i < m.nx); i++)
1092 fprintf(out, "%c", m.map[m.matrix[i][j]].code.c1);
1097 for (i = 0; (i < m.nx); i++)
1099 c = m.map[m.matrix[i][j]].code;
1100 fprintf(out, "%c%c", c.c1, c.c2);
1105 fprintf(out, "\",\n");
1109 fprintf(out, "\"\n");
1114 void write_xpm3(FILE *out, unsigned int flags,
1115 const char *title, const char *legend,
1116 const char *label_x, const char *label_y,
1117 int n_x, int n_y, real axis_x[], real axis_y[],
1118 real *matrix[], real lo, real mid, real hi,
1119 t_rgb rlo, t_rgb rmid, t_rgb rhi, int *nlevels)
1122 * Writes a colormap varying as rlo -> rmid -> rhi.
1127 gmx_fatal(FARGS, "hi (%g) <= lo (%g)", hi, lo);
1130 write_xpm_header(out, title, legend, label_x, label_y, FALSE);
1131 write_xpm_map3(out, n_x, n_y, nlevels, lo, mid, hi, rlo, rmid, rhi);
1132 write_xpm_axis(out, "x", flags & MAT_SPATIAL_X, n_x, axis_x);
1133 write_xpm_axis(out, "y", flags & MAT_SPATIAL_Y, n_y, axis_y);
1134 write_xpm_data3(out, n_x, n_y, matrix, lo, mid, hi, *nlevels);
1137 void write_xpm_split(FILE *out, unsigned int flags,
1138 const char *title, const char *legend,
1139 const char *label_x, const char *label_y,
1140 int n_x, int n_y, real axis_x[], real axis_y[],
1142 real lo_top, real hi_top, int *nlevel_top,
1143 t_rgb rlo_top, t_rgb rhi_top,
1144 real lo_bot, real hi_bot, int *nlevel_bot,
1145 gmx_bool bDiscreteColor,
1146 t_rgb rlo_bot, t_rgb rhi_bot)
1149 * Writes a colormap varying as rlo -> rmid -> rhi.
1152 if (hi_top <= lo_top)
1154 gmx_fatal(FARGS, "hi_top (%g) <= lo_top (%g)", hi_top, lo_top);
1156 if (hi_bot <= lo_bot)
1158 gmx_fatal(FARGS, "hi_bot (%g) <= lo_bot (%g)", hi_bot, lo_bot);
1160 if (bDiscreteColor && (*nlevel_bot >= 16))
1162 gmx_impl("Can not plot more than 16 discrete colors");
1165 write_xpm_header(out, title, legend, label_x, label_y, FALSE);
1166 write_xpm_map_split(out, n_x, n_y, nlevel_top, lo_top, hi_top, rlo_top, rhi_top,
1167 bDiscreteColor, nlevel_bot, lo_bot, hi_bot, rlo_bot, rhi_bot);
1168 write_xpm_axis(out, "x", flags & MAT_SPATIAL_X, n_x, axis_x);
1169 write_xpm_axis(out, "y", flags & MAT_SPATIAL_Y, n_y, axis_y);
1170 write_xpm_data_split(out, n_x, n_y, matrix, lo_top, hi_top, *nlevel_top,
1171 lo_bot, hi_bot, *nlevel_bot);
1174 void write_xpm(FILE *out, unsigned int flags,
1175 const char *title, const char *legend,
1176 const char *label_x, const char *label_y,
1177 int n_x, int n_y, real axis_x[], real axis_y[],
1178 real *matrix[], real lo, real hi,
1179 t_rgb rlo, t_rgb rhi, int *nlevels)
1182 * title matrix title
1183 * legend label for the continuous legend
1184 * label_x label for the x-axis
1185 * label_y label for the y-axis
1186 * n_x, n_y size of the matrix
1187 * axis_x[] the x-ticklabels
1188 * axis_y[] the y-ticklables
1189 * *matrix[] element x,y is matrix[x][y]
1190 * lo output lower than lo is set to lo
1191 * hi output higher than hi is set to hi
1192 * rlo rgb value for level lo
1193 * rhi rgb value for level hi
1194 * nlevels number of color levels for the output
1199 gmx_fatal(FARGS, "hi (%f) <= lo (%f)", hi, lo);
1202 write_xpm_header(out, title, legend, label_x, label_y, FALSE);
1203 write_xpm_map(out, n_x, n_y, nlevels, lo, hi, rlo, rhi);
1204 write_xpm_axis(out, "x", flags & MAT_SPATIAL_X, n_x, axis_x);
1205 write_xpm_axis(out, "y", flags & MAT_SPATIAL_Y, n_y, axis_y);
1206 write_xpm_data(out, n_x, n_y, matrix, lo, hi, *nlevels);