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