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