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