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