Use relative rpath to support relocatable binary packages
[gromacs.git] / include / mdrun.h
blobdb7923a571f9b12b5e3d588d02c2afc933c889ed
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 _mdrun_h
37 #define _mdrun_h
39 #include <stdio.h>
40 #include <time.h>
41 #include "typedefs.h"
42 #include "network.h"
43 #include "tgroup.h"
44 #include "filenm.h"
45 #include "mshift.h"
46 #include "force.h"
47 #include "edsam.h"
48 #include "mdebin.h"
49 #include "vcm.h"
50 #include "vsite.h"
51 #include "pull.h"
52 #include "update.h"
55 #ifdef GMX_THREAD_MPI
56 #include "thread_mpi/threads.h"
57 #endif
59 #ifdef __cplusplus
60 extern "C" {
61 #endif
63 #define MD_POLARISE (1<<2)
64 #define MD_IONIZE (1<<3)
65 #define MD_RERUN (1<<4)
66 #define MD_RERUN_VSITE (1<<5)
67 #define MD_FFSCAN (1<<6)
68 #define MD_SEPPOT (1<<7)
69 #define MD_PARTDEC (1<<9)
70 #define MD_DDBONDCHECK (1<<10)
71 #define MD_DDBONDCOMM (1<<11)
72 #define MD_CONFOUT (1<<12)
73 #define MD_REPRODUCIBLE (1<<13)
74 #define MD_READ_RNG (1<<14)
75 #define MD_APPENDFILES (1<<15)
76 #define MD_APPENDFILESSET (1<<21)
77 #define MD_KEEPANDNUMCPT (1<<16)
78 #define MD_READ_EKIN (1<<17)
79 #define MD_STARTFROMCPT (1<<18)
80 #define MD_RESETCOUNTERSHALFWAY (1<<19)
82 /* Define a number of flags to better control the information
83 * passed to compute_globals in md.c and global_stat.
86 /* We are rerunning the simulation */
87 #define CGLO_RERUNMD (1<<1)
88 /* we are computing the kinetic energy from average velocities */
89 #define CGLO_EKINAVEVEL (1<<2)
90 /* we are removing the center of mass momenta */
91 #define CGLO_STOPCM (1<<3)
92 /* bGStat is defined in do_md */
93 #define CGLO_GSTAT (1<<4)
94 /* Sum the energy terms in global computation */
95 #define CGLO_ENERGY (1<<6)
96 /* Sum the kinetic energy terms in global computation */
97 #define CGLO_TEMPERATURE (1<<7)
98 /* Sum the kinetic energy terms in global computation */
99 #define CGLO_PRESSURE (1<<8)
100 /* Sum the constraint term in global computation */
101 #define CGLO_CONSTRAINT (1<<9)
102 /* we are using an integrator that requires iteration over some steps - currently not used*/
103 #define CGLO_ITERATE (1<<10)
104 /* it is the first time we are iterating (or, only once through is required */
105 #define CGLO_FIRSTITERATE (1<<11)
106 /* Reading ekin from the trajectory */
107 #define CGLO_READEKIN (1<<12)
108 /* we need to reset the ekin rescaling factor here */
109 #define CGLO_SCALEEKIN (1<<13)
111 enum {
112 ddnoSEL, ddnoINTERLEAVE, ddnoPP_PME, ddnoCARTESIAN, ddnoNR
115 typedef struct {
116 double real;
117 #ifdef GMX_CRAY_XT3
118 double proc;
119 #else
120 clock_t proc;
121 #endif
122 double realtime;
123 double proctime;
124 double time_per_step;
125 double last;
126 gmx_large_int_t nsteps_done;
127 } gmx_runtime_t;
129 typedef struct {
130 t_fileio *fp_trn;
131 t_fileio *fp_xtc;
132 int xtc_prec;
133 ener_file_t fp_ene;
134 const char *fn_cpt;
135 gmx_bool bKeepAndNumCPT;
136 int eIntegrator;
137 gmx_bool bExpanded;
138 int elamstats;
139 int simulation_part;
140 FILE *fp_dhdl;
141 FILE *fp_field;
142 } gmx_mdoutf_t;
144 /* Variables for temporary use with the deform option,
145 * used in runner.c and md.c.
146 * (These variables should be stored in the tpx file.)
148 extern gmx_large_int_t deform_init_init_step_tpx;
149 extern matrix deform_init_box_tpx;
150 #ifdef GMX_THREAD_MPI
151 extern tMPI_Thread_mutex_t deform_init_box_mutex;
153 /* The minimum number of atoms per thread. With fewer atoms than this,
154 * the number of threads will get lowered.
156 #define MIN_ATOMS_PER_THREAD 90
157 #endif
160 typedef double gmx_integrator_t(FILE *log,t_commrec *cr,
161 int nfile,const t_filenm fnm[],
162 const output_env_t oenv, gmx_bool bVerbose,
163 gmx_bool bCompact, int nstglobalcomm,
164 gmx_vsite_t *vsite,gmx_constr_t constr,
165 int stepout,
166 t_inputrec *inputrec,
167 gmx_mtop_t *mtop,t_fcdata *fcd,
168 t_state *state,
169 t_mdatoms *mdatoms,
170 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
171 gmx_edsam_t ed,
172 t_forcerec *fr,
173 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
174 gmx_membed_t membed,
175 real cpt_period,real max_hours,
176 const char *deviceOptions,
177 unsigned long Flags,
178 gmx_runtime_t *runtime);
180 typedef struct gmx_global_stat *gmx_global_stat_t;
182 /* ROUTINES from md.c */
184 gmx_integrator_t do_md;
186 gmx_integrator_t do_md_openmm;
190 /* ROUTINES from minimize.c */
192 gmx_integrator_t do_steep;
193 /* Do steepest descents EM */
195 gmx_integrator_t do_cg;
196 /* Do conjugate gradient EM */
198 gmx_integrator_t do_lbfgs;
199 /* Do conjugate gradient L-BFGS */
201 gmx_integrator_t do_nm;
202 /* Do normal mode analysis */
204 /* ROUTINES from tpi.c */
206 gmx_integrator_t do_tpi;
207 /* Do test particle insertion */
210 /* ROUTINES from md_support.c */
212 /* return the number of steps between global communcations */
213 int check_nstglobalcomm(FILE *fplog,t_commrec *cr,
214 int nstglobalcomm,t_inputrec *ir);
216 /* check whether an 'nst'-style parameter p is a multiple of nst, and
217 set it to be one if not, with a warning. */
218 void check_nst_param(FILE *fplog,t_commrec *cr,
219 const char *desc_nst,int nst,
220 const char *desc_p,int *p);
222 /* check which of the multisim simulations has the shortest number of
223 steps and return that number of nsteps */
224 gmx_large_int_t get_multisim_nsteps(const t_commrec *cr,
225 gmx_large_int_t nsteps);
227 void rerun_parallel_comm(t_commrec *cr,t_trxframe *fr,
228 gmx_bool *bNotLastFrame);
230 /* get the conserved energy associated with the ensemble type*/
231 real compute_conserved_from_auxiliary(t_inputrec *ir, t_state *state,
232 t_extmass *MassQ);
234 /* set the lambda values at each step of mdrun when they change */
235 void set_current_lambdas(gmx_large_int_t step, t_lambda *fepvals, gmx_bool bRerunMD,
236 t_trxframe *rerun_fr, t_state *state_global, t_state *state, double lam0[]);
238 /* reset all cycle and time counters. */
239 void reset_all_counters(FILE *fplog,t_commrec *cr,
240 gmx_large_int_t step,
241 gmx_large_int_t *step_rel,t_inputrec *ir,
242 gmx_wallcycle_t wcycle,t_nrnb *nrnb,
243 gmx_runtime_t *runtime);
247 /* ROUTINES from sim_util.c */
248 void do_pbc_first(FILE *log,matrix box,t_forcerec *fr,
249 t_graph *graph,rvec x[]);
251 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
252 gmx_mtop_t *mtop,rvec x[]);
254 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
255 gmx_mtop_t *mtop,rvec x[]);
259 /* ROUTINES from stat.c */
260 gmx_global_stat_t global_stat_init(t_inputrec *ir);
262 void global_stat_destroy(gmx_global_stat_t gs);
264 void global_stat(FILE *log,gmx_global_stat_t gs,
265 t_commrec *cr,gmx_enerdata_t *enerd,
266 tensor fvir,tensor svir,rvec mu_tot,
267 t_inputrec *inputrec,
268 gmx_ekindata_t *ekind,
269 gmx_constr_t constr,t_vcm *vcm,
270 int nsig,real *sig,
271 gmx_mtop_t *top_global, t_state *state_local,
272 gmx_bool bSumEkinhOld, int flags);
273 /* Communicate statistics over cr->mpi_comm_mysim */
275 gmx_mdoutf_t *init_mdoutf(int nfile,const t_filenm fnm[],
276 int mdrun_flags,
277 const t_commrec *cr,const t_inputrec *ir,
278 const output_env_t oenv);
279 /* Returns a pointer to a data structure with all output file pointers
280 * and names required by mdrun.
283 void done_mdoutf(gmx_mdoutf_t *of);
284 /* Close all open output files and free the of pointer */
286 #define MDOF_X (1<<0)
287 #define MDOF_V (1<<1)
288 #define MDOF_F (1<<2)
289 #define MDOF_XTC (1<<3)
290 #define MDOF_CPT (1<<4)
292 void write_traj(FILE *fplog,t_commrec *cr,
293 gmx_mdoutf_t *of,
294 int mdof_flags,
295 gmx_mtop_t *top_global,
296 gmx_large_int_t step,double t,
297 t_state *state_local,t_state *state_global,
298 rvec *f_local,rvec *f_global,
299 int *n_xtc,rvec **x_xtc);
300 /* Routine that writes frames to trn, xtc and/or checkpoint.
301 * What is written is determined by the mdof_flags defined above.
302 * Data is collected to the master node only when necessary.
305 int do_per_step(gmx_large_int_t step,gmx_large_int_t nstep);
306 /* Return TRUE if io should be done */
308 int do_any_io(int step, t_inputrec *ir);
310 /* ROUTINES from sim_util.c */
312 double gmx_gettime();
314 void print_time(FILE *out, gmx_runtime_t *runtime,
315 gmx_large_int_t step,t_inputrec *ir, t_commrec *cr);
317 void runtime_start(gmx_runtime_t *runtime);
319 void runtime_end(gmx_runtime_t *runtime);
321 void runtime_upd_proc(gmx_runtime_t *runtime);
322 /* The processor time should be updated every once in a while,
323 * since on 32-bit manchines it loops after 72 minutes.
326 void print_date_and_time(FILE *log,int pid,const char *title,
327 const gmx_runtime_t *runtime);
329 void nstop_cm(FILE *log,t_commrec *cr,
330 int start,int nr_atoms,real mass[],rvec x[],rvec v[]);
332 void finish_run(FILE *log,t_commrec *cr,const char *confout,
333 t_inputrec *inputrec,
334 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
335 gmx_runtime_t *runtime,
336 gmx_bool bWriteStat);
338 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr);
340 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
341 gmx_large_int_t step, int natoms,
342 matrix box,real lambda,tensor pres,tensor virial,
343 real *prescorr, real *enercorr, real *dvdlcorr);
345 void initialize_lambdas(FILE *fplog,t_inputrec *ir,int *fep_state,real *lambda,double *lam0);
347 void init_npt_masses(t_inputrec *ir, t_state *state, t_extmass *MassQ, gmx_bool bInit);
349 int ExpandedEnsembleDynamics(FILE *log,t_inputrec *ir, gmx_enerdata_t *enerd,
350 t_state *state, t_extmass *MassQ, df_history_t *dfhist,
351 gmx_large_int_t step, gmx_rng_t mcrng,
352 rvec *v, t_mdatoms *mdatoms);
354 void PrintFreeEnergyInfoToFile(FILE *outfile, t_lambda *fep, t_expanded *expand, t_simtemp *simtemp, df_history_t *dfhist,
355 int nlam, int frequency, gmx_large_int_t step);
357 void get_mc_state(gmx_rng_t rng,t_state *state);
359 void set_mc_state(gmx_rng_t rng,t_state *state);
362 typedef enum
364 LIST_SCALARS =0001,
365 LIST_INPUTREC =0002,
366 LIST_TOP =0004,
367 LIST_X =0010,
368 LIST_V =0020,
369 LIST_F =0040,
370 LIST_LOAD =0100
371 } t_listitem;
373 void check_nnodes_top(char *fn,t_topology *top);
374 /* Reset the tpr file to work with one node if necessary */
377 /* check the version */
378 void check_ir_old_tpx_versions(t_commrec *cr,FILE *fplog,
379 t_inputrec *ir,gmx_mtop_t *mtop);
381 /* Allocate and initialize node-local state entries. */
382 void set_state_entries(t_state *state,const t_inputrec *ir,int nnodes);
384 /* Broadcast the data for a simulation, and allocate node-specific settings
385 such as rng generators. */
386 void init_parallel(FILE *log, t_commrec *cr, t_inputrec *inputrec,
387 gmx_mtop_t *mtop);
390 void do_constrain_first(FILE *log,gmx_constr_t constr,
391 t_inputrec *inputrec,t_mdatoms *md,
392 t_state *state,rvec *f,
393 t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
394 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir);
396 void dynamic_load_balancing(gmx_bool bVerbose,t_commrec *cr,real capacity[],
397 int dimension,t_mdatoms *md,t_topology *top,
398 rvec x[],rvec v[],matrix box);
399 /* Perform load balancing, i.e. split the particles over processors
400 * based on their coordinates in the "dimension" direction.
403 int multisim_min(const gmx_multisim_t *ms,int nmin,int n);
404 /* Set an appropriate value for n across the whole multi-simulation */
406 int multisim_nstsimsync(const t_commrec *cr,
407 const t_inputrec *ir,int repl_ex_nst);
408 /* Determine the interval for inter-simulation communication */
410 void init_global_signals(globsig_t *gs,const t_commrec *cr,
411 const t_inputrec *ir,int repl_ex_nst);
412 /* Constructor for globsig_t */
414 void copy_coupling_state(t_state *statea,t_state *stateb,
415 gmx_ekindata_t *ekinda,gmx_ekindata_t *ekindb, t_grpopts* opts);
416 /* Copy stuff from state A to state B */
418 void compute_globals(FILE *fplog, gmx_global_stat_t gstat, t_commrec *cr, t_inputrec *ir,
419 t_forcerec *fr, gmx_ekindata_t *ekind,
420 t_state *state, t_state *state_global, t_mdatoms *mdatoms,
421 t_nrnb *nrnb, t_vcm *vcm, gmx_wallcycle_t wcycle,
422 gmx_enerdata_t *enerd,tensor force_vir, tensor shake_vir, tensor total_vir,
423 tensor pres, rvec mu_tot, gmx_constr_t constr,
424 globsig_t *gs,gmx_bool bInterSimGS,
425 matrix box, gmx_mtop_t *top_global, real *pcurr,
426 int natoms, gmx_bool *bSumEkinhOld, int flags);
427 /* Compute global variables during integration */
429 int mdrunner(int nthreads_requested, FILE *fplog,t_commrec *cr,int nfile,
430 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
431 gmx_bool bCompact, int nstglobalcomm, ivec ddxyz,int dd_node_order,
432 real rdd, real rconstr, const char *dddlb_opt,real dlb_scale,
433 const char *ddcsx,const char *ddcsy,const char *ddcsz,
434 int nstepout, int resetstep, int nmultisim, int repl_ex_nst, int repl_ex_nex,
435 int repl_ex_seed, real pforce,real cpt_period,real max_hours,
436 const char *deviceOptions, unsigned long Flags);
437 /* Driver routine, that calls the different methods */
439 void md_print_warning(const t_commrec *cr,FILE *fplog,const char *buf);
440 /* Print a warning message to stderr on the master node
441 * and to fplog if fplog!=NULL.
444 void init_md(FILE *fplog,
445 t_commrec *cr,t_inputrec *ir, const output_env_t oenv,
446 double *t,double *t0,
447 real *lambda,int *fep_state, double *lam0,
448 t_nrnb *nrnb,gmx_mtop_t *mtop,
449 gmx_update_t *upd,
450 int nfile,const t_filenm fnm[],
451 gmx_mdoutf_t **outf,t_mdebin **mdebin,
452 tensor force_vir,tensor shake_vir,
453 rvec mu_tot,
454 gmx_bool *bSimAnn,t_vcm **vcm,
455 t_state *state, unsigned long Flags);
456 /* Routine in sim_util.c */
458 #ifdef __cplusplus
460 #endif
462 #endif /* _mdrun_h */