Cleaned up high-level boolean variable naming
[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,2017, 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/ewald/pme.h"
61 #include "gromacs/fileio/checkpoint.h"
62 #include "gromacs/fileio/oenv.h"
63 #include "gromacs/fileio/tpxio.h"
64 #include "gromacs/gmxlib/network.h"
65 #include "gromacs/gpu_utils/gpu_utils.h"
66 #include "gromacs/hardware/cpuinfo.h"
67 #include "gromacs/hardware/detecthardware.h"
68 #include "gromacs/hardware/hardwareassign.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/mdmodules.h"
92 #include "gromacs/mdrunutility/threadaffinity.h"
93 #include "gromacs/mdtypes/commrec.h"
94 #include "gromacs/mdtypes/edsamhistory.h"
95 #include "gromacs/mdtypes/energyhistory.h"
96 #include "gromacs/mdtypes/inputrec.h"
97 #include "gromacs/mdtypes/md_enums.h"
98 #include "gromacs/mdtypes/observableshistory.h"
99 #include "gromacs/mdtypes/state.h"
100 #include "gromacs/mdtypes/swaphistory.h"
101 #include "gromacs/pbcutil/pbc.h"
102 #include "gromacs/pulling/pull.h"
103 #include "gromacs/pulling/pull_rotation.h"
104 #include "gromacs/timing/wallcycle.h"
105 #include "gromacs/topology/mtop_util.h"
106 #include "gromacs/trajectory/trajectoryframe.h"
107 #include "gromacs/utility/cstringutil.h"
108 #include "gromacs/utility/exceptions.h"
109 #include "gromacs/utility/fatalerror.h"
110 #include "gromacs/utility/filestream.h"
111 #include "gromacs/utility/gmxassert.h"
112 #include "gromacs/utility/gmxmpi.h"
113 #include "gromacs/utility/logger.h"
114 #include "gromacs/utility/loggerbuilder.h"
115 #include "gromacs/utility/pleasecite.h"
116 #include "gromacs/utility/programcontext.h"
117 #include "gromacs/utility/smalloc.h"
118
119 #include "deform.h"
120 #include "md.h"
121 #include "membed.h"
122 #include "repl_ex.h"
123 #include "resource-division.h"
124
125 #ifdef GMX_FAHCORE
126 #include "corewrap.h"
127 #endif
128
129 //! First step used in pressure scaling
130 gmx_int64_t         deform_init_init_step_tpx;
131 //! Initial box for pressure scaling
132 matrix              deform_init_box_tpx;
133 //! MPI variable for use in pressure scaling
134 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
135
136 #if GMX_THREAD_MPI
137 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
138  * the number of threads will get lowered.
139  */
140 #define MIN_ATOMS_PER_MPI_THREAD    90
141 #define MIN_ATOMS_PER_GPU           900
142
143 struct mdrunner_arglist
144 {
145     gmx_hw_opt_t                     hw_opt;
146     FILE                            *fplog;
147     t_commrec                       *cr;
148     int                              nfile;
149     const t_filenm                  *fnm;
150     const gmx_output_env_t          *oenv;
151     gmx_bool                         bVerbose;
152     int                              nstglobalcomm;
153     ivec                             ddxyz;
154     int                              dd_rank_order;
155     int                              npme;
156     real                             rdd;
157     real                             rconstr;
158     const char                      *dddlb_opt;
159     real                             dlb_scale;
160     const char                      *ddcsx;
161     const char                      *ddcsy;
162     const char                      *ddcsz;
163     const char                      *nbpu_opt;
164     int                              nstlist_cmdline;
165     gmx_int64_t                      nsteps_cmdline;
166     int                              nstepout;
167     int                              resetstep;
168     int                              nmultisim;
169     const ReplicaExchangeParameters *replExParams;
170     real                             pforce;
171     real                             cpt_period;
172     real                             max_hours;
173     int                              imdport;
174     unsigned long                    Flags;
175 };
176
177
178 /* The function used for spawning threads. Extracts the mdrunner()
179    arguments from its one argument and calls mdrunner(), after making
180    a commrec. */
181 static void mdrunner_start_fn(void *arg)
182 {
183     try
184     {
185         struct mdrunner_arglist *mda = (struct mdrunner_arglist*)arg;
186         struct mdrunner_arglist  mc  = *mda; /* copy the arg list to make sure
187                                                 that it's thread-local. This doesn't
188                                                 copy pointed-to items, of course,
189                                                 but those are all const. */
190         t_commrec *cr;                       /* we need a local version of this */
191         FILE      *fplog = nullptr;
192         t_filenm  *fnm;
193
194         fnm = dup_tfn(mc.nfile, mc.fnm);
195
196         cr = reinitialize_commrec_for_this_thread(mc.cr);
197
198         if (MASTER(cr))
199         {
200             fplog = mc.fplog;
201         }
202
203         gmx::mdrunner(&mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
204                       mc.bVerbose, mc.nstglobalcomm,
205                       mc.ddxyz, mc.dd_rank_order, mc.npme, mc.rdd,
206                       mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
207                       mc.ddcsx, mc.ddcsy, mc.ddcsz,
208                       mc.nbpu_opt, mc.nstlist_cmdline,
209                       mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
210                       mc.nmultisim, *mc.replExParams, mc.pforce,
211                       mc.cpt_period, mc.max_hours, mc.imdport, mc.Flags);
212     }
213     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
214 }
215
216
217 /* called by mdrunner() to start a specific number of threads (including
218    the main thread) for thread-parallel runs. This in turn calls mdrunner()
219    for each thread.
220    All options besides nthreads are the same as for mdrunner(). */
221 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
222                                          FILE *fplog, t_commrec *cr, int nfile,
223                                          const t_filenm fnm[], const gmx_output_env_t *oenv, gmx_bool bVerbose,
224                                          int nstglobalcomm,
225                                          ivec ddxyz, int dd_rank_order, int npme,
226                                          real rdd, real rconstr,
227                                          const char *dddlb_opt, real dlb_scale,
228                                          const char *ddcsx, const char *ddcsy, const char *ddcsz,
229                                          const char *nbpu_opt, int nstlist_cmdline,
230                                          gmx_int64_t nsteps_cmdline,
231                                          int nstepout, int resetstep,
232                                          int nmultisim,
233                                          const ReplicaExchangeParameters &replExParams,
234                                          real pforce, real cpt_period, real max_hours,
235                                          unsigned long Flags)
236 {
237     int                      ret;
238     struct mdrunner_arglist *mda;
239     t_commrec               *crn; /* the new commrec */
240     t_filenm                *fnmn;
241
242     /* first check whether we even need to start tMPI */
243     if (hw_opt->nthreads_tmpi < 2)
244     {
245         return cr;
246     }
247
248     /* a few small, one-time, almost unavoidable memory leaks: */
249     snew(mda, 1);
250     fnmn = dup_tfn(nfile, fnm);
251
252     /* fill the data structure to pass as void pointer to thread start fn */
253     /* hw_opt contains pointers, which should all be NULL at this stage */
254     mda->hw_opt          = *hw_opt;
255     mda->fplog           = fplog;
256     mda->cr              = cr;
257     mda->nfile           = nfile;
258     mda->fnm             = fnmn;
259     mda->oenv            = oenv;
260     mda->bVerbose        = bVerbose;
261     mda->nstglobalcomm   = nstglobalcomm;
262     mda->ddxyz[XX]       = ddxyz[XX];
263     mda->ddxyz[YY]       = ddxyz[YY];
264     mda->ddxyz[ZZ]       = ddxyz[ZZ];
265     mda->dd_rank_order   = dd_rank_order;
266     mda->npme            = npme;
267     mda->rdd             = rdd;
268     mda->rconstr         = rconstr;
269     mda->dddlb_opt       = dddlb_opt;
270     mda->dlb_scale       = dlb_scale;
271     mda->ddcsx           = ddcsx;
272     mda->ddcsy           = ddcsy;
273     mda->ddcsz           = ddcsz;
274     mda->nbpu_opt        = nbpu_opt;
275     mda->nstlist_cmdline = nstlist_cmdline;
276     mda->nsteps_cmdline  = nsteps_cmdline;
277     mda->nstepout        = nstepout;
278     mda->resetstep       = resetstep;
279     mda->nmultisim       = nmultisim;
280     mda->replExParams    = &replExParams;
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 nullptr;
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                              bool makeGpuPairList, 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 != nullptr && makeGpuPairList && 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 != nullptr)
390         {
391             fprintf(fp, "%s\n", nve_err);
392         }
393
394         return;
395     }
396
397     if (ir->verletbuf_tol == 0 && makeGpuPairList)
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 != nullptr)
409         {
410             fprintf(fp, "%s\n", vbd_err);
411         }
412
413         return;
414     }
415
416     if (makeGpuPairList)
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, makeGpuPairList, &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, nullptr,
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, nullptr, &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 != nullptr)
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 != nullptr)
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                                   bool                            makeGpuPairList,
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, makeGpuPairList, &ls);
573
574         calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, nullptr, &rlist_new);
575
576         if (rlist_new != ir->rlist)
577         {
578             if (fplog != nullptr)
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, makeGpuPairList, 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 //! Halt the run if there are inconsistences between user choices to run with GPUs and/or hardware detection.
644 static void exitIfCannotForceGpuRun(bool requirePhysicalGpu,
645                                     bool useVerletScheme,
646                                     bool compatibleGpusFound)
647 {
648     /* Was GPU acceleration either explicitly (-nb gpu) or implicitly
649      * (gpu ID passed) requested? */
650     if (!requirePhysicalGpu)
651     {
652         return;
653     }
654
655     if (GMX_GPU == GMX_GPU_NONE)
656     {
657         gmx_fatal(FARGS, "GPU acceleration requested, but %s was compiled without GPU support!",
658                   gmx::getProgramContext().displayName());
659     }
660
661     if (!useVerletScheme)
662     {
663         gmx_fatal(FARGS, "GPU acceleration requested, but can't be used without cutoff-scheme=Verlet");
664     }
665
666     if (!compatibleGpusFound)
667     {
668         gmx_fatal(FARGS, "GPU acceleration requested, but no compatible GPUs were detected.");
669     }
670 }
671
672 /*! \brief Return whether GPU acceleration is useful with the given settings.
673  *
674  * If not, logs a message about falling back to CPU code. */
675 static bool gpuAccelerationIsUseful(const MDLogger   &mdlog,
676                                     const t_inputrec *ir,
677                                     bool              doRerun)
678 {
679     if (doRerun && ir->opts.ngener > 1)
680     {
681         /* Rerun execution time is dominated by I/O and pair search,
682          * so GPUs are not very useful, plus they do not support more
683          * than one energy group. If the user requested GPUs
684          * explicitly, a fatal error is given later.  With non-reruns,
685          * we fall back to a single whole-of system energy group
686          * (which runs much faster than a multiple-energy-groups
687          * implementation would), and issue a note in the .log
688          * file. Users can re-run if they want the information. */
689         GMX_LOG(mdlog.warning).asParagraph().appendText("Multiple energy groups is not implemented for GPUs, so is not useful for this rerun, so falling back to the CPU");
690         return false;
691     }
692
693     return true;
694 }
695
696 //! \brief Return the correct integrator function.
697 static integrator_t *my_integrator(unsigned int ei)
698 {
699     switch (ei)
700     {
701         case eiMD:
702         case eiBD:
703         case eiSD1:
704         case eiVV:
705         case eiVVAK:
706             if (!EI_DYNAMICS(ei))
707             {
708                 GMX_THROW(APIError("do_md integrator would be called for a non-dynamical integrator"));
709             }
710             return do_md;
711         case eiSteep:
712             return do_steep;
713         case eiCG:
714             return do_cg;
715         case eiNM:
716             return do_nm;
717         case eiLBFGS:
718             return do_lbfgs;
719         case eiTPI:
720         case eiTPIC:
721             if (!EI_TPI(ei))
722             {
723                 GMX_THROW(APIError("do_tpi integrator would be called for a non-TPI integrator"));
724             }
725             return do_tpi;
726         case eiSD2_REMOVED:
727             GMX_THROW(NotImplementedError("SD2 integrator has been removed"));
728         default:
729             GMX_THROW(APIError("Non existing integrator selected"));
730     }
731 }
732
733 //! Initializes the logger for mdrun.
734 static gmx::LoggerOwner buildLogger(FILE *fplog, const t_commrec *cr)
735 {
736     gmx::LoggerBuilder builder;
737     if (fplog != nullptr)
738     {
739         builder.addTargetFile(gmx::MDLogger::LogLevel::Info, fplog);
740     }
741     if (cr == nullptr || SIMMASTER(cr))
742     {
743         builder.addTargetStream(gmx::MDLogger::LogLevel::Warning,
744                                 &gmx::TextOutputFile::standardError());
745     }
746     return builder.build();
747 }
748
749 int mdrunner(gmx_hw_opt_t *hw_opt,
750              FILE *fplog, t_commrec *cr, int nfile,
751              const t_filenm fnm[], const gmx_output_env_t *oenv, gmx_bool bVerbose,
752              int nstglobalcomm,
753              ivec ddxyz, int dd_rank_order, int npme, real rdd, real rconstr,
754              const char *dddlb_opt, real dlb_scale,
755              const char *ddcsx, const char *ddcsy, const char *ddcsz,
756              const char *nbpu_opt, int nstlist_cmdline,
757              gmx_int64_t nsteps_cmdline, int nstepout, int resetstep,
758              int gmx_unused nmultisim,
759              const ReplicaExchangeParameters &replExParams,
760              real pforce, real cpt_period, real max_hours,
761              int imdport, unsigned long Flags)
762 {
763     matrix                    box;
764     gmx_ddbox_t               ddbox = {0};
765     int                       npme_major, npme_minor;
766     t_nrnb                   *nrnb;
767     gmx_mtop_t               *mtop          = nullptr;
768     t_mdatoms                *mdatoms       = nullptr;
769     t_forcerec               *fr            = nullptr;
770     t_fcdata                 *fcd           = nullptr;
771     real                      ewaldcoeff_q  = 0;
772     real                      ewaldcoeff_lj = 0;
773     struct gmx_pme_t        **pmedata       = nullptr;
774     gmx_vsite_t              *vsite         = nullptr;
775     gmx_constr_t              constr;
776     int                       nChargePerturbed = -1, nTypePerturbed = 0, status;
777     gmx_wallcycle_t           wcycle;
778     gmx_walltime_accounting_t walltime_accounting = nullptr;
779     int                       rc;
780     gmx_int64_t               reset_counters;
781     int                       nthreads_pme = 1;
782     gmx_membed_t *            membed       = nullptr;
783     gmx_hw_info_t            *hwinfo       = nullptr;
784     /* The master rank decides early on bUseGPU and broadcasts this later */
785     gmx_bool                  bUseGPU            = FALSE;
786
787     /* CAUTION: threads may be started later on in this function, so
788        cr doesn't reflect the final parallel state right now */
789     gmx::MDModules mdModules;
790     t_inputrec     inputrecInstance;
791     t_inputrec    *inputrec = &inputrecInstance;
792     snew(mtop, 1);
793
794     if (Flags & MD_APPENDFILES)
795     {
796         fplog = nullptr;
797     }
798
799     bool doMembed = opt2bSet("-membed", nfile, fnm);
800     bool doRerun  = (Flags & MD_RERUN);
801
802     /* Handle GPU-related user options. Later, we check consistency
803      * with things like whether support is compiled, or tMPI thread
804      * count. */
805     gmx_parse_gpu_ids(&hw_opt->gpu_opt);
806     bool userSetGpuIds = hasUserSetGpuIds(&hw_opt->gpu_opt);
807     bool forceUseCpu   = (strncmp(nbpu_opt, "cpu", 3) == 0);
808     if (userSetGpuIds && forceUseCpu)
809     {
810         gmx_fatal(FARGS, "GPU IDs were specified, and short-ranged interactions were assigned to the CPU. Make no more than one of these choices.");
811     }
812     bool forceUsePhysicalGpu = (strncmp(nbpu_opt, "gpu", 3) == 0) || userSetGpuIds;
813     bool tryUsePhysicalGpu   = (strncmp(nbpu_opt, "auto", 4) == 0) || forceUsePhysicalGpu;
814
815     // Here we assume that SIMMASTER(cr) does not change even after the
816     // threads are started.
817     gmx::LoggerOwner logOwner(buildLogger(fplog, cr));
818     gmx::MDLogger    mdlog(logOwner.logger());
819
820     /* Detect hardware, gather information. This is an operation that is
821      * global for this process (MPI rank). */
822     hwinfo = gmx_detect_hardware(mdlog, cr, tryUsePhysicalGpu);
823
824     gmx_print_detected_hardware(fplog, cr, mdlog, hwinfo);
825
826     if (fplog != nullptr)
827     {
828         /* Print references after all software/hardware printing */
829         please_cite(fplog, "Abraham2015");
830         please_cite(fplog, "Pall2015");
831         please_cite(fplog, "Pronk2013");
832         please_cite(fplog, "Hess2008b");
833         please_cite(fplog, "Spoel2005a");
834         please_cite(fplog, "Lindahl2001a");
835         please_cite(fplog, "Berendsen95a");
836     }
837
838     std::unique_ptr<t_state> stateInstance = std::unique_ptr<t_state>(new t_state);
839     t_state *                state         = stateInstance.get();
840
841     if (SIMMASTER(cr))
842     {
843         /* Read (nearly) all data required for the simulation */
844         read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, state, mtop);
845
846         exitIfCannotForceGpuRun(forceUsePhysicalGpu,
847                                 inputrec->cutoff_scheme == ecutsVERLET,
848                                 compatibleGpusFound(hwinfo->gpu_info));
849
850         if (inputrec->cutoff_scheme == ecutsVERLET)
851         {
852             /* Here the master rank decides if all ranks will use GPUs */
853             bUseGPU = (compatibleGpusFound(hwinfo->gpu_info) ||
854                        getenv("GMX_EMULATE_GPU") != nullptr);
855
856             if (bUseGPU &&
857                 !gpuAccelerationIsUseful(mdlog, inputrec, doRerun))
858             {
859                 /* Fallback message printed by nbnxn_acceleration_supported */
860                 if (forceUsePhysicalGpu)
861                 {
862                     gmx_fatal(FARGS, "GPU acceleration requested, but not supported with the given input settings");
863                 }
864                 bUseGPU = FALSE;
865             }
866
867             prepare_verlet_scheme(fplog, cr,
868                                   inputrec, nstlist_cmdline, mtop, state->box,
869                                   bUseGPU, *hwinfo->cpuInfo);
870         }
871         else
872         {
873             if (nstlist_cmdline > 0)
874             {
875                 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
876             }
877
878             if (compatibleGpusFound(hwinfo->gpu_info))
879             {
880                 GMX_LOG(mdlog.warning).asParagraph().appendText(
881                         "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
882                         "      To use a GPU, set the mdp option: cutoff-scheme = Verlet");
883             }
884
885 #if GMX_TARGET_BGQ
886             md_print_warn(cr, fplog,
887                           "NOTE: There is no SIMD implementation of the group scheme kernels on\n"
888                           "      BlueGene/Q. You will observe better performance from using the\n"
889                           "      Verlet cut-off scheme.\n");
890 #endif
891         }
892     }
893
894     /* Check and update the hardware options for internal consistency */
895     check_and_update_hw_opt_1(hw_opt, cr, npme);
896
897     /* Early check for externally set process affinity. */
898     gmx_check_thread_affinity_set(mdlog, cr,
899                                   hw_opt, hwinfo->nthreads_hw_avail, FALSE);
900
901 #if GMX_THREAD_MPI
902     if (SIMMASTER(cr))
903     {
904         if (npme > 0 && hw_opt->nthreads_tmpi <= 0)
905         {
906             gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
907         }
908
909         /* Since the master knows the cut-off scheme, update hw_opt for this.
910          * This is done later for normal MPI and also once more with tMPI
911          * for all tMPI ranks.
912          */
913         check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
914
915         /* NOW the threads will be started: */
916         hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
917                                                  hw_opt,
918                                                  inputrec, mtop,
919                                                  mdlog, bUseGPU,
920                                                  doMembed);
921
922         if (hw_opt->nthreads_tmpi > 1)
923         {
924             t_commrec *cr_old       = cr;
925             /* now start the threads. */
926             cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
927                                         oenv, bVerbose, nstglobalcomm,
928                                         ddxyz, dd_rank_order, npme, rdd, rconstr,
929                                         dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
930                                         nbpu_opt, nstlist_cmdline,
931                                         nsteps_cmdline, nstepout, resetstep, nmultisim,
932                                         replExParams, pforce,
933                                         cpt_period, max_hours,
934                                         Flags);
935             /* the main thread continues here with a new cr. We don't deallocate
936                the old cr because other threads may still be reading it. */
937             if (cr == nullptr)
938             {
939                 gmx_comm("Failed to spawn threads");
940             }
941         }
942     }
943 #endif
944     /* END OF CAUTION: cr is now reliable */
945
946     if (PAR(cr))
947     {
948         /* now broadcast everything to the non-master nodes/threads: */
949         init_parallel(cr, inputrec, mtop);
950
951         /* The master rank decided on the use of GPUs,
952          * broadcast this information to all ranks.
953          */
954         gmx_bcast_sim(sizeof(bUseGPU), &bUseGPU, cr);
955     }
956     // TODO: Error handling
957     mdModules.assignOptionsToModules(*inputrec->params, nullptr);
958
959     if (fplog != nullptr)
960     {
961         pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
962         fprintf(fplog, "\n");
963     }
964
965     /* now make sure the state is initialized and propagated */
966     set_state_entries(state, inputrec);
967
968     /* A parallel command line option consistency check that we can
969        only do after any threads have started. */
970     if (!PAR(cr) &&
971         (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || npme > 0))
972     {
973         gmx_fatal(FARGS,
974                   "The -dd or -npme option request a parallel simulation, "
975 #if !GMX_MPI
976                   "but %s was compiled without threads or MPI enabled"
977 #else
978 #if GMX_THREAD_MPI
979                   "but the number of MPI-threads (option -ntmpi) is not set or is 1"
980 #else
981                   "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
982 #endif
983 #endif
984                   , output_env_get_program_display_name(oenv)
985                   );
986     }
987
988     if (doRerun &&
989         (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
990     {
991         gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
992     }
993
994     if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
995     {
996         gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
997     }
998
999     if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
1000     {
1001         if (npme > 0)
1002         {
1003             gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
1004                                  "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
1005         }
1006
1007         npme = 0;
1008     }
1009
1010     if (bUseGPU && npme < 0)
1011     {
1012         /* With GPUs we don't automatically use PME-only ranks. PME ranks can
1013          * improve performance with many threads per GPU, since our OpenMP
1014          * scaling is bad, but it's difficult to automate the setup.
1015          */
1016         npme = 0;
1017     }
1018
1019 #ifdef GMX_FAHCORE
1020     if (MASTER(cr))
1021     {
1022         fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1023     }
1024 #endif
1025
1026     /* NMR restraints must be initialized before load_checkpoint,
1027      * since with time averaging the history is added to t_state.
1028      * For proper consistency check we therefore need to extend
1029      * t_state here.
1030      * So the PME-only nodes (if present) will also initialize
1031      * the distance restraints.
1032      */
1033     snew(fcd, 1);
1034
1035     /* This needs to be called before read_checkpoint to extend the state */
1036     init_disres(fplog, mtop, inputrec, cr, fcd, state, replExParams.exchangeInterval > 0);
1037
1038     init_orires(fplog, mtop, as_rvec_array(state->x.data()), inputrec, cr, &(fcd->orires),
1039                 state);
1040
1041     if (inputrecDeform(inputrec))
1042     {
1043         /* Store the deform reference box before reading the checkpoint */
1044         if (SIMMASTER(cr))
1045         {
1046             copy_mat(state->box, box);
1047         }
1048         if (PAR(cr))
1049         {
1050             gmx_bcast(sizeof(box), box, cr);
1051         }
1052         /* Because we do not have the update struct available yet
1053          * in which the reference values should be stored,
1054          * we store them temporarily in static variables.
1055          * This should be thread safe, since they are only written once
1056          * and with identical values.
1057          */
1058         tMPI_Thread_mutex_lock(&deform_init_box_mutex);
1059         deform_init_init_step_tpx = inputrec->init_step;
1060         copy_mat(box, deform_init_box_tpx);
1061         tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
1062     }
1063
1064     ObservablesHistory observablesHistory = {};
1065
1066     if (Flags & MD_STARTFROMCPT)
1067     {
1068         /* Check if checkpoint file exists before doing continuation.
1069          * This way we can use identical input options for the first and subsequent runs...
1070          */
1071         gmx_bool bReadEkin;
1072
1073         load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
1074                         cr, ddxyz, &npme,
1075                         inputrec, state, &bReadEkin, &observablesHistory,
1076                         (Flags & MD_APPENDFILES),
1077                         (Flags & MD_APPENDFILESSET),
1078                         (Flags & MD_REPRODUCIBLE));
1079
1080         if (bReadEkin)
1081         {
1082             Flags |= MD_READ_EKIN;
1083         }
1084     }
1085
1086     if (SIMMASTER(cr) && (Flags & MD_APPENDFILES))
1087     {
1088         gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
1089                      Flags, &fplog);
1090         logOwner = buildLogger(fplog, nullptr);
1091         mdlog    = logOwner.logger();
1092     }
1093
1094     /* override nsteps with value from cmdline */
1095     override_nsteps_cmdline(mdlog, nsteps_cmdline, inputrec);
1096
1097     if (SIMMASTER(cr))
1098     {
1099         copy_mat(state->box, box);
1100     }
1101
1102     if (PAR(cr))
1103     {
1104         gmx_bcast(sizeof(box), box, cr);
1105     }
1106
1107     if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
1108                      inputrec->eI == eiNM))
1109     {
1110         cr->dd = init_domain_decomposition(fplog, cr, Flags, ddxyz, npme,
1111                                            dd_rank_order,
1112                                            rdd, rconstr,
1113                                            dddlb_opt, dlb_scale,
1114                                            ddcsx, ddcsy, ddcsz,
1115                                            mtop, inputrec,
1116                                            box, as_rvec_array(state->x.data()),
1117                                            &ddbox, &npme_major, &npme_minor);
1118     }
1119     else
1120     {
1121         /* PME, if used, is done on all nodes with 1D decomposition */
1122         cr->npmenodes = 0;
1123         cr->duty      = (DUTY_PP | DUTY_PME);
1124         npme_major    = 1;
1125         npme_minor    = 1;
1126
1127         if (inputrec->ePBC == epbcSCREW)
1128         {
1129             gmx_fatal(FARGS,
1130                       "pbc=%s is only implemented with domain decomposition",
1131                       epbc_names[inputrec->ePBC]);
1132         }
1133     }
1134
1135     if (PAR(cr))
1136     {
1137         /* After possible communicator splitting in make_dd_communicators.
1138          * we can set up the intra/inter node communication.
1139          */
1140         gmx_setup_nodecomm(fplog, cr);
1141     }
1142
1143     /* Initialize per-physical-node MPI process/thread ID and counters. */
1144     gmx_init_intranode_counters(cr);
1145 #if GMX_MPI
1146     if (MULTISIM(cr))
1147     {
1148         GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
1149                 "This is simulation %d out of %d running as a composite GROMACS\n"
1150                 "multi-simulation job. Setup for this simulation:\n",
1151                 cr->ms->sim, cr->ms->nsim);
1152     }
1153     GMX_LOG(mdlog.warning).appendTextFormatted(
1154             "Using %d MPI %s\n",
1155             cr->nnodes,
1156 #if GMX_THREAD_MPI
1157             cr->nnodes == 1 ? "thread" : "threads"
1158 #else
1159             cr->nnodes == 1 ? "process" : "processes"
1160 #endif
1161             );
1162     fflush(stderr);
1163 #endif
1164
1165     /* Check and update hw_opt for the cut-off scheme */
1166     check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1167
1168     /* Check and update hw_opt for the number of MPI ranks */
1169     check_and_update_hw_opt_3(hw_opt);
1170
1171     gmx_omp_nthreads_init(mdlog, cr,
1172                           hwinfo->nthreads_hw_avail,
1173                           hw_opt->nthreads_omp,
1174                           hw_opt->nthreads_omp_pme,
1175                           (cr->duty & DUTY_PP) == 0,
1176                           inputrec->cutoff_scheme == ecutsVERLET);
1177
1178 #ifndef NDEBUG
1179     if (EI_TPI(inputrec->eI) &&
1180         inputrec->cutoff_scheme == ecutsVERLET)
1181     {
1182         gmx_feenableexcept();
1183     }
1184 #endif
1185
1186     if (bUseGPU)
1187     {
1188         /* Select GPU id's to use */
1189         gmx_select_rank_gpu_ids(mdlog, cr, &hwinfo->gpu_info,
1190                                 userSetGpuIds, &hw_opt->gpu_opt);
1191     }
1192     else
1193     {
1194         /* Ignore (potentially) manually selected GPUs */
1195         hw_opt->gpu_opt.n_dev_use = 0;
1196     }
1197
1198     /* check consistency across ranks of things like SIMD
1199      * support and number of GPUs selected */
1200     gmx_check_hw_runconf_consistency(mdlog, hwinfo, cr, hw_opt, userSetGpuIds, bUseGPU);
1201
1202     /* Now that we know the setup is consistent, check for efficiency */
1203     check_resource_division_efficiency(hwinfo, hw_opt, hw_opt->gpu_opt.n_dev_use, Flags & MD_NTOMPSET,
1204                                        cr, mdlog);
1205
1206     if (DOMAINDECOMP(cr))
1207     {
1208         /* When we share GPUs over ranks, we need to know this for the DLB */
1209         dd_setup_dlb_resource_sharing(cr, hwinfo, hw_opt);
1210     }
1211
1212     /* getting number of PP/PME threads
1213        PME: env variable should be read only on one node to make sure it is
1214        identical everywhere;
1215      */
1216     nthreads_pme = gmx_omp_nthreads_get(emntPME);
1217
1218     wcycle = wallcycle_init(fplog, resetstep, cr);
1219
1220     if (PAR(cr))
1221     {
1222         /* Master synchronizes its value of reset_counters with all nodes
1223          * including PME only nodes */
1224         reset_counters = wcycle_get_reset_counters(wcycle);
1225         gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1226         wcycle_set_reset_counters(wcycle, reset_counters);
1227     }
1228
1229     // Membrane embedding must be initialized before we call init_forcerec()
1230     if (doMembed)
1231     {
1232         if (MASTER(cr))
1233         {
1234             fprintf(stderr, "Initializing membed");
1235         }
1236         /* Note that membed cannot work in parallel because mtop is
1237          * changed here. Fix this if we ever want to make it run with
1238          * multiple ranks. */
1239         membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
1240     }
1241
1242     snew(nrnb, 1);
1243     if (cr->duty & DUTY_PP)
1244     {
1245         bcast_state(cr, state);
1246
1247         /* Initiate forcerecord */
1248         fr          = mk_forcerec();
1249         fr->hwinfo  = hwinfo;
1250         fr->gpu_opt = &hw_opt->gpu_opt;
1251         init_forcerec(fplog, mdlog, fr, fcd, mdModules.forceProvider(),
1252                       inputrec, mtop, cr, box,
1253                       opt2fn("-table", nfile, fnm),
1254                       opt2fn("-tablep", nfile, fnm),
1255                       getFilenm("-tableb", nfile, fnm),
1256                       nbpu_opt,
1257                       FALSE,
1258                       pforce);
1259
1260         /* Initialize QM-MM */
1261         if (fr->bQMMM)
1262         {
1263             init_QMMMrec(cr, mtop, inputrec, fr);
1264         }
1265
1266         /* Initialize the mdatoms structure.
1267          * mdatoms is not filled with atom data,
1268          * as this can not be done now with domain decomposition.
1269          */
1270         mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1271
1272         /* Initialize the virtual site communication */
1273         vsite = init_vsite(mtop, cr, FALSE);
1274
1275         calc_shifts(box, fr->shift_vec);
1276
1277         /* With periodic molecules the charge groups should be whole at start up
1278          * and the virtual sites should not be far from their proper positions.
1279          */
1280         if (!inputrec->bContinuation && MASTER(cr) &&
1281             !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1282         {
1283             /* Make molecules whole at start of run */
1284             if (fr->ePBC != epbcNONE)
1285             {
1286                 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, as_rvec_array(state->x.data()));
1287             }
1288             if (vsite)
1289             {
1290                 /* Correct initial vsite positions are required
1291                  * for the initial distribution in the domain decomposition
1292                  * and for the initial shell prediction.
1293                  */
1294                 construct_vsites_mtop(vsite, mtop, as_rvec_array(state->x.data()));
1295             }
1296         }
1297
1298         if (EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype))
1299         {
1300             ewaldcoeff_q  = fr->ewaldcoeff_q;
1301             ewaldcoeff_lj = fr->ewaldcoeff_lj;
1302             pmedata       = &fr->pmedata;
1303         }
1304         else
1305         {
1306             pmedata = nullptr;
1307         }
1308     }
1309     else
1310     {
1311         /* This is a PME only node */
1312
1313         /* We don't need the state */
1314         stateInstance.reset();
1315         state         = nullptr;
1316
1317         ewaldcoeff_q  = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1318         ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1319         snew(pmedata, 1);
1320     }
1321
1322     if (hw_opt->thread_affinity != threadaffOFF)
1323     {
1324         /* Before setting affinity, check whether the affinity has changed
1325          * - which indicates that probably the OpenMP library has changed it
1326          * since we first checked).
1327          */
1328         gmx_check_thread_affinity_set(mdlog, cr,
1329                                       hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1330
1331         int nthread_local;
1332         /* threads on this MPI process or TMPI thread */
1333         if (cr->duty & DUTY_PP)
1334         {
1335             nthread_local = gmx_omp_nthreads_get(emntNonbonded);
1336         }
1337         else
1338         {
1339             nthread_local = gmx_omp_nthreads_get(emntPME);
1340         }
1341
1342         /* Set the CPU affinity */
1343         gmx_set_thread_affinity(mdlog, cr, hw_opt, *hwinfo->hardwareTopology,
1344                                 nthread_local, nullptr);
1345     }
1346
1347     /* Initiate PME if necessary,
1348      * either on all nodes or on dedicated PME nodes only. */
1349     if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1350     {
1351         if (mdatoms)
1352         {
1353             nChargePerturbed = mdatoms->nChargePerturbed;
1354             if (EVDW_PME(inputrec->vdwtype))
1355             {
1356                 nTypePerturbed   = mdatoms->nTypePerturbed;
1357             }
1358         }
1359         if (cr->npmenodes > 0)
1360         {
1361             /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1362             gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1363             gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1364         }
1365
1366         if (cr->duty & DUTY_PME)
1367         {
1368             try
1369             {
1370                 status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1371                                       mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1372                                       (Flags & MD_REPRODUCIBLE),
1373                                       ewaldcoeff_q, ewaldcoeff_lj,
1374                                       nthreads_pme);
1375             }
1376             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1377             if (status != 0)
1378             {
1379                 gmx_fatal(FARGS, "Error %d initializing PME", status);
1380             }
1381         }
1382     }
1383
1384
1385     if (EI_DYNAMICS(inputrec->eI))
1386     {
1387         /* Turn on signal handling on all nodes */
1388         /*
1389          * (A user signal from the PME nodes (if any)
1390          * is communicated to the PP nodes.
1391          */
1392         signal_handler_install();
1393     }
1394
1395     if (cr->duty & DUTY_PP)
1396     {
1397         /* Assumes uniform use of the number of OpenMP threads */
1398         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1399
1400         if (inputrec->bPull)
1401         {
1402             /* Initialize pull code */
1403             inputrec->pull_work =
1404                 init_pull(fplog, inputrec->pull, inputrec, nfile, fnm,
1405                           mtop, cr, oenv, inputrec->fepvals->init_lambda,
1406                           EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1407         }
1408
1409         if (inputrec->bRot)
1410         {
1411             /* Initialize enforced rotation code */
1412             init_rot(fplog, inputrec, nfile, fnm, cr, as_rvec_array(state->x.data()), state->box, mtop, oenv,
1413                      bVerbose, Flags);
1414         }
1415
1416         /* Let init_constraints know whether we have essential dynamics constraints.
1417          * TODO: inputrec should tell us whether we use an algorithm, not a file option or the checkpoint
1418          */
1419         bool doEdsam = (opt2fn_null("-ei", nfile, fnm) != nullptr || observablesHistory.edsamHistory);
1420
1421         constr = init_constraints(fplog, mtop, inputrec, doEdsam, cr);
1422
1423         if (DOMAINDECOMP(cr))
1424         {
1425             GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1426             /* This call is not included in init_domain_decomposition mainly
1427              * because fr->cginfo_mb is set later.
1428              */
1429             dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1430                             Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1431         }
1432
1433         /* Now do whatever the user wants us to do (how flexible...) */
1434         my_integrator(inputrec->eI) (fplog, cr, mdlog, nfile, fnm,
1435                                      oenv, bVerbose,
1436                                      nstglobalcomm,
1437                                      vsite, constr,
1438                                      nstepout, mdModules.outputProvider(),
1439                                      inputrec, mtop,
1440                                      fcd, state, &observablesHistory,
1441                                      mdatoms, nrnb, wcycle, fr,
1442                                      replExParams,
1443                                      membed,
1444                                      cpt_period, max_hours,
1445                                      imdport,
1446                                      Flags,
1447                                      walltime_accounting);
1448
1449         if (inputrec->bRot)
1450         {
1451             finish_rot(inputrec->rot);
1452         }
1453
1454         if (inputrec->bPull)
1455         {
1456             finish_pull(inputrec->pull_work);
1457         }
1458
1459     }
1460     else
1461     {
1462         GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1463         /* do PME only */
1464         walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1465         gmx_pmeonly(*pmedata, cr, nrnb, wcycle, walltime_accounting, ewaldcoeff_q, ewaldcoeff_lj, inputrec);
1466     }
1467
1468     wallcycle_stop(wcycle, ewcRUN);
1469
1470     /* Finish up, write some stuff
1471      * if rerunMD, don't write last frame again
1472      */
1473     finish_run(fplog, mdlog, cr,
1474                inputrec, nrnb, wcycle, walltime_accounting,
1475                fr ? fr->nbv : nullptr,
1476                EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1477
1478     // Free PME data
1479     if (pmedata)
1480     {
1481         gmx_pme_destroy(*pmedata); // TODO: pmedata is always a single element list, refactor
1482         pmedata = nullptr;
1483     }
1484
1485     /* Free GPU memory and context */
1486     free_gpu_resources(fr, cr, &hwinfo->gpu_info, fr ? fr->gpu_opt : nullptr);
1487
1488     if (doMembed)
1489     {
1490         free_membed(membed);
1491     }
1492
1493     gmx_hardware_info_free(hwinfo);
1494
1495     /* Does what it says */
1496     print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1497     walltime_accounting_destroy(walltime_accounting);
1498
1499     /* Close logfile already here if we were appending to it */
1500     if (MASTER(cr) && (Flags & MD_APPENDFILES))
1501     {
1502         gmx_log_close(fplog);
1503     }
1504
1505     rc = (int)gmx_get_stop_condition();
1506
1507 #if GMX_THREAD_MPI
1508     /* we need to join all threads. The sub-threads join when they
1509        exit this function, but the master thread needs to be told to
1510        wait for that. */
1511     if (PAR(cr) && MASTER(cr))
1512     {
1513         tMPI_Finalize();
1514     }
1515 #endif
1516
1517     return rc;
1518 }
1519
1520 } // namespace gmx