Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / gmxana / gmx_helix.cpp
blob7f659c1e1530aa72b7269a009a45518bc4e31422
1 /*
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,2016,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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstring>
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/commandline/viewit.h"
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/fileio/tpxio.h"
46 #include "gromacs/fileio/trxio.h"
47 #include "gromacs/fileio/xvgr.h"
48 #include "gromacs/gmxana/fitahx.h"
49 #include "gromacs/gmxana/gmx_ana.h"
50 #include "gromacs/gmxana/gstat.h"
51 #include "gromacs/gmxana/hxprops.h"
52 #include "gromacs/math/utilities.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/pbcutil/rmpbc.h"
55 #include "gromacs/topology/topology.h"
56 #include "gromacs/utility/arraysize.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/utility/smalloc.h"
61 int gmx_helix(int argc, char *argv[])
63 const char *desc[] = {
64 "[THISMODULE] computes all kinds of helix properties. First, the peptide",
65 "is checked to find the longest helical part, as determined by",
66 "hydrogen bonds and [GRK]phi[grk]/[GRK]psi[grk] angles.",
67 "That bit is fitted",
68 "to an ideal helix around the [IT]z[it]-axis and centered around the origin.",
69 "Then the following properties are computed:",
70 "",
71 " * Helix radius (file [TT]radius.xvg[tt]). This is merely the",
72 " RMS deviation in two dimensions for all C[GRK]alpha[grk] atoms.",
73 " it is calculated as [SQRT]([SUM][sum][SUB]i[sub] (x^2(i)+y^2(i)))/N[sqrt] where N is the number",
74 " of backbone atoms. For an ideal helix the radius is 0.23 nm.",
75 " * Twist (file [TT]twist.xvg[tt]). The average helical angle per",
76 " residue is calculated. For an [GRK]alpha[grk]-helix it is 100 degrees,",
77 " for 3-10 helices it will be smaller, and ",
78 " for 5-helices it will be larger.",
79 " * Rise per residue (file [TT]rise.xvg[tt]). The helical rise per",
80 " residue is plotted as the difference in [IT]z[it]-coordinate between C[GRK]alpha[grk]",
81 " atoms. For an ideal helix, this is 0.15 nm.",
82 " * Total helix length (file [TT]len-ahx.xvg[tt]). The total length",
83 " of the",
84 " helix in nm. This is simply the average rise (see above) times the",
85 " number of helical residues (see below).",
86 " * Helix dipole, backbone only (file [TT]dip-ahx.xvg[tt]).",
87 " * RMS deviation from ideal helix, calculated for the C[GRK]alpha[grk]",
88 " atoms only (file [TT]rms-ahx.xvg[tt]).",
89 " * Average C[GRK]alpha[grk] - C[GRK]alpha[grk] dihedral angle (file [TT]phi-ahx.xvg[tt]).",
90 " * Average [GRK]phi[grk] and [GRK]psi[grk] angles (file [TT]phipsi.xvg[tt]).",
91 " * Ellipticity at 222 nm according to Hirst and Brooks."
93 static gmx_bool bCheck = FALSE, bFit = TRUE, bDBG = FALSE, bEV = FALSE;
94 static int rStart = 0, rEnd = 0, r0 = 1;
95 t_pargs pa [] = {
96 { "-r0", FALSE, etINT, {&r0},
97 "The first residue number in the sequence" },
98 { "-q", FALSE, etBOOL, {&bCheck},
99 "Check at every step which part of the sequence is helical" },
100 { "-F", FALSE, etBOOL, {&bFit},
101 "Toggle fit to a perfect helix" },
102 { "-db", FALSE, etBOOL, {&bDBG},
103 "Print debug info" },
104 { "-ev", FALSE, etBOOL, {&bEV},
105 "Write a new 'trajectory' file for ED" },
106 { "-ahxstart", FALSE, etINT, {&rStart},
107 "First residue in helix" },
108 { "-ahxend", FALSE, etINT, {&rEnd},
109 "Last residue in helix" }
112 typedef struct { //NOLINT(clang-analyzer-optin.performance.Padding)
113 FILE *fp, *fp2;
114 gmx_bool bfp2;
115 const char *filenm;
116 const char *title;
117 const char *xaxis;
118 const char *yaxis;
119 real val;
120 } t_xvgrfile;
122 t_xvgrfile xf[efhNR] = {
123 { nullptr, nullptr, TRUE, "radius", "Helix radius", nullptr, "r (nm)", 0.0 },
124 { nullptr, nullptr, TRUE, "twist", "Twist per residue", nullptr, "Angle (deg)", 0.0 },
125 { nullptr, nullptr, TRUE, "rise", "Rise per residue", nullptr, "Rise (nm)", 0.0 },
126 { nullptr, nullptr, FALSE, "len-ahx", "Length of the Helix", nullptr, "Length (nm)", 0.0 },
127 { nullptr, nullptr, FALSE, "dip-ahx", "Helix Backbone Dipole", nullptr, "rq (nm e)", 0.0 },
128 { nullptr, nullptr, TRUE, "rms-ahx", "RMS Deviation from Ideal Helix", nullptr, "RMS (nm)", 0.0 },
129 { nullptr, nullptr, FALSE, "rmsa-ahx", "Average RMSD per Residue", "Residue", "RMS (nm)", 0.0 },
130 { nullptr, nullptr, FALSE, "cd222", "Ellipticity at 222 nm", nullptr, "nm", 0.0 },
131 { nullptr, nullptr, TRUE, "pprms", "RMS Distance from \\8a\\4-helix", nullptr, "deg", 0.0 },
132 { nullptr, nullptr, TRUE, "caphi", "Average Ca-Ca Dihedral", nullptr, "\\8F\\4(deg)", 0.0 },
133 { nullptr, nullptr, TRUE, "phi", "Average \\8F\\4 angles", nullptr, "deg", 0.0 },
134 { nullptr, nullptr, TRUE, "psi", "Average \\8Y\\4 angles", nullptr, "deg", 0.0 },
135 { nullptr, nullptr, TRUE, "hb3", "Average n-n+3 hbond length", nullptr, "nm", 0.0 },
136 { nullptr, nullptr, TRUE, "hb4", "Average n-n+4 hbond length", nullptr, "nm", 0.0 },
137 { nullptr, nullptr, TRUE, "hb5", "Average n-n+5 hbond length", nullptr, "nm", 0.0 },
138 { nullptr, nullptr, FALSE, "JCaHa", "J-Coupling Values", "Residue", "Hz", 0.0 },
139 { nullptr, nullptr, FALSE, "helicity", "Helicity per Residue", "Residue", "% of time", 0.0 }
142 gmx_output_env_t *oenv;
143 char buf[54];
144 t_trxstatus *status;
145 int natoms, nres;
146 t_bb *bb;
147 int i, j, nall, nbb, nca, teller;
148 int *bbindex, *caindex, *allindex;
149 t_topology *top;
150 int ePBC;
151 rvec *x, *xref;
152 real t;
153 real rms;
154 matrix box;
155 gmx_rmpbc_t gpbc = nullptr;
156 gmx_bool bRange;
157 t_filenm fnm[] = {
158 { efTPR, nullptr, nullptr, ffREAD },
159 { efNDX, nullptr, nullptr, ffREAD },
160 { efTRX, "-f", nullptr, ffREAD },
161 { efSTO, "-cz", "zconf", ffWRITE },
163 #define NFILE asize(fnm)
165 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME,
166 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
168 return 0;
171 bRange = (opt2parg_bSet("-ahxstart", asize(pa), pa) &&
172 opt2parg_bSet("-ahxend", asize(pa), pa));
174 top = read_top(ftp2fn(efTPR, NFILE, fnm), &ePBC);
176 natoms = read_first_x(oenv, &status, opt2fn("-f", NFILE, fnm), &t, &x, box);
178 if (natoms != top->atoms.nr)
180 gmx_fatal(FARGS, "Sorry can only run when the number of atoms in the run input file (%d) is equal to the number in the trajectory (%d)",
181 top->atoms.nr, natoms);
184 bb = mkbbind(ftp2fn(efNDX, NFILE, fnm), &nres, &nbb, r0, &nall, &allindex,
185 top->atoms.atomname, top->atoms.atom, top->atoms.resinfo);
186 snew(bbindex, natoms);
187 snew(caindex, nres);
189 fprintf(stderr, "nall=%d\n", nall);
191 /* Open output files, default x-axis is time */
192 for (i = 0; (i < efhNR); i++)
194 sprintf(buf, "%s.xvg", xf[i].filenm);
195 remove(buf);
196 xf[i].fp = xvgropen(buf, xf[i].title,
197 xf[i].xaxis ? xf[i].xaxis : "Time (ps)",
198 xf[i].yaxis, oenv);
199 if (xf[i].bfp2)
201 sprintf(buf, "%s.out", xf[i].filenm);
202 remove(buf);
203 xf[i].fp2 = gmx_ffopen(buf, "w");
207 /* Read reference frame from tpx file to compute helix length */
208 snew(xref, top->atoms.nr);
209 read_tpx(ftp2fn(efTPR, NFILE, fnm),
210 nullptr, nullptr, nullptr, xref, nullptr, nullptr);
211 calc_hxprops(nres, bb, xref);
212 do_start_end(nres, bb, &nbb, bbindex, &nca, caindex, bRange, rStart, rEnd);
213 sfree(xref);
214 if (bDBG)
216 fprintf(stderr, "nca=%d, nbb=%d\n", nca, nbb);
217 pr_bb(stdout, nres, bb);
220 gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
222 teller = 0;
225 if ((teller++ % 10) == 0)
227 fprintf(stderr, "\rt=%.2f", t);
228 fflush(stderr);
230 gmx_rmpbc(gpbc, natoms, box, x);
233 calc_hxprops(nres, bb, x);
234 if (bCheck)
236 do_start_end(nres, bb, &nbb, bbindex, &nca, caindex, FALSE, 0, 0);
239 if (nca >= 5)
241 rms = fit_ahx(nres, bb, natoms, nall, allindex, x, nca, caindex, bFit);
243 if (teller == 1)
245 write_sto_conf(opt2fn("-cz", NFILE, fnm), "Helix fitted to Z-Axis",
246 &(top->atoms), x, nullptr, ePBC, box);
249 xf[efhRAD].val = radius(xf[efhRAD].fp2, nca, caindex, x);
250 xf[efhTWIST].val = twist(nca, caindex, x);
251 xf[efhRISE].val = rise(nca, caindex, x);
252 xf[efhLEN].val = ahx_len(nca, caindex, x);
253 xf[efhCD222].val = ellipticity(nres, bb);
254 xf[efhDIP].val = dip(nbb, bbindex, x, top->atoms.atom);
255 xf[efhRMS].val = rms;
256 xf[efhCPHI].val = ca_phi(nca, caindex, x);
257 xf[efhPPRMS].val = pprms(xf[efhPPRMS].fp2, nres, bb);
259 for (j = 0; (j <= efhCPHI); j++)
261 fprintf(xf[j].fp, "%10g %10g\n", t, xf[j].val);
264 av_phipsi(xf[efhPHI].fp, xf[efhPSI].fp, xf[efhPHI].fp2, xf[efhPSI].fp2,
265 t, nres, bb);
266 av_hblen(xf[efhHB3].fp, xf[efhHB3].fp2,
267 xf[efhHB4].fp, xf[efhHB4].fp2,
268 xf[efhHB5].fp, xf[efhHB5].fp2,
269 t, nres, bb);
272 while (read_next_x(oenv, status, &t, x, box));
273 fprintf(stderr, "\n");
275 gmx_rmpbc_done(gpbc);
277 close_trx(status);
279 for (i = 0; (i < nres); i++)
281 if (bb[i].nrms > 0)
283 fprintf(xf[efhRMSA].fp, "%10d %10g\n", r0+i, bb[i].rmsa/bb[i].nrms);
285 fprintf(xf[efhAHX].fp, "%10d %10g\n", r0+i, (bb[i].nhx*100.0)/static_cast<real>(teller));
286 fprintf(xf[efhJCA].fp, "%10d %10g\n",
287 r0+i, 140.3+(bb[i].jcaha/static_cast<double>(teller)));
290 for (i = 0; (i < efhNR); i++)
292 xvgrclose(xf[i].fp);
293 if (xf[i].bfp2)
295 xvgrclose(xf[i].fp2);
297 do_view(oenv, xf[i].filenm, "-nxy");
300 return 0;