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<<13)
69 #define MD_READ_RNG (1<<14)
70 #define MD_APPENDFILES (1<<15)
71 #define MD_KEEPANDNUMCPT (1<<16)
72 #define MD_READ_EKIN (1<<17)
73 #define MD_STARTFROMCPT (1<<18)
74 #define MD_RESETCOUNTERSHALFWAY (1<<19)
76 /* Define a number of flags to better control the information
77 * passed to compute_globals in md.c and global_stat.
80 /* We are rerunning the simulation */
81 #define CGLO_RERUNMD (1<<1)
82 /* we are computing the kinetic energy from average velocities */
83 #define CGLO_EKINAVEVEL (1<<2)
84 /* we are removing the center of mass momenta */
85 #define CGLO_STOPCM (1<<3)
86 /* bGStat is defined in do_md */
87 #define CGLO_GSTAT (1<<4)
88 /* Sum the energy terms in global computation */
89 #define CGLO_ENERGY (1<<6)
90 /* Sum the kinetic energy terms in global computation */
91 #define CGLO_TEMPERATURE (1<<7)
92 /* Sum the kinetic energy terms in global computation */
93 #define CGLO_PRESSURE (1<<8)
94 /* Sum the constraint term in global computation */
95 #define CGLO_CONSTRAINT (1<<9)
96 /* we are using an integrator that requires iteration over some steps - currently not used*/
97 #define CGLO_ITERATE (1<<10)
98 /* it is the first time we are iterating (or, only once through is required */
99 #define CGLO_FIRSTITERATE (1<<11)
100 /* Reading ekin from the trajectory */
101 #define CGLO_READEKIN (1<<12)
102 /* we need to reset the ekin rescaling factor here */
103 #define CGLO_SCALEEKIN (1<<13)
106 ddnoSEL
, ddnoINTERLEAVE
, ddnoPP_PME
, ddnoCARTESIAN
, ddnoNR
118 double time_per_step
;
120 gmx_large_int_t nsteps_done
;
129 gmx_bool bKeepAndNumCPT
;
136 /* Variables for temporary use with the deform option,
137 * used in runner.c and md.c.
138 * (These variables should be stored in the tpx file.)
140 extern gmx_large_int_t deform_init_init_step_tpx
;
141 extern matrix deform_init_box_tpx
;
142 #ifdef GMX_THREAD_MPI
143 extern tMPI_Thread_mutex_t deform_init_box_mutex
;
145 /* The minimum number of atoms per thread. With fewer atoms than this,
146 * the number of threads will get lowered.
148 #define MIN_ATOMS_PER_THREAD 90
152 typedef double gmx_integrator_t(FILE *log
,t_commrec
*cr
,
153 int nfile
,const t_filenm fnm
[],
154 const output_env_t oenv
, gmx_bool bVerbose
,
155 gmx_bool bCompact
, int nstglobalcomm
,
156 gmx_vsite_t
*vsite
,gmx_constr_t constr
,
158 t_inputrec
*inputrec
,
159 gmx_mtop_t
*mtop
,t_fcdata
*fcd
,
162 t_nrnb
*nrnb
,gmx_wallcycle_t wcycle
,
165 int repl_ex_nst
,int repl_ex_seed
,
166 real cpt_period
,real max_hours
,
167 const char *deviceOptions
,
169 gmx_runtime_t
*runtime
);
171 typedef struct gmx_global_stat
*gmx_global_stat_t
;
173 /* ROUTINES from md.c */
175 gmx_integrator_t do_md
;
177 gmx_integrator_t do_md_openmm
;
181 /* ROUTINES from minimize.c */
183 gmx_integrator_t do_steep
;
184 /* Do steepest descents EM */
186 gmx_integrator_t do_cg
;
187 /* Do conjugate gradient EM */
189 gmx_integrator_t do_lbfgs
;
190 /* Do conjugate gradient L-BFGS */
192 gmx_integrator_t do_nm
;
193 /* Do normal mode analysis */
195 /* ROUTINES from tpi.c */
197 gmx_integrator_t do_tpi
;
198 /* Do test particle insertion */
201 /* ROUTINES from md_support.c */
203 /* return the number of steps between global communcations */
204 int check_nstglobalcomm(FILE *fplog
,t_commrec
*cr
,
205 int nstglobalcomm
,t_inputrec
*ir
);
207 /* check whether an 'nst'-style parameter p is a multiple of nst, and
208 set it to be one if not, with a warning. */
209 void check_nst_param(FILE *fplog
,t_commrec
*cr
,
210 const char *desc_nst
,int nst
,
211 const char *desc_p
,int *p
);
213 /* check which of the multisim simulations has the shortest number of
214 steps and return that number of nsteps */
215 gmx_large_int_t
get_multisim_nsteps(const t_commrec
*cr
,
216 gmx_large_int_t nsteps
);
218 void rerun_parallel_comm(t_commrec
*cr
,t_trxframe
*fr
,
219 gmx_bool
*bNotLastFrame
);
221 /* get the conserved energy associated with the ensemble type*/
222 real
compute_conserved_from_auxiliary(t_inputrec
*ir
, t_state
*state
,
225 /* reset all cycle and time counters. */
226 void reset_all_counters(FILE *fplog
,t_commrec
*cr
,
227 gmx_large_int_t step
,
228 gmx_large_int_t
*step_rel
,t_inputrec
*ir
,
229 gmx_wallcycle_t wcycle
,t_nrnb
*nrnb
,
230 gmx_runtime_t
*runtime
);
234 /* ROUTINES from sim_util.c */
235 void do_pbc_first(FILE *log
,matrix box
,t_forcerec
*fr
,
236 t_graph
*graph
,rvec x
[]);
238 void do_pbc_first_mtop(FILE *fplog
,int ePBC
,matrix box
,
239 gmx_mtop_t
*mtop
,rvec x
[]);
241 void do_pbc_mtop(FILE *fplog
,int ePBC
,matrix box
,
242 gmx_mtop_t
*mtop
,rvec x
[]);
246 /* ROUTINES from stat.c */
247 gmx_global_stat_t
global_stat_init(t_inputrec
*ir
);
249 void global_stat_destroy(gmx_global_stat_t gs
);
251 void global_stat(FILE *log
,gmx_global_stat_t gs
,
252 t_commrec
*cr
,gmx_enerdata_t
*enerd
,
253 tensor fvir
,tensor svir
,rvec mu_tot
,
254 t_inputrec
*inputrec
,
255 gmx_ekindata_t
*ekind
,
256 gmx_constr_t constr
,t_vcm
*vcm
,
258 gmx_mtop_t
*top_global
, t_state
*state_local
,
259 gmx_bool bSumEkinhOld
, int flags
);
260 /* Communicate statistics over cr->mpi_comm_mysim */
262 gmx_mdoutf_t
*init_mdoutf(int nfile
,const t_filenm fnm
[],
264 const t_commrec
*cr
,const t_inputrec
*ir
,
265 const output_env_t oenv
);
266 /* Returns a pointer to a data structure with all output file pointers
267 * and names required by mdrun.
270 void done_mdoutf(gmx_mdoutf_t
*of
);
271 /* Close all open output files and free the of pointer */
273 #define MDOF_X (1<<0)
274 #define MDOF_V (1<<1)
275 #define MDOF_F (1<<2)
276 #define MDOF_XTC (1<<3)
277 #define MDOF_CPT (1<<4)
279 void write_traj(FILE *fplog
,t_commrec
*cr
,
282 gmx_mtop_t
*top_global
,
283 gmx_large_int_t step
,double t
,
284 t_state
*state_local
,t_state
*state_global
,
285 rvec
*f_local
,rvec
*f_global
,
286 int *n_xtc
,rvec
**x_xtc
);
287 /* Routine that writes frames to trn, xtc and/or checkpoint.
288 * What is written is determined by the mdof_flags defined above.
289 * Data is collected to the master node only when necessary.
292 int do_per_step(gmx_large_int_t step
,gmx_large_int_t nstep
);
293 /* Return TRUE if io should be done */
295 int do_any_io(int step
, t_inputrec
*ir
);
297 /* ROUTINES from sim_util.c */
299 double gmx_gettime();
301 void print_time(FILE *out
, gmx_runtime_t
*runtime
,
302 gmx_large_int_t step
,t_inputrec
*ir
, t_commrec
*cr
);
304 void runtime_start(gmx_runtime_t
*runtime
);
306 void runtime_end(gmx_runtime_t
*runtime
);
308 void runtime_upd_proc(gmx_runtime_t
*runtime
);
309 /* The processor time should be updated every once in a while,
310 * since on 32-bit manchines it loops after 72 minutes.
313 void print_date_and_time(FILE *log
,int pid
,const char *title
,
314 const gmx_runtime_t
*runtime
);
316 void nstop_cm(FILE *log
,t_commrec
*cr
,
317 int start
,int nr_atoms
,real mass
[],rvec x
[],rvec v
[]);
319 void finish_run(FILE *log
,t_commrec
*cr
,const char *confout
,
320 t_inputrec
*inputrec
,
321 t_nrnb nrnb
[],gmx_wallcycle_t wcycle
,
322 gmx_runtime_t
*runtime
,
323 gmx_bool bWriteStat
);
325 void calc_enervirdiff(FILE *fplog
,int eDispCorr
,t_forcerec
*fr
);
327 void calc_dispcorr(FILE *fplog
,t_inputrec
*ir
,t_forcerec
*fr
,
328 gmx_large_int_t step
, int natoms
,
329 matrix box
,real lambda
,tensor pres
,tensor virial
,
330 real
*prescorr
, real
*enercorr
, real
*dvdlcorr
);
343 void check_nnodes_top(char *fn
,t_topology
*top
);
344 /* Reset the tpr file to work with one node if necessary */
347 /* check the version */
348 void check_ir_old_tpx_versions(t_commrec
*cr
,FILE *fplog
,
349 t_inputrec
*ir
,gmx_mtop_t
*mtop
);
351 /* Allocate and initialize node-local state entries. */
352 void set_state_entries(t_state
*state
,const t_inputrec
*ir
,int nnodes
);
354 /* Broadcast the data for a simulation, and allocate node-specific settings
355 such as rng generators. */
356 void init_parallel(FILE *log
, t_commrec
*cr
, t_inputrec
*inputrec
,
360 void do_constrain_first(FILE *log
,gmx_constr_t constr
,
361 t_inputrec
*inputrec
,t_mdatoms
*md
,
362 t_state
*state
,rvec
*f
,
363 t_graph
*graph
,t_commrec
*cr
,t_nrnb
*nrnb
,
364 t_forcerec
*fr
, gmx_localtop_t
*top
, tensor shake_vir
);
366 void dynamic_load_balancing(gmx_bool bVerbose
,t_commrec
*cr
,real capacity
[],
367 int dimension
,t_mdatoms
*md
,t_topology
*top
,
368 rvec x
[],rvec v
[],matrix box
);
369 /* Perform load balancing, i.e. split the particles over processors
370 * based on their coordinates in the "dimension" direction.
373 int multisim_min(const gmx_multisim_t
*ms
,int nmin
,int n
);
374 /* Set an appropriate value for n across the whole multi-simulation */
376 int multisim_nstsimsync(const t_commrec
*cr
,
377 const t_inputrec
*ir
,int repl_ex_nst
);
378 /* Determine the interval for inter-simulation communication */
380 void init_global_signals(globsig_t
*gs
,const t_commrec
*cr
,
381 const t_inputrec
*ir
,int repl_ex_nst
);
382 /* Constructor for globsig_t */
384 void copy_coupling_state(t_state
*statea
,t_state
*stateb
,
385 gmx_ekindata_t
*ekinda
,gmx_ekindata_t
*ekindb
, t_grpopts
* opts
);
386 /* Copy stuff from state A to state B */
388 void compute_globals(FILE *fplog
, gmx_global_stat_t gstat
, t_commrec
*cr
, t_inputrec
*ir
,
389 t_forcerec
*fr
, gmx_ekindata_t
*ekind
,
390 t_state
*state
, t_state
*state_global
, t_mdatoms
*mdatoms
,
391 t_nrnb
*nrnb
, t_vcm
*vcm
, gmx_wallcycle_t wcycle
,
392 gmx_enerdata_t
*enerd
,tensor force_vir
, tensor shake_vir
, tensor total_vir
,
393 tensor pres
, rvec mu_tot
, gmx_constr_t constr
,
394 globsig_t
*gs
,gmx_bool bInterSimGS
,
395 matrix box
, gmx_mtop_t
*top_global
, real
*pcurr
,
396 int natoms
, gmx_bool
*bSumEkinhOld
, int flags
);
397 /* Compute global variables during integration */
399 int mdrunner(int nthreads_requested
, FILE *fplog
,t_commrec
*cr
,int nfile
,
400 const t_filenm fnm
[], const output_env_t oenv
, gmx_bool bVerbose
,
401 gmx_bool bCompact
, int nstglobalcomm
, ivec ddxyz
,int dd_node_order
,
402 real rdd
, real rconstr
, const char *dddlb_opt
,real dlb_scale
,
403 const char *ddcsx
,const char *ddcsy
,const char *ddcsz
,
404 int nstepout
, int resetstep
, int nmultisim
, int repl_ex_nst
,
405 int repl_ex_seed
, real pforce
,real cpt_period
,real max_hours
,
406 const char *deviceOptions
, unsigned long Flags
);
407 /* Driver routine, that calls the different methods */
409 void md_print_warning(const t_commrec
*cr
,FILE *fplog
,const char *buf
);
410 /* Print a warning message to stderr on the master node
411 * and to fplog if fplog!=NULL.
414 void init_md(FILE *fplog
,
415 t_commrec
*cr
,t_inputrec
*ir
, const output_env_t oenv
,
416 double *t
,double *t0
,
417 real
*lambda
,double *lam0
,
418 t_nrnb
*nrnb
,gmx_mtop_t
*mtop
,
420 int nfile
,const t_filenm fnm
[],
421 gmx_mdoutf_t
**outf
,t_mdebin
**mdebin
,
422 tensor force_vir
,tensor shake_vir
,
424 gmx_bool
*bSimAnn
,t_vcm
**vcm
,
425 t_state
*state
, unsigned long Flags
);
426 /* Routine in sim_util.c */
432 #endif /* _mdrun_h */