3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
57 #include "gmx_fatal.h"
62 static gmx_bool *bPhobics(int nres, char *resnm[])
68 nb = get_strings("phbres.dat", &cb);
71 for (i = 0; (i < nres); i++)
73 if (search_str(nb, cb, resnm[i]) != -1)
81 void wheel(const char *fn, int nres, char *resnm[], int r0, real rot0, char *title)
83 const real fontsize = 16;
84 const real gray = 0.9;
85 const real fontasp = 0.6;
86 const real fontwidth = fontsize*fontasp;
90 real ring, inner, outer;
99 for (i = 0; (i < nres); i++)
102 sl = strlen(resnm[i]);
103 sign = resnm[i][sl-1];
104 if ((sign == '+') || (sign == '-'))
106 resnm[i][sl-1] = '\0';
108 sprintf(rnms[i], "%s-%d", resnm[i], i+r0);
109 if ((sign == '+') || (sign == '-'))
111 sl = strlen(rnms[i]);
113 rnms[i][sl+1] = '\0';
116 slen = max(slen, (int)strlen(rnms[i]));
118 ring = (2+slen)*fontwidth;
120 box = inner*1.5+(1+(nres / 18))*ring;
122 bPh = bPhobics(nres, resnm);
124 out = ps_open(fn, 0, 0, 2.0*box, 2.0*box);
128 ps_font(out, efontHELV, 1.5*fontsize);
129 ps_translate(out, xc, yc);
132 ps_ctext(out, 0, -fontsize*1.5/2.0, title, eXCenter);
134 ps_font(out, efontHELV, fontsize);
135 ps_rotate(out, rot0);
136 for (i = 0; (i < nres); )
140 ps_color(out, gray, gray, gray);
141 ps_fillarcslice(out, 0, 0, inner, outer, -10, 10);
142 ps_color(out, 0, 0, 0);
144 ps_arcslice(out, 0, 0, inner, outer, -10, 10);
146 ps_ctext(out, inner+fontwidth, -fontsize/2.0, rnms[i], eXLeft);
147 ps_rotate(out, -100);
159 void wheel2(const char *fn, int nres, char *resnm[], real rot0, char *title)
161 const real fontsize = 14;
162 const real gray = 0.9;
163 const real fontasp = 0.45;
165 const real fontwidth = fontsize*fontasp;
169 real ring, inner, outer;
174 for (i = 0; (i < nres); i++)
176 slen = max(slen, (int)strlen(resnm[i]));
178 fprintf(stderr, "slen = %d\n", slen);
179 ring = (slen)*fontwidth;
181 box = (1+(nres / (2*angle)))*outer;
183 out = ps_open(fn, 0, 0, 2.0*box, 2.0*box);
187 ps_font(out, efontHELV, 1.5*fontsize);
188 ps_translate(out, xc, yc);
189 ps_color(out, 0, 0, 0);
192 ps_ctext(out, 0, -fontsize*1.5/2.0, title, eXCenter);
194 ps_font(out, efontHELV, fontsize);
196 ps_rotate(out, rot0);
197 for (i = 0; (i < nres); )
201 ps_color(out, gray, gray, 1.0);
202 ps_fillarcslice(out, 0, 0, inner, outer, -angle, angle);
203 ps_color(out, 0, 0, 0);
205 ps_arcslice(out, 0, 0, inner, outer, -angle, angle);
207 ps_ctext(out, inner+fontwidth, -fontsize/2.0, resnm[i], eXLeft);
208 ps_rotate(out, -2*angle);
211 if ((i % (2*angle)) == 0)
220 int gmx_wheel(int argc, char *argv[])
222 const char *desc[] = {
223 "[TT]g_wheel[tt] plots a helical wheel representation of your sequence.",
224 "The input sequence is in the [TT].dat[tt] file where the first line contains",
225 "the number of residues and each consecutive line contains a residue "
229 static real rot0 = 0;
230 static gmx_bool bNum = TRUE;
231 static char *title = NULL;
234 { "-r0", FALSE, etINT, {&r0},
235 "The first residue number in the sequence" },
236 { "-rot0", FALSE, etREAL, {&rot0},
237 "Rotate around an angle initially (90 degrees makes sense)" },
238 { "-T", FALSE, etSTR, {&title},
239 "Plot a title in the center of the wheel (must be shorter than 10 characters, or it will overwrite the wheel)" },
240 { "-nn", FALSE, etBOOL, {&bNum},
244 { efDAT, "-f", NULL, ffREAD },
245 { efEPS, "-o", NULL, ffWRITE }
247 #define NFILE asize(fnm)
252 parse_common_args(&argc, argv, PCA_BE_NICE, NFILE, fnm, asize(pa), pa,
253 asize(desc), desc, 0, NULL, &oenv);
255 for (i = 1; (i < argc); i++)
257 if (strcmp(argv[i], "-r0") == 0)
259 r0 = strtol(argv[++i], NULL, 10);
260 fprintf(stderr, "First residue is %d\n", r0);
262 else if (strcmp(argv[i], "-rot0") == 0)
264 rot0 = strtod(argv[++i], NULL);
265 fprintf(stderr, "Initial rotation is %g\n", rot0);
267 else if (strcmp(argv[i], "-T") == 0)
269 title = strdup(argv[++i]);
270 fprintf(stderr, "Title will be '%s'\n", title);
272 else if (strcmp(argv[i], "-nn") == 0)
275 fprintf(stderr, "No residue numbers\n");
279 gmx_fatal(FARGS, "Incorrect usage of option %s", argv[i]);
283 nres = get_lines(ftp2fn(efDAT, NFILE, fnm), &resnm);
286 wheel(ftp2fn(efEPS, NFILE, fnm), nres, resnm, r0, rot0, title);
290 wheel2(ftp2fn(efEPS, NFILE, fnm), nres, resnm, rot0, title);