Revertiing my fake merge to make history consistent
[gromacs/adressmacs.git] / include / update.h
blob382ce7237f857187af5c720e381336505f7381c5
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
44 #include "typedefs.h"
45 #include "mshift.h"
46 #include "tgroup.h"
47 #include "network.h"
48 #include "force.h"
49 #include "pull.h"
50 #include "gmx_random.h"
51 #include "maths.h"
53 #ifdef __cplusplus
54 extern "C" {
55 #endif
57 /* Abstract type for stochastic dynamics */
58 typedef struct gmx_update *gmx_update_t;
60 /* Initialize the stochastic dynamics struct */
61 extern gmx_update_t init_update(FILE *fplog,t_inputrec *ir);
63 /* Store the random state from sd in state */
64 extern void get_stochd_state(gmx_update_t sd,t_state *state);
66 /* Set the random in sd from state */
67 extern void set_stochd_state(gmx_update_t sd,t_state *state);
69 /* Store the box at step step
70 * as a reference state for simulations with box deformation.
72 extern void set_deform_reference_box(gmx_update_t upd,
73 gmx_large_int_t step,matrix box);
75 extern void update_tcouple(FILE *fplog,
76 gmx_large_int_t step,
77 t_inputrec *inputrec,
78 t_state *state,
79 gmx_ekindata_t *ekind,
80 gmx_wallcycle_t wcycle,
81 gmx_update_t upd,
82 t_extmass *MassQ,
83 t_mdatoms *md
86 extern void update_pcouple(FILE *fplog,
87 gmx_large_int_t step,
88 t_inputrec *inputrec,
89 t_state *state,
90 matrix pcoupl_mu,
91 matrix M,
92 gmx_wallcycle_t wcycle,
93 gmx_update_t upd,
94 bool bInitStep);
96 extern void update_coords(FILE *fplog,
97 gmx_large_int_t step,
98 t_inputrec *inputrec, /* input record and box stuff */
99 t_mdatoms *md,
100 t_state *state,
101 rvec *f, /* forces on home particles */
102 bool bDoLR,
103 rvec *f_lr,
104 t_fcdata *fcd,
105 gmx_ekindata_t *ekind,
106 matrix M,
107 gmx_wallcycle_t wcycle,
108 gmx_update_t upd,
109 bool bInitStep,
110 int bUpdatePart,
111 t_commrec *cr, /* these shouldn't be here -- need to think about it */
112 t_nrnb *nrnb,
113 gmx_constr_t constr,
114 t_idef *idef);
116 /* Return TRUE if OK, FALSE in case of Shake Error */
118 extern void update_constraints(FILE *fplog,
119 gmx_large_int_t step,
120 real *dvdlambda, /* FEP stuff */
121 t_inputrec *inputrec, /* input record and box stuff */
122 gmx_ekindata_t *ekind,
123 t_mdatoms *md,
124 t_state *state,
125 t_graph *graph,
126 rvec force[], /* forces on home particles */
127 t_idef *idef,
128 tensor vir_part,
129 tensor vir,
130 t_commrec *cr,
131 t_nrnb *nrnb,
132 gmx_wallcycle_t wcycle,
133 gmx_update_t upd,
134 gmx_constr_t constr,
135 bool bInitStep,
136 bool bFirstHalf,
137 bool bCalcVir,
138 real vetanew);
140 /* Return TRUE if OK, FALSE in case of Shake Error */
142 extern void update_box(FILE *fplog,
143 gmx_large_int_t step,
144 t_inputrec *inputrec, /* input record and box stuff */
145 t_mdatoms *md,
146 t_state *state,
147 t_graph *graph,
148 rvec force[], /* forces on home particles */
149 matrix *scale_tot,
150 matrix pcoupl_mu,
151 t_nrnb *nrnb,
152 gmx_wallcycle_t wcycle,
153 gmx_update_t upd,
154 bool bInitStep,
155 bool bFirstHalf);
156 /* Return TRUE if OK, FALSE in case of Shake Error */
158 extern void calc_ke_part(t_state *state,t_grpopts *opts,t_mdatoms *md,
159 gmx_ekindata_t *ekind,t_nrnb *nrnb,bool bEkinAveVel, bool bSaveOld);
161 * Compute the partial kinetic energy for home particles;
162 * will be accumulated in the calling routine.
163 * The tensor is
165 * Ekin = SUM(i) 0.5 m[i] v[i] (x) v[i]
167 * use v[i] = v[i] - u[i] when calculating temperature
169 * u must be accumulated already.
171 * Now also computes the contribution of the kinetic energy to the
172 * free energy
176 extern void
177 init_ekinstate(ekinstate_t *ekinstate,const t_inputrec *ir);
179 extern void
180 update_ekinstate(ekinstate_t *ekinstate,gmx_ekindata_t *ekind);
182 extern void
183 restore_ekinstate_from_state(t_commrec *cr,
184 gmx_ekindata_t *ekind,ekinstate_t *ekinstate);
186 extern void berendsen_tcoupl(t_inputrec *ir,gmx_ekindata_t *ekind,real dt);
188 extern void nosehoover_tcoupl(t_grpopts *opts,gmx_ekindata_t *ekind,real dt,
189 double xi[],double vxi[],t_extmass *MassQ);
191 extern t_state *init_bufstate(const t_state *template_state);
193 extern void destroy_bufstate(t_state *state);
195 extern void trotter_update(t_inputrec *ir, gmx_large_int_t step, gmx_ekindata_t *ekind,
196 gmx_enerdata_t *enerd, t_state *state, tensor vir, t_mdatoms *md,
197 t_extmass *MassQ, int **trotter_seqlist, int trotter_seqno);
199 extern int **init_npt_vars(t_inputrec *ir, t_state *state, t_extmass *Mass, bool bTrotter);
201 extern real NPT_energy(t_inputrec *ir, t_state *state, t_extmass *MassQ);
202 /* computes all the pressure/tempertature control energy terms to get a conserved energy */
204 extern void NBaroT_trotter(t_grpopts *opts, real dt,
205 double xi[],double vxi[],real *veta, t_extmass *MassQ);
207 extern void vrescale_tcoupl(t_inputrec *ir,gmx_ekindata_t *ekind,real dt,
208 double therm_integral[],
209 gmx_rng_t rng);
210 /* Compute temperature scaling. For V-rescale it is done in update. */
212 extern real vrescale_energy(t_grpopts *opts,double therm_integral[]);
213 /* Returns the V-rescale contribution to the conserved energy */
215 extern void rescale_velocities(gmx_ekindata_t *ekind,t_mdatoms *mdatoms,
216 int start,int end,rvec v[]);
217 /* Rescale the velocities with the scaling factor in ekind */
219 extern void update_annealing_target_temp(t_grpopts *opts,real t);
220 /* Set reference temp for simulated annealing at time t*/
222 extern real calc_temp(real ekin,real nrdf);
223 /* Calculate the temperature */
225 extern real calc_pres(int ePBC,int nwall,matrix box,
226 tensor ekin,tensor vir,tensor pres,real Elr);
227 /* Calculate the pressure tensor, returns the scalar pressure.
228 * The unit of pressure is bar, If Elr != 0
229 * a long range correction based on Ewald/PPPM is made (see c-code)
232 extern void parrinellorahman_pcoupl(FILE *fplog,gmx_large_int_t step,
233 t_inputrec *ir,real dt,tensor pres,
234 tensor box,tensor box_rel,tensor boxv,
235 tensor M,matrix mu,
236 bool bFirstStep);
238 extern void berendsen_pcoupl(FILE *fplog,gmx_large_int_t step,
239 t_inputrec *ir,real dt,tensor pres,matrix box,
240 matrix mu);
243 extern void berendsen_pscale(t_inputrec *ir,matrix mu,
244 matrix box,matrix box_rel,
245 int start,int nr_atoms,
246 rvec x[],unsigned short cFREEZE[],
247 t_nrnb *nrnb);
249 extern void correct_ekin(FILE *log,int start,int end,rvec v[],
250 rvec vcm,real mass[],real tmass,tensor ekin);
251 /* Correct ekin for vcm */
254 #ifdef __cplusplus
256 #endif
258 #endif /* _update_h */