Fixes for Windows/ICC--typecasting
[gromacs.git] / include / mdrun.h
blobc4b4e783173fb4fa77245fae15e45a7cc815dcf5
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 "typedefs.h"
45 #include "network.h"
46 #include "tgroup.h"
47 #include "filenm.h"
48 #include "mshift.h"
49 #include "force.h"
50 #include "time.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<<14)
73 #define MD_READ_RNG (1<<15)
74 #define MD_APPENDFILES (1<<16)
75 #define MD_READ_EKIN (1<<17)
76 #define MD_STARTFROMCPT (1<<18)
77 #define MD_RESETCOUNTERSHALFWAY (1<<19)
79 /* Define a number of flags to better control the information
80 * passed to compute_globals in md.c and global_stat.
83 /* We are rerunning the simulation */
84 #define CGLO_RERUNMD (1<<1)
85 /* we are computing the kinetic energy from average velocities */
86 #define CGLO_EKINAVEVEL (1<<2)
87 /* we are removing the center of mass momenta */
88 #define CGLO_STOPCM (1<<3)
89 /* bGStat is defined in do_md */
90 #define CGLO_GSTAT (1<<4)
91 /* bNEMD is defined in do_md */
92 #define CGLO_NEMD (1<<5)
93 /* Sum the energy terms in global computation */
94 #define CGLO_ENERGY (1<<6)
95 /* Sum the kinetic energy terms in global computation */
96 #define CGLO_TEMPERATURE (1<<7)
97 /* Sum the kinetic energy terms in global computation */
98 #define CGLO_PRESSURE (1<<8)
99 /* Sum the constraint term in global computation */
100 #define CGLO_CONSTRAINT (1<<9)
101 /* we are using an integrator that requires iteration over some steps - currently not used*/
102 #define CGLO_ITERATE (1<<10)
103 /* it is the first time we are iterating (or, only once through is required */
104 #define CGLO_FIRSTITERATE (1<<11)
105 /* Reading ekin from the trajectory */
106 #define CGLO_READEKIN (1<<12)
107 /* we need to reset the ekin rescaling factor here */
108 #define CGLO_SCALEEKIN (1<<13)
110 enum {
111 ddnoSEL, ddnoINTERLEAVE, ddnoPP_PME, ddnoCARTESIAN, ddnoNR
114 typedef struct {
115 double real;
116 #ifdef GMX_CRAY_XT3
117 double proc;
118 #else
119 clock_t proc;
120 #endif
121 double realtime;
122 double proctime;
123 double time_per_step;
124 double last;
125 gmx_large_int_t nsteps_done;
126 } gmx_runtime_t;
128 typedef struct {
129 int fp_trn;
130 int fp_xtc;
131 int xtc_prec;
132 ener_file_t fp_ene;
133 const char *fn_cpt;
134 int eIntegrator;
135 int simulation_part;
136 FILE *fp_dhdl;
137 FILE *fp_field;
138 } gmx_mdoutf_t;
140 /* Variables for temporary use with the deform option,
141 * used in runner.c and md.c.
142 * (These variables should be stored in the tpx file.)
144 extern gmx_large_int_t deform_init_init_step_tpx;
145 extern matrix deform_init_box_tpx;
146 #ifdef GMX_THREADS
147 extern tMPI_Thread_mutex_t deform_init_box_mutex;
148 #endif
151 typedef double gmx_integrator_t(FILE *log,t_commrec *cr,
152 int nfile,const t_filenm fnm[],
153 const output_env_t oenv, bool bVerbose,
154 bool bCompact, int nstglobalcomm,
155 gmx_vsite_t *vsite,gmx_constr_t constr,
156 int stepout,
157 t_inputrec *inputrec,
158 gmx_mtop_t *mtop,t_fcdata *fcd,
159 t_state *state,
160 t_mdatoms *mdatoms,
161 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
162 gmx_edsam_t ed,
163 t_forcerec *fr,
164 int repl_ex_nst,int repl_ex_seed,
165 real cpt_period,real max_hours,
166 const char *deviceOptions,
167 unsigned long Flags,
168 gmx_runtime_t *runtime);
170 typedef struct gmx_global_stat *gmx_global_stat_t;
172 /* ROUTINES from md.c */
174 extern gmx_integrator_t do_md;
176 extern gmx_integrator_t do_md_openmm;
178 /* ROUTINES from minimize.c */
180 extern gmx_integrator_t do_steep;
181 /* Do steepest descents EM */
183 extern gmx_integrator_t do_cg;
184 /* Do conjugate gradient EM */
186 extern gmx_integrator_t do_lbfgs;
187 /* Do conjugate gradient L-BFGS */
189 extern gmx_integrator_t do_nm;
190 /* Do normal mode analysis */
192 /* ROUTINES from tpi.c */
194 extern gmx_integrator_t do_tpi;
195 /* Do test particle insertion */
198 /* ROUTINES from sim_util.c */
199 extern void do_pbc_first(FILE *log,matrix box,t_forcerec *fr,
200 t_graph *graph,rvec x[]);
202 extern void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
203 gmx_mtop_t *mtop,rvec x[]);
205 extern void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
206 gmx_mtop_t *mtop,rvec x[]);
209 /* ROUTINES from stat.c */
210 extern gmx_global_stat_t global_stat_init(t_inputrec *ir);
212 extern void global_stat_destroy(gmx_global_stat_t gs);
214 extern void global_stat(FILE *log,gmx_global_stat_t gs,
215 t_commrec *cr,gmx_enerdata_t *enerd,
216 tensor fvir,tensor svir,rvec mu_tot,
217 t_inputrec *inputrec,
218 gmx_ekindata_t *ekind,
219 gmx_constr_t constr,t_vcm *vcm,
220 int nsig,real *sig,
221 gmx_mtop_t *top_global, t_state *state_local,
222 bool bSumEkinhOld, int flags);
223 /* Communicate statistics over cr->mpi_comm_mysim */
225 extern gmx_mdoutf_t *init_mdoutf(int nfile,const t_filenm fnm[],
226 bool bAppendFiles,
227 const t_commrec *cr,const t_inputrec *ir,
228 const output_env_t oenv);
229 /* Returns a pointer to a data structure with all output file pointers
230 * and names required by mdrun.
233 extern void done_mdoutf(gmx_mdoutf_t *of);
234 /* Close all open output files and free the of pointer */
236 #define MDOF_X (1<<0)
237 #define MDOF_V (1<<1)
238 #define MDOF_F (1<<2)
239 #define MDOF_XTC (1<<3)
240 #define MDOF_CPT (1<<4)
242 extern void write_traj(FILE *fplog,t_commrec *cr,
243 gmx_mdoutf_t *of,
244 int mdof_flags,
245 gmx_mtop_t *top_global,
246 gmx_large_int_t step,double t,
247 t_state *state_local,t_state *state_global,
248 rvec *f_local,rvec *f_global,
249 int *n_xtc,rvec **x_xtc);
250 /* Routine that writes frames to trn, xtc and/or checkpoint.
251 * What is written is determined by the mdof_flags defined above.
252 * Data is collected to the master node only when necessary.
255 extern int do_per_step(gmx_large_int_t step,gmx_large_int_t nstep);
256 /* Return TRUE if io should be done */
258 extern int do_any_io(int step, t_inputrec *ir);
260 /* ROUTINES from sim_util.c */
262 extern double gmx_gettime();
264 extern void print_time(FILE *out, gmx_runtime_t *runtime,
265 gmx_large_int_t step,t_inputrec *ir, t_commrec *cr);
267 extern void runtime_start(gmx_runtime_t *runtime);
269 extern void runtime_end(gmx_runtime_t *runtime);
271 extern void runtime_upd_proc(gmx_runtime_t *runtime);
272 /* The processor time should be updated every once in a while,
273 * since on 32-bit manchines it loops after 72 minutes.
276 extern void print_date_and_time(FILE *log,int pid,const char *title,
277 const gmx_runtime_t *runtime);
279 extern void nstop_cm(FILE *log,t_commrec *cr,
280 int start,int nr_atoms,real mass[],rvec x[],rvec v[]);
282 extern void finish_run(FILE *log,t_commrec *cr,const char *confout,
283 t_inputrec *inputrec,
284 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
285 gmx_runtime_t *runtime,
286 bool bWriteStat);
288 extern void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr);
290 extern void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
291 gmx_large_int_t step, int natoms,
292 matrix box,real lambda,tensor pres,tensor virial,
293 real *prescorr, real *enercorr, real *dvdlcorr);
295 typedef enum
297 LIST_SCALARS =0001,
298 LIST_INPUTREC =0002,
299 LIST_TOP =0004,
300 LIST_X =0010,
301 LIST_V =0020,
302 LIST_F =0040,
303 LIST_LOAD =0100
304 } t_listitem;
306 extern void check_nnodes_top(char *fn,t_topology *top);
307 /* Reset the tpr file to work with one node if necessary */
309 extern void init_single(FILE *log, t_inputrec *inputrec, const char *tpbfile,
310 gmx_mtop_t *mtop, t_state *state);
312 * Allocates space for the topology (top), the coordinates x, the
313 * velocities v, masses mass. Reads the parameters, topology,
314 * coordinates and velocities from the file specified in tpbfile
317 /* check the version */
318 void check_ir_old_tpx_versions(t_commrec *cr,FILE *fplog,
319 t_inputrec *ir,gmx_mtop_t *mtop);
322 extern void init_parallel(FILE *log,const char *tpxfile, t_commrec *cr,
323 t_inputrec *inputrec,gmx_mtop_t *mtop,
324 t_state *state, int list);
326 * Loads the data for a simulation from the ring. Parameters, topology
327 * coordinates, velocities, and masses are initialised equal to using
328 * init_single() in the single processor version. The extra argument
329 * f_add is allocated to use for the update of the forces, the load
330 * array specifies in which part of the x and f array the subsystems
331 * of the other processors are located. Homenr0, homenr1, nparts0 and
332 * nparts1 are necessary to calculate the non bonded interaction using
333 * the symmetry and thus calculating every force only once. List is a
334 * facility for logging (and debugging). One can decide to print none or a
335 * set of * selected parameters to the file specified by log. Parameters
336 * are printed by or-ing the corresponding items from t_listitem. A 0
337 * (zero) specifies that nothing is to be printed on the file. The function
338 * returns the number of shifts over the ring to perform to calculate
339 * all interactions.
341 * NOTE: for threaded simulations that don't support parallel runs (at
342 * the moment that's only the LBGFS integrator), this function may
343 * cancel them and re-write the commrec.
346 extern void do_constrain_first(FILE *log,gmx_constr_t constr,
347 t_inputrec *inputrec,t_mdatoms *md,
348 t_state *state,rvec *f,
349 t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
350 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir);
352 extern void dynamic_load_balancing(bool bVerbose,t_commrec *cr,real capacity[],
353 int dimension,t_mdatoms *md,t_topology *top,
354 rvec x[],rvec v[],matrix box);
355 /* Perform load balancing, i.e. split the particles over processors
356 * based on their coordinates in the "dimension" direction.
359 int mdrunner(FILE *fplog,t_commrec *cr,int nfile,const t_filenm fnm[],
360 const output_env_t oenv, bool bVerbose,bool bCompact,
361 int nstglobalcomm, ivec ddxyz,int dd_node_order,real rdd,
362 real rconstr, const char *dddlb_opt,real dlb_scale,
363 const char *ddcsx,const char *ddcsy,const char *ddcsz,
364 int nstepout, int resetstep, int nmultisim, int repl_ex_nst,int repl_ex_seed,
365 real pforce,real cpt_period,real max_hours,
366 const char *deviceOptions, unsigned long Flags);
367 /* Driver routine, that calls the different methods */
369 int mdrunner_threads(int nthreads,
370 FILE *fplog,t_commrec *cr,int nfile,const t_filenm fnm[],
371 const output_env_t oenv, bool bVerbose,bool bCompact,
372 int nstglobalcomm,
373 ivec ddxyz,int dd_node_order,real rdd,real rconstr,
374 const char *dddlb_opt,real dlb_scale,
375 const char *ddcsx,const char *ddcsy,const char *ddcsz,
376 int nstepout,int resetstep,int nmultisim, int repl_ex_nst,
377 int repl_ex_seed, real pforce,real cpt_period,
378 real max_hours, const char *deviceOptions, unsigned long Flags);
379 /* initializes nthread threads before running mdrunner: is the preferred
380 way to start a simulation (even if nthreads=1 and no threads are started) */
382 extern void md_print_warning(const t_commrec *cr,FILE *fplog,const char *buf);
383 /* Print a warning message to stderr on the master node
384 * and to fplog if fplog!=NULL.
387 extern void init_md(FILE *fplog,
388 t_commrec *cr,t_inputrec *ir, const output_env_t oenv,
389 double *t,double *t0,
390 real *lambda,double *lam0,
391 t_nrnb *nrnb,gmx_mtop_t *mtop,
392 gmx_update_t *upd,
393 int nfile,const t_filenm fnm[],
394 gmx_mdoutf_t **outf,t_mdebin **mdebin,
395 tensor force_vir,tensor shake_vir,
396 rvec mu_tot,
397 bool *bNEMD,bool *bSimAnn,t_vcm **vcm,
398 t_state *state, unsigned long Flags);
399 /* Routine in sim_util.c */
401 #ifdef __cplusplus
403 #endif
405 #endif /* _mdrun_h */