Split filenm.*
[gromacs/AngularHB.git] / src / programs / mdrun / runner.cpp
blob0e707ef38ff1a24b1febef45418c4f04989e2b9a
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/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/essentialdynamics/edsam.h"
60 #include "gromacs/ewald/pme.h"
61 #include "gromacs/fileio/checkpoint.h"
62 #include "gromacs/fileio/copyrite.h"
63 #include "gromacs/fileio/filenm.h"
64 #include "gromacs/fileio/oenv.h"
65 #include "gromacs/fileio/tpxio.h"
66 #include "gromacs/fileio/trx.h"
67 #include "gromacs/fileio/txtdump.h"
68 #include "gromacs/gmxlib/disre.h"
69 #include "gromacs/gmxlib/gmx_detect_hardware.h"
70 #include "gromacs/gmxlib/gmx_omp_nthreads.h"
71 #include "gromacs/gmxlib/main.h"
72 #include "gromacs/gmxlib/md_logging.h"
73 #include "gromacs/gmxlib/network.h"
74 #include "gromacs/gmxlib/orires.h"
75 #include "gromacs/gmxlib/sighandler.h"
76 #include "gromacs/gmxlib/thread_affinity.h"
77 #include "gromacs/gmxlib/gpu_utils/gpu_utils.h"
78 #include "gromacs/math/calculate-ewald-splitting-coefficient.h"
79 #include "gromacs/math/vec.h"
80 #include "gromacs/mdlib/calc_verletbuf.h"
81 #include "gromacs/mdlib/constr.h"
82 #include "gromacs/mdlib/force.h"
83 #include "gromacs/mdlib/forcerec.h"
84 #include "gromacs/mdlib/integrator.h"
85 #include "gromacs/mdlib/md_support.h"
86 #include "gromacs/mdlib/mdatoms.h"
87 #include "gromacs/mdlib/mdrun.h"
88 #include "gromacs/mdlib/minimize.h"
89 #include "gromacs/mdlib/nbnxn_search.h"
90 #include "gromacs/mdlib/qmmm.h"
91 #include "gromacs/mdlib/tpi.h"
92 #include "gromacs/mdtypes/inputrec.h"
93 #include "gromacs/mdtypes/md_enums.h"
94 #include "gromacs/mdtypes/state.h"
95 #include "gromacs/pbcutil/pbc.h"
96 #include "gromacs/pulling/pull.h"
97 #include "gromacs/pulling/pull_rotation.h"
98 #include "gromacs/swap/swapcoords.h"
99 #include "gromacs/timing/wallcycle.h"
100 #include "gromacs/topology/mtop_util.h"
101 #include "gromacs/utility/cstringutil.h"
102 #include "gromacs/utility/exceptions.h"
103 #include "gromacs/utility/fatalerror.h"
104 #include "gromacs/utility/gmxassert.h"
105 #include "gromacs/utility/gmxmpi.h"
106 #include "gromacs/utility/smalloc.h"
108 #include "deform.h"
109 #include "md.h"
110 #include "membed.h"
111 #include "repl_ex.h"
112 #include "resource-division.h"
114 #ifdef GMX_FAHCORE
115 #include "corewrap.h"
116 #endif
118 //! First step used in pressure scaling
119 gmx_int64_t deform_init_init_step_tpx;
120 //! Initial box for pressure scaling
121 matrix deform_init_box_tpx;
122 //! MPI variable for use in pressure scaling
123 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
125 #ifdef GMX_THREAD_MPI
126 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
127 * the number of threads will get lowered.
129 #define MIN_ATOMS_PER_MPI_THREAD 90
130 #define MIN_ATOMS_PER_GPU 900
132 struct mdrunner_arglist
134 gmx_hw_opt_t hw_opt;
135 FILE *fplog;
136 t_commrec *cr;
137 int nfile;
138 const t_filenm *fnm;
139 const gmx_output_env_t *oenv;
140 gmx_bool bVerbose;
141 gmx_bool bCompact;
142 int nstglobalcomm;
143 ivec ddxyz;
144 int dd_node_order;
145 real rdd;
146 real rconstr;
147 const char *dddlb_opt;
148 real dlb_scale;
149 const char *ddcsx;
150 const char *ddcsy;
151 const char *ddcsz;
152 const char *nbpu_opt;
153 int nstlist_cmdline;
154 gmx_int64_t nsteps_cmdline;
155 int nstepout;
156 int resetstep;
157 int nmultisim;
158 int repl_ex_nst;
159 int repl_ex_nex;
160 int repl_ex_seed;
161 real pforce;
162 real cpt_period;
163 real max_hours;
164 int imdport;
165 unsigned long Flags;
169 /* The function used for spawning threads. Extracts the mdrunner()
170 arguments from its one argument and calls mdrunner(), after making
171 a commrec. */
172 static void mdrunner_start_fn(void *arg)
176 struct mdrunner_arglist *mda = (struct mdrunner_arglist*)arg;
177 struct mdrunner_arglist mc = *mda; /* copy the arg list to make sure
178 that it's thread-local. This doesn't
179 copy pointed-to items, of course,
180 but those are all const. */
181 t_commrec *cr; /* we need a local version of this */
182 FILE *fplog = NULL;
183 t_filenm *fnm;
185 fnm = dup_tfn(mc.nfile, mc.fnm);
187 cr = reinitialize_commrec_for_this_thread(mc.cr);
189 if (MASTER(cr))
191 fplog = mc.fplog;
194 gmx::mdrunner(&mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
195 mc.bVerbose, mc.bCompact, mc.nstglobalcomm,
196 mc.ddxyz, mc.dd_node_order, mc.rdd,
197 mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
198 mc.ddcsx, mc.ddcsy, mc.ddcsz,
199 mc.nbpu_opt, mc.nstlist_cmdline,
200 mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
201 mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce,
202 mc.cpt_period, mc.max_hours, mc.imdport, mc.Flags);
204 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
208 /* called by mdrunner() to start a specific number of threads (including
209 the main thread) for thread-parallel runs. This in turn calls mdrunner()
210 for each thread.
211 All options besides nthreads are the same as for mdrunner(). */
212 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
213 FILE *fplog, t_commrec *cr, int nfile,
214 const t_filenm fnm[], const gmx_output_env_t *oenv, gmx_bool bVerbose,
215 gmx_bool bCompact, int nstglobalcomm,
216 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
217 const char *dddlb_opt, real dlb_scale,
218 const char *ddcsx, const char *ddcsy, const char *ddcsz,
219 const char *nbpu_opt, int nstlist_cmdline,
220 gmx_int64_t nsteps_cmdline,
221 int nstepout, int resetstep,
222 int nmultisim, int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
223 real pforce, real cpt_period, real max_hours,
224 unsigned long Flags)
226 int ret;
227 struct mdrunner_arglist *mda;
228 t_commrec *crn; /* the new commrec */
229 t_filenm *fnmn;
231 /* first check whether we even need to start tMPI */
232 if (hw_opt->nthreads_tmpi < 2)
234 return cr;
237 /* a few small, one-time, almost unavoidable memory leaks: */
238 snew(mda, 1);
239 fnmn = dup_tfn(nfile, fnm);
241 /* fill the data structure to pass as void pointer to thread start fn */
242 /* hw_opt contains pointers, which should all be NULL at this stage */
243 mda->hw_opt = *hw_opt;
244 mda->fplog = fplog;
245 mda->cr = cr;
246 mda->nfile = nfile;
247 mda->fnm = fnmn;
248 mda->oenv = oenv;
249 mda->bVerbose = bVerbose;
250 mda->bCompact = bCompact;
251 mda->nstglobalcomm = nstglobalcomm;
252 mda->ddxyz[XX] = ddxyz[XX];
253 mda->ddxyz[YY] = ddxyz[YY];
254 mda->ddxyz[ZZ] = ddxyz[ZZ];
255 mda->dd_node_order = dd_node_order;
256 mda->rdd = rdd;
257 mda->rconstr = rconstr;
258 mda->dddlb_opt = dddlb_opt;
259 mda->dlb_scale = dlb_scale;
260 mda->ddcsx = ddcsx;
261 mda->ddcsy = ddcsy;
262 mda->ddcsz = ddcsz;
263 mda->nbpu_opt = nbpu_opt;
264 mda->nstlist_cmdline = nstlist_cmdline;
265 mda->nsteps_cmdline = nsteps_cmdline;
266 mda->nstepout = nstepout;
267 mda->resetstep = resetstep;
268 mda->nmultisim = nmultisim;
269 mda->repl_ex_nst = repl_ex_nst;
270 mda->repl_ex_nex = repl_ex_nex;
271 mda->repl_ex_seed = repl_ex_seed;
272 mda->pforce = pforce;
273 mda->cpt_period = cpt_period;
274 mda->max_hours = max_hours;
275 mda->Flags = Flags;
277 /* now spawn new threads that start mdrunner_start_fn(), while
278 the main thread returns, we set thread affinity later */
279 ret = tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi, TMPI_AFFINITY_NONE,
280 mdrunner_start_fn, (void*)(mda) );
281 if (ret != TMPI_SUCCESS)
283 return NULL;
286 crn = reinitialize_commrec_for_this_thread(cr);
287 return crn;
290 #endif /* GMX_THREAD_MPI */
293 /*! \brief Cost of non-bonded kernels
295 * We determine the extra cost of the non-bonded kernels compared to
296 * a reference nstlist value of 10 (which is the default in grompp).
298 static const int nbnxnReferenceNstlist = 10;
299 //! The values to try when switching
300 const int nstlist_try[] = { 20, 25, 40 };
301 //! Number of elements in the neighborsearch list trials.
302 #define NNSTL sizeof(nstlist_try)/sizeof(nstlist_try[0])
303 /* Increase nstlist until the non-bonded cost increases more than listfac_ok,
304 * but never more than listfac_max.
305 * A standard (protein+)water system at 300K with PME ewald_rtol=1e-5
306 * needs 1.28 at rcoulomb=0.9 and 1.24 at rcoulomb=1.0 to get to nstlist=40.
307 * Note that both CPU and GPU factors are conservative. Performance should
308 * not go down due to this tuning, except with a relatively slow GPU.
309 * On the other hand, at medium/high parallelization or with fast GPUs
310 * nstlist will not be increased enough to reach optimal performance.
312 /* CPU: pair-search is about a factor 1.5 slower than the non-bonded kernel */
313 //! Max OK performance ratio beween force calc and neighbor searching
314 static const float nbnxn_cpu_listfac_ok = 1.05;
315 //! Too high performance ratio beween force calc and neighbor searching
316 static const float nbnxn_cpu_listfac_max = 1.09;
317 /* GPU: pair-search is a factor 1.5-3 slower than the non-bonded kernel */
318 //! Max OK performance ratio beween force calc and neighbor searching
319 static const float nbnxn_gpu_listfac_ok = 1.20;
320 //! Too high performance ratio beween force calc and neighbor searching
321 static const float nbnxn_gpu_listfac_max = 1.30;
323 /*! \brief Try to increase nstlist when using the Verlet cut-off scheme */
324 static void increase_nstlist(FILE *fp, t_commrec *cr,
325 t_inputrec *ir, int nstlist_cmdline,
326 const gmx_mtop_t *mtop, matrix box,
327 gmx_bool bGPU)
329 float listfac_ok, listfac_max;
330 int nstlist_orig, nstlist_prev;
331 verletbuf_list_setup_t ls;
332 real rlistWithReferenceNstlist, rlist_inc, rlist_ok, rlist_max;
333 real rlist_new, rlist_prev;
334 size_t nstlist_ind = 0;
335 t_state state_tmp;
336 gmx_bool bBox, bDD, bCont;
337 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";
338 const char *nve_err = "Can not increase nstlist because an NVE ensemble is used";
339 const char *vbd_err = "Can not increase nstlist because verlet-buffer-tolerance is not set or used";
340 const char *box_err = "Can not increase nstlist because the box is too small";
341 const char *dd_err = "Can not increase nstlist because of domain decomposition limitations";
342 char buf[STRLEN];
343 const float oneThird = 1.0f / 3.0f;
345 if (nstlist_cmdline <= 0)
347 if (ir->nstlist == 1)
349 /* The user probably set nstlist=1 for a reason,
350 * don't mess with the settings.
352 return;
355 if (fp != NULL && bGPU && ir->nstlist < nstlist_try[0])
357 fprintf(fp, nstl_gpu, ir->nstlist);
359 nstlist_ind = 0;
360 while (nstlist_ind < NNSTL && ir->nstlist >= nstlist_try[nstlist_ind])
362 nstlist_ind++;
364 if (nstlist_ind == NNSTL)
366 /* There are no larger nstlist value to try */
367 return;
371 if (EI_MD(ir->eI) && ir->etc == etcNO)
373 if (MASTER(cr))
375 fprintf(stderr, "%s\n", nve_err);
377 if (fp != NULL)
379 fprintf(fp, "%s\n", nve_err);
382 return;
385 if (ir->verletbuf_tol == 0 && bGPU)
387 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");
390 if (ir->verletbuf_tol < 0)
392 if (MASTER(cr))
394 fprintf(stderr, "%s\n", vbd_err);
396 if (fp != NULL)
398 fprintf(fp, "%s\n", vbd_err);
401 return;
404 if (bGPU)
406 listfac_ok = nbnxn_gpu_listfac_ok;
407 listfac_max = nbnxn_gpu_listfac_max;
409 else
411 listfac_ok = nbnxn_cpu_listfac_ok;
412 listfac_max = nbnxn_cpu_listfac_max;
415 nstlist_orig = ir->nstlist;
416 if (nstlist_cmdline > 0)
418 if (fp)
420 sprintf(buf, "Getting nstlist=%d from command line option",
421 nstlist_cmdline);
423 ir->nstlist = nstlist_cmdline;
426 verletbuf_get_list_setup(TRUE, bGPU, &ls);
428 /* Allow rlist to make the list a given factor larger than the list
429 * would be with the reference value for nstlist (10).
431 nstlist_prev = ir->nstlist;
432 ir->nstlist = nbnxnReferenceNstlist;
433 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL,
434 &rlistWithReferenceNstlist);
435 ir->nstlist = nstlist_prev;
437 /* Determine the pair list size increase due to zero interactions */
438 rlist_inc = nbnxn_get_rlist_effective_inc(ls.cluster_size_j,
439 mtop->natoms/det(box));
440 rlist_ok = (rlistWithReferenceNstlist + rlist_inc)*pow(listfac_ok, oneThird) - rlist_inc;
441 rlist_max = (rlistWithReferenceNstlist + rlist_inc)*pow(listfac_max, oneThird) - rlist_inc;
442 if (debug)
444 fprintf(debug, "nstlist tuning: rlist_inc %.3f rlist_ok %.3f rlist_max %.3f\n",
445 rlist_inc, rlist_ok, rlist_max);
448 nstlist_prev = nstlist_orig;
449 rlist_prev = ir->rlist;
452 if (nstlist_cmdline <= 0)
454 ir->nstlist = nstlist_try[nstlist_ind];
457 /* Set the pair-list buffer size in ir */
458 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
460 /* Does rlist fit in the box? */
461 bBox = (sqr(rlist_new) < max_cutoff2(ir->ePBC, box));
462 bDD = TRUE;
463 if (bBox && DOMAINDECOMP(cr))
465 /* Check if rlist fits in the domain decomposition */
466 if (inputrec2nboundeddim(ir) < DIM)
468 gmx_incons("Changing nstlist with domain decomposition and unbounded dimensions is not implemented yet");
470 copy_mat(box, state_tmp.box);
471 bDD = change_dd_cutoff(cr, &state_tmp, ir, rlist_new);
474 if (debug)
476 fprintf(debug, "nstlist %d rlist %.3f bBox %d bDD %d\n",
477 ir->nstlist, rlist_new, bBox, bDD);
480 bCont = FALSE;
482 if (nstlist_cmdline <= 0)
484 if (bBox && bDD && rlist_new <= rlist_max)
486 /* Increase nstlist */
487 nstlist_prev = ir->nstlist;
488 rlist_prev = rlist_new;
489 bCont = (nstlist_ind+1 < NNSTL && rlist_new < rlist_ok);
491 else
493 /* Stick with the previous nstlist */
494 ir->nstlist = nstlist_prev;
495 rlist_new = rlist_prev;
496 bBox = TRUE;
497 bDD = TRUE;
501 nstlist_ind++;
503 while (bCont);
505 if (!bBox || !bDD)
507 gmx_warning(!bBox ? box_err : dd_err);
508 if (fp != NULL)
510 fprintf(fp, "\n%s\n", bBox ? box_err : dd_err);
512 ir->nstlist = nstlist_orig;
514 else if (ir->nstlist != nstlist_orig || rlist_new != ir->rlist)
516 sprintf(buf, "Changing nstlist from %d to %d, rlist from %g to %g",
517 nstlist_orig, ir->nstlist,
518 ir->rlist, rlist_new);
519 if (MASTER(cr))
521 fprintf(stderr, "%s\n\n", buf);
523 if (fp != NULL)
525 fprintf(fp, "%s\n\n", buf);
527 ir->rlist = rlist_new;
531 /*! \brief Initialize variables for Verlet scheme simulation */
532 static void prepare_verlet_scheme(FILE *fplog,
533 t_commrec *cr,
534 t_inputrec *ir,
535 int nstlist_cmdline,
536 const gmx_mtop_t *mtop,
537 matrix box,
538 gmx_bool bUseGPU)
540 /* For NVE simulations, we will retain the initial list buffer */
541 if (ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
543 /* Update the Verlet buffer size for the current run setup */
544 verletbuf_list_setup_t ls;
545 real rlist_new;
547 /* Here we assume SIMD-enabled kernels are being used. But as currently
548 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
549 * and 4x2 gives a larger buffer than 4x4, this is ok.
551 verletbuf_get_list_setup(TRUE, bUseGPU, &ls);
553 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
555 if (rlist_new != ir->rlist)
557 if (fplog != NULL)
559 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
560 ir->rlist, rlist_new,
561 ls.cluster_size_i, ls.cluster_size_j);
563 ir->rlist = rlist_new;
567 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
569 gmx_fatal(FARGS, "Can not set nstlist without %s",
570 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
573 if (EI_DYNAMICS(ir->eI))
575 /* Set or try nstlist values */
576 increase_nstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, bUseGPU);
580 /*! \brief Override the nslist value in inputrec
582 * with value passed on the command line (if any)
584 static void override_nsteps_cmdline(FILE *fplog,
585 gmx_int64_t nsteps_cmdline,
586 t_inputrec *ir,
587 const t_commrec *cr)
589 assert(ir);
590 assert(cr);
592 /* override with anything else than the default -2 */
593 if (nsteps_cmdline > -2)
595 char sbuf_steps[STEPSTRSIZE];
596 char sbuf_msg[STRLEN];
598 ir->nsteps = nsteps_cmdline;
599 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
601 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
602 gmx_step_str(nsteps_cmdline, sbuf_steps),
603 fabs(nsteps_cmdline*ir->delta_t));
605 else
607 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
608 gmx_step_str(nsteps_cmdline, sbuf_steps));
611 md_print_warn(cr, fplog, "%s\n", sbuf_msg);
613 else if (nsteps_cmdline < -2)
615 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %d",
616 nsteps_cmdline);
618 /* Do nothing if nsteps_cmdline == -2 */
621 namespace gmx
624 //! \brief Return the correct integrator function.
625 static integrator_t *my_integrator(unsigned int ei)
627 switch (ei)
629 case eiMD:
630 case eiBD:
631 case eiSD2:
632 case eiSD1:
633 case eiVV:
634 case eiVVAK:
635 if (!EI_DYNAMICS(ei))
637 GMX_THROW(APIError("do_md integrator would be called for a non-dynamical integrator"));
639 return do_md;
640 case eiSteep:
641 return do_steep;
642 case eiCG:
643 return do_cg;
644 case eiNM:
645 return do_nm;
646 case eiLBFGS:
647 return do_lbfgs;
648 case eiTPI:
649 case eiTPIC:
650 if (!EI_TPI(ei))
652 GMX_THROW(APIError("do_tpi integrator would be called for a non-TPI integrator"));
654 return do_tpi;
655 default:
656 GMX_THROW(APIError("Non existing integrator selected"));
660 int mdrunner(gmx_hw_opt_t *hw_opt,
661 FILE *fplog, t_commrec *cr, int nfile,
662 const t_filenm fnm[], const gmx_output_env_t *oenv, gmx_bool bVerbose,
663 gmx_bool bCompact, int nstglobalcomm,
664 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
665 const char *dddlb_opt, real dlb_scale,
666 const char *ddcsx, const char *ddcsy, const char *ddcsz,
667 const char *nbpu_opt, int nstlist_cmdline,
668 gmx_int64_t nsteps_cmdline, int nstepout, int resetstep,
669 int gmx_unused nmultisim, int repl_ex_nst, int repl_ex_nex,
670 int repl_ex_seed, real pforce, real cpt_period, real max_hours,
671 int imdport, unsigned long Flags)
673 gmx_bool bForceUseGPU, bTryUseGPU, bRerunMD;
674 t_inputrec *inputrec;
675 t_state *state = NULL;
676 matrix box;
677 gmx_ddbox_t ddbox = {0};
678 int npme_major, npme_minor;
679 t_nrnb *nrnb;
680 gmx_mtop_t *mtop = NULL;
681 t_mdatoms *mdatoms = NULL;
682 t_forcerec *fr = NULL;
683 t_fcdata *fcd = NULL;
684 real ewaldcoeff_q = 0;
685 real ewaldcoeff_lj = 0;
686 struct gmx_pme_t **pmedata = NULL;
687 gmx_vsite_t *vsite = NULL;
688 gmx_constr_t constr;
689 int nChargePerturbed = -1, nTypePerturbed = 0, status;
690 gmx_wallcycle_t wcycle;
691 gmx_walltime_accounting_t walltime_accounting = NULL;
692 int rc;
693 gmx_int64_t reset_counters;
694 gmx_edsam_t ed = NULL;
695 int nthreads_pme = 1;
696 gmx_membed_t *membed = NULL;
697 gmx_hw_info_t *hwinfo = NULL;
698 /* The master rank decides early on bUseGPU and broadcasts this later */
699 gmx_bool bUseGPU = FALSE;
701 /* CAUTION: threads may be started later on in this function, so
702 cr doesn't reflect the final parallel state right now */
703 snew(inputrec, 1);
704 snew(mtop, 1);
706 if (Flags & MD_APPENDFILES)
708 fplog = NULL;
711 bRerunMD = (Flags & MD_RERUN);
712 bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
713 bTryUseGPU = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
715 /* Detect hardware, gather information. This is an operation that is
716 * global for this process (MPI rank). */
717 hwinfo = gmx_detect_hardware(fplog, cr, bTryUseGPU);
719 gmx_print_detected_hardware(fplog, cr, hwinfo);
721 if (fplog != NULL)
723 /* Print references after all software/hardware printing */
724 please_cite(fplog, "Abraham2015");
725 please_cite(fplog, "Pall2015");
726 please_cite(fplog, "Pronk2013");
727 please_cite(fplog, "Hess2008b");
728 please_cite(fplog, "Spoel2005a");
729 please_cite(fplog, "Lindahl2001a");
730 please_cite(fplog, "Berendsen95a");
733 snew(state, 1);
734 if (SIMMASTER(cr))
736 /* Read (nearly) all data required for the simulation */
737 read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, state, mtop);
739 if (inputrec->cutoff_scheme == ecutsVERLET)
741 /* Here the master rank decides if all ranks will use GPUs */
742 bUseGPU = (hwinfo->gpu_info.n_dev_compatible > 0 ||
743 getenv("GMX_EMULATE_GPU") != NULL);
745 /* TODO add GPU kernels for this and replace this check by:
746 * (bUseGPU && (ir->vdwtype == evdwPME &&
747 * ir->ljpme_combination_rule == eljpmeLB))
748 * update the message text and the content of nbnxn_acceleration_supported.
750 if (bUseGPU &&
751 !nbnxn_gpu_acceleration_supported(fplog, cr, inputrec, bRerunMD))
753 /* Fallback message printed by nbnxn_acceleration_supported */
754 if (bForceUseGPU)
756 gmx_fatal(FARGS, "GPU acceleration requested, but not supported with the given input settings");
758 bUseGPU = FALSE;
761 prepare_verlet_scheme(fplog, cr,
762 inputrec, nstlist_cmdline, mtop, state->box,
763 bUseGPU);
765 else
767 if (nstlist_cmdline > 0)
769 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
772 if (hwinfo->gpu_info.n_dev_compatible > 0)
774 md_print_warn(cr, fplog,
775 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
776 " To use a GPU, set the mdp option: cutoff-scheme = Verlet\n");
779 if (bForceUseGPU)
781 gmx_fatal(FARGS, "GPU requested, but can't be used without cutoff-scheme=Verlet");
784 #ifdef GMX_TARGET_BGQ
785 md_print_warn(cr, fplog,
786 "NOTE: There is no SIMD implementation of the group scheme kernels on\n"
787 " BlueGene/Q. You will observe better performance from using the\n"
788 " Verlet cut-off scheme.\n");
789 #endif
792 if (inputrec->eI == eiSD2)
794 md_print_warn(cr, fplog, "The stochastic dynamics integrator %s is deprecated, since\n"
795 "it is slower than integrator %s and is slightly less accurate\n"
796 "with constraints. Use the %s integrator.",
797 ei_names[inputrec->eI], ei_names[eiSD1], ei_names[eiSD1]);
801 /* Check and update the hardware options for internal consistency */
802 check_and_update_hw_opt_1(hw_opt, cr);
804 /* Early check for externally set process affinity. */
805 gmx_check_thread_affinity_set(fplog, cr,
806 hw_opt, hwinfo->nthreads_hw_avail, FALSE);
808 #ifdef GMX_THREAD_MPI
809 if (SIMMASTER(cr))
811 if (cr->npmenodes > 0 && hw_opt->nthreads_tmpi <= 0)
813 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
816 /* Since the master knows the cut-off scheme, update hw_opt for this.
817 * This is done later for normal MPI and also once more with tMPI
818 * for all tMPI ranks.
820 check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
822 /* NOW the threads will be started: */
823 hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
824 hw_opt,
825 inputrec, mtop,
826 cr, fplog, bUseGPU);
828 if (hw_opt->nthreads_tmpi > 1)
830 t_commrec *cr_old = cr;
831 /* now start the threads. */
832 cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
833 oenv, bVerbose, bCompact, nstglobalcomm,
834 ddxyz, dd_node_order, rdd, rconstr,
835 dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
836 nbpu_opt, nstlist_cmdline,
837 nsteps_cmdline, nstepout, resetstep, nmultisim,
838 repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
839 cpt_period, max_hours,
840 Flags);
841 /* the main thread continues here with a new cr. We don't deallocate
842 the old cr because other threads may still be reading it. */
843 if (cr == NULL)
845 gmx_comm("Failed to spawn threads");
849 #endif
850 /* END OF CAUTION: cr is now reliable */
852 /* g_membed initialisation *
853 * Because we change the mtop, init_membed is called before the init_parallel *
854 * (in case we ever want to make it run in parallel) */
855 if (opt2bSet("-membed", nfile, fnm))
857 if (MASTER(cr))
859 fprintf(stderr, "Initializing membed");
861 membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
864 if (PAR(cr))
866 /* now broadcast everything to the non-master nodes/threads: */
867 init_parallel(cr, inputrec, mtop);
869 /* The master rank decided on the use of GPUs,
870 * broadcast this information to all ranks.
872 gmx_bcast_sim(sizeof(bUseGPU), &bUseGPU, cr);
875 if (fplog != NULL)
877 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
878 fprintf(fplog, "\n");
881 /* now make sure the state is initialized and propagated */
882 set_state_entries(state, inputrec);
884 /* A parallel command line option consistency check that we can
885 only do after any threads have started. */
886 if (!PAR(cr) &&
887 (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
889 gmx_fatal(FARGS,
890 "The -dd or -npme option request a parallel simulation, "
891 #ifndef GMX_MPI
892 "but %s was compiled without threads or MPI enabled"
893 #else
894 #ifdef GMX_THREAD_MPI
895 "but the number of threads (option -nt) is 1"
896 #else
897 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
898 #endif
899 #endif
900 , output_env_get_program_display_name(oenv)
904 if (bRerunMD &&
905 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
907 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
910 if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
912 gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
915 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
917 if (cr->npmenodes > 0)
919 gmx_fatal_collective(FARGS, cr->mpi_comm_mysim, MASTER(cr),
920 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
923 cr->npmenodes = 0;
926 if (bUseGPU && cr->npmenodes < 0)
928 /* With GPUs we don't automatically use PME-only ranks. PME ranks can
929 * improve performance with many threads per GPU, since our OpenMP
930 * scaling is bad, but it's difficult to automate the setup.
932 cr->npmenodes = 0;
935 #ifdef GMX_FAHCORE
936 if (MASTER(cr))
938 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
940 #endif
942 /* NMR restraints must be initialized before load_checkpoint,
943 * since with time averaging the history is added to t_state.
944 * For proper consistency check we therefore need to extend
945 * t_state here.
946 * So the PME-only nodes (if present) will also initialize
947 * the distance restraints.
949 snew(fcd, 1);
951 /* This needs to be called before read_checkpoint to extend the state */
952 init_disres(fplog, mtop, inputrec, cr, fcd, state, repl_ex_nst > 0);
954 init_orires(fplog, mtop, state->x, inputrec, cr, &(fcd->orires),
955 state);
957 if (inputrecDeform(inputrec))
959 /* Store the deform reference box before reading the checkpoint */
960 if (SIMMASTER(cr))
962 copy_mat(state->box, box);
964 if (PAR(cr))
966 gmx_bcast(sizeof(box), box, cr);
968 /* Because we do not have the update struct available yet
969 * in which the reference values should be stored,
970 * we store them temporarily in static variables.
971 * This should be thread safe, since they are only written once
972 * and with identical values.
974 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
975 deform_init_init_step_tpx = inputrec->init_step;
976 copy_mat(box, deform_init_box_tpx);
977 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
980 if (Flags & MD_STARTFROMCPT)
982 /* Check if checkpoint file exists before doing continuation.
983 * This way we can use identical input options for the first and subsequent runs...
985 gmx_bool bReadEkin;
987 load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
988 cr, ddxyz,
989 inputrec, state, &bReadEkin,
990 (Flags & MD_APPENDFILES),
991 (Flags & MD_APPENDFILESSET));
993 if (bReadEkin)
995 Flags |= MD_READ_EKIN;
999 if (MASTER(cr) && (Flags & MD_APPENDFILES))
1001 gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
1002 Flags, &fplog);
1005 /* override nsteps with value from cmdline */
1006 override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
1008 if (SIMMASTER(cr))
1010 copy_mat(state->box, box);
1013 if (PAR(cr))
1015 gmx_bcast(sizeof(box), box, cr);
1018 /* Essential dynamics */
1019 if (opt2bSet("-ei", nfile, fnm))
1021 /* Open input and output files, allocate space for ED data structure */
1022 ed = ed_open(mtop->natoms, &state->edsamstate, nfile, fnm, Flags, oenv, cr);
1025 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
1026 inputrec->eI == eiNM))
1028 cr->dd = init_domain_decomposition(fplog, cr, Flags, ddxyz, rdd, rconstr,
1029 dddlb_opt, dlb_scale,
1030 ddcsx, ddcsy, ddcsz,
1031 mtop, inputrec,
1032 box, state->x,
1033 &ddbox, &npme_major, &npme_minor);
1035 make_dd_communicators(fplog, cr, dd_node_order);
1037 /* Set overallocation to avoid frequent reallocation of arrays */
1038 set_over_alloc_dd(TRUE);
1040 else
1042 /* PME, if used, is done on all nodes with 1D decomposition */
1043 cr->npmenodes = 0;
1044 cr->duty = (DUTY_PP | DUTY_PME);
1045 npme_major = 1;
1046 npme_minor = 1;
1048 if (inputrec->ePBC == epbcSCREW)
1050 gmx_fatal(FARGS,
1051 "pbc=%s is only implemented with domain decomposition",
1052 epbc_names[inputrec->ePBC]);
1056 if (PAR(cr))
1058 /* After possible communicator splitting in make_dd_communicators.
1059 * we can set up the intra/inter node communication.
1061 gmx_setup_nodecomm(fplog, cr);
1064 /* Initialize per-physical-node MPI process/thread ID and counters. */
1065 gmx_init_intranode_counters(cr);
1066 #ifdef GMX_MPI
1067 if (MULTISIM(cr))
1069 md_print_info(cr, fplog,
1070 "This is simulation %d out of %d running as a composite GROMACS\n"
1071 "multi-simulation job. Setup for this simulation:\n\n",
1072 cr->ms->sim, cr->ms->nsim);
1074 md_print_info(cr, fplog, "Using %d MPI %s\n",
1075 cr->nnodes,
1076 #ifdef GMX_THREAD_MPI
1077 cr->nnodes == 1 ? "thread" : "threads"
1078 #else
1079 cr->nnodes == 1 ? "process" : "processes"
1080 #endif
1082 fflush(stderr);
1083 #endif
1085 /* Check and update hw_opt for the cut-off scheme */
1086 check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1088 /* Check and update hw_opt for the number of MPI ranks */
1089 check_and_update_hw_opt_3(hw_opt);
1091 gmx_omp_nthreads_init(fplog, cr,
1092 hwinfo->nthreads_hw_avail,
1093 hw_opt->nthreads_omp,
1094 hw_opt->nthreads_omp_pme,
1095 (cr->duty & DUTY_PP) == 0,
1096 inputrec->cutoff_scheme == ecutsVERLET);
1098 #ifndef NDEBUG
1099 if (EI_TPI(inputrec->eI) &&
1100 inputrec->cutoff_scheme == ecutsVERLET)
1102 gmx_feenableexcept();
1104 #endif
1106 if (bUseGPU)
1108 /* Select GPU id's to use */
1109 gmx_select_gpu_ids(fplog, cr, &hwinfo->gpu_info, bForceUseGPU,
1110 &hw_opt->gpu_opt);
1112 else
1114 /* Ignore (potentially) manually selected GPUs */
1115 hw_opt->gpu_opt.n_dev_use = 0;
1118 /* check consistency across ranks of things like SIMD
1119 * support and number of GPUs selected */
1120 gmx_check_hw_runconf_consistency(fplog, hwinfo, cr, hw_opt, bUseGPU);
1122 /* Now that we know the setup is consistent, check for efficiency */
1123 check_resource_division_efficiency(hwinfo, hw_opt, Flags & MD_NTOMPSET,
1124 cr, fplog);
1126 if (DOMAINDECOMP(cr))
1128 /* When we share GPUs over ranks, we need to know this for the DLB */
1129 dd_setup_dlb_resource_sharing(cr, hwinfo, hw_opt);
1132 /* getting number of PP/PME threads
1133 PME: env variable should be read only on one node to make sure it is
1134 identical everywhere;
1136 nthreads_pme = gmx_omp_nthreads_get(emntPME);
1138 wcycle = wallcycle_init(fplog, resetstep, cr);
1140 if (PAR(cr))
1142 /* Master synchronizes its value of reset_counters with all nodes
1143 * including PME only nodes */
1144 reset_counters = wcycle_get_reset_counters(wcycle);
1145 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1146 wcycle_set_reset_counters(wcycle, reset_counters);
1149 snew(nrnb, 1);
1150 if (cr->duty & DUTY_PP)
1152 bcast_state(cr, state);
1154 /* Initiate forcerecord */
1155 fr = mk_forcerec();
1156 fr->hwinfo = hwinfo;
1157 fr->gpu_opt = &hw_opt->gpu_opt;
1158 init_forcerec(fplog, fr, fcd, inputrec, mtop, cr, box,
1159 opt2fn("-table", nfile, fnm),
1160 opt2fn("-tablep", nfile, fnm),
1161 opt2fn("-tableb", nfile, fnm),
1162 nbpu_opt,
1163 FALSE,
1164 pforce);
1166 /* Initialize QM-MM */
1167 if (fr->bQMMM)
1169 init_QMMMrec(cr, mtop, inputrec, fr);
1172 /* Initialize the mdatoms structure.
1173 * mdatoms is not filled with atom data,
1174 * as this can not be done now with domain decomposition.
1176 mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1178 /* Initialize the virtual site communication */
1179 vsite = init_vsite(mtop, cr, FALSE);
1181 calc_shifts(box, fr->shift_vec);
1183 /* With periodic molecules the charge groups should be whole at start up
1184 * and the virtual sites should not be far from their proper positions.
1186 if (!inputrec->bContinuation && MASTER(cr) &&
1187 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1189 /* Make molecules whole at start of run */
1190 if (fr->ePBC != epbcNONE)
1192 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, state->x);
1194 if (vsite)
1196 /* Correct initial vsite positions are required
1197 * for the initial distribution in the domain decomposition
1198 * and for the initial shell prediction.
1200 construct_vsites_mtop(vsite, mtop, state->x);
1204 if (EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype))
1206 ewaldcoeff_q = fr->ewaldcoeff_q;
1207 ewaldcoeff_lj = fr->ewaldcoeff_lj;
1208 pmedata = &fr->pmedata;
1210 else
1212 pmedata = NULL;
1215 else
1217 /* This is a PME only node */
1219 /* We don't need the state */
1220 done_state(state);
1222 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1223 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1224 snew(pmedata, 1);
1227 if (hw_opt->thread_affinity != threadaffOFF)
1229 /* Before setting affinity, check whether the affinity has changed
1230 * - which indicates that probably the OpenMP library has changed it
1231 * since we first checked).
1233 gmx_check_thread_affinity_set(fplog, cr,
1234 hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1236 /* Set the CPU affinity */
1237 gmx_set_thread_affinity(fplog, cr, hw_opt, hwinfo);
1240 /* Initiate PME if necessary,
1241 * either on all nodes or on dedicated PME nodes only. */
1242 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1244 if (mdatoms)
1246 nChargePerturbed = mdatoms->nChargePerturbed;
1247 if (EVDW_PME(inputrec->vdwtype))
1249 nTypePerturbed = mdatoms->nTypePerturbed;
1252 if (cr->npmenodes > 0)
1254 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1255 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1256 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1259 if (cr->duty & DUTY_PME)
1261 status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1262 mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1263 (Flags & MD_REPRODUCIBLE), nthreads_pme);
1264 if (status != 0)
1266 gmx_fatal(FARGS, "Error %d initializing PME", status);
1272 if (EI_DYNAMICS(inputrec->eI))
1274 /* Turn on signal handling on all nodes */
1276 * (A user signal from the PME nodes (if any)
1277 * is communicated to the PP nodes.
1279 signal_handler_install();
1282 if (cr->duty & DUTY_PP)
1284 /* Assumes uniform use of the number of OpenMP threads */
1285 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1287 if (inputrec->bPull)
1289 /* Initialize pull code */
1290 inputrec->pull_work =
1291 init_pull(fplog, inputrec->pull, inputrec, nfile, fnm,
1292 mtop, cr, oenv, inputrec->fepvals->init_lambda,
1293 EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1296 if (inputrec->bRot)
1298 /* Initialize enforced rotation code */
1299 init_rot(fplog, inputrec, nfile, fnm, cr, state->x, box, mtop, oenv,
1300 bVerbose, Flags);
1303 if (inputrec->eSwapCoords != eswapNO)
1305 /* Initialize ion swapping code */
1306 init_swapcoords(fplog, bVerbose, inputrec, opt2fn_master("-swap", nfile, fnm, cr),
1307 mtop, state->x, state->box, &state->swapstate, cr, oenv, Flags);
1310 constr = init_constraints(fplog, mtop, inputrec, ed, state, cr);
1312 if (DOMAINDECOMP(cr))
1314 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1315 dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1316 Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1318 set_dd_parameters(fplog, cr->dd, dlb_scale, inputrec, &ddbox);
1320 setup_dd_grid(fplog, cr->dd);
1323 /* Now do whatever the user wants us to do (how flexible...) */
1324 my_integrator(inputrec->eI) (fplog, cr, nfile, fnm,
1325 oenv, bVerbose, bCompact,
1326 nstglobalcomm,
1327 vsite, constr,
1328 nstepout, inputrec, mtop,
1329 fcd, state,
1330 mdatoms, nrnb, wcycle, ed, fr,
1331 repl_ex_nst, repl_ex_nex, repl_ex_seed,
1332 membed,
1333 cpt_period, max_hours,
1334 imdport,
1335 Flags,
1336 walltime_accounting);
1338 if (inputrec->bPull)
1340 finish_pull(inputrec->pull_work);
1343 if (inputrec->bRot)
1345 finish_rot(inputrec->rot);
1349 else
1351 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1352 /* do PME only */
1353 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1354 gmx_pmeonly(*pmedata, cr, nrnb, wcycle, walltime_accounting, ewaldcoeff_q, ewaldcoeff_lj, inputrec);
1357 wallcycle_stop(wcycle, ewcRUN);
1359 /* Finish up, write some stuff
1360 * if rerunMD, don't write last frame again
1362 finish_run(fplog, cr,
1363 inputrec, nrnb, wcycle, walltime_accounting,
1364 fr ? fr->nbv : NULL,
1365 EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1368 /* Free GPU memory and context */
1369 free_gpu_resources(fr, cr, &hwinfo->gpu_info, fr ? fr->gpu_opt : NULL);
1371 if (membed != nullptr)
1373 free_membed(membed);
1376 gmx_hardware_info_free(hwinfo);
1378 /* Does what it says */
1379 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1380 walltime_accounting_destroy(walltime_accounting);
1382 /* Close logfile already here if we were appending to it */
1383 if (MASTER(cr) && (Flags & MD_APPENDFILES))
1385 gmx_log_close(fplog);
1388 rc = (int)gmx_get_stop_condition();
1390 done_ed(&ed);
1392 #ifdef GMX_THREAD_MPI
1393 /* we need to join all threads. The sub-threads join when they
1394 exit this function, but the master thread needs to be told to
1395 wait for that. */
1396 if (PAR(cr) && MASTER(cr))
1398 tMPI_Finalize();
1400 #endif
1402 return rc;
1405 } // namespace gmx