Merge branch 'master' of git.gromacs.org:gromacs
[gromacs/rigid-bodies.git] / include / mdrun.h
blob13fc174dfd5583746375b7d3329cd11449998948
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 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
43 #include <stdio.h>
44 #include <time.h>
45 #include "typedefs.h"
46 #include "network.h"
47 #include "tgroup.h"
48 #include "filenm.h"
49 #include "mshift.h"
50 #include "force.h"
51 #include "edsam.h"
52 #include "mdebin.h"
53 #include "vcm.h"
54 #include "vsite.h"
55 #include "pull.h"
56 #include "update.h"
58 #ifdef __cplusplus
59 extern "C" {
60 #endif
62 #define MD_POLARISE (1<<2)
63 #define MD_IONIZE (1<<3)
64 #define MD_RERUN (1<<4)
65 #define MD_RERUN_VSITE (1<<5)
66 #define MD_FFSCAN (1<<6)
67 #define MD_SEPPOT (1<<7)
68 #define MD_PARTDEC (1<<9)
69 #define MD_DDBONDCHECK (1<<10)
70 #define MD_DDBONDCOMM (1<<11)
71 #define MD_CONFOUT (1<<12)
72 #define MD_REPRODUCIBLE (1<<13)
73 #define MD_READ_RNG (1<<14)
74 #define MD_APPENDFILES (1<<15)
75 #define MD_KEEPANDNUMCPT (1<<16)
76 #define MD_READ_EKIN (1<<17)
77 #define MD_STARTFROMCPT (1<<18)
78 #define MD_RESETCOUNTERSHALFWAY (1<<19)
80 /* Define a number of flags to better control the information
81 * passed to compute_globals in md.c and global_stat.
84 /* We are rerunning the simulation */
85 #define CGLO_RERUNMD (1<<1)
86 /* we are computing the kinetic energy from average velocities */
87 #define CGLO_EKINAVEVEL (1<<2)
88 /* we are removing the center of mass momenta */
89 #define CGLO_STOPCM (1<<3)
90 /* bGStat is defined in do_md */
91 #define CGLO_GSTAT (1<<4)
92 /* bNEMD is defined in do_md */
93 #define CGLO_NEMD (1<<5)
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 bool bKeepAndNumCPT;
136 int eIntegrator;
137 int simulation_part;
138 FILE *fp_dhdl;
139 FILE *fp_field;
140 } gmx_mdoutf_t;
142 /* Variables for temporary use with the deform option,
143 * used in runner.c and md.c.
144 * (These variables should be stored in the tpx file.)
146 extern gmx_large_int_t deform_init_init_step_tpx;
147 extern matrix deform_init_box_tpx;
148 #ifdef GMX_THREADS
149 extern tMPI_Thread_mutex_t deform_init_box_mutex;
151 /* The minimum number of atoms per thread. With fewer atoms than this,
152 * the number of threads will get lowered.
154 #define MIN_ATOMS_PER_THREAD 90
155 #endif
158 typedef double gmx_integrator_t(FILE *log,t_commrec *cr,
159 int nfile,const t_filenm fnm[],
160 const output_env_t oenv, bool bVerbose,
161 bool bCompact, int nstglobalcomm,
162 gmx_vsite_t *vsite,gmx_constr_t constr,
163 int stepout,
164 t_inputrec *inputrec,
165 gmx_mtop_t *mtop,t_fcdata *fcd,
166 t_state *state,
167 t_mdatoms *mdatoms,
168 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
169 gmx_edsam_t ed,
170 t_forcerec *fr,
171 int repl_ex_nst,int repl_ex_seed,
172 real cpt_period,real max_hours,
173 const char *deviceOptions,
174 unsigned long Flags,
175 gmx_runtime_t *runtime);
177 typedef struct gmx_global_stat *gmx_global_stat_t;
179 /* ROUTINES from md.c */
181 extern gmx_integrator_t do_md;
183 extern gmx_integrator_t do_md_openmm;
185 /* ROUTINES from minimize.c */
187 extern gmx_integrator_t do_steep;
188 /* Do steepest descents EM */
190 extern gmx_integrator_t do_cg;
191 /* Do conjugate gradient EM */
193 extern gmx_integrator_t do_lbfgs;
194 /* Do conjugate gradient L-BFGS */
196 extern gmx_integrator_t do_nm;
197 /* Do normal mode analysis */
199 /* ROUTINES from tpi.c */
201 extern gmx_integrator_t do_tpi;
202 /* Do test particle insertion */
205 /* ROUTINES from sim_util.c */
206 extern void do_pbc_first(FILE *log,matrix box,t_forcerec *fr,
207 t_graph *graph,rvec x[]);
209 extern void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
210 gmx_mtop_t *mtop,rvec x[]);
212 extern void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
213 gmx_mtop_t *mtop,rvec x[]);
216 /* ROUTINES from stat.c */
217 extern gmx_global_stat_t global_stat_init(t_inputrec *ir);
219 extern void global_stat_destroy(gmx_global_stat_t gs);
221 extern void global_stat(FILE *log,gmx_global_stat_t gs,
222 t_commrec *cr,gmx_enerdata_t *enerd,
223 tensor fvir,tensor svir,rvec mu_tot,
224 t_inputrec *inputrec,
225 gmx_ekindata_t *ekind,
226 gmx_constr_t constr,t_vcm *vcm,
227 int nsig,real *sig,
228 gmx_mtop_t *top_global, t_state *state_local,
229 bool bSumEkinhOld, int flags);
230 /* Communicate statistics over cr->mpi_comm_mysim */
232 extern gmx_mdoutf_t *init_mdoutf(int nfile,const t_filenm fnm[],
233 int mdrun_flags,
234 const t_commrec *cr,const t_inputrec *ir,
235 const output_env_t oenv);
236 /* Returns a pointer to a data structure with all output file pointers
237 * and names required by mdrun.
240 extern void done_mdoutf(gmx_mdoutf_t *of);
241 /* Close all open output files and free the of pointer */
243 #define MDOF_X (1<<0)
244 #define MDOF_V (1<<1)
245 #define MDOF_F (1<<2)
246 #define MDOF_XTC (1<<3)
247 #define MDOF_CPT (1<<4)
249 extern void write_traj(FILE *fplog,t_commrec *cr,
250 gmx_mdoutf_t *of,
251 int mdof_flags,
252 gmx_mtop_t *top_global,
253 gmx_large_int_t step,double t,
254 t_state *state_local,t_state *state_global,
255 rvec *f_local,rvec *f_global,
256 int *n_xtc,rvec **x_xtc);
257 /* Routine that writes frames to trn, xtc and/or checkpoint.
258 * What is written is determined by the mdof_flags defined above.
259 * Data is collected to the master node only when necessary.
262 extern int do_per_step(gmx_large_int_t step,gmx_large_int_t nstep);
263 /* Return TRUE if io should be done */
265 extern int do_any_io(int step, t_inputrec *ir);
267 /* ROUTINES from sim_util.c */
269 extern double gmx_gettime();
271 extern void print_time(FILE *out, gmx_runtime_t *runtime,
272 gmx_large_int_t step,t_inputrec *ir, t_commrec *cr);
274 extern void runtime_start(gmx_runtime_t *runtime);
276 extern void runtime_end(gmx_runtime_t *runtime);
278 extern void runtime_upd_proc(gmx_runtime_t *runtime);
279 /* The processor time should be updated every once in a while,
280 * since on 32-bit manchines it loops after 72 minutes.
283 extern void print_date_and_time(FILE *log,int pid,const char *title,
284 const gmx_runtime_t *runtime);
286 extern void nstop_cm(FILE *log,t_commrec *cr,
287 int start,int nr_atoms,real mass[],rvec x[],rvec v[]);
289 extern void finish_run(FILE *log,t_commrec *cr,const char *confout,
290 t_inputrec *inputrec,
291 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
292 gmx_runtime_t *runtime,
293 bool bWriteStat);
295 extern void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr);
297 extern void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
298 gmx_large_int_t step, int natoms,
299 matrix box,real lambda,tensor pres,tensor virial,
300 real *prescorr, real *enercorr, real *dvdlcorr);
302 typedef enum
304 LIST_SCALARS =0001,
305 LIST_INPUTREC =0002,
306 LIST_TOP =0004,
307 LIST_X =0010,
308 LIST_V =0020,
309 LIST_F =0040,
310 LIST_LOAD =0100
311 } t_listitem;
313 extern void check_nnodes_top(char *fn,t_topology *top);
314 /* Reset the tpr file to work with one node if necessary */
317 /* check the version */
318 void check_ir_old_tpx_versions(t_commrec *cr,FILE *fplog,
319 t_inputrec *ir,gmx_mtop_t *mtop);
321 /* Allocate and initialize node-local state entries. */
322 void set_state_entries(t_state *state,t_inputrec *ir,int nnodes);
324 /* Broadcast the data for a simulation, and allocate node-specific settings
325 such as rng generators. */
326 extern void init_parallel(FILE *log, t_commrec *cr, t_inputrec *inputrec,
327 gmx_mtop_t *mtop, t_state *state);
330 extern void do_constrain_first(FILE *log,gmx_constr_t constr,
331 t_inputrec *inputrec,t_mdatoms *md,
332 t_state *state,rvec *f,
333 t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
334 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir);
336 extern void dynamic_load_balancing(bool bVerbose,t_commrec *cr,real capacity[],
337 int dimension,t_mdatoms *md,t_topology *top,
338 rvec x[],rvec v[],matrix box);
339 /* Perform load balancing, i.e. split the particles over processors
340 * based on their coordinates in the "dimension" direction.
343 int mdrunner(int nthreads_requested, FILE *fplog,t_commrec *cr,int nfile,
344 const t_filenm fnm[], const output_env_t oenv, bool bVerbose,
345 bool bCompact, int nstglobalcomm, ivec ddxyz,int dd_node_order,
346 real rdd, real rconstr, const char *dddlb_opt,real dlb_scale,
347 const char *ddcsx,const char *ddcsy,const char *ddcsz,
348 int nstepout, int resetstep, int nmultisim, int repl_ex_nst,
349 int repl_ex_seed, real pforce,real cpt_period,real max_hours,
350 const char *deviceOptions, unsigned long Flags);
351 /* Driver routine, that calls the different methods */
353 extern void md_print_warning(const t_commrec *cr,FILE *fplog,const char *buf);
354 /* Print a warning message to stderr on the master node
355 * and to fplog if fplog!=NULL.
358 extern void init_md(FILE *fplog,
359 t_commrec *cr,t_inputrec *ir, const output_env_t oenv,
360 double *t,double *t0,
361 real *lambda,double *lam0,
362 t_nrnb *nrnb,gmx_mtop_t *mtop,
363 gmx_update_t *upd,
364 int nfile,const t_filenm fnm[],
365 gmx_mdoutf_t **outf,t_mdebin **mdebin,
366 tensor force_vir,tensor shake_vir,
367 rvec mu_tot,
368 bool *bNEMD,bool *bSimAnn,t_vcm **vcm,
369 t_state *state, unsigned long Flags);
370 /* Routine in sim_util.c */
372 #ifdef __cplusplus
374 #endif
376 #endif /* _mdrun_h */