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,2014,2015,2017,2018, 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.
49 #include "gromacs/fileio/gmxfio.h"
50 #include "gromacs/math/utilities.h"
51 #include "gromacs/utility/binaryinformation.h"
52 #include "gromacs/utility/cstringutil.h"
53 #include "gromacs/utility/exceptions.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/futil.h"
56 #include "gromacs/utility/programcontext.h"
57 #include "gromacs/utility/smalloc.h"
59 static const char mapper[] = "ABCDEFGHIJKLMNOPQRSTUVWXYZabcdefghijklmnopqrstuvwxyz0123456789!@#$%^&*()-_=+{}|;:',<.>/?";
60 #define NMAP (long int)strlen(mapper)
62 #define MAX_XPM_LINELENGTH 4096
64 real **mk_matrix(int nx, int ny, gmx_bool b1D)
75 for (i = 0; (i < nx); i++)
90 void done_matrix(int nx, real ***m)
94 for (i = 0; (i < nx); i++)
102 gmx_bool matelmt_cmp(t_xpmelmt e1, t_xpmelmt e2)
104 return (e1.c1 == e2.c1) && (e1.c2 == e2.c2);
107 t_matelmt searchcmap(int n, t_mapping map[], t_xpmelmt c)
111 for (i = 0; (i < n); i++)
113 if (matelmt_cmp(map[i].code, c))
122 int getcmap(FILE *in, const char *fn, t_mapping **map)
126 char code[STRLEN], desc[STRLEN];
130 if (fgets2(line, STRLEN-1, in) == nullptr)
132 gmx_fatal(FARGS, "Not enough lines in colormap file %s"
133 "(just wanted to read number of entries)", fn);
135 sscanf(line, "%d", &n);
137 for (i = 0; (i < n); i++)
139 if (fgets2(line, STRLEN-1, in) == nullptr)
141 gmx_fatal(FARGS, "Not enough lines in colormap file %s"
142 "(should be %d, found only %d)", fn, n+1, i);
144 sscanf(line, "%s%s%lf%lf%lf", code, desc, &r, &g, &b);
145 m[i].code.c1 = code[0];
147 m[i].desc = gmx_strdup(desc);
157 int readcmap(const char *fn, t_mapping **map)
163 n = getcmap(in, fn, map);
169 void printcmap(FILE *out, int n, t_mapping map[])
173 fprintf(out, "%d\n", n);
174 for (i = 0; (i < n); i++)
176 fprintf(out, "%c%c %20s %10g %10g %10g\n",
177 map[i].code.c1 ? map[i].code.c1 : ' ',
178 map[i].code.c2 ? map[i].code.c2 : ' ',
179 map[i].desc, map[i].rgb.r, map[i].rgb.g, map[i].rgb.b);
183 void writecmap(const char *fn, int n, t_mapping map[])
187 out = gmx_fio_fopen(fn, "w");
188 printcmap(out, n, map);
192 static char *fgetline(char **line, int llmax, int *llalloc, FILE *in)
196 if (llmax > *llalloc)
198 srenew(*line, llmax+1);
201 fg = fgets(*line, llmax, in);
207 static void skipstr(char *line)
213 while ((line[c] != ' ') && (line[c] != '\0'))
218 while (line[c] != '\0')
226 static char *line2string(char **line)
230 if (*line != nullptr)
232 while (((*line)[0] != '\"' ) && ( (*line)[0] != '\0' ))
237 if ((*line)[0] != '\"')
244 while (( (*line)[i] != '\"' ) && ( (*line)[i] != '\0' ))
249 if ((*line)[i] != '\"')
262 static void parsestring(char *line, const char *label, char *string)
264 if (strstr(line, label))
266 if (strstr(line, label) < strchr(line, '\"'))
269 strcpy(string, line);
274 static void read_xpm_entry(FILE *in, t_matrix *mm)
277 char *line_buf = nullptr, *line = nullptr, *str, buf[256] = {0};
278 int i, m, col_len, nch = 0, n_axis_x, n_axis_y, llmax;
280 unsigned int r, g, b;
282 gmx_bool bGetOnWithIt, bSetLine;
290 mm->matrix = nullptr;
291 mm->axis_x = nullptr;
292 mm->axis_y = nullptr;
293 mm->bDiscrete = FALSE;
297 while ((nullptr != fgetline(&line_buf, llmax, &llalloc, in)) &&
298 (std::strncmp(line_buf, "static", 6) != 0))
301 parsestring(line, "title", (mm->title));
302 parsestring(line, "legend", (mm->legend));
303 parsestring(line, "x-label", (mm->label_x));
304 parsestring(line, "y-label", (mm->label_y));
305 parsestring(line, "type", buf);
308 if (!line || strncmp(line, "static", 6) != 0)
310 gmx_input("Invalid XPixMap");
313 if (buf[0] && (gmx_strcasecmp(buf, "Discrete") == 0))
315 mm->bDiscrete = TRUE;
320 fprintf(debug, "%s %s %s %s\n",
321 mm->title, mm->legend, mm->label_x, mm->label_y);
325 bGetOnWithIt = FALSE;
326 while (!bGetOnWithIt && (nullptr != fgetline(&line_buf, llmax, &llalloc, in)))
329 while (( line[0] != '\"' ) && ( line[0] != '\0' ))
337 sscanf(line, "%d %d %d %d", &(mm->nx), &(mm->ny), &(mm->nmap), &nch);
340 gmx_fatal(FARGS, "Sorry can only read xpm's with at most 2 caracters per pixel\n");
342 if (mm->nx <= 0 || mm->ny <= 0)
344 gmx_fatal(FARGS, "Dimensions of xpm-file have to be larger than 0\n");
346 llmax = std::max(STRLEN, mm->nx+10);
352 fprintf(debug, "mm->nx %d mm->ny %d mm->nmap %d nch %d\n",
353 mm->nx, mm->ny, mm->nmap, nch);
357 gmx_fatal(FARGS, "Number of characters per pixel not found in xpm\n");
363 while ((m < mm->nmap) && (nullptr != fgetline(&line_buf, llmax, &llalloc, in)))
365 line = std::strchr(line_buf, '\"');
369 /* Read xpm color map entry */
370 map[m].code.c1 = line[0];
377 map[m].code.c2 = line[1];
380 str = std::strchr(line, '#');
385 while (std::isxdigit(str[col_len]))
391 sscanf(line, "%*s #%2x%2x%2x", &r, &g, &b);
392 map[m].rgb.r = r/255.0;
393 map[m].rgb.g = g/255.0;
394 map[m].rgb.b = b/255.0;
396 else if (col_len == 12)
398 sscanf(line, "%*s #%4x%4x%4x", &r, &g, &b);
399 map[m].rgb.r = r/65535.0;
400 map[m].rgb.g = g/65535.0;
401 map[m].rgb.b = b/65535.0;
405 gmx_file("Unsupported or invalid colormap in X PixMap");
410 str = std::strchr(line, 'c');
417 gmx_file("Unsupported or invalid colormap in X PixMap");
419 fprintf(stderr, "Using white for color \"%s", str);
424 line = std::strchr(line, '\"');
427 map[m].desc = gmx_strdup(line);
433 gmx_fatal(FARGS, "Number of read colors map entries (%d) does not match the number in the header (%d)", m, mm->nmap);
437 /* Read axes, if there are any */
448 if (strstr(line, "x-axis"))
450 line = std::strstr(line, "x-axis");
452 if (mm->axis_x == nullptr)
454 snew(mm->axis_x, mm->nx + 1);
456 while (sscanf(line, "%lf", &u) == 1)
458 if (n_axis_x > mm->nx)
460 gmx_fatal(FARGS, "Too many x-axis labels in xpm (max %d)", mm->nx);
462 else if (n_axis_x == mm->nx)
464 mm->flags |= MAT_SPATIAL_X;
466 mm->axis_x[n_axis_x] = u;
471 else if (std::strstr(line, "y-axis"))
473 line = std::strstr(line, "y-axis");
475 if (mm->axis_y == nullptr)
477 snew(mm->axis_y, mm->ny + 1);
479 while (sscanf(line, "%lf", &u) == 1)
481 if (n_axis_y > mm->ny)
483 gmx_file("Too many y-axis labels in xpm");
485 else if (n_axis_y == mm->ny)
487 mm->flags |= MAT_SPATIAL_Y;
489 mm->axis_y[n_axis_y] = u;
495 while ((line[0] != '\"') && (nullptr != fgetline(&line_buf, llmax, &llalloc, in)));
498 snew(mm->matrix, mm->nx);
499 for (i = 0; i < mm->nx; i++)
501 snew(mm->matrix[i], mm->ny);
512 if (m%(1+mm->ny/100) == 0)
514 fprintf(stderr, "%3d%%\b\b\b\b", (100*(mm->ny-m))/mm->ny);
516 while ((line[0] != '\"') && (line[0] != '\0'))
522 gmx_fatal(FARGS, "Not enough caracters in row %d of the matrix\n", m+1);
527 for (i = 0; i < mm->nx; i++)
536 c.c2 = line[nch*i+1];
538 mm->matrix[i][m] = searchcmap(mm->nmap, mm->map, c);
543 while ((m >= 0) && (nullptr != fgetline(&line_buf, llmax, &llalloc, in)));
546 gmx_incons("Not enough rows in the matrix");
552 int read_xpm_matrix(const char *fnm, t_matrix **mat)
555 char *line = nullptr;
559 in = gmx_fio_fopen(fnm, "r");
562 while (nullptr != fgetline(&line, STRLEN, &llalloc, in))
564 if (std::strstr(line, "/* XPM */"))
566 srenew(*mat, nmat+1);
567 read_xpm_entry(in, &(*mat)[nmat]);
575 gmx_file("Invalid XPixMap");
583 real **matrix2real(t_matrix *in, real **out)
594 for (i = 0; i < nmap; i++)
596 if ((map[i].desc == nullptr) || (sscanf(map[i].desc, "%lf", &tmp) != 1))
598 fprintf(stderr, "Could not convert matrix to reals,\n"
599 "color map entry %d has a non-real description: \"%s\"\n",
610 for (i = 0; i < in->nx; i++)
612 snew(out[i], in->ny);
615 for (i = 0; i < in->nx; i++)
617 for (j = 0; j < in->ny; j++)
619 out[i][j] = rmap[in->matrix[i][j]];
625 fprintf(stderr, "Converted a %dx%d matrix with %d levels to reals\n",
626 in->nx, in->ny, nmap);
631 static void write_xpm_header(FILE *out,
632 const char *title, const char *legend,
633 const char *label_x, const char *label_y,
636 fprintf(out, "/* XPM */\n");
637 fprintf(out, "/* This file can be converted to EPS by the GROMACS program xpm2ps */\n");
638 fprintf(out, "/* title: \"%s\" */\n", title);
639 fprintf(out, "/* legend: \"%s\" */\n", legend);
640 fprintf(out, "/* x-label: \"%s\" */\n", label_x);
641 fprintf(out, "/* y-label: \"%s\" */\n", label_y);
644 fprintf(out, "/* type: \"Discrete\" */\n");
648 fprintf(out, "/* type: \"Continuous\" */\n");
652 static int calc_nmid(int nlevels, real lo, real mid, real hi)
654 /* Take care that we have at least 1 entry in the mid to hi range
656 return std::min(std::max(0, static_cast<int>(((mid-lo)/(hi-lo))*(nlevels-1))),
660 static void write_xpm_map3(FILE *out, int n_x, int n_y, int *nlevels,
661 real lo, real mid, real hi,
662 t_rgb rlo, t_rgb rmid, t_rgb rhi)
665 double r, g, b, clev_lo, clev_hi;
667 if (*nlevels > NMAP*NMAP)
669 fprintf(stderr, "Warning, too many levels (%d) in matrix, using %d only\n",
670 *nlevels, static_cast<int>(NMAP*NMAP));
671 *nlevels = NMAP*NMAP;
673 else if (*nlevels < 2)
675 fprintf(stderr, "Warning, too few levels (%d) in matrix, using 2 instead\n",
679 if (!((mid >= lo) && (mid < hi)))
681 gmx_fatal(FARGS, "Lo: %f, Mid: %f, Hi: %f\n", lo, mid, hi);
684 fprintf(out, "static char *gromacs_xpm[] = {\n");
685 fprintf(out, "\"%d %d %d %d\",\n",
686 n_x, n_y, *nlevels, (*nlevels <= NMAP) ? 1 : 2);
688 nmid = calc_nmid(*nlevels, lo, mid, hi);
690 clev_hi = (*nlevels - 1 - nmid);
691 for (i = 0; (i < nmid); i++)
693 r = rlo.r+(i*(rmid.r-rlo.r)/clev_lo);
694 g = rlo.g+(i*(rmid.g-rlo.g)/clev_lo);
695 b = rlo.b+(i*(rmid.b-rlo.b)/clev_lo);
696 fprintf(out, "\"%c%c c #%02X%02X%02X \" /* \"%.3g\" */,\n",
698 (*nlevels <= NMAP) ? ' ' : mapper[i/NMAP],
699 static_cast<unsigned int>(std::round(255*r)),
700 static_cast<unsigned int>(std::round(255*g)),
701 static_cast<unsigned int>(std::round(255*b)),
702 ((nmid - i)*lo + i*mid)/clev_lo);
704 for (i = 0; (i < (*nlevels-nmid)); i++)
706 r = rmid.r+(i*(rhi.r-rmid.r)/clev_hi);
707 g = rmid.g+(i*(rhi.g-rmid.g)/clev_hi);
708 b = rmid.b+(i*(rhi.b-rmid.b)/clev_hi);
709 fprintf(out, "\"%c%c c #%02X%02X%02X \" /* \"%.3g\" */,\n",
710 mapper[(i+nmid) % NMAP],
711 (*nlevels <= NMAP) ? ' ' : mapper[(i+nmid)/NMAP],
712 static_cast<unsigned int>(std::round(255*r)),
713 static_cast<unsigned int>(std::round(255*g)),
714 static_cast<unsigned int>(std::round(255*b)),
715 ((*nlevels - 1 - nmid - i)*mid + i*hi)/clev_hi);
719 static void pr_simple_cmap(FILE *out, real lo, real hi, int nlevel, t_rgb rlo,
725 for (i = 0; (i < nlevel); i++)
727 fac = (i+1.0)/(nlevel);
728 r = rlo.r+fac*(rhi.r-rlo.r);
729 g = rlo.g+fac*(rhi.g-rlo.g);
730 b = rlo.b+fac*(rhi.b-rlo.b);
731 fprintf(out, "\"%c%c c #%02X%02X%02X \" /* \"%.3g\" */,\n",
732 mapper[(i+i0) % NMAP],
733 (nlevel <= NMAP) ? ' ' : mapper[(i+i0)/NMAP],
734 static_cast<unsigned int>(std::round(255*r)),
735 static_cast<unsigned int>(std::round(255*g)),
736 static_cast<unsigned int>(std::round(255*b)),
741 static void pr_discrete_cmap(FILE *out, int *nlevel, int i0)
744 { 1.0, 1.0, 1.0 }, /* white */
745 { 1.0, 0.0, 0.0 }, /* red */
746 { 1.0, 1.0, 0.0 }, /* yellow */
747 { 0.0, 0.0, 1.0 }, /* blue */
748 { 0.0, 1.0, 0.0 }, /* green */
749 { 1.0, 0.0, 1.0 }, /* purple */
750 { 1.0, 0.4, 0.0 }, /* orange */
751 { 0.0, 1.0, 1.0 }, /* cyan */
752 { 1.0, 0.4, 0.4 }, /* pink */
753 { 1.0, 1.0, 0.0 }, /* yellow */
754 { 0.4, 0.4, 1.0 }, /* lightblue */
755 { 0.4, 1.0, 0.4 }, /* lightgreen */
756 { 1.0, 0.4, 1.0 }, /* lightpurple */
757 { 1.0, 0.7, 0.4 }, /* lightorange */
758 { 0.4, 1.0, 1.0 }, /* lightcyan */
759 { 0.0, 0.0, 0.0 } /* black */
764 *nlevel = std::min(16, *nlevel);
766 for (i = 0; (i < n); i++)
768 fprintf(out, "\"%c%c c #%02X%02X%02X \" /* \"%3d\" */,\n",
769 mapper[(i+i0) % NMAP],
770 (n <= NMAP) ? ' ' : mapper[(i+i0)/NMAP],
771 static_cast<unsigned int>(round(255*rgbd[i].r)),
772 static_cast<unsigned int>(round(255*rgbd[i].g)),
773 static_cast<unsigned int>(round(255*rgbd[i].b)),
780 static void write_xpm_map_split(FILE *out, int n_x, int n_y,
781 const int *nlevel_top, real lo_top, real hi_top,
782 t_rgb rlo_top, t_rgb rhi_top,
783 gmx_bool bDiscreteColor,
784 int *nlevel_bot, real lo_bot, real hi_bot,
785 t_rgb rlo_bot, t_rgb rhi_bot)
789 ntot = *nlevel_top + *nlevel_bot;
792 gmx_fatal(FARGS, "Warning, too many levels (%d) in matrix", ntot);
795 fprintf(out, "static char *gromacs_xpm[] = {\n");
796 fprintf(out, "\"%d %d %d %d\",\n", n_x, n_y, ntot, 1);
800 pr_discrete_cmap(out, nlevel_bot, 0);
804 pr_simple_cmap(out, lo_bot, hi_bot, *nlevel_bot, rlo_bot, rhi_bot, 0);
807 pr_simple_cmap(out, lo_top, hi_top, *nlevel_top, rlo_top, rhi_top, *nlevel_bot);
811 static void write_xpm_map(FILE *out, int n_x, int n_y, int *nlevels,
812 real lo, real hi, t_rgb rlo, t_rgb rhi)
815 real invlevel, r, g, b;
817 if (*nlevels > NMAP*NMAP)
819 fprintf(stderr, "Warning, too many levels (%d) in matrix, using %d only\n",
820 *nlevels, static_cast<int>(NMAP*NMAP));
821 *nlevels = NMAP*NMAP;
823 else if (*nlevels < 2)
825 fprintf(stderr, "Warning, too few levels (%d) in matrix, using 2 instead\n", *nlevels);
829 fprintf(out, "static char *gromacs_xpm[] = {\n");
830 fprintf(out, "\"%d %d %d %d\",\n",
831 n_x, n_y, *nlevels, (*nlevels <= NMAP) ? 1 : 2);
833 invlevel = 1.0/(*nlevels-1);
834 for (i = 0; (i < *nlevels); i++)
837 r = (nlo*rlo.r+i*rhi.r)*invlevel;
838 g = (nlo*rlo.g+i*rhi.g)*invlevel;
839 b = (nlo*rlo.b+i*rhi.b)*invlevel;
840 fprintf(out, "\"%c%c c #%02X%02X%02X \" /* \"%.3g\" */,\n",
841 mapper[i % NMAP], (*nlevels <= NMAP) ? ' ' : mapper[i/NMAP],
842 static_cast<unsigned int>(std::round(255*r)),
843 static_cast<unsigned int>(std::round(255*g)),
844 static_cast<unsigned int>(std::round(255*b)),
845 (nlo*lo+i*hi)*invlevel);
849 static void write_xpm_axis(FILE *out, const char *axis, gmx_bool bSpatial,
856 for (i = 0; i < (bSpatial ? n+1 : n); i++)
862 fprintf(out, "*/\n");
864 fprintf(out, "/* %s-axis: ", axis);
866 fprintf(out, "%g ", label[i]);
868 fprintf(out, "*/\n");
872 static void write_xpm_data(FILE *out, int n_x, int n_y, real **mat,
873 real lo, real hi, int nlevels)
878 invlevel = (nlevels-1)/(hi-lo);
879 for (j = n_y-1; (j >= 0); j--)
881 if (j%(1+n_y/100) == 0)
883 fprintf(stderr, "%3d%%\b\b\b\b", (100*(n_y-j))/n_y);
886 for (i = 0; (i < n_x); i++)
888 c = std::round((mat[i][j]-lo)*invlevel);
899 fprintf(out, "%c", mapper[c]);
903 fprintf(out, "%c%c", mapper[c % NMAP], mapper[c / NMAP]);
908 fprintf(out, "\",\n");
912 fprintf(out, "\"\n");
917 static void write_xpm_data3(FILE *out, int n_x, int n_y, real **mat,
918 real lo, real mid, real hi, int nlevels)
920 int i, j, c = 0, nmid;
921 real invlev_lo, invlev_hi;
923 nmid = calc_nmid(nlevels, lo, mid, hi);
924 invlev_hi = (nlevels-1-nmid)/(hi-mid);
925 invlev_lo = (nmid)/(mid-lo);
927 for (j = n_y-1; (j >= 0); j--)
929 if (j%(1+n_y/100) == 0)
931 fprintf(stderr, "%3d%%\b\b\b\b", (100*(n_y-j))/n_y);
934 for (i = 0; (i < n_x); i++)
936 if (mat[i][j] >= mid)
938 c = nmid+std::round((mat[i][j]-mid)*invlev_hi);
940 else if (mat[i][j] >= lo)
942 c = std::round((mat[i][j]-lo)*invlev_lo);
959 fprintf(out, "%c", mapper[c]);
963 fprintf(out, "%c%c", mapper[c % NMAP], mapper[c / NMAP]);
968 fprintf(out, "\",\n");
972 fprintf(out, "\"\n");
977 static void write_xpm_data_split(FILE *out, int n_x, int n_y, real **mat,
978 real lo_top, real hi_top, int nlevel_top,
979 real lo_bot, real hi_bot, int nlevel_bot)
982 real invlev_top, invlev_bot;
984 invlev_top = (nlevel_top-1)/(hi_top-lo_top);
985 invlev_bot = (nlevel_bot-1)/(hi_bot-lo_bot);
987 for (j = n_y-1; (j >= 0); j--)
989 if (j % (1+n_y/100) == 0)
991 fprintf(stderr, "%3d%%\b\b\b\b", (100*(n_y-j))/n_y);
994 for (i = 0; (i < n_x); i++)
998 c = nlevel_bot+std::round((mat[i][j]-lo_top)*invlev_top);
999 if ((c < nlevel_bot) || (c >= nlevel_bot+nlevel_top))
1001 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, mat[i][j]);
1006 c = std::round((mat[i][j]-lo_bot)*invlev_bot);
1007 if ((c < 0) || (c >= nlevel_bot+nlevel_bot))
1009 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, mat[i][j]);
1017 fprintf(out, "%c", mapper[c]);
1021 fprintf(out, "\",\n");
1025 fprintf(out, "\"\n");
1030 void write_xpm_m(FILE *out, t_matrix m)
1032 /* Writes a t_matrix struct to .xpm file */
1038 bOneChar = (m.map[0].code.c2 == 0);
1039 write_xpm_header(out, m.title, m.legend, m.label_x, m.label_y,
1041 fprintf(out, "static char *gromacs_xpm[] = {\n");
1042 fprintf(out, "\"%d %d %d %d\",\n", m.nx, m.ny, m.nmap, bOneChar ? 1 : 2);
1043 for (i = 0; (i < m.nmap); i++)
1045 fprintf(out, "\"%c%c c #%02X%02X%02X \" /* \"%s\" */,\n",
1047 bOneChar ? ' ' : m.map[i].code.c2,
1048 static_cast<unsigned int>(round(m.map[i].rgb.r*255)),
1049 static_cast<unsigned int>(round(m.map[i].rgb.g*255)),
1050 static_cast<unsigned int>(round(m.map[i].rgb.b*255)), m.map[i].desc);
1052 write_xpm_axis(out, "x", m.flags & MAT_SPATIAL_X, m.nx, m.axis_x);
1053 write_xpm_axis(out, "y", m.flags & MAT_SPATIAL_Y, m.ny, m.axis_y);
1054 for (j = m.ny-1; (j >= 0); j--)
1056 if (j%(1+m.ny/100) == 0)
1058 fprintf(stderr, "%3d%%\b\b\b\b", (100*(m.ny-j))/m.ny);
1063 for (i = 0; (i < m.nx); i++)
1065 fprintf(out, "%c", m.map[m.matrix[i][j]].code.c1);
1070 for (i = 0; (i < m.nx); i++)
1072 c = m.map[m.matrix[i][j]].code;
1073 fprintf(out, "%c%c", c.c1, c.c2);
1078 fprintf(out, "\",\n");
1082 fprintf(out, "\"\n");
1087 void write_xpm3(FILE *out, unsigned int flags,
1088 const std::string &title, const std::string &legend,
1089 const std::string &label_x, const std::string &label_y,
1090 int n_x, int n_y, real axis_x[], real axis_y[],
1091 real *mat[], real lo, real mid, real hi,
1092 t_rgb rlo, t_rgb rmid, t_rgb rhi, int *nlevels)
1095 * Writes a colormap varying as rlo -> rmid -> rhi.
1100 gmx_fatal(FARGS, "hi (%g) <= lo (%g)", hi, lo);
1103 write_xpm_header(out, title.c_str(), legend.c_str(), label_x.c_str(), label_y.c_str(), FALSE);
1104 write_xpm_map3(out, n_x, n_y, nlevels, lo, mid, hi, rlo, rmid, rhi);
1105 write_xpm_axis(out, "x", flags & MAT_SPATIAL_X, n_x, axis_x);
1106 write_xpm_axis(out, "y", flags & MAT_SPATIAL_Y, n_y, axis_y);
1107 write_xpm_data3(out, n_x, n_y, mat, lo, mid, hi, *nlevels);
1110 void write_xpm_split(FILE *out, unsigned int flags,
1111 const std::string &title, const std::string &legend,
1112 const std::string &label_x, const std::string &label_y,
1113 int n_x, int n_y, real axis_x[], real axis_y[],
1115 real lo_top, real hi_top, int *nlevel_top,
1116 t_rgb rlo_top, t_rgb rhi_top,
1117 real lo_bot, real hi_bot, int *nlevel_bot,
1118 gmx_bool bDiscreteColor,
1119 t_rgb rlo_bot, t_rgb rhi_bot)
1122 * Writes a colormap varying as rlo -> rmid -> rhi.
1125 if (hi_top <= lo_top)
1127 gmx_fatal(FARGS, "hi_top (%g) <= lo_top (%g)", hi_top, lo_top);
1129 if (hi_bot <= lo_bot)
1131 gmx_fatal(FARGS, "hi_bot (%g) <= lo_bot (%g)", hi_bot, lo_bot);
1133 if (bDiscreteColor && (*nlevel_bot >= 16))
1135 gmx_impl("Can not plot more than 16 discrete colors");
1138 write_xpm_header(out, title.c_str(), legend.c_str(), label_x.c_str(), label_y.c_str(), FALSE);
1139 write_xpm_map_split(out, n_x, n_y, nlevel_top, lo_top, hi_top, rlo_top, rhi_top,
1140 bDiscreteColor, nlevel_bot, lo_bot, hi_bot, rlo_bot, rhi_bot);
1141 write_xpm_axis(out, "x", flags & MAT_SPATIAL_X, n_x, axis_x);
1142 write_xpm_axis(out, "y", flags & MAT_SPATIAL_Y, n_y, axis_y);
1143 write_xpm_data_split(out, n_x, n_y, mat, lo_top, hi_top, *nlevel_top,
1144 lo_bot, hi_bot, *nlevel_bot);
1147 void write_xpm(FILE *out, unsigned int flags,
1148 const std::string &title, const std::string &legend,
1149 const std::string &label_x, const std::string &label_y,
1150 int n_x, int n_y, real axis_x[], real axis_y[],
1151 real *mat[], real lo, real hi,
1152 t_rgb rlo, t_rgb rhi, int *nlevels)
1155 * title matrix title
1156 * legend label for the continuous legend
1157 * label_x label for the x-axis
1158 * label_y label for the y-axis
1159 * n_x, n_y size of the matrix
1160 * axis_x[] the x-ticklabels
1161 * axis_y[] the y-ticklables
1162 * *matrix[] element x,y is matrix[x][y]
1163 * lo output lower than lo is set to lo
1164 * hi output higher than hi is set to hi
1165 * rlo rgb value for level lo
1166 * rhi rgb value for level hi
1167 * nlevels number of color levels for the output
1172 gmx_fatal(FARGS, "hi (%f) <= lo (%f)", hi, lo);
1175 write_xpm_header(out, title.c_str(), legend.c_str(), label_x.c_str(), label_y.c_str(), FALSE);
1176 write_xpm_map(out, n_x, n_y, nlevels, lo, hi, rlo, rhi);
1177 write_xpm_axis(out, "x", flags & MAT_SPATIAL_X, n_x, axis_x);
1178 write_xpm_axis(out, "y", flags & MAT_SPATIAL_Y, n_y, axis_y);
1179 write_xpm_data(out, n_x, n_y, mat, lo, hi, *nlevels);