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