Fixed segfaults with empty selections.
[gromacs/qmmm-gamess-us.git] / include / force.h
blob256859066fe982182912fd30d686c66f5a2298fc
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
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
32 * And Hey:
33 * Gromacs Runs On Most of All Computer Systems
36 #ifndef _force_h
37 #define _force_h
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
44 #include "typedefs.h"
45 #include "pbc.h"
46 #include "network.h"
47 #include "tgroup.h"
48 #include "vsite.h"
49 #include "genborn.h"
51 static const char *sepdvdlformat=" %-30s V %12.5e dVdl %12.5e\n";
53 extern void calc_vir(FILE *fplog,int nxf,rvec x[],rvec f[],tensor vir,
54 bool bScrewPBC,matrix box);
55 /* Calculate virial for nxf atoms, and add it to vir */
57 extern 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 extern 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 extern void calc_rffac(FILE *fplog,int eel,real eps_r,real eps_rf,
72 real Rc,real Temp,
73 real zsq,matrix box,
74 real *kappa,real *krf,real *crf);
75 /* Determine the reaction-field constants */
77 extern void init_generalized_rf(FILE *fplog,
78 const gmx_mtop_t *mtop,const t_inputrec *ir,
79 t_forcerec *fr);
80 /* Initialize the generalized reaction field parameters */
83 /* In wall.c */
84 extern 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,
87 t_forcerec *fr);
89 extern 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 extern t_forcerec *mk_forcerec(void);
96 #define GMX_MAKETABLES_FORCEUSER (1<<0)
97 #define GMX_MAKETABLES_14ONLY (1<<1)
99 extern t_forcetable make_tables(FILE *fp,const output_env_t oenv,
100 const t_forcerec *fr, bool bVerbose,
101 const char *fn, real rtab,int flags);
102 /* Return tables for inner loops. When bVerbose the tables are printed
103 * to .xvg files
106 extern 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 extern t_forcetable make_gb_table(FILE *out,const output_env_t oenv,
113 const t_forcerec *fr,
114 const char *fn,
115 real rtab);
117 extern void pr_forcerec(FILE *fplog,t_forcerec *fr,t_commrec *cr);
119 extern void
120 forcerec_set_ranges(t_forcerec *fr,
121 int ncg_home,int ncg_force,
122 int natoms_force,int natoms_f_novirsum);
123 /* Set the number of cg's and atoms for the force calculation */
125 extern void init_forcerec(FILE *fplog,
126 const output_env_t oenv,
127 t_forcerec *fr,
128 t_fcdata *fcd,
129 const t_inputrec *ir,
130 const gmx_mtop_t *mtop,
131 const t_commrec *cr,
132 matrix box,
133 bool bMolEpot,
134 const char *tabfn,
135 const char *tabpfn,
136 const char *tabbfn,
137 bool bNoSolvOpt,
138 real print_force);
139 /* The Force rec struct must be created with mk_forcerec
140 * The booleans have the following meaning:
141 * bSetQ: Copy the charges [ only necessary when they change ]
142 * bMolEpot: Use the free energy stuff per molecule
143 * print_force >= 0: print forces for atoms with force >= print_force
146 extern void init_enerdata(int ngener,int n_flambda,gmx_enerdata_t *enerd);
147 /* Intializes the energy storage struct */
149 extern void destroy_enerdata(gmx_enerdata_t *enerd);
150 /* Free all memory associated with enerd */
152 extern void reset_enerdata(t_grpopts *opts,
153 t_forcerec *fr,bool bNS,
154 gmx_enerdata_t *enerd,
155 bool bMaster);
156 /* Resets the energy data, if bNS=TRUE also zeros the long-range part */
158 extern void sum_epot(t_grpopts *opts,gmx_enerdata_t *enerd);
159 /* Locally sum the non-bonded potential energy terms */
161 extern void sum_dhdl(gmx_enerdata_t *enerd,double lambda,t_inputrec *ir);
162 /* Sum the free energy contributions */
164 extern void update_forcerec(FILE *fplog,t_forcerec *fr,matrix box);
165 /* Updates parameters in the forcerec that are time dependent */
167 /* Compute the average C6 and C12 params for LJ corrections */
168 extern void set_avcsixtwelve(FILE *fplog,t_forcerec *fr,
169 const gmx_mtop_t *mtop);
171 /* The state has changed */
172 #define GMX_FORCE_STATECHANGED (1<<0)
173 /* The box might have changed */
174 #define GMX_FORCE_DYNAMICBOX (1<<1)
175 /* Do neighbor searching */
176 #define GMX_FORCE_NS (1<<2)
177 /* Calculate bonded energies/forces */
178 #define GMX_FORCE_DOLR (1<<3)
179 /* Calculate long-range energies/forces */
180 #define GMX_FORCE_BONDED (1<<4)
181 /* Store long-range forces in a separate array */
182 #define GMX_FORCE_SEPLRF (1<<5)
183 /* Calculate non-bonded energies/forces */
184 #define GMX_FORCE_NONBONDED (1<<6)
185 /* Calculate forces (not only energies) */
186 #define GMX_FORCE_FORCES (1<<7)
187 /* Calculate the virial */
188 #define GMX_FORCE_VIRIAL (1<<8)
189 /* Calculate dHdl */
190 #define GMX_FORCE_DHDL (1<<9)
191 /* Normally one want all energy terms and forces */
192 #define GMX_FORCE_ALLFORCES (GMX_FORCE_BONDED | GMX_FORCE_NONBONDED | GMX_FORCE_FORCES)
194 extern void do_force(FILE *log,t_commrec *cr,
195 t_inputrec *inputrec,
196 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
197 gmx_localtop_t *top,
198 gmx_mtop_t *mtop,
199 gmx_groups_t *groups,
200 matrix box,rvec x[],history_t *hist,
201 rvec f[],
202 tensor vir_force,
203 t_mdatoms *mdatoms,
204 gmx_enerdata_t *enerd,t_fcdata *fcd,
205 real lambda,t_graph *graph,
206 t_forcerec *fr,gmx_vsite_t *vsite,rvec mu_tot,
207 double t,FILE *field,gmx_edsam_t ed,
208 bool bBornRadii,
209 int flags);
210 /* Communicate coordinates (if parallel).
211 * Do neighbor searching (if necessary).
212 * Calculate forces.
213 * Communicate forces (if parallel).
214 * Spread forces for vsites (if present).
216 * f is always required.
219 extern void ns(FILE *fplog,
220 t_forcerec *fr,
221 rvec x[],
222 matrix box,
223 gmx_groups_t *groups,
224 t_grpopts *opts,
225 gmx_localtop_t *top,
226 t_mdatoms *md,
227 t_commrec *cr,
228 t_nrnb *nrnb,
229 real lambda,
230 real *dvdlambda,
231 gmx_grppairener_t *grppener,
232 bool bFillGrid,
233 bool bDoLongRange,
234 bool bDoForces,
235 rvec *f);
236 /* Call the neighborsearcher */
238 extern void do_force_lowlevel(FILE *fplog,
239 gmx_large_int_t step,
240 t_forcerec *fr,
241 t_inputrec *ir,
242 t_idef *idef,
243 t_commrec *cr,
244 t_nrnb *nrnb,
245 gmx_wallcycle_t wcycle,
246 t_mdatoms *md,
247 t_grpopts *opts,
248 rvec x[],
249 history_t *hist,
250 rvec f[],
251 gmx_enerdata_t *enerd,
252 t_fcdata *fcd,
253 gmx_mtop_t *mtop,
254 gmx_localtop_t *top,
255 gmx_genborn_t *born,
256 t_atomtypes *atype,
257 bool bBornRadii,
258 matrix box,
259 real lambda,
260 t_graph *graph,
261 t_blocka *excl,
262 rvec mu_tot[2],
263 int flags,
264 float *cycles_pme);
265 /* Call all the force routines */
268 #endif /* _force_h */