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