Make PBC type enumeration into PbcType enum class
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_helix.cpp
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,2015,2016,2017 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include <cmath>
41 #include <cstring>
42
43 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/commandline/viewit.h"
45 #include "gromacs/fileio/confio.h"
46 #include "gromacs/fileio/tpxio.h"
47 #include "gromacs/fileio/trxio.h"
48 #include "gromacs/fileio/xvgr.h"
49 #include "gromacs/gmxana/fitahx.h"
50 #include "gromacs/gmxana/gmx_ana.h"
51 #include "gromacs/gmxana/gstat.h"
52 #include "gromacs/gmxana/hxprops.h"
53 #include "gromacs/math/utilities.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/pbcutil/rmpbc.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/utility/arraysize.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/futil.h"
60 #include "gromacs/utility/smalloc.h"
61
62 int gmx_helix(int argc, char* argv[])
63 {
64     const char* desc[] = {
65         "[THISMODULE] computes all kinds of helix properties. First, the peptide",
66         "is checked to find the longest helical part, as determined by",
67         "hydrogen bonds and [GRK]phi[grk]/[GRK]psi[grk] angles.",
68         "That bit is fitted",
69         "to an ideal helix around the [IT]z[it]-axis and centered around the origin.",
70         "Then the following properties are computed:",
71         "",
72         " * Helix radius (file [TT]radius.xvg[tt]). This is merely the",
73         "   RMS deviation in two dimensions for all C[GRK]alpha[grk] atoms.",
74         "   it is calculated as [SQRT]([SUM][sum][SUB]i[sub] (x^2(i)+y^2(i)))/N[sqrt] where N is ",
75         "   the number of backbone atoms. For an ideal helix the radius is 0.23 nm.",
76         " * Twist (file [TT]twist.xvg[tt]). The average helical angle per",
77         "   residue is calculated. For an [GRK]alpha[grk]-helix it is 100 degrees,",
78         "   for 3-10 helices it will be smaller, and ",
79         "   for 5-helices it will be larger.",
80         " * Rise per residue (file [TT]rise.xvg[tt]). The helical rise per",
81         "   residue is plotted as the difference in [IT]z[it]-coordinate between C[GRK]alpha[grk]",
82         "   atoms. For an ideal helix, this is 0.15 nm.",
83         " * Total helix length (file [TT]len-ahx.xvg[tt]). The total length",
84         "   of the",
85         "   helix in nm. This is simply the average rise (see above) times the",
86         "   number of helical residues (see below).",
87         " * Helix dipole, backbone only (file [TT]dip-ahx.xvg[tt]).",
88         " * RMS deviation from ideal helix, calculated for the C[GRK]alpha[grk]",
89         "   atoms only (file [TT]rms-ahx.xvg[tt]).",
90         " * Average C[GRK]alpha[grk] - C[GRK]alpha[grk] dihedral angle (file [TT]phi-ahx.xvg[tt]).",
91         " * Average [GRK]phi[grk] and [GRK]psi[grk] angles (file [TT]phipsi.xvg[tt]).",
92         " * Ellipticity at 222 nm according to Hirst and Brooks."
93     };
94     static gmx_bool bCheck = FALSE, bFit = TRUE, bDBG = FALSE, bEV = FALSE;
95     static int      rStart = 0, rEnd = 0, r0 = 1;
96     t_pargs pa[] = { { "-r0", FALSE, etINT, { &r0 }, "The first residue number in the sequence" },
97                      { "-q",
98                        FALSE,
99                        etBOOL,
100                        { &bCheck },
101                        "Check at every step which part of the sequence is helical" },
102                      { "-F", FALSE, etBOOL, { &bFit }, "Toggle fit to a perfect helix" },
103                      { "-db", FALSE, etBOOL, { &bDBG }, "Print debug info" },
104                      { "-ev", FALSE, etBOOL, { &bEV }, "Write a new 'trajectory' file for ED" },
105                      { "-ahxstart", FALSE, etINT, { &rStart }, "First residue in helix" },
106                      { "-ahxend", FALSE, etINT, { &rEnd }, "Last residue in helix" } };
107
108     typedef struct
109     { //NOLINT(clang-analyzer-optin.performance.Padding)
110         FILE *      fp, *fp2;
111         gmx_bool    bfp2;
112         const char* filenm;
113         const char* title;
114         const char* xaxis;
115         const char* yaxis;
116         real        val;
117     } t_xvgrfile;
118
119     t_xvgrfile xf[efhNR] = {
120         { nullptr, nullptr, TRUE, "radius", "Helix radius", nullptr, "r (nm)", 0.0 },
121         { nullptr, nullptr, TRUE, "twist", "Twist per residue", nullptr, "Angle (deg)", 0.0 },
122         { nullptr, nullptr, TRUE, "rise", "Rise per residue", nullptr, "Rise (nm)", 0.0 },
123         { nullptr, nullptr, FALSE, "len-ahx", "Length of the Helix", nullptr, "Length (nm)", 0.0 },
124         { nullptr, nullptr, FALSE, "dip-ahx", "Helix Backbone Dipole", nullptr, "rq (nm e)", 0.0 },
125         { nullptr, nullptr, TRUE, "rms-ahx", "RMS Deviation from Ideal Helix", nullptr, "RMS (nm)", 0.0 },
126         { nullptr, nullptr, FALSE, "rmsa-ahx", "Average RMSD per Residue", "Residue", "RMS (nm)", 0.0 },
127         { nullptr, nullptr, FALSE, "cd222", "Ellipticity at 222 nm", nullptr, "nm", 0.0 },
128         { nullptr, nullptr, TRUE, "pprms", "RMS Distance from \\8a\\4-helix", nullptr, "deg", 0.0 },
129         { nullptr, nullptr, TRUE, "caphi", "Average Ca-Ca Dihedral", nullptr, "\\8F\\4(deg)", 0.0 },
130         { nullptr, nullptr, TRUE, "phi", "Average \\8F\\4 angles", nullptr, "deg", 0.0 },
131         { nullptr, nullptr, TRUE, "psi", "Average \\8Y\\4 angles", nullptr, "deg", 0.0 },
132         { nullptr, nullptr, TRUE, "hb3", "Average n-n+3 hbond length", nullptr, "nm", 0.0 },
133         { nullptr, nullptr, TRUE, "hb4", "Average n-n+4 hbond length", nullptr, "nm", 0.0 },
134         { nullptr, nullptr, TRUE, "hb5", "Average n-n+5 hbond length", nullptr, "nm", 0.0 },
135         { nullptr, nullptr, FALSE, "JCaHa", "J-Coupling Values", "Residue", "Hz", 0.0 },
136         { nullptr, nullptr, FALSE, "helicity", "Helicity per Residue", "Residue", "% of time", 0.0 }
137     };
138
139     gmx_output_env_t* oenv;
140     char              buf[54];
141     t_trxstatus*      status;
142     int               natoms, nres;
143     t_bb*             bb;
144     int               i, j, nall, nbb, nca, teller;
145     int *             bbindex, *caindex, *allindex;
146     t_topology*       top;
147     PbcType           pbcType;
148     rvec *            x, *xref;
149     real              t;
150     real              rms;
151     matrix            box;
152     gmx_rmpbc_t       gpbc = nullptr;
153     gmx_bool          bRange;
154     t_filenm          fnm[] = {
155         { efTPR, nullptr, nullptr, ffREAD },
156         { efNDX, nullptr, nullptr, ffREAD },
157         { efTRX, "-f", nullptr, ffREAD },
158         { efSTO, "-cz", "zconf", ffWRITE },
159     };
160 #define NFILE asize(fnm)
161
162     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME, NFILE, fnm, asize(pa), pa,
163                            asize(desc), desc, 0, nullptr, &oenv))
164     {
165         return 0;
166     }
167
168     bRange = (opt2parg_bSet("-ahxstart", asize(pa), pa) && opt2parg_bSet("-ahxend", asize(pa), pa));
169
170     top = read_top(ftp2fn(efTPR, NFILE, fnm), &pbcType);
171
172     natoms = read_first_x(oenv, &status, opt2fn("-f", NFILE, fnm), &t, &x, box);
173
174     if (natoms != top->atoms.nr)
175     {
176         gmx_fatal(FARGS,
177                   "Sorry can only run when the number of atoms in the run input file (%d) is equal "
178                   "to the number in the trajectory (%d)",
179                   top->atoms.nr, natoms);
180     }
181
182     bb = mkbbind(ftp2fn(efNDX, NFILE, fnm), &nres, &nbb, r0, &nall, &allindex, top->atoms.atomname,
183                  top->atoms.atom, top->atoms.resinfo);
184     snew(bbindex, natoms);
185     snew(caindex, nres);
186
187     fprintf(stderr, "nall=%d\n", nall);
188
189     /* Open output files, default x-axis is time */
190     for (i = 0; (i < efhNR); i++)
191     {
192         sprintf(buf, "%s.xvg", xf[i].filenm);
193         remove(buf);
194         xf[i].fp = xvgropen(buf, xf[i].title, xf[i].xaxis ? xf[i].xaxis : "Time (ps)", xf[i].yaxis, oenv);
195         if (xf[i].bfp2)
196         {
197             sprintf(buf, "%s.out", xf[i].filenm);
198             remove(buf);
199             xf[i].fp2 = gmx_ffopen(buf, "w");
200         }
201     }
202
203     /* Read reference frame from tpx file to compute helix length */
204     snew(xref, top->atoms.nr);
205     read_tpx(ftp2fn(efTPR, NFILE, fnm), nullptr, nullptr, nullptr, xref, nullptr, nullptr);
206     calc_hxprops(nres, bb, xref);
207     do_start_end(nres, bb, &nbb, bbindex, &nca, caindex, bRange, rStart, rEnd);
208     sfree(xref);
209     if (bDBG)
210     {
211         fprintf(stderr, "nca=%d, nbb=%d\n", nca, nbb);
212         pr_bb(stdout, nres, bb);
213     }
214
215     gpbc = gmx_rmpbc_init(&top->idef, pbcType, natoms);
216
217     teller = 0;
218     do
219     {
220         if ((teller++ % 10) == 0)
221         {
222             fprintf(stderr, "\rt=%.2f", t);
223             fflush(stderr);
224         }
225         gmx_rmpbc(gpbc, natoms, box, x);
226
227
228         calc_hxprops(nres, bb, x);
229         if (bCheck)
230         {
231             do_start_end(nres, bb, &nbb, bbindex, &nca, caindex, FALSE, 0, 0);
232         }
233
234         if (nca >= 5)
235         {
236             rms = fit_ahx(nres, bb, natoms, nall, allindex, x, nca, caindex, bFit);
237
238             if (teller == 1)
239             {
240                 write_sto_conf(opt2fn("-cz", NFILE, fnm), "Helix fitted to Z-Axis", &(top->atoms),
241                                x, nullptr, pbcType, box);
242             }
243
244             xf[efhRAD].val   = radius(xf[efhRAD].fp2, nca, caindex, x);
245             xf[efhTWIST].val = twist(nca, caindex, x);
246             xf[efhRISE].val  = rise(nca, caindex, x);
247             xf[efhLEN].val   = ahx_len(nca, caindex, x);
248             xf[efhCD222].val = ellipticity(nres, bb);
249             xf[efhDIP].val   = dip(nbb, bbindex, x, top->atoms.atom);
250             xf[efhRMS].val   = rms;
251             xf[efhCPHI].val  = ca_phi(nca, caindex, x);
252             xf[efhPPRMS].val = pprms(xf[efhPPRMS].fp2, nres, bb);
253
254             for (j = 0; (j <= efhCPHI); j++)
255             {
256                 fprintf(xf[j].fp, "%10g  %10g\n", t, xf[j].val);
257             }
258
259             av_phipsi(xf[efhPHI].fp, xf[efhPSI].fp, xf[efhPHI].fp2, xf[efhPSI].fp2, t, nres, bb);
260             av_hblen(xf[efhHB3].fp, xf[efhHB3].fp2, xf[efhHB4].fp, xf[efhHB4].fp2, xf[efhHB5].fp,
261                      xf[efhHB5].fp2, t, nres, bb);
262         }
263     } while (read_next_x(oenv, status, &t, x, box));
264     fprintf(stderr, "\n");
265
266     gmx_rmpbc_done(gpbc);
267
268     close_trx(status);
269
270     for (i = 0; (i < nres); i++)
271     {
272         if (bb[i].nrms > 0)
273         {
274             fprintf(xf[efhRMSA].fp, "%10d  %10g\n", r0 + i, bb[i].rmsa / bb[i].nrms);
275         }
276         fprintf(xf[efhAHX].fp, "%10d  %10g\n", r0 + i, (bb[i].nhx * 100.0) / static_cast<real>(teller));
277         fprintf(xf[efhJCA].fp, "%10d  %10g\n", r0 + i,
278                 140.3 + (bb[i].jcaha / static_cast<double>(teller)));
279     }
280
281     for (i = 0; (i < efhNR); i++)
282     {
283         xvgrclose(xf[i].fp);
284         if (xf[i].bfp2)
285         {
286             xvgrclose(xf[i].fp2);
287         }
288         do_view(oenv, xf[i].filenm, "-nxy");
289     }
290
291     return 0;
292 }