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, 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 #ifndef GMX_GMXANA_GSTAT_H
38 #define GMX_GMXANA_GSTAT_H
40 #include "gromacs/commandline/pargs.h"
41 #include "gromacs/legacyheaders/oenv.h"
42 #include "gromacs/legacyheaders/typedefs.h"
43 #include "gromacs/topology/index.h"
49 struct gmx_residuetype_t;
51 /* must correspond with 'leg' g_chi.c:727 */
53 edPhi = 0, edPsi, edOmega, edChi1, edChi2, edChi3, edChi4, edChi5, edChi6, edMax
57 edPrintST = 0, edPrintRO
62 #define MAXCHI edMax-NONCHI
63 #define NROT 4 /* number of rotamers: 1=g(-), 2=t, 3=g(+), 0=other */
66 int minCalpha, minC, H, N, C, O, Cn[MAXCHI+3];
67 } t_dihatms; /* Cn[0]=N, Cn[1]=Ca, Cn[2]=Cb etc. */
72 int index; /* Index for amino acids (histograms) */
73 int j0[edMax]; /* Index in dih array (phi angle is first...) */
78 real rot_occ[edMax][NROT];
83 const char *name; /* Description of the J coupling constant */
84 real A, B, C; /* Karplus coefficients */
85 real offset; /* Offset for dihedral angle in histogram (e.g. -M_PI/3) */
86 real Jc; /* Resulting Jcoupling */
87 real Jcsig; /* Standard deviation in Jc */
90 void calc_distribution_props(int nh, int histo[],
91 real start, int nkkk, t_karplus kkk[],
93 /* This routine takes a dihedral distribution and calculates
94 * coupling constants and dihedral order parameters of it.
96 * nh is the number of points
97 * histo is the array of datapoints which is assumed to span
99 * start is the starting angle of the histogram, this can be either 0
101 * nkkk is the number of karplus sets (multiple coupling constants may be
102 * derived from a single angle)
103 * kkk are the constants for calculating J coupling constants using a
104 * Karplus equation according to
107 * J = A cos theta + B cos theta + C
109 * where theta is phi - offset (phi is the angle in the histogram)
110 * offset is subtracted from phi before substitution in the Karplus
112 * S2 is the resulting dihedral order parameter
116 void ana_dih_trans(const char *fn_trans, const char *fn_histo,
117 real **dih, int nframes, int nangles,
118 const char *grpname, real *time, gmx_bool bRb,
119 const output_env_t oenv);
121 * Analyse dihedral transitions, by counting transitions per dihedral
122 * and per frame. The total number of transitions is printed to
123 * stderr, as well as the average time between transitions.
125 * is wrapper to low_ana_dih_trans, which also passes in and out the
126 number of transitions per dihedral per residue. that uses struc dlist
127 which is not external, so pp2shift.h must be included.
129 * Dihedrals are supposed to be in either of three minima,
130 * (trans, gauche+, gauche-)
132 * fn_trans output file name for #transitions per timeframe
133 * fn_histo output file name for transition time histogram
134 * dih the actual dihedral angles
135 * nframes number of times frames
136 * nangles number of angles
137 * grpname a string for the header of plots
138 * time array (size nframes) of times of trajectory frames
139 * bRb determines whether the polymer convention is used
143 void low_ana_dih_trans(gmx_bool bTrans, const char *fn_trans,
144 gmx_bool bHisto, const char *fn_histo, int maxchi,
145 real **dih, int nlist, t_dlist dlist[],
146 int nframes, int nangles, const char *grpname,
147 int multiplicity[], real *time, gmx_bool bRb,
148 real core_frac, const output_env_t oenv);
149 /* as above but passes dlist so can copy occupancies into it, and multiplicity[]
150 * (1..nangles, corresp to dih[this][], so can have non-3 multiplicity of
151 * rotamers. Also production of xvg output files is conditional
152 * and the fractional width of each rotamer can be set ie for a 3 fold
153 * dihedral with core_frac = 0.5 only the central 60 degrees is assigned
154 * to each rotamer, the rest goes to rotamer zero */
158 void read_ang_dih(const char *trj_fn,
159 gmx_bool bAngles, gmx_bool bSaveAll, gmx_bool bRb, gmx_bool bPBC,
160 int maxangstat, int angstat[],
161 int *nframes, real **time,
162 int isize, atom_id index[],
166 const output_env_t oenv);
168 * Read a trajectory and calculate angles and dihedrals.
170 * trj_fn file name of trajectory
171 * bAngles do we have to read angles or dihedrals
172 * bSaveAll do we have to store all in the dih array
173 * bRb do we have Ryckaert-Bellemans dihedrals (trans = 0)
174 * bPBC compute angles module 2 Pi
175 * maxangstat number of entries in distribution array
176 * angstat angle distribution
177 * *nframes number of frames read
178 * time simulation time at each time frame
179 * isize number of entries in the index, when angles 3*number of angles
180 * else 4*number of angles
181 * index atom numbers that define the angles or dihedrals
182 * (i,j,k) resp (i,j,k,l)
183 * trans_frac number of dihedrals in trans
184 * aver_angle average angle at each time frame
185 * dih all angles at each time frame
188 void make_histo(FILE *log,
189 int ndata, real data[], int npoints, int histo[],
190 real minx, real maxx);
192 * Make a histogram from data. The min and max of the data array can
193 * be determined (if minx == 0 and maxx == 0)
194 * and the index in the histogram is computed from
195 * ind = npoints/(max(data) - min(data))
197 * log write error output to this file
198 * ndata number of points in data
200 * npoints number of points in histogram
201 * histo histogram array. This is NOT set to zero, to allow you
202 * to add multiple histograms
203 * minx start of the histogram
204 * maxx end of the histogram
205 * if both are 0, these values are computed by the routine itself
208 void normalize_histo(int npoints, int histo[], real dx, real normhisto[]);
210 * Normalize a histogram so that the integral over the histo is 1
212 * npoints number of points in the histo array
213 * histo input histogram
214 * dx distance between points on the X-axis
215 * normhisto normalized output histogram
218 /* Routines from pp2shift (anadih.c etc.) */
220 void do_pp2shifts(FILE *fp, int nframes,
221 int nlist, t_dlist dlist[], real **dih);
223 gmx_bool has_dihedral(int Dih, t_dlist *dl);
225 t_dlist *mk_dlist(FILE *log,
226 t_atoms *atoms, int *nlist,
227 gmx_bool bPhi, gmx_bool bPsi, gmx_bool bChi, gmx_bool bHChi,
228 int maxchi, int r0, struct gmx_residuetype_t *rt);
230 void pr_dlist(FILE *fp, int nl, t_dlist dl[], real dt, int printtype,
231 gmx_bool bPhi, gmx_bool bPsi, gmx_bool bChi, gmx_bool bOmega, int maxchi);
233 int pr_trans(FILE *fp, int nl, t_dlist dl[], real dt, int Xi);
235 void mk_chi_lookup (int **lookup, int maxchi,
236 int nlist, t_dlist dlist[]);
238 void mk_multiplicity_lookup (int *multiplicity, int maxchi,
239 int nlist, t_dlist dlist[], int nangle);
241 void get_chi_product_traj (real **dih, int nframes,
242 int nlist, int maxchi, t_dlist dlist[],
243 real time[], int **lookup, int *multiplicity,
244 gmx_bool bRb, gmx_bool bNormalize,
245 real core_frac, gmx_bool bAll, const char *fnall,
246 const output_env_t oenv);
248 void print_one (const output_env_t oenv, const char *base,
250 const char *title, const char *ylabel, int nf,
251 real time[], real data[]);
253 /* Routines from g_hbond */
254 void analyse_corr(int n, real t[], real ct[], real nt[], real kt[],
255 real sigma_ct[], real sigma_nt[], real sigma_kt[],
256 real fit_start, real temp);
258 void compute_derivative(int nn, real x[], real y[], real dydx[]);