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