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