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