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