3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gromacs Runs On Most of All Computer Systems
52 /***********************************************
54 * A U T O C O R R E L A T I O N
56 ***********************************************/
58 extern real
LegendreP(real x
,unsigned long m
);
60 #define eacNormal (1<<0)
62 #define eacVector (1<<2)
63 #define eacRcross (1<<3 | eacVector)
64 #define eacP0 (1<<4 | eacVector)
65 #define eacP1 (1<<5 | eacVector)
66 #define eacP2 (1<<6 | eacVector)
67 #define eacP3 (1<<7 | eacVector)
68 #define eacP4 (1<<8 | eacVector)
69 #define eacIden (1<<9)
72 effnNONE
, effnEXP1
, effnEXP2
, effnEXP3
, effnVAC
,
73 effnEXP5
, effnEXP7
, effnEXP9
, effnERREST
, effnNR
76 /* must correspond with 'leg' g_chi.c:727 */
77 enum { edPhi
=0, edPsi
, edOmega
, edChi1
, edChi2
, edChi3
, edChi4
, edChi5
, edChi6
, edMax
};
79 enum { edPrintST
=0,edPrintRO
} ;
83 #define MAXCHI edMax-NONCHI
84 #define NROT 4 /* number of rotamers: 1=g(-), 2=t, 3=g(+), 0=other */
87 int minO
,minC
,H
,N
,C
,O
,Cn
[MAXCHI
+3];
88 } t_dihatms
; /* Cn[0]=N, Cn[1]=Ca, Cn[2]=Cb etc. */
93 int index
; /* Index for amino acids (histograms) */
94 int j0
[edMax
]; /* Index in dih array (phi angle is first...) */
99 real rot_occ
[edMax
][NROT
];
103 extern int nfp_ffn
[effnNR
];
105 extern const char *s_ffn
[effnNR
+2];
107 extern const char *longs_ffn
[effnNR
];
109 int sffn2effn(const char **sffn
);
110 /* Returns the ffn enum corresponding to the selected enum option in sffn */
112 extern t_pargs
*add_acf_pargs(int *npargs
,t_pargs
*pa
);
113 /* Add options for autocorr to the current set of options.
114 * *npargs must be initialised to the number of elements in pa,
115 * it will be incremented appropriately.
118 extern void cross_corr(int n
,real f
[],real g
[],real corr
[]);
119 /* Simple minded cross correlation algorithm */
121 extern real
fit_acf(int ncorr
,int fitfn
,const output_env_t oenv
,bool bVerbose
,
122 real tbeginfit
,real tendfit
,real dt
,real c1
[],real
*fit
);
123 /* Fit an ACF to a given function */
125 extern void do_autocorr(const char *fn
,const output_env_t oenv
,
127 int nframes
,int nitem
,real
**c1
,
128 real dt
,unsigned long mode
,bool bAver
);
129 /* Calls low_do_autocorr (see below). After calling add_acf_pargs */
131 extern void low_do_autocorr(const char *fn
,const output_env_t oenv
,
132 const char *title
, int nframes
,int nitem
,
133 int nout
,real
**c1
, real dt
,unsigned long mode
,
134 int nrestart
, bool bAver
,bool bNormalize
,
135 bool bVerbose
,real tbeginfit
,real tendfit
,
136 int nfitparm
,int nskip
);
138 * do_autocorr calculates autocorrelation functions for many things.
139 * It takes a 2 d array containing nitem arrays of length nframes
140 * for each item the ACF is calculated.
142 * A number of "modes" exist for computation of the ACF
144 * if (mode == eacNormal) {
145 * C(t) = < X (tau) * X (tau+t) >
147 * else if (mode == eacCos) {
148 * C(t) = < cos (X(tau) - X(tau+t)) >
150 * else if (mode == eacIden) { **not fully supported yet**
151 * C(t) = < (X(tau) == X(tau+t)) >
153 * else if (mode == eacVector) {
154 * C(t) = < X(tau) * X(tau+t)
156 * else if (mode == eacP1) {
157 * C(t) = < cos (X(tau) * X(tau+t) >
159 * else if (mode == eacP2) {
160 * C(t) = 1/2 * < 3 cos (X(tau) * X(tau+t) - 1 >
162 * else if (mode == eacRcross) {
163 * C(t) = < ( X(tau) * X(tau+t) )^2 >
166 * For modes eacVector, eacP1, eacP2 and eacRcross the input should be
167 * 3 x nframes long, where each triplet is taken as a 3D vector
169 * For mode eacCos inputdata must be in radians, not degrees!
171 * Other parameters are:
173 * fn is output filename (.xvg) where the correlation function(s) are printed
174 * title is the title in the output file
175 * nframes is the number of frames in the time series
176 * nitem is the number of items
177 * c1 is an array of dimension [ 0 .. nitem-1 ] [ 0 .. nframes-1 ]
178 * on output, this array is filled with the correlation function
180 * nrestart is the number of steps between restarts for direct ACFs
181 * (i.e. without FFT) When set to 1 all points are used as
182 * time origin for averaging
183 * dt is the time between frames
184 * bAver If set, all ndih C(t) functions are averaged into a single
186 * (bFour If set, will use fast fourier transform (FFT) for evaluating
187 * the ACF: removed option, now on the command line only)
188 * bNormalize If set, all ACFs will be normalized to start at 0
189 * nskip Determines whether steps a re skipped in the output
193 const char *name
; /* Description of the J coupling constant */
194 real A
,B
,C
; /* Karplus coefficients */
195 real offset
; /* Offset for dihedral angle in histogram (e.g. -M_PI/3) */
196 real Jc
; /* Resulting Jcoupling */
197 real Jcsig
; /* Standard deviation in Jc */
200 extern void calc_distribution_props(int nh
,int histo
[],
201 real start
,int nkkk
, t_karplus kkk
[],
203 /* This routine takes a dihedral distribution and calculates
204 * coupling constants and dihedral order parameters of it.
206 * nh is the number of points
207 * histo is the array of datapoints which is assumed to span
209 * start is the starting angle of the histogram, this can be either 0
211 * nkkk is the number of karplus sets (multiple coupling constants may be
212 * derived from a single angle)
213 * kkk are the constants for calculating J coupling constants using a
214 * Karplus equation according to
217 * J = A cos theta + B cos theta + C
219 * where theta is phi - offset (phi is the angle in the histogram)
220 * offset is subtracted from phi before substitution in the Karplus
222 * S2 is the resulting dihedral order parameter
227 /***********************************************
229 * F I T R O U T I N E S
231 ***********************************************/
232 extern void do_expfit(int ndata
,real c1
[],real dt
,
233 real begintimefit
,real endtimefit
);
235 extern void expfit(int n
, real x
[], real y
[], real Dy
[],
238 /* This procedure fits y=exp(a+bx) for n (x,y) pairs to determine a and b.
239 * The uncertainties in the y values must be in the vector Dy.
240 * The standard deviations of a and b, sa and sb, are also calculated.
242 * Routine from Computers in physics, 7(3) (1993), p. 280-285.
245 extern void ana_dih_trans(const char *fn_trans
,const char *fn_histo
,
246 real
**dih
,int nframes
,int nangles
,
247 const char *grpname
,real t0
,real dt
,bool bRb
,
248 const output_env_t oenv
);
250 * Analyse dihedral transitions, by counting transitions per dihedral
251 * and per frame. The total number of transitions is printed to
252 * stderr, as well as the average time between transitions.
254 * is wrapper to low_ana_dih_trans, which also passes in and out the
255 number of transitions per dihedral per residue. that uses struc dlist
256 which is not external, so pp2shift.h must be included.
258 * Dihedrals are supposed to be in either of three minima,
259 * (trans, gauche+, gauche-)
261 * fn_trans output file name for #transitions per timeframe
262 * fn_histo output file name for transition time histogram
263 * dih the actual dihedral angles
264 * nframes number of times frames
265 * nangles number of angles
266 * grpname a string for the header of plots
267 * t0,dt starting time resp. time between time frames
268 * bRb determines whether the polymer convention is used
272 extern void low_ana_dih_trans(bool bTrans
, const char *fn_trans
,
273 bool bHisto
, const char *fn_histo
, int maxchi
,
274 real
**dih
, int nlist
, t_dlist dlist
[],
275 int nframes
, int nangles
, const char *grpname
,
276 int xity
[], real t0
, real dt
, bool bRb
,
277 real core_frac
, const output_env_t oenv
);
278 /* as above but passes dlist so can copy occupancies into it, and xity[]
279 * (1..nangles, corresp to dih[this][], so can have non-3 multiplicity of
280 * rotamers. Also production of xvg output files is conditional
281 * and the fractional width of each rotamer can be set ie for a 3 fold
282 * dihedral with core_frac = 0.5 only the central 60 degrees is assigned
283 * to each rotamer, the rest goes to rotamer zero */
287 extern void read_ang_dih(const char *trj_fn
,
288 bool bAngles
,bool bSaveAll
,bool bRb
,bool bPBC
,
289 int maxangstat
,int angstat
[],
290 int *nframes
,real
**time
,
291 int isize
,atom_id index
[],
295 const output_env_t oenv
);
297 * Read a trajectory and calculate angles and dihedrals.
299 * trj_fn file name of trajectory
300 * tpb_fn file name of tpb file
301 * bAngles do we have to read angles or dihedrals
302 * bSaveAll do we have to store all in the dih array
303 * bRb do we have Ryckaert-Bellemans dihedrals (trans = 0)
304 * bPBC compute angles module 2 Pi
305 * maxangstat number of entries in distribution array
306 * angstat angle distribution
307 * *nframes number of frames read
308 * time simulation time at each time frame
309 * isize number of entries in the index, when angles 3*number of angles
310 * else 4*number of angles
311 * index atom numbers that define the angles or dihedrals
312 * (i,j,k) resp (i,j,k,l)
313 * trans_frac number of dihedrals in trans
314 * aver_angle average angle at each time frame
315 * dih all angles at each time frame
318 extern void make_histo(FILE *log
,
319 int ndata
,real data
[],int npoints
,int histo
[],
320 real minx
,real maxx
);
322 * Make a histogram from data. The min and max of the data array can
323 * be determined (if minx == 0 and maxx == 0)
324 * and the index in the histogram is computed from
325 * ind = npoints/(max(data) - min(data))
327 * log write error output to this file
328 * ndata number of points in data
330 * npoints number of points in histogram
331 * histo histogram array. This is NOT set to zero, to allow you
332 * to add multiple histograms
333 * minx start of the histogram
334 * maxx end of the histogram
335 * if both are 0, these values are computed by the routine itself
338 extern void normalize_histo(int npoints
,int histo
[],real dx
,real normhisto
[]);
340 * Normalize a histogram so that the integral over the histo is 1
342 * npoints number of points in the histo array
343 * histo input histogram
344 * dx distance between points on the X-axis
345 * normhisto normalized output histogram
348 extern real
fit_function(int eFitFn
,real
*parm
,real x
);
349 /* Returns the value of fit function eFitFn at x */
351 /* Use Levenberg-Marquardt method to fit to a nfitparm parameter exponential */
352 /* or to a transverse current autocorrelation function */
353 /* Or: "There is no KILL like OVERKILL", Dr. Ir. D. van der Spoel */
354 extern real
do_lmfit(int ndata
,real c1
[],real sig
[],real dt
,real
*x
,
355 real begintimefit
,real endtimefit
,const output_env_t oenv
,
356 bool bVerbose
, int eFitFn
,real fitparms
[],int fix
);
358 * If x == NULL, the timestep dt will be used to create a time axis.
359 * fix fixes fit parameter i at it's starting value, when the i'th bit
363 extern real
evaluate_integral(int n
,real x
[],real y
[],real dy
[],
364 real aver_start
,real
*stddev
);
365 /* Integrate data in y, and, if given, use dy as weighting
366 * aver_start should be set to a value where the function has
370 extern real
print_and_integrate(FILE *fp
,int n
,real dt
,
371 real c
[],real
*fit
,int nskip
);
372 /* Integrate the data in c[] from 0 to n using trapezium rule.
373 * If fp != NULL output is written to it
374 * nskip determines whether all elements are written to the output file
375 * (written when i % nskip == 0)
376 * If fit != NULL the fit is also written.
379 extern int get_acfnout(void);
380 /* Return the output length for the correlation function
381 * Works only AFTER do_auto_corr has been called!
384 extern int get_acffitfn(void);
385 /* Return the fit function type.
386 * Works only AFTER do_auto_corr has been called!
389 /* Routines from pp2shift (anadih.c etc.) */
391 extern void do_pp2shifts(FILE *fp
,int nframes
,
392 int nlist
,t_dlist dlist
[],real
**dih
);
394 extern bool has_dihedral(int Dih
,t_dlist
*dl
);
396 extern t_dlist
*mk_dlist(FILE *log
,
397 t_atoms
*atoms
, int *nlist
,
398 bool bPhi
, bool bPsi
, bool bChi
, bool bHChi
,
399 int maxchi
,int r0
,int naa
,char **aa
);
401 extern void pr_dlist(FILE *fp
,int nl
,t_dlist dl
[],real dt
, int printtype
,
402 bool bPhi
, bool bPsi
,bool bChi
,bool bOmega
, int maxchi
);
404 extern int pr_trans(FILE *fp
,int nl
,t_dlist dl
[],real dt
,int Xi
);
406 extern void mk_chi_lookup (int **lookup
, int maxchi
, real
**dih
,
407 int nlist
, t_dlist dlist
[]) ;
409 extern void mk_multiplicity_lookup (int *xity
, int maxchi
, real
**dih
,
410 int nlist
, t_dlist dlist
[],int nangle
) ;
412 extern void get_chi_product_traj (real
**dih
,int nframes
,int nangles
,
413 int nlist
,int maxchi
, t_dlist dlist
[],
414 real time
[], int **lookup
,int *xity
,
415 bool bRb
,bool bNormalize
,
416 real core_frac
,bool bAll
,const char *fnall
,
417 const output_env_t oenv
);
419 extern void print_one (const output_env_t oenv
, const char *base
,
421 const char *title
, const char *ylabel
,int nf
,
422 real time
[],real data
[]);
424 /* Routines from g_hbond */
425 extern void analyse_corr(int n
,real t
[],real ct
[],real nt
[],real kt
[],
426 real sigma_ct
[],real sigma_nt
[],real sigma_kt
[],
427 real fit_start
,real temp
,real smooth_tail_start
,
428 const output_env_t oenv
);
430 extern void compute_derivative(int nn
,real x
[],real y
[],real dydx
[]);