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 * GRoups of Organic Molecules in ACtion for Science
43 typedef real rvec5
[5];
45 /* Distance restraining stuff */
47 int dr_weighting
; /* Weighting of pairs in one restraint */
48 bool dr_bMixed
; /* Use sqrt of the instantaneous times *
49 * the time averaged violation */
50 real dr_fc
; /* Force constant for disres, *
51 * which is multiplied by a (possibly) *
52 * different factor for each restraint */
53 real dr_tau
; /* Time constant for disres */
54 real ETerm
; /* multiplication factor for time averaging */
55 real ETerm1
; /* 1 - ETerm1 */
56 real exp_min_t_tau
; /* Factor for slowly switching on the force */
57 int nres
; /* The number of distance restraints */
58 int npair
; /* The number of distance restraint pairs */
59 real sumviol
; /* The sum of violations */
60 real
*rt
; /* The calculated instantaneous distance (npr) */
61 real
*rm3tav
; /* The calculated time averaged distance (npr) */
62 real
*Rtl_6
; /* The calculated instantaneous r^-6 (nr) */
63 real
*Rt_6
; /* The calculated inst. ens. averaged r^-6 (nr) */
64 real
*Rtav_6
; /* The calculated time and ens. averaged r^-6 (nr) */
65 int nsystems
; /* The number of systems for ensemble averaging */
67 MPI_Comm mpi_comm_ensemble
; /* For ensemble averaging */
72 /* Orientation restraining stuff */
74 real fc
; /* Force constant for the restraints */
75 real edt
; /* Multiplication factor for time averaging */
76 real edt_1
; /* 1 - edt */
77 real exp_min_t_tau
; /* Factor for slowly switching on the force */
78 int nr
; /* The number of orientation restraints */
79 int nex
; /* The number of experiments */
80 int nref
; /* The number of atoms for the fit */
81 real
*mref
; /* The masses of the reference atoms */
82 rvec
*xref
; /* The reference coordinates for the fit (nref) */
83 rvec
*xtmp
; /* Temporary array for fitting (nref) */
84 matrix R
; /* Rotation matrix to rotate to the reference coor. */
85 tensor
*S
; /* Array of order tensors for each experiment (nexp) */
86 rvec5
*Dinsl
; /* The order matrix D for all restraints (nr x 5) */
87 rvec5
*Dins
; /* The ensemble averaged D (nr x 5) */
88 rvec5
*Dtav
; /* The time and ensemble averaged D (nr x 5) */
89 real
*oinsl
; /* The calculated instantaneous orientations */
90 real
*oins
; /* The calculated emsemble averaged orientations */
91 real
*otav
; /* The calculated time and ensemble averaged orient. */
92 real rmsdev
; /* The weighted (using kfac) RMS deviation */
93 rvec5
*tmp
; /* An array of temporary 5-vectors (nex); */
94 real
***TMP
; /* An array of temporary 5x5 matrices (nex); */
95 real
*eig
; /* Eigenvalues/vectors, for output only (nex x 12) */
99 * Data struct used in the force calculation routines
100 * for storing the tables for bonded interactions and
101 * for storing information which is needed in following steps
102 * (for instance for time averaging in distance retraints)
103 * or for storing output, since force routines only return the potential.
106 bondedtable_t
*bondtab
;
107 bondedtable_t
*angletab
;
108 bondedtable_t
*dihtab
;