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