2 * This file is part of the GROMACS molecular simulation package.
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, 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.
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.
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.
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.
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.
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.
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"
58 static gmx_bool
*bPhobics(int nres
, char *resnm
[])
64 nb
= get_lines("phbres.dat", &cb
);
67 for (i
= 0; (i
< nres
); i
++)
69 if (search_str(nb
, cb
, resnm
[i
]) != -1)
77 static void wheel(const char *fn
, int nres
, char *resnm
[], int r0
, real rot0
, char *title
)
79 const real fontsize
= 16;
80 const real gray
= 0.9;
81 const real fontasp
= 0.6;
82 const real fontwidth
= fontsize
*fontasp
;
86 real ring
, inner
, outer
;
95 for (i
= 0; (i
< nres
); i
++)
98 sl
= std::strlen(resnm
[i
]);
99 sign
= resnm
[i
][sl
-1];
100 if ((sign
== '+') || (sign
== '-'))
102 resnm
[i
][sl
-1] = '\0';
104 sprintf(rnms
[i
], "%s-%d", resnm
[i
], i
+r0
);
105 if ((sign
== '+') || (sign
== '-'))
107 sl
= std::strlen(rnms
[i
]);
109 rnms
[i
][sl
+1] = '\0';
112 slen
= std::max(slen
, static_cast<int>(std::strlen(rnms
[i
])));
114 ring
= (2+slen
)*fontwidth
;
116 box
= inner
*1.5+(1+int{nres
/18})*ring
;
118 bPh
= bPhobics(nres
, resnm
);
120 out
= ps_open(fn
, 0, 0, 2.0*box
, 2.0*box
);
124 ps_font(out
, efontHELV
, 1.5*fontsize
);
125 ps_translate(out
, xc
, yc
);
128 ps_ctext(out
, 0, -fontsize
*1.5/2.0, title
, eXCenter
);
130 ps_font(out
, efontHELV
, fontsize
);
131 ps_rotate(out
, rot0
);
132 for (i
= 0; (i
< nres
); )
136 ps_color(out
, gray
, gray
, gray
);
137 ps_fillarcslice(out
, 0, 0, inner
, outer
, -10, 10);
138 ps_color(out
, 0, 0, 0);
140 ps_arcslice(out
, 0, 0, inner
, outer
, -10, 10);
142 ps_ctext(out
, inner
+fontwidth
, -fontsize
/2.0, rnms
[i
], eXLeft
);
143 ps_rotate(out
, -100);
155 static void wheel2(const char *fn
, int nres
, char *resnm
[], real rot0
, char *title
)
157 const real fontsize
= 14;
158 const real gray
= 0.9;
159 const real fontasp
= 0.45;
161 const real fontwidth
= fontsize
*fontasp
;
165 real ring
, inner
, outer
;
170 for (i
= 0; (i
< nres
); i
++)
172 slen
= std::max(slen
, static_cast<int>(strlen(resnm
[i
])));
174 fprintf(stderr
, "slen = %d\n", slen
);
175 ring
= slen
*fontwidth
;
177 box
= (1+gmx::exactDiv(nres
, 2*angle
))*outer
;
179 out
= ps_open(fn
, 0, 0, 2.0*box
, 2.0*box
);
183 ps_font(out
, efontHELV
, 1.5*fontsize
);
184 ps_translate(out
, xc
, yc
);
185 ps_color(out
, 0, 0, 0);
188 ps_ctext(out
, 0, -fontsize
*1.5/2.0, title
, eXCenter
);
190 ps_font(out
, efontHELV
, fontsize
);
192 ps_rotate(out
, rot0
);
193 for (i
= 0; (i
< nres
); )
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);
201 ps_arcslice(out
, 0, 0, inner
, outer
, -angle
, angle
);
203 ps_ctext(out
, inner
+fontwidth
, -fontsize
/2.0, resnm
[i
], eXLeft
);
204 ps_rotate(out
, -2*angle
);
207 if ((i
% (2*angle
)) == 0)
216 int gmx_wheel(int argc
, char *argv
[])
218 const char *desc
[] = {
219 "[THISMODULE] plots a helical wheel representation of your sequence.",
220 "The input sequence is in the [REF].dat[ref] file where the first line contains",
221 "the number of residues and each consecutive line contains a residue "
224 gmx_output_env_t
*oenv
;
225 static real rot0
= 0;
226 static gmx_bool bNum
= TRUE
;
227 static char *title
= nullptr;
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
},
240 { efDAT
, "-f", nullptr, ffREAD
},
241 { efEPS
, "-o", nullptr, ffWRITE
}
243 #define NFILE asize(fnm)
248 if (!parse_common_args(&argc
, argv
, 0, NFILE
, fnm
, asize(pa
), pa
,
249 asize(desc
), desc
, 0, nullptr, &oenv
))
254 for (i
= 1; (i
< argc
); i
++)
256 if (std::strcmp(argv
[i
], "-r0") == 0)
258 r0
= std::strtol(argv
[++i
], nullptr, 10);
259 fprintf(stderr
, "First residue is %d\n", r0
);
261 else if (std::strcmp(argv
[i
], "-rot0") == 0)
263 rot0
= strtod(argv
[++i
], nullptr);
264 fprintf(stderr
, "Initial rotation is %g\n", rot0
);
266 else if (std::strcmp(argv
[i
], "-T") == 0)
268 title
= gmx_strdup(argv
[++i
]);
269 fprintf(stderr
, "Title will be '%s'\n", title
);
271 else if (std::strcmp(argv
[i
], "-nn") == 0)
274 fprintf(stderr
, "No residue numbers\n");
278 gmx_fatal(FARGS
, "Incorrect usage of option %s", argv
[i
]);
282 nres
= get_lines(ftp2fn(efDAT
, NFILE
, fnm
), &resnm
);
285 wheel(ftp2fn(efEPS
, NFILE
, fnm
), nres
, resnm
, r0
, rot0
, title
);
289 wheel2(ftp2fn(efEPS
, NFILE
, fnm
), nres
, resnm
, rot0
, title
);