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
58 #define MD_POLARISE (1<<2)
59 #define MD_IONIZE (1<<3)
60 #define MD_RERUN (1<<4)
61 #define MD_RERUN_VSITE (1<<5)
62 #define MD_FFSCAN (1<<6)
63 #define MD_SEPPOT (1<<7)
64 #define MD_PARTDEC (1<<9)
65 #define MD_DDBONDCHECK (1<<10)
66 #define MD_DDBONDCOMM (1<<11)
67 #define MD_CONFOUT (1<<12)
68 #define MD_REPRODUCIBLE (1<<14)
69 #define MD_READ_RNG (1<<15)
70 #define MD_APPENDFILES (1<<16)
71 #define MD_READ_EKIN (1<<17)
72 #define MD_STARTFROMCPT (1<<18)
76 ddnoSEL
, ddnoINTERLEAVE
, ddnoPP_PME
, ddnoCARTESIAN
, ddnoNR
90 gmx_large_int_t nsteps_done
;
93 typedef double gmx_integrator_t(FILE *log
,t_commrec
*cr
,
94 int nfile
,const t_filenm fnm
[],
95 const output_env_t oenv
, bool bVerbose
,
96 bool bCompact
, int nstglobalcomm
,
97 gmx_vsite_t
*vsite
,gmx_constr_t constr
,
100 gmx_mtop_t
*mtop
,t_fcdata
*fcd
,
103 t_nrnb
*nrnb
,gmx_wallcycle_t wcycle
,
106 int repl_ex_nst
,int repl_ex_seed
,
107 real cpt_period
,real max_hours
,
109 gmx_runtime_t
*runtime
);
111 typedef struct gmx_global_stat
*gmx_global_stat_t
;
113 /* ROUTINES from md.c */
115 extern gmx_integrator_t do_md
;
117 /* ROUTINES from minimize.c */
119 extern gmx_integrator_t do_steep
;
120 /* Do steepest descents EM */
122 extern gmx_integrator_t do_cg
;
123 /* Do conjugate gradient EM */
125 extern gmx_integrator_t do_lbfgs
;
126 /* Do conjugate gradient L-BFGS */
128 extern gmx_integrator_t do_nm
;
129 /* Do normal mode analysis */
131 extern gmx_integrator_t do_tpi
;
132 /* Do test particle insertion */
135 /* ROUTINES from sim_util.c */
136 extern void do_pbc_first(FILE *log
,matrix box
,t_forcerec
*fr
,
137 t_graph
*graph
,rvec x
[]);
139 extern void do_pbc_first_mtop(FILE *fplog
,int ePBC
,matrix box
,
140 gmx_mtop_t
*mtop
,rvec x
[]);
142 extern void do_pbc_mtop(FILE *fplog
,int ePBC
,matrix box
,
143 gmx_mtop_t
*mtop
,rvec x
[]);
146 /* ROUTINES from stat.c */
147 extern gmx_global_stat_t
global_stat_init(t_inputrec
*ir
);
149 extern void global_stat_destroy(gmx_global_stat_t gs
);
151 extern void global_stat(FILE *log
,gmx_global_stat_t gs
,
152 t_commrec
*cr
,gmx_enerdata_t
*enerd
,
153 tensor fvir
,tensor svir
,rvec mu_tot
,
154 t_inputrec
*inputrec
,
155 gmx_ekindata_t
*ekind
,bool bSumEkinhOld
,
156 gmx_constr_t constr
,t_vcm
*vcm
,
157 int *nabnsb
,real
*chkpt
,real
*terminate
,
158 gmx_mtop_t
*top_global
, t_state
*state_local
);
159 /* Communicate statistics over cr->mpi_comm_mysim */
161 void write_traj(FILE *fplog
,t_commrec
*cr
,
162 int fp_trn
,bool bX
,bool bV
,bool bF
,
163 int fp_xtc
,bool bXTC
,int xtc_prec
,
164 const char *fn_cpt
,bool bCPT
,
165 gmx_mtop_t
*top_global
,
166 int eIntegrator
,int simulation_part
,gmx_large_int_t step
,double t
,
167 t_state
*state_local
,t_state
*state_global
,
168 rvec
*f_local
,rvec
*f_global
,
169 int *n_xtc
,rvec
**x_xtc
);
170 /* Routine that writes frames to trn, xtc and/or checkpoint.
171 * Data is collected to the master node only when necessary.
174 extern int do_per_step(gmx_large_int_t step
,gmx_large_int_t nstep
);
175 /* Return TRUE if io should be done */
177 extern int do_any_io(int step
, t_inputrec
*ir
);
179 /* ROUTINES from sim_util.c */
181 extern void print_time(FILE *out
, gmx_runtime_t
*runtime
,
182 gmx_large_int_t step
,t_inputrec
*ir
, t_commrec
*cr
);
184 extern void runtime_start(gmx_runtime_t
*runtime
);
186 extern void runtime_end(gmx_runtime_t
*runtime
);
188 extern void runtime_upd_proc(gmx_runtime_t
*runtime
);
189 /* The processor time should be updated every once in a while,
190 * since on 32-bit manchines it loops after 72 minutes.
193 extern void print_date_and_time(FILE *log
,int pid
,const char *title
,
194 const gmx_runtime_t
*runtime
);
196 extern void nstop_cm(FILE *log
,t_commrec
*cr
,
197 int start
,int nr_atoms
,real mass
[],rvec x
[],rvec v
[]);
199 extern void finish_run(FILE *log
,t_commrec
*cr
,const char *confout
,
200 t_inputrec
*inputrec
,
201 t_nrnb nrnb
[],gmx_wallcycle_t wcycle
,
202 gmx_runtime_t
*runtime
,
205 extern void calc_enervirdiff(FILE *fplog
,int eDispCorr
,t_forcerec
*fr
);
207 extern void calc_dispcorr(FILE *fplog
,t_inputrec
*ir
,t_forcerec
*fr
,
208 gmx_large_int_t step
,int natoms
,matrix box
,real lambda
,
209 tensor pres
,tensor virial
,gmx_enerdata_t
*enerd
);
223 extern void check_nnodes_top(char *fn
,t_topology
*top
);
224 /* Reset the tpr file to work with one node if necessary */
226 extern void init_single(FILE *log
, t_inputrec
*inputrec
, const char *tpbfile
,
227 gmx_mtop_t
*mtop
, t_state
*state
);
229 * Allocates space for the topology (top), the coordinates x, the
230 * velocities v, masses mass. Reads the parameters, topology,
231 * coordinates and velocities from the file specified in tpbfile
234 extern void init_parallel(FILE *log
,const char *tpxfile
, t_commrec
*cr
,
235 t_inputrec
*inputrec
,gmx_mtop_t
*mtop
,
236 t_state
*state
, int list
);
238 * Loads the data for a simulation from the ring. Parameters, topology
239 * coordinates, velocities, and masses are initialised equal to using
240 * init_single() in the single processor version. The extra argument
241 * f_add is allocated to use for the update of the forces, the load
242 * array specifies in which part of the x and f array the subsystems
243 * of the other processors are located. Homenr0, homenr1, nparts0 and
244 * nparts1 are necessary to calculate the non bonded interaction using
245 * the symmetry and thus calculating every force only once. List is a
246 * facility for logging (and debugging). One can decide to print none or a
247 * set of * selected parameters to the file specified by log. Parameters are
248 * printed by or-ing the corresponding items from t_listitem. A 0 (zero)
249 * specifies that nothing is to be printed on the file. The function
250 * returns the number of shifts over the ring to perform to calculate
254 extern void do_constrain_first(FILE *log
,gmx_constr_t constr
,
255 t_inputrec
*inputrec
,t_mdatoms
*md
,
257 t_graph
*graph
,t_commrec
*cr
,t_nrnb
*nrnb
,
258 t_forcerec
*fr
,t_idef
*idef
);
260 extern void dynamic_load_balancing(bool bVerbose
,t_commrec
*cr
,real capacity
[],
261 int dimension
,t_mdatoms
*md
,t_topology
*top
,
262 rvec x
[],rvec v
[],matrix box
);
263 /* Perform load balancing, i.e. split the particles over processors
264 * based on their coordinates in the "dimension" direction.
267 int mdrunner(FILE *fplog
,t_commrec
*cr
,int nfile
,const t_filenm fnm
[],
268 const output_env_t oenv
, bool bVerbose
,bool bCompact
,
269 int nstglobalcomm
, ivec ddxyz
,int dd_node_order
,real rdd
,
270 real rconstr
, const char *dddlb_opt
,real dlb_scale
,
271 const char *ddcsx
,const char *ddcsy
,const char *ddcsz
,
272 int nstepout
, int resetstep
, int nmultisim
, int repl_ex_nst
,int repl_ex_seed
,
273 real pforce
,real cpt_period
,real max_hours
,
274 unsigned long Flags
);
275 /* Driver routine, that calls the different methods */
277 int mdrunner_threads(int nthreads
,
278 FILE *fplog
,t_commrec
*cr
,int nfile
,const t_filenm fnm
[],
279 const output_env_t oenv
, bool bVerbose
,bool bCompact
,
281 ivec ddxyz
,int dd_node_order
,real rdd
,real rconstr
,
282 const char *dddlb_opt
,real dlb_scale
,
283 const char *ddcsx
,const char *ddcsy
,const char *ddcsz
,
284 int nstepout
,int resetstep
,int nmultisim
, int repl_ex_nst
,
285 int repl_ex_seed
, real pforce
,real cpt_period
,
286 real max_hours
, unsigned long Flags
);
287 /* initializes nthread threads before running mdrunner: is the preferred
288 way to start a simulation (even if nthreads=1 and no threads are started) */
291 extern void init_md(FILE *fplog
,
292 t_commrec
*cr
,t_inputrec
*ir
, const output_env_t oenv
,
293 double *t
,double *t0
,
294 real
*lambda
,double *lam0
,
295 t_nrnb
*nrnb
,gmx_mtop_t
*mtop
,
297 int nfile
,const t_filenm fnm
[],
298 int *fp_trn
,int *fp_xtc
,ener_file_t
*fp_ene
,
300 FILE **fp_dhdl
,FILE **fp_field
,
302 tensor force_vir
,tensor shake_vir
,
304 bool *bNEMD
,bool *bSimAnn
,t_vcm
**vcm
, unsigned long Flags
);
305 /* Routine in sim_util.c */
307 #endif /* _mdrun_h */