Merge branch release-5-1 into release-2016
[alexxy/gromacs.git] / src / programs / mdrun / runner.cpp
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,2014,2015,2016, 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 /*! \internal \file
38  *
39  * \brief Implements the MD runner routine calling all integrators.
40  *
41  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42  * \ingroup module_mdlib
43  */
44 #include "gmxpre.h"
45
46 #include "runner.h"
47
48 #include "config.h"
49
50 #include <assert.h>
51 #include <signal.h>
52 #include <stdlib.h>
53 #include <string.h>
54
55 #include <algorithm>
56
57 #include "gromacs/commandline/filenm.h"
58 #include "gromacs/domdec/domdec.h"
59 #include "gromacs/domdec/domdec_struct.h"
60 #include "gromacs/essentialdynamics/edsam.h"
61 #include "gromacs/ewald/pme.h"
62 #include "gromacs/fileio/checkpoint.h"
63 #include "gromacs/fileio/oenv.h"
64 #include "gromacs/fileio/tpxio.h"
65 #include "gromacs/gmxlib/md_logging.h"
66 #include "gromacs/gmxlib/network.h"
67 #include "gromacs/gpu_utils/gpu_utils.h"
68 #include "gromacs/hardware/cpuinfo.h"
69 #include "gromacs/hardware/detecthardware.h"
70 #include "gromacs/listed-forces/disre.h"
71 #include "gromacs/listed-forces/orires.h"
72 #include "gromacs/math/calculate-ewald-splitting-coefficient.h"
73 #include "gromacs/math/functions.h"
74 #include "gromacs/math/utilities.h"
75 #include "gromacs/math/vec.h"
76 #include "gromacs/mdlib/calc_verletbuf.h"
77 #include "gromacs/mdlib/constr.h"
78 #include "gromacs/mdlib/force.h"
79 #include "gromacs/mdlib/forcerec.h"
80 #include "gromacs/mdlib/gmx_omp_nthreads.h"
81 #include "gromacs/mdlib/integrator.h"
82 #include "gromacs/mdlib/main.h"
83 #include "gromacs/mdlib/md_support.h"
84 #include "gromacs/mdlib/mdatoms.h"
85 #include "gromacs/mdlib/mdrun.h"
86 #include "gromacs/mdlib/minimize.h"
87 #include "gromacs/mdlib/nbnxn_search.h"
88 #include "gromacs/mdlib/qmmm.h"
89 #include "gromacs/mdlib/sighandler.h"
90 #include "gromacs/mdlib/sim_util.h"
91 #include "gromacs/mdlib/tpi.h"
92 #include "gromacs/mdrunutility/threadaffinity.h"
93 #include "gromacs/mdtypes/commrec.h"
94 #include "gromacs/mdtypes/inputrec.h"
95 #include "gromacs/mdtypes/md_enums.h"
96 #include "gromacs/mdtypes/state.h"
97 #include "gromacs/pbcutil/pbc.h"
98 #include "gromacs/pulling/pull.h"
99 #include "gromacs/pulling/pull_rotation.h"
100 #include "gromacs/timing/wallcycle.h"
101 #include "gromacs/topology/mtop_util.h"
102 #include "gromacs/trajectory/trajectoryframe.h"
103 #include "gromacs/utility/cstringutil.h"
104 #include "gromacs/utility/exceptions.h"
105 #include "gromacs/utility/fatalerror.h"
106 #include "gromacs/utility/gmxassert.h"
107 #include "gromacs/utility/gmxmpi.h"
108 #include "gromacs/utility/pleasecite.h"
109 #include "gromacs/utility/smalloc.h"
110
111 #include "deform.h"
112 #include "md.h"
113 #include "repl_ex.h"
114 #include "resource-division.h"
115
116 #ifdef GMX_FAHCORE
117 #include "corewrap.h"
118 #endif
119
120 //! First step used in pressure scaling
121 gmx_int64_t         deform_init_init_step_tpx;
122 //! Initial box for pressure scaling
123 matrix              deform_init_box_tpx;
124 //! MPI variable for use in pressure scaling
125 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
126
127 #if GMX_THREAD_MPI
128 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
129  * the number of threads will get lowered.
130  */
131 #define MIN_ATOMS_PER_MPI_THREAD    90
132 #define MIN_ATOMS_PER_GPU           900
133
134 struct mdrunner_arglist
135 {
136     gmx_hw_opt_t            hw_opt;
137     FILE                   *fplog;
138     t_commrec              *cr;
139     int                     nfile;
140     const t_filenm         *fnm;
141     const gmx_output_env_t *oenv;
142     gmx_bool                bVerbose;
143     int                     nstglobalcomm;
144     ivec                    ddxyz;
145     int                     dd_rank_order;
146     int                     npme;
147     real                    rdd;
148     real                    rconstr;
149     const char             *dddlb_opt;
150     real                    dlb_scale;
151     const char             *ddcsx;
152     const char             *ddcsy;
153     const char             *ddcsz;
154     const char             *nbpu_opt;
155     int                     nstlist_cmdline;
156     gmx_int64_t             nsteps_cmdline;
157     int                     nstepout;
158     int                     resetstep;
159     int                     nmultisim;
160     int                     repl_ex_nst;
161     int                     repl_ex_nex;
162     int                     repl_ex_seed;
163     real                    pforce;
164     real                    cpt_period;
165     real                    max_hours;
166     int                     imdport;
167     unsigned long           Flags;
168 };
169
170
171 /* The function used for spawning threads. Extracts the mdrunner()
172    arguments from its one argument and calls mdrunner(), after making
173    a commrec. */
174 static void mdrunner_start_fn(void *arg)
175 {
176     try
177     {
178         struct mdrunner_arglist *mda = (struct mdrunner_arglist*)arg;
179         struct mdrunner_arglist  mc  = *mda; /* copy the arg list to make sure
180                                                 that it's thread-local. This doesn't
181                                                 copy pointed-to items, of course,
182                                                 but those are all const. */
183         t_commrec *cr;                       /* we need a local version of this */
184         FILE      *fplog = NULL;
185         t_filenm  *fnm;
186
187         fnm = dup_tfn(mc.nfile, mc.fnm);
188
189         cr = reinitialize_commrec_for_this_thread(mc.cr);
190
191         if (MASTER(cr))
192         {
193             fplog = mc.fplog;
194         }
195
196         gmx::mdrunner(&mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
197                       mc.bVerbose, mc.nstglobalcomm,
198                       mc.ddxyz, mc.dd_rank_order, mc.npme, mc.rdd,
199                       mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
200                       mc.ddcsx, mc.ddcsy, mc.ddcsz,
201                       mc.nbpu_opt, mc.nstlist_cmdline,
202                       mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
203                       mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce,
204                       mc.cpt_period, mc.max_hours, mc.imdport, mc.Flags);
205     }
206     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
207 }
208
209
210 /* called by mdrunner() to start a specific number of threads (including
211    the main thread) for thread-parallel runs. This in turn calls mdrunner()
212    for each thread.
213    All options besides nthreads are the same as for mdrunner(). */
214 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
215                                          FILE *fplog, t_commrec *cr, int nfile,
216                                          const t_filenm fnm[], const gmx_output_env_t *oenv, gmx_bool bVerbose,
217                                          int nstglobalcomm,
218                                          ivec ddxyz, int dd_rank_order, int npme,
219                                          real rdd, real rconstr,
220                                          const char *dddlb_opt, real dlb_scale,
221                                          const char *ddcsx, const char *ddcsy, const char *ddcsz,
222                                          const char *nbpu_opt, int nstlist_cmdline,
223                                          gmx_int64_t nsteps_cmdline,
224                                          int nstepout, int resetstep,
225                                          int nmultisim, int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
226                                          real pforce, real cpt_period, real max_hours,
227                                          unsigned long Flags)
228 {
229     int                      ret;
230     struct mdrunner_arglist *mda;
231     t_commrec               *crn; /* the new commrec */
232     t_filenm                *fnmn;
233
234     /* first check whether we even need to start tMPI */
235     if (hw_opt->nthreads_tmpi < 2)
236     {
237         return cr;
238     }
239
240     /* a few small, one-time, almost unavoidable memory leaks: */
241     snew(mda, 1);
242     fnmn = dup_tfn(nfile, fnm);
243
244     /* fill the data structure to pass as void pointer to thread start fn */
245     /* hw_opt contains pointers, which should all be NULL at this stage */
246     mda->hw_opt          = *hw_opt;
247     mda->fplog           = fplog;
248     mda->cr              = cr;
249     mda->nfile           = nfile;
250     mda->fnm             = fnmn;
251     mda->oenv            = oenv;
252     mda->bVerbose        = bVerbose;
253     mda->nstglobalcomm   = nstglobalcomm;
254     mda->ddxyz[XX]       = ddxyz[XX];
255     mda->ddxyz[YY]       = ddxyz[YY];
256     mda->ddxyz[ZZ]       = ddxyz[ZZ];
257     mda->dd_rank_order   = dd_rank_order;
258     mda->npme            = npme;
259     mda->rdd             = rdd;
260     mda->rconstr         = rconstr;
261     mda->dddlb_opt       = dddlb_opt;
262     mda->dlb_scale       = dlb_scale;
263     mda->ddcsx           = ddcsx;
264     mda->ddcsy           = ddcsy;
265     mda->ddcsz           = ddcsz;
266     mda->nbpu_opt        = nbpu_opt;
267     mda->nstlist_cmdline = nstlist_cmdline;
268     mda->nsteps_cmdline  = nsteps_cmdline;
269     mda->nstepout        = nstepout;
270     mda->resetstep       = resetstep;
271     mda->nmultisim       = nmultisim;
272     mda->repl_ex_nst     = repl_ex_nst;
273     mda->repl_ex_nex     = repl_ex_nex;
274     mda->repl_ex_seed    = repl_ex_seed;
275     mda->pforce          = pforce;
276     mda->cpt_period      = cpt_period;
277     mda->max_hours       = max_hours;
278     mda->Flags           = Flags;
279
280     /* now spawn new threads that start mdrunner_start_fn(), while
281        the main thread returns, we set thread affinity later */
282     ret = tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi, TMPI_AFFINITY_NONE,
283                        mdrunner_start_fn, (void*)(mda) );
284     if (ret != TMPI_SUCCESS)
285     {
286         return NULL;
287     }
288
289     crn = reinitialize_commrec_for_this_thread(cr);
290     return crn;
291 }
292
293 #endif /* GMX_THREAD_MPI */
294
295
296 /*! \brief Cost of non-bonded kernels
297  *
298  * We determine the extra cost of the non-bonded kernels compared to
299  * a reference nstlist value of 10 (which is the default in grompp).
300  */
301 static const int    nbnxnReferenceNstlist = 10;
302 //! The values to try when switching
303 const int           nstlist_try[] = { 20, 25, 40 };
304 //! Number of elements in the neighborsearch list trials.
305 #define NNSTL  sizeof(nstlist_try)/sizeof(nstlist_try[0])
306 /* Increase nstlist until the non-bonded cost increases more than listfac_ok,
307  * but never more than listfac_max.
308  * A standard (protein+)water system at 300K with PME ewald_rtol=1e-5
309  * needs 1.28 at rcoulomb=0.9 and 1.24 at rcoulomb=1.0 to get to nstlist=40.
310  * Note that both CPU and GPU factors are conservative. Performance should
311  * not go down due to this tuning, except with a relatively slow GPU.
312  * On the other hand, at medium/high parallelization or with fast GPUs
313  * nstlist will not be increased enough to reach optimal performance.
314  */
315 /* CPU: pair-search is about a factor 1.5 slower than the non-bonded kernel */
316 //! Max OK performance ratio beween force calc and neighbor searching
317 static const float  nbnxn_cpu_listfac_ok    = 1.05;
318 //! Too high performance ratio beween force calc and neighbor searching
319 static const float  nbnxn_cpu_listfac_max   = 1.09;
320 /* CPU: pair-search is about a factor 2-3 slower than the non-bonded kernel */
321 //! Max OK performance ratio beween force calc and neighbor searching
322 static const float  nbnxn_knl_listfac_ok    = 1.22;
323 //! Too high performance ratio beween force calc and neighbor searching
324 static const float  nbnxn_knl_listfac_max   = 1.3;
325 /* GPU: pair-search is a factor 1.5-3 slower than the non-bonded kernel */
326 //! Max OK performance ratio beween force calc and neighbor searching
327 static const float  nbnxn_gpu_listfac_ok    = 1.20;
328 //! Too high performance ratio beween force calc and neighbor searching
329 static const float  nbnxn_gpu_listfac_max   = 1.30;
330
331 /*! \brief Try to increase nstlist when using the Verlet cut-off scheme */
332 static void increase_nstlist(FILE *fp, t_commrec *cr,
333                              t_inputrec *ir, int nstlist_cmdline,
334                              const gmx_mtop_t *mtop, matrix box,
335                              gmx_bool bGPU, const gmx::CpuInfo &cpuinfo)
336 {
337     float                  listfac_ok, listfac_max;
338     int                    nstlist_orig, nstlist_prev;
339     verletbuf_list_setup_t ls;
340     real                   rlistWithReferenceNstlist, rlist_inc, rlist_ok, rlist_max;
341     real                   rlist_new, rlist_prev;
342     size_t                 nstlist_ind = 0;
343     t_state                state_tmp;
344     gmx_bool               bBox, bDD, bCont;
345     const char            *nstl_gpu = "\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";
346     const char            *nve_err  = "Can not increase nstlist because an NVE ensemble is used";
347     const char            *vbd_err  = "Can not increase nstlist because verlet-buffer-tolerance is not set or used";
348     const char            *box_err  = "Can not increase nstlist because the box is too small";
349     const char            *dd_err   = "Can not increase nstlist because of domain decomposition limitations";
350     char                   buf[STRLEN];
351
352     if (nstlist_cmdline <= 0)
353     {
354         if (ir->nstlist == 1)
355         {
356             /* The user probably set nstlist=1 for a reason,
357              * don't mess with the settings.
358              */
359             return;
360         }
361
362         if (fp != NULL && bGPU && ir->nstlist < nstlist_try[0])
363         {
364             fprintf(fp, nstl_gpu, ir->nstlist);
365         }
366         nstlist_ind = 0;
367         while (nstlist_ind < NNSTL && ir->nstlist >= nstlist_try[nstlist_ind])
368         {
369             nstlist_ind++;
370         }
371         if (nstlist_ind == NNSTL)
372         {
373             /* There are no larger nstlist value to try */
374             return;
375         }
376     }
377
378     if (EI_MD(ir->eI) && ir->etc == etcNO)
379     {
380         if (MASTER(cr))
381         {
382             fprintf(stderr, "%s\n", nve_err);
383         }
384         if (fp != NULL)
385         {
386             fprintf(fp, "%s\n", nve_err);
387         }
388
389         return;
390     }
391
392     if (ir->verletbuf_tol == 0 && bGPU)
393     {
394         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");
395     }
396
397     if (ir->verletbuf_tol < 0)
398     {
399         if (MASTER(cr))
400         {
401             fprintf(stderr, "%s\n", vbd_err);
402         }
403         if (fp != NULL)
404         {
405             fprintf(fp, "%s\n", vbd_err);
406         }
407
408         return;
409     }
410
411     if (bGPU)
412     {
413         listfac_ok  = nbnxn_gpu_listfac_ok;
414         listfac_max = nbnxn_gpu_listfac_max;
415     }
416     else if (cpuinfo.feature(gmx::CpuInfo::Feature::X86_Avx512ER))
417     {
418         listfac_ok  = nbnxn_knl_listfac_ok;
419         listfac_max = nbnxn_knl_listfac_max;
420     }
421     else
422     {
423         listfac_ok  = nbnxn_cpu_listfac_ok;
424         listfac_max = nbnxn_cpu_listfac_max;
425     }
426
427     nstlist_orig = ir->nstlist;
428     if (nstlist_cmdline > 0)
429     {
430         if (fp)
431         {
432             sprintf(buf, "Getting nstlist=%d from command line option",
433                     nstlist_cmdline);
434         }
435         ir->nstlist = nstlist_cmdline;
436     }
437
438     verletbuf_get_list_setup(TRUE, bGPU, &ls);
439
440     /* Allow rlist to make the list a given factor larger than the list
441      * would be with the reference value for nstlist (10).
442      */
443     nstlist_prev = ir->nstlist;
444     ir->nstlist  = nbnxnReferenceNstlist;
445     calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL,
446                             &rlistWithReferenceNstlist);
447     ir->nstlist  = nstlist_prev;
448
449     /* Determine the pair list size increase due to zero interactions */
450     rlist_inc = nbnxn_get_rlist_effective_inc(ls.cluster_size_j,
451                                               mtop->natoms/det(box));
452     rlist_ok  = (rlistWithReferenceNstlist + rlist_inc)*std::cbrt(listfac_ok) - rlist_inc;
453     rlist_max = (rlistWithReferenceNstlist + rlist_inc)*std::cbrt(listfac_max) - rlist_inc;
454     if (debug)
455     {
456         fprintf(debug, "nstlist tuning: rlist_inc %.3f rlist_ok %.3f rlist_max %.3f\n",
457                 rlist_inc, rlist_ok, rlist_max);
458     }
459
460     nstlist_prev = nstlist_orig;
461     rlist_prev   = ir->rlist;
462     do
463     {
464         if (nstlist_cmdline <= 0)
465         {
466             ir->nstlist = nstlist_try[nstlist_ind];
467         }
468
469         /* Set the pair-list buffer size in ir */
470         calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
471
472         /* Does rlist fit in the box? */
473         bBox = (gmx::square(rlist_new) < max_cutoff2(ir->ePBC, box));
474         bDD  = TRUE;
475         if (bBox && DOMAINDECOMP(cr))
476         {
477             /* Check if rlist fits in the domain decomposition */
478             if (inputrec2nboundeddim(ir) < DIM)
479             {
480                 gmx_incons("Changing nstlist with domain decomposition and unbounded dimensions is not implemented yet");
481             }
482             copy_mat(box, state_tmp.box);
483             bDD = change_dd_cutoff(cr, &state_tmp, ir, rlist_new);
484         }
485
486         if (debug)
487         {
488             fprintf(debug, "nstlist %d rlist %.3f bBox %d bDD %d\n",
489                     ir->nstlist, rlist_new, bBox, bDD);
490         }
491
492         bCont = FALSE;
493
494         if (nstlist_cmdline <= 0)
495         {
496             if (bBox && bDD && rlist_new <= rlist_max)
497             {
498                 /* Increase nstlist */
499                 nstlist_prev = ir->nstlist;
500                 rlist_prev   = rlist_new;
501                 bCont        = (nstlist_ind+1 < NNSTL && rlist_new < rlist_ok);
502             }
503             else
504             {
505                 /* Stick with the previous nstlist */
506                 ir->nstlist = nstlist_prev;
507                 rlist_new   = rlist_prev;
508                 bBox        = TRUE;
509                 bDD         = TRUE;
510             }
511         }
512
513         nstlist_ind++;
514     }
515     while (bCont);
516
517     if (!bBox || !bDD)
518     {
519         gmx_warning(!bBox ? box_err : dd_err);
520         if (fp != NULL)
521         {
522             fprintf(fp, "\n%s\n", bBox ? box_err : dd_err);
523         }
524         ir->nstlist = nstlist_orig;
525     }
526     else if (ir->nstlist != nstlist_orig || rlist_new != ir->rlist)
527     {
528         sprintf(buf, "Changing nstlist from %d to %d, rlist from %g to %g",
529                 nstlist_orig, ir->nstlist,
530                 ir->rlist, rlist_new);
531         if (MASTER(cr))
532         {
533             fprintf(stderr, "%s\n\n", buf);
534         }
535         if (fp != NULL)
536         {
537             fprintf(fp, "%s\n\n", buf);
538         }
539         ir->rlist     = rlist_new;
540     }
541 }
542
543 /*! \brief Initialize variables for Verlet scheme simulation */
544 static void prepare_verlet_scheme(FILE                           *fplog,
545                                   t_commrec                      *cr,
546                                   t_inputrec                     *ir,
547                                   int                             nstlist_cmdline,
548                                   const gmx_mtop_t               *mtop,
549                                   matrix                          box,
550                                   gmx_bool                        bUseGPU,
551                                   const gmx::CpuInfo             &cpuinfo)
552 {
553     /* For NVE simulations, we will retain the initial list buffer */
554     if (EI_DYNAMICS(ir->eI) &&
555         ir->verletbuf_tol > 0 &&
556         !(EI_MD(ir->eI) && ir->etc == etcNO))
557     {
558         /* Update the Verlet buffer size for the current run setup */
559         verletbuf_list_setup_t ls;
560         real                   rlist_new;
561
562         /* Here we assume SIMD-enabled kernels are being used. But as currently
563          * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
564          * and 4x2 gives a larger buffer than 4x4, this is ok.
565          */
566         verletbuf_get_list_setup(TRUE, bUseGPU, &ls);
567
568         calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
569
570         if (rlist_new != ir->rlist)
571         {
572             if (fplog != NULL)
573             {
574                 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
575                         ir->rlist, rlist_new,
576                         ls.cluster_size_i, ls.cluster_size_j);
577             }
578             ir->rlist     = rlist_new;
579         }
580     }
581
582     if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
583     {
584         gmx_fatal(FARGS, "Can not set nstlist without %s",
585                   !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
586     }
587
588     if (EI_DYNAMICS(ir->eI))
589     {
590         /* Set or try nstlist values */
591         increase_nstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, bUseGPU, cpuinfo);
592     }
593 }
594
595 /*! \brief Override the nslist value in inputrec
596  *
597  * with value passed on the command line (if any)
598  */
599 static void override_nsteps_cmdline(FILE            *fplog,
600                                     gmx_int64_t      nsteps_cmdline,
601                                     t_inputrec      *ir,
602                                     const t_commrec *cr)
603 {
604     assert(ir);
605     assert(cr);
606
607     /* override with anything else than the default -2 */
608     if (nsteps_cmdline > -2)
609     {
610         char sbuf_steps[STEPSTRSIZE];
611         char sbuf_msg[STRLEN];
612
613         ir->nsteps = nsteps_cmdline;
614         if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
615         {
616             sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
617                     gmx_step_str(nsteps_cmdline, sbuf_steps),
618                     fabs(nsteps_cmdline*ir->delta_t));
619         }
620         else
621         {
622             sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
623                     gmx_step_str(nsteps_cmdline, sbuf_steps));
624         }
625
626         md_print_warn(cr, fplog, "%s\n", sbuf_msg);
627     }
628     else if (nsteps_cmdline < -2)
629     {
630         gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %d",
631                   nsteps_cmdline);
632     }
633     /* Do nothing if nsteps_cmdline == -2 */
634 }
635
636 namespace gmx
637 {
638
639 //! \brief Return the correct integrator function.
640 static integrator_t *my_integrator(unsigned int ei)
641 {
642     switch (ei)
643     {
644         case eiMD:
645         case eiBD:
646         case eiSD1:
647         case eiVV:
648         case eiVVAK:
649             if (!EI_DYNAMICS(ei))
650             {
651                 GMX_THROW(APIError("do_md integrator would be called for a non-dynamical integrator"));
652             }
653             return do_md;
654         case eiSteep:
655             return do_steep;
656         case eiCG:
657             return do_cg;
658         case eiNM:
659             return do_nm;
660         case eiLBFGS:
661             return do_lbfgs;
662         case eiTPI:
663         case eiTPIC:
664             if (!EI_TPI(ei))
665             {
666                 GMX_THROW(APIError("do_tpi integrator would be called for a non-TPI integrator"));
667             }
668             return do_tpi;
669         case eiSD2_REMOVED:
670             GMX_THROW(NotImplementedError("SD2 integrator has been removed"));
671         default:
672             GMX_THROW(APIError("Non existing integrator selected"));
673     }
674 }
675
676 int mdrunner(gmx_hw_opt_t *hw_opt,
677              FILE *fplog, t_commrec *cr, int nfile,
678              const t_filenm fnm[], const gmx_output_env_t *oenv, gmx_bool bVerbose,
679              int nstglobalcomm,
680              ivec ddxyz, int dd_rank_order, int npme, real rdd, real rconstr,
681              const char *dddlb_opt, real dlb_scale,
682              const char *ddcsx, const char *ddcsy, const char *ddcsz,
683              const char *nbpu_opt, int nstlist_cmdline,
684              gmx_int64_t nsteps_cmdline, int nstepout, int resetstep,
685              int gmx_unused nmultisim, int repl_ex_nst, int repl_ex_nex,
686              int repl_ex_seed, real pforce, real cpt_period, real max_hours,
687              int imdport, unsigned long Flags)
688 {
689     gmx_bool                  bForceUseGPU, bTryUseGPU, bRerunMD;
690     t_inputrec               *inputrec;
691     t_state                  *state = NULL;
692     matrix                    box;
693     gmx_ddbox_t               ddbox = {0};
694     int                       npme_major, npme_minor;
695     t_nrnb                   *nrnb;
696     gmx_mtop_t               *mtop          = NULL;
697     t_mdatoms                *mdatoms       = NULL;
698     t_forcerec               *fr            = NULL;
699     t_fcdata                 *fcd           = NULL;
700     real                      ewaldcoeff_q  = 0;
701     real                      ewaldcoeff_lj = 0;
702     struct gmx_pme_t        **pmedata       = NULL;
703     gmx_vsite_t              *vsite         = NULL;
704     gmx_constr_t              constr;
705     int                       nChargePerturbed = -1, nTypePerturbed = 0, status;
706     gmx_wallcycle_t           wcycle;
707     gmx_walltime_accounting_t walltime_accounting = NULL;
708     int                       rc;
709     gmx_int64_t               reset_counters;
710     gmx_edsam_t               ed           = NULL;
711     int                       nthreads_pme = 1;
712     gmx_hw_info_t            *hwinfo       = NULL;
713     /* The master rank decides early on bUseGPU and broadcasts this later */
714     gmx_bool                  bUseGPU            = FALSE;
715
716     /* CAUTION: threads may be started later on in this function, so
717        cr doesn't reflect the final parallel state right now */
718     snew(inputrec, 1);
719     snew(mtop, 1);
720
721     if (Flags & MD_APPENDFILES)
722     {
723         fplog = NULL;
724     }
725
726     bRerunMD     = (Flags & MD_RERUN);
727     bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
728     bTryUseGPU   = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
729
730     /* Detect hardware, gather information. This is an operation that is
731      * global for this process (MPI rank). */
732     hwinfo = gmx_detect_hardware(fplog, cr, bTryUseGPU);
733
734     gmx_print_detected_hardware(fplog, cr, hwinfo);
735
736     if (fplog != NULL)
737     {
738         /* Print references after all software/hardware printing */
739         please_cite(fplog, "Abraham2015");
740         please_cite(fplog, "Pall2015");
741         please_cite(fplog, "Pronk2013");
742         please_cite(fplog, "Hess2008b");
743         please_cite(fplog, "Spoel2005a");
744         please_cite(fplog, "Lindahl2001a");
745         please_cite(fplog, "Berendsen95a");
746     }
747
748     snew(state, 1);
749     if (SIMMASTER(cr))
750     {
751         /* Read (nearly) all data required for the simulation */
752         read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, state, mtop);
753
754         if (inputrec->cutoff_scheme == ecutsVERLET)
755         {
756             /* Here the master rank decides if all ranks will use GPUs */
757             bUseGPU = (hwinfo->gpu_info.n_dev_compatible > 0 ||
758                        getenv("GMX_EMULATE_GPU") != NULL);
759
760             /* TODO add GPU kernels for this and replace this check by:
761              * (bUseGPU && (ir->vdwtype == evdwPME &&
762              *               ir->ljpme_combination_rule == eljpmeLB))
763              * update the message text and the content of nbnxn_acceleration_supported.
764              */
765             if (bUseGPU &&
766                 !nbnxn_gpu_acceleration_supported(fplog, cr, inputrec, bRerunMD))
767             {
768                 /* Fallback message printed by nbnxn_acceleration_supported */
769                 if (bForceUseGPU)
770                 {
771                     gmx_fatal(FARGS, "GPU acceleration requested, but not supported with the given input settings");
772                 }
773                 bUseGPU = FALSE;
774             }
775
776             prepare_verlet_scheme(fplog, cr,
777                                   inputrec, nstlist_cmdline, mtop, state->box,
778                                   bUseGPU, *hwinfo->cpuInfo);
779         }
780         else
781         {
782             if (nstlist_cmdline > 0)
783             {
784                 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
785             }
786
787             if (hwinfo->gpu_info.n_dev_compatible > 0)
788             {
789                 md_print_warn(cr, fplog,
790                               "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
791                               "      To use a GPU, set the mdp option: cutoff-scheme = Verlet\n");
792             }
793
794             if (bForceUseGPU)
795             {
796                 gmx_fatal(FARGS, "GPU requested, but can't be used without cutoff-scheme=Verlet");
797             }
798
799 #if GMX_TARGET_BGQ
800             md_print_warn(cr, fplog,
801                           "NOTE: There is no SIMD implementation of the group scheme kernels on\n"
802                           "      BlueGene/Q. You will observe better performance from using the\n"
803                           "      Verlet cut-off scheme.\n");
804 #endif
805         }
806     }
807
808     /* Check and update the hardware options for internal consistency */
809     check_and_update_hw_opt_1(hw_opt, cr, npme);
810
811     /* Early check for externally set process affinity. */
812     gmx_check_thread_affinity_set(fplog, cr,
813                                   hw_opt, hwinfo->nthreads_hw_avail, FALSE);
814
815 #if GMX_THREAD_MPI
816     if (SIMMASTER(cr))
817     {
818         if (npme > 0 && hw_opt->nthreads_tmpi <= 0)
819         {
820             gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
821         }
822
823         /* Since the master knows the cut-off scheme, update hw_opt for this.
824          * This is done later for normal MPI and also once more with tMPI
825          * for all tMPI ranks.
826          */
827         check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
828
829         /* NOW the threads will be started: */
830         hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
831                                                  hw_opt,
832                                                  inputrec, mtop,
833                                                  cr, fplog, bUseGPU);
834
835         if (hw_opt->nthreads_tmpi > 1)
836         {
837             t_commrec *cr_old       = cr;
838             /* now start the threads. */
839             cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
840                                         oenv, bVerbose, nstglobalcomm,
841                                         ddxyz, dd_rank_order, npme, rdd, rconstr,
842                                         dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
843                                         nbpu_opt, nstlist_cmdline,
844                                         nsteps_cmdline, nstepout, resetstep, nmultisim,
845                                         repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
846                                         cpt_period, max_hours,
847                                         Flags);
848             /* the main thread continues here with a new cr. We don't deallocate
849                the old cr because other threads may still be reading it. */
850             if (cr == NULL)
851             {
852                 gmx_comm("Failed to spawn threads");
853             }
854         }
855     }
856 #endif
857     /* END OF CAUTION: cr is now reliable */
858
859     if (PAR(cr))
860     {
861         /* now broadcast everything to the non-master nodes/threads: */
862         init_parallel(cr, inputrec, mtop);
863
864         /* The master rank decided on the use of GPUs,
865          * broadcast this information to all ranks.
866          */
867         gmx_bcast_sim(sizeof(bUseGPU), &bUseGPU, cr);
868     }
869
870     if (fplog != NULL)
871     {
872         pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
873         fprintf(fplog, "\n");
874     }
875
876     /* now make sure the state is initialized and propagated */
877     set_state_entries(state, inputrec);
878
879     /* A parallel command line option consistency check that we can
880        only do after any threads have started. */
881     if (!PAR(cr) &&
882         (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || npme > 0))
883     {
884         gmx_fatal(FARGS,
885                   "The -dd or -npme option request a parallel simulation, "
886 #if !GMX_MPI
887                   "but %s was compiled without threads or MPI enabled"
888 #else
889 #if GMX_THREAD_MPI
890                   "but the number of MPI-threads (option -ntmpi) is not set or is 1"
891 #else
892                   "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
893 #endif
894 #endif
895                   , output_env_get_program_display_name(oenv)
896                   );
897     }
898
899     if (bRerunMD &&
900         (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
901     {
902         gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
903     }
904
905     if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
906     {
907         gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
908     }
909
910     if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
911     {
912         if (npme > 0)
913         {
914             gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
915                                  "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
916         }
917
918         npme = 0;
919     }
920
921     if (bUseGPU && npme < 0)
922     {
923         /* With GPUs we don't automatically use PME-only ranks. PME ranks can
924          * improve performance with many threads per GPU, since our OpenMP
925          * scaling is bad, but it's difficult to automate the setup.
926          */
927         npme = 0;
928     }
929
930 #ifdef GMX_FAHCORE
931     if (MASTER(cr))
932     {
933         fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
934     }
935 #endif
936
937     /* NMR restraints must be initialized before load_checkpoint,
938      * since with time averaging the history is added to t_state.
939      * For proper consistency check we therefore need to extend
940      * t_state here.
941      * So the PME-only nodes (if present) will also initialize
942      * the distance restraints.
943      */
944     snew(fcd, 1);
945
946     /* This needs to be called before read_checkpoint to extend the state */
947     init_disres(fplog, mtop, inputrec, cr, fcd, state, repl_ex_nst > 0);
948
949     init_orires(fplog, mtop, state->x, inputrec, cr, &(fcd->orires),
950                 state);
951
952     if (inputrecDeform(inputrec))
953     {
954         /* Store the deform reference box before reading the checkpoint */
955         if (SIMMASTER(cr))
956         {
957             copy_mat(state->box, box);
958         }
959         if (PAR(cr))
960         {
961             gmx_bcast(sizeof(box), box, cr);
962         }
963         /* Because we do not have the update struct available yet
964          * in which the reference values should be stored,
965          * we store them temporarily in static variables.
966          * This should be thread safe, since they are only written once
967          * and with identical values.
968          */
969         tMPI_Thread_mutex_lock(&deform_init_box_mutex);
970         deform_init_init_step_tpx = inputrec->init_step;
971         copy_mat(box, deform_init_box_tpx);
972         tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
973     }
974
975     if (Flags & MD_STARTFROMCPT)
976     {
977         /* Check if checkpoint file exists before doing continuation.
978          * This way we can use identical input options for the first and subsequent runs...
979          */
980         gmx_bool bReadEkin;
981
982         load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
983                         cr, ddxyz, &npme,
984                         inputrec, state, &bReadEkin,
985                         (Flags & MD_APPENDFILES),
986                         (Flags & MD_APPENDFILESSET));
987
988         if (bReadEkin)
989         {
990             Flags |= MD_READ_EKIN;
991         }
992     }
993
994     if (MASTER(cr) && (Flags & MD_APPENDFILES))
995     {
996         gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
997                      Flags, &fplog);
998     }
999
1000     /* override nsteps with value from cmdline */
1001     override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
1002
1003     if (SIMMASTER(cr))
1004     {
1005         copy_mat(state->box, box);
1006     }
1007
1008     if (PAR(cr))
1009     {
1010         gmx_bcast(sizeof(box), box, cr);
1011     }
1012
1013     // TODO This should move to do_md(), because it only makes sense
1014     // with dynamical integrators, but there is no test coverage and
1015     // it interacts with constraints, somehow.
1016     /* Essential dynamics */
1017     if (opt2bSet("-ei", nfile, fnm))
1018     {
1019         /* Open input and output files, allocate space for ED data structure */
1020         ed = ed_open(mtop->natoms, &state->edsamstate, nfile, fnm, Flags, oenv, cr);
1021     }
1022
1023     if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
1024                      inputrec->eI == eiNM))
1025     {
1026         cr->dd = init_domain_decomposition(fplog, cr, Flags, ddxyz, npme,
1027                                            dd_rank_order,
1028                                            rdd, rconstr,
1029                                            dddlb_opt, dlb_scale,
1030                                            ddcsx, ddcsy, ddcsz,
1031                                            mtop, inputrec,
1032                                            box, state->x,
1033                                            &ddbox, &npme_major, &npme_minor);
1034     }
1035     else
1036     {
1037         /* PME, if used, is done on all nodes with 1D decomposition */
1038         cr->npmenodes = 0;
1039         cr->duty      = (DUTY_PP | DUTY_PME);
1040         npme_major    = 1;
1041         npme_minor    = 1;
1042
1043         if (inputrec->ePBC == epbcSCREW)
1044         {
1045             gmx_fatal(FARGS,
1046                       "pbc=%s is only implemented with domain decomposition",
1047                       epbc_names[inputrec->ePBC]);
1048         }
1049     }
1050
1051     if (PAR(cr))
1052     {
1053         /* After possible communicator splitting in make_dd_communicators.
1054          * we can set up the intra/inter node communication.
1055          */
1056         gmx_setup_nodecomm(fplog, cr);
1057     }
1058
1059     /* Initialize per-physical-node MPI process/thread ID and counters. */
1060     gmx_init_intranode_counters(cr);
1061 #if GMX_MPI
1062     if (MULTISIM(cr))
1063     {
1064         md_print_info(cr, fplog,
1065                       "This is simulation %d out of %d running as a composite GROMACS\n"
1066                       "multi-simulation job. Setup for this simulation:\n\n",
1067                       cr->ms->sim, cr->ms->nsim);
1068     }
1069     md_print_info(cr, fplog, "Using %d MPI %s\n",
1070                   cr->nnodes,
1071 #if GMX_THREAD_MPI
1072                   cr->nnodes == 1 ? "thread" : "threads"
1073 #else
1074                   cr->nnodes == 1 ? "process" : "processes"
1075 #endif
1076                   );
1077     fflush(stderr);
1078 #endif
1079
1080     /* Check and update hw_opt for the cut-off scheme */
1081     check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1082
1083     /* Check and update hw_opt for the number of MPI ranks */
1084     check_and_update_hw_opt_3(hw_opt);
1085
1086     gmx_omp_nthreads_init(fplog, cr,
1087                           hwinfo->nthreads_hw_avail,
1088                           hw_opt->nthreads_omp,
1089                           hw_opt->nthreads_omp_pme,
1090                           (cr->duty & DUTY_PP) == 0,
1091                           inputrec->cutoff_scheme == ecutsVERLET);
1092
1093 #ifndef NDEBUG
1094     if (EI_TPI(inputrec->eI) &&
1095         inputrec->cutoff_scheme == ecutsVERLET)
1096     {
1097         gmx_feenableexcept();
1098     }
1099 #endif
1100
1101     if (bUseGPU)
1102     {
1103         /* Select GPU id's to use */
1104         gmx_select_gpu_ids(fplog, cr, &hwinfo->gpu_info, bForceUseGPU,
1105                            &hw_opt->gpu_opt);
1106     }
1107     else
1108     {
1109         /* Ignore (potentially) manually selected GPUs */
1110         hw_opt->gpu_opt.n_dev_use = 0;
1111     }
1112
1113     /* check consistency across ranks of things like SIMD
1114      * support and number of GPUs selected */
1115     gmx_check_hw_runconf_consistency(fplog, hwinfo, cr, hw_opt, bUseGPU);
1116
1117     /* Now that we know the setup is consistent, check for efficiency */
1118     check_resource_division_efficiency(hwinfo, hw_opt, Flags & MD_NTOMPSET,
1119                                        cr, fplog);
1120
1121     if (DOMAINDECOMP(cr))
1122     {
1123         /* When we share GPUs over ranks, we need to know this for the DLB */
1124         dd_setup_dlb_resource_sharing(cr, hwinfo, hw_opt);
1125     }
1126
1127     /* getting number of PP/PME threads
1128        PME: env variable should be read only on one node to make sure it is
1129        identical everywhere;
1130      */
1131     nthreads_pme = gmx_omp_nthreads_get(emntPME);
1132
1133     wcycle = wallcycle_init(fplog, resetstep, cr);
1134
1135     if (PAR(cr))
1136     {
1137         /* Master synchronizes its value of reset_counters with all nodes
1138          * including PME only nodes */
1139         reset_counters = wcycle_get_reset_counters(wcycle);
1140         gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1141         wcycle_set_reset_counters(wcycle, reset_counters);
1142     }
1143
1144     snew(nrnb, 1);
1145     if (cr->duty & DUTY_PP)
1146     {
1147         bcast_state(cr, state);
1148
1149         /* Initiate forcerecord */
1150         fr          = mk_forcerec();
1151         fr->hwinfo  = hwinfo;
1152         fr->gpu_opt = &hw_opt->gpu_opt;
1153         init_forcerec(fplog, fr, fcd, inputrec, mtop, cr, box,
1154                       opt2fn("-table", nfile, fnm),
1155                       opt2fn("-tablep", nfile, fnm),
1156                       getFilenm("-tableb", nfile, fnm),
1157                       nbpu_opt,
1158                       FALSE,
1159                       pforce);
1160
1161         /* Initialize QM-MM */
1162         if (fr->bQMMM)
1163         {
1164             init_QMMMrec(cr, mtop, inputrec, fr);
1165         }
1166
1167         /* Initialize the mdatoms structure.
1168          * mdatoms is not filled with atom data,
1169          * as this can not be done now with domain decomposition.
1170          */
1171         mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1172
1173         /* Initialize the virtual site communication */
1174         vsite = init_vsite(mtop, cr, FALSE);
1175
1176         calc_shifts(box, fr->shift_vec);
1177
1178         /* With periodic molecules the charge groups should be whole at start up
1179          * and the virtual sites should not be far from their proper positions.
1180          */
1181         if (!inputrec->bContinuation && MASTER(cr) &&
1182             !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1183         {
1184             /* Make molecules whole at start of run */
1185             if (fr->ePBC != epbcNONE)
1186             {
1187                 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, state->x);
1188             }
1189             if (vsite)
1190             {
1191                 /* Correct initial vsite positions are required
1192                  * for the initial distribution in the domain decomposition
1193                  * and for the initial shell prediction.
1194                  */
1195                 construct_vsites_mtop(vsite, mtop, state->x);
1196             }
1197         }
1198
1199         if (EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype))
1200         {
1201             ewaldcoeff_q  = fr->ewaldcoeff_q;
1202             ewaldcoeff_lj = fr->ewaldcoeff_lj;
1203             pmedata       = &fr->pmedata;
1204         }
1205         else
1206         {
1207             pmedata = NULL;
1208         }
1209     }
1210     else
1211     {
1212         /* This is a PME only node */
1213
1214         /* We don't need the state */
1215         done_state(state);
1216
1217         ewaldcoeff_q  = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1218         ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1219         snew(pmedata, 1);
1220     }
1221
1222     if (hw_opt->thread_affinity != threadaffOFF)
1223     {
1224         /* Before setting affinity, check whether the affinity has changed
1225          * - which indicates that probably the OpenMP library has changed it
1226          * since we first checked).
1227          */
1228         gmx_check_thread_affinity_set(fplog, cr,
1229                                       hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1230
1231         /* Set the CPU affinity */
1232         gmx_set_thread_affinity(fplog, cr, hw_opt, hwinfo);
1233     }
1234
1235     /* Initiate PME if necessary,
1236      * either on all nodes or on dedicated PME nodes only. */
1237     if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1238     {
1239         if (mdatoms)
1240         {
1241             nChargePerturbed = mdatoms->nChargePerturbed;
1242             if (EVDW_PME(inputrec->vdwtype))
1243             {
1244                 nTypePerturbed   = mdatoms->nTypePerturbed;
1245             }
1246         }
1247         if (cr->npmenodes > 0)
1248         {
1249             /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1250             gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1251             gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1252         }
1253
1254         if (cr->duty & DUTY_PME)
1255         {
1256             status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1257                                   mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1258                                   (Flags & MD_REPRODUCIBLE), nthreads_pme);
1259             if (status != 0)
1260             {
1261                 gmx_fatal(FARGS, "Error %d initializing PME", status);
1262             }
1263         }
1264     }
1265
1266
1267     if (EI_DYNAMICS(inputrec->eI))
1268     {
1269         /* Turn on signal handling on all nodes */
1270         /*
1271          * (A user signal from the PME nodes (if any)
1272          * is communicated to the PP nodes.
1273          */
1274         signal_handler_install();
1275     }
1276
1277     if (cr->duty & DUTY_PP)
1278     {
1279         /* Assumes uniform use of the number of OpenMP threads */
1280         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1281
1282         if (inputrec->bPull)
1283         {
1284             /* Initialize pull code */
1285             inputrec->pull_work =
1286                 init_pull(fplog, inputrec->pull, inputrec, nfile, fnm,
1287                           mtop, cr, oenv, inputrec->fepvals->init_lambda,
1288                           EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1289         }
1290
1291         if (inputrec->bRot)
1292         {
1293             /* Initialize enforced rotation code */
1294             init_rot(fplog, inputrec, nfile, fnm, cr, state->x, state->box, mtop, oenv,
1295                      bVerbose, Flags);
1296         }
1297
1298         constr = init_constraints(fplog, mtop, inputrec, ed, state, cr);
1299
1300         if (DOMAINDECOMP(cr))
1301         {
1302             GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1303             /* This call is not included in init_domain_decomposition mainly
1304              * because fr->cginfo_mb is set later.
1305              */
1306             dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1307                             Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1308         }
1309
1310         /* Now do whatever the user wants us to do (how flexible...) */
1311         my_integrator(inputrec->eI) (fplog, cr, nfile, fnm,
1312                                      oenv, bVerbose,
1313                                      nstglobalcomm,
1314                                      vsite, constr,
1315                                      nstepout, inputrec, mtop,
1316                                      fcd, state,
1317                                      mdatoms, nrnb, wcycle, ed, fr,
1318                                      repl_ex_nst, repl_ex_nex, repl_ex_seed,
1319                                      cpt_period, max_hours,
1320                                      imdport,
1321                                      Flags,
1322                                      walltime_accounting);
1323
1324         if (inputrec->bRot)
1325         {
1326             finish_rot(inputrec->rot);
1327         }
1328
1329         if (inputrec->bPull)
1330         {
1331             finish_pull(inputrec->pull_work);
1332         }
1333
1334     }
1335     else
1336     {
1337         GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1338         /* do PME only */
1339         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1340         gmx_pmeonly(*pmedata, cr, nrnb, wcycle, walltime_accounting, ewaldcoeff_q, ewaldcoeff_lj, inputrec);
1341     }
1342
1343     wallcycle_stop(wcycle, ewcRUN);
1344
1345     /* Finish up, write some stuff
1346      * if rerunMD, don't write last frame again
1347      */
1348     finish_run(fplog, cr,
1349                inputrec, nrnb, wcycle, walltime_accounting,
1350                fr ? fr->nbv : NULL,
1351                EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1352
1353
1354     /* Free GPU memory and context */
1355     free_gpu_resources(fr, cr, &hwinfo->gpu_info, fr ? fr->gpu_opt : NULL);
1356
1357     gmx_hardware_info_free(hwinfo);
1358
1359     /* Does what it says */
1360     print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1361     walltime_accounting_destroy(walltime_accounting);
1362
1363     /* Close logfile already here if we were appending to it */
1364     if (MASTER(cr) && (Flags & MD_APPENDFILES))
1365     {
1366         gmx_log_close(fplog);
1367     }
1368
1369     rc = (int)gmx_get_stop_condition();
1370
1371     done_ed(&ed);
1372
1373 #if GMX_THREAD_MPI
1374     /* we need to join all threads. The sub-threads join when they
1375        exit this function, but the master thread needs to be told to
1376        wait for that. */
1377     if (PAR(cr) && MASTER(cr))
1378     {
1379         tMPI_Finalize();
1380     }
1381 #endif
1382
1383     return rc;
1384 }
1385
1386 } // namespace gmx