d8f8a383f57a41053695bff1c9ae5f43ed3bbbfb
[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 #if defined(HAVE_SCHED_H) && (defined(HAVE_SCHED_GETAFFINITY) || defined(HAVE_SCHED_SETAFFINITY))
40 #define _GNU_SOURCE
41 #include <sched.h>
42 #include <sys/syscall.h>
43 #endif
44 #include <signal.h>
45 #include <stdlib.h>
46 #include <string.h>
47 #include <assert.h>
48
49 #include "typedefs.h"
50 #include "smalloc.h"
51 #include "sysstuff.h"
52 #include "statutil.h"
53 #include "mdrun.h"
54 #include "md_logging.h"
55 #include "md_support.h"
56 #include "network.h"
57 #include "pull.h"
58 #include "names.h"
59 #include "disre.h"
60 #include "orires.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 "gmx_detect_hardware.h"
77 #include "gmx_omp_nthreads.h"
78 #include "pull_rotation.h"
79 #include "calc_verletbuf.h"
80 #include "../mdlib/nbnxn_search.h"
81 #include "../mdlib/nbnxn_consts.h"
82 #include "gmx_fatal_collective.h"
83 #include "membed.h"
84 #include "md_openmm.h"
85 #include "gmx_omp.h"
86
87 #ifdef GMX_LIB_MPI
88 #include <mpi.h>
89 #endif
90 #ifdef GMX_THREAD_MPI
91 #include "tmpi.h"
92 #endif
93
94 #ifdef GMX_FAHCORE
95 #include "corewrap.h"
96 #endif
97
98 #ifdef GMX_OPENMM
99 #include "md_openmm.h"
100 #endif
101
102 #include "gpu_utils.h"
103 #include "nbnxn_cuda_data_mgmt.h"
104
105 typedef struct { 
106     gmx_integrator_t *func;
107 } gmx_intp_t;
108
109 /* The array should match the eI array in include/types/enums.h */
110 #ifdef GMX_OPENMM  /* FIXME do_md_openmm needs fixing */
111 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}};
112 #else
113 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}};
114 #endif
115
116 gmx_large_int_t     deform_init_init_step_tpx;
117 matrix              deform_init_box_tpx;
118 #ifdef GMX_THREAD_MPI
119 tMPI_Thread_mutex_t deform_init_box_mutex=TMPI_THREAD_MUTEX_INITIALIZER;
120 #endif
121
122
123 #ifdef GMX_THREAD_MPI
124 struct mdrunner_arglist
125 {
126     gmx_hw_opt_t *hw_opt;
127     FILE *fplog;
128     t_commrec *cr;
129     int nfile;
130     const t_filenm *fnm;
131     output_env_t oenv;
132     gmx_bool bVerbose;
133     gmx_bool bCompact;
134     int nstglobalcomm;
135     ivec ddxyz;
136     int dd_node_order;
137     real rdd;
138     real rconstr;
139     const char *dddlb_opt;
140     real dlb_scale;
141     const char *ddcsx;
142     const char *ddcsy;
143     const char *ddcsz;
144     const char *nbpu_opt;
145     int nsteps_cmdline;
146     int nstepout;
147     int resetstep;
148     int nmultisim;
149     int repl_ex_nst;
150     int repl_ex_nex;
151     int repl_ex_seed;
152     real pforce;
153     real cpt_period;
154     real max_hours;
155     const char *deviceOptions;
156     unsigned long Flags;
157     int ret; /* return value */
158 };
159
160
161 /* The function used for spawning threads. Extracts the mdrunner() 
162    arguments from its one argument and calls mdrunner(), after making
163    a commrec. */
164 static void mdrunner_start_fn(void *arg)
165 {
166     struct mdrunner_arglist *mda=(struct mdrunner_arglist*)arg;
167     struct mdrunner_arglist mc=*mda; /* copy the arg list to make sure 
168                                         that it's thread-local. This doesn't
169                                         copy pointed-to items, of course,
170                                         but those are all const. */
171     t_commrec *cr;  /* we need a local version of this */
172     FILE *fplog=NULL;
173     t_filenm *fnm;
174
175     fnm = dup_tfn(mc.nfile, mc.fnm);
176
177     cr = init_par_threads(mc.cr);
178
179     if (MASTER(cr))
180     {
181         fplog=mc.fplog;
182     }
183
184     mda->ret=mdrunner(mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv, 
185                       mc.bVerbose, mc.bCompact, mc.nstglobalcomm, 
186                       mc.ddxyz, mc.dd_node_order, mc.rdd,
187                       mc.rconstr, mc.dddlb_opt, mc.dlb_scale, 
188                       mc.ddcsx, mc.ddcsy, mc.ddcsz,
189                       mc.nbpu_opt,
190                       mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
191                       mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce, 
192                       mc.cpt_period, mc.max_hours, mc.deviceOptions, mc.Flags);
193 }
194
195 /* called by mdrunner() to start a specific number of threads (including 
196    the main thread) for thread-parallel runs. This in turn calls mdrunner()
197    for each thread. 
198    All options besides nthreads are the same as for mdrunner(). */
199 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt, 
200               FILE *fplog,t_commrec *cr,int nfile, 
201               const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
202               gmx_bool bCompact, int nstglobalcomm,
203               ivec ddxyz,int dd_node_order,real rdd,real rconstr,
204               const char *dddlb_opt,real dlb_scale,
205               const char *ddcsx,const char *ddcsy,const char *ddcsz,
206               const char *nbpu_opt,
207               int nsteps_cmdline, int nstepout,int resetstep,
208               int nmultisim,int repl_ex_nst,int repl_ex_nex, int repl_ex_seed,
209               real pforce,real cpt_period, real max_hours, 
210               const char *deviceOptions, unsigned long Flags)
211 {
212     int ret;
213     struct mdrunner_arglist *mda;
214     t_commrec *crn; /* the new commrec */
215     t_filenm *fnmn;
216
217     /* first check whether we even need to start tMPI */
218     if (hw_opt->nthreads_tmpi < 2)
219     {
220         return cr;
221     }
222
223     /* a few small, one-time, almost unavoidable memory leaks: */
224     snew(mda,1);
225     fnmn=dup_tfn(nfile, fnm);
226
227     /* fill the data structure to pass as void pointer to thread start fn */
228     mda->hw_opt=hw_opt;
229     mda->fplog=fplog;
230     mda->cr=cr;
231     mda->nfile=nfile;
232     mda->fnm=fnmn;
233     mda->oenv=oenv;
234     mda->bVerbose=bVerbose;
235     mda->bCompact=bCompact;
236     mda->nstglobalcomm=nstglobalcomm;
237     mda->ddxyz[XX]=ddxyz[XX];
238     mda->ddxyz[YY]=ddxyz[YY];
239     mda->ddxyz[ZZ]=ddxyz[ZZ];
240     mda->dd_node_order=dd_node_order;
241     mda->rdd=rdd;
242     mda->rconstr=rconstr;
243     mda->dddlb_opt=dddlb_opt;
244     mda->dlb_scale=dlb_scale;
245     mda->ddcsx=ddcsx;
246     mda->ddcsy=ddcsy;
247     mda->ddcsz=ddcsz;
248     mda->nbpu_opt=nbpu_opt;
249     mda->nsteps_cmdline=nsteps_cmdline;
250     mda->nstepout=nstepout;
251     mda->resetstep=resetstep;
252     mda->nmultisim=nmultisim;
253     mda->repl_ex_nst=repl_ex_nst;
254     mda->repl_ex_nex=repl_ex_nex;
255     mda->repl_ex_seed=repl_ex_seed;
256     mda->pforce=pforce;
257     mda->cpt_period=cpt_period;
258     mda->max_hours=max_hours;
259     mda->deviceOptions=deviceOptions;
260     mda->Flags=Flags;
261
262     fprintf(stderr, "Starting %d tMPI threads\n",hw_opt->nthreads_tmpi);
263     fflush(stderr);
264     /* now spawn new threads that start mdrunner_start_fn(), while 
265        the main thread returns */
266     ret=tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi,
267                      (hw_opt->bThreadPinning ? TMPI_AFFINITY_ALL_CORES : TMPI_AFFINITY_NONE),
268                      mdrunner_start_fn, (void*)(mda) );
269     if (ret!=TMPI_SUCCESS)
270         return NULL;
271
272     /* make a new comm_rec to reflect the new situation */
273     crn=init_par_threads(cr);
274     return crn;
275 }
276
277
278 static int get_tmpi_omp_thread_division(const gmx_hw_info_t *hwinfo,
279                                         const gmx_hw_opt_t *hw_opt,
280                                         int nthreads_tot,
281                                         int ngpu)
282 {
283     int nthreads_tmpi;
284
285     /* There are no separate PME nodes here, as we ensured in
286      * check_and_update_hw_opt that nthreads_tmpi>0 with PME nodes
287      * and a conditional ensures we would not have ended up here.
288      * Note that separate PME nodes might be switched on later.
289      */
290     if (ngpu > 0)
291     {
292         nthreads_tmpi = ngpu;
293         if (nthreads_tot > 0 && nthreads_tot < nthreads_tmpi)
294         {
295             nthreads_tmpi = nthreads_tot;
296         }
297     }
298     else if (hw_opt->nthreads_omp > 0)
299     {
300         /* Here we could oversubscribe, when we do, we issue a warning later */
301         nthreads_tmpi = max(1,nthreads_tot/hw_opt->nthreads_omp);
302     }
303     else
304     {
305         /* TODO choose nthreads_omp based on hardware topology
306            when we have a hardware topology detection library */
307         /* In general, when running up to 4 threads, OpenMP should be faster.
308          * Note: on AMD Bulldozer we should avoid running OpenMP over two dies.
309          * On Intel>=Nehalem running OpenMP on a single CPU is always faster,
310          * even on two CPUs it's usually faster (but with many OpenMP threads
311          * it could be faster not to use HT, currently we always use HT).
312          * On Nehalem/Westmere we want to avoid running 16 threads over
313          * two CPUs with HT, so we need a limit<16; thus we use 12.
314          * A reasonable limit for Intel Sandy and Ivy bridge,
315          * not knowing the topology, is 16 threads.
316          */
317         const int nthreads_omp_always_faster             =  4;
318         const int nthreads_omp_always_faster_Nehalem     = 12;
319         const int nthreads_omp_always_faster_SandyBridge = 16;
320         const int first_model_Nehalem     = 0x1A;
321         const int first_model_SandyBridge = 0x2A;
322         gmx_bool bIntel_Family6;
323
324         bIntel_Family6 =
325             (gmx_cpuid_vendor(hwinfo->cpuid_info) == GMX_CPUID_VENDOR_INTEL &&
326              gmx_cpuid_family(hwinfo->cpuid_info) == 6);
327
328         if (nthreads_tot <= nthreads_omp_always_faster ||
329             (bIntel_Family6 &&
330              ((gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_Nehalem && nthreads_tot <= nthreads_omp_always_faster_Nehalem) ||
331               (gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_SandyBridge && nthreads_tot <= nthreads_omp_always_faster_SandyBridge))))
332         {
333             /* Use pure OpenMP parallelization */
334             nthreads_tmpi = 1;
335         }
336         else
337         {
338             /* Don't use OpenMP parallelization */
339             nthreads_tmpi = nthreads_tot;
340         }
341     }
342
343     return nthreads_tmpi;
344 }
345
346
347 /* Get the number of threads to use for thread-MPI based on how many
348  * were requested, which algorithms we're using,
349  * and how many particles there are.
350  * At the point we have already called check_and_update_hw_opt.
351  * Thus all options should be internally consistent and consistent
352  * with the hardware, except that ntmpi could be larger than #GPU.
353  */
354 static int get_nthreads_mpi(gmx_hw_info_t *hwinfo,
355                             gmx_hw_opt_t *hw_opt,
356                             t_inputrec *inputrec, gmx_mtop_t *mtop,
357                             const t_commrec *cr,
358                             FILE *fplog)
359 {
360     int nthreads_hw,nthreads_tot_max,nthreads_tmpi,nthreads_new,ngpu;
361     int min_atoms_per_mpi_thread;
362     char *env;
363     char sbuf[STRLEN];
364     gmx_bool bCanUseGPU;
365
366     if (hw_opt->nthreads_tmpi > 0)
367     {
368         /* Trivial, return right away */
369         return hw_opt->nthreads_tmpi;
370     }
371
372     nthreads_hw = hwinfo->nthreads_hw_avail;
373
374     /* How many total (#tMPI*#OpenMP) threads can we start? */ 
375     if (hw_opt->nthreads_tot > 0)
376     {
377         nthreads_tot_max = hw_opt->nthreads_tot;
378     }
379     else
380     {
381         nthreads_tot_max = nthreads_hw;
382     }
383
384     bCanUseGPU = (inputrec->cutoff_scheme == ecutsVERLET && hwinfo->bCanUseGPU);
385     if (bCanUseGPU)
386     {
387         ngpu = hwinfo->gpu_info.ncuda_dev_use;
388     }
389     else
390     {
391         ngpu = 0;
392     }
393
394     nthreads_tmpi =
395         get_tmpi_omp_thread_division(hwinfo,hw_opt,nthreads_tot_max,ngpu);
396
397     if (inputrec->eI == eiNM || EI_TPI(inputrec->eI))
398     {
399         /* Steps are divided over the nodes iso splitting the atoms */
400         min_atoms_per_mpi_thread = 0;
401     }
402     else
403     {
404         if (bCanUseGPU)
405         {
406             min_atoms_per_mpi_thread = MIN_ATOMS_PER_GPU;
407         }
408         else
409         {
410             min_atoms_per_mpi_thread = MIN_ATOMS_PER_MPI_THREAD;
411         }
412     }
413
414     /* Check if an algorithm does not support parallel simulation.  */
415     if (nthreads_tmpi != 1 &&
416         ( inputrec->eI == eiLBFGS ||
417           inputrec->coulombtype == eelEWALD ) )
418     {
419         nthreads_tmpi = 1;
420
421         md_print_warn(cr,fplog,"The integration or electrostatics algorithm doesn't support parallel runs. Using a single thread-MPI thread.\n");
422         if (hw_opt->nthreads_tmpi > nthreads_tmpi)
423         {
424             gmx_fatal(FARGS,"You asked for more than 1 thread-MPI thread, but an algorithm doesn't support that");
425         }
426     }
427     else if (mtop->natoms/nthreads_tmpi < min_atoms_per_mpi_thread)
428     {
429         /* the thread number was chosen automatically, but there are too many
430            threads (too few atoms per thread) */
431         nthreads_new = max(1,mtop->natoms/min_atoms_per_mpi_thread);
432
433         /* Avoid partial use of Hyper-Threading */
434         if (gmx_cpuid_x86_smt(hwinfo->cpuid_info) == GMX_CPUID_X86_SMT_ENABLED &&
435             nthreads_new > nthreads_hw/2 && nthreads_new < nthreads_hw)
436         {
437             nthreads_new = nthreads_hw/2;
438         }
439
440         /* Avoid large prime numbers in the thread count */
441         if (nthreads_new >= 6)
442         {
443             /* Use only 6,8,10 with additional factors of 2 */
444             int fac;
445
446             fac = 2;
447             while (3*fac*2 <= nthreads_new)
448             {
449                 fac *= 2;
450             }
451
452             nthreads_new = (nthreads_new/fac)*fac;
453         }
454         else
455         {
456             /* Avoid 5 */
457             if (nthreads_new == 5)
458             {
459                 nthreads_new = 4;
460             }
461         }
462
463         nthreads_tmpi = nthreads_new;
464
465         fprintf(stderr,"\n");
466         fprintf(stderr,"NOTE: Parallelization is limited by the small number of atoms,\n");
467         fprintf(stderr,"      only starting %d thread-MPI threads.\n",nthreads_tmpi);
468         fprintf(stderr,"      You can use the -nt and/or -ntmpi option to optimize the number of threads.\n\n");
469     }
470
471     return nthreads_tmpi;
472 }
473 #endif /* GMX_THREAD_MPI */
474
475
476 /* Environment variable for setting nstlist */
477 static const char*  NSTLIST_ENVVAR          =  "GMX_NSTLIST";
478 /* Try to increase nstlist when using a GPU with nstlist less than this */
479 static const int    NSTLIST_GPU_ENOUGH      = 20;
480 /* Increase nstlist until the non-bonded cost increases more than this factor */
481 static const float  NBNXN_GPU_LIST_OK_FAC   = 1.25;
482 /* Don't increase nstlist beyond a non-bonded cost increases of this factor */
483 static const float  NBNXN_GPU_LIST_MAX_FAC  = 1.40;
484
485 /* Try to increase nstlist when running on a GPU */
486 static void increase_nstlist(FILE *fp,t_commrec *cr,
487                              t_inputrec *ir,const gmx_mtop_t *mtop,matrix box)
488 {
489     char *env;
490     int  nstlist_orig,nstlist_prev;
491     verletbuf_list_setup_t ls;
492     real rlist_inc,rlist_ok,rlist_max,rlist_new,rlist_prev;
493     int  i;
494     t_state state_tmp;
495     gmx_bool bBox,bDD,bCont;
496     const char *nstl_fmt="\nFor optimal performance with a GPU nstlist (now %d) should be larger.\nThe optimum depends on your CPU and GPU resources.\nYou might want to try several nstlist values.\n";
497     const char *vbd_err="Can not increase nstlist for GPU run because verlet-buffer-drift is not set or used";
498     const char *box_err="Can not increase nstlist for GPU run because the box is too small";
499     const char *dd_err ="Can not increase nstlist for GPU run because of domain decomposition limitations";
500     char buf[STRLEN];
501
502     /* Number of + nstlist alternative values to try when switching  */
503     const int nstl[]={ 20, 25, 40, 50 };
504 #define NNSTL  sizeof(nstl)/sizeof(nstl[0])
505
506     env = getenv(NSTLIST_ENVVAR);
507     if (env == NULL)
508     {
509         if (fp != NULL)
510         {
511             fprintf(fp,nstl_fmt,ir->nstlist);
512         }
513     }
514
515     if (ir->verletbuf_drift == 0)
516     {
517         gmx_fatal(FARGS,"You are using an old tpr file with a GPU, please generate a new tpr file with an up to date version of grompp");
518     }
519
520     if (ir->verletbuf_drift < 0)
521     {
522         if (MASTER(cr))
523         {
524             fprintf(stderr,"%s\n",vbd_err);
525         }
526         if (fp != NULL)
527         {
528             fprintf(fp,"%s\n",vbd_err);
529         }
530
531         return;
532     }
533
534     nstlist_orig = ir->nstlist;
535     if (env != NULL)
536     {
537         sprintf(buf,"Getting nstlist from environment variable GMX_NSTLIST=%s",env);
538         if (MASTER(cr))
539         {
540             fprintf(stderr,"%s\n",buf);
541         }
542         if (fp != NULL)
543         {
544             fprintf(fp,"%s\n",buf);
545         }
546         sscanf(env,"%d",&ir->nstlist);
547     }
548
549     verletbuf_get_list_setup(TRUE,&ls);
550
551     /* Allow rlist to make the list double the size of the cut-off sphere */
552     rlist_inc = nbnxn_get_rlist_effective_inc(NBNXN_GPU_CLUSTER_SIZE,mtop->natoms/det(box));
553     rlist_ok  = (max(ir->rvdw,ir->rcoulomb) + rlist_inc)*pow(NBNXN_GPU_LIST_OK_FAC,1.0/3.0) - rlist_inc;
554     rlist_max = (max(ir->rvdw,ir->rcoulomb) + rlist_inc)*pow(NBNXN_GPU_LIST_MAX_FAC,1.0/3.0) - rlist_inc;
555     if (debug)
556     {
557         fprintf(debug,"GPU nstlist tuning: rlist_inc %.3f rlist_max %.3f\n",
558                 rlist_inc,rlist_max);
559     }
560
561     i = 0;
562     nstlist_prev = nstlist_orig;
563     rlist_prev   = ir->rlist;
564     do
565     {
566         if (env == NULL)
567         {
568             ir->nstlist = nstl[i];
569         }
570
571         /* Set the pair-list buffer size in ir */
572         calc_verlet_buffer_size(mtop,det(box),ir,ir->verletbuf_drift,&ls,
573                                 NULL,&rlist_new);
574
575         /* Does rlist fit in the box? */
576         bBox = (sqr(rlist_new) < max_cutoff2(ir->ePBC,box));
577         bDD  = TRUE;
578         if (bBox && DOMAINDECOMP(cr))
579         {
580             /* Check if rlist fits in the domain decomposition */
581             if (inputrec2nboundeddim(ir) < DIM)
582             {
583                 gmx_incons("Changing nstlist with domain decomposition and unbounded dimensions is not implemented yet");
584             }
585             copy_mat(box,state_tmp.box);
586             bDD = change_dd_cutoff(cr,&state_tmp,ir,rlist_new);
587         }
588
589         bCont = FALSE;
590
591         if (env == NULL)
592         {
593             if (bBox && bDD && rlist_new <= rlist_max)
594             {
595                 /* Increase nstlist */
596                 nstlist_prev = ir->nstlist;
597                 rlist_prev   = rlist_new;
598                 bCont = (i+1 < NNSTL && rlist_new < rlist_ok);
599             }
600             else
601             {
602                 /* Stick with the previous nstlist */
603                 ir->nstlist = nstlist_prev;
604                 rlist_new   = rlist_prev;
605                 bBox = TRUE;
606                 bDD  = TRUE;
607             }
608         }
609
610         i++;
611     }
612     while (bCont);
613
614     if (!bBox || !bDD)
615     {
616         gmx_warning(!bBox ? box_err : dd_err);
617         if (fp != NULL)
618         {
619             fprintf(fp,"\n%s\n",bBox ? box_err : dd_err);
620         }
621         ir->nstlist = nstlist_orig;
622     }
623     else if (ir->nstlist != nstlist_orig || rlist_new != ir->rlist)
624     {
625         sprintf(buf,"Changing nstlist from %d to %d, rlist from %g to %g",
626                 nstlist_orig,ir->nstlist,
627                 ir->rlist,rlist_new);
628         if (MASTER(cr))
629         {
630             fprintf(stderr,"%s\n\n",buf);
631         }
632         if (fp != NULL)
633         {
634             fprintf(fp,"%s\n\n",buf);
635         }
636         ir->rlist     = rlist_new;
637         ir->rlistlong = rlist_new;
638     }
639 }
640
641 static void prepare_verlet_scheme(FILE *fplog,
642                                   gmx_hw_info_t *hwinfo,
643                                   t_commrec *cr,
644                                   gmx_hw_opt_t *hw_opt,
645                                   const char *nbpu_opt,
646                                   t_inputrec *ir,
647                                   const gmx_mtop_t *mtop,
648                                   matrix box,
649                                   gmx_bool *bUseGPU)
650 {
651     /* Here we only check for GPU usage on the MPI master process,
652      * as here we don't know how many GPUs we will use yet.
653      * We check for a GPU on all processes later.
654      */
655     *bUseGPU = hwinfo->bCanUseGPU || (getenv("GMX_EMULATE_GPU") != NULL);
656
657     if (ir->verletbuf_drift > 0)
658     {
659         /* Update the Verlet buffer size for the current run setup */
660         verletbuf_list_setup_t ls;
661         real rlist_new;
662
663         /* Here we assume CPU acceleration is on. But as currently
664          * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
665          * and 4x2 gives a larger buffer than 4x4, this is ok.
666          */
667         verletbuf_get_list_setup(*bUseGPU,&ls);
668
669         calc_verlet_buffer_size(mtop,det(box),ir,
670                                 ir->verletbuf_drift,&ls,
671                                 NULL,&rlist_new);
672         if (rlist_new != ir->rlist)
673         {
674             if (fplog != NULL)
675             {
676                 fprintf(fplog,"\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
677                         ir->rlist,rlist_new,
678                         ls.cluster_size_i,ls.cluster_size_j);
679             }
680             ir->rlist     = rlist_new;
681             ir->rlistlong = rlist_new;
682         }
683     }
684
685     /* With GPU or emulation we should check nstlist for performance */
686     if ((EI_DYNAMICS(ir->eI) &&
687          *bUseGPU &&
688          ir->nstlist < NSTLIST_GPU_ENOUGH) ||
689         getenv(NSTLIST_ENVVAR) != NULL)
690     {
691         /* Choose a better nstlist */
692         increase_nstlist(fplog,cr,ir,mtop,box);
693     }
694 }
695
696 static void convert_to_verlet_scheme(FILE *fplog,
697                                      t_inputrec *ir,
698                                      gmx_mtop_t *mtop,real box_vol)
699 {
700     char *conv_mesg="Converting input file with group cut-off scheme to the Verlet cut-off scheme";
701
702     md_print_warn(NULL,fplog,"%s\n",conv_mesg);
703
704     ir->cutoff_scheme   = ecutsVERLET;
705     ir->verletbuf_drift = 0.005;
706
707     if (ir->rcoulomb != ir->rvdw)
708     {
709         gmx_fatal(FARGS,"The VdW and Coulomb cut-offs are different, whereas the Verlet scheme only supports equal cut-offs");
710     }
711
712     if (ir->vdwtype == evdwUSER || EEL_USER(ir->coulombtype))
713     {
714         gmx_fatal(FARGS,"User non-bonded potentials are not (yet) supported with the Verlet scheme");
715     }
716     else if (EVDW_SWITCHED(ir->vdwtype) || EEL_SWITCHED(ir->coulombtype))
717     {
718         md_print_warn(NULL,fplog,"Converting switched or shifted interactions to a shifted potential (without force shift), this will lead to slightly different interaction potentials");
719
720         if (EVDW_SWITCHED(ir->vdwtype))
721         {
722             ir->vdwtype = evdwCUT;
723         }
724         if (EEL_SWITCHED(ir->coulombtype))
725         {
726             if (EEL_FULL(ir->coulombtype))
727             {
728                 /* With full electrostatic only PME can be switched */
729                 ir->coulombtype = eelPME;
730             }
731             else
732             {
733                 md_print_warn(NULL,fplog,"NOTE: Replacing %s electrostatics with reaction-field with epsilon-rf=inf\n",eel_names[ir->coulombtype]);
734                 ir->coulombtype = eelRF;
735                 ir->epsilon_rf  = 0.0;
736             }
737         }
738
739         /* We set the target energy drift to a small number.
740          * Note that this is only for testing. For production the user
741          * should think about this and set the mdp options.
742          */
743         ir->verletbuf_drift = 1e-4;
744     }
745
746     if (inputrec2nboundeddim(ir) != 3)
747     {
748         gmx_fatal(FARGS,"Can only convert old tpr files to the Verlet cut-off scheme with 3D pbc");
749     }
750
751     if (ir->efep != efepNO || ir->implicit_solvent != eisNO)
752     {
753         gmx_fatal(FARGS,"Will not convert old tpr files to the Verlet cut-off scheme with free-energy calculations or implicit solvent");
754     }
755
756     if (EI_DYNAMICS(ir->eI) && !(EI_MD(ir->eI) && ir->etc == etcNO))
757     {
758         verletbuf_list_setup_t ls;
759
760         verletbuf_get_list_setup(FALSE,&ls);
761         calc_verlet_buffer_size(mtop,box_vol,ir,ir->verletbuf_drift,&ls,
762                                 NULL,&ir->rlist);
763     }
764     else
765     {
766         ir->verletbuf_drift = -1;
767         ir->rlist           = 1.05*max(ir->rvdw,ir->rcoulomb);
768     }
769
770     gmx_mtop_remove_chargegroups(mtop);
771 }
772
773 /* Check the process affinity mask and if it is found to be non-zero,
774  * will honor it and disable mdrun internal affinity setting.
775  * This function should be called first before the OpenMP library gets
776  * initialized with the last argument FALSE (which will detect affinity
777  * set by external tools like taskset), and later, after the OpenMP
778  * initialization, with the last argument TRUE to detect affinity changes
779  * made by the OpenMP library.
780  *
781  * Note that this will only work on Linux as we use a GNU feature. */
782 static void check_cpu_affinity_set(FILE *fplog, const t_commrec *cr,
783                                    gmx_hw_opt_t *hw_opt, int ncpus,
784                                    gmx_bool bAfterOpenmpInit)
785 {
786 #ifdef HAVE_SCHED_GETAFFINITY
787     cpu_set_t mask_current;
788     int       i, ret, cpu_count, cpu_set;
789     gmx_bool  bAllSet;
790
791     assert(hw_opt);
792     if (!hw_opt->bThreadPinning)
793     {
794         /* internal affinity setting is off, don't bother checking process affinity */
795         return;
796     }
797
798     CPU_ZERO(&mask_current);
799     if ((ret = sched_getaffinity(0, sizeof(cpu_set_t), &mask_current)) != 0)
800     {
801         /* failed to query affinity mask, will just return */
802         if (debug)
803         {
804             fprintf(debug, "Failed to query affinity mask (error %d)", ret);
805         }
806         return;
807     }
808
809     /* Before proceeding with the actual check, make sure that the number of
810      * detected CPUs is >= the CPUs in the current set. */
811     if (ncpus < CPU_COUNT(&mask_current))
812     {
813         if (debug)
814         {
815             fprintf(debug, "%d CPUs detected, but %d was returned by CPU_COUNT",
816                     ncpus, CPU_COUNT(&mask_current));
817         }
818         return;
819     }
820
821     bAllSet = TRUE;
822     for (i = 0; (i < ncpus && i < CPU_SETSIZE); i++)
823     {
824         bAllSet = bAllSet && (CPU_ISSET(i, &mask_current) != 0);
825     }
826
827     if (!bAllSet)
828     {
829         if (!bAfterOpenmpInit)
830         {
831             md_print_warn(cr, fplog,
832                           "Non-default process affinity set, disabling internal affinity");
833         }
834         else
835         {
836             md_print_warn(cr, fplog,
837                           "Non-default process affinity set probably by the OpenMP library, "
838                           "disabling internal affinity");
839         }
840         hw_opt->bThreadPinning = FALSE;
841
842         if (debug)
843         {
844             fprintf(debug, "Non-default affinity mask found\n");
845         }
846     }
847     else
848     {
849         if (debug)
850         {
851             fprintf(debug, "Default affinity mask found\n");
852         }
853     }
854 #endif /* HAVE_SCHED_GETAFFINITY */
855 }
856
857 /* Set CPU affinity. Can be important for performance.
858    On some systems (e.g. Cray) CPU Affinity is set by default.
859    But default assigning doesn't work (well) with only some ranks
860    having threads. This causes very low performance.
861    External tools have cumbersome syntax for setting affinity
862    in the case that only some ranks have threads.
863    Thus it is important that GROMACS sets the affinity internally
864    if only PME is using threads.
865 */
866 static void set_cpu_affinity(FILE *fplog,
867                              const t_commrec *cr,
868                              gmx_hw_opt_t *hw_opt,
869                              int nthreads_pme,
870                              const gmx_hw_info_t *hwinfo,
871                              const t_inputrec *inputrec)
872 {
873 #if defined GMX_THREAD_MPI
874     /* With the number of TMPI threads equal to the number of cores
875      * we already pinned in thread-MPI, so don't pin again here.
876      */
877     if (hw_opt->nthreads_tmpi == tMPI_Thread_get_hw_number())
878     {
879         return;
880     }
881 #endif
882
883 #ifdef GMX_OPENMP /* TODO: actually we could do this even without OpenMP?! */
884 #ifdef HAVE_SCHED_SETAFFINITY
885     if (hw_opt->bThreadPinning)
886     {
887         int thread, nthread_local, nthread_node, nthread_hw_max, nphyscore;
888         int offset;
889         char *env;
890
891         /* threads on this MPI process or TMPI thread */
892         if (cr->duty & DUTY_PP)
893         {
894             nthread_local = gmx_omp_nthreads_get(emntNonbonded);
895         }
896         else
897         {
898             nthread_local = gmx_omp_nthreads_get(emntPME);
899         }
900
901         /* map the current process to cores */
902         thread = 0;
903         nthread_node = nthread_local;
904 #ifdef GMX_MPI
905         if (PAR(cr) || MULTISIM(cr))
906         {
907             /* We need to determine a scan of the thread counts in this
908              * compute node.
909              */
910             MPI_Comm comm_intra;
911
912             MPI_Comm_split(MPI_COMM_WORLD,gmx_hostname_num(),cr->nodeid_intra,
913                            &comm_intra);
914             MPI_Scan(&nthread_local,&thread,1,MPI_INT,MPI_SUM,comm_intra);
915             /* MPI_Scan is inclusive, but here we need exclusive */
916             thread -= nthread_local;
917             /* Get the total number of threads on this physical node */
918             MPI_Allreduce(&nthread_local,&nthread_node,1,MPI_INT,MPI_SUM,comm_intra);
919             MPI_Comm_free(&comm_intra);
920         }
921 #endif
922
923         offset = 0;
924         if (hw_opt->core_pinning_offset > 0)
925         {
926             offset = hw_opt->core_pinning_offset;
927             if (SIMMASTER(cr))
928             {
929                 fprintf(stderr, "Applying core pinning offset %d\n", offset);
930             }
931             if (fplog)
932             {
933                 fprintf(fplog, "Applying core pinning offset %d\n", offset);
934             }
935         }
936
937         /* With Intel Hyper-Threading enabled, we want to pin consecutive
938          * threads to physical cores when using more threads than physical
939          * cores or when the user requests so.
940          */
941         nthread_hw_max = hwinfo->nthreads_hw_avail;
942         nphyscore = -1;
943         if (hw_opt->bPinHyperthreading ||
944             (gmx_cpuid_x86_smt(hwinfo->cpuid_info) == GMX_CPUID_X86_SMT_ENABLED &&
945              nthread_node > nthread_hw_max/2 && getenv("GMX_DISABLE_PINHT") == NULL))
946         {
947             if (gmx_cpuid_x86_smt(hwinfo->cpuid_info) != GMX_CPUID_X86_SMT_ENABLED)
948             {
949                 /* We print to stderr on all processes, as we might have
950                  * different settings on different physical nodes.
951                  */
952                 if (gmx_cpuid_vendor(hwinfo->cpuid_info) != GMX_CPUID_VENDOR_INTEL)
953                 {
954                     md_print_warn(NULL, fplog, "Pinning for Hyper-Threading layout requested, "
955                                   "but non-Intel CPU detected (vendor: %s)\n",
956                                   gmx_cpuid_vendor_string[gmx_cpuid_vendor(hwinfo->cpuid_info)]);
957                 }
958                 else
959                 {
960                     md_print_warn(NULL, fplog, "Pinning for Hyper-Threading layout requested, "
961                                   "but the CPU detected does not have Intel Hyper-Threading support "
962                                   "(or it is turned off)\n");
963                 }
964             }
965             nphyscore = nthread_hw_max/2;
966
967             if (SIMMASTER(cr))
968             {
969                 fprintf(stderr, "Pinning to Hyper-Threading cores with %d physical cores in a compute node\n",
970                         nphyscore);
971             }
972             if (fplog)
973             {
974                 fprintf(fplog, "Pinning to Hyper-Threading cores with %d physical cores in a compute node\n",
975                         nphyscore);
976             }
977         }
978
979         /* set the per-thread affinity */
980 #pragma omp parallel firstprivate(thread) num_threads(nthread_local)
981         {
982             cpu_set_t mask;
983             int core;
984
985             CPU_ZERO(&mask);
986             thread += gmx_omp_get_thread_num();
987             if (nphyscore <= 0)
988             {
989                 core = offset + thread;
990             }
991             else
992             {
993                 /* Lock pairs of threads to the same hyperthreaded core */
994                 core = offset + thread/2 + (thread % 2)*nphyscore;
995             }
996             CPU_SET(core, &mask);
997             sched_setaffinity((pid_t) syscall (SYS_gettid), sizeof(cpu_set_t), &mask);
998         }
999     }
1000 #endif /* HAVE_SCHED_SETAFFINITY */
1001 #endif /* GMX_OPENMP */
1002 }
1003
1004
1005 static void check_and_update_hw_opt(gmx_hw_opt_t *hw_opt,
1006                                     int cutoff_scheme)
1007 {
1008     gmx_omp_nthreads_read_env(&hw_opt->nthreads_omp);
1009
1010 #ifndef GMX_THREAD_MPI
1011     if (hw_opt->nthreads_tot > 0)
1012     {
1013         gmx_fatal(FARGS,"Setting the total number of threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
1014     }
1015     if (hw_opt->nthreads_tmpi > 0)
1016     {
1017         gmx_fatal(FARGS,"Setting the number of thread-MPI threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
1018     }
1019 #endif
1020
1021     if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp_pme <= 0)
1022     {
1023         /* We have the same number of OpenMP threads for PP and PME processes,
1024          * thus we can perform several consistency checks.
1025          */
1026         if (hw_opt->nthreads_tmpi > 0 &&
1027             hw_opt->nthreads_omp > 0 &&
1028             hw_opt->nthreads_tot != hw_opt->nthreads_tmpi*hw_opt->nthreads_omp)
1029         {
1030             gmx_fatal(FARGS,"The total number of threads requested (%d) does not match the thread-MPI threads (%d) times the OpenMP threads (%d) requested",
1031                       hw_opt->nthreads_tot,hw_opt->nthreads_tmpi,hw_opt->nthreads_omp);
1032         }
1033
1034         if (hw_opt->nthreads_tmpi > 0 &&
1035             hw_opt->nthreads_tot % hw_opt->nthreads_tmpi != 0)
1036         {
1037             gmx_fatal(FARGS,"The total number of threads requested (%d) is not divisible by the number of thread-MPI threads requested (%d)",
1038                       hw_opt->nthreads_tot,hw_opt->nthreads_tmpi);
1039         }
1040
1041         if (hw_opt->nthreads_omp > 0 &&
1042             hw_opt->nthreads_tot % hw_opt->nthreads_omp != 0)
1043         {
1044             gmx_fatal(FARGS,"The total number of threads requested (%d) is not divisible by the number of OpenMP threads requested (%d)",
1045                       hw_opt->nthreads_tot,hw_opt->nthreads_omp);
1046         }
1047
1048         if (hw_opt->nthreads_tmpi > 0 &&
1049             hw_opt->nthreads_omp <= 0)
1050         {
1051             hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
1052         }
1053     }
1054
1055 #ifndef GMX_OPENMP
1056     if (hw_opt->nthreads_omp > 1)
1057     {
1058         gmx_fatal(FARGS,"OpenMP threads are requested, but Gromacs was compiled without OpenMP support");
1059     }
1060 #endif
1061
1062     if (cutoff_scheme == ecutsGROUP)
1063     {
1064         /* We only have OpenMP support for PME only nodes */
1065         if (hw_opt->nthreads_omp > 1)
1066         {
1067             gmx_fatal(FARGS,"OpenMP threads have been requested with cut-off scheme %s, but these are only supported with cut-off scheme %s",
1068                       ecutscheme_names[cutoff_scheme],
1069                       ecutscheme_names[ecutsVERLET]);
1070         }
1071         hw_opt->nthreads_omp = 1;
1072     }
1073
1074     if (hw_opt->nthreads_omp_pme > 0 && hw_opt->nthreads_omp <= 0)
1075     {
1076         gmx_fatal(FARGS,"You need to specify -ntomp in addition to -ntomp_pme");
1077     }
1078
1079     if (hw_opt->nthreads_tot == 1)
1080     {
1081         hw_opt->nthreads_tmpi = 1;
1082
1083         if (hw_opt->nthreads_omp > 1)
1084         {
1085             gmx_fatal(FARGS,"You requested %d OpenMP threads with %d total threads",
1086                       hw_opt->nthreads_tmpi,hw_opt->nthreads_tot);
1087         }
1088         hw_opt->nthreads_omp = 1;
1089     }
1090
1091     if (hw_opt->nthreads_omp_pme <= 0 && hw_opt->nthreads_omp > 0)
1092     {
1093         hw_opt->nthreads_omp_pme = hw_opt->nthreads_omp;
1094     }
1095
1096     if (debug)
1097     {
1098         fprintf(debug,"hw_opt: nt %d ntmpi %d ntomp %d ntomp_pme %d gpu_id '%s'\n",
1099                 hw_opt->nthreads_tot,
1100                 hw_opt->nthreads_tmpi,
1101                 hw_opt->nthreads_omp,
1102                 hw_opt->nthreads_omp_pme,
1103                 hw_opt->gpu_id!=NULL ? hw_opt->gpu_id : "");
1104                 
1105     }
1106 }
1107
1108
1109 /* Override the value in inputrec with value passed on the command line (if any) */
1110 static void override_nsteps_cmdline(FILE *fplog,
1111                                     int nsteps_cmdline,
1112                                     t_inputrec *ir,
1113                                     const t_commrec *cr)
1114 {
1115     assert(ir);
1116     assert(cr);
1117
1118     /* override with anything else than the default -2 */
1119     if (nsteps_cmdline > -2)
1120     {
1121         char stmp[STRLEN];
1122
1123         ir->nsteps = nsteps_cmdline;
1124         if (EI_DYNAMICS(ir->eI))
1125         {
1126             sprintf(stmp, "Overriding nsteps with value passed on the command line: %d steps, %.3f ps",
1127                     nsteps_cmdline, nsteps_cmdline*ir->delta_t);
1128         }
1129         else
1130         {
1131             sprintf(stmp, "Overriding nsteps with value passed on the command line: %d steps",
1132                     nsteps_cmdline);
1133         }
1134
1135         md_print_warn(cr, fplog, "%s\n", stmp);
1136     }
1137 }
1138
1139 /* Data structure set by SIMMASTER which needs to be passed to all nodes
1140  * before the other nodes have read the tpx file and called gmx_detect_hardware.
1141  */
1142 typedef struct {
1143     int cutoff_scheme; /* The cutoff scheme from inputrec_t */
1144     gmx_bool bUseGPU;       /* Use GPU or GPU emulation          */
1145 } master_inf_t;
1146
1147 int mdrunner(gmx_hw_opt_t *hw_opt,
1148              FILE *fplog,t_commrec *cr,int nfile,
1149              const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
1150              gmx_bool bCompact, int nstglobalcomm,
1151              ivec ddxyz,int dd_node_order,real rdd,real rconstr,
1152              const char *dddlb_opt,real dlb_scale,
1153              const char *ddcsx,const char *ddcsy,const char *ddcsz,
1154              const char *nbpu_opt,
1155              int nsteps_cmdline, int nstepout,int resetstep,
1156              int nmultisim,int repl_ex_nst,int repl_ex_nex,
1157              int repl_ex_seed, real pforce,real cpt_period,real max_hours,
1158              const char *deviceOptions, unsigned long Flags)
1159 {
1160     gmx_bool   bForceUseGPU,bTryUseGPU;
1161     double     nodetime=0,realtime;
1162     t_inputrec *inputrec;
1163     t_state    *state=NULL;
1164     matrix     box;
1165     gmx_ddbox_t ddbox={0};
1166     int        npme_major,npme_minor;
1167     real       tmpr1,tmpr2;
1168     t_nrnb     *nrnb;
1169     gmx_mtop_t *mtop=NULL;
1170     t_mdatoms  *mdatoms=NULL;
1171     t_forcerec *fr=NULL;
1172     t_fcdata   *fcd=NULL;
1173     real       ewaldcoeff=0;
1174     gmx_pme_t  *pmedata=NULL;
1175     gmx_vsite_t *vsite=NULL;
1176     gmx_constr_t constr;
1177     int        i,m,nChargePerturbed=-1,status,nalloc;
1178     char       *gro;
1179     gmx_wallcycle_t wcycle;
1180     gmx_bool       bReadRNG,bReadEkin;
1181     int        list;
1182     gmx_runtime_t runtime;
1183     int        rc;
1184     gmx_large_int_t reset_counters;
1185     gmx_edsam_t ed=NULL;
1186     t_commrec   *cr_old=cr; 
1187     int         nthreads_pme=1;
1188     int         nthreads_pp=1;
1189     gmx_membed_t membed=NULL;
1190     gmx_hw_info_t *hwinfo=NULL;
1191     master_inf_t minf={-1,FALSE};
1192
1193     /* CAUTION: threads may be started later on in this function, so
1194        cr doesn't reflect the final parallel state right now */
1195     snew(inputrec,1);
1196     snew(mtop,1);
1197     
1198     if (Flags & MD_APPENDFILES) 
1199     {
1200         fplog = NULL;
1201     }
1202
1203     bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
1204     bTryUseGPU   = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
1205
1206     snew(state,1);
1207     if (SIMMASTER(cr)) 
1208     {
1209         /* Read (nearly) all data required for the simulation */
1210         read_tpx_state(ftp2fn(efTPX,nfile,fnm),inputrec,state,NULL,mtop);
1211
1212         if (inputrec->cutoff_scheme != ecutsVERLET &&
1213             ((Flags & MD_TESTVERLET) || getenv("GMX_VERLET_SCHEME") != NULL))
1214         {
1215             convert_to_verlet_scheme(fplog,inputrec,mtop,det(state->box));
1216         }
1217
1218         /* Detect hardware, gather information. With tMPI only thread 0 does it
1219          * and after threads are started broadcasts hwinfo around. */
1220         snew(hwinfo, 1);
1221         gmx_detect_hardware(fplog, hwinfo, cr,
1222                             bForceUseGPU, bTryUseGPU, hw_opt->gpu_id);
1223
1224         minf.cutoff_scheme = inputrec->cutoff_scheme;
1225         minf.bUseGPU       = FALSE;
1226
1227         if (inputrec->cutoff_scheme == ecutsVERLET)
1228         {
1229             prepare_verlet_scheme(fplog,hwinfo,cr,hw_opt,nbpu_opt,
1230                                   inputrec,mtop,state->box,
1231                                   &minf.bUseGPU);
1232         }
1233         else if (hwinfo->bCanUseGPU)
1234         {
1235             md_print_warn(cr,fplog,
1236                           "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
1237                           "      To use a GPU, set the mdp option: cutoff-scheme = Verlet\n"
1238                           "      (for quick performance testing you can use the -testverlet option)\n");
1239
1240             if (bForceUseGPU)
1241             {
1242                 gmx_fatal(FARGS,"GPU requested, but can't be used without cutoff-scheme=Verlet");
1243             }
1244         }
1245     }
1246 #ifndef GMX_THREAD_MPI
1247     if (PAR(cr))
1248     {
1249         gmx_bcast_sim(sizeof(minf),&minf,cr);
1250     }
1251 #endif
1252     if (minf.bUseGPU && cr->npmenodes == -1)
1253     {
1254         /* Don't automatically use PME-only nodes with GPUs */
1255         cr->npmenodes = 0;
1256     }
1257
1258     /* Check for externally set OpenMP affinity and turn off internal
1259      * pinning if any is found. We need to do this check early to tell
1260      * thread-MPI whether it should do pinning when spawning threads.
1261      */
1262     gmx_omp_check_thread_affinity(fplog, cr, hw_opt);
1263
1264 #ifdef GMX_THREAD_MPI
1265     /* With thread-MPI inputrec is only set here on the master thread */
1266     if (SIMMASTER(cr))
1267 #endif
1268     {
1269         check_and_update_hw_opt(hw_opt,minf.cutoff_scheme);
1270
1271 #ifdef GMX_THREAD_MPI
1272         /* Early check for externally set process affinity. Can't do over all
1273          * MPI processes because hwinfo is not available everywhere, but with
1274          * thread-MPI it's needed as pinning might get turned off which needs
1275          * to be known before starting thread-MPI. */
1276         check_cpu_affinity_set(fplog,
1277                                NULL,
1278                                hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1279 #endif
1280
1281 #ifdef GMX_THREAD_MPI
1282         if (cr->npmenodes > 0 && hw_opt->nthreads_tmpi <= 0)
1283         {
1284             gmx_fatal(FARGS,"You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME nodes");
1285         }
1286 #endif
1287
1288         if (hw_opt->nthreads_omp_pme != hw_opt->nthreads_omp &&
1289             cr->npmenodes <= 0)
1290         {
1291             gmx_fatal(FARGS,"You need to explicitly specify the number of PME nodes (-npme) when using different number of OpenMP threads for PP and PME nodes");
1292         }
1293     }
1294
1295 #ifdef GMX_THREAD_MPI
1296     if (SIMMASTER(cr))
1297     {
1298         /* NOW the threads will be started: */
1299         hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
1300                                                  hw_opt,
1301                                                  inputrec, mtop,
1302                                                  cr, fplog);
1303         if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp <= 0)
1304         {
1305             hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
1306         }
1307
1308         if (hw_opt->nthreads_tmpi > 1)
1309         {
1310             /* now start the threads. */
1311             cr=mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm, 
1312                                       oenv, bVerbose, bCompact, nstglobalcomm, 
1313                                       ddxyz, dd_node_order, rdd, rconstr, 
1314                                       dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
1315                                       nbpu_opt,
1316                                       nsteps_cmdline, nstepout, resetstep, nmultisim, 
1317                                       repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
1318                                       cpt_period, max_hours, deviceOptions, 
1319                                       Flags);
1320             /* the main thread continues here with a new cr. We don't deallocate
1321                the old cr because other threads may still be reading it. */
1322             if (cr == NULL)
1323             {
1324                 gmx_comm("Failed to spawn threads");
1325             }
1326         }
1327     }
1328 #endif
1329     /* END OF CAUTION: cr is now reliable */
1330
1331     /* g_membed initialisation *
1332      * Because we change the mtop, init_membed is called before the init_parallel *
1333      * (in case we ever want to make it run in parallel) */
1334     if (opt2bSet("-membed",nfile,fnm))
1335     {
1336         if (MASTER(cr))
1337         {
1338             fprintf(stderr,"Initializing membed");
1339         }
1340         membed = init_membed(fplog,nfile,fnm,mtop,inputrec,state,cr,&cpt_period);
1341     }
1342
1343     if (PAR(cr))
1344     {
1345         /* now broadcast everything to the non-master nodes/threads: */
1346         init_parallel(fplog, cr, inputrec, mtop);
1347
1348         /* This check needs to happen after get_nthreads_mpi() */
1349         if (inputrec->cutoff_scheme == ecutsVERLET && (Flags & MD_PARTDEC))
1350         {
1351             gmx_fatal_collective(FARGS,cr,NULL,
1352                                  "The Verlet cut-off scheme is not supported with particle decomposition.\n"
1353                                  "You can achieve the same effect as particle decomposition by running in parallel using only OpenMP threads.");
1354         }
1355     }
1356     if (fplog != NULL)
1357     {
1358         pr_inputrec(fplog,0,"Input Parameters",inputrec,FALSE);
1359     }
1360
1361 #if defined GMX_THREAD_MPI
1362     /* With tMPI we detected on thread 0 and we'll just pass the hwinfo pointer
1363      * to the other threads  -- slightly uncool, but works fine, just need to
1364      * make sure that the data doesn't get freed twice. */
1365     if (cr->nnodes > 1)
1366     {
1367         if (!SIMMASTER(cr))
1368         {
1369             snew(hwinfo, 1);
1370         }
1371         gmx_bcast(sizeof(&hwinfo), &hwinfo, cr);
1372     }
1373 #else
1374     if (PAR(cr) && !SIMMASTER(cr))
1375     {
1376         /* now we have inputrec on all nodes, can run the detection */
1377         /* TODO: perhaps it's better to propagate within a node instead? */
1378         snew(hwinfo, 1);
1379         gmx_detect_hardware(fplog, hwinfo, cr,
1380                                  bForceUseGPU, bTryUseGPU, hw_opt->gpu_id);
1381     }
1382
1383     /* Now do the affinity check with MPI/no-MPI (done earlier with thread-MPI). */
1384     check_cpu_affinity_set(fplog, cr,
1385                            hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1386 #endif
1387
1388     /* now make sure the state is initialized and propagated */
1389     set_state_entries(state,inputrec,cr->nnodes);
1390
1391     /* remove when vv and rerun works correctly! */
1392     if (PAR(cr) && EI_VV(inputrec->eI) && ((Flags & MD_RERUN) || (Flags & MD_RERUN_VSITE)))
1393     {
1394         gmx_fatal(FARGS,
1395                   "Currently can't do velocity verlet with rerun in parallel.");
1396     }
1397
1398     /* A parallel command line option consistency check that we can
1399        only do after any threads have started. */
1400     if (!PAR(cr) &&
1401         (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
1402     {
1403         gmx_fatal(FARGS,
1404                   "The -dd or -npme option request a parallel simulation, "
1405 #ifndef GMX_MPI
1406                   "but %s was compiled without threads or MPI enabled"
1407 #else
1408 #ifdef GMX_THREAD_MPI
1409                   "but the number of threads (option -nt) is 1"
1410 #else
1411                   "but %s was not started through mpirun/mpiexec or only one process was requested through mpirun/mpiexec"
1412 #endif
1413 #endif
1414                   , ShortProgram()
1415             );
1416     }
1417
1418     if ((Flags & MD_RERUN) &&
1419         (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1420     {
1421         gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
1422     }
1423
1424     if (can_use_allvsall(inputrec,mtop,TRUE,cr,fplog) && PAR(cr))
1425     {
1426         /* All-vs-all loops do not work with domain decomposition */
1427         Flags |= MD_PARTDEC;
1428     }
1429
1430     if (!EEL_PME(inputrec->coulombtype) || (Flags & MD_PARTDEC))
1431     {
1432         if (cr->npmenodes > 0)
1433         {
1434             if (!EEL_PME(inputrec->coulombtype))
1435             {
1436                 gmx_fatal_collective(FARGS,cr,NULL,
1437                                      "PME nodes are requested, but the system does not use PME electrostatics");
1438             }
1439             if (Flags & MD_PARTDEC)
1440             {
1441                 gmx_fatal_collective(FARGS,cr,NULL,
1442                                      "PME nodes are requested, but particle decomposition does not support separate PME nodes");
1443             }
1444         }
1445
1446         cr->npmenodes = 0;
1447     }
1448
1449 #ifdef GMX_FAHCORE
1450     fcRegisterSteps(inputrec->nsteps,inputrec->init_step);
1451 #endif
1452
1453     /* NMR restraints must be initialized before load_checkpoint,
1454      * since with time averaging the history is added to t_state.
1455      * For proper consistency check we therefore need to extend
1456      * t_state here.
1457      * So the PME-only nodes (if present) will also initialize
1458      * the distance restraints.
1459      */
1460     snew(fcd,1);
1461
1462     /* This needs to be called before read_checkpoint to extend the state */
1463     init_disres(fplog,mtop,inputrec,cr,Flags & MD_PARTDEC,fcd,state);
1464
1465     if (gmx_mtop_ftype_count(mtop,F_ORIRES) > 0)
1466     {
1467         if (PAR(cr) && !(Flags & MD_PARTDEC))
1468         {
1469             gmx_fatal(FARGS,"Orientation restraints do not work (yet) with domain decomposition, use particle decomposition (mdrun option -pd)");
1470         }
1471         /* Orientation restraints */
1472         if (MASTER(cr))
1473         {
1474             init_orires(fplog,mtop,state->x,inputrec,cr->ms,&(fcd->orires),
1475                         state);
1476         }
1477     }
1478
1479     if (DEFORM(*inputrec))
1480     {
1481         /* Store the deform reference box before reading the checkpoint */
1482         if (SIMMASTER(cr))
1483         {
1484             copy_mat(state->box,box);
1485         }
1486         if (PAR(cr))
1487         {
1488             gmx_bcast(sizeof(box),box,cr);
1489         }
1490         /* Because we do not have the update struct available yet
1491          * in which the reference values should be stored,
1492          * we store them temporarily in static variables.
1493          * This should be thread safe, since they are only written once
1494          * and with identical values.
1495          */
1496 #ifdef GMX_THREAD_MPI
1497         tMPI_Thread_mutex_lock(&deform_init_box_mutex);
1498 #endif
1499         deform_init_init_step_tpx = inputrec->init_step;
1500         copy_mat(box,deform_init_box_tpx);
1501 #ifdef GMX_THREAD_MPI
1502         tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
1503 #endif
1504     }
1505
1506     if (opt2bSet("-cpi",nfile,fnm)) 
1507     {
1508         /* Check if checkpoint file exists before doing continuation.
1509          * This way we can use identical input options for the first and subsequent runs...
1510          */
1511         if( gmx_fexist_master(opt2fn_master("-cpi",nfile,fnm,cr),cr) )
1512         {
1513             load_checkpoint(opt2fn_master("-cpi",nfile,fnm,cr),&fplog,
1514                             cr,Flags & MD_PARTDEC,ddxyz,
1515                             inputrec,state,&bReadRNG,&bReadEkin,
1516                             (Flags & MD_APPENDFILES),
1517                             (Flags & MD_APPENDFILESSET));
1518             
1519             if (bReadRNG)
1520             {
1521                 Flags |= MD_READ_RNG;
1522             }
1523             if (bReadEkin)
1524             {
1525                 Flags |= MD_READ_EKIN;
1526             }
1527         }
1528     }
1529
1530     if (((MASTER(cr) || (Flags & MD_SEPPOT)) && (Flags & MD_APPENDFILES))
1531 #ifdef GMX_THREAD_MPI
1532         /* With thread MPI only the master node/thread exists in mdrun.c,
1533          * therefore non-master nodes need to open the "seppot" log file here.
1534          */
1535         || (!MASTER(cr) && (Flags & MD_SEPPOT))
1536 #endif
1537         )
1538     {
1539         gmx_log_open(ftp2fn(efLOG,nfile,fnm),cr,!(Flags & MD_SEPPOT),
1540                              Flags,&fplog);
1541     }
1542
1543     /* override nsteps with value from cmdline */
1544     override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
1545
1546     if (SIMMASTER(cr)) 
1547     {
1548         copy_mat(state->box,box);
1549     }
1550
1551     if (PAR(cr)) 
1552     {
1553         gmx_bcast(sizeof(box),box,cr);
1554     }
1555
1556     /* Essential dynamics */
1557     if (opt2bSet("-ei",nfile,fnm))
1558     {
1559         /* Open input and output files, allocate space for ED data structure */
1560         ed = ed_open(nfile,fnm,Flags,cr);
1561     }
1562
1563     if (PAR(cr) && !((Flags & MD_PARTDEC) ||
1564                      EI_TPI(inputrec->eI) ||
1565                      inputrec->eI == eiNM))
1566     {
1567         cr->dd = init_domain_decomposition(fplog,cr,Flags,ddxyz,rdd,rconstr,
1568                                            dddlb_opt,dlb_scale,
1569                                            ddcsx,ddcsy,ddcsz,
1570                                            mtop,inputrec,
1571                                            box,state->x,
1572                                            &ddbox,&npme_major,&npme_minor);
1573
1574         make_dd_communicators(fplog,cr,dd_node_order);
1575
1576         /* Set overallocation to avoid frequent reallocation of arrays */
1577         set_over_alloc_dd(TRUE);
1578     }
1579     else
1580     {
1581         /* PME, if used, is done on all nodes with 1D decomposition */
1582         cr->npmenodes = 0;
1583         cr->duty = (DUTY_PP | DUTY_PME);
1584         npme_major = 1;
1585         npme_minor = 1;
1586         if (!EI_TPI(inputrec->eI))
1587         {
1588             npme_major = cr->nnodes;
1589         }
1590         
1591         if (inputrec->ePBC == epbcSCREW)
1592         {
1593             gmx_fatal(FARGS,
1594                       "pbc=%s is only implemented with domain decomposition",
1595                       epbc_names[inputrec->ePBC]);
1596         }
1597     }
1598
1599     if (PAR(cr))
1600     {
1601         /* After possible communicator splitting in make_dd_communicators.
1602          * we can set up the intra/inter node communication.
1603          */
1604         gmx_setup_nodecomm(fplog,cr);
1605     }
1606
1607     /* Initialize per-node process ID and counters. */
1608     gmx_init_intra_counters(cr);
1609
1610 #ifdef GMX_MPI
1611     md_print_info(cr,fplog,"Using %d MPI %s\n",
1612                   cr->nnodes,
1613 #ifdef GMX_THREAD_MPI
1614                   cr->nnodes==1 ? "thread" : "threads"
1615 #else
1616                   cr->nnodes==1 ? "process" : "processes"
1617 #endif
1618                   );
1619 #endif
1620
1621     gmx_omp_nthreads_init(fplog, cr,
1622                           hwinfo->nthreads_hw_avail,
1623                           hw_opt->nthreads_omp,
1624                           hw_opt->nthreads_omp_pme,
1625                           (cr->duty & DUTY_PP) == 0,
1626                           inputrec->cutoff_scheme == ecutsVERLET);
1627
1628     gmx_check_hw_runconf_consistency(fplog, hwinfo, cr, hw_opt->nthreads_tmpi, minf.bUseGPU);
1629
1630     /* getting number of PP/PME threads
1631        PME: env variable should be read only on one node to make sure it is 
1632        identical everywhere;
1633      */
1634     /* TODO nthreads_pp is only used for pinning threads.
1635      * This is a temporary solution until we have a hw topology library.
1636      */
1637     nthreads_pp  = gmx_omp_nthreads_get(emntNonbonded);
1638     nthreads_pme = gmx_omp_nthreads_get(emntPME);
1639
1640     wcycle = wallcycle_init(fplog,resetstep,cr,nthreads_pp,nthreads_pme);
1641
1642     if (PAR(cr))
1643     {
1644         /* Master synchronizes its value of reset_counters with all nodes 
1645          * including PME only nodes */
1646         reset_counters = wcycle_get_reset_counters(wcycle);
1647         gmx_bcast_sim(sizeof(reset_counters),&reset_counters,cr);
1648         wcycle_set_reset_counters(wcycle, reset_counters);
1649     }
1650
1651     snew(nrnb,1);
1652     if (cr->duty & DUTY_PP)
1653     {
1654         /* For domain decomposition we allocate dynamically
1655          * in dd_partition_system.
1656          */
1657         if (DOMAINDECOMP(cr))
1658         {
1659             bcast_state_setup(cr,state);
1660         }
1661         else
1662         {
1663             if (PAR(cr))
1664             {
1665                 bcast_state(cr,state,TRUE);
1666             }
1667         }
1668
1669         /* Initiate forcerecord */
1670         fr = mk_forcerec();
1671         fr->hwinfo = hwinfo;
1672         init_forcerec(fplog,oenv,fr,fcd,inputrec,mtop,cr,box,FALSE,
1673                       opt2fn("-table",nfile,fnm),
1674                       opt2fn("-tabletf",nfile,fnm),
1675                       opt2fn("-tablep",nfile,fnm),
1676                       opt2fn("-tableb",nfile,fnm),
1677                       nbpu_opt,
1678                       FALSE,pforce);
1679
1680         /* version for PCA_NOT_READ_NODE (see md.c) */
1681         /*init_forcerec(fplog,fr,fcd,inputrec,mtop,cr,box,FALSE,
1682           "nofile","nofile","nofile","nofile",FALSE,pforce);
1683           */        
1684         fr->bSepDVDL = ((Flags & MD_SEPPOT) == MD_SEPPOT);
1685
1686         /* Initialize QM-MM */
1687         if(fr->bQMMM)
1688         {
1689             init_QMMMrec(cr,box,mtop,inputrec,fr);
1690         }
1691
1692         /* Initialize the mdatoms structure.
1693          * mdatoms is not filled with atom data,
1694          * as this can not be done now with domain decomposition.
1695          */
1696         mdatoms = init_mdatoms(fplog,mtop,inputrec->efep!=efepNO);
1697
1698         /* Initialize the virtual site communication */
1699         vsite = init_vsite(mtop,cr,FALSE);
1700
1701         calc_shifts(box,fr->shift_vec);
1702
1703         /* With periodic molecules the charge groups should be whole at start up
1704          * and the virtual sites should not be far from their proper positions.
1705          */
1706         if (!inputrec->bContinuation && MASTER(cr) &&
1707             !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1708         {
1709             /* Make molecules whole at start of run */
1710             if (fr->ePBC != epbcNONE)
1711             {
1712                 do_pbc_first_mtop(fplog,inputrec->ePBC,box,mtop,state->x);
1713             }
1714             if (vsite)
1715             {
1716                 /* Correct initial vsite positions are required
1717                  * for the initial distribution in the domain decomposition
1718                  * and for the initial shell prediction.
1719                  */
1720                 construct_vsites_mtop(fplog,vsite,mtop,state->x);
1721             }
1722         }
1723
1724         if (EEL_PME(fr->eeltype))
1725         {
1726             ewaldcoeff = fr->ewaldcoeff;
1727             pmedata = &fr->pmedata;
1728         }
1729         else
1730         {
1731             pmedata = NULL;
1732         }
1733     }
1734     else
1735     {
1736         /* This is a PME only node */
1737
1738         /* We don't need the state */
1739         done_state(state);
1740
1741         ewaldcoeff = calc_ewaldcoeff(inputrec->rcoulomb, inputrec->ewald_rtol);
1742         snew(pmedata,1);
1743     }
1744
1745     /* Before setting affinity, check whether the affinity has changed
1746      * - which indicates that probably the OpenMP library has changed it since
1747      * we first checked). */
1748     check_cpu_affinity_set(fplog, cr, hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1749
1750     /* Set the CPU affinity */
1751     set_cpu_affinity(fplog,cr,hw_opt,nthreads_pme,hwinfo,inputrec);
1752
1753     /* Initiate PME if necessary,
1754      * either on all nodes or on dedicated PME nodes only. */
1755     if (EEL_PME(inputrec->coulombtype))
1756     {
1757         if (mdatoms)
1758         {
1759             nChargePerturbed = mdatoms->nChargePerturbed;
1760         }
1761         if (cr->npmenodes > 0)
1762         {
1763             /* The PME only nodes need to know nChargePerturbed */
1764             gmx_bcast_sim(sizeof(nChargePerturbed),&nChargePerturbed,cr);
1765         }
1766
1767         if (cr->duty & DUTY_PME)
1768         {
1769             status = gmx_pme_init(pmedata,cr,npme_major,npme_minor,inputrec,
1770                                   mtop ? mtop->natoms : 0,nChargePerturbed,
1771                                   (Flags & MD_REPRODUCIBLE),nthreads_pme);
1772             if (status != 0) 
1773             {
1774                 gmx_fatal(FARGS,"Error %d initializing PME",status);
1775             }
1776         }
1777     }
1778
1779
1780     if (integrator[inputrec->eI].func == do_md
1781 #ifdef GMX_OPENMM
1782         ||
1783         integrator[inputrec->eI].func == do_md_openmm
1784 #endif
1785         )
1786     {
1787         /* Turn on signal handling on all nodes */
1788         /*
1789          * (A user signal from the PME nodes (if any)
1790          * is communicated to the PP nodes.
1791          */
1792         signal_handler_install();
1793     }
1794
1795     if (cr->duty & DUTY_PP)
1796     {
1797         if (inputrec->ePull != epullNO)
1798         {
1799             /* Initialize pull code */
1800             init_pull(fplog,inputrec,nfile,fnm,mtop,cr,oenv, inputrec->fepvals->init_lambda,
1801                       EI_DYNAMICS(inputrec->eI) && MASTER(cr),Flags);
1802         }
1803         
1804         if (inputrec->bRot)
1805         {
1806            /* Initialize enforced rotation code */
1807            init_rot(fplog,inputrec,nfile,fnm,cr,state->x,box,mtop,oenv,
1808                     bVerbose,Flags);
1809         }
1810
1811         constr = init_constraints(fplog,mtop,inputrec,ed,state,cr);
1812
1813         if (DOMAINDECOMP(cr))
1814         {
1815             dd_init_bondeds(fplog,cr->dd,mtop,vsite,constr,inputrec,
1816                             Flags & MD_DDBONDCHECK,fr->cginfo_mb);
1817
1818             set_dd_parameters(fplog,cr->dd,dlb_scale,inputrec,fr,&ddbox);
1819
1820             setup_dd_grid(fplog,cr->dd);
1821         }
1822
1823         /* Now do whatever the user wants us to do (how flexible...) */
1824         integrator[inputrec->eI].func(fplog,cr,nfile,fnm,
1825                                       oenv,bVerbose,bCompact,
1826                                       nstglobalcomm,
1827                                       vsite,constr,
1828                                       nstepout,inputrec,mtop,
1829                                       fcd,state,
1830                                       mdatoms,nrnb,wcycle,ed,fr,
1831                                       repl_ex_nst,repl_ex_nex,repl_ex_seed,
1832                                       membed,
1833                                       cpt_period,max_hours,
1834                                       deviceOptions,
1835                                       Flags,
1836                                       &runtime);
1837
1838         if (inputrec->ePull != epullNO)
1839         {
1840             finish_pull(fplog,inputrec->pull);
1841         }
1842         
1843         if (inputrec->bRot)
1844         {
1845             finish_rot(fplog,inputrec->rot);
1846         }
1847
1848     } 
1849     else 
1850     {
1851         /* do PME only */
1852         gmx_pmeonly(*pmedata,cr,nrnb,wcycle,ewaldcoeff,FALSE,inputrec);
1853     }
1854
1855     if (EI_DYNAMICS(inputrec->eI) || EI_TPI(inputrec->eI))
1856     {
1857         /* Some timing stats */  
1858         if (SIMMASTER(cr))
1859         {
1860             if (runtime.proc == 0)
1861             {
1862                 runtime.proc = runtime.real;
1863             }
1864         }
1865         else
1866         {
1867             runtime.real = 0;
1868         }
1869     }
1870
1871     wallcycle_stop(wcycle,ewcRUN);
1872
1873     /* Finish up, write some stuff
1874      * if rerunMD, don't write last frame again 
1875      */
1876     finish_run(fplog,cr,ftp2fn(efSTO,nfile,fnm),
1877                inputrec,nrnb,wcycle,&runtime,
1878                fr != NULL && fr->nbv != NULL && fr->nbv->bUseGPU ?
1879                  nbnxn_cuda_get_timings(fr->nbv->cu_nbv) : NULL,
1880                nthreads_pp, 
1881                EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1882
1883     if ((cr->duty & DUTY_PP) && fr->nbv != NULL && fr->nbv->bUseGPU)
1884     {
1885         char gpu_err_str[STRLEN];
1886
1887         /* free GPU memory and uninitialize GPU (by destroying the context) */
1888         nbnxn_cuda_free(fplog, fr->nbv->cu_nbv);
1889
1890         if (!free_gpu(gpu_err_str))
1891         {
1892             gmx_warning("On node %d failed to free GPU #%d: %s",
1893                         cr->nodeid, get_current_gpu_device_id(), gpu_err_str);
1894         }
1895     }
1896
1897     if (opt2bSet("-membed",nfile,fnm))
1898     {
1899         sfree(membed);
1900     }
1901
1902 #ifdef GMX_THREAD_MPI
1903     if (PAR(cr) && SIMMASTER(cr))
1904 #endif
1905     {
1906         gmx_hardware_info_free(hwinfo);
1907     }
1908
1909     /* Does what it says */  
1910     print_date_and_time(fplog,cr->nodeid,"Finished mdrun",&runtime);
1911
1912     /* Close logfile already here if we were appending to it */
1913     if (MASTER(cr) && (Flags & MD_APPENDFILES))
1914     {
1915         gmx_log_close(fplog);
1916     }   
1917
1918     rc=(int)gmx_get_stop_condition();
1919
1920 #ifdef GMX_THREAD_MPI
1921     /* we need to join all threads. The sub-threads join when they
1922        exit this function, but the master thread needs to be told to 
1923        wait for that. */
1924     if (PAR(cr) && MASTER(cr))
1925     {
1926         tMPI_Finalize();
1927     }
1928 #endif
1929
1930     return rc;
1931 }