Fix membed with partial revert of 29943f
[gromacs/AngularHB.git] / src / programs / mdrun / runner.cpp
blob6fb9d5f5707f52af4c545f5a144880af27d4053b
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 "membed.h"
114 #include "repl_ex.h"
115 #include "resource-division.h"
117 #ifdef GMX_FAHCORE
118 #include "corewrap.h"
119 #endif
121 //! First step used in pressure scaling
122 gmx_int64_t deform_init_init_step_tpx;
123 //! Initial box for pressure scaling
124 matrix deform_init_box_tpx;
125 //! MPI variable for use in pressure scaling
126 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
128 #if GMX_THREAD_MPI
129 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
130 * the number of threads will get lowered.
132 #define MIN_ATOMS_PER_MPI_THREAD 90
133 #define MIN_ATOMS_PER_GPU 900
135 struct mdrunner_arglist
137 gmx_hw_opt_t hw_opt;
138 FILE *fplog;
139 t_commrec *cr;
140 int nfile;
141 const t_filenm *fnm;
142 const gmx_output_env_t *oenv;
143 gmx_bool bVerbose;
144 int nstglobalcomm;
145 ivec ddxyz;
146 int dd_rank_order;
147 int npme;
148 real rdd;
149 real rconstr;
150 const char *dddlb_opt;
151 real dlb_scale;
152 const char *ddcsx;
153 const char *ddcsy;
154 const char *ddcsz;
155 const char *nbpu_opt;
156 int nstlist_cmdline;
157 gmx_int64_t nsteps_cmdline;
158 int nstepout;
159 int resetstep;
160 int nmultisim;
161 int repl_ex_nst;
162 int repl_ex_nex;
163 int repl_ex_seed;
164 real pforce;
165 real cpt_period;
166 real max_hours;
167 int imdport;
168 unsigned long Flags;
172 /* The function used for spawning threads. Extracts the mdrunner()
173 arguments from its one argument and calls mdrunner(), after making
174 a commrec. */
175 static void mdrunner_start_fn(void *arg)
179 struct mdrunner_arglist *mda = (struct mdrunner_arglist*)arg;
180 struct mdrunner_arglist mc = *mda; /* copy the arg list to make sure
181 that it's thread-local. This doesn't
182 copy pointed-to items, of course,
183 but those are all const. */
184 t_commrec *cr; /* we need a local version of this */
185 FILE *fplog = NULL;
186 t_filenm *fnm;
188 fnm = dup_tfn(mc.nfile, mc.fnm);
190 cr = reinitialize_commrec_for_this_thread(mc.cr);
192 if (MASTER(cr))
194 fplog = mc.fplog;
197 gmx::mdrunner(&mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
198 mc.bVerbose, mc.nstglobalcomm,
199 mc.ddxyz, mc.dd_rank_order, mc.npme, mc.rdd,
200 mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
201 mc.ddcsx, mc.ddcsy, mc.ddcsz,
202 mc.nbpu_opt, mc.nstlist_cmdline,
203 mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
204 mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce,
205 mc.cpt_period, mc.max_hours, mc.imdport, mc.Flags);
207 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
211 /* called by mdrunner() to start a specific number of threads (including
212 the main thread) for thread-parallel runs. This in turn calls mdrunner()
213 for each thread.
214 All options besides nthreads are the same as for mdrunner(). */
215 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
216 FILE *fplog, t_commrec *cr, int nfile,
217 const t_filenm fnm[], const gmx_output_env_t *oenv, gmx_bool bVerbose,
218 int nstglobalcomm,
219 ivec ddxyz, int dd_rank_order, int npme,
220 real rdd, real rconstr,
221 const char *dddlb_opt, real dlb_scale,
222 const char *ddcsx, const char *ddcsy, const char *ddcsz,
223 const char *nbpu_opt, int nstlist_cmdline,
224 gmx_int64_t nsteps_cmdline,
225 int nstepout, int resetstep,
226 int nmultisim, int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
227 real pforce, real cpt_period, real max_hours,
228 unsigned long Flags)
230 int ret;
231 struct mdrunner_arglist *mda;
232 t_commrec *crn; /* the new commrec */
233 t_filenm *fnmn;
235 /* first check whether we even need to start tMPI */
236 if (hw_opt->nthreads_tmpi < 2)
238 return cr;
241 /* a few small, one-time, almost unavoidable memory leaks: */
242 snew(mda, 1);
243 fnmn = dup_tfn(nfile, fnm);
245 /* fill the data structure to pass as void pointer to thread start fn */
246 /* hw_opt contains pointers, which should all be NULL at this stage */
247 mda->hw_opt = *hw_opt;
248 mda->fplog = fplog;
249 mda->cr = cr;
250 mda->nfile = nfile;
251 mda->fnm = fnmn;
252 mda->oenv = oenv;
253 mda->bVerbose = bVerbose;
254 mda->nstglobalcomm = nstglobalcomm;
255 mda->ddxyz[XX] = ddxyz[XX];
256 mda->ddxyz[YY] = ddxyz[YY];
257 mda->ddxyz[ZZ] = ddxyz[ZZ];
258 mda->dd_rank_order = dd_rank_order;
259 mda->npme = npme;
260 mda->rdd = rdd;
261 mda->rconstr = rconstr;
262 mda->dddlb_opt = dddlb_opt;
263 mda->dlb_scale = dlb_scale;
264 mda->ddcsx = ddcsx;
265 mda->ddcsy = ddcsy;
266 mda->ddcsz = ddcsz;
267 mda->nbpu_opt = nbpu_opt;
268 mda->nstlist_cmdline = nstlist_cmdline;
269 mda->nsteps_cmdline = nsteps_cmdline;
270 mda->nstepout = nstepout;
271 mda->resetstep = resetstep;
272 mda->nmultisim = nmultisim;
273 mda->repl_ex_nst = repl_ex_nst;
274 mda->repl_ex_nex = repl_ex_nex;
275 mda->repl_ex_seed = repl_ex_seed;
276 mda->pforce = pforce;
277 mda->cpt_period = cpt_period;
278 mda->max_hours = max_hours;
279 mda->Flags = Flags;
281 /* now spawn new threads that start mdrunner_start_fn(), while
282 the main thread returns, we set thread affinity later */
283 ret = tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi, TMPI_AFFINITY_NONE,
284 mdrunner_start_fn, (void*)(mda) );
285 if (ret != TMPI_SUCCESS)
287 return NULL;
290 crn = reinitialize_commrec_for_this_thread(cr);
291 return crn;
294 #endif /* GMX_THREAD_MPI */
297 /*! \brief Cost of non-bonded kernels
299 * We determine the extra cost of the non-bonded kernels compared to
300 * a reference nstlist value of 10 (which is the default in grompp).
302 static const int nbnxnReferenceNstlist = 10;
303 //! The values to try when switching
304 const int nstlist_try[] = { 20, 25, 40 };
305 //! Number of elements in the neighborsearch list trials.
306 #define NNSTL sizeof(nstlist_try)/sizeof(nstlist_try[0])
307 /* Increase nstlist until the non-bonded cost increases more than listfac_ok,
308 * but never more than listfac_max.
309 * A standard (protein+)water system at 300K with PME ewald_rtol=1e-5
310 * needs 1.28 at rcoulomb=0.9 and 1.24 at rcoulomb=1.0 to get to nstlist=40.
311 * Note that both CPU and GPU factors are conservative. Performance should
312 * not go down due to this tuning, except with a relatively slow GPU.
313 * On the other hand, at medium/high parallelization or with fast GPUs
314 * nstlist will not be increased enough to reach optimal performance.
316 /* CPU: pair-search is about a factor 1.5 slower than the non-bonded kernel */
317 //! Max OK performance ratio beween force calc and neighbor searching
318 static const float nbnxn_cpu_listfac_ok = 1.05;
319 //! Too high performance ratio beween force calc and neighbor searching
320 static const float nbnxn_cpu_listfac_max = 1.09;
321 /* CPU: pair-search is about a factor 2-3 slower than the non-bonded kernel */
322 //! Max OK performance ratio beween force calc and neighbor searching
323 static const float nbnxn_knl_listfac_ok = 1.22;
324 //! Too high performance ratio beween force calc and neighbor searching
325 static const float nbnxn_knl_listfac_max = 1.3;
326 /* GPU: pair-search is a factor 1.5-3 slower than the non-bonded kernel */
327 //! Max OK performance ratio beween force calc and neighbor searching
328 static const float nbnxn_gpu_listfac_ok = 1.20;
329 //! Too high performance ratio beween force calc and neighbor searching
330 static const float nbnxn_gpu_listfac_max = 1.30;
332 /*! \brief Try to increase nstlist when using the Verlet cut-off scheme */
333 static void increase_nstlist(FILE *fp, t_commrec *cr,
334 t_inputrec *ir, int nstlist_cmdline,
335 const gmx_mtop_t *mtop, matrix box,
336 gmx_bool bGPU, const gmx::CpuInfo &cpuinfo)
338 float listfac_ok, listfac_max;
339 int nstlist_orig, nstlist_prev;
340 verletbuf_list_setup_t ls;
341 real rlistWithReferenceNstlist, rlist_inc, rlist_ok, rlist_max;
342 real rlist_new, rlist_prev;
343 size_t nstlist_ind = 0;
344 t_state state_tmp;
345 gmx_bool bBox, bDD, bCont;
346 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";
347 const char *nve_err = "Can not increase nstlist because an NVE ensemble is used";
348 const char *vbd_err = "Can not increase nstlist because verlet-buffer-tolerance is not set or used";
349 const char *box_err = "Can not increase nstlist because the box is too small";
350 const char *dd_err = "Can not increase nstlist because of domain decomposition limitations";
351 char buf[STRLEN];
353 if (nstlist_cmdline <= 0)
355 if (ir->nstlist == 1)
357 /* The user probably set nstlist=1 for a reason,
358 * don't mess with the settings.
360 return;
363 if (fp != NULL && bGPU && ir->nstlist < nstlist_try[0])
365 fprintf(fp, nstl_gpu, ir->nstlist);
367 nstlist_ind = 0;
368 while (nstlist_ind < NNSTL && ir->nstlist >= nstlist_try[nstlist_ind])
370 nstlist_ind++;
372 if (nstlist_ind == NNSTL)
374 /* There are no larger nstlist value to try */
375 return;
379 if (EI_MD(ir->eI) && ir->etc == etcNO)
381 if (MASTER(cr))
383 fprintf(stderr, "%s\n", nve_err);
385 if (fp != NULL)
387 fprintf(fp, "%s\n", nve_err);
390 return;
393 if (ir->verletbuf_tol == 0 && bGPU)
395 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");
398 if (ir->verletbuf_tol < 0)
400 if (MASTER(cr))
402 fprintf(stderr, "%s\n", vbd_err);
404 if (fp != NULL)
406 fprintf(fp, "%s\n", vbd_err);
409 return;
412 if (bGPU)
414 listfac_ok = nbnxn_gpu_listfac_ok;
415 listfac_max = nbnxn_gpu_listfac_max;
417 else if (cpuinfo.feature(gmx::CpuInfo::Feature::X86_Avx512ER))
419 listfac_ok = nbnxn_knl_listfac_ok;
420 listfac_max = nbnxn_knl_listfac_max;
422 else
424 listfac_ok = nbnxn_cpu_listfac_ok;
425 listfac_max = nbnxn_cpu_listfac_max;
428 nstlist_orig = ir->nstlist;
429 if (nstlist_cmdline > 0)
431 if (fp)
433 sprintf(buf, "Getting nstlist=%d from command line option",
434 nstlist_cmdline);
436 ir->nstlist = nstlist_cmdline;
439 verletbuf_get_list_setup(TRUE, bGPU, &ls);
441 /* Allow rlist to make the list a given factor larger than the list
442 * would be with the reference value for nstlist (10).
444 nstlist_prev = ir->nstlist;
445 ir->nstlist = nbnxnReferenceNstlist;
446 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL,
447 &rlistWithReferenceNstlist);
448 ir->nstlist = nstlist_prev;
450 /* Determine the pair list size increase due to zero interactions */
451 rlist_inc = nbnxn_get_rlist_effective_inc(ls.cluster_size_j,
452 mtop->natoms/det(box));
453 rlist_ok = (rlistWithReferenceNstlist + rlist_inc)*std::cbrt(listfac_ok) - rlist_inc;
454 rlist_max = (rlistWithReferenceNstlist + rlist_inc)*std::cbrt(listfac_max) - rlist_inc;
455 if (debug)
457 fprintf(debug, "nstlist tuning: rlist_inc %.3f rlist_ok %.3f rlist_max %.3f\n",
458 rlist_inc, rlist_ok, rlist_max);
461 nstlist_prev = nstlist_orig;
462 rlist_prev = ir->rlist;
465 if (nstlist_cmdline <= 0)
467 ir->nstlist = nstlist_try[nstlist_ind];
470 /* Set the pair-list buffer size in ir */
471 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
473 /* Does rlist fit in the box? */
474 bBox = (gmx::square(rlist_new) < max_cutoff2(ir->ePBC, box));
475 bDD = TRUE;
476 if (bBox && DOMAINDECOMP(cr))
478 /* Check if rlist fits in the domain decomposition */
479 if (inputrec2nboundeddim(ir) < DIM)
481 gmx_incons("Changing nstlist with domain decomposition and unbounded dimensions is not implemented yet");
483 copy_mat(box, state_tmp.box);
484 bDD = change_dd_cutoff(cr, &state_tmp, ir, rlist_new);
487 if (debug)
489 fprintf(debug, "nstlist %d rlist %.3f bBox %d bDD %d\n",
490 ir->nstlist, rlist_new, bBox, bDD);
493 bCont = FALSE;
495 if (nstlist_cmdline <= 0)
497 if (bBox && bDD && rlist_new <= rlist_max)
499 /* Increase nstlist */
500 nstlist_prev = ir->nstlist;
501 rlist_prev = rlist_new;
502 bCont = (nstlist_ind+1 < NNSTL && rlist_new < rlist_ok);
504 else
506 /* Stick with the previous nstlist */
507 ir->nstlist = nstlist_prev;
508 rlist_new = rlist_prev;
509 bBox = TRUE;
510 bDD = TRUE;
514 nstlist_ind++;
516 while (bCont);
518 if (!bBox || !bDD)
520 gmx_warning(!bBox ? box_err : dd_err);
521 if (fp != NULL)
523 fprintf(fp, "\n%s\n", bBox ? box_err : dd_err);
525 ir->nstlist = nstlist_orig;
527 else if (ir->nstlist != nstlist_orig || rlist_new != ir->rlist)
529 sprintf(buf, "Changing nstlist from %d to %d, rlist from %g to %g",
530 nstlist_orig, ir->nstlist,
531 ir->rlist, rlist_new);
532 if (MASTER(cr))
534 fprintf(stderr, "%s\n\n", buf);
536 if (fp != NULL)
538 fprintf(fp, "%s\n\n", buf);
540 ir->rlist = rlist_new;
544 /*! \brief Initialize variables for Verlet scheme simulation */
545 static void prepare_verlet_scheme(FILE *fplog,
546 t_commrec *cr,
547 t_inputrec *ir,
548 int nstlist_cmdline,
549 const gmx_mtop_t *mtop,
550 matrix box,
551 gmx_bool bUseGPU,
552 const gmx::CpuInfo &cpuinfo)
554 /* For NVE simulations, we will retain the initial list buffer */
555 if (EI_DYNAMICS(ir->eI) &&
556 ir->verletbuf_tol > 0 &&
557 !(EI_MD(ir->eI) && ir->etc == etcNO))
559 /* Update the Verlet buffer size for the current run setup */
560 verletbuf_list_setup_t ls;
561 real rlist_new;
563 /* Here we assume SIMD-enabled kernels are being used. But as currently
564 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
565 * and 4x2 gives a larger buffer than 4x4, this is ok.
567 verletbuf_get_list_setup(TRUE, bUseGPU, &ls);
569 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
571 if (rlist_new != ir->rlist)
573 if (fplog != NULL)
575 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
576 ir->rlist, rlist_new,
577 ls.cluster_size_i, ls.cluster_size_j);
579 ir->rlist = rlist_new;
583 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
585 gmx_fatal(FARGS, "Can not set nstlist without %s",
586 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
589 if (EI_DYNAMICS(ir->eI))
591 /* Set or try nstlist values */
592 increase_nstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, bUseGPU, cpuinfo);
596 /*! \brief Override the nslist value in inputrec
598 * with value passed on the command line (if any)
600 static void override_nsteps_cmdline(FILE *fplog,
601 gmx_int64_t nsteps_cmdline,
602 t_inputrec *ir,
603 const t_commrec *cr)
605 assert(ir);
606 assert(cr);
608 /* override with anything else than the default -2 */
609 if (nsteps_cmdline > -2)
611 char sbuf_steps[STEPSTRSIZE];
612 char sbuf_msg[STRLEN];
614 ir->nsteps = nsteps_cmdline;
615 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
617 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
618 gmx_step_str(nsteps_cmdline, sbuf_steps),
619 fabs(nsteps_cmdline*ir->delta_t));
621 else
623 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
624 gmx_step_str(nsteps_cmdline, sbuf_steps));
627 md_print_warn(cr, fplog, "%s\n", sbuf_msg);
629 else if (nsteps_cmdline < -2)
631 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %d",
632 nsteps_cmdline);
634 /* Do nothing if nsteps_cmdline == -2 */
637 namespace gmx
640 //! \brief Return the correct integrator function.
641 static integrator_t *my_integrator(unsigned int ei)
643 switch (ei)
645 case eiMD:
646 case eiBD:
647 case eiSD1:
648 case eiVV:
649 case eiVVAK:
650 if (!EI_DYNAMICS(ei))
652 GMX_THROW(APIError("do_md integrator would be called for a non-dynamical integrator"));
654 return do_md;
655 case eiSteep:
656 return do_steep;
657 case eiCG:
658 return do_cg;
659 case eiNM:
660 return do_nm;
661 case eiLBFGS:
662 return do_lbfgs;
663 case eiTPI:
664 case eiTPIC:
665 if (!EI_TPI(ei))
667 GMX_THROW(APIError("do_tpi integrator would be called for a non-TPI integrator"));
669 return do_tpi;
670 case eiSD2_REMOVED:
671 GMX_THROW(NotImplementedError("SD2 integrator has been removed"));
672 default:
673 GMX_THROW(APIError("Non existing integrator selected"));
677 int mdrunner(gmx_hw_opt_t *hw_opt,
678 FILE *fplog, t_commrec *cr, int nfile,
679 const t_filenm fnm[], const gmx_output_env_t *oenv, gmx_bool bVerbose,
680 int nstglobalcomm,
681 ivec ddxyz, int dd_rank_order, int npme, real rdd, real rconstr,
682 const char *dddlb_opt, real dlb_scale,
683 const char *ddcsx, const char *ddcsy, const char *ddcsz,
684 const char *nbpu_opt, int nstlist_cmdline,
685 gmx_int64_t nsteps_cmdline, int nstepout, int resetstep,
686 int gmx_unused nmultisim, int repl_ex_nst, int repl_ex_nex,
687 int repl_ex_seed, real pforce, real cpt_period, real max_hours,
688 int imdport, unsigned long Flags)
690 gmx_bool bForceUseGPU, bTryUseGPU, bRerunMD;
691 t_inputrec *inputrec;
692 t_state *state = NULL;
693 matrix box;
694 gmx_ddbox_t ddbox = {0};
695 int npme_major, npme_minor;
696 t_nrnb *nrnb;
697 gmx_mtop_t *mtop = NULL;
698 t_mdatoms *mdatoms = NULL;
699 t_forcerec *fr = NULL;
700 t_fcdata *fcd = NULL;
701 real ewaldcoeff_q = 0;
702 real ewaldcoeff_lj = 0;
703 struct gmx_pme_t **pmedata = NULL;
704 gmx_vsite_t *vsite = NULL;
705 gmx_constr_t constr;
706 int nChargePerturbed = -1, nTypePerturbed = 0, status;
707 gmx_wallcycle_t wcycle;
708 gmx_walltime_accounting_t walltime_accounting = NULL;
709 int rc;
710 gmx_int64_t reset_counters;
711 gmx_edsam_t ed = NULL;
712 int nthreads_pme = 1;
713 gmx_membed_t * membed = NULL;
714 gmx_hw_info_t *hwinfo = NULL;
715 /* The master rank decides early on bUseGPU and broadcasts this later */
716 gmx_bool bUseGPU = FALSE;
718 /* CAUTION: threads may be started later on in this function, so
719 cr doesn't reflect the final parallel state right now */
720 snew(inputrec, 1);
721 snew(mtop, 1);
723 if (Flags & MD_APPENDFILES)
725 fplog = NULL;
728 bRerunMD = (Flags & MD_RERUN);
729 bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
730 bTryUseGPU = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
732 /* Detect hardware, gather information. This is an operation that is
733 * global for this process (MPI rank). */
734 hwinfo = gmx_detect_hardware(fplog, cr, bTryUseGPU);
736 gmx_print_detected_hardware(fplog, cr, hwinfo);
738 if (fplog != NULL)
740 /* Print references after all software/hardware printing */
741 please_cite(fplog, "Abraham2015");
742 please_cite(fplog, "Pall2015");
743 please_cite(fplog, "Pronk2013");
744 please_cite(fplog, "Hess2008b");
745 please_cite(fplog, "Spoel2005a");
746 please_cite(fplog, "Lindahl2001a");
747 please_cite(fplog, "Berendsen95a");
750 snew(state, 1);
751 if (SIMMASTER(cr))
753 /* Read (nearly) all data required for the simulation */
754 read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, state, mtop);
756 if (inputrec->cutoff_scheme == ecutsVERLET)
758 /* Here the master rank decides if all ranks will use GPUs */
759 bUseGPU = (hwinfo->gpu_info.n_dev_compatible > 0 ||
760 getenv("GMX_EMULATE_GPU") != NULL);
762 /* TODO add GPU kernels for this and replace this check by:
763 * (bUseGPU && (ir->vdwtype == evdwPME &&
764 * ir->ljpme_combination_rule == eljpmeLB))
765 * update the message text and the content of nbnxn_acceleration_supported.
767 if (bUseGPU &&
768 !nbnxn_gpu_acceleration_supported(fplog, cr, inputrec, bRerunMD))
770 /* Fallback message printed by nbnxn_acceleration_supported */
771 if (bForceUseGPU)
773 gmx_fatal(FARGS, "GPU acceleration requested, but not supported with the given input settings");
775 bUseGPU = FALSE;
778 prepare_verlet_scheme(fplog, cr,
779 inputrec, nstlist_cmdline, mtop, state->box,
780 bUseGPU, *hwinfo->cpuInfo);
782 else
784 if (nstlist_cmdline > 0)
786 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
789 if (hwinfo->gpu_info.n_dev_compatible > 0)
791 md_print_warn(cr, fplog,
792 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
793 " To use a GPU, set the mdp option: cutoff-scheme = Verlet\n");
796 if (bForceUseGPU)
798 gmx_fatal(FARGS, "GPU requested, but can't be used without cutoff-scheme=Verlet");
801 #if GMX_TARGET_BGQ
802 md_print_warn(cr, fplog,
803 "NOTE: There is no SIMD implementation of the group scheme kernels on\n"
804 " BlueGene/Q. You will observe better performance from using the\n"
805 " Verlet cut-off scheme.\n");
806 #endif
810 /* Check and update the hardware options for internal consistency */
811 check_and_update_hw_opt_1(hw_opt, cr, npme);
813 /* Early check for externally set process affinity. */
814 gmx_check_thread_affinity_set(fplog, cr,
815 hw_opt, hwinfo->nthreads_hw_avail, FALSE);
817 #if GMX_THREAD_MPI
818 if (SIMMASTER(cr))
820 if (npme > 0 && hw_opt->nthreads_tmpi <= 0)
822 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
825 /* Since the master knows the cut-off scheme, update hw_opt for this.
826 * This is done later for normal MPI and also once more with tMPI
827 * for all tMPI ranks.
829 check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
831 /* NOW the threads will be started: */
832 hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
833 hw_opt,
834 inputrec, mtop,
835 cr, fplog, bUseGPU);
837 if (hw_opt->nthreads_tmpi > 1)
839 t_commrec *cr_old = cr;
840 /* now start the threads. */
841 cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
842 oenv, bVerbose, nstglobalcomm,
843 ddxyz, dd_rank_order, npme, rdd, rconstr,
844 dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
845 nbpu_opt, nstlist_cmdline,
846 nsteps_cmdline, nstepout, resetstep, nmultisim,
847 repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
848 cpt_period, max_hours,
849 Flags);
850 /* the main thread continues here with a new cr. We don't deallocate
851 the old cr because other threads may still be reading it. */
852 if (cr == NULL)
854 gmx_comm("Failed to spawn threads");
858 #endif
859 /* END OF CAUTION: cr is now reliable */
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 || npme > 0))
886 gmx_fatal(FARGS,
887 "The -dd or -npme option request a parallel simulation, "
888 #if !GMX_MPI
889 "but %s was compiled without threads or MPI enabled"
890 #else
891 #if GMX_THREAD_MPI
892 "but the number of MPI-threads (option -ntmpi) is not set or 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 (npme > 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 npme = 0;
923 if (bUseGPU && npme < 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 npme = 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, &npme,
986 inputrec, state, &bReadEkin,
987 (Flags & MD_APPENDFILES),
988 (Flags & MD_APPENDFILESSET),
989 (Flags & MD_REPRODUCIBLE));
991 if (bReadEkin)
993 Flags |= MD_READ_EKIN;
997 if (MASTER(cr) && (Flags & MD_APPENDFILES))
999 gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
1000 Flags, &fplog);
1003 /* override nsteps with value from cmdline */
1004 override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
1006 if (SIMMASTER(cr))
1008 copy_mat(state->box, box);
1011 if (PAR(cr))
1013 gmx_bcast(sizeof(box), box, cr);
1016 // TODO This should move to do_md(), because it only makes sense
1017 // with dynamical integrators, but there is no test coverage and
1018 // it interacts with constraints, somehow.
1019 /* Essential dynamics */
1020 if (opt2bSet("-ei", nfile, fnm))
1022 /* Open input and output files, allocate space for ED data structure */
1023 ed = ed_open(mtop->natoms, &state->edsamstate, nfile, fnm, Flags, oenv, cr);
1026 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
1027 inputrec->eI == eiNM))
1029 cr->dd = init_domain_decomposition(fplog, cr, Flags, ddxyz, npme,
1030 dd_rank_order,
1031 rdd, rconstr,
1032 dddlb_opt, dlb_scale,
1033 ddcsx, ddcsy, ddcsz,
1034 mtop, inputrec,
1035 box, state->x,
1036 &ddbox, &npme_major, &npme_minor);
1038 else
1040 /* PME, if used, is done on all nodes with 1D decomposition */
1041 cr->npmenodes = 0;
1042 cr->duty = (DUTY_PP | DUTY_PME);
1043 npme_major = 1;
1044 npme_minor = 1;
1046 if (inputrec->ePBC == epbcSCREW)
1048 gmx_fatal(FARGS,
1049 "pbc=%s is only implemented with domain decomposition",
1050 epbc_names[inputrec->ePBC]);
1054 if (PAR(cr))
1056 /* After possible communicator splitting in make_dd_communicators.
1057 * we can set up the intra/inter node communication.
1059 gmx_setup_nodecomm(fplog, cr);
1062 /* Initialize per-physical-node MPI process/thread ID and counters. */
1063 gmx_init_intranode_counters(cr);
1064 #if GMX_MPI
1065 if (MULTISIM(cr))
1067 md_print_info(cr, fplog,
1068 "This is simulation %d out of %d running as a composite GROMACS\n"
1069 "multi-simulation job. Setup for this simulation:\n\n",
1070 cr->ms->sim, cr->ms->nsim);
1072 md_print_info(cr, fplog, "Using %d MPI %s\n",
1073 cr->nnodes,
1074 #if GMX_THREAD_MPI
1075 cr->nnodes == 1 ? "thread" : "threads"
1076 #else
1077 cr->nnodes == 1 ? "process" : "processes"
1078 #endif
1080 fflush(stderr);
1081 #endif
1083 /* Check and update hw_opt for the cut-off scheme */
1084 check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1086 /* Check and update hw_opt for the number of MPI ranks */
1087 check_and_update_hw_opt_3(hw_opt);
1089 gmx_omp_nthreads_init(fplog, cr,
1090 hwinfo->nthreads_hw_avail,
1091 hw_opt->nthreads_omp,
1092 hw_opt->nthreads_omp_pme,
1093 (cr->duty & DUTY_PP) == 0,
1094 inputrec->cutoff_scheme == ecutsVERLET);
1096 #ifndef NDEBUG
1097 if (EI_TPI(inputrec->eI) &&
1098 inputrec->cutoff_scheme == ecutsVERLET)
1100 gmx_feenableexcept();
1102 #endif
1104 if (bUseGPU)
1106 /* Select GPU id's to use */
1107 gmx_select_gpu_ids(fplog, cr, &hwinfo->gpu_info, bForceUseGPU,
1108 &hw_opt->gpu_opt);
1110 else
1112 /* Ignore (potentially) manually selected GPUs */
1113 hw_opt->gpu_opt.n_dev_use = 0;
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 // Membrane embedding must be initialized before we call init_forcerec()
1148 if (opt2bSet("-membed", nfile, fnm))
1150 if (MASTER(cr))
1152 fprintf(stderr, "Initializing membed");
1154 /* Note that membed cannot work in parallel because mtop is
1155 * changed here. Fix this if we ever want to make it run with
1156 * multiple ranks. */
1157 membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
1160 snew(nrnb, 1);
1161 if (cr->duty & DUTY_PP)
1163 bcast_state(cr, state);
1165 /* Initiate forcerecord */
1166 fr = mk_forcerec();
1167 fr->hwinfo = hwinfo;
1168 fr->gpu_opt = &hw_opt->gpu_opt;
1169 init_forcerec(fplog, fr, fcd, inputrec, mtop, cr, box,
1170 opt2fn("-table", nfile, fnm),
1171 opt2fn("-tablep", nfile, fnm),
1172 getFilenm("-tableb", nfile, fnm),
1173 nbpu_opt,
1174 FALSE,
1175 pforce);
1177 /* Initialize QM-MM */
1178 if (fr->bQMMM)
1180 init_QMMMrec(cr, mtop, inputrec, fr);
1183 /* Initialize the mdatoms structure.
1184 * mdatoms is not filled with atom data,
1185 * as this can not be done now with domain decomposition.
1187 mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1189 /* Initialize the virtual site communication */
1190 vsite = init_vsite(mtop, cr, FALSE);
1192 calc_shifts(box, fr->shift_vec);
1194 /* With periodic molecules the charge groups should be whole at start up
1195 * and the virtual sites should not be far from their proper positions.
1197 if (!inputrec->bContinuation && MASTER(cr) &&
1198 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1200 /* Make molecules whole at start of run */
1201 if (fr->ePBC != epbcNONE)
1203 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, state->x);
1205 if (vsite)
1207 /* Correct initial vsite positions are required
1208 * for the initial distribution in the domain decomposition
1209 * and for the initial shell prediction.
1211 construct_vsites_mtop(vsite, mtop, state->x);
1215 if (EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype))
1217 ewaldcoeff_q = fr->ewaldcoeff_q;
1218 ewaldcoeff_lj = fr->ewaldcoeff_lj;
1219 pmedata = &fr->pmedata;
1221 else
1223 pmedata = NULL;
1226 else
1228 /* This is a PME only node */
1230 /* We don't need the state */
1231 done_state(state);
1233 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1234 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1235 snew(pmedata, 1);
1238 if (hw_opt->thread_affinity != threadaffOFF)
1240 /* Before setting affinity, check whether the affinity has changed
1241 * - which indicates that probably the OpenMP library has changed it
1242 * since we first checked).
1244 gmx_check_thread_affinity_set(fplog, cr,
1245 hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1247 /* Set the CPU affinity */
1248 gmx_set_thread_affinity(fplog, cr, hw_opt, hwinfo);
1251 /* Initiate PME if necessary,
1252 * either on all nodes or on dedicated PME nodes only. */
1253 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1255 if (mdatoms)
1257 nChargePerturbed = mdatoms->nChargePerturbed;
1258 if (EVDW_PME(inputrec->vdwtype))
1260 nTypePerturbed = mdatoms->nTypePerturbed;
1263 if (cr->npmenodes > 0)
1265 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1266 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1267 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1270 if (cr->duty & DUTY_PME)
1272 status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1273 mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1274 (Flags & MD_REPRODUCIBLE), nthreads_pme);
1275 if (status != 0)
1277 gmx_fatal(FARGS, "Error %d initializing PME", status);
1283 if (EI_DYNAMICS(inputrec->eI))
1285 /* Turn on signal handling on all nodes */
1287 * (A user signal from the PME nodes (if any)
1288 * is communicated to the PP nodes.
1290 signal_handler_install();
1293 if (cr->duty & DUTY_PP)
1295 /* Assumes uniform use of the number of OpenMP threads */
1296 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1298 if (inputrec->bPull)
1300 /* Initialize pull code */
1301 inputrec->pull_work =
1302 init_pull(fplog, inputrec->pull, inputrec, nfile, fnm,
1303 mtop, cr, oenv, inputrec->fepvals->init_lambda,
1304 EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1307 if (inputrec->bRot)
1309 /* Initialize enforced rotation code */
1310 init_rot(fplog, inputrec, nfile, fnm, cr, state->x, state->box, mtop, oenv,
1311 bVerbose, Flags);
1314 constr = init_constraints(fplog, mtop, inputrec, ed, state, cr);
1316 if (DOMAINDECOMP(cr))
1318 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1319 /* This call is not included in init_domain_decomposition mainly
1320 * because fr->cginfo_mb is set later.
1322 dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1323 Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1326 /* Now do whatever the user wants us to do (how flexible...) */
1327 my_integrator(inputrec->eI) (fplog, cr, nfile, fnm,
1328 oenv, bVerbose,
1329 nstglobalcomm,
1330 vsite, constr,
1331 nstepout, inputrec, mtop,
1332 fcd, state,
1333 mdatoms, nrnb, wcycle, ed, fr,
1334 repl_ex_nst, repl_ex_nex, repl_ex_seed,
1335 membed,
1336 cpt_period, max_hours,
1337 imdport,
1338 Flags,
1339 walltime_accounting);
1341 if (inputrec->bRot)
1343 finish_rot(inputrec->rot);
1346 if (inputrec->bPull)
1348 finish_pull(inputrec->pull_work);
1352 else
1354 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1355 /* do PME only */
1356 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1357 gmx_pmeonly(*pmedata, cr, nrnb, wcycle, walltime_accounting, ewaldcoeff_q, ewaldcoeff_lj, inputrec);
1360 wallcycle_stop(wcycle, ewcRUN);
1362 /* Finish up, write some stuff
1363 * if rerunMD, don't write last frame again
1365 finish_run(fplog, cr,
1366 inputrec, nrnb, wcycle, walltime_accounting,
1367 fr ? fr->nbv : NULL,
1368 EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1371 /* Free GPU memory and context */
1372 free_gpu_resources(fr, cr, &hwinfo->gpu_info, fr ? fr->gpu_opt : NULL);
1374 if (membed != nullptr)
1376 free_membed(membed);
1379 gmx_hardware_info_free(hwinfo);
1381 /* Does what it says */
1382 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1383 walltime_accounting_destroy(walltime_accounting);
1385 /* Close logfile already here if we were appending to it */
1386 if (MASTER(cr) && (Flags & MD_APPENDFILES))
1388 gmx_log_close(fplog);
1391 rc = (int)gmx_get_stop_condition();
1393 done_ed(&ed);
1395 #if GMX_THREAD_MPI
1396 /* we need to join all threads. The sub-threads join when they
1397 exit this function, but the master thread needs to be told to
1398 wait for that. */
1399 if (PAR(cr) && MASTER(cr))
1401 tMPI_Finalize();
1403 #endif
1405 return rc;
1408 } // namespace gmx