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