Merge branch 'master' of git://git.gromacs.org/gromacs
[gromacs/adressmacs.git] / include / mdrun.h
blobdb0f7d7c05922c8a8d25b0e88896230bd3fae2a4
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"
53 #include "membed.h"
55 #ifdef __cplusplus
56 extern "C" {
57 #endif
59 #define MD_POLARISE (1<<2)
60 #define MD_IONIZE (1<<3)
61 #define MD_RERUN (1<<4)
62 #define MD_RERUN_VSITE (1<<5)
63 #define MD_FFSCAN (1<<6)
64 #define MD_SEPPOT (1<<7)
65 #define MD_PARTDEC (1<<9)
66 #define MD_DDBONDCHECK (1<<10)
67 #define MD_DDBONDCOMM (1<<11)
68 #define MD_CONFOUT (1<<12)
69 #define MD_REPRODUCIBLE (1<<13)
70 #define MD_READ_RNG (1<<14)
71 #define MD_APPENDFILES (1<<15)
72 #define MD_KEEPANDNUMCPT (1<<16)
73 #define MD_READ_EKIN (1<<17)
74 #define MD_STARTFROMCPT (1<<18)
75 #define MD_RESETCOUNTERSHALFWAY (1<<19)
77 /* Define a number of flags to better control the information
78 * passed to compute_globals in md.c and global_stat.
81 /* We are rerunning the simulation */
82 #define CGLO_RERUNMD (1<<1)
83 /* we are computing the kinetic energy from average velocities */
84 #define CGLO_EKINAVEVEL (1<<2)
85 /* we are removing the center of mass momenta */
86 #define CGLO_STOPCM (1<<3)
87 /* bGStat is defined in do_md */
88 #define CGLO_GSTAT (1<<4)
89 /* Sum the energy terms in global computation */
90 #define CGLO_ENERGY (1<<6)
91 /* Sum the kinetic energy terms in global computation */
92 #define CGLO_TEMPERATURE (1<<7)
93 /* Sum the kinetic energy terms in global computation */
94 #define CGLO_PRESSURE (1<<8)
95 /* Sum the constraint term in global computation */
96 #define CGLO_CONSTRAINT (1<<9)
97 /* we are using an integrator that requires iteration over some steps - currently not used*/
98 #define CGLO_ITERATE (1<<10)
99 /* it is the first time we are iterating (or, only once through is required */
100 #define CGLO_FIRSTITERATE (1<<11)
101 /* Reading ekin from the trajectory */
102 #define CGLO_READEKIN (1<<12)
103 /* we need to reset the ekin rescaling factor here */
104 #define CGLO_SCALEEKIN (1<<13)
106 enum {
107 ddnoSEL, ddnoINTERLEAVE, ddnoPP_PME, ddnoCARTESIAN, ddnoNR
110 typedef struct {
111 double real;
112 #ifdef GMX_CRAY_XT3
113 double proc;
114 #else
115 clock_t proc;
116 #endif
117 double realtime;
118 double proctime;
119 double time_per_step;
120 double last;
121 gmx_large_int_t nsteps_done;
122 } gmx_runtime_t;
124 typedef struct {
125 t_fileio *fp_trn;
126 t_fileio *fp_xtc;
127 int xtc_prec;
128 ener_file_t fp_ene;
129 const char *fn_cpt;
130 gmx_bool bKeepAndNumCPT;
131 int eIntegrator;
132 int simulation_part;
133 FILE *fp_dhdl;
134 FILE *fp_field;
135 } gmx_mdoutf_t;
137 /* Variables for temporary use with the deform option,
138 * used in runner.c and md.c.
139 * (These variables should be stored in the tpx file.)
141 extern gmx_large_int_t deform_init_init_step_tpx;
142 extern matrix deform_init_box_tpx;
143 #ifdef GMX_THREADS
144 extern tMPI_Thread_mutex_t deform_init_box_mutex;
146 /* The minimum number of atoms per thread. With fewer atoms than this,
147 * the number of threads will get lowered.
149 #define MIN_ATOMS_PER_THREAD 90
150 #endif
153 typedef double gmx_integrator_t(FILE *log,t_commrec *cr,
154 int nfile,const t_filenm fnm[],
155 const output_env_t oenv, gmx_bool bVerbose,
156 gmx_bool bCompact, int nstglobalcomm,
157 gmx_vsite_t *vsite,gmx_constr_t constr,
158 int stepout,
159 t_inputrec *inputrec,
160 gmx_mtop_t *mtop,t_fcdata *fcd,
161 t_state *state,
162 t_mdatoms *mdatoms,
163 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
164 gmx_edsam_t ed,
165 t_forcerec *fr,
166 int repl_ex_nst,int repl_ex_seed,
167 gmx_membed_t *membed,
168 real cpt_period,real max_hours,
169 const char *deviceOptions,
170 unsigned long Flags,
171 gmx_runtime_t *runtime);
173 typedef struct gmx_global_stat *gmx_global_stat_t;
175 /* ROUTINES from md.c */
177 gmx_integrator_t do_md;
179 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 sim_util.c */
202 void do_pbc_first(FILE *log,matrix box,t_forcerec *fr,
203 t_graph *graph,rvec x[]);
205 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
206 gmx_mtop_t *mtop,rvec x[]);
208 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
209 gmx_mtop_t *mtop,rvec x[]);
212 /* ROUTINES from stat.c */
213 gmx_global_stat_t global_stat_init(t_inputrec *ir);
215 void global_stat_destroy(gmx_global_stat_t gs);
217 void global_stat(FILE *log,gmx_global_stat_t gs,
218 t_commrec *cr,gmx_enerdata_t *enerd,
219 tensor fvir,tensor svir,rvec mu_tot,
220 t_inputrec *inputrec,
221 gmx_ekindata_t *ekind,
222 gmx_constr_t constr,t_vcm *vcm,
223 int nsig,real *sig,
224 gmx_mtop_t *top_global, t_state *state_local,
225 gmx_bool bSumEkinhOld, int flags);
226 /* Communicate statistics over cr->mpi_comm_mysim */
228 gmx_mdoutf_t *init_mdoutf(int nfile,const t_filenm fnm[],
229 int mdrun_flags,
230 const t_commrec *cr,const t_inputrec *ir,
231 const output_env_t oenv);
232 /* Returns a pointer to a data structure with all output file pointers
233 * and names required by mdrun.
236 void done_mdoutf(gmx_mdoutf_t *of);
237 /* Close all open output files and free the of pointer */
239 #define MDOF_X (1<<0)
240 #define MDOF_V (1<<1)
241 #define MDOF_F (1<<2)
242 #define MDOF_XTC (1<<3)
243 #define MDOF_CPT (1<<4)
245 void write_traj(FILE *fplog,t_commrec *cr,
246 gmx_mdoutf_t *of,
247 int mdof_flags,
248 gmx_mtop_t *top_global,
249 gmx_large_int_t step,double t,
250 t_state *state_local,t_state *state_global,
251 rvec *f_local,rvec *f_global,
252 int *n_xtc,rvec **x_xtc);
253 /* Routine that writes frames to trn, xtc and/or checkpoint.
254 * What is written is determined by the mdof_flags defined above.
255 * Data is collected to the master node only when necessary.
258 int do_per_step(gmx_large_int_t step,gmx_large_int_t nstep);
259 /* Return TRUE if io should be done */
261 int do_any_io(int step, t_inputrec *ir);
263 /* ROUTINES from sim_util.c */
265 double gmx_gettime();
267 void print_time(FILE *out, gmx_runtime_t *runtime,
268 gmx_large_int_t step,t_inputrec *ir, t_commrec *cr);
270 void runtime_start(gmx_runtime_t *runtime);
272 void runtime_end(gmx_runtime_t *runtime);
274 void runtime_upd_proc(gmx_runtime_t *runtime);
275 /* The processor time should be updated every once in a while,
276 * since on 32-bit manchines it loops after 72 minutes.
279 void print_date_and_time(FILE *log,int pid,const char *title,
280 const gmx_runtime_t *runtime);
282 void nstop_cm(FILE *log,t_commrec *cr,
283 int start,int nr_atoms,real mass[],rvec x[],rvec v[]);
285 void finish_run(FILE *log,t_commrec *cr,const char *confout,
286 t_inputrec *inputrec,
287 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
288 gmx_runtime_t *runtime,
289 gmx_bool bWriteStat);
291 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr);
293 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
294 gmx_large_int_t step, int natoms,
295 matrix box,real lambda,tensor pres,tensor virial,
296 real *prescorr, real *enercorr, real *dvdlcorr);
298 typedef enum
300 LIST_SCALARS =0001,
301 LIST_INPUTREC =0002,
302 LIST_TOP =0004,
303 LIST_X =0010,
304 LIST_V =0020,
305 LIST_F =0040,
306 LIST_LOAD =0100
307 } t_listitem;
309 void check_nnodes_top(char *fn,t_topology *top);
310 /* Reset the tpr file to work with one node if necessary */
313 /* check the version */
314 void check_ir_old_tpx_versions(t_commrec *cr,FILE *fplog,
315 t_inputrec *ir,gmx_mtop_t *mtop);
317 /* Allocate and initialize node-local state entries. */
318 void set_state_entries(t_state *state,const t_inputrec *ir,int nnodes);
320 /* Broadcast the data for a simulation, and allocate node-specific settings
321 such as rng generators. */
322 void init_parallel(FILE *log, t_commrec *cr, t_inputrec *inputrec,
323 gmx_mtop_t *mtop);
326 void do_constrain_first(FILE *log,gmx_constr_t constr,
327 t_inputrec *inputrec,t_mdatoms *md,
328 t_state *state,rvec *f,
329 t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
330 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir);
332 void dynamic_load_balancing(gmx_bool bVerbose,t_commrec *cr,real capacity[],
333 int dimension,t_mdatoms *md,t_topology *top,
334 rvec x[],rvec v[],matrix box);
335 /* Perform load balancing, i.e. split the particles over processors
336 * based on their coordinates in the "dimension" direction.
339 int mdrunner(int nthreads_requested, FILE *fplog,t_commrec *cr,int nfile,
340 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
341 gmx_bool bCompact, int nstglobalcomm, ivec ddxyz,int dd_node_order,
342 real rdd, real rconstr, const char *dddlb_opt,real dlb_scale,
343 const char *ddcsx,const char *ddcsy,const char *ddcsz,
344 int nstepout, int resetstep, int nmultisim, int repl_ex_nst,
345 int repl_ex_seed, real pforce,real cpt_period,real max_hours,
346 const char *deviceOptions, unsigned long Flags);
347 /* Driver routine, that calls the different methods */
349 void md_print_warning(const t_commrec *cr,FILE *fplog,const char *buf);
350 /* Print a warning message to stderr on the master node
351 * and to fplog if fplog!=NULL.
354 void init_md(FILE *fplog,
355 t_commrec *cr,t_inputrec *ir, const output_env_t oenv,
356 double *t,double *t0,
357 real *lambda,double *lam0,
358 t_nrnb *nrnb,gmx_mtop_t *mtop,
359 gmx_update_t *upd,
360 int nfile,const t_filenm fnm[],
361 gmx_mdoutf_t **outf,t_mdebin **mdebin,
362 tensor force_vir,tensor shake_vir,
363 rvec mu_tot,
364 gmx_bool *bSimAnn,t_vcm **vcm,
365 t_state *state, unsigned long Flags);
366 /* Routine in sim_util.c */
368 #ifdef __cplusplus
370 #endif
372 #endif /* _mdrun_h */