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