Velocities and forces for selections.
[gromacs/qmmm-gamess-us.git] / include / update.h
blobd060fe85634a2c9aebbd6429301364b79d26fd19
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 _update_h
37 #define _update_h
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
43 #include "typedefs.h"
44 #include "mshift.h"
45 #include "tgroup.h"
46 #include "network.h"
47 #include "force.h"
48 #include "pull.h"
49 #include "gmx_random.h"
51 #ifdef __cplusplus
52 extern "C" {
53 #endif
55 /* Abstract type for stochastic dynamics */
56 typedef struct gmx_update *gmx_update_t;
58 /* Initialize the stochastic dynamics struct */
59 extern gmx_update_t init_update(FILE *fplog,t_inputrec *ir);
61 /* Store the random state from sd in state */
62 extern void get_stochd_state(gmx_update_t sd,t_state *state);
64 /* Set the random in sd from state */
65 extern void set_stochd_state(gmx_update_t sd,t_state *state);
67 /* Store the box at step step
68 * as a reference state for simulations with box deformation.
70 extern void set_deform_reference_box(gmx_update_t upd,
71 gmx_large_int_t step,matrix box);
73 extern void update(FILE *fplog,
74 gmx_large_int_t step,
75 real *dvdlambda, /* FEP stuff */
76 t_inputrec *inputrec, /* input record and box stuff */
77 t_mdatoms *md,
78 t_state *state,
79 t_graph *graph,
80 rvec *f, /* forces on home particles */
81 bool bDoLR,
82 rvec *f_lr, /* LR forces, only with dDoLR */
83 t_fcdata *fcd,
84 t_idef *idef,
85 gmx_ekindata_t *ekind,
86 matrix *scale_tot,
87 t_commrec *cr,
88 t_nrnb *nrnb,
89 gmx_wallcycle_t wcycle,
90 gmx_update_t upd,
91 gmx_constr_t constr,
92 bool bCalcVir,
93 tensor vir_part,
94 bool bNEMD,
95 bool bInitStep);
96 /* Return TRUE if OK, FALSE in case of Shake Error */
98 extern void calc_ke_part(t_state *state,t_grpopts *opts,t_mdatoms *md,
99 gmx_ekindata_t *ekind,t_nrnb *nrnb);
101 * Compute the partial kinetic energy for home particles;
102 * will be accumulated in the calling routine.
103 * The tensor is
105 * Ekin = SUM(i) 0.5 m[i] v[i] (x) v[i]
107 * use v[i] = v[i] - u[i] when calculating temperature
109 * u must be accumulated already.
111 * Now also computes the contribution of the kinetic energy to the
112 * free energy
116 extern void
117 init_ekinstate(ekinstate_t *ekinstate,t_inputrec *ir);
119 extern void
120 update_ekinstate(ekinstate_t *ekinstate,gmx_ekindata_t *ekind);
122 extern void
123 restore_ekinstate_from_state(t_commrec *cr,
124 gmx_ekindata_t *ekind,ekinstate_t *ekinstate);
126 extern void berendsen_tcoupl(t_grpopts *opts,gmx_ekindata_t *ekind,real dt);
128 extern void nosehoover_tcoupl(t_grpopts *opts,gmx_ekindata_t *ekind,real dt,
129 real xi[],double ixi[]);
130 /* Compute temperature scaling. For Nose-Hoover it is done in update. */
132 extern real nosehoover_energy(t_grpopts *opts,gmx_ekindata_t *ekind,
133 real *xi,double *ixi);
134 /* Returns the Nose-Hoover contribution to the conserved energy */
136 extern void vrescale_tcoupl(t_grpopts *opts,gmx_ekindata_t *ekind,real dt,
137 double therm_integral[],
138 gmx_rng_t rng);
139 /* Compute temperature scaling. For V-rescale it is done in update. */
141 extern real vrescale_energy(t_grpopts *opts,double therm_integral[]);
142 /* Returns the V-rescale contribution to the conserved energy */
144 /* Set reference temp for simulated annealing at time t*/
145 extern void update_annealing_target_temp(t_grpopts *opts,real t);
147 extern real calc_temp(real ekin,real nrdf);
148 /* Calculate the temperature */
150 extern real calc_pres(int ePBC,int nwall,matrix box,
151 tensor ekin,tensor vir,tensor pres,real Elr);
152 /* Calculate the pressure tensor, returns the scalar pressure.
153 * The unit of pressure is bar, If Elr != 0
154 * a long range correction based on Ewald/PPPM is made (see c-code)
157 extern void parrinellorahman_pcoupl(FILE *fplog,gmx_large_int_t step,
158 t_inputrec *ir,real dt,tensor pres,
159 tensor box,tensor box_rel,tensor boxv,
160 tensor M,matrix mu,
161 bool bFirstStep);
163 extern void berendsen_pcoupl(FILE *fplog,gmx_large_int_t step,
164 t_inputrec *ir,real dt,tensor pres,matrix box,
165 matrix mu);
168 extern void berendsen_pscale(t_inputrec *ir,matrix mu,
169 matrix box,matrix box_rel,
170 int start,int nr_atoms,
171 rvec x[],unsigned short cFREEZE[],
172 t_nrnb *nrnb);
174 extern void correct_ekin(FILE *log,int start,int end,rvec v[],
175 rvec vcm,real mass[],real tmass,tensor ekin);
176 /* Correct ekin for vcm */
178 #ifdef __cplusplus
180 #endif
182 #endif /* _update_h */