924eda09b184b05a8a8bb7f822ea7e948d9e5cba
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_wheel.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,2019, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include <cmath>
40 #include <cstdio>
41 #include <cstdlib>
42 #include <cstring>
43
44 #include <algorithm>
45
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/fileio/writeps.h"
48 #include "gromacs/gmxana/gmx_ana.h"
49 #include "gromacs/gmxana/gstat.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/utility/arraysize.h"
52 #include "gromacs/utility/cstringutil.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/smalloc.h"
56 #include "gromacs/utility/strdb.h"
57
58 static gmx_bool* bPhobics(int nres, char* resnm[])
59 {
60     int       i, nb;
61     char**    cb;
62     gmx_bool* bb;
63
64     nb = get_lines("phbres.dat", &cb);
65     snew(bb, nres);
66
67     for (i = 0; (i < nres); i++)
68     {
69         if (search_str(nb, cb, resnm[i]) != -1)
70         {
71             bb[i] = TRUE;
72         }
73     }
74     return bb;
75 }
76
77 static void wheel(const char* fn, int nres, char* resnm[], int r0, real rot0, char* title)
78 {
79     const real fontsize  = 16;
80     const real gray      = 0.9;
81     const real fontasp   = 0.6;
82     const real fontwidth = fontsize * fontasp;
83
84     int       i, sl, slen;
85     real      ring, inner, outer;
86     real      xc, yc, box;
87     gmx_bool* bPh;
88     char**    rnms;
89     char      sign;
90
91     inner = 75.0;
92     slen  = 0;
93     snew(rnms, nres);
94     for (i = 0; (i < nres); i++)
95     {
96         snew(rnms[i], 256);
97         sl   = std::strlen(resnm[i]);
98         sign = resnm[i][sl - 1];
99         if ((sign == '+') || (sign == '-'))
100         {
101             resnm[i][sl - 1] = '\0';
102         }
103         sprintf(rnms[i], "%s-%d", resnm[i], i + r0);
104         if ((sign == '+') || (sign == '-'))
105         {
106             sl              = std::strlen(rnms[i]);
107             rnms[i][sl]     = sign;
108             rnms[i][sl + 1] = '\0';
109         }
110
111         slen = std::max(slen, static_cast<int>(std::strlen(rnms[i])));
112     }
113     ring  = (2 + slen) * fontwidth;
114     outer = inner + ring;
115     box   = inner * 1.5 + (1 + int{ nres / 18 }) * ring;
116
117     bPh = bPhobics(nres, resnm);
118
119     t_psdata out = ps_open(fn, 0, 0, 2.0 * box, 2.0 * box);
120     xc           = box;
121     yc           = box;
122
123     ps_font(&out, efontHELV, 1.5 * fontsize);
124     ps_translate(&out, xc, yc);
125     if (title)
126     {
127         ps_ctext(&out, 0, -fontsize * 1.5 / 2.0, title, eXCenter);
128     }
129     ps_font(&out, efontHELV, fontsize);
130     ps_rotate(&out, rot0);
131     for (i = 0; (i < nres);)
132     {
133         if (bPh[i])
134         {
135             ps_color(&out, gray, gray, gray);
136             ps_fillarcslice(&out, 0, 0, inner, outer, -10, 10);
137             ps_color(&out, 0, 0, 0);
138         }
139         ps_arcslice(&out, 0, 0, inner, outer, -10, 10);
140
141         ps_ctext(&out, inner + fontwidth, -fontsize / 2.0, rnms[i], eXLeft);
142         ps_rotate(&out, -100);
143         i++;
144
145         if ((i % 18) == 0)
146         {
147             inner = outer;
148             outer += ring;
149         }
150     }
151     ps_close(&out);
152 }
153
154 static void wheel2(const char* fn, int nres, char* resnm[], real rot0, char* title)
155 {
156     const real fontsize  = 14;
157     const real gray      = 0.9;
158     const real fontasp   = 0.45;
159     const int  angle     = 9;
160     const real fontwidth = fontsize * fontasp;
161
162     int  i, slen;
163     real ring, inner, outer;
164     real xc, yc, box;
165
166     inner = 60.0;
167     slen  = 0;
168     for (i = 0; (i < nres); i++)
169     {
170         slen = std::max(slen, static_cast<int>(strlen(resnm[i])));
171     }
172     fprintf(stderr, "slen = %d\n", slen);
173     ring  = slen * fontwidth;
174     outer = inner + ring;
175     box   = (1 + gmx::exactDiv(nres, 2 * angle)) * outer;
176
177     t_psdata out = ps_open(fn, 0, 0, 2.0 * box, 2.0 * box);
178     xc           = box;
179     yc           = box;
180
181     ps_font(&out, efontHELV, 1.5 * fontsize);
182     ps_translate(&out, xc, yc);
183     ps_color(&out, 0, 0, 0);
184     if (title)
185     {
186         ps_ctext(&out, 0, -fontsize * 1.5 / 2.0, title, eXCenter);
187     }
188     ps_font(&out, efontHELV, fontsize);
189
190     ps_rotate(&out, rot0);
191     for (i = 0; (i < nres);)
192     {
193         if ((i % 5) == 4)
194         {
195             ps_color(&out, gray, gray, 1.0);
196             ps_fillarcslice(&out, 0, 0, inner, outer, -angle, angle);
197             ps_color(&out, 0, 0, 0);
198         }
199         ps_arcslice(&out, 0, 0, inner, outer, -angle, angle);
200
201         ps_ctext(&out, inner + fontwidth, -fontsize / 2.0, resnm[i], eXLeft);
202         ps_rotate(&out, -2 * angle);
203         i++;
204
205         if ((i % (2 * angle)) == 0)
206         {
207             inner = outer;
208             outer += ring;
209         }
210     }
211     ps_close(&out);
212 }
213
214 int gmx_wheel(int argc, char* argv[])
215 {
216     const char* desc[] = {
217         "[THISMODULE] plots a helical wheel representation of your sequence.",
218         "The input sequence is in the [REF].dat[ref] file where the first line contains",
219         "the number of residues and each consecutive line contains a residue "
220         "name."
221     };
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)
242
243     int    i, nres;
244     char** resnm;
245
246     if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
247     {
248         return 0;
249     }
250
251     for (i = 1; (i < argc); i++)
252     {
253         if (std::strcmp(argv[i], "-r0") == 0)
254         {
255             r0 = std::strtol(argv[++i], nullptr, 10);
256             fprintf(stderr, "First residue is %d\n", r0);
257         }
258         else if (std::strcmp(argv[i], "-rot0") == 0)
259         {
260             rot0 = strtod(argv[++i], nullptr);
261             fprintf(stderr, "Initial rotation is %g\n", rot0);
262         }
263         else if (std::strcmp(argv[i], "-T") == 0)
264         {
265             title = gmx_strdup(argv[++i]);
266             fprintf(stderr, "Title will be '%s'\n", title);
267         }
268         else if (std::strcmp(argv[i], "-nn") == 0)
269         {
270             bNum = FALSE;
271             fprintf(stderr, "No residue numbers\n");
272         }
273         else
274         {
275             gmx_fatal(FARGS, "Incorrect usage of option %s", argv[i]);
276         }
277     }
278
279     nres = get_lines(ftp2fn(efDAT, NFILE, fnm), &resnm);
280     if (bNum)
281     {
282         wheel(ftp2fn(efEPS, NFILE, fnm), nres, resnm, r0, rot0, title);
283     }
284     else
285     {
286         wheel2(ftp2fn(efEPS, NFILE, fnm), nres, resnm, rot0, title);
287     }
288
289     return 0;
290 }