Sort all includes in src/gromacs
[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/commandline/pargs.h"
45 #include "gromacs/fileio/strdb.h"
46 #include "gromacs/fileio/writeps.h"
47 #include "gromacs/gmxana/gmx_ana.h"
48 #include "gromacs/gmxana/gstat.h"
49 #include "gromacs/legacyheaders/macros.h"
50 #include "gromacs/legacyheaders/typedefs.h"
51 #include "gromacs/math/vec.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
57 static gmx_bool *bPhobics(int nres, char *resnm[])
58 {
59     int       i, nb;
60     char    **cb;
61     gmx_bool *bb;
62
63     nb = get_strings("phbres.dat", &cb);
64     snew(bb, nres);
65
66     for (i = 0; (i < nres); i++)
67     {
68         if (search_str(nb, cb, resnm[i]) != -1)
69         {
70             bb[i] = TRUE;
71         }
72     }
73     return bb;
74 }
75
76 void wheel(const char *fn, int nres, char *resnm[], int r0, real rot0, char *title)
77 {
78     const real fontsize  = 16;
79     const real gray      = 0.9;
80     const real fontasp   = 0.6;
81     const real fontwidth = fontsize*fontasp;
82
83     t_psdata   out;
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   = 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            = strlen(rnms[i]);
107             rnms[i][sl]   = sign;
108             rnms[i][sl+1] = '\0';
109         }
110
111         slen = max(slen, (int)strlen(rnms[i]));
112     }
113     ring  = (2+slen)*fontwidth;
114     outer = inner+ring;
115     box   = inner*1.5+(1+(nres / 18))*ring;
116
117     bPh = bPhobics(nres, resnm);
118
119     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 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     t_psdata   out;
163     int        i, slen;
164     real       ring, inner, outer;
165     real       xc, yc, box;
166
167     inner = 60.0;
168     slen  = 0;
169     for (i = 0; (i < nres); i++)
170     {
171         slen = max(slen, (int)strlen(resnm[i]));
172     }
173     fprintf(stderr, "slen = %d\n", slen);
174     ring  = (slen)*fontwidth;
175     outer = inner+ring;
176     box   = (1+(nres / (2*angle)))*outer;
177
178     out = ps_open(fn, 0, 0, 2.0*box, 2.0*box);
179     xc  = box;
180     yc  = box;
181
182     ps_font(out, efontHELV, 1.5*fontsize);
183     ps_translate(out, xc, yc);
184     ps_color(out, 0, 0, 0);
185     if (title)
186     {
187         ps_ctext(out, 0, -fontsize*1.5/2.0, title, eXCenter);
188     }
189     ps_font(out, efontHELV, fontsize);
190
191     ps_rotate(out, rot0);
192     for (i = 0; (i < nres); )
193     {
194         if ((i % 5) == 4)
195         {
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);
199         }
200         ps_arcslice(out, 0, 0, inner, outer, -angle, angle);
201
202         ps_ctext(out, inner+fontwidth, -fontsize/2.0, resnm[i], eXLeft);
203         ps_rotate(out, -2*angle);
204         i++;
205
206         if ((i % (2*angle)) == 0)
207         {
208             inner  = outer;
209             outer += ring;
210         }
211     }
212     ps_close(out);
213 }
214
215 int gmx_wheel(int argc, char *argv[])
216 {
217     const char     *desc[] = {
218         "[THISMODULE] plots a helical wheel representation of your sequence.",
219         "The input sequence is in the [TT].dat[tt] file where the first line contains",
220         "the number of residues and each consecutive line contains a residue "
221         "name."
222     };
223     output_env_t    oenv;
224     static real     rot0  = 0;
225     static gmx_bool bNum  = TRUE;
226     static char    *title = NULL;
227     static int      r0    = 1;
228     t_pargs         pa [] = {
229         { "-r0",  FALSE, etINT, {&r0},
230           "The first residue number in the sequence" },
231         { "-rot0", FALSE, etREAL, {&rot0},
232           "Rotate around an angle initially (90 degrees makes sense)" },
233         { "-T",   FALSE, etSTR, {&title},
234           "Plot a title in the center of the wheel (must be shorter than 10 characters, or it will overwrite the wheel)" },
235         { "-nn",  FALSE, etBOOL, {&bNum},
236           "Toggle numbers" }
237     };
238     t_filenm        fnm[] = {
239         { efDAT, "-f", NULL,  ffREAD  },
240         { efEPS, "-o", NULL,  ffWRITE }
241     };
242 #define NFILE asize(fnm)
243
244     int    i, nres;
245     char **resnm;
246
247     if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
248                            asize(desc), desc, 0, NULL, &oenv))
249     {
250         return 0;
251     }
252
253     for (i = 1; (i < argc); i++)
254     {
255         if (strcmp(argv[i], "-r0") == 0)
256         {
257             r0 = strtol(argv[++i], NULL, 10);
258             fprintf(stderr, "First residue is %d\n", r0);
259         }
260         else if (strcmp(argv[i], "-rot0") == 0)
261         {
262             rot0 = strtod(argv[++i], NULL);
263             fprintf(stderr, "Initial rotation is %g\n", rot0);
264         }
265         else if (strcmp(argv[i], "-T") == 0)
266         {
267             title = gmx_strdup(argv[++i]);
268             fprintf(stderr, "Title will be '%s'\n", title);
269         }
270         else if (strcmp(argv[i], "-nn") == 0)
271         {
272             bNum = FALSE;
273             fprintf(stderr, "No residue numbers\n");
274         }
275         else
276         {
277             gmx_fatal(FARGS, "Incorrect usage of option %s", argv[i]);
278         }
279     }
280
281     nres = get_lines(ftp2fn(efDAT, NFILE, fnm), &resnm);
282     if (bNum)
283     {
284         wheel(ftp2fn(efEPS, NFILE, fnm), nres, resnm, r0, rot0, title);
285     }
286     else
287     {
288         wheel2(ftp2fn(efEPS, NFILE, fnm), nres, resnm, rot0, title);
289     }
290
291     return 0;
292 }