Enforced rotation
[alexxy/gromacs.git] / src / kernel / runner.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
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.
15
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.
20  * 
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.
27  * 
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.
30  * 
31  * For more info, check our website at http://www.gromacs.org
32  * 
33  * And Hey:
34  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <signal.h>
41 #include <stdlib.h>
42
43 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
44 /* _isnan() */
45 #include <float.h>
46 #endif
47
48 #include "typedefs.h"
49 #include "smalloc.h"
50 #include "sysstuff.h"
51 #include "statutil.h"
52 #include "mdrun.h"
53 #include "network.h"
54 #include "pull.h"
55 #include "names.h"
56 #include "disre.h"
57 #include "orires.h"
58 #include "dihre.h"
59 #include "pppm.h"
60 #include "pme.h"
61 #include "mdatoms.h"
62 #include "repl_ex.h"
63 #include "qmmm.h"
64 #include "mpelogging.h"
65 #include "domdec.h"
66 #include "partdec.h"
67 #include "coulomb.h"
68 #include "constr.h"
69 #include "mvdata.h"
70 #include "checkpoint.h"
71 #include "mtop_util.h"
72 #include "sighandler.h"
73 #include "tpxio.h"
74 #include "txtdump.h"
75 #include "pull_rotation.h"
76 #include "md_openmm.h"
77
78 #ifdef GMX_LIB_MPI
79 #include <mpi.h>
80 #endif
81 #ifdef GMX_THREADS
82 #include "tmpi.h"
83 #endif
84
85 #ifdef GMX_FAHCORE
86 #include "corewrap.h"
87 #endif
88
89 #ifdef GMX_OPENMM
90 #include "md_openmm.h"
91 #endif
92
93
94 typedef struct { 
95     gmx_integrator_t *func;
96 } gmx_intp_t;
97
98 /* The array should match the eI array in include/types/enums.h */
99 #ifdef GMX_OPENMM  /* FIXME do_md_openmm needs fixing */
100 const gmx_intp_t integrator[eiNR] = { {do_md_openmm}, {do_md_openmm}, {do_md_openmm}, {do_md_openmm}, {do_md_openmm}, {do_md_openmm}, {do_md_openmm}, {do_md_openmm}, {do_md_openmm}, {do_md_openmm}, {do_md_openmm},{do_md_openmm}};
101 #else
102 const gmx_intp_t integrator[eiNR] = { {do_md}, {do_steep}, {do_cg}, {do_md}, {do_md}, {do_nm}, {do_lbfgs}, {do_tpi}, {do_tpi}, {do_md}, {do_md},{do_md}};
103 #endif
104
105 gmx_large_int_t     deform_init_init_step_tpx;
106 matrix              deform_init_box_tpx;
107 #ifdef GMX_THREADS
108 tMPI_Thread_mutex_t deform_init_box_mutex=TMPI_THREAD_MUTEX_INITIALIZER;
109 #endif
110
111
112 #ifdef GMX_THREADS
113 struct mdrunner_arglist
114 {
115     FILE *fplog;
116     t_commrec *cr;
117     int nfile;
118     const t_filenm *fnm;
119     output_env_t oenv;
120     gmx_bool bVerbose;
121     gmx_bool bCompact;
122     int nstglobalcomm;
123     ivec ddxyz;
124     int dd_node_order;
125     real rdd;
126     real rconstr;
127     const char *dddlb_opt;
128     real dlb_scale;
129     const char *ddcsx;
130     const char *ddcsy;
131     const char *ddcsz;
132     int nstepout;
133     int resetstep;
134     int nmultisim;
135     int repl_ex_nst;
136     int repl_ex_seed;
137     real pforce;
138     real cpt_period;
139     real max_hours;
140     const char *deviceOptions;
141     unsigned long Flags;
142     int ret; /* return value */
143 };
144
145
146 /* The function used for spawning threads. Extracts the mdrunner() 
147    arguments from its one argument and calls mdrunner(), after making
148    a commrec. */
149 static void mdrunner_start_fn(void *arg)
150 {
151     struct mdrunner_arglist *mda=(struct mdrunner_arglist*)arg;
152     struct mdrunner_arglist mc=*mda; /* copy the arg list to make sure 
153                                         that it's thread-local. This doesn't
154                                         copy pointed-to items, of course,
155                                         but those are all const. */
156     t_commrec *cr;  /* we need a local version of this */
157     FILE *fplog=NULL;
158     t_filenm *fnm;
159
160     fnm = dup_tfn(mc.nfile, mc.fnm);
161
162     cr = init_par_threads(mc.cr);
163
164     if (MASTER(cr))
165     {
166         fplog=mc.fplog;
167     }
168
169     mda->ret=mdrunner(cr->nnodes, fplog, cr, mc.nfile, fnm, mc.oenv, 
170                       mc.bVerbose, mc.bCompact, mc.nstglobalcomm, 
171                       mc.ddxyz, mc.dd_node_order, mc.rdd,
172                       mc.rconstr, mc.dddlb_opt, mc.dlb_scale, 
173                       mc.ddcsx, mc.ddcsy, mc.ddcsz, mc.nstepout, mc.resetstep, 
174                       mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_seed, mc.pforce, 
175                       mc.cpt_period, mc.max_hours, mc.deviceOptions, mc.Flags);
176 }
177
178 /* called by mdrunner() to start a specific number of threads (including 
179    the main thread) for thread-parallel runs. This in turn calls mdrunner()
180    for each thread. 
181    All options besides nthreads are the same as for mdrunner(). */
182 static t_commrec *mdrunner_start_threads(int nthreads, 
183               FILE *fplog,t_commrec *cr,int nfile, 
184               const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
185               gmx_bool bCompact, int nstglobalcomm,
186               ivec ddxyz,int dd_node_order,real rdd,real rconstr,
187               const char *dddlb_opt,real dlb_scale,
188               const char *ddcsx,const char *ddcsy,const char *ddcsz,
189               int nstepout,int resetstep,int nmultisim,int repl_ex_nst,
190               int repl_ex_seed, real pforce,real cpt_period, real max_hours, 
191               const char *deviceOptions, unsigned long Flags)
192 {
193     int ret;
194     struct mdrunner_arglist *mda;
195     t_commrec *crn; /* the new commrec */
196     t_filenm *fnmn;
197
198     /* first check whether we even need to start tMPI */
199     if (nthreads<2)
200         return cr;
201
202     /* a few small, one-time, almost unavoidable memory leaks: */
203     snew(mda,1);
204     fnmn=dup_tfn(nfile, fnm);
205
206     /* fill the data structure to pass as void pointer to thread start fn */
207     mda->fplog=fplog;
208     mda->cr=cr;
209     mda->nfile=nfile;
210     mda->fnm=fnmn;
211     mda->oenv=oenv;
212     mda->bVerbose=bVerbose;
213     mda->bCompact=bCompact;
214     mda->nstglobalcomm=nstglobalcomm;
215     mda->ddxyz[XX]=ddxyz[XX];
216     mda->ddxyz[YY]=ddxyz[YY];
217     mda->ddxyz[ZZ]=ddxyz[ZZ];
218     mda->dd_node_order=dd_node_order;
219     mda->rdd=rdd;
220     mda->rconstr=rconstr;
221     mda->dddlb_opt=dddlb_opt;
222     mda->dlb_scale=dlb_scale;
223     mda->ddcsx=ddcsx;
224     mda->ddcsy=ddcsy;
225     mda->ddcsz=ddcsz;
226     mda->nstepout=nstepout;
227     mda->resetstep=resetstep;
228     mda->nmultisim=nmultisim;
229     mda->repl_ex_nst=repl_ex_nst;
230     mda->repl_ex_seed=repl_ex_seed;
231     mda->pforce=pforce;
232     mda->cpt_period=cpt_period;
233     mda->max_hours=max_hours;
234     mda->deviceOptions=deviceOptions;
235     mda->Flags=Flags;
236
237     fprintf(stderr, "Starting %d threads\n",nthreads);
238     fflush(stderr);
239     /* now spawn new threads that start mdrunner_start_fn(), while 
240        the main thread returns */
241     ret=tMPI_Init_fn(TRUE, nthreads, mdrunner_start_fn, (void*)(mda) );
242     if (ret!=TMPI_SUCCESS)
243         return NULL;
244
245     /* make a new comm_rec to reflect the new situation */
246     crn=init_par_threads(cr);
247     return crn;
248 }
249
250
251 /* get the number of threads based on how many there were requested, 
252    which algorithms we're using, and how many particles there are. */
253 static int get_nthreads(int nthreads_requested, t_inputrec *inputrec,
254                         gmx_mtop_t *mtop)
255 {
256     int nthreads,nthreads_new;
257     int min_atoms_per_thread;
258     char *env;
259
260     nthreads = nthreads_requested;
261
262     /* determine # of hardware threads. */
263     if (nthreads_requested < 1)
264     {
265         if ((env = getenv("GMX_MAX_THREADS")) != NULL)
266         {
267             nthreads = 0;
268             sscanf(env,"%d",&nthreads);
269             if (nthreads < 1)
270             {
271                 gmx_fatal(FARGS,"GMX_MAX_THREADS (%d) should be larger than 0",
272                           nthreads);
273             }
274         }
275         else
276         {
277             nthreads = tMPI_Thread_get_hw_number();
278         }
279     }
280
281     if (inputrec->eI == eiNM || EI_TPI(inputrec->eI))
282     {
283         /* Steps are divided over the nodes iso splitting the atoms */
284         min_atoms_per_thread = 0;
285     }
286     else
287     {
288         min_atoms_per_thread = MIN_ATOMS_PER_THREAD;
289     }
290
291     /* Check if an algorithm does not support parallel simulation.  */
292     if (nthreads != 1 && 
293         ( inputrec->eI == eiLBFGS ||
294           inputrec->coulombtype == eelEWALD ) )
295     {
296         fprintf(stderr,"\nThe integration or electrostatics algorithm doesn't support parallel runs. Not starting any threads.\n");
297         nthreads = 1;
298     }
299     else if (nthreads_requested < 1 &&
300              mtop->natoms/nthreads < min_atoms_per_thread)
301     {
302         /* the thread number was chosen automatically, but there are too many
303            threads (too few atoms per thread) */
304         nthreads_new = max(1,mtop->natoms/min_atoms_per_thread);
305
306         if (nthreads_new > 8 || (nthreads == 8 && nthreads_new > 4))
307         {
308             /* Use only multiples of 4 above 8 threads
309              * or with an 8-core processor
310              * (to avoid 6 threads on 8 core processors with 4 real cores).
311              */
312             nthreads_new = (nthreads_new/4)*4;
313         }
314         else if (nthreads_new > 4)
315         {
316             /* Avoid 5 or 7 threads */
317             nthreads_new = (nthreads_new/2)*2;
318         }
319
320         nthreads = nthreads_new;
321
322         fprintf(stderr,"\n");
323         fprintf(stderr,"NOTE: Parallelization is limited by the small number of atoms,\n");
324         fprintf(stderr,"      only starting %d threads.\n",nthreads);
325         fprintf(stderr,"      You can use the -nt option to optimize the number of threads.\n\n");
326     }
327     return nthreads;
328 }
329 #endif
330
331
332 int mdrunner(int nthreads_requested, FILE *fplog,t_commrec *cr,int nfile,
333              const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
334              gmx_bool bCompact, int nstglobalcomm,
335              ivec ddxyz,int dd_node_order,real rdd,real rconstr,
336              const char *dddlb_opt,real dlb_scale,
337              const char *ddcsx,const char *ddcsy,const char *ddcsz,
338              int nstepout,int resetstep,int nmultisim,int repl_ex_nst,
339              int repl_ex_seed, real pforce,real cpt_period,real max_hours,
340              const char *deviceOptions, unsigned long Flags)
341 {
342     double     nodetime=0,realtime;
343     t_inputrec *inputrec;
344     t_state    *state=NULL;
345     matrix     box;
346     gmx_ddbox_t ddbox={0};
347     int        npme_major,npme_minor;
348     real       tmpr1,tmpr2;
349     t_nrnb     *nrnb;
350     gmx_mtop_t *mtop=NULL;
351     t_mdatoms  *mdatoms=NULL;
352     t_forcerec *fr=NULL;
353     t_fcdata   *fcd=NULL;
354     real       ewaldcoeff=0;
355     gmx_pme_t  *pmedata=NULL;
356     gmx_vsite_t *vsite=NULL;
357     gmx_constr_t constr;
358     int        i,m,nChargePerturbed=-1,status,nalloc;
359     char       *gro;
360     gmx_wallcycle_t wcycle;
361     gmx_bool       bReadRNG,bReadEkin;
362     int        list;
363     gmx_runtime_t runtime;
364     int        rc;
365     gmx_large_int_t reset_counters;
366     gmx_edsam_t ed=NULL;
367     t_commrec   *cr_old=cr; 
368     int         nthreads=1;
369
370     /* CAUTION: threads may be started later on in this function, so
371        cr doesn't reflect the final parallel state right now */
372     snew(inputrec,1);
373     snew(mtop,1);
374
375     if (bVerbose && SIMMASTER(cr))
376     {
377         fprintf(stderr,"Getting Loaded...\n");
378     }
379     
380     if (Flags & MD_APPENDFILES) 
381     {
382         fplog = NULL;
383     }
384
385     snew(state,1);
386     if (MASTER(cr)) 
387     {
388         /* Read (nearly) all data required for the simulation */
389         read_tpx_state(ftp2fn(efTPX,nfile,fnm),inputrec,state,NULL,mtop);
390
391         /* NOW the threads will be started: */
392 #ifdef GMX_THREADS
393         nthreads = get_nthreads(nthreads_requested, inputrec, mtop);
394
395         if (nthreads > 1)
396         {
397             /* now start the threads. */
398             cr=mdrunner_start_threads(nthreads, fplog, cr_old, nfile, fnm, 
399                                       oenv, bVerbose, bCompact, nstglobalcomm, 
400                                       ddxyz, dd_node_order, rdd, rconstr, 
401                                       dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
402                                       nstepout, resetstep, nmultisim, 
403                                       repl_ex_nst, repl_ex_seed, pforce, 
404                                       cpt_period, max_hours, deviceOptions, 
405                                       Flags);
406             /* the main thread continues here with a new cr. We don't deallocate
407                the old cr because other threads may still be reading it. */
408             if (cr == NULL)
409             {
410                 gmx_comm("Failed to spawn threads");
411             }
412         }
413 #endif
414     }
415     /* END OF CAUTION: cr is now reliable */
416
417     if (PAR(cr))
418     {
419         /* now broadcast everything to the non-master nodes/threads: */
420         init_parallel(fplog, cr, inputrec, mtop);
421     }
422     if (fplog != NULL)
423     {
424         pr_inputrec(fplog,0,"Input Parameters",inputrec,FALSE);
425     }
426
427     /* now make sure the state is initialized and propagated */
428     set_state_entries(state,inputrec,cr->nnodes);
429
430     /* A parallel command line option consistency check that we can
431        only do after any threads have started. */
432     if (!PAR(cr) &&
433         (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
434     {
435         gmx_fatal(FARGS,
436                   "The -dd or -npme option request a parallel simulation, "
437 #ifndef GMX_MPI
438                   "but mdrun was compiled without threads or MPI enabled"
439 #else
440 #ifdef GMX_THREADS
441                   "but the number of threads (option -nt) is 1"
442 #else
443                   "but mdrun was not started through mpirun/mpiexec or only one process was requested through mpirun/mpiexec" 
444 #endif
445 #endif
446             );
447     }
448
449     if ((Flags & MD_RERUN) &&
450         (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
451     {
452         gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
453     }
454
455     if (can_use_allvsall(inputrec,mtop,TRUE,cr,fplog))
456     {
457         /* All-vs-all loops do not work with domain decomposition */
458         Flags |= MD_PARTDEC;
459     }
460
461     if (!EEL_PME(inputrec->coulombtype) || (Flags & MD_PARTDEC))
462     {
463         cr->npmenodes = 0;
464     }
465
466 #ifdef GMX_FAHCORE
467     fcRegisterSteps(inputrec->nsteps,inputrec->init_step);
468 #endif
469
470     /* NMR restraints must be initialized before load_checkpoint,
471      * since with time averaging the history is added to t_state.
472      * For proper consistency check we therefore need to extend
473      * t_state here.
474      * So the PME-only nodes (if present) will also initialize
475      * the distance restraints.
476      */
477     snew(fcd,1);
478
479     /* This needs to be called before read_checkpoint to extend the state */
480     init_disres(fplog,mtop,inputrec,cr,Flags & MD_PARTDEC,fcd,state);
481
482     if (gmx_mtop_ftype_count(mtop,F_ORIRES) > 0)
483     {
484         if (PAR(cr) && !(Flags & MD_PARTDEC))
485         {
486             gmx_fatal(FARGS,"Orientation restraints do not work (yet) with domain decomposition, use particle decomposition (mdrun option -pd)");
487         }
488         /* Orientation restraints */
489         if (MASTER(cr))
490         {
491             init_orires(fplog,mtop,state->x,inputrec,cr->ms,&(fcd->orires),
492                         state);
493         }
494     }
495
496     if (DEFORM(*inputrec))
497     {
498         /* Store the deform reference box before reading the checkpoint */
499         if (SIMMASTER(cr))
500         {
501             copy_mat(state->box,box);
502         }
503         if (PAR(cr))
504         {
505             gmx_bcast(sizeof(box),box,cr);
506         }
507         /* Because we do not have the update struct available yet
508          * in which the reference values should be stored,
509          * we store them temporarily in static variables.
510          * This should be thread safe, since they are only written once
511          * and with identical values.
512          */
513 #ifdef GMX_THREADS
514         tMPI_Thread_mutex_lock(&deform_init_box_mutex);
515 #endif
516         deform_init_init_step_tpx = inputrec->init_step;
517         copy_mat(box,deform_init_box_tpx);
518 #ifdef GMX_THREADS
519         tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
520 #endif
521     }
522
523     if (opt2bSet("-cpi",nfile,fnm)) 
524     {
525         /* Check if checkpoint file exists before doing continuation.
526          * This way we can use identical input options for the first and subsequent runs...
527          */
528         if( gmx_fexist_master(opt2fn_master("-cpi",nfile,fnm,cr),cr) )
529         {
530             load_checkpoint(opt2fn_master("-cpi",nfile,fnm,cr),&fplog,
531                             cr,Flags & MD_PARTDEC,ddxyz,
532                             inputrec,state,&bReadRNG,&bReadEkin,
533                             (Flags & MD_APPENDFILES));
534             
535             if (bReadRNG)
536             {
537                 Flags |= MD_READ_RNG;
538             }
539             if (bReadEkin)
540             {
541                 Flags |= MD_READ_EKIN;
542             }
543         }
544     }
545
546     if (((MASTER(cr) || (Flags & MD_SEPPOT)) && (Flags & MD_APPENDFILES))
547 #ifdef GMX_THREADS
548         /* With thread MPI only the master node/thread exists in mdrun.c,
549          * therefore non-master nodes need to open the "seppot" log file here.
550          */
551         || (!MASTER(cr) && (Flags & MD_SEPPOT))
552 #endif
553         )
554     {
555         gmx_log_open(ftp2fn(efLOG,nfile,fnm),cr,!(Flags & MD_SEPPOT),
556                              Flags,&fplog);
557     }
558
559     if (SIMMASTER(cr)) 
560     {
561         copy_mat(state->box,box);
562     }
563
564     if (PAR(cr)) 
565     {
566         gmx_bcast(sizeof(box),box,cr);
567     }
568
569     /* Essential dynamics */
570     if (opt2bSet("-ei",nfile,fnm))
571     {
572         /* Open input and output files, allocate space for ED data structure */
573         ed = ed_open(nfile,fnm,Flags,cr);
574     }
575
576     if (bVerbose && SIMMASTER(cr))
577     {
578         fprintf(stderr,"Loaded with Money\n\n");
579     }
580
581     if (PAR(cr) && !((Flags & MD_PARTDEC) ||
582                      EI_TPI(inputrec->eI) ||
583                      inputrec->eI == eiNM))
584     {
585         cr->dd = init_domain_decomposition(fplog,cr,Flags,ddxyz,rdd,rconstr,
586                                            dddlb_opt,dlb_scale,
587                                            ddcsx,ddcsy,ddcsz,
588                                            mtop,inputrec,
589                                            box,state->x,
590                                            &ddbox,&npme_major,&npme_minor);
591
592         make_dd_communicators(fplog,cr,dd_node_order);
593
594         /* Set overallocation to avoid frequent reallocation of arrays */
595         set_over_alloc_dd(TRUE);
596     }
597     else
598     {
599         /* PME, if used, is done on all nodes with 1D decomposition */
600         cr->npmenodes = 0;
601         cr->duty = (DUTY_PP | DUTY_PME);
602         npme_major = 1;
603         npme_minor = 1;
604         if (!EI_TPI(inputrec->eI))
605         {
606             npme_major = cr->nnodes;
607         }
608         
609         if (inputrec->ePBC == epbcSCREW)
610         {
611             gmx_fatal(FARGS,
612                       "pbc=%s is only implemented with domain decomposition",
613                       epbc_names[inputrec->ePBC]);
614         }
615     }
616
617     if (PAR(cr))
618     {
619         /* After possible communicator splitting in make_dd_communicators.
620          * we can set up the intra/inter node communication.
621          */
622         gmx_setup_nodecomm(fplog,cr);
623     }
624
625     wcycle = wallcycle_init(fplog,resetstep,cr);
626     if (PAR(cr))
627     {
628         /* Master synchronizes its value of reset_counters with all nodes 
629          * including PME only nodes */
630         reset_counters = wcycle_get_reset_counters(wcycle);
631         gmx_bcast_sim(sizeof(reset_counters),&reset_counters,cr);
632         wcycle_set_reset_counters(wcycle, reset_counters);
633     }
634
635
636     snew(nrnb,1);
637     if (cr->duty & DUTY_PP)
638     {
639         /* For domain decomposition we allocate dynamically
640          * in dd_partition_system.
641          */
642         if (DOMAINDECOMP(cr))
643         {
644             bcast_state_setup(cr,state);
645         }
646         else
647         {
648             if (PAR(cr))
649             {
650                 bcast_state(cr,state,TRUE);
651             }
652         }
653
654         /* Dihedral Restraints */
655         if (gmx_mtop_ftype_count(mtop,F_DIHRES) > 0)
656         {
657             init_dihres(fplog,mtop,inputrec,fcd);
658         }
659
660         /* Initiate forcerecord */
661         fr = mk_forcerec();
662         init_forcerec(fplog,oenv,fr,fcd,inputrec,mtop,cr,box,FALSE,
663                       opt2fn("-table",nfile,fnm),
664                       opt2fn("-tablep",nfile,fnm),
665                       opt2fn("-tableb",nfile,fnm),FALSE,pforce);
666
667         /* version for PCA_NOT_READ_NODE (see md.c) */
668         /*init_forcerec(fplog,fr,fcd,inputrec,mtop,cr,box,FALSE,
669           "nofile","nofile","nofile",FALSE,pforce);
670           */        
671         fr->bSepDVDL = ((Flags & MD_SEPPOT) == MD_SEPPOT);
672
673         /* Initialize QM-MM */
674         if(fr->bQMMM)
675         {
676             init_QMMMrec(cr,box,mtop,inputrec,fr);
677         }
678
679         /* Initialize the mdatoms structure.
680          * mdatoms is not filled with atom data,
681          * as this can not be done now with domain decomposition.
682          */
683         mdatoms = init_mdatoms(fplog,mtop,inputrec->efep!=efepNO);
684
685         /* Initialize the virtual site communication */
686         vsite = init_vsite(mtop,cr);
687
688         calc_shifts(box,fr->shift_vec);
689
690         /* With periodic molecules the charge groups should be whole at start up
691          * and the virtual sites should not be far from their proper positions.
692          */
693         if (!inputrec->bContinuation && MASTER(cr) &&
694             !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
695         {
696             /* Make molecules whole at start of run */
697             if (fr->ePBC != epbcNONE)
698             {
699                 do_pbc_first_mtop(fplog,inputrec->ePBC,box,mtop,state->x);
700             }
701             if (vsite)
702             {
703                 /* Correct initial vsite positions are required
704                  * for the initial distribution in the domain decomposition
705                  * and for the initial shell prediction.
706                  */
707                 construct_vsites_mtop(fplog,vsite,mtop,state->x);
708             }
709         }
710
711         /* Initiate PPPM if necessary */
712         if (fr->eeltype == eelPPPM)
713         {
714             if (mdatoms->nChargePerturbed)
715             {
716                 gmx_fatal(FARGS,"Free energy with %s is not implemented",
717                           eel_names[fr->eeltype]);
718             }
719             status = gmx_pppm_init(fplog,cr,oenv,FALSE,TRUE,box,
720                                    getenv("GMXGHAT"),inputrec, (Flags & MD_REPRODUCIBLE));
721             if (status != 0)
722             {
723                 gmx_fatal(FARGS,"Error %d initializing PPPM",status);
724             }
725         }
726
727         if (EEL_PME(fr->eeltype))
728         {
729             ewaldcoeff = fr->ewaldcoeff;
730             pmedata = &fr->pmedata;
731         }
732         else
733         {
734             pmedata = NULL;
735         }
736     }
737     else
738     {
739         /* This is a PME only node */
740
741         /* We don't need the state */
742         done_state(state);
743
744         ewaldcoeff = calc_ewaldcoeff(inputrec->rcoulomb, inputrec->ewald_rtol);
745         snew(pmedata,1);
746     }
747
748     /* Initiate PME if necessary,
749      * either on all nodes or on dedicated PME nodes only. */
750     if (EEL_PME(inputrec->coulombtype))
751     {
752         if (mdatoms)
753         {
754             nChargePerturbed = mdatoms->nChargePerturbed;
755         }
756         if (cr->npmenodes > 0)
757         {
758             /* The PME only nodes need to know nChargePerturbed */
759             gmx_bcast_sim(sizeof(nChargePerturbed),&nChargePerturbed,cr);
760         }
761         if (cr->duty & DUTY_PME)
762         {
763             status = gmx_pme_init(pmedata,cr,npme_major,npme_minor,inputrec,
764                                   mtop ? mtop->natoms : 0,nChargePerturbed,
765                                   (Flags & MD_REPRODUCIBLE));
766             if (status != 0) 
767             {
768                 gmx_fatal(FARGS,"Error %d initializing PME",status);
769             }
770         }
771     }
772
773
774     if (integrator[inputrec->eI].func == do_md
775 #ifdef GMX_OPENMM
776         ||
777         integrator[inputrec->eI].func == do_md_openmm
778 #endif
779         )
780     {
781         /* Turn on signal handling on all nodes */
782         /*
783          * (A user signal from the PME nodes (if any)
784          * is communicated to the PP nodes.
785          */
786         signal_handler_install();
787     }
788
789     if (cr->duty & DUTY_PP)
790     {
791         if (inputrec->ePull != epullNO)
792         {
793             /* Initialize pull code */
794             init_pull(fplog,inputrec,nfile,fnm,mtop,cr,oenv,
795                       EI_DYNAMICS(inputrec->eI) && MASTER(cr),Flags);
796         }
797         
798         if (inputrec->bRot)
799         {
800            /* Initialize enforced rotation code */
801            init_rot(fplog,inputrec,nfile,fnm,cr,state->x,state->box,mtop,oenv,
802                     bVerbose,Flags);
803         }
804
805         constr = init_constraints(fplog,mtop,inputrec,ed,state,cr);
806
807         if (DOMAINDECOMP(cr))
808         {
809             dd_init_bondeds(fplog,cr->dd,mtop,vsite,constr,inputrec,
810                             Flags & MD_DDBONDCHECK,fr->cginfo_mb);
811
812             set_dd_parameters(fplog,cr->dd,dlb_scale,inputrec,fr,&ddbox);
813
814             setup_dd_grid(fplog,cr->dd);
815         }
816
817         /* Now do whatever the user wants us to do (how flexible...) */
818         integrator[inputrec->eI].func(fplog,cr,nfile,fnm,
819                                       oenv,bVerbose,bCompact,
820                                       nstglobalcomm,
821                                       vsite,constr,
822                                       nstepout,inputrec,mtop,
823                                       fcd,state,
824                                       mdatoms,nrnb,wcycle,ed,fr,
825                                       repl_ex_nst,repl_ex_seed,
826                                       cpt_period,max_hours,
827                                       deviceOptions,
828                                       Flags,
829                                       &runtime);
830
831         if (inputrec->ePull != epullNO)
832         {
833             finish_pull(fplog,inputrec->pull);
834         }
835         
836         if (inputrec->bRot)
837         {
838             finish_rot(fplog,inputrec->rot);
839         }
840
841     } 
842     else 
843     {
844         /* do PME only */
845         gmx_pmeonly(*pmedata,cr,nrnb,wcycle,ewaldcoeff,FALSE,inputrec);
846     }
847
848     if (EI_DYNAMICS(inputrec->eI) || EI_TPI(inputrec->eI))
849     {
850         /* Some timing stats */  
851         if (SIMMASTER(cr))
852         {
853             if (runtime.proc == 0)
854             {
855                 runtime.proc = runtime.real;
856             }
857         }
858         else
859         {
860             runtime.real = 0;
861         }
862     }
863
864     wallcycle_stop(wcycle,ewcRUN);
865
866     /* Finish up, write some stuff
867      * if rerunMD, don't write last frame again 
868      */
869     finish_run(fplog,cr,ftp2fn(efSTO,nfile,fnm),
870                inputrec,nrnb,wcycle,&runtime,
871                EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
872
873     /* Does what it says */  
874     print_date_and_time(fplog,cr->nodeid,"Finished mdrun",&runtime);
875
876     /* Close logfile already here if we were appending to it */
877     if (MASTER(cr) && (Flags & MD_APPENDFILES))
878     {
879         gmx_log_close(fplog);
880     }   
881
882     rc=(int)gmx_get_stop_condition();
883
884 #ifdef GMX_THREADS
885     /* we need to join all threads. The sub-threads join when they
886        exit this function, but the master thread needs to be told to 
887        wait for that. */
888     if (PAR(cr) && MASTER(cr))
889     {
890         tMPI_Finalize();
891     }
892 #endif
893
894     return rc;
895 }