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
51 static const char *sepdvdlformat
=" %-30s V %12.5e dVdl %12.5e\n";
53 void calc_vir(FILE *fplog
,int nxf
,rvec x
[],rvec f
[],tensor vir
,
54 gmx_bool bScrewPBC
,matrix box
);
55 /* Calculate virial for nxf atoms, and add it to vir */
57 void f_calc_vir(FILE *fplog
,int i0
,int i1
,rvec x
[],rvec f
[],tensor vir
,
58 t_graph
*g
,rvec shift_vec
[]);
59 /* Calculate virial taking periodicity into account */
61 real
RF_excl_correction(FILE *fplog
,
62 const t_forcerec
*fr
,t_graph
*g
,
63 const t_mdatoms
*mdatoms
,const t_blocka
*excl
,
64 rvec x
[],rvec f
[],rvec
*fshift
,const t_pbc
*pbc
,
65 real lambda
,real
*dvdlambda
);
66 /* Calculate the reaction-field energy correction for this node:
67 * epsfac q_i q_j (k_rf r_ij^2 - c_rf)
68 * and force correction for all excluded pairs, including self pairs.
71 void calc_rffac(FILE *fplog
,int eel
,real eps_r
,real eps_rf
,
74 real
*kappa
,real
*krf
,real
*crf
);
75 /* Determine the reaction-field constants */
77 void init_generalized_rf(FILE *fplog
,
78 const gmx_mtop_t
*mtop
,const t_inputrec
*ir
,
80 /* Initialize the generalized reaction field parameters */
84 void make_wall_tables(FILE *fplog
,const output_env_t oenv
,
85 const t_inputrec
*ir
,const char *tabfn
,
86 const gmx_groups_t
*groups
,
89 real
do_walls(t_inputrec
*ir
,t_forcerec
*fr
,matrix box
,t_mdatoms
*md
,
90 rvec x
[],rvec f
[],real lambda
,real Vlj
[],t_nrnb
*nrnb
);
94 t_forcerec
*mk_forcerec(void);
96 #define GMX_MAKETABLES_FORCEUSER (1<<0)
97 #define GMX_MAKETABLES_14ONLY (1<<1)
99 t_forcetable
make_tables(FILE *fp
,const output_env_t oenv
,
100 const t_forcerec
*fr
, gmx_bool bVerbose
,
101 const char *fn
, real rtab
,int flags
);
102 /* Return tables for inner loops. When bVerbose the tables are printed
106 bondedtable_t
make_bonded_table(FILE *fplog
,char *fn
,int angle
);
107 /* Return a table for bonded interactions,
108 * angle should be: bonds 0, angles 1, dihedrals 2
111 /* Return a table for GB calculations */
112 t_forcetable
make_gb_table(FILE *out
,const output_env_t oenv
,
113 const t_forcerec
*fr
,
117 /* Read a table for AdResS Thermo Force calculations */
118 extern t_forcetable
make_atf_table(FILE *out
,const output_env_t oenv
,
119 const t_forcerec
*fr
,
123 void pr_forcerec(FILE *fplog
,t_forcerec
*fr
,t_commrec
*cr
);
126 forcerec_set_ranges(t_forcerec
*fr
,
127 int ncg_home
,int ncg_force
,
129 int natoms_force_constr
,int natoms_f_novirsum
);
130 /* Set the number of cg's and atoms for the force calculation */
132 gmx_bool
can_use_allvsall(const t_inputrec
*ir
, const gmx_mtop_t
*mtop
,
133 gmx_bool bPrintNote
,t_commrec
*cr
,FILE *fp
);
134 /* Returns if we can use all-vs-all loops.
135 * If bPrintNote==TRUE, prints a note, if necessary, to stderr
136 * and fp (if !=NULL) on the master node.
139 void init_forcerec(FILE *fplog
,
140 const output_env_t oenv
,
143 const t_inputrec
*ir
,
144 const gmx_mtop_t
*mtop
,
154 /* The Force rec struct must be created with mk_forcerec
155 * The gmx_booleans have the following meaning:
156 * bSetQ: Copy the charges [ only necessary when they change ]
157 * bMolEpot: Use the free energy stuff per molecule
158 * print_force >= 0: print forces for atoms with force >= print_force
161 void init_enerdata(int ngener
,int n_flambda
,gmx_enerdata_t
*enerd
);
162 /* Intializes the energy storage struct */
164 void destroy_enerdata(gmx_enerdata_t
*enerd
);
165 /* Free all memory associated with enerd */
167 void reset_enerdata(t_grpopts
*opts
,
168 t_forcerec
*fr
,gmx_bool bNS
,
169 gmx_enerdata_t
*enerd
,
171 /* Resets the energy data, if bNS=TRUE also zeros the long-range part */
173 void sum_epot(t_grpopts
*opts
,gmx_enerdata_t
*enerd
);
174 /* Locally sum the non-bonded potential energy terms */
176 void sum_dhdl(gmx_enerdata_t
*enerd
,double lambda
,t_inputrec
*ir
);
177 /* Sum the free energy contributions */
179 void update_forcerec(FILE *fplog
,t_forcerec
*fr
,matrix box
);
180 /* Updates parameters in the forcerec that are time dependent */
182 /* Compute the average C6 and C12 params for LJ corrections */
183 void set_avcsixtwelve(FILE *fplog
,t_forcerec
*fr
,
184 const gmx_mtop_t
*mtop
);
186 /* The state has changed */
187 #define GMX_FORCE_STATECHANGED (1<<0)
188 /* The box might have changed */
189 #define GMX_FORCE_DYNAMICBOX (1<<1)
190 /* Do neighbor searching */
191 #define GMX_FORCE_NS (1<<2)
192 /* Calculate bonded energies/forces */
193 #define GMX_FORCE_DOLR (1<<3)
194 /* Calculate long-range energies/forces */
195 #define GMX_FORCE_BONDED (1<<4)
196 /* Store long-range forces in a separate array */
197 #define GMX_FORCE_SEPLRF (1<<5)
198 /* Calculate non-bonded energies/forces */
199 #define GMX_FORCE_NONBONDED (1<<6)
200 /* Calculate forces (not only energies) */
201 #define GMX_FORCE_FORCES (1<<7)
202 /* Calculate the virial */
203 #define GMX_FORCE_VIRIAL (1<<8)
205 #define GMX_FORCE_DHDL (1<<9)
206 /* Normally one want all energy terms and forces */
207 #define GMX_FORCE_ALLFORCES (GMX_FORCE_BONDED | GMX_FORCE_NONBONDED | GMX_FORCE_FORCES)
209 void do_force(FILE *log
,t_commrec
*cr
,
210 t_inputrec
*inputrec
,
211 gmx_large_int_t step
,t_nrnb
*nrnb
,gmx_wallcycle_t wcycle
,
214 gmx_groups_t
*groups
,
215 matrix box
,rvec x
[],history_t
*hist
,
219 gmx_enerdata_t
*enerd
,t_fcdata
*fcd
,
220 real lambda
,t_graph
*graph
,
221 t_forcerec
*fr
,gmx_vsite_t
*vsite
,rvec mu_tot
,
222 double t
,FILE *field
,gmx_edsam_t ed
,
225 /* Communicate coordinates (if parallel).
226 * Do neighbor searching (if necessary).
228 * Communicate forces (if parallel).
229 * Spread forces for vsites (if present).
231 * f is always required.
238 gmx_groups_t
*groups
,
246 gmx_grppairener_t
*grppener
,
248 gmx_bool bDoLongRange
,
251 /* Call the neighborsearcher */
253 void do_force_lowlevel(FILE *fplog
,
254 gmx_large_int_t step
,
260 gmx_wallcycle_t wcycle
,
266 gmx_enerdata_t
*enerd
,
280 /* Call all the force routines */
286 #endif /* _force_h */