Move some mdrun code out of gmxlib/
[gromacs.git] / src / programs / mdrun / runner.cpp
blob7d240d248c6586bef53ff790d0c7f7d59dc64eed
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, 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/gmx_detect_hardware.h"
66 #include "gromacs/gmxlib/gmx_omp_nthreads.h"
67 #include "gromacs/gmxlib/main.h"
68 #include "gromacs/gmxlib/md_logging.h"
69 #include "gromacs/gmxlib/network.h"
70 #include "gromacs/gpu_utils/gpu_utils.h"
71 #include "gromacs/listed-forces/disre.h"
72 #include "gromacs/listed-forces/orires.h"
73 #include "gromacs/math/calculate-ewald-splitting-coefficient.h"
74 #include "gromacs/math/vec.h"
75 #include "gromacs/mdlib/calc_verletbuf.h"
76 #include "gromacs/mdlib/constr.h"
77 #include "gromacs/mdlib/force.h"
78 #include "gromacs/mdlib/forcerec.h"
79 #include "gromacs/mdlib/integrator.h"
80 #include "gromacs/mdlib/md_support.h"
81 #include "gromacs/mdlib/mdatoms.h"
82 #include "gromacs/mdlib/mdrun.h"
83 #include "gromacs/mdlib/minimize.h"
84 #include "gromacs/mdlib/nbnxn_search.h"
85 #include "gromacs/mdlib/qmmm.h"
86 #include "gromacs/mdlib/sighandler.h"
87 #include "gromacs/mdlib/tpi.h"
88 #include "gromacs/mdrunutility/threadaffinity.h"
89 #include "gromacs/mdtypes/inputrec.h"
90 #include "gromacs/mdtypes/md_enums.h"
91 #include "gromacs/mdtypes/state.h"
92 #include "gromacs/pbcutil/pbc.h"
93 #include "gromacs/pulling/pull.h"
94 #include "gromacs/pulling/pull_rotation.h"
95 #include "gromacs/swap/swapcoords.h"
96 #include "gromacs/timing/wallcycle.h"
97 #include "gromacs/topology/mtop_util.h"
98 #include "gromacs/trajectory/trajectoryframe.h"
99 #include "gromacs/utility/cstringutil.h"
100 #include "gromacs/utility/exceptions.h"
101 #include "gromacs/utility/fatalerror.h"
102 #include "gromacs/utility/gmxassert.h"
103 #include "gromacs/utility/gmxmpi.h"
104 #include "gromacs/utility/pleasecite.h"
105 #include "gromacs/utility/smalloc.h"
107 #include "deform.h"
108 #include "md.h"
109 #include "membed.h"
110 #include "repl_ex.h"
111 #include "resource-division.h"
113 #ifdef GMX_FAHCORE
114 #include "corewrap.h"
115 #endif
117 //! First step used in pressure scaling
118 gmx_int64_t deform_init_init_step_tpx;
119 //! Initial box for pressure scaling
120 matrix deform_init_box_tpx;
121 //! MPI variable for use in pressure scaling
122 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
124 #ifdef GMX_THREAD_MPI
125 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
126 * the number of threads will get lowered.
128 #define MIN_ATOMS_PER_MPI_THREAD 90
129 #define MIN_ATOMS_PER_GPU 900
131 struct mdrunner_arglist
133 gmx_hw_opt_t hw_opt;
134 FILE *fplog;
135 t_commrec *cr;
136 int nfile;
137 const t_filenm *fnm;
138 const gmx_output_env_t *oenv;
139 gmx_bool bVerbose;
140 int nstglobalcomm;
141 ivec ddxyz;
142 int dd_node_order;
143 real rdd;
144 real rconstr;
145 const char *dddlb_opt;
146 real dlb_scale;
147 const char *ddcsx;
148 const char *ddcsy;
149 const char *ddcsz;
150 const char *nbpu_opt;
151 int nstlist_cmdline;
152 gmx_int64_t nsteps_cmdline;
153 int nstepout;
154 int resetstep;
155 int nmultisim;
156 int repl_ex_nst;
157 int repl_ex_nex;
158 int repl_ex_seed;
159 real pforce;
160 real cpt_period;
161 real max_hours;
162 int imdport;
163 unsigned long Flags;
167 /* The function used for spawning threads. Extracts the mdrunner()
168 arguments from its one argument and calls mdrunner(), after making
169 a commrec. */
170 static void mdrunner_start_fn(void *arg)
174 struct mdrunner_arglist *mda = (struct mdrunner_arglist*)arg;
175 struct mdrunner_arglist mc = *mda; /* copy the arg list to make sure
176 that it's thread-local. This doesn't
177 copy pointed-to items, of course,
178 but those are all const. */
179 t_commrec *cr; /* we need a local version of this */
180 FILE *fplog = NULL;
181 t_filenm *fnm;
183 fnm = dup_tfn(mc.nfile, mc.fnm);
185 cr = reinitialize_commrec_for_this_thread(mc.cr);
187 if (MASTER(cr))
189 fplog = mc.fplog;
192 gmx::mdrunner(&mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
193 mc.bVerbose, mc.nstglobalcomm,
194 mc.ddxyz, mc.dd_node_order, mc.rdd,
195 mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
196 mc.ddcsx, mc.ddcsy, mc.ddcsz,
197 mc.nbpu_opt, mc.nstlist_cmdline,
198 mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
199 mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce,
200 mc.cpt_period, mc.max_hours, mc.imdport, mc.Flags);
202 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
206 /* called by mdrunner() to start a specific number of threads (including
207 the main thread) for thread-parallel runs. This in turn calls mdrunner()
208 for each thread.
209 All options besides nthreads are the same as for mdrunner(). */
210 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
211 FILE *fplog, t_commrec *cr, int nfile,
212 const t_filenm fnm[], const gmx_output_env_t *oenv, gmx_bool bVerbose,
213 int nstglobalcomm,
214 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
215 const char *dddlb_opt, real dlb_scale,
216 const char *ddcsx, const char *ddcsy, const char *ddcsz,
217 const char *nbpu_opt, int nstlist_cmdline,
218 gmx_int64_t nsteps_cmdline,
219 int nstepout, int resetstep,
220 int nmultisim, int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
221 real pforce, real cpt_period, real max_hours,
222 unsigned long Flags)
224 int ret;
225 struct mdrunner_arglist *mda;
226 t_commrec *crn; /* the new commrec */
227 t_filenm *fnmn;
229 /* first check whether we even need to start tMPI */
230 if (hw_opt->nthreads_tmpi < 2)
232 return cr;
235 /* a few small, one-time, almost unavoidable memory leaks: */
236 snew(mda, 1);
237 fnmn = dup_tfn(nfile, fnm);
239 /* fill the data structure to pass as void pointer to thread start fn */
240 /* hw_opt contains pointers, which should all be NULL at this stage */
241 mda->hw_opt = *hw_opt;
242 mda->fplog = fplog;
243 mda->cr = cr;
244 mda->nfile = nfile;
245 mda->fnm = fnmn;
246 mda->oenv = oenv;
247 mda->bVerbose = bVerbose;
248 mda->nstglobalcomm = nstglobalcomm;
249 mda->ddxyz[XX] = ddxyz[XX];
250 mda->ddxyz[YY] = ddxyz[YY];
251 mda->ddxyz[ZZ] = ddxyz[ZZ];
252 mda->dd_node_order = dd_node_order;
253 mda->rdd = rdd;
254 mda->rconstr = rconstr;
255 mda->dddlb_opt = dddlb_opt;
256 mda->dlb_scale = dlb_scale;
257 mda->ddcsx = ddcsx;
258 mda->ddcsy = ddcsy;
259 mda->ddcsz = ddcsz;
260 mda->nbpu_opt = nbpu_opt;
261 mda->nstlist_cmdline = nstlist_cmdline;
262 mda->nsteps_cmdline = nsteps_cmdline;
263 mda->nstepout = nstepout;
264 mda->resetstep = resetstep;
265 mda->nmultisim = nmultisim;
266 mda->repl_ex_nst = repl_ex_nst;
267 mda->repl_ex_nex = repl_ex_nex;
268 mda->repl_ex_seed = repl_ex_seed;
269 mda->pforce = pforce;
270 mda->cpt_period = cpt_period;
271 mda->max_hours = max_hours;
272 mda->Flags = Flags;
274 /* now spawn new threads that start mdrunner_start_fn(), while
275 the main thread returns, we set thread affinity later */
276 ret = tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi, TMPI_AFFINITY_NONE,
277 mdrunner_start_fn, (void*)(mda) );
278 if (ret != TMPI_SUCCESS)
280 return NULL;
283 crn = reinitialize_commrec_for_this_thread(cr);
284 return crn;
287 #endif /* GMX_THREAD_MPI */
290 /*! \brief Cost of non-bonded kernels
292 * We determine the extra cost of the non-bonded kernels compared to
293 * a reference nstlist value of 10 (which is the default in grompp).
295 static const int nbnxnReferenceNstlist = 10;
296 //! The values to try when switching
297 const int nstlist_try[] = { 20, 25, 40 };
298 //! Number of elements in the neighborsearch list trials.
299 #define NNSTL sizeof(nstlist_try)/sizeof(nstlist_try[0])
300 /* Increase nstlist until the non-bonded cost increases more than listfac_ok,
301 * but never more than listfac_max.
302 * A standard (protein+)water system at 300K with PME ewald_rtol=1e-5
303 * needs 1.28 at rcoulomb=0.9 and 1.24 at rcoulomb=1.0 to get to nstlist=40.
304 * Note that both CPU and GPU factors are conservative. Performance should
305 * not go down due to this tuning, except with a relatively slow GPU.
306 * On the other hand, at medium/high parallelization or with fast GPUs
307 * nstlist will not be increased enough to reach optimal performance.
309 /* CPU: pair-search is about a factor 1.5 slower than the non-bonded kernel */
310 //! Max OK performance ratio beween force calc and neighbor searching
311 static const float nbnxn_cpu_listfac_ok = 1.05;
312 //! Too high performance ratio beween force calc and neighbor searching
313 static const float nbnxn_cpu_listfac_max = 1.09;
314 /* GPU: pair-search is a factor 1.5-3 slower than the non-bonded kernel */
315 //! Max OK performance ratio beween force calc and neighbor searching
316 static const float nbnxn_gpu_listfac_ok = 1.20;
317 //! Too high performance ratio beween force calc and neighbor searching
318 static const float nbnxn_gpu_listfac_max = 1.30;
320 /*! \brief Try to increase nstlist when using the Verlet cut-off scheme */
321 static void increase_nstlist(FILE *fp, t_commrec *cr,
322 t_inputrec *ir, int nstlist_cmdline,
323 const gmx_mtop_t *mtop, matrix box,
324 gmx_bool bGPU)
326 float listfac_ok, listfac_max;
327 int nstlist_orig, nstlist_prev;
328 verletbuf_list_setup_t ls;
329 real rlistWithReferenceNstlist, rlist_inc, rlist_ok, rlist_max;
330 real rlist_new, rlist_prev;
331 size_t nstlist_ind = 0;
332 t_state state_tmp;
333 gmx_bool bBox, bDD, bCont;
334 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";
335 const char *nve_err = "Can not increase nstlist because an NVE ensemble is used";
336 const char *vbd_err = "Can not increase nstlist because verlet-buffer-tolerance is not set or used";
337 const char *box_err = "Can not increase nstlist because the box is too small";
338 const char *dd_err = "Can not increase nstlist because of domain decomposition limitations";
339 char buf[STRLEN];
340 const float oneThird = 1.0f / 3.0f;
342 if (nstlist_cmdline <= 0)
344 if (ir->nstlist == 1)
346 /* The user probably set nstlist=1 for a reason,
347 * don't mess with the settings.
349 return;
352 if (fp != NULL && bGPU && ir->nstlist < nstlist_try[0])
354 fprintf(fp, nstl_gpu, ir->nstlist);
356 nstlist_ind = 0;
357 while (nstlist_ind < NNSTL && ir->nstlist >= nstlist_try[nstlist_ind])
359 nstlist_ind++;
361 if (nstlist_ind == NNSTL)
363 /* There are no larger nstlist value to try */
364 return;
368 if (EI_MD(ir->eI) && ir->etc == etcNO)
370 if (MASTER(cr))
372 fprintf(stderr, "%s\n", nve_err);
374 if (fp != NULL)
376 fprintf(fp, "%s\n", nve_err);
379 return;
382 if (ir->verletbuf_tol == 0 && bGPU)
384 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");
387 if (ir->verletbuf_tol < 0)
389 if (MASTER(cr))
391 fprintf(stderr, "%s\n", vbd_err);
393 if (fp != NULL)
395 fprintf(fp, "%s\n", vbd_err);
398 return;
401 if (bGPU)
403 listfac_ok = nbnxn_gpu_listfac_ok;
404 listfac_max = nbnxn_gpu_listfac_max;
406 else
408 listfac_ok = nbnxn_cpu_listfac_ok;
409 listfac_max = nbnxn_cpu_listfac_max;
412 nstlist_orig = ir->nstlist;
413 if (nstlist_cmdline > 0)
415 if (fp)
417 sprintf(buf, "Getting nstlist=%d from command line option",
418 nstlist_cmdline);
420 ir->nstlist = nstlist_cmdline;
423 verletbuf_get_list_setup(TRUE, bGPU, &ls);
425 /* Allow rlist to make the list a given factor larger than the list
426 * would be with the reference value for nstlist (10).
428 nstlist_prev = ir->nstlist;
429 ir->nstlist = nbnxnReferenceNstlist;
430 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL,
431 &rlistWithReferenceNstlist);
432 ir->nstlist = nstlist_prev;
434 /* Determine the pair list size increase due to zero interactions */
435 rlist_inc = nbnxn_get_rlist_effective_inc(ls.cluster_size_j,
436 mtop->natoms/det(box));
437 rlist_ok = (rlistWithReferenceNstlist + rlist_inc)*pow(listfac_ok, oneThird) - rlist_inc;
438 rlist_max = (rlistWithReferenceNstlist + rlist_inc)*pow(listfac_max, oneThird) - rlist_inc;
439 if (debug)
441 fprintf(debug, "nstlist tuning: rlist_inc %.3f rlist_ok %.3f rlist_max %.3f\n",
442 rlist_inc, rlist_ok, rlist_max);
445 nstlist_prev = nstlist_orig;
446 rlist_prev = ir->rlist;
449 if (nstlist_cmdline <= 0)
451 ir->nstlist = nstlist_try[nstlist_ind];
454 /* Set the pair-list buffer size in ir */
455 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
457 /* Does rlist fit in the box? */
458 bBox = (sqr(rlist_new) < max_cutoff2(ir->ePBC, box));
459 bDD = TRUE;
460 if (bBox && DOMAINDECOMP(cr))
462 /* Check if rlist fits in the domain decomposition */
463 if (inputrec2nboundeddim(ir) < DIM)
465 gmx_incons("Changing nstlist with domain decomposition and unbounded dimensions is not implemented yet");
467 copy_mat(box, state_tmp.box);
468 bDD = change_dd_cutoff(cr, &state_tmp, ir, rlist_new);
471 if (debug)
473 fprintf(debug, "nstlist %d rlist %.3f bBox %d bDD %d\n",
474 ir->nstlist, rlist_new, bBox, bDD);
477 bCont = FALSE;
479 if (nstlist_cmdline <= 0)
481 if (bBox && bDD && rlist_new <= rlist_max)
483 /* Increase nstlist */
484 nstlist_prev = ir->nstlist;
485 rlist_prev = rlist_new;
486 bCont = (nstlist_ind+1 < NNSTL && rlist_new < rlist_ok);
488 else
490 /* Stick with the previous nstlist */
491 ir->nstlist = nstlist_prev;
492 rlist_new = rlist_prev;
493 bBox = TRUE;
494 bDD = TRUE;
498 nstlist_ind++;
500 while (bCont);
502 if (!bBox || !bDD)
504 gmx_warning(!bBox ? box_err : dd_err);
505 if (fp != NULL)
507 fprintf(fp, "\n%s\n", bBox ? box_err : dd_err);
509 ir->nstlist = nstlist_orig;
511 else if (ir->nstlist != nstlist_orig || rlist_new != ir->rlist)
513 sprintf(buf, "Changing nstlist from %d to %d, rlist from %g to %g",
514 nstlist_orig, ir->nstlist,
515 ir->rlist, rlist_new);
516 if (MASTER(cr))
518 fprintf(stderr, "%s\n\n", buf);
520 if (fp != NULL)
522 fprintf(fp, "%s\n\n", buf);
524 ir->rlist = rlist_new;
528 /*! \brief Initialize variables for Verlet scheme simulation */
529 static void prepare_verlet_scheme(FILE *fplog,
530 t_commrec *cr,
531 t_inputrec *ir,
532 int nstlist_cmdline,
533 const gmx_mtop_t *mtop,
534 matrix box,
535 gmx_bool bUseGPU)
537 /* For NVE simulations, we will retain the initial list buffer */
538 if (ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
540 /* Update the Verlet buffer size for the current run setup */
541 verletbuf_list_setup_t ls;
542 real rlist_new;
544 /* Here we assume SIMD-enabled kernels are being used. But as currently
545 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
546 * and 4x2 gives a larger buffer than 4x4, this is ok.
548 verletbuf_get_list_setup(TRUE, bUseGPU, &ls);
550 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
552 if (rlist_new != ir->rlist)
554 if (fplog != NULL)
556 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
557 ir->rlist, rlist_new,
558 ls.cluster_size_i, ls.cluster_size_j);
560 ir->rlist = rlist_new;
564 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
566 gmx_fatal(FARGS, "Can not set nstlist without %s",
567 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
570 if (EI_DYNAMICS(ir->eI))
572 /* Set or try nstlist values */
573 increase_nstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, bUseGPU);
577 /*! \brief Override the nslist value in inputrec
579 * with value passed on the command line (if any)
581 static void override_nsteps_cmdline(FILE *fplog,
582 gmx_int64_t nsteps_cmdline,
583 t_inputrec *ir,
584 const t_commrec *cr)
586 assert(ir);
587 assert(cr);
589 /* override with anything else than the default -2 */
590 if (nsteps_cmdline > -2)
592 char sbuf_steps[STEPSTRSIZE];
593 char sbuf_msg[STRLEN];
595 ir->nsteps = nsteps_cmdline;
596 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
598 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
599 gmx_step_str(nsteps_cmdline, sbuf_steps),
600 fabs(nsteps_cmdline*ir->delta_t));
602 else
604 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
605 gmx_step_str(nsteps_cmdline, sbuf_steps));
608 md_print_warn(cr, fplog, "%s\n", sbuf_msg);
610 else if (nsteps_cmdline < -2)
612 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %d",
613 nsteps_cmdline);
615 /* Do nothing if nsteps_cmdline == -2 */
618 namespace gmx
621 //! \brief Return the correct integrator function.
622 static integrator_t *my_integrator(unsigned int ei)
624 switch (ei)
626 case eiMD:
627 case eiBD:
628 case eiSD2:
629 case eiSD1:
630 case eiVV:
631 case eiVVAK:
632 if (!EI_DYNAMICS(ei))
634 GMX_THROW(APIError("do_md integrator would be called for a non-dynamical integrator"));
636 return do_md;
637 case eiSteep:
638 return do_steep;
639 case eiCG:
640 return do_cg;
641 case eiNM:
642 return do_nm;
643 case eiLBFGS:
644 return do_lbfgs;
645 case eiTPI:
646 case eiTPIC:
647 if (!EI_TPI(ei))
649 GMX_THROW(APIError("do_tpi integrator would be called for a non-TPI integrator"));
651 return do_tpi;
652 default:
653 GMX_THROW(APIError("Non existing integrator selected"));
657 int mdrunner(gmx_hw_opt_t *hw_opt,
658 FILE *fplog, t_commrec *cr, int nfile,
659 const t_filenm fnm[], const gmx_output_env_t *oenv, gmx_bool bVerbose,
660 int nstglobalcomm,
661 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
662 const char *dddlb_opt, real dlb_scale,
663 const char *ddcsx, const char *ddcsy, const char *ddcsz,
664 const char *nbpu_opt, int nstlist_cmdline,
665 gmx_int64_t nsteps_cmdline, int nstepout, int resetstep,
666 int gmx_unused nmultisim, int repl_ex_nst, int repl_ex_nex,
667 int repl_ex_seed, real pforce, real cpt_period, real max_hours,
668 int imdport, unsigned long Flags)
670 gmx_bool bForceUseGPU, bTryUseGPU, bRerunMD;
671 t_inputrec *inputrec;
672 t_state *state = NULL;
673 matrix box;
674 gmx_ddbox_t ddbox = {0};
675 int npme_major, npme_minor;
676 t_nrnb *nrnb;
677 gmx_mtop_t *mtop = NULL;
678 t_mdatoms *mdatoms = NULL;
679 t_forcerec *fr = NULL;
680 t_fcdata *fcd = NULL;
681 real ewaldcoeff_q = 0;
682 real ewaldcoeff_lj = 0;
683 struct gmx_pme_t **pmedata = NULL;
684 gmx_vsite_t *vsite = NULL;
685 gmx_constr_t constr;
686 int nChargePerturbed = -1, nTypePerturbed = 0, status;
687 gmx_wallcycle_t wcycle;
688 gmx_walltime_accounting_t walltime_accounting = NULL;
689 int rc;
690 gmx_int64_t reset_counters;
691 gmx_edsam_t ed = NULL;
692 int nthreads_pme = 1;
693 gmx_membed_t *membed = NULL;
694 gmx_hw_info_t *hwinfo = NULL;
695 /* The master rank decides early on bUseGPU and broadcasts this later */
696 gmx_bool bUseGPU = FALSE;
698 /* CAUTION: threads may be started later on in this function, so
699 cr doesn't reflect the final parallel state right now */
700 snew(inputrec, 1);
701 snew(mtop, 1);
703 if (Flags & MD_APPENDFILES)
705 fplog = NULL;
708 bRerunMD = (Flags & MD_RERUN);
709 bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
710 bTryUseGPU = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
712 /* Detect hardware, gather information. This is an operation that is
713 * global for this process (MPI rank). */
714 hwinfo = gmx_detect_hardware(fplog, cr, bTryUseGPU);
716 gmx_print_detected_hardware(fplog, cr, hwinfo);
718 if (fplog != NULL)
720 /* Print references after all software/hardware printing */
721 please_cite(fplog, "Abraham2015");
722 please_cite(fplog, "Pall2015");
723 please_cite(fplog, "Pronk2013");
724 please_cite(fplog, "Hess2008b");
725 please_cite(fplog, "Spoel2005a");
726 please_cite(fplog, "Lindahl2001a");
727 please_cite(fplog, "Berendsen95a");
730 snew(state, 1);
731 if (SIMMASTER(cr))
733 /* Read (nearly) all data required for the simulation */
734 read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, state, mtop);
736 if (inputrec->cutoff_scheme == ecutsVERLET)
738 /* Here the master rank decides if all ranks will use GPUs */
739 bUseGPU = (hwinfo->gpu_info.n_dev_compatible > 0 ||
740 getenv("GMX_EMULATE_GPU") != NULL);
742 /* TODO add GPU kernels for this and replace this check by:
743 * (bUseGPU && (ir->vdwtype == evdwPME &&
744 * ir->ljpme_combination_rule == eljpmeLB))
745 * update the message text and the content of nbnxn_acceleration_supported.
747 if (bUseGPU &&
748 !nbnxn_gpu_acceleration_supported(fplog, cr, inputrec, bRerunMD))
750 /* Fallback message printed by nbnxn_acceleration_supported */
751 if (bForceUseGPU)
753 gmx_fatal(FARGS, "GPU acceleration requested, but not supported with the given input settings");
755 bUseGPU = FALSE;
758 prepare_verlet_scheme(fplog, cr,
759 inputrec, nstlist_cmdline, mtop, state->box,
760 bUseGPU);
762 else
764 if (nstlist_cmdline > 0)
766 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
769 if (hwinfo->gpu_info.n_dev_compatible > 0)
771 md_print_warn(cr, fplog,
772 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
773 " To use a GPU, set the mdp option: cutoff-scheme = Verlet\n");
776 if (bForceUseGPU)
778 gmx_fatal(FARGS, "GPU requested, but can't be used without cutoff-scheme=Verlet");
781 #ifdef GMX_TARGET_BGQ
782 md_print_warn(cr, fplog,
783 "NOTE: There is no SIMD implementation of the group scheme kernels on\n"
784 " BlueGene/Q. You will observe better performance from using the\n"
785 " Verlet cut-off scheme.\n");
786 #endif
789 if (inputrec->eI == eiSD2)
791 md_print_warn(cr, fplog, "The stochastic dynamics integrator %s is deprecated, since\n"
792 "it is slower than integrator %s and is slightly less accurate\n"
793 "with constraints. Use the %s integrator.",
794 ei_names[inputrec->eI], ei_names[eiSD1], ei_names[eiSD1]);
798 /* Check and update the hardware options for internal consistency */
799 check_and_update_hw_opt_1(hw_opt, cr);
801 /* Early check for externally set process affinity. */
802 gmx_check_thread_affinity_set(fplog, cr,
803 hw_opt, hwinfo->nthreads_hw_avail, FALSE);
805 #ifdef GMX_THREAD_MPI
806 if (SIMMASTER(cr))
808 if (cr->npmenodes > 0 && hw_opt->nthreads_tmpi <= 0)
810 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
813 /* Since the master knows the cut-off scheme, update hw_opt for this.
814 * This is done later for normal MPI and also once more with tMPI
815 * for all tMPI ranks.
817 check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
819 /* NOW the threads will be started: */
820 hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
821 hw_opt,
822 inputrec, mtop,
823 cr, fplog, bUseGPU);
825 if (hw_opt->nthreads_tmpi > 1)
827 t_commrec *cr_old = cr;
828 /* now start the threads. */
829 cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
830 oenv, bVerbose, nstglobalcomm,
831 ddxyz, dd_node_order, rdd, rconstr,
832 dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
833 nbpu_opt, nstlist_cmdline,
834 nsteps_cmdline, nstepout, resetstep, nmultisim,
835 repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
836 cpt_period, max_hours,
837 Flags);
838 /* the main thread continues here with a new cr. We don't deallocate
839 the old cr because other threads may still be reading it. */
840 if (cr == NULL)
842 gmx_comm("Failed to spawn threads");
846 #endif
847 /* END OF CAUTION: cr is now reliable */
849 /* g_membed initialisation *
850 * Because we change the mtop, init_membed is called before the init_parallel *
851 * (in case we ever want to make it run in parallel) */
852 if (opt2bSet("-membed", nfile, fnm))
854 if (MASTER(cr))
856 fprintf(stderr, "Initializing membed");
858 membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
861 if (PAR(cr))
863 /* now broadcast everything to the non-master nodes/threads: */
864 init_parallel(cr, inputrec, mtop);
866 /* The master rank decided on the use of GPUs,
867 * broadcast this information to all ranks.
869 gmx_bcast_sim(sizeof(bUseGPU), &bUseGPU, cr);
872 if (fplog != NULL)
874 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
875 fprintf(fplog, "\n");
878 /* now make sure the state is initialized and propagated */
879 set_state_entries(state, inputrec);
881 /* A parallel command line option consistency check that we can
882 only do after any threads have started. */
883 if (!PAR(cr) &&
884 (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
886 gmx_fatal(FARGS,
887 "The -dd or -npme option request a parallel simulation, "
888 #ifndef GMX_MPI
889 "but %s was compiled without threads or MPI enabled"
890 #else
891 #ifdef GMX_THREAD_MPI
892 "but the number of threads (option -nt) is 1"
893 #else
894 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
895 #endif
896 #endif
897 , output_env_get_program_display_name(oenv)
901 if (bRerunMD &&
902 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
904 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
907 if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
909 gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
912 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
914 if (cr->npmenodes > 0)
916 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
917 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
920 cr->npmenodes = 0;
923 if (bUseGPU && cr->npmenodes < 0)
925 /* With GPUs we don't automatically use PME-only ranks. PME ranks can
926 * improve performance with many threads per GPU, since our OpenMP
927 * scaling is bad, but it's difficult to automate the setup.
929 cr->npmenodes = 0;
932 #ifdef GMX_FAHCORE
933 if (MASTER(cr))
935 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
937 #endif
939 /* NMR restraints must be initialized before load_checkpoint,
940 * since with time averaging the history is added to t_state.
941 * For proper consistency check we therefore need to extend
942 * t_state here.
943 * So the PME-only nodes (if present) will also initialize
944 * the distance restraints.
946 snew(fcd, 1);
948 /* This needs to be called before read_checkpoint to extend the state */
949 init_disres(fplog, mtop, inputrec, cr, fcd, state, repl_ex_nst > 0);
951 init_orires(fplog, mtop, state->x, inputrec, cr, &(fcd->orires),
952 state);
954 if (inputrecDeform(inputrec))
956 /* Store the deform reference box before reading the checkpoint */
957 if (SIMMASTER(cr))
959 copy_mat(state->box, box);
961 if (PAR(cr))
963 gmx_bcast(sizeof(box), box, cr);
965 /* Because we do not have the update struct available yet
966 * in which the reference values should be stored,
967 * we store them temporarily in static variables.
968 * This should be thread safe, since they are only written once
969 * and with identical values.
971 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
972 deform_init_init_step_tpx = inputrec->init_step;
973 copy_mat(box, deform_init_box_tpx);
974 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
977 if (Flags & MD_STARTFROMCPT)
979 /* Check if checkpoint file exists before doing continuation.
980 * This way we can use identical input options for the first and subsequent runs...
982 gmx_bool bReadEkin;
984 load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
985 cr, ddxyz,
986 inputrec, state, &bReadEkin,
987 (Flags & MD_APPENDFILES),
988 (Flags & MD_APPENDFILESSET));
990 if (bReadEkin)
992 Flags |= MD_READ_EKIN;
996 if (MASTER(cr) && (Flags & MD_APPENDFILES))
998 gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
999 Flags, &fplog);
1002 /* override nsteps with value from cmdline */
1003 override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
1005 if (SIMMASTER(cr))
1007 copy_mat(state->box, box);
1010 if (PAR(cr))
1012 gmx_bcast(sizeof(box), box, cr);
1015 /* Essential dynamics */
1016 if (opt2bSet("-ei", nfile, fnm))
1018 /* Open input and output files, allocate space for ED data structure */
1019 ed = ed_open(mtop->natoms, &state->edsamstate, nfile, fnm, Flags, oenv, cr);
1022 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
1023 inputrec->eI == eiNM))
1025 cr->dd = init_domain_decomposition(fplog, cr, Flags, ddxyz, rdd, rconstr,
1026 dddlb_opt, dlb_scale,
1027 ddcsx, ddcsy, ddcsz,
1028 mtop, inputrec,
1029 box, state->x,
1030 &ddbox, &npme_major, &npme_minor);
1032 make_dd_communicators(fplog, cr, dd_node_order);
1034 /* Set overallocation to avoid frequent reallocation of arrays */
1035 set_over_alloc_dd(TRUE);
1037 else
1039 /* PME, if used, is done on all nodes with 1D decomposition */
1040 cr->npmenodes = 0;
1041 cr->duty = (DUTY_PP | DUTY_PME);
1042 npme_major = 1;
1043 npme_minor = 1;
1045 if (inputrec->ePBC == epbcSCREW)
1047 gmx_fatal(FARGS,
1048 "pbc=%s is only implemented with domain decomposition",
1049 epbc_names[inputrec->ePBC]);
1053 if (PAR(cr))
1055 /* After possible communicator splitting in make_dd_communicators.
1056 * we can set up the intra/inter node communication.
1058 gmx_setup_nodecomm(fplog, cr);
1061 /* Initialize per-physical-node MPI process/thread ID and counters. */
1062 gmx_init_intranode_counters(cr);
1063 #ifdef GMX_MPI
1064 if (MULTISIM(cr))
1066 md_print_info(cr, fplog,
1067 "This is simulation %d out of %d running as a composite GROMACS\n"
1068 "multi-simulation job. Setup for this simulation:\n\n",
1069 cr->ms->sim, cr->ms->nsim);
1071 md_print_info(cr, fplog, "Using %d MPI %s\n",
1072 cr->nnodes,
1073 #ifdef GMX_THREAD_MPI
1074 cr->nnodes == 1 ? "thread" : "threads"
1075 #else
1076 cr->nnodes == 1 ? "process" : "processes"
1077 #endif
1079 fflush(stderr);
1080 #endif
1082 /* Check and update hw_opt for the cut-off scheme */
1083 check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1085 /* Check and update hw_opt for the number of MPI ranks */
1086 check_and_update_hw_opt_3(hw_opt);
1088 gmx_omp_nthreads_init(fplog, cr,
1089 hwinfo->nthreads_hw_avail,
1090 hw_opt->nthreads_omp,
1091 hw_opt->nthreads_omp_pme,
1092 (cr->duty & DUTY_PP) == 0,
1093 inputrec->cutoff_scheme == ecutsVERLET);
1095 #ifndef NDEBUG
1096 if (EI_TPI(inputrec->eI) &&
1097 inputrec->cutoff_scheme == ecutsVERLET)
1099 gmx_feenableexcept();
1101 #endif
1103 if (bUseGPU)
1105 /* Select GPU id's to use */
1106 gmx_select_gpu_ids(fplog, cr, &hwinfo->gpu_info, bForceUseGPU,
1107 &hw_opt->gpu_opt);
1109 else
1111 /* Ignore (potentially) manually selected GPUs */
1112 hw_opt->gpu_opt.n_dev_use = 0;
1115 /* check consistency across ranks of things like SIMD
1116 * support and number of GPUs selected */
1117 gmx_check_hw_runconf_consistency(fplog, hwinfo, cr, hw_opt, bUseGPU);
1119 /* Now that we know the setup is consistent, check for efficiency */
1120 check_resource_division_efficiency(hwinfo, hw_opt, Flags & MD_NTOMPSET,
1121 cr, fplog);
1123 if (DOMAINDECOMP(cr))
1125 /* When we share GPUs over ranks, we need to know this for the DLB */
1126 dd_setup_dlb_resource_sharing(cr, hwinfo, hw_opt);
1129 /* getting number of PP/PME threads
1130 PME: env variable should be read only on one node to make sure it is
1131 identical everywhere;
1133 nthreads_pme = gmx_omp_nthreads_get(emntPME);
1135 wcycle = wallcycle_init(fplog, resetstep, cr);
1137 if (PAR(cr))
1139 /* Master synchronizes its value of reset_counters with all nodes
1140 * including PME only nodes */
1141 reset_counters = wcycle_get_reset_counters(wcycle);
1142 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1143 wcycle_set_reset_counters(wcycle, reset_counters);
1146 snew(nrnb, 1);
1147 if (cr->duty & DUTY_PP)
1149 bcast_state(cr, state);
1151 /* Initiate forcerecord */
1152 fr = mk_forcerec();
1153 fr->hwinfo = hwinfo;
1154 fr->gpu_opt = &hw_opt->gpu_opt;
1155 init_forcerec(fplog, fr, fcd, inputrec, mtop, cr, box,
1156 opt2fn("-table", nfile, fnm),
1157 opt2fn("-tablep", nfile, fnm),
1158 opt2fn("-tableb", nfile, fnm),
1159 nbpu_opt,
1160 FALSE,
1161 pforce);
1163 /* Initialize QM-MM */
1164 if (fr->bQMMM)
1166 init_QMMMrec(cr, mtop, inputrec, fr);
1169 /* Initialize the mdatoms structure.
1170 * mdatoms is not filled with atom data,
1171 * as this can not be done now with domain decomposition.
1173 mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1175 /* Initialize the virtual site communication */
1176 vsite = init_vsite(mtop, cr, FALSE);
1178 calc_shifts(box, fr->shift_vec);
1180 /* With periodic molecules the charge groups should be whole at start up
1181 * and the virtual sites should not be far from their proper positions.
1183 if (!inputrec->bContinuation && MASTER(cr) &&
1184 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1186 /* Make molecules whole at start of run */
1187 if (fr->ePBC != epbcNONE)
1189 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, state->x);
1191 if (vsite)
1193 /* Correct initial vsite positions are required
1194 * for the initial distribution in the domain decomposition
1195 * and for the initial shell prediction.
1197 construct_vsites_mtop(vsite, mtop, state->x);
1201 if (EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype))
1203 ewaldcoeff_q = fr->ewaldcoeff_q;
1204 ewaldcoeff_lj = fr->ewaldcoeff_lj;
1205 pmedata = &fr->pmedata;
1207 else
1209 pmedata = NULL;
1212 else
1214 /* This is a PME only node */
1216 /* We don't need the state */
1217 done_state(state);
1219 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1220 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1221 snew(pmedata, 1);
1224 if (hw_opt->thread_affinity != threadaffOFF)
1226 /* Before setting affinity, check whether the affinity has changed
1227 * - which indicates that probably the OpenMP library has changed it
1228 * since we first checked).
1230 gmx_check_thread_affinity_set(fplog, cr,
1231 hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1233 /* Set the CPU affinity */
1234 gmx_set_thread_affinity(fplog, cr, hw_opt, hwinfo);
1237 /* Initiate PME if necessary,
1238 * either on all nodes or on dedicated PME nodes only. */
1239 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1241 if (mdatoms)
1243 nChargePerturbed = mdatoms->nChargePerturbed;
1244 if (EVDW_PME(inputrec->vdwtype))
1246 nTypePerturbed = mdatoms->nTypePerturbed;
1249 if (cr->npmenodes > 0)
1251 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1252 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1253 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1256 if (cr->duty & DUTY_PME)
1258 status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1259 mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1260 (Flags & MD_REPRODUCIBLE), nthreads_pme);
1261 if (status != 0)
1263 gmx_fatal(FARGS, "Error %d initializing PME", status);
1269 if (EI_DYNAMICS(inputrec->eI))
1271 /* Turn on signal handling on all nodes */
1273 * (A user signal from the PME nodes (if any)
1274 * is communicated to the PP nodes.
1276 signal_handler_install();
1279 if (cr->duty & DUTY_PP)
1281 /* Assumes uniform use of the number of OpenMP threads */
1282 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1284 if (inputrec->bPull)
1286 /* Initialize pull code */
1287 inputrec->pull_work =
1288 init_pull(fplog, inputrec->pull, inputrec, nfile, fnm,
1289 mtop, cr, oenv, inputrec->fepvals->init_lambda,
1290 EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1293 if (inputrec->bRot)
1295 /* Initialize enforced rotation code */
1296 init_rot(fplog, inputrec, nfile, fnm, cr, state->x, box, mtop, oenv,
1297 bVerbose, Flags);
1300 if (inputrec->eSwapCoords != eswapNO)
1302 /* Initialize ion swapping code */
1303 init_swapcoords(fplog, bVerbose, inputrec, opt2fn_master("-swap", nfile, fnm, cr),
1304 mtop, state->x, state->box, &state->swapstate, cr, oenv, Flags);
1307 constr = init_constraints(fplog, mtop, inputrec, ed, state, cr);
1309 if (DOMAINDECOMP(cr))
1311 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1312 dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1313 Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1315 set_dd_parameters(fplog, cr->dd, dlb_scale, inputrec, &ddbox);
1317 setup_dd_grid(fplog, cr->dd);
1320 /* Now do whatever the user wants us to do (how flexible...) */
1321 my_integrator(inputrec->eI) (fplog, cr, nfile, fnm,
1322 oenv, bVerbose,
1323 nstglobalcomm,
1324 vsite, constr,
1325 nstepout, inputrec, mtop,
1326 fcd, state,
1327 mdatoms, nrnb, wcycle, ed, fr,
1328 repl_ex_nst, repl_ex_nex, repl_ex_seed,
1329 membed,
1330 cpt_period, max_hours,
1331 imdport,
1332 Flags,
1333 walltime_accounting);
1335 if (inputrec->bPull)
1337 finish_pull(inputrec->pull_work);
1340 if (inputrec->bRot)
1342 finish_rot(inputrec->rot);
1346 else
1348 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1349 /* do PME only */
1350 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1351 gmx_pmeonly(*pmedata, cr, nrnb, wcycle, walltime_accounting, ewaldcoeff_q, ewaldcoeff_lj, inputrec);
1354 wallcycle_stop(wcycle, ewcRUN);
1356 /* Finish up, write some stuff
1357 * if rerunMD, don't write last frame again
1359 finish_run(fplog, cr,
1360 inputrec, nrnb, wcycle, walltime_accounting,
1361 fr ? fr->nbv : NULL,
1362 EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1365 /* Free GPU memory and context */
1366 free_gpu_resources(fr, cr, &hwinfo->gpu_info, fr ? fr->gpu_opt : NULL);
1368 if (membed != nullptr)
1370 free_membed(membed);
1373 gmx_hardware_info_free(hwinfo);
1375 /* Does what it says */
1376 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1377 walltime_accounting_destroy(walltime_accounting);
1379 /* Close logfile already here if we were appending to it */
1380 if (MASTER(cr) && (Flags & MD_APPENDFILES))
1382 gmx_log_close(fplog);
1385 rc = (int)gmx_get_stop_condition();
1387 done_ed(&ed);
1389 #ifdef GMX_THREAD_MPI
1390 /* we need to join all threads. The sub-threads join when they
1391 exit this function, but the master thread needs to be told to
1392 wait for that. */
1393 if (PAR(cr) && MASTER(cr))
1395 tMPI_Finalize();
1397 #endif
1399 return rc;
1402 } // namespace gmx