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