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);
92 t_forcerec *mk_forcerec(void);
94 #define GMX_MAKETABLES_FORCEUSER (1<<0)
95 #define GMX_MAKETABLES_14ONLY (1<<1)
97 t_forcetable make_tables(FILE *fp,const output_env_t oenv,
98 const t_forcerec *fr, gmx_bool bVerbose,
99 const char *fn, real rtab,int flags);
100 /* Return tables for inner loops. When bVerbose the tables are printed
104 bondedtable_t make_bonded_table(FILE *fplog,char *fn,int angle);
105 /* Return a table for bonded interactions,
106 * angle should be: bonds 0, angles 1, dihedrals 2
109 /* Return a table for GB calculations */
110 t_forcetable make_gb_table(FILE *out,const output_env_t oenv,
111 const t_forcerec *fr,
115 /* Read a table for AdResS Thermo Force calculations */
116 extern t_forcetable make_atf_table(FILE *out,const output_env_t oenv,
117 const t_forcerec *fr,
121 void pr_forcerec(FILE *fplog,t_forcerec *fr,t_commrec *cr);
124 forcerec_set_ranges(t_forcerec *fr,
125 int ncg_home,int ncg_force,
127 int natoms_force_constr,int natoms_f_novirsum);
128 /* Set the number of cg's and atoms for the force calculation */
130 gmx_bool can_use_allvsall(const t_inputrec *ir, const gmx_mtop_t *mtop,
131 gmx_bool bPrintNote,t_commrec *cr,FILE *fp);
132 /* Returns if we can use all-vs-all loops.
133 * If bPrintNote==TRUE, prints a note, if necessary, to stderr
134 * and fp (if !=NULL) on the master node.
137 void init_forcerec(FILE *fplog,
138 const output_env_t oenv,
141 const t_inputrec *ir,
142 const gmx_mtop_t *mtop,
152 /* The Force rec struct must be created with mk_forcerec
153 * The gmx_booleans have the following meaning:
154 * bSetQ: Copy the charges [ only necessary when they change ]
155 * bMolEpot: Use the free energy stuff per molecule
156 * print_force >= 0: print forces for atoms with force >= print_force
159 void init_enerdata(int ngener,int n_lambda,gmx_enerdata_t *enerd);
160 /* Intializes the energy storage struct */
162 void destroy_enerdata(gmx_enerdata_t *enerd);
163 /* Free all memory associated with enerd */
165 void reset_enerdata(t_grpopts *opts,
166 t_forcerec *fr,gmx_bool bNS,
167 gmx_enerdata_t *enerd,
169 /* Resets the energy data, if bNS=TRUE also zeros the long-range part */
171 void sum_epot(t_grpopts *opts,gmx_enerdata_t *enerd);
172 /* Locally sum the non-bonded potential energy terms */
174 void sum_dhdl(gmx_enerdata_t *enerd,real *lambda,t_lambda *fepvals);
175 /* Sum the free energy contributions */
177 void update_forcerec(FILE *fplog,t_forcerec *fr,matrix box);
178 /* Updates parameters in the forcerec that are time dependent */
180 /* Compute the average C6 and C12 params for LJ corrections */
181 void set_avcsixtwelve(FILE *fplog,t_forcerec *fr,
182 const gmx_mtop_t *mtop);
184 /* The state has changed */
185 #define GMX_FORCE_STATECHANGED (1<<0)
186 /* The box might have changed */
187 #define GMX_FORCE_DYNAMICBOX (1<<1)
188 /* Do neighbor searching */
189 #define GMX_FORCE_NS (1<<2)
190 /* Calculate bonded energies/forces */
191 #define GMX_FORCE_DOLR (1<<3)
192 /* Calculate long-range energies/forces */
193 #define GMX_FORCE_BONDED (1<<4)
194 /* Store long-range forces in a separate array */
195 #define GMX_FORCE_SEPLRF (1<<5)
196 /* Calculate non-bonded energies/forces */
197 #define GMX_FORCE_NONBONDED (1<<6)
198 /* Calculate forces (not only energies) */
199 #define GMX_FORCE_FORCES (1<<7)
200 /* Calculate the virial */
201 #define GMX_FORCE_VIRIAL (1<<8)
203 #define GMX_FORCE_DHDL (1<<9)
204 /* Normally one want all energy terms and forces */
205 #define GMX_FORCE_ALLFORCES (GMX_FORCE_BONDED | GMX_FORCE_NONBONDED | GMX_FORCE_FORCES)
207 void do_force(FILE *log,t_commrec *cr,
208 t_inputrec *inputrec,
209 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
212 gmx_groups_t *groups,
213 matrix box,rvec x[],history_t *hist,
217 gmx_enerdata_t *enerd,t_fcdata *fcd,
218 real *lambda,t_graph *graph,
219 t_forcerec *fr,gmx_vsite_t *vsite,rvec mu_tot,
220 double t,FILE *field,gmx_edsam_t ed,
223 /* Communicate coordinates (if parallel).
224 * Do neighbor searching (if necessary).
226 * Communicate forces (if parallel).
227 * Spread forces for vsites (if present).
229 * f is always required.
236 gmx_groups_t *groups,
244 gmx_grppairener_t *grppener,
246 gmx_bool bDoLongRange,
249 /* Call the neighborsearcher */
251 void do_force_lowlevel(FILE *fplog,
252 gmx_large_int_t step,
258 gmx_wallcycle_t wcycle,
264 gmx_enerdata_t *enerd,
279 /* Call all the force routines */
285 #endif /* _force_h */