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.
7 * Copyright (c) 2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
51 #include "gromacs/fileio/gmxfio.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/utilities.h"
54 #include "gromacs/utility/binaryinformation.h"
55 #include "gromacs/utility/cstringutil.h"
56 #include "gromacs/utility/exceptions.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/utility/programcontext.h"
60 #include "gromacs/utility/smalloc.h"
64 static const char mapper
[] =
65 "ABCDEFGHIJKLMNOPQRSTUVWXYZabcdefghijklmnopqrstuvwxyz0123456789!@#$%^&*()-_=+{}|;:',<.>/"
67 #define NMAP static_cast<long int>(sizeof(mapper) / sizeof(mapper[0]))
69 #define MAX_XPM_LINELENGTH 4096
71 real
** mk_matrix(int nx
, int ny
, gmx_bool b1D
)
82 for (i
= 0; (i
< nx
); i
++)
86 m
[i
] = &(m
[0][i
* ny
]);
97 void done_matrix(int nx
, real
*** m
)
101 for (i
= 0; (i
< nx
); i
++)
109 static bool operator==(t_xpmelmt e1
, t_xpmelmt e2
)
111 return (e1
.c1
== e2
.c1
) && (e1
.c2
== e2
.c2
);
114 //! Return the index of the first element that matches \c c, or -1 if not found.
115 t_matelmt
searchcmap(ArrayRef
<const t_mapping
> map
, t_xpmelmt c
)
117 auto findIt
= std::find_if(map
.begin(), map
.end(), [&c
](const auto& m
) { return (m
.code
== c
); });
118 return (findIt
== map
.end()) ? -1 : (findIt
- map
.begin());
121 //! Read the mapping table from in, return number of entries
122 static std::vector
<t_mapping
> getcmap(FILE* in
, const char* fn
)
126 char code
[STRLEN
], desc
[STRLEN
];
128 std::vector
<t_mapping
> m
;
130 if (fgets2(line
, STRLEN
- 1, in
) == nullptr)
133 "Not enough lines in colormap file %s"
134 "(just wanted to read number of entries)",
137 sscanf(line
, "%d", &n
);
139 for (i
= 0; (i
< n
); i
++)
141 if (fgets2(line
, STRLEN
- 1, in
) == nullptr)
144 "Not enough lines in colormap file %s"
145 "(should be %d, found only %d)",
148 sscanf(line
, "%s%s%lf%lf%lf", code
, desc
, &r
, &g
, &b
);
149 m
[i
].code
.c1
= code
[0];
151 m
[i
].desc
= gmx_strdup(desc
);
160 std::vector
<t_mapping
> readcmap(const char* fn
)
162 FilePtr in
= openLibraryFile(fn
);
163 return getcmap(in
.get(), fn
);
166 void printcmap(FILE* out
, int n
, t_mapping map
[])
170 fprintf(out
, "%d\n", n
);
171 for (i
= 0; (i
< n
); i
++)
173 fprintf(out
, "%c%c %20s %10g %10g %10g\n", map
[i
].code
.c1
? map
[i
].code
.c1
: ' ',
174 map
[i
].code
.c2
? map
[i
].code
.c2
: ' ', map
[i
].desc
, map
[i
].rgb
.r
, map
[i
].rgb
.g
,
179 void writecmap(const char* fn
, int n
, t_mapping map
[])
183 out
= gmx_fio_fopen(fn
, "w");
184 printcmap(out
, n
, map
);
188 static char* fgetline(char** line
, int llmax
, int* llalloc
, FILE* in
)
192 if (llmax
> *llalloc
)
194 srenew(*line
, llmax
+ 1);
197 fg
= fgets(*line
, llmax
, in
);
203 static void skipstr(char* line
)
209 while ((line
[c
] != ' ') && (line
[c
] != '\0'))
214 while (line
[c
] != '\0')
216 line
[c
- i
] = line
[c
];
222 static char* line2string(char** line
)
226 if (*line
!= nullptr)
228 while (((*line
)[0] != '\"') && ((*line
)[0] != '\0'))
233 if ((*line
)[0] != '\"')
240 while (((*line
)[i
] != '\"') && ((*line
)[i
] != '\0'))
245 if ((*line
)[i
] != '\"')
258 //! If a label named \c label is found in \c line, return it. Otherwise return empty string.
259 static std::string
findLabelInLine(const std::string
& line
, const std::string
& label
)
261 std::regex
re(".*" + label
+ "\"(.*)\"");
263 if (std::regex_search(line
, match
, re
) && match
.size() > 1)
267 return std::string();
270 //! Read and return a matrix from \c in
271 static t_matrix
read_xpm_entry(FILE* in
)
273 char * line_buf
= nullptr, *line
= nullptr, *str
, buf
[256] = { 0 };
274 int i
, m
, col_len
, nch
= 0, llmax
;
276 unsigned int r
, g
, b
;
278 gmx_bool bGetOnWithIt
, bSetLine
;
285 while ((nullptr != fgetline(&line_buf
, llmax
, &llalloc
, in
))
286 && (std::strncmp(line_buf
, "static", 6) != 0))
288 std::string lineString
= line_buf
;
289 mm
.title
= findLabelInLine(lineString
, "title");
290 mm
.legend
= findLabelInLine(lineString
, "legend");
291 mm
.label_x
= findLabelInLine(lineString
, "x-label");
292 mm
.label_y
= findLabelInLine(lineString
, "y-label");
293 findLabelInLine(lineString
, "type"); // discard the returned string
296 if (!line_buf
|| strncmp(line_buf
, "static", 6) != 0)
298 gmx_input("Invalid XPixMap");
301 if (buf
[0] && (gmx_strcasecmp(buf
, "Discrete") == 0))
308 fprintf(debug
, "%s %s %s %s\n", mm
.title
.c_str(), mm
.legend
.c_str(), mm
.label_x
.c_str(),
314 bGetOnWithIt
= FALSE
;
315 while (!bGetOnWithIt
&& (nullptr != fgetline(&line_buf
, llmax
, &llalloc
, in
)))
318 while ((line
[0] != '\"') && (line
[0] != '\0'))
326 sscanf(line
, "%d %d %d %d", &(mm
.nx
), &(mm
.ny
), &nmap
, &nch
);
329 gmx_fatal(FARGS
, "Sorry can only read xpm's with at most 2 characters per pixel\n");
331 if (mm
.nx
<= 0 || mm
.ny
<= 0)
333 gmx_fatal(FARGS
, "Dimensions of xpm-file have to be larger than 0\n");
335 llmax
= std::max(STRLEN
, mm
.nx
+ 10);
341 fprintf(debug
, "mm.nx %d mm.ny %d nmap %d nch %d\n", mm
.nx
, mm
.ny
, nmap
, nch
);
345 gmx_fatal(FARGS
, "Number of characters per pixel not found in xpm\n");
351 while ((m
< nmap
) && (nullptr != fgetline(&line_buf
, llmax
, &llalloc
, in
)))
353 line
= std::strchr(line_buf
, '\"');
357 /* Read xpm color map entry */
358 mm
.map
[m
].code
.c1
= line
[0];
361 mm
.map
[m
].code
.c2
= 0;
365 mm
.map
[m
].code
.c2
= line
[1];
368 str
= std::strchr(line
, '#');
373 while (std::isxdigit(str
[col_len
]))
379 sscanf(line
, "%*s #%2x%2x%2x", &r
, &g
, &b
);
380 mm
.map
[m
].rgb
.r
= r
/ 255.0;
381 mm
.map
[m
].rgb
.g
= g
/ 255.0;
382 mm
.map
[m
].rgb
.b
= b
/ 255.0;
384 else if (col_len
== 12)
386 sscanf(line
, "%*s #%4x%4x%4x", &r
, &g
, &b
);
387 mm
.map
[m
].rgb
.r
= r
/ 65535.0;
388 mm
.map
[m
].rgb
.g
= g
/ 65535.0;
389 mm
.map
[m
].rgb
.b
= b
/ 65535.0;
393 gmx_file("Unsupported or invalid colormap in X PixMap");
398 str
= std::strchr(line
, 'c');
405 gmx_file("Unsupported or invalid colormap in X PixMap");
407 fprintf(stderr
, "Using white for color \"%s", str
);
412 line
= std::strchr(line
, '\"');
415 mm
.map
[m
].desc
= gmx_strdup(line
);
422 "Number of read colors map entries (%d) does not match the number in the header "
427 /* Read axes, if there are any */
436 GMX_RELEASE_ASSERT(line
, "Need to have valid line to parse");
437 if (strstr(line
, "x-axis"))
439 line
= std::strstr(line
, "x-axis");
442 mm
.axis_x
.reserve(mm
.nx
+ 1);
443 while (sscanf(line
, "%lf", &u
) == 1)
445 if (ssize(mm
.axis_x
) > mm
.nx
)
447 gmx_fatal(FARGS
, "Too many x-axis labels in xpm (max %d)", mm
.nx
);
449 else if (ssize(mm
.axis_x
) == mm
.nx
)
451 mm
.flags
|= MAT_SPATIAL_X
;
453 mm
.axis_x
.push_back(u
);
457 else if (std::strstr(line
, "y-axis"))
459 line
= std::strstr(line
, "y-axis");
462 mm
.axis_y
.reserve(mm
.ny
+ 1);
463 while (sscanf(line
, "%lf", &u
) == 1)
465 if (ssize(mm
.axis_y
) > mm
.ny
)
467 gmx_fatal(FARGS
, "Too many y-axis labels in xpm (max %d)", mm
.ny
);
469 else if (ssize(mm
.axis_y
) == mm
.ny
)
471 mm
.flags
|= MAT_SPATIAL_Y
;
473 mm
.axis_y
.push_back(u
);
477 } while ((line
[0] != '\"') && (nullptr != fgetline(&line_buf
, llmax
, &llalloc
, in
)));
480 mm
.matrix
.resize(mm
.nx
, mm
.ny
);
481 int rowIndex
= mm
.ny
- 1;
490 if (rowIndex
% (1 + mm
.ny
/ 100) == 0)
492 fprintf(stderr
, "%3d%%\b\b\b\b", (100 * (mm
.ny
- rowIndex
)) / mm
.ny
);
494 while ((line
[0] != '\"') && (line
[0] != '\0'))
500 gmx_fatal(FARGS
, "Not enough characters in row %d of the matrix\n", rowIndex
+ 1);
505 for (i
= 0; i
< mm
.nx
; i
++)
507 c
.c1
= line
[nch
* i
];
514 c
.c2
= line
[nch
* i
+ 1];
516 mm
.matrix(i
, rowIndex
) = searchcmap(mm
.map
, c
);
520 } while ((rowIndex
>= 0) && (nullptr != fgetline(&line_buf
, llmax
, &llalloc
, in
)));
523 gmx_incons("Not enough rows in the matrix");
530 std::vector
<t_matrix
> read_xpm_matrix(const char* fnm
)
533 char* line
= nullptr;
536 in
= gmx_fio_fopen(fnm
, "r");
538 std::vector
<t_matrix
> mat
;
539 while (nullptr != fgetline(&line
, STRLEN
, &llalloc
, in
))
541 if (std::strstr(line
, "/* XPM */"))
543 mat
.emplace_back(read_xpm_entry(in
));
550 gmx_file("Invalid XPixMap");
558 real
** matrix2real(t_matrix
* in
, real
** out
)
562 std::vector
<real
> rmap(in
->map
.size());
564 for (gmx::index i
= 0; i
!= ssize(in
->map
); ++i
)
566 if ((in
->map
[i
].desc
== nullptr) || (sscanf(in
->map
[i
].desc
, "%lf", &tmp
) != 1))
569 "Could not convert matrix to reals,\n"
570 "color map entry %zd has a non-real description: \"%s\"\n",
580 for (int i
= 0; i
< in
->nx
; i
++)
582 snew(out
[i
], in
->ny
);
585 for (int i
= 0; i
< in
->nx
; i
++)
587 for (int j
= 0; j
< in
->ny
; j
++)
589 out
[i
][j
] = rmap
[in
->matrix(i
, j
)];
593 fprintf(stderr
, "Converted a %dx%d matrix with %zu levels to reals\n", in
->nx
, in
->ny
, in
->map
.size());
598 static void write_xpm_header(FILE* out
,
599 const std::string
& title
,
600 const std::string
& legend
,
601 const std::string
& label_x
,
602 const std::string
& label_y
,
605 fprintf(out
, "/* XPM */\n");
606 fprintf(out
, "/* This file can be converted to EPS by the GROMACS program xpm2ps */\n");
607 fprintf(out
, "/* title: \"%s\" */\n", title
.c_str());
608 fprintf(out
, "/* legend: \"%s\" */\n", legend
.c_str());
609 fprintf(out
, "/* x-label: \"%s\" */\n", label_x
.c_str());
610 fprintf(out
, "/* y-label: \"%s\" */\n", label_y
.c_str());
613 fprintf(out
, "/* type: \"Discrete\" */\n");
617 fprintf(out
, "/* type: \"Continuous\" */\n");
621 static int calc_nmid(int nlevels
, real lo
, real mid
, real hi
)
623 /* Take care that we have at least 1 entry in the mid to hi range
625 return std::min(std::max(0, static_cast<int>(((mid
- lo
) / (hi
- lo
)) * (nlevels
- 1))), nlevels
- 1);
629 write_xpm_map3(FILE* out
, int n_x
, int n_y
, int* nlevels
, real lo
, real mid
, real hi
, t_rgb rlo
, t_rgb rmid
, t_rgb rhi
)
632 double r
, g
, b
, clev_lo
, clev_hi
;
634 if (*nlevels
> NMAP
* NMAP
)
636 fprintf(stderr
, "Warning, too many levels (%d) in matrix, using %d only\n", *nlevels
,
637 static_cast<int>(NMAP
* NMAP
));
638 *nlevels
= NMAP
* NMAP
;
640 else if (*nlevels
< 2)
642 fprintf(stderr
, "Warning, too few levels (%d) in matrix, using 2 instead\n", *nlevels
);
645 if (!((mid
>= lo
) && (mid
< hi
)))
647 gmx_fatal(FARGS
, "Lo: %f, Mid: %f, Hi: %f\n", lo
, mid
, hi
);
650 fprintf(out
, "static char *gromacs_xpm[] = {\n");
651 fprintf(out
, "\"%d %d %d %d\",\n", n_x
, n_y
, *nlevels
, (*nlevels
<= NMAP
) ? 1 : 2);
653 nmid
= calc_nmid(*nlevels
, lo
, mid
, hi
);
655 clev_hi
= (*nlevels
- 1 - nmid
);
656 for (i
= 0; (i
< nmid
); i
++)
658 r
= rlo
.r
+ (i
* (rmid
.r
- rlo
.r
) / clev_lo
);
659 g
= rlo
.g
+ (i
* (rmid
.g
- rlo
.g
) / clev_lo
);
660 b
= rlo
.b
+ (i
* (rmid
.b
- rlo
.b
) / clev_lo
);
661 fprintf(out
, "\"%c%c c #%02X%02X%02X \" /* \"%.3g\" */,\n", mapper
[i
% NMAP
],
662 (*nlevels
<= NMAP
) ? ' ' : mapper
[i
/ NMAP
],
663 static_cast<unsigned int>(std::round(255 * r
)),
664 static_cast<unsigned int>(std::round(255 * g
)),
665 static_cast<unsigned int>(std::round(255 * b
)), ((nmid
- i
) * lo
+ i
* mid
) / clev_lo
);
667 for (i
= 0; (i
< (*nlevels
- nmid
)); i
++)
669 r
= rmid
.r
+ (i
* (rhi
.r
- rmid
.r
) / clev_hi
);
670 g
= rmid
.g
+ (i
* (rhi
.g
- rmid
.g
) / clev_hi
);
671 b
= rmid
.b
+ (i
* (rhi
.b
- rmid
.b
) / clev_hi
);
672 fprintf(out
, "\"%c%c c #%02X%02X%02X \" /* \"%.3g\" */,\n", mapper
[(i
+ nmid
) % NMAP
],
673 (*nlevels
<= NMAP
) ? ' ' : mapper
[(i
+ nmid
) / NMAP
],
674 static_cast<unsigned int>(std::round(255 * r
)),
675 static_cast<unsigned int>(std::round(255 * g
)),
676 static_cast<unsigned int>(std::round(255 * b
)),
677 ((*nlevels
- 1 - nmid
- i
) * mid
+ i
* hi
) / clev_hi
);
681 static void pr_simple_cmap(FILE* out
, real lo
, real hi
, int nlevel
, t_rgb rlo
, t_rgb rhi
, int i0
)
686 for (i
= 0; (i
< nlevel
); i
++)
688 fac
= (i
+ 1.0) / (nlevel
);
689 r
= rlo
.r
+ fac
* (rhi
.r
- rlo
.r
);
690 g
= rlo
.g
+ fac
* (rhi
.g
- rlo
.g
);
691 b
= rlo
.b
+ fac
* (rhi
.b
- rlo
.b
);
692 fprintf(out
, "\"%c%c c #%02X%02X%02X \" /* \"%.3g\" */,\n", mapper
[(i
+ i0
) % NMAP
],
693 (nlevel
<= NMAP
) ? ' ' : mapper
[(i
+ i0
) / NMAP
],
694 static_cast<unsigned int>(std::round(255 * r
)),
695 static_cast<unsigned int>(std::round(255 * g
)),
696 static_cast<unsigned int>(std::round(255 * b
)), lo
+ fac
* (hi
- lo
));
700 static void pr_discrete_cmap(FILE* out
, int* nlevel
, int i0
)
703 { 1.0, 1.0, 1.0 }, /* white */
704 { 1.0, 0.0, 0.0 }, /* red */
705 { 1.0, 1.0, 0.0 }, /* yellow */
706 { 0.0, 0.0, 1.0 }, /* blue */
707 { 0.0, 1.0, 0.0 }, /* green */
708 { 1.0, 0.0, 1.0 }, /* purple */
709 { 1.0, 0.4, 0.0 }, /* orange */
710 { 0.0, 1.0, 1.0 }, /* cyan */
711 { 1.0, 0.4, 0.4 }, /* pink */
712 { 1.0, 1.0, 0.0 }, /* yellow */
713 { 0.4, 0.4, 1.0 }, /* lightblue */
714 { 0.4, 1.0, 0.4 }, /* lightgreen */
715 { 1.0, 0.4, 1.0 }, /* lightpurple */
716 { 1.0, 0.7, 0.4 }, /* lightorange */
717 { 0.4, 1.0, 1.0 }, /* lightcyan */
718 { 0.0, 0.0, 0.0 } /* black */
723 *nlevel
= std::min(16, *nlevel
);
725 for (i
= 0; (i
< n
); i
++)
727 fprintf(out
, "\"%c%c c #%02X%02X%02X \" /* \"%3d\" */,\n", mapper
[(i
+ i0
) % NMAP
],
728 (n
<= NMAP
) ? ' ' : mapper
[(i
+ i0
) / NMAP
],
729 static_cast<unsigned int>(round(255 * rgbd
[i
].r
)),
730 static_cast<unsigned int>(round(255 * rgbd
[i
].g
)),
731 static_cast<unsigned int>(round(255 * rgbd
[i
].b
)), i
);
736 static void write_xpm_map_split(FILE* out
,
739 const int* nlevel_top
,
744 gmx_bool bDiscreteColor
,
753 ntot
= *nlevel_top
+ *nlevel_bot
;
756 gmx_fatal(FARGS
, "Warning, too many levels (%d) in matrix", ntot
);
759 fprintf(out
, "static char *gromacs_xpm[] = {\n");
760 fprintf(out
, "\"%d %d %d %d\",\n", n_x
, n_y
, ntot
, 1);
764 pr_discrete_cmap(out
, nlevel_bot
, 0);
768 pr_simple_cmap(out
, lo_bot
, hi_bot
, *nlevel_bot
, rlo_bot
, rhi_bot
, 0);
771 pr_simple_cmap(out
, lo_top
, hi_top
, *nlevel_top
, rlo_top
, rhi_top
, *nlevel_bot
);
775 static void write_xpm_map(FILE* out
, int n_x
, int n_y
, int* nlevels
, real lo
, real hi
, t_rgb rlo
, t_rgb rhi
)
778 real invlevel
, r
, g
, b
;
780 if (*nlevels
> NMAP
* NMAP
)
782 fprintf(stderr
, "Warning, too many levels (%d) in matrix, using %d only\n", *nlevels
,
783 static_cast<int>(NMAP
* NMAP
));
784 *nlevels
= NMAP
* NMAP
;
786 else if (*nlevels
< 2)
788 fprintf(stderr
, "Warning, too few levels (%d) in matrix, using 2 instead\n", *nlevels
);
792 fprintf(out
, "static char *gromacs_xpm[] = {\n");
793 fprintf(out
, "\"%d %d %d %d\",\n", n_x
, n_y
, *nlevels
, (*nlevels
<= NMAP
) ? 1 : 2);
795 invlevel
= 1.0 / (*nlevels
- 1);
796 for (i
= 0; (i
< *nlevels
); i
++)
798 nlo
= *nlevels
- 1 - i
;
799 r
= (nlo
* rlo
.r
+ i
* rhi
.r
) * invlevel
;
800 g
= (nlo
* rlo
.g
+ i
* rhi
.g
) * invlevel
;
801 b
= (nlo
* rlo
.b
+ i
* rhi
.b
) * invlevel
;
802 fprintf(out
, "\"%c%c c #%02X%02X%02X \" /* \"%.3g\" */,\n", mapper
[i
% NMAP
],
803 (*nlevels
<= NMAP
) ? ' ' : mapper
[i
/ NMAP
],
804 static_cast<unsigned int>(std::round(255 * r
)),
805 static_cast<unsigned int>(std::round(255 * g
)),
806 static_cast<unsigned int>(std::round(255 * b
)), (nlo
* lo
+ i
* hi
) * invlevel
);
810 static void writeXpmAxis(FILE* out
, const char* axis
, ArrayRef
<const real
> label
)
816 for (gmx::index i
= 0; i
!= ssize(label
); ++i
)
822 fprintf(out
, "*/\n");
824 fprintf(out
, "/* %s-axis: ", axis
);
826 fprintf(out
, "%g ", label
[i
]);
828 fprintf(out
, "*/\n");
831 static void write_xpm_data(FILE* out
, int n_x
, int n_y
, real
** mat
, real lo
, real hi
, int nlevels
)
836 invlevel
= (nlevels
- 1) / (hi
- lo
);
837 for (j
= n_y
- 1; (j
>= 0); j
--)
839 if (j
% (1 + n_y
/ 100) == 0)
841 fprintf(stderr
, "%3d%%\b\b\b\b", (100 * (n_y
- j
)) / n_y
);
844 for (i
= 0; (i
< n_x
); i
++)
846 c
= roundToInt((mat
[i
][j
] - lo
) * invlevel
);
857 fprintf(out
, "%c", mapper
[c
]);
861 fprintf(out
, "%c%c", mapper
[c
% NMAP
], mapper
[c
/ NMAP
]);
866 fprintf(out
, "\",\n");
870 fprintf(out
, "\"\n");
875 static void write_xpm_data3(FILE* out
, int n_x
, int n_y
, real
** mat
, real lo
, real mid
, real hi
, int nlevels
)
877 int i
, j
, c
= 0, nmid
;
878 real invlev_lo
, invlev_hi
;
880 nmid
= calc_nmid(nlevels
, lo
, mid
, hi
);
881 invlev_hi
= (nlevels
- 1 - nmid
) / (hi
- mid
);
882 invlev_lo
= (nmid
) / (mid
- lo
);
884 for (j
= n_y
- 1; (j
>= 0); j
--)
886 if (j
% (1 + n_y
/ 100) == 0)
888 fprintf(stderr
, "%3d%%\b\b\b\b", (100 * (n_y
- j
)) / n_y
);
891 for (i
= 0; (i
< n_x
); i
++)
893 if (mat
[i
][j
] >= mid
)
895 c
= nmid
+ roundToInt((mat
[i
][j
] - mid
) * invlev_hi
);
897 else if (mat
[i
][j
] >= lo
)
899 c
= roundToInt((mat
[i
][j
] - lo
) * invlev_lo
);
916 fprintf(out
, "%c", mapper
[c
]);
920 fprintf(out
, "%c%c", mapper
[c
% NMAP
], mapper
[c
/ NMAP
]);
925 fprintf(out
, "\",\n");
929 fprintf(out
, "\"\n");
934 static void write_xpm_data_split(FILE* out
,
946 real invlev_top
, invlev_bot
;
948 invlev_top
= (nlevel_top
- 1) / (hi_top
- lo_top
);
949 invlev_bot
= (nlevel_bot
- 1) / (hi_bot
- lo_bot
);
951 for (j
= n_y
- 1; (j
>= 0); j
--)
953 if (j
% (1 + n_y
/ 100) == 0)
955 fprintf(stderr
, "%3d%%\b\b\b\b", (100 * (n_y
- j
)) / n_y
);
958 for (i
= 0; (i
< n_x
); i
++)
962 c
= nlevel_bot
+ roundToInt((mat
[i
][j
] - lo_top
) * invlev_top
);
963 if ((c
< nlevel_bot
) || (c
>= nlevel_bot
+ nlevel_top
))
966 "Range checking i = %d, j = %d, c = %d, bot = %d, top = %d "
968 i
, j
, c
, nlevel_bot
, nlevel_top
, mat
[i
][j
]);
973 c
= roundToInt((mat
[i
][j
] - lo_bot
) * invlev_bot
);
974 if ((c
< 0) || (c
>= nlevel_bot
+ nlevel_bot
))
977 "Range checking i = %d, j = %d, c = %d, bot = %d, top = %d "
979 i
, j
, c
, nlevel_bot
, nlevel_top
, mat
[i
][j
]);
987 fprintf(out
, "%c", mapper
[c
]);
991 fprintf(out
, "\",\n");
995 fprintf(out
, "\"\n");
1000 void write_xpm_m(FILE* out
, t_matrix m
)
1005 bOneChar
= (m
.map
[0].code
.c2
== 0);
1006 write_xpm_header(out
, m
.title
, m
.legend
, m
.label_x
, m
.label_y
, m
.bDiscrete
);
1007 fprintf(out
, "static char *gromacs_xpm[] = {\n");
1008 fprintf(out
, "\"%d %d %zu %d\",\n", m
.nx
, m
.ny
, m
.map
.size(), bOneChar
? 1 : 2);
1009 for (const auto& map
: m
.map
)
1011 fprintf(out
, "\"%c%c c #%02X%02X%02X \" /* \"%s\" */,\n", map
.code
.c1
,
1012 bOneChar
? ' ' : map
.code
.c2
, static_cast<unsigned int>(round(map
.rgb
.r
* 255)),
1013 static_cast<unsigned int>(round(map
.rgb
.g
* 255)),
1014 static_cast<unsigned int>(round(map
.rgb
.b
* 255)), map
.desc
);
1016 writeXpmAxis(out
, "x", m
.axis_x
);
1017 writeXpmAxis(out
, "y", m
.axis_y
);
1018 for (int j
= m
.ny
- 1; (j
>= 0); j
--)
1020 if (j
% (1 + m
.ny
/ 100) == 0)
1022 fprintf(stderr
, "%3d%%\b\b\b\b", (100 * (m
.ny
- j
)) / m
.ny
);
1027 for (int i
= 0; (i
< m
.nx
); i
++)
1029 fprintf(out
, "%c", m
.map
[m
.matrix(i
, j
)].code
.c1
);
1034 for (int i
= 0; (i
< m
.nx
); i
++)
1036 c
= m
.map
[m
.matrix(i
, j
)].code
;
1037 fprintf(out
, "%c%c", c
.c1
, c
.c2
);
1042 fprintf(out
, "\",\n");
1046 fprintf(out
, "\"\n");
1051 void write_xpm3(FILE* out
,
1053 const std::string
& title
,
1054 const std::string
& legend
,
1055 const std::string
& label_x
,
1056 const std::string
& label_y
,
1071 * Writes a colormap varying as rlo -> rmid -> rhi.
1076 gmx_fatal(FARGS
, "hi (%g) <= lo (%g)", hi
, lo
);
1079 write_xpm_header(out
, title
, legend
, label_x
, label_y
, FALSE
);
1080 write_xpm_map3(out
, n_x
, n_y
, nlevels
, lo
, mid
, hi
, rlo
, rmid
, rhi
);
1081 writeXpmAxis(out
, "x", ArrayRef
<real
>(axis_x
, axis_x
+ n_x
+ ((flags
& MAT_SPATIAL_X
) != 0U ? 1 : 0)));
1082 writeXpmAxis(out
, "y", ArrayRef
<real
>(axis_y
, axis_y
+ n_y
+ ((flags
& MAT_SPATIAL_Y
) != 0U ? 1 : 0)));
1083 write_xpm_data3(out
, n_x
, n_y
, mat
, lo
, mid
, hi
, *nlevels
);
1086 void write_xpm_split(FILE* out
,
1088 const std::string
& title
,
1089 const std::string
& legend
,
1090 const std::string
& label_x
,
1091 const std::string
& label_y
,
1105 gmx_bool bDiscreteColor
,
1110 * Writes a colormap varying as rlo -> rmid -> rhi.
1113 if (hi_top
<= lo_top
)
1115 gmx_fatal(FARGS
, "hi_top (%g) <= lo_top (%g)", hi_top
, lo_top
);
1117 if (hi_bot
<= lo_bot
)
1119 gmx_fatal(FARGS
, "hi_bot (%g) <= lo_bot (%g)", hi_bot
, lo_bot
);
1121 if (bDiscreteColor
&& (*nlevel_bot
>= 16))
1123 gmx_impl("Can not plot more than 16 discrete colors");
1126 write_xpm_header(out
, title
, legend
, label_x
, label_y
, FALSE
);
1127 write_xpm_map_split(out
, n_x
, n_y
, nlevel_top
, lo_top
, hi_top
, rlo_top
, rhi_top
, bDiscreteColor
,
1128 nlevel_bot
, lo_bot
, hi_bot
, rlo_bot
, rhi_bot
);
1129 writeXpmAxis(out
, "x", ArrayRef
<real
>(axis_x
, axis_x
+ n_x
+ ((flags
& MAT_SPATIAL_X
) != 0U ? 1 : 0)));
1130 writeXpmAxis(out
, "y", ArrayRef
<real
>(axis_y
, axis_y
+ n_y
+ ((flags
& MAT_SPATIAL_Y
) != 0U ? 1 : 0)));
1131 write_xpm_data_split(out
, n_x
, n_y
, mat
, lo_top
, hi_top
, *nlevel_top
, lo_bot
, hi_bot
, *nlevel_bot
);
1134 void write_xpm(FILE* out
,
1136 const std::string
& title
,
1137 const std::string
& legend
,
1138 const std::string
& label_x
,
1139 const std::string
& label_y
,
1152 * title matrix title
1153 * legend label for the continuous legend
1154 * label_x label for the x-axis
1155 * label_y label for the y-axis
1156 * n_x, n_y size of the matrix
1157 * axis_x[] the x-ticklabels
1158 * axis_y[] the y-ticklables
1159 * *matrix[] element x,y is matrix[x][y]
1160 * lo output lower than lo is set to lo
1161 * hi output higher than hi is set to hi
1162 * rlo rgb value for level lo
1163 * rhi rgb value for level hi
1164 * nlevels number of color levels for the output
1169 gmx_fatal(FARGS
, "hi (%f) <= lo (%f)", hi
, lo
);
1172 write_xpm_header(out
, title
, legend
, label_x
, label_y
, FALSE
);
1173 write_xpm_map(out
, n_x
, n_y
, nlevels
, lo
, hi
, rlo
, rhi
);
1174 writeXpmAxis(out
, "x", ArrayRef
<real
>(axis_x
, axis_x
+ n_x
+ ((flags
& MAT_SPATIAL_X
) != 0U ? 1 : 0)));
1175 writeXpmAxis(out
, "y", ArrayRef
<real
>(axis_y
, axis_y
+ n_y
+ ((flags
& MAT_SPATIAL_Y
) != 0U ? 1 : 0)));
1176 write_xpm_data(out
, n_x
, n_y
, mat
, lo
, hi
, *nlevels
);