mdrun without OpenMP with thread-MPI now uses all cores
[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(const 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         /* Dims/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                                   const gmx_hw_info_t *hwinfo,
632                                   t_commrec           *cr,
633                                   const char          *nbpu_opt,
634                                   t_inputrec          *ir,
635                                   const gmx_mtop_t    *mtop,
636                                   matrix               box,
637                                   gmx_bool            *bUseGPU)
638 {
639     /* Here we only check for GPU usage on the MPI master process,
640      * as here we don't know how many GPUs we will use yet.
641      * We check for a GPU on all processes later.
642      */
643     *bUseGPU = hwinfo->bCanUseGPU || (getenv("GMX_EMULATE_GPU") != NULL);
644
645     if (ir->verletbuf_drift > 0)
646     {
647         /* Update the Verlet buffer size for the current run setup */
648         verletbuf_list_setup_t ls;
649         real                   rlist_new;
650
651         /* Here we assume CPU acceleration is on. But as currently
652          * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
653          * and 4x2 gives a larger buffer than 4x4, this is ok.
654          */
655         verletbuf_get_list_setup(*bUseGPU, &ls);
656
657         calc_verlet_buffer_size(mtop, det(box), ir,
658                                 ir->verletbuf_drift, &ls,
659                                 NULL, &rlist_new);
660         if (rlist_new != ir->rlist)
661         {
662             if (fplog != NULL)
663             {
664                 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
665                         ir->rlist, rlist_new,
666                         ls.cluster_size_i, ls.cluster_size_j);
667             }
668             ir->rlist     = rlist_new;
669             ir->rlistlong = rlist_new;
670         }
671     }
672
673     /* With GPU or emulation we should check nstlist for performance */
674     if ((EI_DYNAMICS(ir->eI) &&
675          *bUseGPU &&
676          ir->nstlist < NSTLIST_GPU_ENOUGH) ||
677         getenv(NSTLIST_ENVVAR) != NULL)
678     {
679         /* Choose a better nstlist */
680         increase_nstlist(fplog, cr, ir, mtop, box);
681     }
682 }
683
684 static void convert_to_verlet_scheme(FILE *fplog,
685                                      t_inputrec *ir,
686                                      gmx_mtop_t *mtop, real box_vol)
687 {
688     char *conv_mesg = "Converting input file with group cut-off scheme to the Verlet cut-off scheme";
689
690     md_print_warn(NULL, fplog, "%s\n", conv_mesg);
691
692     ir->cutoff_scheme   = ecutsVERLET;
693     ir->verletbuf_drift = 0.005;
694
695     if (ir->rcoulomb != ir->rvdw)
696     {
697         gmx_fatal(FARGS, "The VdW and Coulomb cut-offs are different, whereas the Verlet scheme only supports equal cut-offs");
698     }
699
700     if (ir->vdwtype == evdwUSER || EEL_USER(ir->coulombtype))
701     {
702         gmx_fatal(FARGS, "User non-bonded potentials are not (yet) supported with the Verlet scheme");
703     }
704     else if (EVDW_SWITCHED(ir->vdwtype) || EEL_SWITCHED(ir->coulombtype))
705     {
706         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");
707
708         if (EVDW_SWITCHED(ir->vdwtype))
709         {
710             ir->vdwtype = evdwCUT;
711         }
712         if (EEL_SWITCHED(ir->coulombtype))
713         {
714             if (EEL_FULL(ir->coulombtype))
715             {
716                 /* With full electrostatic only PME can be switched */
717                 ir->coulombtype = eelPME;
718             }
719             else
720             {
721                 md_print_warn(NULL, fplog, "NOTE: Replacing %s electrostatics with reaction-field with epsilon-rf=inf\n", eel_names[ir->coulombtype]);
722                 ir->coulombtype = eelRF;
723                 ir->epsilon_rf  = 0.0;
724             }
725         }
726
727         /* We set the target energy drift to a small number.
728          * Note that this is only for testing. For production the user
729          * should think about this and set the mdp options.
730          */
731         ir->verletbuf_drift = 1e-4;
732     }
733
734     if (inputrec2nboundeddim(ir) != 3)
735     {
736         gmx_fatal(FARGS, "Can only convert old tpr files to the Verlet cut-off scheme with 3D pbc");
737     }
738
739     if (ir->efep != efepNO || ir->implicit_solvent != eisNO)
740     {
741         gmx_fatal(FARGS, "Will not convert old tpr files to the Verlet cut-off scheme with free-energy calculations or implicit solvent");
742     }
743
744     if (EI_DYNAMICS(ir->eI) && !(EI_MD(ir->eI) && ir->etc == etcNO))
745     {
746         verletbuf_list_setup_t ls;
747
748         verletbuf_get_list_setup(FALSE, &ls);
749         calc_verlet_buffer_size(mtop, box_vol, ir, ir->verletbuf_drift, &ls,
750                                 NULL, &ir->rlist);
751     }
752     else
753     {
754         ir->verletbuf_drift = -1;
755         ir->rlist           = 1.05*max(ir->rvdw, ir->rcoulomb);
756     }
757
758     gmx_mtop_remove_chargegroups(mtop);
759 }
760
761 static void check_and_update_hw_opt(gmx_hw_opt_t *hw_opt,
762                                     int           cutoff_scheme,
763                                     gmx_bool      bIsSimMaster)
764 {
765     gmx_omp_nthreads_read_env(&hw_opt->nthreads_omp, bIsSimMaster);
766
767 #ifndef GMX_THREAD_MPI
768     if (hw_opt->nthreads_tot > 0)
769     {
770         gmx_fatal(FARGS, "Setting the total number of threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
771     }
772     if (hw_opt->nthreads_tmpi > 0)
773     {
774         gmx_fatal(FARGS, "Setting the number of thread-MPI threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
775     }
776 #endif
777
778 #ifndef GMX_OPENMP
779     if (hw_opt->nthreads_omp > 1)
780     {
781         gmx_fatal(FARGS, "More than 1 OpenMP thread requested, but Gromacs was compiled without OpenMP support");
782     }
783     hw_opt->nthreads_omp = 1;
784 #endif
785
786     if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp_pme <= 0)
787     {
788         /* We have the same number of OpenMP threads for PP and PME processes,
789          * thus we can perform several consistency checks.
790          */
791         if (hw_opt->nthreads_tmpi > 0 &&
792             hw_opt->nthreads_omp > 0 &&
793             hw_opt->nthreads_tot != hw_opt->nthreads_tmpi*hw_opt->nthreads_omp)
794         {
795             gmx_fatal(FARGS, "The total number of threads requested (%d) does not match the thread-MPI threads (%d) times the OpenMP threads (%d) requested",
796                       hw_opt->nthreads_tot, hw_opt->nthreads_tmpi, hw_opt->nthreads_omp);
797         }
798
799         if (hw_opt->nthreads_tmpi > 0 &&
800             hw_opt->nthreads_tot % hw_opt->nthreads_tmpi != 0)
801         {
802             gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of thread-MPI threads requested (%d)",
803                       hw_opt->nthreads_tot, hw_opt->nthreads_tmpi);
804         }
805
806         if (hw_opt->nthreads_omp > 0 &&
807             hw_opt->nthreads_tot % hw_opt->nthreads_omp != 0)
808         {
809             gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of OpenMP threads requested (%d)",
810                       hw_opt->nthreads_tot, hw_opt->nthreads_omp);
811         }
812
813         if (hw_opt->nthreads_tmpi > 0 &&
814             hw_opt->nthreads_omp <= 0)
815         {
816             hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
817         }
818     }
819
820 #ifndef GMX_OPENMP
821     if (hw_opt->nthreads_omp > 1)
822     {
823         gmx_fatal(FARGS, "OpenMP threads are requested, but Gromacs was compiled without OpenMP support");
824     }
825 #endif
826
827     if (cutoff_scheme == ecutsGROUP)
828     {
829         /* We only have OpenMP support for PME only nodes */
830         if (hw_opt->nthreads_omp > 1)
831         {
832             gmx_fatal(FARGS, "OpenMP threads have been requested with cut-off scheme %s, but these are only supported with cut-off scheme %s",
833                       ecutscheme_names[cutoff_scheme],
834                       ecutscheme_names[ecutsVERLET]);
835         }
836         hw_opt->nthreads_omp = 1;
837     }
838
839     if (hw_opt->nthreads_omp_pme > 0 && hw_opt->nthreads_omp <= 0)
840     {
841         gmx_fatal(FARGS, "You need to specify -ntomp in addition to -ntomp_pme");
842     }
843
844     if (hw_opt->nthreads_tot == 1)
845     {
846         hw_opt->nthreads_tmpi = 1;
847
848         if (hw_opt->nthreads_omp > 1)
849         {
850             gmx_fatal(FARGS, "You requested %d OpenMP threads with %d total threads",
851                       hw_opt->nthreads_tmpi, hw_opt->nthreads_tot);
852         }
853         hw_opt->nthreads_omp = 1;
854     }
855
856     if (hw_opt->nthreads_omp_pme <= 0 && hw_opt->nthreads_omp > 0)
857     {
858         hw_opt->nthreads_omp_pme = hw_opt->nthreads_omp;
859     }
860
861     if (debug)
862     {
863         fprintf(debug, "hw_opt: nt %d ntmpi %d ntomp %d ntomp_pme %d gpu_id '%s'\n",
864                 hw_opt->nthreads_tot,
865                 hw_opt->nthreads_tmpi,
866                 hw_opt->nthreads_omp,
867                 hw_opt->nthreads_omp_pme,
868                 hw_opt->gpu_id != NULL ? hw_opt->gpu_id : "");
869
870     }
871 }
872
873
874 /* Override the value in inputrec with value passed on the command line (if any) */
875 static void override_nsteps_cmdline(FILE            *fplog,
876                                     gmx_large_int_t  nsteps_cmdline,
877                                     t_inputrec      *ir,
878                                     const t_commrec *cr)
879 {
880     char sbuf[STEPSTRSIZE];
881
882     assert(ir);
883     assert(cr);
884
885     /* override with anything else than the default -2 */
886     if (nsteps_cmdline > -2)
887     {
888         char stmp[STRLEN];
889
890         ir->nsteps = nsteps_cmdline;
891         if (EI_DYNAMICS(ir->eI))
892         {
893             sprintf(stmp, "Overriding nsteps with value passed on the command line: %s steps, %.3f ps",
894                     gmx_step_str(nsteps_cmdline, sbuf),
895                     nsteps_cmdline*ir->delta_t);
896         }
897         else
898         {
899             sprintf(stmp, "Overriding nsteps with value passed on the command line: %s steps",
900                     gmx_step_str(nsteps_cmdline, sbuf));
901         }
902
903         md_print_warn(cr, fplog, "%s\n", stmp);
904     }
905 }
906
907 /* Data structure set by SIMMASTER which needs to be passed to all nodes
908  * before the other nodes have read the tpx file and called gmx_detect_hardware.
909  */
910 typedef struct {
911     int      cutoff_scheme; /* The cutoff scheme from inputrec_t */
912     gmx_bool bUseGPU;       /* Use GPU or GPU emulation          */
913 } master_inf_t;
914
915 int mdrunner(gmx_hw_opt_t *hw_opt,
916              FILE *fplog, t_commrec *cr, int nfile,
917              const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
918              gmx_bool bCompact, int nstglobalcomm,
919              ivec ddxyz, int dd_node_order, real rdd, real rconstr,
920              const char *dddlb_opt, real dlb_scale,
921              const char *ddcsx, const char *ddcsy, const char *ddcsz,
922              const char *nbpu_opt,
923              gmx_large_int_t nsteps_cmdline, int nstepout, int resetstep,
924              int nmultisim, int repl_ex_nst, int repl_ex_nex,
925              int repl_ex_seed, real pforce, real cpt_period, real max_hours,
926              const char *deviceOptions, unsigned long Flags)
927 {
928     gmx_bool        bForceUseGPU, bTryUseGPU;
929     double          nodetime = 0, realtime;
930     t_inputrec     *inputrec;
931     t_state        *state = NULL;
932     matrix          box;
933     gmx_ddbox_t     ddbox = {0};
934     int             npme_major, npme_minor;
935     real            tmpr1, tmpr2;
936     t_nrnb         *nrnb;
937     gmx_mtop_t     *mtop       = NULL;
938     t_mdatoms      *mdatoms    = NULL;
939     t_forcerec     *fr         = NULL;
940     t_fcdata       *fcd        = NULL;
941     real            ewaldcoeff = 0;
942     gmx_pme_t      *pmedata    = NULL;
943     gmx_vsite_t    *vsite      = NULL;
944     gmx_constr_t    constr;
945     int             i, m, nChargePerturbed = -1, status, nalloc;
946     char           *gro;
947     gmx_wallcycle_t wcycle;
948     gmx_bool        bReadRNG, bReadEkin;
949     int             list;
950     gmx_runtime_t   runtime;
951     int             rc;
952     gmx_large_int_t reset_counters;
953     gmx_edsam_t     ed           = NULL;
954     t_commrec      *cr_old       = cr;
955     int             nthreads_pme = 1;
956     int             nthreads_pp  = 1;
957     gmx_membed_t    membed       = NULL;
958     gmx_hw_info_t  *hwinfo       = NULL;
959     master_inf_t    minf         = {-1, FALSE};
960
961     /* CAUTION: threads may be started later on in this function, so
962        cr doesn't reflect the final parallel state right now */
963     snew(inputrec, 1);
964     snew(mtop, 1);
965
966     if (Flags & MD_APPENDFILES)
967     {
968         fplog = NULL;
969     }
970
971     bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
972     bTryUseGPU   = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
973
974     /* Detect hardware, gather information. This is an operation that is
975      * global for this process (MPI rank). */
976     hwinfo = gmx_detect_hardware(fplog, cr,
977                                  bForceUseGPU, bTryUseGPU, hw_opt->gpu_id);
978
979
980     snew(state, 1);
981     if (SIMMASTER(cr))
982     {
983         /* Read (nearly) all data required for the simulation */
984         read_tpx_state(ftp2fn(efTPX, nfile, fnm), inputrec, state, NULL, mtop);
985
986         if (inputrec->cutoff_scheme != ecutsVERLET &&
987             ((Flags & MD_TESTVERLET) || getenv("GMX_VERLET_SCHEME") != NULL))
988         {
989             convert_to_verlet_scheme(fplog, inputrec, mtop, det(state->box));
990         }
991
992
993         minf.cutoff_scheme = inputrec->cutoff_scheme;
994         minf.bUseGPU       = FALSE;
995
996         if (inputrec->cutoff_scheme == ecutsVERLET)
997         {
998             prepare_verlet_scheme(fplog, hwinfo, cr, nbpu_opt,
999                                   inputrec, mtop, state->box,
1000                                   &minf.bUseGPU);
1001         }
1002         else if (hwinfo->bCanUseGPU)
1003         {
1004             md_print_warn(cr, fplog,
1005                           "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
1006                           "      To use a GPU, set the mdp option: cutoff-scheme = Verlet\n"
1007                           "      (for quick performance testing you can use the -testverlet option)\n");
1008
1009             if (bForceUseGPU)
1010             {
1011                 gmx_fatal(FARGS, "GPU requested, but can't be used without cutoff-scheme=Verlet");
1012             }
1013         }
1014     }
1015 #ifndef GMX_THREAD_MPI
1016     if (PAR(cr))
1017     {
1018         gmx_bcast_sim(sizeof(minf), &minf, cr);
1019     }
1020 #endif
1021     if (minf.bUseGPU && cr->npmenodes == -1)
1022     {
1023         /* Don't automatically use PME-only nodes with GPUs */
1024         cr->npmenodes = 0;
1025     }
1026
1027     /* Check for externally set OpenMP affinity and turn off internal
1028      * pinning if any is found. We need to do this check early to tell
1029      * thread-MPI whether it should do pinning when spawning threads.
1030      * TODO: the above no longer holds, we should move these checks down
1031      */
1032     gmx_omp_check_thread_affinity(fplog, cr, hw_opt);
1033
1034 #ifdef GMX_THREAD_MPI
1035     /* With thread-MPI inputrec is only set here on the master thread */
1036     if (SIMMASTER(cr))
1037 #endif
1038     {
1039         check_and_update_hw_opt(hw_opt, minf.cutoff_scheme, SIMMASTER(cr));
1040
1041 #ifdef GMX_THREAD_MPI
1042         /* Early check for externally set process affinity. Can't do over all
1043          * MPI processes because hwinfo is not available everywhere, but with
1044          * thread-MPI it's needed as pinning might get turned off which needs
1045          * to be known before starting thread-MPI. */
1046         gmx_check_thread_affinity_set(fplog,
1047                                       NULL,
1048                                       hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1049 #endif
1050
1051 #ifdef GMX_THREAD_MPI
1052         if (cr->npmenodes > 0 && hw_opt->nthreads_tmpi <= 0)
1053         {
1054             gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME nodes");
1055         }
1056 #endif
1057
1058         if (hw_opt->nthreads_omp_pme != hw_opt->nthreads_omp &&
1059             cr->npmenodes <= 0)
1060         {
1061             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");
1062         }
1063     }
1064
1065 #ifdef GMX_THREAD_MPI
1066     if (SIMMASTER(cr))
1067     {
1068         /* NOW the threads will be started: */
1069         hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
1070                                                  hw_opt,
1071                                                  inputrec, mtop,
1072                                                  cr, fplog);
1073         if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp <= 0)
1074         {
1075             hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
1076         }
1077
1078         if (hw_opt->nthreads_tmpi > 1)
1079         {
1080             /* now start the threads. */
1081             cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
1082                                         oenv, bVerbose, bCompact, nstglobalcomm,
1083                                         ddxyz, dd_node_order, rdd, rconstr,
1084                                         dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
1085                                         nbpu_opt,
1086                                         nsteps_cmdline, nstepout, resetstep, nmultisim,
1087                                         repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
1088                                         cpt_period, max_hours, deviceOptions,
1089                                         Flags);
1090             /* the main thread continues here with a new cr. We don't deallocate
1091                the old cr because other threads may still be reading it. */
1092             if (cr == NULL)
1093             {
1094                 gmx_comm("Failed to spawn threads");
1095             }
1096         }
1097     }
1098 #endif
1099     /* END OF CAUTION: cr is now reliable */
1100
1101     /* g_membed initialisation *
1102      * Because we change the mtop, init_membed is called before the init_parallel *
1103      * (in case we ever want to make it run in parallel) */
1104     if (opt2bSet("-membed", nfile, fnm))
1105     {
1106         if (MASTER(cr))
1107         {
1108             fprintf(stderr, "Initializing membed");
1109         }
1110         membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
1111     }
1112
1113     if (PAR(cr))
1114     {
1115         /* now broadcast everything to the non-master nodes/threads: */
1116         init_parallel(fplog, cr, inputrec, mtop);
1117
1118         /* This check needs to happen after get_nthreads_mpi() */
1119         if (inputrec->cutoff_scheme == ecutsVERLET && (Flags & MD_PARTDEC))
1120         {
1121             gmx_fatal_collective(FARGS, cr, NULL,
1122                                  "The Verlet cut-off scheme is not supported with particle decomposition.\n"
1123                                  "You can achieve the same effect as particle decomposition by running in parallel using only OpenMP threads.");
1124         }
1125     }
1126     if (fplog != NULL)
1127     {
1128         pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
1129     }
1130
1131     /* now make sure the state is initialized and propagated */
1132     set_state_entries(state, inputrec, cr->nnodes);
1133
1134     /* A parallel command line option consistency check that we can
1135        only do after any threads have started. */
1136     if (!PAR(cr) &&
1137         (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
1138     {
1139         gmx_fatal(FARGS,
1140                   "The -dd or -npme option request a parallel simulation, "
1141 #ifndef GMX_MPI
1142                   "but %s was compiled without threads or MPI enabled"
1143 #else
1144 #ifdef GMX_THREAD_MPI
1145                   "but the number of threads (option -nt) is 1"
1146 #else
1147                   "but %s was not started through mpirun/mpiexec or only one process was requested through mpirun/mpiexec"
1148 #endif
1149 #endif
1150                   , ShortProgram()
1151                   );
1152     }
1153
1154     if ((Flags & MD_RERUN) &&
1155         (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1156     {
1157         gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
1158     }
1159
1160     if (can_use_allvsall(inputrec, mtop, TRUE, cr, fplog) && PAR(cr))
1161     {
1162         /* Simple neighbour searching and (also?) all-vs-all loops
1163          * do not work with domain decomposition. */
1164         Flags |= MD_PARTDEC;
1165     }
1166
1167     if (!EEL_PME(inputrec->coulombtype) || (Flags & MD_PARTDEC))
1168     {
1169         if (cr->npmenodes > 0)
1170         {
1171             if (!EEL_PME(inputrec->coulombtype))
1172             {
1173                 gmx_fatal_collective(FARGS, cr, NULL,
1174                                      "PME nodes are requested, but the system does not use PME electrostatics");
1175             }
1176             if (Flags & MD_PARTDEC)
1177             {
1178                 gmx_fatal_collective(FARGS, cr, NULL,
1179                                      "PME nodes are requested, but particle decomposition does not support separate PME nodes");
1180             }
1181         }
1182
1183         cr->npmenodes = 0;
1184     }
1185
1186 #ifdef GMX_FAHCORE
1187     fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1188 #endif
1189
1190     /* NMR restraints must be initialized before load_checkpoint,
1191      * since with time averaging the history is added to t_state.
1192      * For proper consistency check we therefore need to extend
1193      * t_state here.
1194      * So the PME-only nodes (if present) will also initialize
1195      * the distance restraints.
1196      */
1197     snew(fcd, 1);
1198
1199     /* This needs to be called before read_checkpoint to extend the state */
1200     init_disres(fplog, mtop, inputrec, cr, Flags & MD_PARTDEC, fcd, state, repl_ex_nst > 0);
1201
1202     if (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0)
1203     {
1204         if (PAR(cr) && !(Flags & MD_PARTDEC))
1205         {
1206             gmx_fatal(FARGS, "Orientation restraints do not work (yet) with domain decomposition, use particle decomposition (mdrun option -pd)");
1207         }
1208         /* Orientation restraints */
1209         if (MASTER(cr))
1210         {
1211             init_orires(fplog, mtop, state->x, inputrec, cr->ms, &(fcd->orires),
1212                         state);
1213         }
1214     }
1215
1216     if (DEFORM(*inputrec))
1217     {
1218         /* Store the deform reference box before reading the checkpoint */
1219         if (SIMMASTER(cr))
1220         {
1221             copy_mat(state->box, box);
1222         }
1223         if (PAR(cr))
1224         {
1225             gmx_bcast(sizeof(box), box, cr);
1226         }
1227         /* Because we do not have the update struct available yet
1228          * in which the reference values should be stored,
1229          * we store them temporarily in static variables.
1230          * This should be thread safe, since they are only written once
1231          * and with identical values.
1232          */
1233 #ifdef GMX_THREAD_MPI
1234         tMPI_Thread_mutex_lock(&deform_init_box_mutex);
1235 #endif
1236         deform_init_init_step_tpx = inputrec->init_step;
1237         copy_mat(box, deform_init_box_tpx);
1238 #ifdef GMX_THREAD_MPI
1239         tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
1240 #endif
1241     }
1242
1243     if (opt2bSet("-cpi", nfile, fnm))
1244     {
1245         /* Check if checkpoint file exists before doing continuation.
1246          * This way we can use identical input options for the first and subsequent runs...
1247          */
1248         if (gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr) )
1249         {
1250             load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
1251                             cr, Flags & MD_PARTDEC, ddxyz,
1252                             inputrec, state, &bReadRNG, &bReadEkin,
1253                             (Flags & MD_APPENDFILES),
1254                             (Flags & MD_APPENDFILESSET));
1255
1256             if (bReadRNG)
1257             {
1258                 Flags |= MD_READ_RNG;
1259             }
1260             if (bReadEkin)
1261             {
1262                 Flags |= MD_READ_EKIN;
1263             }
1264         }
1265     }
1266
1267     if (((MASTER(cr) || (Flags & MD_SEPPOT)) && (Flags & MD_APPENDFILES))
1268 #ifdef GMX_THREAD_MPI
1269         /* With thread MPI only the master node/thread exists in mdrun.c,
1270          * therefore non-master nodes need to open the "seppot" log file here.
1271          */
1272         || (!MASTER(cr) && (Flags & MD_SEPPOT))
1273 #endif
1274         )
1275     {
1276         gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr, !(Flags & MD_SEPPOT),
1277                      Flags, &fplog);
1278     }
1279
1280     /* override nsteps with value from cmdline */
1281     override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
1282
1283     if (SIMMASTER(cr))
1284     {
1285         copy_mat(state->box, box);
1286     }
1287
1288     if (PAR(cr))
1289     {
1290         gmx_bcast(sizeof(box), box, cr);
1291     }
1292
1293     /* Essential dynamics */
1294     if (opt2bSet("-ei", nfile, fnm))
1295     {
1296         /* Open input and output files, allocate space for ED data structure */
1297         ed = ed_open(mtop->natoms, &state->edsamstate, nfile, fnm, Flags, oenv, cr);
1298     }
1299
1300     if (PAR(cr) && !((Flags & MD_PARTDEC) ||
1301                      EI_TPI(inputrec->eI) ||
1302                      inputrec->eI == eiNM))
1303     {
1304         cr->dd = init_domain_decomposition(fplog, cr, Flags, ddxyz, rdd, rconstr,
1305                                            dddlb_opt, dlb_scale,
1306                                            ddcsx, ddcsy, ddcsz,
1307                                            mtop, inputrec,
1308                                            box, state->x,
1309                                            &ddbox, &npme_major, &npme_minor);
1310
1311         make_dd_communicators(fplog, cr, dd_node_order);
1312
1313         /* Set overallocation to avoid frequent reallocation of arrays */
1314         set_over_alloc_dd(TRUE);
1315     }
1316     else
1317     {
1318         /* PME, if used, is done on all nodes with 1D decomposition */
1319         cr->npmenodes = 0;
1320         cr->duty      = (DUTY_PP | DUTY_PME);
1321         npme_major    = 1;
1322         npme_minor    = 1;
1323         /* NM and TPI perform single node energy calculations in parallel */
1324         if (!(inputrec->eI == eiNM || EI_TPI(inputrec->eI)))
1325         {
1326             npme_major = cr->nnodes;
1327         }
1328
1329         if (inputrec->ePBC == epbcSCREW)
1330         {
1331             gmx_fatal(FARGS,
1332                       "pbc=%s is only implemented with domain decomposition",
1333                       epbc_names[inputrec->ePBC]);
1334         }
1335     }
1336
1337     if (PAR(cr))
1338     {
1339         /* After possible communicator splitting in make_dd_communicators.
1340          * we can set up the intra/inter node communication.
1341          */
1342         gmx_setup_nodecomm(fplog, cr);
1343     }
1344
1345     /* Initialize per-physical-node MPI process/thread ID and counters. */
1346     gmx_init_intranode_counters(cr);
1347
1348 #ifdef GMX_MPI
1349     md_print_info(cr, fplog, "Using %d MPI %s\n",
1350                   cr->nnodes,
1351 #ifdef GMX_THREAD_MPI
1352                   cr->nnodes == 1 ? "thread" : "threads"
1353 #else
1354                   cr->nnodes == 1 ? "process" : "processes"
1355 #endif
1356                   );
1357     fflush(stderr);
1358 #endif
1359
1360     gmx_omp_nthreads_init(fplog, cr,
1361                           hwinfo->nthreads_hw_avail,
1362                           hw_opt->nthreads_omp,
1363                           hw_opt->nthreads_omp_pme,
1364                           (cr->duty & DUTY_PP) == 0,
1365                           inputrec->cutoff_scheme == ecutsVERLET);
1366
1367     /* check consistency and decide on the number of gpus to use. */
1368     gmx_check_hw_runconf_consistency(fplog, hwinfo, cr, hw_opt->nthreads_tmpi,
1369                                      minf.bUseGPU);
1370
1371     /* getting number of PP/PME threads
1372        PME: env variable should be read only on one node to make sure it is
1373        identical everywhere;
1374      */
1375     /* TODO nthreads_pp is only used for pinning threads.
1376      * This is a temporary solution until we have a hw topology library.
1377      */
1378     nthreads_pp  = gmx_omp_nthreads_get(emntNonbonded);
1379     nthreads_pme = gmx_omp_nthreads_get(emntPME);
1380
1381     wcycle = wallcycle_init(fplog, resetstep, cr, nthreads_pp, nthreads_pme);
1382
1383     if (PAR(cr))
1384     {
1385         /* Master synchronizes its value of reset_counters with all nodes
1386          * including PME only nodes */
1387         reset_counters = wcycle_get_reset_counters(wcycle);
1388         gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1389         wcycle_set_reset_counters(wcycle, reset_counters);
1390     }
1391
1392     snew(nrnb, 1);
1393     if (cr->duty & DUTY_PP)
1394     {
1395         /* For domain decomposition we allocate dynamically
1396          * in dd_partition_system.
1397          */
1398         if (DOMAINDECOMP(cr))
1399         {
1400             bcast_state_setup(cr, state);
1401         }
1402         else
1403         {
1404             if (PAR(cr))
1405             {
1406                 bcast_state(cr, state, TRUE);
1407             }
1408         }
1409
1410         /* Initiate forcerecord */
1411         fr         = mk_forcerec();
1412         fr->hwinfo = hwinfo;
1413         init_forcerec(fplog, oenv, fr, fcd, inputrec, mtop, cr, box, FALSE,
1414                       opt2fn("-table", nfile, fnm),
1415                       opt2fn("-tabletf", nfile, fnm),
1416                       opt2fn("-tablep", nfile, fnm),
1417                       opt2fn("-tableb", nfile, fnm),
1418                       nbpu_opt,
1419                       FALSE, pforce);
1420
1421         /* version for PCA_NOT_READ_NODE (see md.c) */
1422         /*init_forcerec(fplog,fr,fcd,inputrec,mtop,cr,box,FALSE,
1423            "nofile","nofile","nofile","nofile",FALSE,pforce);
1424          */
1425         fr->bSepDVDL = ((Flags & MD_SEPPOT) == MD_SEPPOT);
1426
1427         /* Initialize QM-MM */
1428         if (fr->bQMMM)
1429         {
1430             init_QMMMrec(cr, box, mtop, inputrec, fr);
1431         }
1432
1433         /* Initialize the mdatoms structure.
1434          * mdatoms is not filled with atom data,
1435          * as this can not be done now with domain decomposition.
1436          */
1437         mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1438
1439         if (mdatoms->nPerturbed > 0 && inputrec->cutoff_scheme == ecutsVERLET)
1440         {
1441             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.");
1442         }
1443
1444         /* Initialize the virtual site communication */
1445         vsite = init_vsite(mtop, cr, FALSE);
1446
1447         calc_shifts(box, fr->shift_vec);
1448
1449         /* With periodic molecules the charge groups should be whole at start up
1450          * and the virtual sites should not be far from their proper positions.
1451          */
1452         if (!inputrec->bContinuation && MASTER(cr) &&
1453             !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1454         {
1455             /* Make molecules whole at start of run */
1456             if (fr->ePBC != epbcNONE)
1457             {
1458                 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, state->x);
1459             }
1460             if (vsite)
1461             {
1462                 /* Correct initial vsite positions are required
1463                  * for the initial distribution in the domain decomposition
1464                  * and for the initial shell prediction.
1465                  */
1466                 construct_vsites_mtop(fplog, vsite, mtop, state->x);
1467             }
1468         }
1469
1470         if (EEL_PME(fr->eeltype))
1471         {
1472             ewaldcoeff = fr->ewaldcoeff;
1473             pmedata    = &fr->pmedata;
1474         }
1475         else
1476         {
1477             pmedata = NULL;
1478         }
1479     }
1480     else
1481     {
1482         /* This is a PME only node */
1483
1484         /* We don't need the state */
1485         done_state(state);
1486
1487         ewaldcoeff = calc_ewaldcoeff(inputrec->rcoulomb, inputrec->ewald_rtol);
1488         snew(pmedata, 1);
1489     }
1490
1491     if (hw_opt->thread_affinity != threadaffOFF)
1492     {
1493         /* Before setting affinity, check whether the affinity has changed
1494          * - which indicates that probably the OpenMP library has changed it
1495          * since we first checked).
1496          */
1497         gmx_check_thread_affinity_set(fplog, cr,
1498                                       hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1499
1500         /* Set the CPU affinity */
1501         gmx_set_thread_affinity(fplog, cr, hw_opt, nthreads_pme, hwinfo,
1502                                 inputrec);
1503     }
1504
1505     /* Initiate PME if necessary,
1506      * either on all nodes or on dedicated PME nodes only. */
1507     if (EEL_PME(inputrec->coulombtype))
1508     {
1509         if (mdatoms)
1510         {
1511             nChargePerturbed = mdatoms->nChargePerturbed;
1512         }
1513         if (cr->npmenodes > 0)
1514         {
1515             /* The PME only nodes need to know nChargePerturbed */
1516             gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1517         }
1518
1519         if (cr->duty & DUTY_PME)
1520         {
1521             status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1522                                   mtop ? mtop->natoms : 0, nChargePerturbed,
1523                                   (Flags & MD_REPRODUCIBLE), nthreads_pme);
1524             if (status != 0)
1525             {
1526                 gmx_fatal(FARGS, "Error %d initializing PME", status);
1527             }
1528         }
1529     }
1530
1531
1532     if (integrator[inputrec->eI].func == do_md)
1533     {
1534         /* Turn on signal handling on all nodes */
1535         /*
1536          * (A user signal from the PME nodes (if any)
1537          * is communicated to the PP nodes.
1538          */
1539         signal_handler_install();
1540     }
1541
1542     if (cr->duty & DUTY_PP)
1543     {
1544         if (inputrec->ePull != epullNO)
1545         {
1546             /* Initialize pull code */
1547             init_pull(fplog, inputrec, nfile, fnm, mtop, cr, oenv, inputrec->fepvals->init_lambda,
1548                       EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1549         }
1550
1551         if (inputrec->bRot)
1552         {
1553             /* Initialize enforced rotation code */
1554             init_rot(fplog, inputrec, nfile, fnm, cr, state->x, box, mtop, oenv,
1555                      bVerbose, Flags);
1556         }
1557
1558         constr = init_constraints(fplog, mtop, inputrec, ed, state, cr);
1559
1560         if (DOMAINDECOMP(cr))
1561         {
1562             dd_init_bondeds(fplog, cr->dd, mtop, vsite, constr, inputrec,
1563                             Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1564
1565             set_dd_parameters(fplog, cr->dd, dlb_scale, inputrec, fr, &ddbox);
1566
1567             setup_dd_grid(fplog, cr->dd);
1568         }
1569
1570         /* Now do whatever the user wants us to do (how flexible...) */
1571         integrator[inputrec->eI].func(fplog, cr, nfile, fnm,
1572                                       oenv, bVerbose, bCompact,
1573                                       nstglobalcomm,
1574                                       vsite, constr,
1575                                       nstepout, inputrec, mtop,
1576                                       fcd, state,
1577                                       mdatoms, nrnb, wcycle, ed, fr,
1578                                       repl_ex_nst, repl_ex_nex, repl_ex_seed,
1579                                       membed,
1580                                       cpt_period, max_hours,
1581                                       deviceOptions,
1582                                       Flags,
1583                                       &runtime);
1584
1585         if (inputrec->ePull != epullNO)
1586         {
1587             finish_pull(fplog, inputrec->pull);
1588         }
1589
1590         if (inputrec->bRot)
1591         {
1592             finish_rot(inputrec->rot);
1593         }
1594
1595     }
1596     else
1597     {
1598         /* do PME only */
1599         gmx_pmeonly(*pmedata, cr, nrnb, wcycle, ewaldcoeff, FALSE, inputrec);
1600     }
1601
1602     if (EI_DYNAMICS(inputrec->eI) || EI_TPI(inputrec->eI))
1603     {
1604         /* Some timing stats */
1605         if (SIMMASTER(cr))
1606         {
1607             if (runtime.proc == 0)
1608             {
1609                 runtime.proc = runtime.real;
1610             }
1611         }
1612         else
1613         {
1614             runtime.real = 0;
1615         }
1616     }
1617
1618     wallcycle_stop(wcycle, ewcRUN);
1619
1620     /* Finish up, write some stuff
1621      * if rerunMD, don't write last frame again
1622      */
1623     finish_run(fplog, cr, ftp2fn(efSTO, nfile, fnm),
1624                inputrec, nrnb, wcycle, &runtime,
1625                fr != NULL && fr->nbv != NULL && fr->nbv->bUseGPU ?
1626                nbnxn_cuda_get_timings(fr->nbv->cu_nbv) : NULL,
1627                nthreads_pp,
1628                EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1629
1630     if ((cr->duty & DUTY_PP) && fr->nbv != NULL && fr->nbv->bUseGPU)
1631     {
1632         char gpu_err_str[STRLEN];
1633
1634         /* free GPU memory and uninitialize GPU (by destroying the context) */
1635         nbnxn_cuda_free(fplog, fr->nbv->cu_nbv);
1636
1637         if (!free_gpu(gpu_err_str))
1638         {
1639             gmx_warning("On node %d failed to free GPU #%d: %s",
1640                         cr->nodeid, get_current_gpu_device_id(), gpu_err_str);
1641         }
1642     }
1643
1644     if (opt2bSet("-membed", nfile, fnm))
1645     {
1646         sfree(membed);
1647     }
1648
1649     gmx_hardware_info_free(hwinfo);
1650
1651     /* Does what it says */
1652     print_date_and_time(fplog, cr->nodeid, "Finished mdrun", &runtime);
1653
1654     /* Close logfile already here if we were appending to it */
1655     if (MASTER(cr) && (Flags & MD_APPENDFILES))
1656     {
1657         gmx_log_close(fplog);
1658     }
1659
1660     rc = (int)gmx_get_stop_condition();
1661
1662 #ifdef GMX_THREAD_MPI
1663     /* we need to join all threads. The sub-threads join when they
1664        exit this function, but the master thread needs to be told to
1665        wait for that. */
1666     if (PAR(cr) && MASTER(cr))
1667     {
1668         tMPI_Finalize();
1669     }
1670 #endif
1671
1672     return rc;
1673 }