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