4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
12 * Copyright (c) 1991-1999
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
17 * GROMACS: A message-passing parallel molecular dynamics implementation
18 * H.J.C. Berendsen, D. van der Spoel and R. van Drunen
19 * Comp. Phys. Comm. 91, 43-56 (1995)
21 * Also check out our WWW page:
22 * http://md.chem.rug.nl/~gmx
27 * Great Red Oystrich Makes All Chemists Sane
33 static char *SRCID_do_gct_h
= "$Id$";
40 enum { eoPres
, eoEpot
, eoVir
, eoDist
, eoMu
, eoForce
, eoFx
, eoFy
, eoFz
,
41 eoPx
, eoPy
, eoPz
, eoObsNR
,
42 eoPolarizability
=eoObsNR
, eoDipole
,
43 eoMemory
, eoInter
, eoUseVirial
, eoNR
};
44 extern char *eoNames
[eoNR
];
47 int at_i
,at_j
; /* Atom type # for i and j */
48 int eObs
; /* Observable to couple to */
49 bool bPrint
; /* Does this struct have to be printed */
50 real c6
,c12
; /* Actual value of params */
51 real xi_6
,xi_12
; /* Constants for coupling C6 and C12 */
55 int at_i
,at_j
; /* Atom type # for i and j */
56 int eObs
; /* Observable to couple to */
57 bool bPrint
; /* Does this struct have to be printed */
58 real a
,b
,c
; /* Actual value of params */
59 real xi_a
,xi_b
,xi_c
; /* Constants for coupling A, B and C */
63 int at_i
; /* Atom type */
64 int eObs
; /* Observable to couple to */
65 bool bPrint
; /* Does this struct have to be printed */
66 real Q
; /* Actual value of charge */
67 real xi_Q
; /* Constant for coupling Q */
71 int type
; /* Type number in the iparams struct */
72 int eObs
; /* Observable to couple to */
73 t_iparams xi
; /* Parameters that need to be changed */
78 real act_value
[eoObsNR
];
79 real av_value
[eoObsNR
];
80 real ref_value
[eoObsNR
];
81 bool bObsUsed
[eoObsNR
];
94 extern void write_gct(char *fn
,t_coupl_rec
*tcr
,t_idef
*idef
);
96 extern void read_gct(char *fn
,t_coupl_rec
*tcr
);
98 extern void comm_tcr(FILE *log
,t_commrec
*cr
,t_coupl_rec
**tcr
);
100 extern void copy_ff(t_coupl_rec
*tcr
,t_forcerec
*fr
,t_mdatoms
*md
,
103 extern t_coupl_rec
*init_coupling(FILE *log
,int nfile
,t_filenm fnm
[],
104 t_commrec
*cr
,t_forcerec
*fr
,t_mdatoms
*md
,
107 extern void calc_force(int natom
,rvec f
[],rvec fff
[]);
109 extern void calc_f_dev(int natoms
,real charge
[],rvec x
[],rvec f
[],
110 t_idef
*idef
,real
*xiH
,real
*xiS
);
112 extern void do_coupling(FILE *log
,int nfile
,t_filenm fnm
[],
113 t_coupl_rec
*tcr
,real t
,int step
,real ener
[],
114 t_forcerec
*fr
,t_inputrec
*ir
,bool bMaster
,
115 t_mdatoms
*md
,t_idef
*idef
,real mu_aver
,int nmols
,
116 t_commrec
*cr
,matrix box
,tensor virial
,
117 tensor pres
,rvec mu_tot
,
118 rvec x
[],rvec f
[],bool bDoIt
);