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