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