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