Merge "Merge branch release-5-1 into master"
[gromacs/AngularHB.git] / src / programs / mdrun / runner.cpp
blob9e395d347fc949d08c64f8dcaaa902430a0ffdf7
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/essentialdynamics/edsam.h"
59 #include "gromacs/ewald/pme.h"
60 #include "gromacs/fileio/oenv.h"
61 #include "gromacs/fileio/tpxio.h"
62 #include "gromacs/fileio/trx.h"
63 #include "gromacs/gmxlib/disre.h"
64 #include "gromacs/gmxlib/main.h"
65 #include "gromacs/gmxlib/md_logging.h"
66 #include "gromacs/gmxlib/orires.h"
67 #include "gromacs/gmxlib/sighandler.h"
68 #include "gromacs/gmxlib/gpu_utils/gpu_utils.h"
69 #include "gromacs/legacyheaders/checkpoint.h"
70 #include "gromacs/legacyheaders/copyrite.h"
71 #include "gromacs/legacyheaders/force.h"
72 #include "gromacs/legacyheaders/gmx_detect_hardware.h"
73 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
74 #include "gromacs/legacyheaders/gmx_thread_affinity.h"
75 #include "gromacs/legacyheaders/inputrec.h"
76 #include "gromacs/legacyheaders/names.h"
77 #include "gromacs/legacyheaders/network.h"
78 #include "gromacs/legacyheaders/txtdump.h"
79 #include "gromacs/legacyheaders/typedefs.h"
80 #include "gromacs/math/calculate-ewald-splitting-coefficient.h"
81 #include "gromacs/math/vec.h"
82 #include "gromacs/mdlib/calc_verletbuf.h"
83 #include "gromacs/mdlib/constr.h"
84 #include "gromacs/mdlib/forcerec.h"
85 #include "gromacs/mdlib/integrator.h"
86 #include "gromacs/mdlib/md_support.h"
87 #include "gromacs/mdlib/mdatoms.h"
88 #include "gromacs/mdlib/mdrun.h"
89 #include "gromacs/mdlib/minimize.h"
90 #include "gromacs/mdlib/nbnxn_search.h"
91 #include "gromacs/mdlib/qmmm.h"
92 #include "gromacs/mdlib/tpi.h"
93 #include "gromacs/pbcutil/pbc.h"
94 #include "gromacs/pulling/pull.h"
95 #include "gromacs/pulling/pull_rotation.h"
96 #include "gromacs/swap/swapcoords.h"
97 #include "gromacs/timing/wallcycle.h"
98 #include "gromacs/topology/mtop_util.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/smalloc.h"
106 #include "deform.h"
107 #include "md.h"
108 #include "membed.h"
109 #include "repl_ex.h"
110 #include "resource-division.h"
112 #ifdef GMX_FAHCORE
113 #include "corewrap.h"
114 #endif
116 //! First step used in pressure scaling
117 gmx_int64_t deform_init_init_step_tpx;
118 //! Initial box for pressure scaling
119 matrix deform_init_box_tpx;
120 //! MPI variable for use in pressure scaling
121 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
123 #ifdef GMX_THREAD_MPI
124 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
125 * the number of threads will get lowered.
127 #define MIN_ATOMS_PER_MPI_THREAD 90
128 #define MIN_ATOMS_PER_GPU 900
130 struct mdrunner_arglist
132 gmx_hw_opt_t hw_opt;
133 FILE *fplog;
134 t_commrec *cr;
135 int nfile;
136 const t_filenm *fnm;
137 const gmx_output_env_t *oenv;
138 gmx_bool bVerbose;
139 gmx_bool bCompact;
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)
172 struct mdrunner_arglist *mda = (struct mdrunner_arglist*)arg;
173 struct mdrunner_arglist mc = *mda; /* copy the arg list to make sure
174 that it's thread-local. This doesn't
175 copy pointed-to items, of course,
176 but those are all const. */
177 t_commrec *cr; /* we need a local version of this */
178 FILE *fplog = NULL;
179 t_filenm *fnm;
181 fnm = dup_tfn(mc.nfile, mc.fnm);
183 cr = reinitialize_commrec_for_this_thread(mc.cr);
185 if (MASTER(cr))
187 fplog = mc.fplog;
190 gmx::mdrunner(&mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
191 mc.bVerbose, mc.bCompact, mc.nstglobalcomm,
192 mc.ddxyz, mc.dd_node_order, mc.rdd,
193 mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
194 mc.ddcsx, mc.ddcsy, mc.ddcsz,
195 mc.nbpu_opt, mc.nstlist_cmdline,
196 mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
197 mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce,
198 mc.cpt_period, mc.max_hours, mc.imdport, mc.Flags);
201 /* called by mdrunner() to start a specific number of threads (including
202 the main thread) for thread-parallel runs. This in turn calls mdrunner()
203 for each thread.
204 All options besides nthreads are the same as for mdrunner(). */
205 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
206 FILE *fplog, t_commrec *cr, int nfile,
207 const t_filenm fnm[], const gmx_output_env_t *oenv, gmx_bool bVerbose,
208 gmx_bool bCompact, int nstglobalcomm,
209 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
210 const char *dddlb_opt, real dlb_scale,
211 const char *ddcsx, const char *ddcsy, const char *ddcsz,
212 const char *nbpu_opt, int nstlist_cmdline,
213 gmx_int64_t nsteps_cmdline,
214 int nstepout, int resetstep,
215 int nmultisim, int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
216 real pforce, real cpt_period, real max_hours,
217 unsigned long Flags)
219 int ret;
220 struct mdrunner_arglist *mda;
221 t_commrec *crn; /* the new commrec */
222 t_filenm *fnmn;
224 /* first check whether we even need to start tMPI */
225 if (hw_opt->nthreads_tmpi < 2)
227 return cr;
230 /* a few small, one-time, almost unavoidable memory leaks: */
231 snew(mda, 1);
232 fnmn = dup_tfn(nfile, fnm);
234 /* fill the data structure to pass as void pointer to thread start fn */
235 /* hw_opt contains pointers, which should all be NULL at this stage */
236 mda->hw_opt = *hw_opt;
237 mda->fplog = fplog;
238 mda->cr = cr;
239 mda->nfile = nfile;
240 mda->fnm = fnmn;
241 mda->oenv = oenv;
242 mda->bVerbose = bVerbose;
243 mda->bCompact = bCompact;
244 mda->nstglobalcomm = nstglobalcomm;
245 mda->ddxyz[XX] = ddxyz[XX];
246 mda->ddxyz[YY] = ddxyz[YY];
247 mda->ddxyz[ZZ] = ddxyz[ZZ];
248 mda->dd_node_order = dd_node_order;
249 mda->rdd = rdd;
250 mda->rconstr = rconstr;
251 mda->dddlb_opt = dddlb_opt;
252 mda->dlb_scale = dlb_scale;
253 mda->ddcsx = ddcsx;
254 mda->ddcsy = ddcsy;
255 mda->ddcsz = ddcsz;
256 mda->nbpu_opt = nbpu_opt;
257 mda->nstlist_cmdline = nstlist_cmdline;
258 mda->nsteps_cmdline = nsteps_cmdline;
259 mda->nstepout = nstepout;
260 mda->resetstep = resetstep;
261 mda->nmultisim = nmultisim;
262 mda->repl_ex_nst = repl_ex_nst;
263 mda->repl_ex_nex = repl_ex_nex;
264 mda->repl_ex_seed = repl_ex_seed;
265 mda->pforce = pforce;
266 mda->cpt_period = cpt_period;
267 mda->max_hours = max_hours;
268 mda->Flags = Flags;
270 /* now spawn new threads that start mdrunner_start_fn(), while
271 the main thread returns, we set thread affinity later */
272 ret = tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi, TMPI_AFFINITY_NONE,
273 mdrunner_start_fn, (void*)(mda) );
274 if (ret != TMPI_SUCCESS)
276 return NULL;
279 crn = reinitialize_commrec_for_this_thread(cr);
280 return crn;
283 #endif /* GMX_THREAD_MPI */
286 /*! \brief Cost of non-bonded kernels
288 * We determine the extra cost of the non-bonded kernels compared to
289 * a reference nstlist value of 10 (which is the default in grompp).
291 static const int nbnxnReferenceNstlist = 10;
292 //! The values to try when switching
293 const int nstlist_try[] = { 20, 25, 40 };
294 //! Number of elements in the neighborsearch list trials.
295 #define NNSTL sizeof(nstlist_try)/sizeof(nstlist_try[0])
296 /* Increase nstlist until the non-bonded cost increases more than listfac_ok,
297 * but never more than listfac_max.
298 * A standard (protein+)water system at 300K with PME ewald_rtol=1e-5
299 * needs 1.28 at rcoulomb=0.9 and 1.24 at rcoulomb=1.0 to get to nstlist=40.
300 * Note that both CPU and GPU factors are conservative. Performance should
301 * not go down due to this tuning, except with a relatively slow GPU.
302 * On the other hand, at medium/high parallelization or with fast GPUs
303 * nstlist will not be increased enough to reach optimal performance.
305 /* CPU: pair-search is about a factor 1.5 slower than the non-bonded kernel */
306 //! Max OK performance ratio beween force calc and neighbor searching
307 static const float nbnxn_cpu_listfac_ok = 1.05;
308 //! Too high performance ratio beween force calc and neighbor searching
309 static const float nbnxn_cpu_listfac_max = 1.09;
310 /* GPU: pair-search is a factor 1.5-3 slower than the non-bonded kernel */
311 //! Max OK performance ratio beween force calc and neighbor searching
312 static const float nbnxn_gpu_listfac_ok = 1.20;
313 //! Too high performance ratio beween force calc and neighbor searching
314 static const float nbnxn_gpu_listfac_max = 1.30;
316 /*! \brief Try to increase nstlist when using the Verlet cut-off scheme */
317 static void increase_nstlist(FILE *fp, t_commrec *cr,
318 t_inputrec *ir, int nstlist_cmdline,
319 const gmx_mtop_t *mtop, matrix box,
320 gmx_bool bGPU)
322 float listfac_ok, listfac_max;
323 int nstlist_orig, nstlist_prev;
324 verletbuf_list_setup_t ls;
325 real rlistWithReferenceNstlist, rlist_inc, rlist_ok, rlist_max;
326 real rlist_new, rlist_prev;
327 size_t nstlist_ind = 0;
328 t_state state_tmp;
329 gmx_bool bBox, bDD, bCont;
330 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";
331 const char *nve_err = "Can not increase nstlist because an NVE ensemble is used";
332 const char *vbd_err = "Can not increase nstlist because verlet-buffer-tolerance is not set or used";
333 const char *box_err = "Can not increase nstlist because the box is too small";
334 const char *dd_err = "Can not increase nstlist because of domain decomposition limitations";
335 char buf[STRLEN];
336 const float oneThird = 1.0f / 3.0f;
338 if (nstlist_cmdline <= 0)
340 if (ir->nstlist == 1)
342 /* The user probably set nstlist=1 for a reason,
343 * don't mess with the settings.
345 return;
348 if (fp != NULL && bGPU && ir->nstlist < nstlist_try[0])
350 fprintf(fp, nstl_gpu, ir->nstlist);
352 nstlist_ind = 0;
353 while (nstlist_ind < NNSTL && ir->nstlist >= nstlist_try[nstlist_ind])
355 nstlist_ind++;
357 if (nstlist_ind == NNSTL)
359 /* There are no larger nstlist value to try */
360 return;
364 if (EI_MD(ir->eI) && ir->etc == etcNO)
366 if (MASTER(cr))
368 fprintf(stderr, "%s\n", nve_err);
370 if (fp != NULL)
372 fprintf(fp, "%s\n", nve_err);
375 return;
378 if (ir->verletbuf_tol == 0 && bGPU)
380 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");
383 if (ir->verletbuf_tol < 0)
385 if (MASTER(cr))
387 fprintf(stderr, "%s\n", vbd_err);
389 if (fp != NULL)
391 fprintf(fp, "%s\n", vbd_err);
394 return;
397 if (bGPU)
399 listfac_ok = nbnxn_gpu_listfac_ok;
400 listfac_max = nbnxn_gpu_listfac_max;
402 else
404 listfac_ok = nbnxn_cpu_listfac_ok;
405 listfac_max = nbnxn_cpu_listfac_max;
408 nstlist_orig = ir->nstlist;
409 if (nstlist_cmdline > 0)
411 if (fp)
413 sprintf(buf, "Getting nstlist=%d from command line option",
414 nstlist_cmdline);
416 ir->nstlist = nstlist_cmdline;
419 verletbuf_get_list_setup(TRUE, bGPU, &ls);
421 /* Allow rlist to make the list a given factor larger than the list
422 * would be with the reference value for nstlist (10).
424 nstlist_prev = ir->nstlist;
425 ir->nstlist = nbnxnReferenceNstlist;
426 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL,
427 &rlistWithReferenceNstlist);
428 ir->nstlist = nstlist_prev;
430 /* Determine the pair list size increase due to zero interactions */
431 rlist_inc = nbnxn_get_rlist_effective_inc(ls.cluster_size_j,
432 mtop->natoms/det(box));
433 rlist_ok = (rlistWithReferenceNstlist + rlist_inc)*pow(listfac_ok, oneThird) - rlist_inc;
434 rlist_max = (rlistWithReferenceNstlist + rlist_inc)*pow(listfac_max, oneThird) - rlist_inc;
435 if (debug)
437 fprintf(debug, "nstlist tuning: rlist_inc %.3f rlist_ok %.3f rlist_max %.3f\n",
438 rlist_inc, rlist_ok, rlist_max);
441 nstlist_prev = nstlist_orig;
442 rlist_prev = ir->rlist;
445 if (nstlist_cmdline <= 0)
447 ir->nstlist = nstlist_try[nstlist_ind];
450 /* Set the pair-list buffer size in ir */
451 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
453 /* Does rlist fit in the box? */
454 bBox = (sqr(rlist_new) < max_cutoff2(ir->ePBC, box));
455 bDD = TRUE;
456 if (bBox && DOMAINDECOMP(cr))
458 /* Check if rlist fits in the domain decomposition */
459 if (inputrec2nboundeddim(ir) < DIM)
461 gmx_incons("Changing nstlist with domain decomposition and unbounded dimensions is not implemented yet");
463 copy_mat(box, state_tmp.box);
464 bDD = change_dd_cutoff(cr, &state_tmp, ir, rlist_new);
467 if (debug)
469 fprintf(debug, "nstlist %d rlist %.3f bBox %d bDD %d\n",
470 ir->nstlist, rlist_new, bBox, bDD);
473 bCont = FALSE;
475 if (nstlist_cmdline <= 0)
477 if (bBox && bDD && rlist_new <= rlist_max)
479 /* Increase nstlist */
480 nstlist_prev = ir->nstlist;
481 rlist_prev = rlist_new;
482 bCont = (nstlist_ind+1 < NNSTL && rlist_new < rlist_ok);
484 else
486 /* Stick with the previous nstlist */
487 ir->nstlist = nstlist_prev;
488 rlist_new = rlist_prev;
489 bBox = TRUE;
490 bDD = TRUE;
494 nstlist_ind++;
496 while (bCont);
498 if (!bBox || !bDD)
500 gmx_warning(!bBox ? box_err : dd_err);
501 if (fp != NULL)
503 fprintf(fp, "\n%s\n", bBox ? box_err : dd_err);
505 ir->nstlist = nstlist_orig;
507 else if (ir->nstlist != nstlist_orig || rlist_new != ir->rlist)
509 sprintf(buf, "Changing nstlist from %d to %d, rlist from %g to %g",
510 nstlist_orig, ir->nstlist,
511 ir->rlist, rlist_new);
512 if (MASTER(cr))
514 fprintf(stderr, "%s\n\n", buf);
516 if (fp != NULL)
518 fprintf(fp, "%s\n\n", buf);
520 ir->rlist = rlist_new;
521 ir->rlistlong = rlist_new;
525 /*! \brief Initialize variables for Verlet scheme simulation */
526 static void prepare_verlet_scheme(FILE *fplog,
527 t_commrec *cr,
528 t_inputrec *ir,
529 int nstlist_cmdline,
530 const gmx_mtop_t *mtop,
531 matrix box,
532 gmx_bool bUseGPU)
534 /* For NVE simulations, we will retain the initial list buffer */
535 if (ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
537 /* Update the Verlet buffer size for the current run setup */
538 verletbuf_list_setup_t ls;
539 real rlist_new;
541 /* Here we assume SIMD-enabled kernels are being used. But as currently
542 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
543 * and 4x2 gives a larger buffer than 4x4, this is ok.
545 verletbuf_get_list_setup(TRUE, bUseGPU, &ls);
547 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
549 if (rlist_new != ir->rlist)
551 if (fplog != NULL)
553 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
554 ir->rlist, rlist_new,
555 ls.cluster_size_i, ls.cluster_size_j);
557 ir->rlist = rlist_new;
558 ir->rlistlong = rlist_new;
562 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
564 gmx_fatal(FARGS, "Can not set nstlist without %s",
565 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
568 if (EI_DYNAMICS(ir->eI))
570 /* Set or try nstlist values */
571 increase_nstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, bUseGPU);
575 /*! \brief Override the nslist value in inputrec
577 * with value passed on the command line (if any)
579 static void override_nsteps_cmdline(FILE *fplog,
580 gmx_int64_t nsteps_cmdline,
581 t_inputrec *ir,
582 const t_commrec *cr)
584 assert(ir);
585 assert(cr);
587 /* override with anything else than the default -2 */
588 if (nsteps_cmdline > -2)
590 char sbuf_steps[STEPSTRSIZE];
591 char sbuf_msg[STRLEN];
593 ir->nsteps = nsteps_cmdline;
594 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
596 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
597 gmx_step_str(nsteps_cmdline, sbuf_steps),
598 fabs(nsteps_cmdline*ir->delta_t));
600 else
602 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
603 gmx_step_str(nsteps_cmdline, sbuf_steps));
606 md_print_warn(cr, fplog, "%s\n", sbuf_msg);
608 else if (nsteps_cmdline < -2)
610 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %d",
611 nsteps_cmdline);
613 /* Do nothing if nsteps_cmdline == -2 */
616 namespace gmx
619 //! \brief Return the correct integrator function.
620 static integrator_t *my_integrator(unsigned int ei)
622 switch (ei)
624 case eiMD:
625 case eiBD:
626 case eiSD2:
627 case eiSD1:
628 case eiVV:
629 case eiVVAK:
630 if (!EI_DYNAMICS(ei))
632 GMX_THROW(APIError("do_md integrator would be called for a non-dynamical integrator"));
634 return do_md;
635 case eiSteep:
636 return do_steep;
637 case eiCG:
638 return do_cg;
639 case eiNM:
640 return do_nm;
641 case eiLBFGS:
642 return do_lbfgs;
643 case eiTPI:
644 case eiTPIC:
645 if (!EI_TPI(ei))
647 GMX_THROW(APIError("do_tpi integrator would be called for a non-TPI integrator"));
649 return do_tpi;
650 default:
651 GMX_THROW(APIError("Non existing integrator selected"));
655 int mdrunner(gmx_hw_opt_t *hw_opt,
656 FILE *fplog, t_commrec *cr, int nfile,
657 const t_filenm fnm[], const gmx_output_env_t *oenv, gmx_bool bVerbose,
658 gmx_bool bCompact, int nstglobalcomm,
659 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
660 const char *dddlb_opt, real dlb_scale,
661 const char *ddcsx, const char *ddcsy, const char *ddcsz,
662 const char *nbpu_opt, int nstlist_cmdline,
663 gmx_int64_t nsteps_cmdline, int nstepout, int resetstep,
664 int gmx_unused nmultisim, int repl_ex_nst, int repl_ex_nex,
665 int repl_ex_seed, real pforce, real cpt_period, real max_hours,
666 int imdport, unsigned long Flags)
668 gmx_bool bForceUseGPU, bTryUseGPU, bRerunMD;
669 t_inputrec *inputrec;
670 t_state *state = NULL;
671 matrix box;
672 gmx_ddbox_t ddbox = {0};
673 int npme_major, npme_minor;
674 t_nrnb *nrnb;
675 gmx_mtop_t *mtop = NULL;
676 t_mdatoms *mdatoms = NULL;
677 t_forcerec *fr = NULL;
678 t_fcdata *fcd = NULL;
679 real ewaldcoeff_q = 0;
680 real ewaldcoeff_lj = 0;
681 struct gmx_pme_t **pmedata = NULL;
682 gmx_vsite_t *vsite = NULL;
683 gmx_constr_t constr;
684 int nChargePerturbed = -1, nTypePerturbed = 0, status;
685 gmx_wallcycle_t wcycle;
686 gmx_bool bReadEkin;
687 gmx_walltime_accounting_t walltime_accounting = NULL;
688 int rc;
689 gmx_int64_t reset_counters;
690 gmx_edsam_t ed = NULL;
691 int nthreads_pme = 1;
692 int nthreads_pp = 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, "Pall2015");
722 please_cite(fplog, "Pronk2013");
723 please_cite(fplog, "Hess2008b");
724 please_cite(fplog, "Spoel2005a");
725 please_cite(fplog, "Lindahl2001a");
726 please_cite(fplog, "Berendsen95a");
729 snew(state, 1);
730 if (SIMMASTER(cr))
732 /* Read (nearly) all data required for the simulation */
733 read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, state, mtop);
735 if (inputrec->cutoff_scheme == ecutsVERLET)
737 /* Here the master rank decides if all ranks will use GPUs */
738 bUseGPU = (hwinfo->gpu_info.n_dev_compatible > 0 ||
739 getenv("GMX_EMULATE_GPU") != NULL);
741 /* TODO add GPU kernels for this and replace this check by:
742 * (bUseGPU && (ir->vdwtype == evdwPME &&
743 * ir->ljpme_combination_rule == eljpmeLB))
744 * update the message text and the content of nbnxn_acceleration_supported.
746 if (bUseGPU &&
747 !nbnxn_gpu_acceleration_supported(fplog, cr, inputrec, bRerunMD))
749 /* Fallback message printed by nbnxn_acceleration_supported */
750 if (bForceUseGPU)
752 gmx_fatal(FARGS, "GPU acceleration requested, but not supported with the given input settings");
754 bUseGPU = FALSE;
757 prepare_verlet_scheme(fplog, cr,
758 inputrec, nstlist_cmdline, mtop, state->box,
759 bUseGPU);
761 else
763 if (nstlist_cmdline > 0)
765 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
768 if (hwinfo->gpu_info.n_dev_compatible > 0)
770 md_print_warn(cr, fplog,
771 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
772 " To use a GPU, set the mdp option: cutoff-scheme = Verlet\n");
775 if (bForceUseGPU)
777 gmx_fatal(FARGS, "GPU requested, but can't be used without cutoff-scheme=Verlet");
780 #ifdef GMX_TARGET_BGQ
781 md_print_warn(cr, fplog,
782 "NOTE: There is no SIMD implementation of the group scheme kernels on\n"
783 " BlueGene/Q. You will observe better performance from using the\n"
784 " Verlet cut-off scheme.\n");
785 #endif
788 if (inputrec->eI == eiSD2)
790 md_print_warn(cr, fplog, "The stochastic dynamics integrator %s is deprecated, since\n"
791 "it is slower than integrator %s and is slightly less accurate\n"
792 "with constraints. Use the %s integrator.",
793 ei_names[inputrec->eI], ei_names[eiSD1], ei_names[eiSD1]);
797 /* Check and update the hardware options for internal consistency */
798 check_and_update_hw_opt_1(hw_opt, cr);
800 /* Early check for externally set process affinity. */
801 gmx_check_thread_affinity_set(fplog, cr,
802 hw_opt, hwinfo->nthreads_hw_avail, FALSE);
804 #ifdef GMX_THREAD_MPI
805 if (SIMMASTER(cr))
807 if (cr->npmenodes > 0 && hw_opt->nthreads_tmpi <= 0)
809 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
812 /* Since the master knows the cut-off scheme, update hw_opt for this.
813 * This is done later for normal MPI and also once more with tMPI
814 * for all tMPI ranks.
816 check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
818 /* NOW the threads will be started: */
819 hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
820 hw_opt,
821 inputrec, mtop,
822 cr, fplog, bUseGPU);
824 if (hw_opt->nthreads_tmpi > 1)
826 t_commrec *cr_old = cr;
827 /* now start the threads. */
828 cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
829 oenv, bVerbose, bCompact, nstglobalcomm,
830 ddxyz, dd_node_order, rdd, rconstr,
831 dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
832 nbpu_opt, nstlist_cmdline,
833 nsteps_cmdline, nstepout, resetstep, nmultisim,
834 repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
835 cpt_period, max_hours,
836 Flags);
837 /* the main thread continues here with a new cr. We don't deallocate
838 the old cr because other threads may still be reading it. */
839 if (cr == NULL)
841 gmx_comm("Failed to spawn threads");
845 #endif
846 /* END OF CAUTION: cr is now reliable */
848 /* g_membed initialisation *
849 * Because we change the mtop, init_membed is called before the init_parallel *
850 * (in case we ever want to make it run in parallel) */
851 if (opt2bSet("-membed", nfile, fnm))
853 if (MASTER(cr))
855 fprintf(stderr, "Initializing membed");
857 membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
860 if (PAR(cr))
862 /* now broadcast everything to the non-master nodes/threads: */
863 init_parallel(cr, inputrec, mtop);
865 /* The master rank decided on the use of GPUs,
866 * broadcast this information to all ranks.
868 gmx_bcast_sim(sizeof(bUseGPU), &bUseGPU, cr);
871 if (fplog != NULL)
873 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
874 fprintf(fplog, "\n");
877 /* now make sure the state is initialized and propagated */
878 set_state_entries(state, inputrec);
880 /* A parallel command line option consistency check that we can
881 only do after any threads have started. */
882 if (!PAR(cr) &&
883 (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
885 gmx_fatal(FARGS,
886 "The -dd or -npme option request a parallel simulation, "
887 #ifndef GMX_MPI
888 "but %s was compiled without threads or MPI enabled"
889 #else
890 #ifdef GMX_THREAD_MPI
891 "but the number of threads (option -nt) is 1"
892 #else
893 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
894 #endif
895 #endif
896 , output_env_get_program_display_name(oenv)
900 if (bRerunMD &&
901 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
903 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
906 if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
908 gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
911 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
913 if (cr->npmenodes > 0)
915 gmx_fatal_collective(FARGS, cr, NULL,
916 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
919 cr->npmenodes = 0;
922 if (bUseGPU && cr->npmenodes < 0)
924 /* With GPUs we don't automatically use PME-only ranks. PME ranks can
925 * improve performance with many threads per GPU, since our OpenMP
926 * scaling is bad, but it's difficult to automate the setup.
928 cr->npmenodes = 0;
931 #ifdef GMX_FAHCORE
932 if (MASTER(cr))
934 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
936 #endif
938 /* NMR restraints must be initialized before load_checkpoint,
939 * since with time averaging the history is added to t_state.
940 * For proper consistency check we therefore need to extend
941 * t_state here.
942 * So the PME-only nodes (if present) will also initialize
943 * the distance restraints.
945 snew(fcd, 1);
947 /* This needs to be called before read_checkpoint to extend the state */
948 init_disres(fplog, mtop, inputrec, cr, fcd, state, repl_ex_nst > 0);
950 init_orires(fplog, mtop, state->x, inputrec, cr, &(fcd->orires),
951 state);
953 if (DEFORM(*inputrec))
955 /* Store the deform reference box before reading the checkpoint */
956 if (SIMMASTER(cr))
958 copy_mat(state->box, box);
960 if (PAR(cr))
962 gmx_bcast(sizeof(box), box, cr);
964 /* Because we do not have the update struct available yet
965 * in which the reference values should be stored,
966 * we store them temporarily in static variables.
967 * This should be thread safe, since they are only written once
968 * and with identical values.
970 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
971 deform_init_init_step_tpx = inputrec->init_step;
972 copy_mat(box, deform_init_box_tpx);
973 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
976 if (opt2bSet("-cpi", nfile, fnm))
978 /* Check if checkpoint file exists before doing continuation.
979 * This way we can use identical input options for the first and subsequent runs...
981 if (gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr) )
983 load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
984 cr, ddxyz,
985 inputrec, state, &bReadEkin,
986 (Flags & MD_APPENDFILES),
987 (Flags & MD_APPENDFILESSET));
989 if (bReadEkin)
991 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 /* TODO nthreads_pp is only used for pinning threads.
1134 * This is a temporary solution until we have a hw topology library.
1136 nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
1137 nthreads_pme = gmx_omp_nthreads_get(emntPME);
1139 wcycle = wallcycle_init(fplog, resetstep, cr, nthreads_pp, nthreads_pme);
1141 if (PAR(cr))
1143 /* Master synchronizes its value of reset_counters with all nodes
1144 * including PME only nodes */
1145 reset_counters = wcycle_get_reset_counters(wcycle);
1146 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1147 wcycle_set_reset_counters(wcycle, reset_counters);
1150 snew(nrnb, 1);
1151 if (cr->duty & DUTY_PP)
1153 bcast_state(cr, state);
1155 /* Initiate forcerecord */
1156 fr = mk_forcerec();
1157 fr->hwinfo = hwinfo;
1158 fr->gpu_opt = &hw_opt->gpu_opt;
1159 init_forcerec(fplog, fr, fcd, inputrec, mtop, cr, box,
1160 opt2fn("-table", nfile, fnm),
1161 opt2fn("-tabletf", nfile, fnm),
1162 opt2fn("-tablep", nfile, fnm),
1163 opt2fn("-tableb", nfile, fnm),
1164 nbpu_opt,
1165 FALSE,
1166 pforce);
1168 /* Initialize QM-MM */
1169 if (fr->bQMMM)
1171 init_QMMMrec(cr, mtop, inputrec, fr);
1174 /* Initialize the mdatoms structure.
1175 * mdatoms is not filled with atom data,
1176 * as this can not be done now with domain decomposition.
1178 mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1180 /* Initialize the virtual site communication */
1181 vsite = init_vsite(mtop, cr, FALSE);
1183 calc_shifts(box, fr->shift_vec);
1185 /* With periodic molecules the charge groups should be whole at start up
1186 * and the virtual sites should not be far from their proper positions.
1188 if (!inputrec->bContinuation && MASTER(cr) &&
1189 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1191 /* Make molecules whole at start of run */
1192 if (fr->ePBC != epbcNONE)
1194 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, state->x);
1196 if (vsite)
1198 /* Correct initial vsite positions are required
1199 * for the initial distribution in the domain decomposition
1200 * and for the initial shell prediction.
1202 construct_vsites_mtop(vsite, mtop, state->x);
1206 if (EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype))
1208 ewaldcoeff_q = fr->ewaldcoeff_q;
1209 ewaldcoeff_lj = fr->ewaldcoeff_lj;
1210 pmedata = &fr->pmedata;
1212 else
1214 pmedata = NULL;
1217 else
1219 /* This is a PME only node */
1221 /* We don't need the state */
1222 done_state(state);
1224 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1225 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1226 snew(pmedata, 1);
1229 if (hw_opt->thread_affinity != threadaffOFF)
1231 /* Before setting affinity, check whether the affinity has changed
1232 * - which indicates that probably the OpenMP library has changed it
1233 * since we first checked).
1235 gmx_check_thread_affinity_set(fplog, cr,
1236 hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1238 /* Set the CPU affinity */
1239 gmx_set_thread_affinity(fplog, cr, hw_opt, hwinfo);
1242 /* Initiate PME if necessary,
1243 * either on all nodes or on dedicated PME nodes only. */
1244 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1246 if (mdatoms)
1248 nChargePerturbed = mdatoms->nChargePerturbed;
1249 if (EVDW_PME(inputrec->vdwtype))
1251 nTypePerturbed = mdatoms->nTypePerturbed;
1254 if (cr->npmenodes > 0)
1256 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1257 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1258 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1261 if (cr->duty & DUTY_PME)
1263 status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1264 mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1265 (Flags & MD_REPRODUCIBLE), nthreads_pme);
1266 if (status != 0)
1268 gmx_fatal(FARGS, "Error %d initializing PME", status);
1274 if (EI_DYNAMICS(inputrec->eI))
1276 /* Turn on signal handling on all nodes */
1278 * (A user signal from the PME nodes (if any)
1279 * is communicated to the PP nodes.
1281 signal_handler_install();
1284 if (cr->duty & DUTY_PP)
1286 /* Assumes uniform use of the number of OpenMP threads */
1287 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1289 if (inputrec->bPull)
1291 /* Initialize pull code */
1292 inputrec->pull_work =
1293 init_pull(fplog, inputrec->pull, inputrec, nfile, fnm,
1294 mtop, cr, oenv, inputrec->fepvals->init_lambda,
1295 EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1298 if (inputrec->bRot)
1300 /* Initialize enforced rotation code */
1301 init_rot(fplog, inputrec, nfile, fnm, cr, state->x, box, mtop, oenv,
1302 bVerbose, Flags);
1305 if (inputrec->eSwapCoords != eswapNO)
1307 /* Initialize ion swapping code */
1308 init_swapcoords(fplog, bVerbose, inputrec, opt2fn_master("-swap", nfile, fnm, cr),
1309 mtop, state->x, state->box, &state->swapstate, cr, oenv, Flags);
1312 constr = init_constraints(fplog, mtop, inputrec, ed, state, cr);
1314 if (DOMAINDECOMP(cr))
1316 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1317 dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1318 Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1320 set_dd_parameters(fplog, cr->dd, dlb_scale, inputrec, &ddbox);
1322 setup_dd_grid(fplog, cr->dd);
1325 /* Now do whatever the user wants us to do (how flexible...) */
1326 my_integrator(inputrec->eI) (fplog, cr, nfile, fnm,
1327 oenv, bVerbose, bCompact,
1328 nstglobalcomm,
1329 vsite, constr,
1330 nstepout, inputrec, mtop,
1331 fcd, state,
1332 mdatoms, nrnb, wcycle, ed, fr,
1333 repl_ex_nst, repl_ex_nex, repl_ex_seed,
1334 membed,
1335 cpt_period, max_hours,
1336 imdport,
1337 Flags,
1338 walltime_accounting);
1340 if (inputrec->bPull)
1342 finish_pull(inputrec->pull_work);
1345 if (inputrec->bRot)
1347 finish_rot(inputrec->rot);
1351 else
1353 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1354 /* do PME only */
1355 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1356 gmx_pmeonly(*pmedata, cr, nrnb, wcycle, walltime_accounting, ewaldcoeff_q, ewaldcoeff_lj, inputrec);
1359 wallcycle_stop(wcycle, ewcRUN);
1361 /* Finish up, write some stuff
1362 * if rerunMD, don't write last frame again
1364 finish_run(fplog, cr,
1365 inputrec, nrnb, wcycle, walltime_accounting,
1366 fr ? fr->nbv : NULL,
1367 EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1370 /* Free GPU memory and context */
1371 free_gpu_resources(fr, cr, &hwinfo->gpu_info, fr ? fr->gpu_opt : NULL);
1373 if (membed != nullptr)
1375 free_membed(membed);
1378 gmx_hardware_info_free(hwinfo);
1380 /* Does what it says */
1381 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1382 walltime_accounting_destroy(walltime_accounting);
1384 /* Close logfile already here if we were appending to it */
1385 if (MASTER(cr) && (Flags & MD_APPENDFILES))
1387 gmx_log_close(fplog);
1390 rc = (int)gmx_get_stop_condition();
1392 done_ed(&ed);
1394 #ifdef GMX_THREAD_MPI
1395 /* we need to join all threads. The sub-threads join when they
1396 exit this function, but the master thread needs to be told to
1397 wait for that. */
1398 if (PAR(cr) && MASTER(cr))
1400 tMPI_Finalize();
1402 #endif
1404 return rc;
1407 } // namespace gmx