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
45 #include "gmx_random.h"
52 /* Abstract type for stochastic dynamics */
53 typedef struct gmx_update
*gmx_update_t
;
55 /* Initialize the stochastic dynamics struct */
56 gmx_update_t
init_update(FILE *fplog
,t_inputrec
*ir
);
58 /* Store the random state from sd in state */
59 void get_stochd_state(gmx_update_t sd
,t_state
*state
);
61 /* Set the random in sd from state */
62 void set_stochd_state(gmx_update_t sd
,t_state
*state
);
64 /* Store the box at step step
65 * as a reference state for simulations with box deformation.
67 void set_deform_reference_box(gmx_update_t upd
,
68 gmx_large_int_t step
,matrix box
);
70 void update_tcouple(FILE *fplog
,
74 gmx_ekindata_t
*ekind
,
75 gmx_wallcycle_t wcycle
,
81 void update_pcouple(FILE *fplog
,
87 gmx_wallcycle_t wcycle
,
91 void update_coords(FILE *fplog
,
93 t_inputrec
*inputrec
, /* input record and box stuff */
96 rvec
*f
, /* forces on home particles */
100 gmx_ekindata_t
*ekind
,
102 gmx_wallcycle_t wcycle
,
106 t_commrec
*cr
, /* these shouldn't be here -- need to think about it */
111 /* Return TRUE if OK, FALSE in case of Shake Error */
113 void update_constraints(FILE *fplog
,
114 gmx_large_int_t step
,
115 real
*dvdlambda
, /* FEP stuff */
116 t_inputrec
*inputrec
, /* input record and box stuff */
117 gmx_ekindata_t
*ekind
,
121 rvec force
[], /* forces on home particles */
127 gmx_wallcycle_t wcycle
,
135 /* Return TRUE if OK, FALSE in case of Shake Error */
137 void update_box(FILE *fplog
,
138 gmx_large_int_t step
,
139 t_inputrec
*inputrec
, /* input record and box stuff */
143 rvec force
[], /* forces on home particles */
147 gmx_wallcycle_t wcycle
,
150 gmx_bool bFirstHalf
);
151 /* Return TRUE if OK, FALSE in case of Shake Error */
153 void calc_ke_part(t_state
*state
,t_grpopts
*opts
,t_mdatoms
*md
,
154 gmx_ekindata_t
*ekind
,t_nrnb
*nrnb
,gmx_bool bEkinAveVel
, gmx_bool bSaveOld
);
156 * Compute the partial kinetic energy for home particles;
157 * will be accumulated in the calling routine.
160 * Ekin = SUM(i) 0.5 m[i] v[i] (x) v[i]
162 * use v[i] = v[i] - u[i] when calculating temperature
164 * u must be accumulated already.
166 * Now also computes the contribution of the kinetic energy to the
172 init_ekinstate(ekinstate_t
*ekinstate
,const t_inputrec
*ir
);
175 update_ekinstate(ekinstate_t
*ekinstate
,gmx_ekindata_t
*ekind
);
178 restore_ekinstate_from_state(t_commrec
*cr
,
179 gmx_ekindata_t
*ekind
,ekinstate_t
*ekinstate
);
181 void berendsen_tcoupl(t_inputrec
*ir
,gmx_ekindata_t
*ekind
,real dt
);
183 void nosehoover_tcoupl(t_grpopts
*opts
,gmx_ekindata_t
*ekind
,real dt
,
184 double xi
[],double vxi
[],t_extmass
*MassQ
);
186 t_state
*init_bufstate(const t_state
*template_state
);
188 void destroy_bufstate(t_state
*state
);
190 void trotter_update(t_inputrec
*ir
, gmx_large_int_t step
, gmx_ekindata_t
*ekind
,
191 gmx_enerdata_t
*enerd
, t_state
*state
, tensor vir
, t_mdatoms
*md
,
192 t_extmass
*MassQ
, int **trotter_seqlist
, int trotter_seqno
);
194 int **init_npt_vars(t_inputrec
*ir
, t_state
*state
, t_extmass
*Mass
, gmx_bool bTrotter
);
196 real
NPT_energy(t_inputrec
*ir
, t_state
*state
, t_extmass
*MassQ
);
197 /* computes all the pressure/tempertature control energy terms to get a conserved energy */
199 void NBaroT_trotter(t_grpopts
*opts
, real dt
,
200 double xi
[],double vxi
[],real
*veta
, t_extmass
*MassQ
);
202 void vrescale_tcoupl(t_inputrec
*ir
,gmx_ekindata_t
*ekind
,real dt
,
203 double therm_integral
[],
205 /* Compute temperature scaling. For V-rescale it is done in update. */
207 real
vrescale_energy(t_grpopts
*opts
,double therm_integral
[]);
208 /* Returns the V-rescale contribution to the conserved energy */
210 void rescale_velocities(gmx_ekindata_t
*ekind
,t_mdatoms
*mdatoms
,
211 int start
,int end
,rvec v
[]);
212 /* Rescale the velocities with the scaling factor in ekind */
214 void update_annealing_target_temp(t_grpopts
*opts
,real t
);
215 /* Set reference temp for simulated annealing at time t*/
217 real
calc_temp(real ekin
,real nrdf
);
218 /* Calculate the temperature */
220 real
calc_pres(int ePBC
,int nwall
,matrix box
,
221 tensor ekin
,tensor vir
,tensor pres
,real Elr
);
222 /* Calculate the pressure tensor, returns the scalar pressure.
223 * The unit of pressure is bar, If Elr != 0
224 * a long range correction based on Ewald/PPPM is made (see c-code)
227 void parrinellorahman_pcoupl(FILE *fplog
,gmx_large_int_t step
,
228 t_inputrec
*ir
,real dt
,tensor pres
,
229 tensor box
,tensor box_rel
,tensor boxv
,
231 gmx_bool bFirstStep
);
233 void berendsen_pcoupl(FILE *fplog
,gmx_large_int_t step
,
234 t_inputrec
*ir
,real dt
,tensor pres
,matrix box
,
238 void berendsen_pscale(t_inputrec
*ir
,matrix mu
,
239 matrix box
,matrix box_rel
,
240 int start
,int nr_atoms
,
241 rvec x
[],unsigned short cFREEZE
[],
244 void correct_ekin(FILE *log
,int start
,int end
,rvec v
[],
245 rvec vcm
,real mass
[],real tmass
,tensor ekin
);
246 /* Correct ekin for vcm */
253 #endif /* _update_h */