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