Update instructions in containers.rst
[gromacs.git] / src / gromacs / gmxana / gmx_wheel.cpp
blobdbb5514dd03df56fca5ad6b5d135e02e369dac08
1 /*
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.
38 #include "gmxpre.h"
40 #include <cmath>
41 #include <cstdio>
42 #include <cstdlib>
43 #include <cstring>
45 #include <algorithm>
47 #include "gromacs/commandline/pargs.h"
48 #include "gromacs/fileio/writeps.h"
49 #include "gromacs/gmxana/gmx_ana.h"
50 #include "gromacs/gmxana/gstat.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/utility/arraysize.h"
53 #include "gromacs/utility/cstringutil.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/futil.h"
56 #include "gromacs/utility/smalloc.h"
57 #include "gromacs/utility/strdb.h"
59 static gmx_bool* bPhobics(int nres, char* resnm[])
61 int i, nb;
62 char** cb;
63 gmx_bool* bb;
65 nb = get_lines("phbres.dat", &cb);
66 snew(bb, nres);
68 for (i = 0; (i < nres); i++)
70 if (search_str(nb, cb, resnm[i]) != -1)
72 bb[i] = TRUE;
75 return bb;
78 static void wheel(const char* fn, int nres, char* resnm[], int r0, real rot0, char* title)
80 const real fontsize = 16;
81 const real gray = 0.9;
82 const real fontasp = 0.6;
83 const real fontwidth = fontsize * fontasp;
85 int i, sl, slen;
86 real ring, inner, outer;
87 real xc, yc, box;
88 gmx_bool* bPh;
89 char** rnms;
90 char sign;
92 inner = 75.0;
93 slen = 0;
94 snew(rnms, nres);
95 for (i = 0; (i < nres); i++)
97 snew(rnms[i], 256);
98 sl = std::strlen(resnm[i]);
99 sign = resnm[i][sl - 1];
100 if ((sign == '+') || (sign == '-'))
102 resnm[i][sl - 1] = '\0';
104 sprintf(rnms[i], "%s-%d", resnm[i], i + r0);
105 if ((sign == '+') || (sign == '-'))
107 sl = std::strlen(rnms[i]);
108 rnms[i][sl] = sign;
109 rnms[i][sl + 1] = '\0';
112 slen = std::max(slen, static_cast<int>(std::strlen(rnms[i])));
114 ring = (2 + slen) * fontwidth;
115 outer = inner + ring;
116 box = inner * 1.5 + (1 + int{ nres / 18 }) * ring;
118 bPh = bPhobics(nres, resnm);
120 t_psdata out = ps_open(fn, 0, 0, 2.0 * box, 2.0 * box);
121 xc = box;
122 yc = box;
124 ps_font(&out, efontHELV, 1.5 * fontsize);
125 ps_translate(&out, xc, yc);
126 if (title)
128 ps_ctext(&out, 0, -fontsize * 1.5 / 2.0, title, eXCenter);
130 ps_font(&out, efontHELV, fontsize);
131 ps_rotate(&out, rot0);
132 for (i = 0; (i < nres);)
134 if (bPh[i])
136 ps_color(&out, gray, gray, gray);
137 ps_fillarcslice(&out, 0, 0, inner, outer, -10, 10);
138 ps_color(&out, 0, 0, 0);
140 ps_arcslice(&out, 0, 0, inner, outer, -10, 10);
142 ps_ctext(&out, inner + fontwidth, -fontsize / 2.0, rnms[i], eXLeft);
143 ps_rotate(&out, -100);
144 i++;
146 if ((i % 18) == 0)
148 inner = outer;
149 outer += ring;
152 ps_close(&out);
155 static void wheel2(const char* fn, int nres, char* resnm[], real rot0, char* title)
157 const real fontsize = 14;
158 const real gray = 0.9;
159 const real fontasp = 0.45;
160 const int angle = 9;
161 const real fontwidth = fontsize * fontasp;
163 int i, slen;
164 real ring, inner, outer;
165 real xc, yc, box;
167 inner = 60.0;
168 slen = 0;
169 for (i = 0; (i < nres); i++)
171 slen = std::max(slen, static_cast<int>(strlen(resnm[i])));
173 fprintf(stderr, "slen = %d\n", slen);
174 ring = slen * fontwidth;
175 outer = inner + ring;
176 box = (1 + gmx::exactDiv(nres, 2 * angle)) * outer;
178 t_psdata out = ps_open(fn, 0, 0, 2.0 * box, 2.0 * box);
179 xc = box;
180 yc = box;
182 ps_font(&out, efontHELV, 1.5 * fontsize);
183 ps_translate(&out, xc, yc);
184 ps_color(&out, 0, 0, 0);
185 if (title)
187 ps_ctext(&out, 0, -fontsize * 1.5 / 2.0, title, eXCenter);
189 ps_font(&out, efontHELV, fontsize);
191 ps_rotate(&out, rot0);
192 for (i = 0; (i < nres);)
194 if ((i % 5) == 4)
196 ps_color(&out, gray, gray, 1.0);
197 ps_fillarcslice(&out, 0, 0, inner, outer, -angle, angle);
198 ps_color(&out, 0, 0, 0);
200 ps_arcslice(&out, 0, 0, inner, outer, -angle, angle);
202 ps_ctext(&out, inner + fontwidth, -fontsize / 2.0, resnm[i], eXLeft);
203 ps_rotate(&out, -2 * angle);
204 i++;
206 if ((i % (2 * angle)) == 0)
208 inner = outer;
209 outer += ring;
212 ps_close(&out);
215 int gmx_wheel(int argc, char* argv[])
217 const char* desc[] = {
218 "[THISMODULE] plots a helical wheel representation of your sequence.",
219 "The input sequence is in the [REF].dat[ref] file where the first line contains",
220 "the number of residues and each consecutive line contains a residue name."
222 gmx_output_env_t* oenv;
223 static real rot0 = 0;
224 static gmx_bool bNum = TRUE;
225 static char* title = nullptr;
226 static int r0 = 1;
227 t_pargs pa[] = { { "-r0", FALSE, etINT, { &r0 }, "The first residue number in the sequence" },
228 { "-rot0",
229 FALSE,
230 etREAL,
231 { &rot0 },
232 "Rotate around an angle initially (90 degrees makes sense)" },
233 { "-T",
234 FALSE,
235 etSTR,
236 { &title },
237 "Plot a title in the center of the wheel (must be shorter than 10 "
238 "characters, or it will overwrite the wheel)" },
239 { "-nn", FALSE, etBOOL, { &bNum }, "Toggle numbers" } };
240 t_filenm fnm[] = { { efDAT, "-f", nullptr, ffREAD }, { efEPS, "-o", nullptr, ffWRITE } };
241 #define NFILE asize(fnm)
243 int i, nres;
244 char** resnm;
246 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
248 return 0;
251 for (i = 1; (i < argc); i++)
253 if (std::strcmp(argv[i], "-r0") == 0)
255 r0 = std::strtol(argv[++i], nullptr, 10);
256 fprintf(stderr, "First residue is %d\n", r0);
258 else if (std::strcmp(argv[i], "-rot0") == 0)
260 rot0 = strtod(argv[++i], nullptr);
261 fprintf(stderr, "Initial rotation is %g\n", rot0);
263 else if (std::strcmp(argv[i], "-T") == 0)
265 title = gmx_strdup(argv[++i]);
266 fprintf(stderr, "Title will be '%s'\n", title);
268 else if (std::strcmp(argv[i], "-nn") == 0)
270 bNum = FALSE;
271 fprintf(stderr, "No residue numbers\n");
273 else
275 gmx_fatal(FARGS, "Incorrect usage of option %s", argv[i]);
279 nres = get_lines(ftp2fn(efDAT, NFILE, fnm), &resnm);
280 if (bNum)
282 wheel(ftp2fn(efEPS, NFILE, fnm), nres, resnm, r0, rot0, title);
284 else
286 wheel2(ftp2fn(efEPS, NFILE, fnm), nres, resnm, rot0, title);
289 return 0;