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