Removing scripts and introducing wrapper binaries instead.
[gromacs.git] / include / mdrun.h
blob3cc7768b106e8e041fa4dac602568110f83d7706
1 /*
2 * $Id$
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * Gromacs Runs On Most of All Computer Systems
37 #ifndef _mdrun_h
38 #define _mdrun_h
40 #ifdef HAVE_CONFIG_H
41 #include <config.h>
42 #endif
44 #include <stdio.h>
45 #include "typedefs.h"
46 #include "network.h"
47 #include "tgroup.h"
48 #include "filenm.h"
49 #include "nsb.h"
50 #include "mshift.h"
51 #include "force.h"
52 #include "time.h"
53 #include "edsam.h"
54 #include "mdebin.h"
55 #include "vcm.h"
56 #include "dummies.h"
57 #include "pull.h"
59 #define MD_MULTISIM (1<<0)
60 #define MD_GLAS (1<<1)
61 #define MD_POLARISE (1<<2)
62 #define MD_IONIZE (1<<3)
63 #define MD_RERUN (1<<4)
64 #define MD_FFSCAN (1<<6)
65 #define MD_SEPDVDL (1<<7)
67 /* ROUTINES from md.c */
68 extern time_t do_md(FILE *log,t_commrec *cr,t_commrec *mcr,
69 int nfile,t_filenm fnm[],
70 bool bVerbose,bool bCompact,bool bDummies,
71 t_comm_dummies *dummycomm,int stepout,
72 t_parm *parm,t_groups *grps,
73 t_topology *top,real ener[],t_fcdata *fcd,
74 t_state *state,rvec vold[],rvec vt[],rvec f[],
75 rvec buf[],t_mdatoms *mdatoms,
76 t_nsborder *nsb,t_nrnb nrnb[],
77 t_graph *graph,t_edsamyn *edyn,
78 t_forcerec *fr,rvec box_size,
79 unsigned long Flags);
81 /* ROUTINES from minimize.c */
82 extern time_t do_steep(FILE *log,int nfile,t_filenm fnm[],
83 t_parm *parm,t_topology *top,
84 t_groups *grps,t_nsborder *nsb,
85 t_state *state,rvec grad[],rvec buf[],t_mdatoms *mdatoms,
86 tensor ekin,real ener[],t_fcdata *fcd,t_nrnb nrnb[],
87 bool bVerbose,bool bDummies,t_comm_dummies *dummycomm,
88 t_commrec *cr,t_commrec *mcr,
89 t_graph *graph,t_forcerec *fr,rvec box_size);
90 /* Do steepest descents EM or something like that! */
92 extern time_t do_cg(FILE *log,int nfile,t_filenm fnm[],
93 t_parm *parm,t_topology *top,
94 t_groups *grps,t_nsborder *nsb,
95 t_state *state,rvec grad[],rvec buf[],t_mdatoms *mdatoms,
96 tensor ekin,real ener[],t_fcdata *fcd,t_nrnb nrnb[],
97 bool bVerbose,bool bDummies,t_comm_dummies *dummycomm,
98 t_commrec *cr,t_commrec *mcr,
99 t_graph *graph,t_forcerec *fr,rvec box_size);
100 /* Do conjugate gradients EM! */
102 extern time_t do_lbfgs(FILE *log,int nfile,t_filenm fnm[],
103 t_parm *parm,t_topology *top,
104 t_groups *grps,t_nsborder *nsb, t_state *state,
105 rvec grad[],rvec buf[],t_mdatoms *mdatoms,
106 tensor ekin,real ener[],t_fcdata *fcd,t_nrnb nrnb[],
107 bool bVerbose,bool bDummies,t_comm_dummies *dummycomm,
108 t_commrec *cr,t_commrec *mcr,
109 t_graph *graph,t_forcerec *fr,rvec box_size);
110 /* Do conjugate gradients EM! */
113 extern time_t do_nm(FILE *log,t_commrec *cr,int nfile,t_filenm fnm[],
114 bool bVerbose,bool bCompact,int stepout,
115 t_parm *parm,t_groups *grps,
116 t_topology *top,real ener[],t_fcdata *fcd,
117 t_state *state,rvec vold[],rvec vt[],rvec f[],
118 rvec buf[],t_mdatoms *mdatoms,
119 t_nsborder *nsb,t_nrnb nrnb[],
120 t_graph *graph,t_edsamyn *edyn,
121 t_forcerec *fr,rvec box_size);
122 /* Do normal mode analysis */
124 /* ROUTINES from runner.c */
125 extern bool optRerunMDset (int nfile, t_filenm fnm[]);
127 extern void do_pbc_first(FILE *log,matrix box,rvec box_size,t_forcerec *fr,
128 t_graph *graph,rvec x[]);
130 extern void set_pot_bools(t_inputrec *ir,t_topology *top,
131 bool *bLR,bool *bLJLR,bool *bBHAM,bool *b14);
132 /* Initiate some bools for the potential energy calculation */
134 /* ROUTINES from stat.c */
135 extern void global_stat(FILE *log,
136 t_commrec *cr,real ener[],
137 tensor fvir,tensor svir,
138 t_grpopts *opts,t_groups *grps,
139 t_nrnb *mynrnb,t_nrnb nrnb[],
140 t_vcm *vcm,real *terminate);
141 /* Communicate statistics around the ring */
143 extern int write_traj(FILE *log,t_commrec *cr,char *traj,t_nsborder *nsb,
144 int step,real t,real lambda,t_nrnb nr_nb[],
145 int natoms,rvec *xx,rvec *vv,rvec *ff,matrix box);
146 /* Routine to output statusfiles during a run, as specified in
147 * in parm->ir. If any of the pointers xx,vv,ff or ener is not NULL
148 * it is written to the trajectory file.
149 * Also write the energies etc. to the log file.
150 * Returns the file handle (to be closed with close_trn).
153 extern int do_per_step(int step,int nstep);
154 /* Return TRUE if io should be done */
156 extern int do_any_io(int step, t_inputrec *ir);
158 extern void write_xtc_traj(FILE *log,t_commrec *cr,
159 char *xtc_traj,t_nsborder *nsb,t_mdatoms *md,
160 int step,real t,rvec *xx,
161 matrix box,real prec);
163 extern void close_xtc_traj(void);
165 /* ROUTINES from sim_util.c */
166 extern void update_mdatoms(t_mdatoms *md,real lambda, bool bFirst);
167 /* Compute fields from mdatoms struct (invmass etc.) which may change
168 * due to lambda dependent FEP calculations.
169 * If bFirst all values are set, this is necessary once in the
170 * first step.
171 * You only have to call this routine again if lambda changes.
174 extern void print_time(FILE *out,time_t start,int step,t_inputrec *ir);
176 extern time_t print_date_and_time(FILE *log,int pid,char *title);
178 extern void do_force(FILE *log,t_commrec *cr,t_commrec *mcr,
179 t_parm *parm,t_nsborder *nsb,
180 tensor vir_part,tensor pme_vir,
181 int step,t_nrnb *nrnb,t_topology *top,t_groups *grps,
182 matrix box,rvec x[],rvec f[],rvec buf[],
183 t_mdatoms *mdatoms,real ener[],t_fcdata *fcd,
184 bool bVerbose,real lambda,t_graph *graph,
185 bool bNS,bool bNBFonly,t_forcerec *fr, rvec mu_tot,
186 bool bGatherOnly,real t,FILE *field);
187 extern void sum_lrforces(rvec f[],t_forcerec *fr,int start,int homenr);
189 extern void calc_virial(FILE *log,int start,int homenr,rvec x[],rvec f[],
190 tensor vir_part,tensor pme_vir,
191 t_graph *graph,matrix box,
192 t_nrnb *nrnb,t_forcerec *fr,bool bTweak);
194 extern void nstop_cm(FILE *log,t_commrec *cr,
195 int start,int nr_atoms,real mass[],rvec x[],rvec v[]);
197 extern void finish_run(FILE *log,t_commrec *cr,char *confout, t_nsborder *nsb,
198 t_topology *top, t_parm *parm,t_nrnb nrnb[],
199 double nodetime,double realtime,int step,
200 bool bWriteStat);
202 extern void calc_dispcorr(FILE *log,int eDispCorr,t_forcerec *fr,int natoms,
203 matrix box,tensor pres,tensor virial,real ener[]);
206 /* STUFF from init.c */
207 extern void write_parm(FILE *log,char *title,int pid,t_parm *parm);
208 /* Write parm for debugging */
210 typedef enum
212 LIST_SCALARS =0001,
213 LIST_PARM =0002,
214 LIST_TOP =0004,
215 LIST_X =0010,
216 LIST_V =0020,
217 LIST_F =0040,
218 LIST_LOAD =0100
219 } t_listitem;
221 extern void init_single(FILE *log,
222 t_parm *parm, char *tpbfile, t_topology *top,
223 t_state *state,t_mdatoms **mdatoms,
224 t_nsborder *nsb);
226 * Allocates space for the topology (top), the coordinates x, the
227 * velocities v, masses mass. Reads the parameters, topology,
228 * coordinates and velocities from the file specified in tpbfile
231 extern void distribute_parts(int left,int right,int pid,int nprocs,
232 t_parm *parm,char *tpbfile,int nstDlb);
234 * Reads the parameters, topology, coordinates and velocities for the
235 * multi processor version of the program from the file specified in
236 * parm->files[STATUS_NM]. This file should also contain a so called
237 * split descriptor which describes how to distribute particles over
238 * the system. It then selects for all subsystems the appropriate data
239 * and sends this to the processor using the left and right channels.
240 * At last it sends its own subsystem down the ring where it is buffered.
241 * Its own buffers for reading the data from the file are freed, and it
242 * is now possible to reload this processor from the ring by using the
243 * init_parts() routine.
244 * The routine also creates a renum array which can be used for writing
245 * out the x,v and f for analysis purpose.
248 extern void init_parts(FILE *log,t_commrec *cr,
249 t_parm *parm,t_topology *top,
250 t_state *state,t_mdatoms **mdatoms,
251 t_nsborder *nsb,int list,
252 bool *bParallelDummies,
253 t_comm_dummies *dummycomm);
255 * Loads the data for a simulation from the ring. Parameters, topology
256 * coordinates, velocities, and masses are initialised equal to using
257 * init_single() in the single processor version. The extra argument
258 * f_add is allocated to use for the update of the forces, the load
259 * array specifies in which part of the x and f array the subsystems
260 * of the other processors are located. Homenr0, homenr1, nparts0 and
261 * nparts1 are necessary to calculate the non bonded interaction using
262 * the symmetry and thus calculating every force only once. List is a facility
263 * for logging (and debugging). One can decide to print none or a set of
264 * selected parameters to the file specified by log. Parameters are
265 * printed by or-ing the corresponding items from t_listitem. A 0 (zero)
266 * specifies that nothing is to be printed on the file. The function
267 * returns the number of shifts over the ring to perform to calculate
268 * all interactions.
271 extern void start_time(void);
272 /* Start timing routines */
274 extern void update_time(void);
275 /* Update the timer.This must be done at least every INT_MAX microseconds,
276 * or 2400 s, in order to give reliable answers.
279 extern double node_time(void);
280 /* Return the node time so far in seconds. */
282 extern void do_shakefirst(FILE *log,bool bTYZ,real ener[],
283 t_parm *parm,t_nsborder *nsb,t_mdatoms *md,
284 t_state *state,rvec vold[],rvec buf[],rvec f[],
285 t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
286 t_groups *grps,t_forcerec *fr,t_topology *top,
287 t_edsamyn *edyn,t_pull *pulldata);
289 extern void dynamic_load_balancing(bool bVerbose,t_commrec *cr,real capacity[],
290 int dimension,t_mdatoms *md,t_topology *top,
291 rvec x[],rvec v[],matrix box);
292 /* Perform load balancing, i.e. split the particles over processors
293 * based on their coordinates in the "dimension" direction.
296 extern void mdrunner(t_commrec *cr,t_commrec *mcr,int nfile,t_filenm fnm[],
297 bool bVerbose,bool bCompact,
298 int nDlb,int nstepout,t_edsamyn *edyn,
299 unsigned long Flags);
300 /* Driver routine, that calls the different methods */
302 extern void init_md(t_commrec *cr,t_inputrec *ir,tensor box,real *t,real *t0,
303 real *lambda,real *lam0,
304 t_nrnb *mynrnb,bool *bTYZ,t_topology *top,
305 int nfile,t_filenm fnm[],char **traj,
306 char **xtc_traj,int *fp_ene,
307 FILE **fp_dgdl,FILE **fp_field,
308 t_mdebin **mdebin,t_groups *grps,
309 tensor force_vir,tensor pme_vir,
310 tensor shake_vir,t_mdatoms *mdatoms,rvec mu_tot,
311 bool *bNEMD,bool *bSimAnn,t_vcm **vcm,t_nsborder *nsb);
312 /* Routine in sim_util.c */
314 extern void init_em(FILE *log,const char *title,t_parm *parm,
315 real *lambda,t_nrnb *mynrnb,rvec mu_tot,
316 matrix box,rvec box_size,
317 t_forcerec *fr,t_mdatoms *mdatoms,t_topology *top,
318 t_nsborder *nsb,
319 t_commrec *cr,t_vcm **vcm,int *start,int *end);
320 /* Routine in minimize.c */
322 #endif /* _mdrun_h */