Merge branch 'release-5-0'
[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,2015, 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 "gmxpre.h"
39
40 #include "config.h"
41
42 #include <assert.h>
43 #include <signal.h>
44 #include <stdlib.h>
45 #include <string.h>
46
47 #include <algorithm>
48
49 #include "gromacs/domdec/domdec.h"
50 #include "gromacs/essentialdynamics/edsam.h"
51 #include "gromacs/ewald/pme.h"
52 #include "gromacs/fileio/tpxio.h"
53 #include "gromacs/legacyheaders/checkpoint.h"
54 #include "gromacs/legacyheaders/constr.h"
55 #include "gromacs/legacyheaders/disre.h"
56 #include "gromacs/legacyheaders/force.h"
57 #include "gromacs/legacyheaders/gmx_detect_hardware.h"
58 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
59 #include "gromacs/legacyheaders/gmx_thread_affinity.h"
60 #include "gromacs/legacyheaders/inputrec.h"
61 #include "gromacs/legacyheaders/main.h"
62 #include "gromacs/legacyheaders/md_logging.h"
63 #include "gromacs/legacyheaders/md_support.h"
64 #include "gromacs/legacyheaders/mdatoms.h"
65 #include "gromacs/legacyheaders/mdrun.h"
66 #include "gromacs/legacyheaders/names.h"
67 #include "gromacs/legacyheaders/network.h"
68 #include "gromacs/legacyheaders/oenv.h"
69 #include "gromacs/legacyheaders/orires.h"
70 #include "gromacs/legacyheaders/qmmm.h"
71 #include "gromacs/legacyheaders/sighandler.h"
72 #include "gromacs/legacyheaders/txtdump.h"
73 #include "gromacs/legacyheaders/typedefs.h"
74 #include "gromacs/math/calculate-ewald-splitting-coefficient.h"
75 #include "gromacs/math/vec.h"
76 #include "gromacs/mdlib/calc_verletbuf.h"
77 #include "gromacs/mdlib/nbnxn_consts.h"
78 #include "gromacs/mdlib/nbnxn_search.h"
79 #include "gromacs/pbcutil/pbc.h"
80 #include "gromacs/pulling/pull.h"
81 #include "gromacs/pulling/pull_rotation.h"
82 #include "gromacs/swap/swapcoords.h"
83 #include "gromacs/timing/wallcycle.h"
84 #include "gromacs/topology/mtop_util.h"
85 #include "gromacs/utility/gmxassert.h"
86 #include "gromacs/utility/gmxmpi.h"
87 #include "gromacs/utility/smalloc.h"
88
89 #include "deform.h"
90 #include "membed.h"
91 #include "repl_ex.h"
92
93 #ifdef GMX_FAHCORE
94 #include "corewrap.h"
95 #endif
96
97 #include "gromacs/gmxlib/gpu_utils/gpu_utils.h"
98
99 typedef struct {
100     gmx_integrator_t *func;
101 } gmx_intp_t;
102
103 /* The array should match the eI array in include/types/enums.h */
104 const gmx_intp_t    integrator[eiNR] = { {do_md}, {do_steep}, {do_cg}, {do_md}, {do_md}, {do_nm}, {do_lbfgs}, {do_tpi}, {do_tpi}, {do_md}, {do_md}, {do_md}};
105
106 gmx_int64_t         deform_init_init_step_tpx;
107 matrix              deform_init_box_tpx;
108 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
109
110
111 #ifdef GMX_THREAD_MPI
112 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
113  * the number of threads will get lowered.
114  */
115 #define MIN_ATOMS_PER_MPI_THREAD    90
116 #define MIN_ATOMS_PER_GPU           900
117
118 struct mdrunner_arglist
119 {
120     gmx_hw_opt_t    hw_opt;
121     FILE           *fplog;
122     t_commrec      *cr;
123     int             nfile;
124     const t_filenm *fnm;
125     output_env_t    oenv;
126     gmx_bool        bVerbose;
127     gmx_bool        bCompact;
128     int             nstglobalcomm;
129     ivec            ddxyz;
130     int             dd_node_order;
131     real            rdd;
132     real            rconstr;
133     const char     *dddlb_opt;
134     real            dlb_scale;
135     const char     *ddcsx;
136     const char     *ddcsy;
137     const char     *ddcsz;
138     const char     *nbpu_opt;
139     int             nstlist_cmdline;
140     gmx_int64_t     nsteps_cmdline;
141     int             nstepout;
142     int             resetstep;
143     int             nmultisim;
144     int             repl_ex_nst;
145     int             repl_ex_nex;
146     int             repl_ex_seed;
147     real            pforce;
148     real            cpt_period;
149     real            max_hours;
150     const char     *deviceOptions;
151     int             imdport;
152     unsigned long   Flags;
153 };
154
155
156 /* The function used for spawning threads. Extracts the mdrunner()
157    arguments from its one argument and calls mdrunner(), after making
158    a commrec. */
159 static void mdrunner_start_fn(void *arg)
160 {
161     struct mdrunner_arglist *mda = (struct mdrunner_arglist*)arg;
162     struct mdrunner_arglist  mc  = *mda; /* copy the arg list to make sure
163                                             that it's thread-local. This doesn't
164                                             copy pointed-to items, of course,
165                                             but those are all const. */
166     t_commrec *cr;                       /* we need a local version of this */
167     FILE      *fplog = NULL;
168     t_filenm  *fnm;
169
170     fnm = dup_tfn(mc.nfile, mc.fnm);
171
172     cr = reinitialize_commrec_for_this_thread(mc.cr);
173
174     if (MASTER(cr))
175     {
176         fplog = mc.fplog;
177     }
178
179     mdrunner(&mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
180              mc.bVerbose, mc.bCompact, mc.nstglobalcomm,
181              mc.ddxyz, mc.dd_node_order, mc.rdd,
182              mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
183              mc.ddcsx, mc.ddcsy, mc.ddcsz,
184              mc.nbpu_opt, mc.nstlist_cmdline,
185              mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
186              mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce,
187              mc.cpt_period, mc.max_hours, mc.deviceOptions, mc.imdport, mc.Flags);
188 }
189
190 /* called by mdrunner() to start a specific number of threads (including
191    the main thread) for thread-parallel runs. This in turn calls mdrunner()
192    for each thread.
193    All options besides nthreads are the same as for mdrunner(). */
194 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
195                                          FILE *fplog, t_commrec *cr, int nfile,
196                                          const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
197                                          gmx_bool bCompact, int nstglobalcomm,
198                                          ivec ddxyz, int dd_node_order, real rdd, real rconstr,
199                                          const char *dddlb_opt, real dlb_scale,
200                                          const char *ddcsx, const char *ddcsy, const char *ddcsz,
201                                          const char *nbpu_opt, int nstlist_cmdline,
202                                          gmx_int64_t nsteps_cmdline,
203                                          int nstepout, int resetstep,
204                                          int nmultisim, int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
205                                          real pforce, real cpt_period, real max_hours,
206                                          const char *deviceOptions, unsigned long Flags)
207 {
208     int                      ret;
209     struct mdrunner_arglist *mda;
210     t_commrec               *crn; /* the new commrec */
211     t_filenm                *fnmn;
212
213     /* first check whether we even need to start tMPI */
214     if (hw_opt->nthreads_tmpi < 2)
215     {
216         return cr;
217     }
218
219     /* a few small, one-time, almost unavoidable memory leaks: */
220     snew(mda, 1);
221     fnmn = dup_tfn(nfile, fnm);
222
223     /* fill the data structure to pass as void pointer to thread start fn */
224     /* hw_opt contains pointers, which should all be NULL at this stage */
225     mda->hw_opt          = *hw_opt;
226     mda->fplog           = fplog;
227     mda->cr              = cr;
228     mda->nfile           = nfile;
229     mda->fnm             = fnmn;
230     mda->oenv            = oenv;
231     mda->bVerbose        = bVerbose;
232     mda->bCompact        = bCompact;
233     mda->nstglobalcomm   = nstglobalcomm;
234     mda->ddxyz[XX]       = ddxyz[XX];
235     mda->ddxyz[YY]       = ddxyz[YY];
236     mda->ddxyz[ZZ]       = ddxyz[ZZ];
237     mda->dd_node_order   = dd_node_order;
238     mda->rdd             = rdd;
239     mda->rconstr         = rconstr;
240     mda->dddlb_opt       = dddlb_opt;
241     mda->dlb_scale       = dlb_scale;
242     mda->ddcsx           = ddcsx;
243     mda->ddcsy           = ddcsy;
244     mda->ddcsz           = ddcsz;
245     mda->nbpu_opt        = nbpu_opt;
246     mda->nstlist_cmdline = nstlist_cmdline;
247     mda->nsteps_cmdline  = nsteps_cmdline;
248     mda->nstepout        = nstepout;
249     mda->resetstep       = resetstep;
250     mda->nmultisim       = nmultisim;
251     mda->repl_ex_nst     = repl_ex_nst;
252     mda->repl_ex_nex     = repl_ex_nex;
253     mda->repl_ex_seed    = repl_ex_seed;
254     mda->pforce          = pforce;
255     mda->cpt_period      = cpt_period;
256     mda->max_hours       = max_hours;
257     mda->deviceOptions   = deviceOptions;
258     mda->Flags           = Flags;
259
260     /* now spawn new threads that start mdrunner_start_fn(), while
261        the main thread returns, we set thread affinity later */
262     ret = tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi, TMPI_AFFINITY_NONE,
263                        mdrunner_start_fn, (void*)(mda) );
264     if (ret != TMPI_SUCCESS)
265     {
266         return NULL;
267     }
268
269     crn = reinitialize_commrec_for_this_thread(cr);
270     return crn;
271 }
272
273
274 static int get_tmpi_omp_thread_division(const gmx_hw_info_t *hwinfo,
275                                         const gmx_hw_opt_t  *hw_opt,
276                                         int                  nthreads_tot,
277                                         int                  ngpu)
278 {
279     int nthreads_tmpi;
280
281     /* There are no separate PME nodes here, as we ensured in
282      * check_and_update_hw_opt that nthreads_tmpi>0 with PME nodes
283      * and a conditional ensures we would not have ended up here.
284      * Note that separate PME nodes might be switched on later.
285      */
286     if (ngpu > 0)
287     {
288         nthreads_tmpi = ngpu;
289         if (nthreads_tot > 0 && nthreads_tot < nthreads_tmpi)
290         {
291             nthreads_tmpi = nthreads_tot;
292         }
293     }
294     else if (hw_opt->nthreads_omp > 0)
295     {
296         /* Here we could oversubscribe, when we do, we issue a warning later */
297         nthreads_tmpi = std::max(1, nthreads_tot/hw_opt->nthreads_omp);
298     }
299     else
300     {
301         /* TODO choose nthreads_omp based on hardware topology
302            when we have a hardware topology detection library */
303         /* In general, when running up to 4 threads, OpenMP should be faster.
304          * Note: on AMD Bulldozer we should avoid running OpenMP over two dies.
305          * On Intel>=Nehalem running OpenMP on a single CPU is always faster,
306          * even on two CPUs it's usually faster (but with many OpenMP threads
307          * it could be faster not to use HT, currently we always use HT).
308          * On Nehalem/Westmere we want to avoid running 16 threads over
309          * two CPUs with HT, so we need a limit<16; thus we use 12.
310          * A reasonable limit for Intel Sandy and Ivy bridge,
311          * not knowing the topology, is 16 threads.
312          * Below we check for Intel and AVX, which for now includes
313          * Sandy/Ivy Bridge, Has/Broadwell. By checking for AVX instead of
314          * model numbers we ensure also future Intel CPUs are covered.
315          */
316         const int nthreads_omp_always_faster             =  4;
317         const int nthreads_omp_always_faster_Nehalem     = 12;
318         const int nthreads_omp_always_faster_Intel_AVX   = 16;
319         gmx_bool  bIntelAVX;
320
321         bIntelAVX =
322             (gmx_cpuid_vendor(hwinfo->cpuid_info) == GMX_CPUID_VENDOR_INTEL &&
323              gmx_cpuid_feature(hwinfo->cpuid_info, GMX_CPUID_FEATURE_X86_AVX));
324
325         if (nthreads_tot <= nthreads_omp_always_faster ||
326             ((gmx_cpuid_is_intel_nehalem(hwinfo->cpuid_info) && nthreads_tot <= nthreads_omp_always_faster_Nehalem) ||
327              (bIntelAVX && nthreads_tot <= nthreads_omp_always_faster_Intel_AVX)))
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 print_hw_opt(FILE *fp, const gmx_hw_opt_t *hw_opt)
759 {
760     fprintf(fp, "hw_opt: nt %d ntmpi %d ntomp %d ntomp_pme %d gpu_id '%s'\n",
761             hw_opt->nthreads_tot,
762             hw_opt->nthreads_tmpi,
763             hw_opt->nthreads_omp,
764             hw_opt->nthreads_omp_pme,
765             hw_opt->gpu_opt.gpu_id != NULL ? hw_opt->gpu_opt.gpu_id : "");
766 }
767
768 /* Checks we can do when we don't (yet) know the cut-off scheme */
769 static void check_and_update_hw_opt_1(gmx_hw_opt_t *hw_opt,
770                                       gmx_bool      bIsSimMaster)
771 {
772     gmx_omp_nthreads_read_env(&hw_opt->nthreads_omp, bIsSimMaster);
773
774 #ifndef GMX_THREAD_MPI
775     if (hw_opt->nthreads_tot > 0)
776     {
777         gmx_fatal(FARGS, "Setting the total number of threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
778     }
779     if (hw_opt->nthreads_tmpi > 0)
780     {
781         gmx_fatal(FARGS, "Setting the number of thread-MPI threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
782     }
783 #endif
784
785 #ifndef GMX_OPENMP
786     if (hw_opt->nthreads_omp > 1)
787     {
788         gmx_fatal(FARGS, "More than 1 OpenMP thread requested, but Gromacs was compiled without OpenMP support");
789     }
790     hw_opt->nthreads_omp = 1;
791 #endif
792
793     if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp_pme <= 0)
794     {
795         /* We have the same number of OpenMP threads for PP and PME processes,
796          * thus we can perform several consistency checks.
797          */
798         if (hw_opt->nthreads_tmpi > 0 &&
799             hw_opt->nthreads_omp > 0 &&
800             hw_opt->nthreads_tot != hw_opt->nthreads_tmpi*hw_opt->nthreads_omp)
801         {
802             gmx_fatal(FARGS, "The total number of threads requested (%d) does not match the thread-MPI threads (%d) times the OpenMP threads (%d) requested",
803                       hw_opt->nthreads_tot, hw_opt->nthreads_tmpi, hw_opt->nthreads_omp);
804         }
805
806         if (hw_opt->nthreads_tmpi > 0 &&
807             hw_opt->nthreads_tot % hw_opt->nthreads_tmpi != 0)
808         {
809             gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of thread-MPI threads requested (%d)",
810                       hw_opt->nthreads_tot, hw_opt->nthreads_tmpi);
811         }
812
813         if (hw_opt->nthreads_omp > 0 &&
814             hw_opt->nthreads_tot % hw_opt->nthreads_omp != 0)
815         {
816             gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of OpenMP threads requested (%d)",
817                       hw_opt->nthreads_tot, hw_opt->nthreads_omp);
818         }
819
820         if (hw_opt->nthreads_tmpi > 0 &&
821             hw_opt->nthreads_omp <= 0)
822         {
823             hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
824         }
825     }
826
827 #ifndef GMX_OPENMP
828     if (hw_opt->nthreads_omp > 1)
829     {
830         gmx_fatal(FARGS, "OpenMP threads are requested, but Gromacs was compiled without OpenMP support");
831     }
832 #endif
833
834     if (hw_opt->nthreads_omp_pme > 0 && hw_opt->nthreads_omp <= 0)
835     {
836         gmx_fatal(FARGS, "You need to specify -ntomp in addition to -ntomp_pme");
837     }
838
839     if (hw_opt->nthreads_tot == 1)
840     {
841         hw_opt->nthreads_tmpi = 1;
842
843         if (hw_opt->nthreads_omp > 1)
844         {
845             gmx_fatal(FARGS, "You requested %d OpenMP threads with %d total threads",
846                       hw_opt->nthreads_tmpi, hw_opt->nthreads_tot);
847         }
848         hw_opt->nthreads_omp = 1;
849     }
850
851     if (hw_opt->nthreads_omp_pme <= 0 && hw_opt->nthreads_omp > 0)
852     {
853         hw_opt->nthreads_omp_pme = hw_opt->nthreads_omp;
854     }
855
856     /* Parse GPU IDs, if provided.
857      * We check consistency with the tMPI thread count later.
858      */
859     gmx_parse_gpu_ids(&hw_opt->gpu_opt);
860
861 #ifdef GMX_THREAD_MPI
862     if (hw_opt->gpu_opt.ncuda_dev_use > 0 && hw_opt->nthreads_tmpi == 0)
863     {
864         /* Set the number of MPI threads equal to the number of GPUs */
865         hw_opt->nthreads_tmpi = hw_opt->gpu_opt.ncuda_dev_use;
866
867         if (hw_opt->nthreads_tot > 0 &&
868             hw_opt->nthreads_tmpi > hw_opt->nthreads_tot)
869         {
870             /* We have more GPUs than total threads requested.
871              * We choose to (later) generate a mismatch error,
872              * instead of launching more threads than requested.
873              */
874             hw_opt->nthreads_tmpi = hw_opt->nthreads_tot;
875         }
876     }
877 #endif
878
879     if (debug)
880     {
881         print_hw_opt(debug, hw_opt);
882     }
883 }
884
885 /* Checks we can do when we know the cut-off scheme */
886 static void check_and_update_hw_opt_2(gmx_hw_opt_t *hw_opt,
887                                       int           cutoff_scheme)
888 {
889     if (cutoff_scheme == ecutsGROUP)
890     {
891         /* We only have OpenMP support for PME only nodes */
892         if (hw_opt->nthreads_omp > 1)
893         {
894             gmx_fatal(FARGS, "OpenMP threads have been requested with cut-off scheme %s, but these are only supported with cut-off scheme %s",
895                       ecutscheme_names[cutoff_scheme],
896                       ecutscheme_names[ecutsVERLET]);
897         }
898         hw_opt->nthreads_omp = 1;
899     }
900
901     if (hw_opt->nthreads_omp_pme <= 0 && hw_opt->nthreads_omp > 0)
902     {
903         hw_opt->nthreads_omp_pme = hw_opt->nthreads_omp;
904     }
905
906     if (debug)
907     {
908         print_hw_opt(debug, hw_opt);
909     }
910 }
911
912
913 /* Override the value in inputrec with value passed on the command line (if any) */
914 static void override_nsteps_cmdline(FILE            *fplog,
915                                     gmx_int64_t      nsteps_cmdline,
916                                     t_inputrec      *ir,
917                                     const t_commrec *cr)
918 {
919     assert(ir);
920     assert(cr);
921
922     /* override with anything else than the default -2 */
923     if (nsteps_cmdline > -2)
924     {
925         char sbuf_steps[STEPSTRSIZE];
926         char sbuf_msg[STRLEN];
927
928         ir->nsteps = nsteps_cmdline;
929         if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
930         {
931             sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
932                     gmx_step_str(nsteps_cmdline, sbuf_steps),
933                     fabs(nsteps_cmdline*ir->delta_t));
934         }
935         else
936         {
937             sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
938                     gmx_step_str(nsteps_cmdline, sbuf_steps));
939         }
940
941         md_print_warn(cr, fplog, "%s\n", sbuf_msg);
942     }
943     else if (nsteps_cmdline < -2)
944     {
945         gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %d",
946                   nsteps_cmdline);
947     }
948     /* Do nothing if nsteps_cmdline == -2 */
949 }
950
951 int mdrunner(gmx_hw_opt_t *hw_opt,
952              FILE *fplog, t_commrec *cr, int nfile,
953              const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
954              gmx_bool bCompact, int nstglobalcomm,
955              ivec ddxyz, int dd_node_order, real rdd, real rconstr,
956              const char *dddlb_opt, real dlb_scale,
957              const char *ddcsx, const char *ddcsy, const char *ddcsz,
958              const char *nbpu_opt, int nstlist_cmdline,
959              gmx_int64_t nsteps_cmdline, int nstepout, int resetstep,
960              int gmx_unused nmultisim, int repl_ex_nst, int repl_ex_nex,
961              int repl_ex_seed, real pforce, real cpt_period, real max_hours,
962              const char *deviceOptions, int imdport, unsigned long Flags)
963 {
964     gmx_bool                  bForceUseGPU, bTryUseGPU, bRerunMD, bCantUseGPU;
965     t_inputrec               *inputrec;
966     t_state                  *state = NULL;
967     matrix                    box;
968     gmx_ddbox_t               ddbox = {0};
969     int                       npme_major, npme_minor;
970     t_nrnb                   *nrnb;
971     gmx_mtop_t               *mtop          = NULL;
972     t_mdatoms                *mdatoms       = NULL;
973     t_forcerec               *fr            = NULL;
974     t_fcdata                 *fcd           = NULL;
975     real                      ewaldcoeff_q  = 0;
976     real                      ewaldcoeff_lj = 0;
977     struct gmx_pme_t        **pmedata       = NULL;
978     gmx_vsite_t              *vsite         = NULL;
979     gmx_constr_t              constr;
980     int                       nChargePerturbed = -1, nTypePerturbed = 0, status;
981     gmx_wallcycle_t           wcycle;
982     gmx_bool                  bReadEkin;
983     gmx_walltime_accounting_t walltime_accounting = NULL;
984     int                       rc;
985     gmx_int64_t               reset_counters;
986     gmx_edsam_t               ed           = NULL;
987     int                       nthreads_pme = 1;
988     int                       nthreads_pp  = 1;
989     gmx_membed_t              membed       = NULL;
990     gmx_hw_info_t            *hwinfo       = NULL;
991     /* The master rank decides early on bUseGPU and broadcasts this later */
992     gmx_bool                  bUseGPU      = FALSE;
993
994     /* CAUTION: threads may be started later on in this function, so
995        cr doesn't reflect the final parallel state right now */
996     snew(inputrec, 1);
997     snew(mtop, 1);
998
999     if (Flags & MD_APPENDFILES)
1000     {
1001         fplog = NULL;
1002     }
1003
1004     bRerunMD     = (Flags & MD_RERUN);
1005     bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
1006     bTryUseGPU   = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
1007     /* Rerun execution time is dominated by I/O and pair search, so
1008      * GPUs are not very useful, plus they do not support more than
1009      * one energy group. Don't select them when they can't be used,
1010      * unless the user requested it, then fatal_error is called later.
1011      *
1012      * TODO it would be nice to notify the user that if this check
1013      * causes GPUs not to be used that this is what is happening, and
1014      * why, but that will be easier to do after some future
1015      * cleanup. */
1016     bCantUseGPU = bRerunMD && (inputrec->opts.ngener > 1);
1017     bTryUseGPU  = bTryUseGPU && !(bCantUseGPU && !bForceUseGPU);
1018
1019     /* Detect hardware, gather information. This is an operation that is
1020      * global for this process (MPI rank). */
1021     hwinfo = gmx_detect_hardware(fplog, cr, bTryUseGPU);
1022
1023
1024     snew(state, 1);
1025     if (SIMMASTER(cr))
1026     {
1027         /* Read (nearly) all data required for the simulation */
1028         read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, state, NULL, mtop);
1029
1030         if (inputrec->cutoff_scheme == ecutsVERLET)
1031         {
1032             /* Here the master rank decides if all ranks will use GPUs */
1033             bUseGPU = (hwinfo->gpu_info.ncuda_dev_compatible > 0 ||
1034                        getenv("GMX_EMULATE_GPU") != NULL);
1035
1036             /* TODO add GPU kernels for this and replace this check by:
1037              * (bUseGPU && (ir->vdwtype == evdwPME &&
1038              *               ir->ljpme_combination_rule == eljpmeLB))
1039              * update the message text and the content of nbnxn_acceleration_supported.
1040              */
1041             if (bUseGPU &&
1042                 !nbnxn_acceleration_supported(fplog, cr, inputrec, bUseGPU))
1043             {
1044                 /* Fallback message printed by nbnxn_acceleration_supported */
1045                 if (bForceUseGPU)
1046                 {
1047                     gmx_fatal(FARGS, "GPU acceleration requested, but not supported with the given input settings");
1048                 }
1049                 bUseGPU = FALSE;
1050             }
1051
1052             prepare_verlet_scheme(fplog, cr,
1053                                   inputrec, nstlist_cmdline, mtop, state->box,
1054                                   bUseGPU);
1055         }
1056         else
1057         {
1058             if (nstlist_cmdline > 0)
1059             {
1060                 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
1061             }
1062
1063             if (hwinfo->gpu_info.ncuda_dev_compatible > 0)
1064             {
1065                 md_print_warn(cr, fplog,
1066                               "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
1067                               "      To use a GPU, set the mdp option: cutoff-scheme = Verlet\n");
1068             }
1069
1070             if (bForceUseGPU)
1071             {
1072                 gmx_fatal(FARGS, "GPU requested, but can't be used without cutoff-scheme=Verlet");
1073             }
1074
1075 #ifdef GMX_TARGET_BGQ
1076             md_print_warn(cr, fplog,
1077                           "NOTE: There is no SIMD implementation of the group scheme kernels on\n"
1078                           "      BlueGene/Q. You will observe better performance from using the\n"
1079                           "      Verlet cut-off scheme.\n");
1080 #endif
1081         }
1082
1083         if (inputrec->eI == eiSD2)
1084         {
1085             md_print_warn(cr, fplog, "The stochastic dynamics integrator %s is deprecated, since\n"
1086                           "it is slower than integrator %s and is slightly less accurate\n"
1087                           "with constraints. Use the %s integrator.",
1088                           ei_names[inputrec->eI], ei_names[eiSD1], ei_names[eiSD1]);
1089         }
1090     }
1091
1092     /* Check and update the hardware options for internal consistency */
1093     check_and_update_hw_opt_1(hw_opt, SIMMASTER(cr));
1094
1095     /* Early check for externally set process affinity. */
1096     gmx_check_thread_affinity_set(fplog, cr,
1097                                   hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1098     if (SIMMASTER(cr))
1099     {
1100
1101 #ifdef GMX_THREAD_MPI
1102         if (cr->npmenodes > 0 && hw_opt->nthreads_tmpi <= 0)
1103         {
1104             gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
1105         }
1106 #endif
1107
1108         if (hw_opt->nthreads_omp_pme != hw_opt->nthreads_omp &&
1109             cr->npmenodes <= 0)
1110         {
1111             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");
1112         }
1113     }
1114
1115 #ifdef GMX_THREAD_MPI
1116     if (SIMMASTER(cr))
1117     {
1118         /* Since the master knows the cut-off scheme, update hw_opt for this.
1119          * This is done later for normal MPI and also once more with tMPI
1120          * for all tMPI ranks.
1121          */
1122         check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1123
1124         /* NOW the threads will be started: */
1125         hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
1126                                                  hw_opt,
1127                                                  inputrec, mtop,
1128                                                  cr, fplog);
1129         if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp <= 0)
1130         {
1131             hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
1132         }
1133
1134         if (hw_opt->nthreads_tmpi > 1)
1135         {
1136             t_commrec *cr_old       = cr;
1137             /* now start the threads. */
1138             cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
1139                                         oenv, bVerbose, bCompact, nstglobalcomm,
1140                                         ddxyz, dd_node_order, rdd, rconstr,
1141                                         dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
1142                                         nbpu_opt, nstlist_cmdline,
1143                                         nsteps_cmdline, nstepout, resetstep, nmultisim,
1144                                         repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
1145                                         cpt_period, max_hours, deviceOptions,
1146                                         Flags);
1147             /* the main thread continues here with a new cr. We don't deallocate
1148                the old cr because other threads may still be reading it. */
1149             if (cr == NULL)
1150             {
1151                 gmx_comm("Failed to spawn threads");
1152             }
1153         }
1154     }
1155 #endif
1156     /* END OF CAUTION: cr is now reliable */
1157
1158     /* g_membed initialisation *
1159      * Because we change the mtop, init_membed is called before the init_parallel *
1160      * (in case we ever want to make it run in parallel) */
1161     if (opt2bSet("-membed", nfile, fnm))
1162     {
1163         if (MASTER(cr))
1164         {
1165             fprintf(stderr, "Initializing membed");
1166         }
1167         membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
1168     }
1169
1170     if (PAR(cr))
1171     {
1172         /* now broadcast everything to the non-master nodes/threads: */
1173         init_parallel(cr, inputrec, mtop);
1174     }
1175     if (fplog != NULL)
1176     {
1177         pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
1178     }
1179
1180     /* now make sure the state is initialized and propagated */
1181     set_state_entries(state, inputrec);
1182
1183     /* A parallel command line option consistency check that we can
1184        only do after any threads have started. */
1185     if (!PAR(cr) &&
1186         (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
1187     {
1188         gmx_fatal(FARGS,
1189                   "The -dd or -npme option request a parallel simulation, "
1190 #ifndef GMX_MPI
1191                   "but %s was compiled without threads or MPI enabled"
1192 #else
1193 #ifdef GMX_THREAD_MPI
1194                   "but the number of threads (option -nt) is 1"
1195 #else
1196                   "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
1197 #endif
1198 #endif
1199                   , output_env_get_program_display_name(oenv)
1200                   );
1201     }
1202
1203     if (bRerunMD &&
1204         (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1205     {
1206         gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
1207     }
1208
1209     if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
1210     {
1211         gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
1212     }
1213
1214     if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
1215     {
1216         if (cr->npmenodes > 0)
1217         {
1218             gmx_fatal_collective(FARGS, cr, NULL,
1219                                  "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
1220         }
1221
1222         cr->npmenodes = 0;
1223     }
1224
1225 #ifdef GMX_FAHCORE
1226     if (MASTER(cr))
1227     {
1228         fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1229     }
1230 #endif
1231
1232     /* NMR restraints must be initialized before load_checkpoint,
1233      * since with time averaging the history is added to t_state.
1234      * For proper consistency check we therefore need to extend
1235      * t_state here.
1236      * So the PME-only nodes (if present) will also initialize
1237      * the distance restraints.
1238      */
1239     snew(fcd, 1);
1240
1241     /* This needs to be called before read_checkpoint to extend the state */
1242     init_disres(fplog, mtop, inputrec, cr, fcd, state, repl_ex_nst > 0);
1243
1244     init_orires(fplog, mtop, state->x, inputrec, cr, &(fcd->orires),
1245                 state);
1246
1247     if (DEFORM(*inputrec))
1248     {
1249         /* Store the deform reference box before reading the checkpoint */
1250         if (SIMMASTER(cr))
1251         {
1252             copy_mat(state->box, box);
1253         }
1254         if (PAR(cr))
1255         {
1256             gmx_bcast(sizeof(box), box, cr);
1257         }
1258         /* Because we do not have the update struct available yet
1259          * in which the reference values should be stored,
1260          * we store them temporarily in static variables.
1261          * This should be thread safe, since they are only written once
1262          * and with identical values.
1263          */
1264         tMPI_Thread_mutex_lock(&deform_init_box_mutex);
1265         deform_init_init_step_tpx = inputrec->init_step;
1266         copy_mat(box, deform_init_box_tpx);
1267         tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
1268     }
1269
1270     if (opt2bSet("-cpi", nfile, fnm))
1271     {
1272         /* Check if checkpoint file exists before doing continuation.
1273          * This way we can use identical input options for the first and subsequent runs...
1274          */
1275         if (gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr) )
1276         {
1277             load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
1278                             cr, ddxyz,
1279                             inputrec, state, &bReadEkin,
1280                             (Flags & MD_APPENDFILES),
1281                             (Flags & MD_APPENDFILESSET));
1282
1283             if (bReadEkin)
1284             {
1285                 Flags |= MD_READ_EKIN;
1286             }
1287         }
1288     }
1289
1290     if (MASTER(cr) && (Flags & MD_APPENDFILES))
1291     {
1292         gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
1293                      Flags, &fplog);
1294     }
1295
1296     /* override nsteps with value from cmdline */
1297     override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
1298
1299     if (SIMMASTER(cr))
1300     {
1301         copy_mat(state->box, box);
1302     }
1303
1304     if (PAR(cr))
1305     {
1306         gmx_bcast(sizeof(box), box, cr);
1307     }
1308
1309     /* Essential dynamics */
1310     if (opt2bSet("-ei", nfile, fnm))
1311     {
1312         /* Open input and output files, allocate space for ED data structure */
1313         ed = ed_open(mtop->natoms, &state->edsamstate, nfile, fnm, Flags, oenv, cr);
1314     }
1315
1316     if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
1317                      inputrec->eI == eiNM))
1318     {
1319         cr->dd = init_domain_decomposition(fplog, cr, Flags, ddxyz, rdd, rconstr,
1320                                            dddlb_opt, dlb_scale,
1321                                            ddcsx, ddcsy, ddcsz,
1322                                            mtop, inputrec,
1323                                            box, state->x,
1324                                            &ddbox, &npme_major, &npme_minor);
1325
1326         make_dd_communicators(fplog, cr, dd_node_order);
1327
1328         /* Set overallocation to avoid frequent reallocation of arrays */
1329         set_over_alloc_dd(TRUE);
1330     }
1331     else
1332     {
1333         /* PME, if used, is done on all nodes with 1D decomposition */
1334         cr->npmenodes = 0;
1335         cr->duty      = (DUTY_PP | DUTY_PME);
1336         npme_major    = 1;
1337         npme_minor    = 1;
1338
1339         if (inputrec->ePBC == epbcSCREW)
1340         {
1341             gmx_fatal(FARGS,
1342                       "pbc=%s is only implemented with domain decomposition",
1343                       epbc_names[inputrec->ePBC]);
1344         }
1345     }
1346
1347     if (PAR(cr))
1348     {
1349         /* After possible communicator splitting in make_dd_communicators.
1350          * we can set up the intra/inter node communication.
1351          */
1352         gmx_setup_nodecomm(fplog, cr);
1353     }
1354
1355     /* Initialize per-physical-node MPI process/thread ID and counters. */
1356     gmx_init_intranode_counters(cr);
1357 #ifdef GMX_MPI
1358     if (MULTISIM(cr))
1359     {
1360         md_print_info(cr, fplog,
1361                       "This is simulation %d out of %d running as a composite Gromacs\n"
1362                       "multi-simulation job. Setup for this simulation:\n\n",
1363                       cr->ms->sim, cr->ms->nsim);
1364     }
1365     md_print_info(cr, fplog, "Using %d MPI %s\n",
1366                   cr->nnodes,
1367 #ifdef GMX_THREAD_MPI
1368                   cr->nnodes == 1 ? "thread" : "threads"
1369 #else
1370                   cr->nnodes == 1 ? "process" : "processes"
1371 #endif
1372                   );
1373     fflush(stderr);
1374 #endif
1375
1376     /* Check and update hw_opt for the cut-off scheme */
1377     check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1378
1379     gmx_omp_nthreads_init(fplog, cr,
1380                           hwinfo->nthreads_hw_avail,
1381                           hw_opt->nthreads_omp,
1382                           hw_opt->nthreads_omp_pme,
1383                           (cr->duty & DUTY_PP) == 0,
1384                           inputrec->cutoff_scheme == ecutsVERLET);
1385
1386 #ifndef NDEBUG
1387     if (integrator[inputrec->eI].func != do_tpi &&
1388         inputrec->cutoff_scheme == ecutsVERLET)
1389     {
1390         gmx_feenableexcept();
1391     }
1392 #endif
1393     if (PAR(cr))
1394     {
1395         /* The master rank decided on the use of GPUs,
1396          * broadcast this information to all ranks.
1397          */
1398         gmx_bcast_sim(sizeof(bUseGPU), &bUseGPU, cr);
1399     }
1400
1401     if (bUseGPU)
1402     {
1403         if (cr->npmenodes == -1)
1404         {
1405             /* Don't automatically use PME-only nodes with GPUs */
1406             cr->npmenodes = 0;
1407         }
1408
1409         /* Select GPU id's to use */
1410         gmx_select_gpu_ids(fplog, cr, &hwinfo->gpu_info, bForceUseGPU,
1411                            &hw_opt->gpu_opt);
1412     }
1413     else
1414     {
1415         /* Ignore (potentially) manually selected GPUs */
1416         hw_opt->gpu_opt.ncuda_dev_use = 0;
1417     }
1418
1419     /* check consistency across ranks of things like SIMD
1420      * support and number of GPUs selected */
1421     gmx_check_hw_runconf_consistency(fplog, hwinfo, cr, hw_opt, bUseGPU);
1422
1423     if (DOMAINDECOMP(cr))
1424     {
1425         /* When we share GPUs over ranks, we need to know this for the DLB */
1426         dd_setup_dlb_resource_sharing(cr, hwinfo, hw_opt);
1427     }
1428
1429     /* getting number of PP/PME threads
1430        PME: env variable should be read only on one node to make sure it is
1431        identical everywhere;
1432      */
1433     /* TODO nthreads_pp is only used for pinning threads.
1434      * This is a temporary solution until we have a hw topology library.
1435      */
1436     nthreads_pp  = gmx_omp_nthreads_get(emntNonbonded);
1437     nthreads_pme = gmx_omp_nthreads_get(emntPME);
1438
1439     wcycle = wallcycle_init(fplog, resetstep, cr, nthreads_pp, nthreads_pme);
1440
1441     if (PAR(cr))
1442     {
1443         /* Master synchronizes its value of reset_counters with all nodes
1444          * including PME only nodes */
1445         reset_counters = wcycle_get_reset_counters(wcycle);
1446         gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1447         wcycle_set_reset_counters(wcycle, reset_counters);
1448     }
1449
1450     snew(nrnb, 1);
1451     if (cr->duty & DUTY_PP)
1452     {
1453         bcast_state(cr, state);
1454
1455         /* Initiate forcerecord */
1456         fr          = mk_forcerec();
1457         fr->hwinfo  = hwinfo;
1458         fr->gpu_opt = &hw_opt->gpu_opt;
1459         init_forcerec(fplog, oenv, fr, fcd, inputrec, mtop, cr, box,
1460                       opt2fn("-table", nfile, fnm),
1461                       opt2fn("-tabletf", nfile, fnm),
1462                       opt2fn("-tablep", nfile, fnm),
1463                       opt2fn("-tableb", nfile, fnm),
1464                       nbpu_opt,
1465                       FALSE,
1466                       pforce);
1467
1468         /* version for PCA_NOT_READ_NODE (see md.c) */
1469         /*init_forcerec(fplog,fr,fcd,inputrec,mtop,cr,box,FALSE,
1470            "nofile","nofile","nofile","nofile",FALSE,pforce);
1471          */
1472
1473         /* Initialize QM-MM */
1474         if (fr->bQMMM)
1475         {
1476             init_QMMMrec(cr, mtop, inputrec, fr);
1477         }
1478
1479         /* Initialize the mdatoms structure.
1480          * mdatoms is not filled with atom data,
1481          * as this can not be done now with domain decomposition.
1482          */
1483         mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1484
1485         /* Initialize the virtual site communication */
1486         vsite = init_vsite(mtop, cr, FALSE);
1487
1488         calc_shifts(box, fr->shift_vec);
1489
1490         /* With periodic molecules the charge groups should be whole at start up
1491          * and the virtual sites should not be far from their proper positions.
1492          */
1493         if (!inputrec->bContinuation && MASTER(cr) &&
1494             !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1495         {
1496             /* Make molecules whole at start of run */
1497             if (fr->ePBC != epbcNONE)
1498             {
1499                 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, state->x);
1500             }
1501             if (vsite)
1502             {
1503                 /* Correct initial vsite positions are required
1504                  * for the initial distribution in the domain decomposition
1505                  * and for the initial shell prediction.
1506                  */
1507                 construct_vsites_mtop(vsite, mtop, state->x);
1508             }
1509         }
1510
1511         if (EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype))
1512         {
1513             ewaldcoeff_q  = fr->ewaldcoeff_q;
1514             ewaldcoeff_lj = fr->ewaldcoeff_lj;
1515             pmedata       = &fr->pmedata;
1516         }
1517         else
1518         {
1519             pmedata = NULL;
1520         }
1521     }
1522     else
1523     {
1524         /* This is a PME only node */
1525
1526         /* We don't need the state */
1527         done_state(state);
1528
1529         ewaldcoeff_q  = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1530         ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1531         snew(pmedata, 1);
1532     }
1533
1534     if (hw_opt->thread_affinity != threadaffOFF)
1535     {
1536         /* Before setting affinity, check whether the affinity has changed
1537          * - which indicates that probably the OpenMP library has changed it
1538          * since we first checked).
1539          */
1540         gmx_check_thread_affinity_set(fplog, cr,
1541                                       hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1542
1543         /* Set the CPU affinity */
1544         gmx_set_thread_affinity(fplog, cr, hw_opt, hwinfo);
1545     }
1546
1547     /* Initiate PME if necessary,
1548      * either on all nodes or on dedicated PME nodes only. */
1549     if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1550     {
1551         if (mdatoms)
1552         {
1553             nChargePerturbed = mdatoms->nChargePerturbed;
1554             if (EVDW_PME(inputrec->vdwtype))
1555             {
1556                 nTypePerturbed   = mdatoms->nTypePerturbed;
1557             }
1558         }
1559         if (cr->npmenodes > 0)
1560         {
1561             /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1562             gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1563             gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1564         }
1565
1566         if (cr->duty & DUTY_PME)
1567         {
1568             status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1569                                   mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1570                                   (Flags & MD_REPRODUCIBLE), nthreads_pme);
1571             if (status != 0)
1572             {
1573                 gmx_fatal(FARGS, "Error %d initializing PME", status);
1574             }
1575         }
1576     }
1577
1578
1579     if (integrator[inputrec->eI].func == do_md)
1580     {
1581         /* Turn on signal handling on all nodes */
1582         /*
1583          * (A user signal from the PME nodes (if any)
1584          * is communicated to the PP nodes.
1585          */
1586         signal_handler_install();
1587     }
1588
1589     if (cr->duty & DUTY_PP)
1590     {
1591         /* Assumes uniform use of the number of OpenMP threads */
1592         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1593
1594         if (inputrec->bPull)
1595         {
1596             /* Initialize pull code */
1597             init_pull(fplog, inputrec, nfile, fnm, mtop, cr, oenv, inputrec->fepvals->init_lambda,
1598                       EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1599         }
1600
1601         if (inputrec->bRot)
1602         {
1603             /* Initialize enforced rotation code */
1604             init_rot(fplog, inputrec, nfile, fnm, cr, state->x, box, mtop, oenv,
1605                      bVerbose, Flags);
1606         }
1607
1608         if (inputrec->eSwapCoords != eswapNO)
1609         {
1610             /* Initialize ion swapping code */
1611             init_swapcoords(fplog, bVerbose, inputrec, opt2fn_master("-swap", nfile, fnm, cr),
1612                             mtop, state->x, state->box, &state->swapstate, cr, oenv, Flags);
1613         }
1614
1615         constr = init_constraints(fplog, mtop, inputrec, ed, state, cr);
1616
1617         if (DOMAINDECOMP(cr))
1618         {
1619             GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1620             dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1621                             Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1622
1623             set_dd_parameters(fplog, cr->dd, dlb_scale, inputrec, &ddbox);
1624
1625             setup_dd_grid(fplog, cr->dd);
1626         }
1627
1628         /* Now do whatever the user wants us to do (how flexible...) */
1629         integrator[inputrec->eI].func(fplog, cr, nfile, fnm,
1630                                       oenv, bVerbose, bCompact,
1631                                       nstglobalcomm,
1632                                       vsite, constr,
1633                                       nstepout, inputrec, mtop,
1634                                       fcd, state,
1635                                       mdatoms, nrnb, wcycle, ed, fr,
1636                                       repl_ex_nst, repl_ex_nex, repl_ex_seed,
1637                                       membed,
1638                                       cpt_period, max_hours,
1639                                       deviceOptions,
1640                                       imdport,
1641                                       Flags,
1642                                       walltime_accounting);
1643
1644         if (inputrec->bPull)
1645         {
1646             finish_pull(inputrec->pull);
1647         }
1648
1649         if (inputrec->bRot)
1650         {
1651             finish_rot(inputrec->rot);
1652         }
1653
1654     }
1655     else
1656     {
1657         GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1658         /* do PME only */
1659         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1660         gmx_pmeonly(*pmedata, cr, nrnb, wcycle, walltime_accounting, ewaldcoeff_q, ewaldcoeff_lj, inputrec);
1661     }
1662
1663     wallcycle_stop(wcycle, ewcRUN);
1664
1665     /* Finish up, write some stuff
1666      * if rerunMD, don't write last frame again
1667      */
1668     finish_run(fplog, cr,
1669                inputrec, nrnb, wcycle, walltime_accounting,
1670                fr ? fr->nbv : NULL,
1671                EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1672
1673
1674     /* Free GPU memory and context */
1675     free_gpu_resources(fr, cr, &hwinfo->gpu_info, fr ? fr->gpu_opt : NULL);
1676
1677     if (opt2bSet("-membed", nfile, fnm))
1678     {
1679         sfree(membed);
1680     }
1681
1682     gmx_hardware_info_free(hwinfo);
1683
1684     /* Does what it says */
1685     print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1686     walltime_accounting_destroy(walltime_accounting);
1687
1688     /* Close logfile already here if we were appending to it */
1689     if (MASTER(cr) && (Flags & MD_APPENDFILES))
1690     {
1691         gmx_log_close(fplog);
1692     }
1693
1694     rc = (int)gmx_get_stop_condition();
1695
1696     done_ed(&ed);
1697
1698 #ifdef GMX_THREAD_MPI
1699     /* we need to join all threads. The sub-threads join when they
1700        exit this function, but the master thread needs to be told to
1701        wait for that. */
1702     if (PAR(cr) && MASTER(cr))
1703     {
1704         tMPI_Finalize();
1705     }
1706 #endif
1707
1708     return rc;
1709 }