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