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