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.
49 #include "gromacs/domdec/domdec.h"
50 #include "gromacs/essentialdynamics/edsam.h"
51 #include "gromacs/ewald/pme.h"
52 #include "gromacs/fileio/tpxio.h"
53 #include "gromacs/gmxlib/gpu_utils/gpu_utils.h"
54 #include "gromacs/legacyheaders/checkpoint.h"
55 #include "gromacs/legacyheaders/constr.h"
56 #include "gromacs/legacyheaders/copyrite.h"
57 #include "gromacs/legacyheaders/disre.h"
58 #include "gromacs/legacyheaders/force.h"
59 #include "gromacs/legacyheaders/gmx_detect_hardware.h"
60 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
61 #include "gromacs/legacyheaders/gmx_thread_affinity.h"
62 #include "gromacs/legacyheaders/inputrec.h"
63 #include "gromacs/legacyheaders/main.h"
64 #include "gromacs/legacyheaders/md_logging.h"
65 #include "gromacs/legacyheaders/md_support.h"
66 #include "gromacs/legacyheaders/mdatoms.h"
67 #include "gromacs/legacyheaders/mdrun.h"
68 #include "gromacs/legacyheaders/names.h"
69 #include "gromacs/legacyheaders/network.h"
70 #include "gromacs/legacyheaders/oenv.h"
71 #include "gromacs/legacyheaders/orires.h"
72 #include "gromacs/legacyheaders/qmmm.h"
73 #include "gromacs/legacyheaders/sighandler.h"
74 #include "gromacs/legacyheaders/txtdump.h"
75 #include "gromacs/legacyheaders/typedefs.h"
76 #include "gromacs/math/calculate-ewald-splitting-coefficient.h"
77 #include "gromacs/math/vec.h"
78 #include "gromacs/mdlib/calc_verletbuf.h"
79 #include "gromacs/mdlib/nbnxn_search.h"
80 #include "gromacs/pbcutil/pbc.h"
81 #include "gromacs/pulling/pull.h"
82 #include "gromacs/pulling/pull_rotation.h"
83 #include "gromacs/swap/swapcoords.h"
84 #include "gromacs/timing/wallcycle.h"
85 #include "gromacs/topology/mtop_util.h"
86 #include "gromacs/utility/gmxassert.h"
87 #include "gromacs/utility/gmxmpi.h"
88 #include "gromacs/utility/smalloc.h"
93 #include "resource-division.h"
100 gmx_integrator_t *func;
103 /* The array should match the eI array in include/types/enums.h */
104 const gmx_intp_t integrator[eiNR] = { {do_md}, {do_steep}, {do_cg}, {do_md}, {do_md}, {do_nm}, {do_lbfgs}, {do_tpi}, {do_tpi}, {do_md}, {do_md}, {do_md}};
106 gmx_int64_t deform_init_init_step_tpx;
107 matrix deform_init_box_tpx;
108 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
111 #ifdef GMX_THREAD_MPI
112 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
113 * the number of threads will get lowered.
115 #define MIN_ATOMS_PER_MPI_THREAD 90
116 #define MIN_ATOMS_PER_GPU 900
118 struct mdrunner_arglist
133 const char *dddlb_opt;
138 const char *nbpu_opt;
140 gmx_int64_t nsteps_cmdline;
155 /* The function used for spawning threads. Extracts the mdrunner()
156 arguments from its one argument and calls mdrunner(), after making
158 static void mdrunner_start_fn(void *arg)
160 struct mdrunner_arglist *mda = (struct mdrunner_arglist*)arg;
161 struct mdrunner_arglist mc = *mda; /* copy the arg list to make sure
162 that it's thread-local. This doesn't
163 copy pointed-to items, of course,
164 but those are all const. */
165 t_commrec *cr; /* we need a local version of this */
169 fnm = dup_tfn(mc.nfile, mc.fnm);
171 cr = reinitialize_commrec_for_this_thread(mc.cr);
178 mdrunner(&mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
179 mc.bVerbose, mc.bCompact, mc.nstglobalcomm,
180 mc.ddxyz, mc.dd_node_order, mc.rdd,
181 mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
182 mc.ddcsx, mc.ddcsy, mc.ddcsz,
183 mc.nbpu_opt, mc.nstlist_cmdline,
184 mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
185 mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce,
186 mc.cpt_period, mc.max_hours, mc.imdport, mc.Flags);
189 /* called by mdrunner() to start a specific number of threads (including
190 the main thread) for thread-parallel runs. This in turn calls mdrunner()
192 All options besides nthreads are the same as for mdrunner(). */
193 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
194 FILE *fplog, t_commrec *cr, int nfile,
195 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
196 gmx_bool bCompact, int nstglobalcomm,
197 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
198 const char *dddlb_opt, real dlb_scale,
199 const char *ddcsx, const char *ddcsy, const char *ddcsz,
200 const char *nbpu_opt, int nstlist_cmdline,
201 gmx_int64_t nsteps_cmdline,
202 int nstepout, int resetstep,
203 int nmultisim, int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
204 real pforce, real cpt_period, real max_hours,
208 struct mdrunner_arglist *mda;
209 t_commrec *crn; /* the new commrec */
212 /* first check whether we even need to start tMPI */
213 if (hw_opt->nthreads_tmpi < 2)
218 /* a few small, one-time, almost unavoidable memory leaks: */
220 fnmn = dup_tfn(nfile, fnm);
222 /* fill the data structure to pass as void pointer to thread start fn */
223 /* hw_opt contains pointers, which should all be NULL at this stage */
224 mda->hw_opt = *hw_opt;
230 mda->bVerbose = bVerbose;
231 mda->bCompact = bCompact;
232 mda->nstglobalcomm = nstglobalcomm;
233 mda->ddxyz[XX] = ddxyz[XX];
234 mda->ddxyz[YY] = ddxyz[YY];
235 mda->ddxyz[ZZ] = ddxyz[ZZ];
236 mda->dd_node_order = dd_node_order;
238 mda->rconstr = rconstr;
239 mda->dddlb_opt = dddlb_opt;
240 mda->dlb_scale = dlb_scale;
244 mda->nbpu_opt = nbpu_opt;
245 mda->nstlist_cmdline = nstlist_cmdline;
246 mda->nsteps_cmdline = nsteps_cmdline;
247 mda->nstepout = nstepout;
248 mda->resetstep = resetstep;
249 mda->nmultisim = nmultisim;
250 mda->repl_ex_nst = repl_ex_nst;
251 mda->repl_ex_nex = repl_ex_nex;
252 mda->repl_ex_seed = repl_ex_seed;
253 mda->pforce = pforce;
254 mda->cpt_period = cpt_period;
255 mda->max_hours = max_hours;
258 /* now spawn new threads that start mdrunner_start_fn(), while
259 the main thread returns, we set thread affinity later */
260 ret = tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi, TMPI_AFFINITY_NONE,
261 mdrunner_start_fn, (void*)(mda) );
262 if (ret != TMPI_SUCCESS)
267 crn = reinitialize_commrec_for_this_thread(cr);
271 #endif /* GMX_THREAD_MPI */
274 /* We determine the extra cost of the non-bonded kernels compared to
275 * a reference nstlist value of 10 (which is the default in grompp).
277 static const int nbnxnReferenceNstlist = 10;
278 /* The values to try when switching */
279 const int nstlist_try[] = { 20, 25, 40 };
280 #define NNSTL sizeof(nstlist_try)/sizeof(nstlist_try[0])
281 /* Increase nstlist until the non-bonded cost increases more than listfac_ok,
282 * but never more than listfac_max.
283 * A standard (protein+)water system at 300K with PME ewald_rtol=1e-5
284 * needs 1.28 at rcoulomb=0.9 and 1.24 at rcoulomb=1.0 to get to nstlist=40.
285 * Note that both CPU and GPU factors are conservative. Performance should
286 * not go down due to this tuning, except with a relatively slow GPU.
287 * On the other hand, at medium/high parallelization or with fast GPUs
288 * nstlist will not be increased enough to reach optimal performance.
290 /* CPU: pair-search is about a factor 1.5 slower than the non-bonded kernel */
291 static const float nbnxn_cpu_listfac_ok = 1.05;
292 static const float nbnxn_cpu_listfac_max = 1.09;
293 /* GPU: pair-search is a factor 1.5-3 slower than the non-bonded kernel */
294 static const float nbnxn_gpu_listfac_ok = 1.20;
295 static const float nbnxn_gpu_listfac_max = 1.30;
297 /* Try to increase nstlist when using the Verlet cut-off scheme */
298 static void increase_nstlist(FILE *fp, t_commrec *cr,
299 t_inputrec *ir, int nstlist_cmdline,
300 const gmx_mtop_t *mtop, matrix box,
303 float listfac_ok, listfac_max;
304 int nstlist_orig, nstlist_prev;
305 verletbuf_list_setup_t ls;
306 real rlistWithReferenceNstlist, rlist_inc, rlist_ok, rlist_max;
307 real rlist_new, rlist_prev;
308 size_t nstlist_ind = 0;
310 gmx_bool bBox, bDD, bCont;
311 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";
312 const char *nve_err = "Can not increase nstlist because an NVE ensemble is used";
313 const char *vbd_err = "Can not increase nstlist because verlet-buffer-tolerance is not set or used";
314 const char *box_err = "Can not increase nstlist because the box is too small";
315 const char *dd_err = "Can not increase nstlist because of domain decomposition limitations";
317 const float oneThird = 1.0f / 3.0f;
319 if (nstlist_cmdline <= 0)
321 if (ir->nstlist == 1)
323 /* The user probably set nstlist=1 for a reason,
324 * don't mess with the settings.
329 if (fp != NULL && bGPU && ir->nstlist < nstlist_try[0])
331 fprintf(fp, nstl_gpu, ir->nstlist);
334 while (nstlist_ind < NNSTL && ir->nstlist >= nstlist_try[nstlist_ind])
338 if (nstlist_ind == NNSTL)
340 /* There are no larger nstlist value to try */
345 if (EI_MD(ir->eI) && ir->etc == etcNO)
349 fprintf(stderr, "%s\n", nve_err);
353 fprintf(fp, "%s\n", nve_err);
359 if (ir->verletbuf_tol == 0 && bGPU)
361 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");
364 if (ir->verletbuf_tol < 0)
368 fprintf(stderr, "%s\n", vbd_err);
372 fprintf(fp, "%s\n", vbd_err);
380 listfac_ok = nbnxn_gpu_listfac_ok;
381 listfac_max = nbnxn_gpu_listfac_max;
385 listfac_ok = nbnxn_cpu_listfac_ok;
386 listfac_max = nbnxn_cpu_listfac_max;
389 nstlist_orig = ir->nstlist;
390 if (nstlist_cmdline > 0)
394 sprintf(buf, "Getting nstlist=%d from command line option",
397 ir->nstlist = nstlist_cmdline;
400 verletbuf_get_list_setup(TRUE, bGPU, &ls);
402 /* Allow rlist to make the list a given factor larger than the list
403 * would be with the reference value for nstlist (10).
405 nstlist_prev = ir->nstlist;
406 ir->nstlist = nbnxnReferenceNstlist;
407 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL,
408 &rlistWithReferenceNstlist);
409 ir->nstlist = nstlist_prev;
411 /* Determine the pair list size increase due to zero interactions */
412 rlist_inc = nbnxn_get_rlist_effective_inc(ls.cluster_size_j,
413 mtop->natoms/det(box));
414 rlist_ok = (rlistWithReferenceNstlist + rlist_inc)*pow(listfac_ok, oneThird) - rlist_inc;
415 rlist_max = (rlistWithReferenceNstlist + rlist_inc)*pow(listfac_max, oneThird) - rlist_inc;
418 fprintf(debug, "nstlist tuning: rlist_inc %.3f rlist_ok %.3f rlist_max %.3f\n",
419 rlist_inc, rlist_ok, rlist_max);
422 nstlist_prev = nstlist_orig;
423 rlist_prev = ir->rlist;
426 if (nstlist_cmdline <= 0)
428 ir->nstlist = nstlist_try[nstlist_ind];
431 /* Set the pair-list buffer size in ir */
432 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
434 /* Does rlist fit in the box? */
435 bBox = (sqr(rlist_new) < max_cutoff2(ir->ePBC, box));
437 if (bBox && DOMAINDECOMP(cr))
439 /* Check if rlist fits in the domain decomposition */
440 if (inputrec2nboundeddim(ir) < DIM)
442 gmx_incons("Changing nstlist with domain decomposition and unbounded dimensions is not implemented yet");
444 copy_mat(box, state_tmp.box);
445 bDD = change_dd_cutoff(cr, &state_tmp, ir, rlist_new);
450 fprintf(debug, "nstlist %d rlist %.3f bBox %d bDD %d\n",
451 ir->nstlist, rlist_new, bBox, bDD);
456 if (nstlist_cmdline <= 0)
458 if (bBox && bDD && rlist_new <= rlist_max)
460 /* Increase nstlist */
461 nstlist_prev = ir->nstlist;
462 rlist_prev = rlist_new;
463 bCont = (nstlist_ind+1 < NNSTL && rlist_new < rlist_ok);
467 /* Stick with the previous nstlist */
468 ir->nstlist = nstlist_prev;
469 rlist_new = rlist_prev;
481 gmx_warning(!bBox ? box_err : dd_err);
484 fprintf(fp, "\n%s\n", bBox ? box_err : dd_err);
486 ir->nstlist = nstlist_orig;
488 else if (ir->nstlist != nstlist_orig || rlist_new != ir->rlist)
490 sprintf(buf, "Changing nstlist from %d to %d, rlist from %g to %g",
491 nstlist_orig, ir->nstlist,
492 ir->rlist, rlist_new);
495 fprintf(stderr, "%s\n\n", buf);
499 fprintf(fp, "%s\n\n", buf);
501 ir->rlist = rlist_new;
502 ir->rlistlong = rlist_new;
506 static void prepare_verlet_scheme(FILE *fplog,
510 const gmx_mtop_t *mtop,
514 /* For NVE simulations, we will retain the initial list buffer */
515 if (ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
517 /* Update the Verlet buffer size for the current run setup */
518 verletbuf_list_setup_t ls;
521 /* Here we assume SIMD-enabled kernels are being used. But as currently
522 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
523 * and 4x2 gives a larger buffer than 4x4, this is ok.
525 verletbuf_get_list_setup(TRUE, bUseGPU, &ls);
527 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
529 if (rlist_new != ir->rlist)
533 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
534 ir->rlist, rlist_new,
535 ls.cluster_size_i, ls.cluster_size_j);
537 ir->rlist = rlist_new;
538 ir->rlistlong = rlist_new;
542 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
544 gmx_fatal(FARGS, "Can not set nstlist without %s",
545 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
548 if (EI_DYNAMICS(ir->eI))
550 /* Set or try nstlist values */
551 increase_nstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, bUseGPU);
555 /* Override the value in inputrec with value passed on the command line (if any) */
556 static void override_nsteps_cmdline(FILE *fplog,
557 gmx_int64_t nsteps_cmdline,
564 /* override with anything else than the default -2 */
565 if (nsteps_cmdline > -2)
567 char sbuf_steps[STEPSTRSIZE];
568 char sbuf_msg[STRLEN];
570 ir->nsteps = nsteps_cmdline;
571 if (EI_DYNAMICS(ir->eI) && nsteps_cmdline != -1)
573 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps, %.3g ps",
574 gmx_step_str(nsteps_cmdline, sbuf_steps),
575 fabs(nsteps_cmdline*ir->delta_t));
579 sprintf(sbuf_msg, "Overriding nsteps with value passed on the command line: %s steps",
580 gmx_step_str(nsteps_cmdline, sbuf_steps));
583 md_print_warn(cr, fplog, "%s\n", sbuf_msg);
585 else if (nsteps_cmdline < -2)
587 gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %d",
590 /* Do nothing if nsteps_cmdline == -2 */
593 int mdrunner(gmx_hw_opt_t *hw_opt,
594 FILE *fplog, t_commrec *cr, int nfile,
595 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
596 gmx_bool bCompact, int nstglobalcomm,
597 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
598 const char *dddlb_opt, real dlb_scale,
599 const char *ddcsx, const char *ddcsy, const char *ddcsz,
600 const char *nbpu_opt, int nstlist_cmdline,
601 gmx_int64_t nsteps_cmdline, int nstepout, int resetstep,
602 int gmx_unused nmultisim, int repl_ex_nst, int repl_ex_nex,
603 int repl_ex_seed, real pforce, real cpt_period, real max_hours,
604 int imdport, unsigned long Flags)
606 gmx_bool bForceUseGPU, bTryUseGPU, bRerunMD;
607 t_inputrec *inputrec;
608 t_state *state = NULL;
610 gmx_ddbox_t ddbox = {0};
611 int npme_major, npme_minor;
613 gmx_mtop_t *mtop = NULL;
614 t_mdatoms *mdatoms = NULL;
615 t_forcerec *fr = NULL;
616 t_fcdata *fcd = NULL;
617 real ewaldcoeff_q = 0;
618 real ewaldcoeff_lj = 0;
619 struct gmx_pme_t **pmedata = NULL;
620 gmx_vsite_t *vsite = NULL;
622 int nChargePerturbed = -1, nTypePerturbed = 0, status;
623 gmx_wallcycle_t wcycle;
625 gmx_walltime_accounting_t walltime_accounting = NULL;
627 gmx_int64_t reset_counters;
628 gmx_edsam_t ed = NULL;
629 int nthreads_pme = 1;
631 gmx_membed_t membed = NULL;
632 gmx_hw_info_t *hwinfo = NULL;
633 /* The master rank decides early on bUseGPU and broadcasts this later */
634 gmx_bool bUseGPU = FALSE;
636 /* CAUTION: threads may be started later on in this function, so
637 cr doesn't reflect the final parallel state right now */
641 if (Flags & MD_APPENDFILES)
646 bRerunMD = (Flags & MD_RERUN);
647 bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
648 bTryUseGPU = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
650 /* Detect hardware, gather information. This is an operation that is
651 * global for this process (MPI rank). */
652 hwinfo = gmx_detect_hardware(fplog, cr, bTryUseGPU);
654 gmx_print_detected_hardware(fplog, cr, hwinfo);
658 /* Print references after all software/hardware printing */
659 please_cite(fplog, "Pall2015");
660 please_cite(fplog, "Pronk2013");
661 please_cite(fplog, "Hess2008b");
662 please_cite(fplog, "Spoel2005a");
663 please_cite(fplog, "Lindahl2001a");
664 please_cite(fplog, "Berendsen95a");
670 /* Read (nearly) all data required for the simulation */
671 read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, state, NULL, mtop);
673 if (inputrec->cutoff_scheme == ecutsVERLET)
675 /* Here the master rank decides if all ranks will use GPUs */
676 bUseGPU = (hwinfo->gpu_info.n_dev_compatible > 0 ||
677 getenv("GMX_EMULATE_GPU") != NULL);
679 /* TODO add GPU kernels for this and replace this check by:
680 * (bUseGPU && (ir->vdwtype == evdwPME &&
681 * ir->ljpme_combination_rule == eljpmeLB))
682 * update the message text and the content of nbnxn_acceleration_supported.
685 !nbnxn_gpu_acceleration_supported(fplog, cr, inputrec, bRerunMD))
687 /* Fallback message printed by nbnxn_acceleration_supported */
690 gmx_fatal(FARGS, "GPU acceleration requested, but not supported with the given input settings");
695 prepare_verlet_scheme(fplog, cr,
696 inputrec, nstlist_cmdline, mtop, state->box,
701 if (nstlist_cmdline > 0)
703 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
706 if (hwinfo->gpu_info.n_dev_compatible > 0)
708 md_print_warn(cr, fplog,
709 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
710 " To use a GPU, set the mdp option: cutoff-scheme = Verlet\n");
715 gmx_fatal(FARGS, "GPU requested, but can't be used without cutoff-scheme=Verlet");
718 #ifdef GMX_TARGET_BGQ
719 md_print_warn(cr, fplog,
720 "NOTE: There is no SIMD implementation of the group scheme kernels on\n"
721 " BlueGene/Q. You will observe better performance from using the\n"
722 " Verlet cut-off scheme.\n");
726 if (inputrec->eI == eiSD2)
728 md_print_warn(cr, fplog, "The stochastic dynamics integrator %s is deprecated, since\n"
729 "it is slower than integrator %s and is slightly less accurate\n"
730 "with constraints. Use the %s integrator.",
731 ei_names[inputrec->eI], ei_names[eiSD1], ei_names[eiSD1]);
735 /* Check and update the hardware options for internal consistency */
736 check_and_update_hw_opt_1(hw_opt, cr);
738 /* Early check for externally set process affinity. */
739 gmx_check_thread_affinity_set(fplog, cr,
740 hw_opt, hwinfo->nthreads_hw_avail, FALSE);
742 #ifdef GMX_THREAD_MPI
745 if (cr->npmenodes > 0 && hw_opt->nthreads_tmpi <= 0)
747 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
750 /* Since the master knows the cut-off scheme, update hw_opt for this.
751 * This is done later for normal MPI and also once more with tMPI
752 * for all tMPI ranks.
754 check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
756 /* NOW the threads will be started: */
757 hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
762 if (hw_opt->nthreads_tmpi > 1)
764 t_commrec *cr_old = cr;
765 /* now start the threads. */
766 cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
767 oenv, bVerbose, bCompact, nstglobalcomm,
768 ddxyz, dd_node_order, rdd, rconstr,
769 dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
770 nbpu_opt, nstlist_cmdline,
771 nsteps_cmdline, nstepout, resetstep, nmultisim,
772 repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
773 cpt_period, max_hours,
775 /* the main thread continues here with a new cr. We don't deallocate
776 the old cr because other threads may still be reading it. */
779 gmx_comm("Failed to spawn threads");
784 /* END OF CAUTION: cr is now reliable */
786 /* g_membed initialisation *
787 * Because we change the mtop, init_membed is called before the init_parallel *
788 * (in case we ever want to make it run in parallel) */
789 if (opt2bSet("-membed", nfile, fnm))
793 fprintf(stderr, "Initializing membed");
795 membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
800 /* now broadcast everything to the non-master nodes/threads: */
801 init_parallel(cr, inputrec, mtop);
803 /* The master rank decided on the use of GPUs,
804 * broadcast this information to all ranks.
806 gmx_bcast_sim(sizeof(bUseGPU), &bUseGPU, cr);
811 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
812 fprintf(fplog, "\n");
815 /* now make sure the state is initialized and propagated */
816 set_state_entries(state, inputrec);
818 /* A parallel command line option consistency check that we can
819 only do after any threads have started. */
821 (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
824 "The -dd or -npme option request a parallel simulation, "
826 "but %s was compiled without threads or MPI enabled"
828 #ifdef GMX_THREAD_MPI
829 "but the number of threads (option -nt) is 1"
831 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
834 , output_env_get_program_display_name(oenv)
839 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
841 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
844 if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
846 gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
849 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
851 if (cr->npmenodes > 0)
853 gmx_fatal_collective(FARGS, cr, NULL,
854 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
860 if (bUseGPU && cr->npmenodes < 0)
862 /* With GPUs we don't automatically use PME-only ranks. PME ranks can
863 * improve performance with many threads per GPU, since our OpenMP
864 * scaling is bad, but it's difficult to automate the setup.
872 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
876 /* NMR restraints must be initialized before load_checkpoint,
877 * since with time averaging the history is added to t_state.
878 * For proper consistency check we therefore need to extend
880 * So the PME-only nodes (if present) will also initialize
881 * the distance restraints.
885 /* This needs to be called before read_checkpoint to extend the state */
886 init_disres(fplog, mtop, inputrec, cr, fcd, state, repl_ex_nst > 0);
888 init_orires(fplog, mtop, state->x, inputrec, cr, &(fcd->orires),
891 if (DEFORM(*inputrec))
893 /* Store the deform reference box before reading the checkpoint */
896 copy_mat(state->box, box);
900 gmx_bcast(sizeof(box), box, cr);
902 /* Because we do not have the update struct available yet
903 * in which the reference values should be stored,
904 * we store them temporarily in static variables.
905 * This should be thread safe, since they are only written once
906 * and with identical values.
908 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
909 deform_init_init_step_tpx = inputrec->init_step;
910 copy_mat(box, deform_init_box_tpx);
911 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
914 if (opt2bSet("-cpi", nfile, fnm))
916 /* Check if checkpoint file exists before doing continuation.
917 * This way we can use identical input options for the first and subsequent runs...
919 if (gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr) )
921 load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
923 inputrec, state, &bReadEkin,
924 (Flags & MD_APPENDFILES),
925 (Flags & MD_APPENDFILESSET));
929 Flags |= MD_READ_EKIN;
934 if (MASTER(cr) && (Flags & MD_APPENDFILES))
936 gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
940 /* override nsteps with value from cmdline */
941 override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
945 copy_mat(state->box, box);
950 gmx_bcast(sizeof(box), box, cr);
953 /* Essential dynamics */
954 if (opt2bSet("-ei", nfile, fnm))
956 /* Open input and output files, allocate space for ED data structure */
957 ed = ed_open(mtop->natoms, &state->edsamstate, nfile, fnm, Flags, oenv, cr);
960 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
961 inputrec->eI == eiNM))
963 cr->dd = init_domain_decomposition(fplog, cr, Flags, ddxyz, rdd, rconstr,
964 dddlb_opt, dlb_scale,
968 &ddbox, &npme_major, &npme_minor);
970 make_dd_communicators(fplog, cr, dd_node_order);
972 /* Set overallocation to avoid frequent reallocation of arrays */
973 set_over_alloc_dd(TRUE);
977 /* PME, if used, is done on all nodes with 1D decomposition */
979 cr->duty = (DUTY_PP | DUTY_PME);
983 if (inputrec->ePBC == epbcSCREW)
986 "pbc=%s is only implemented with domain decomposition",
987 epbc_names[inputrec->ePBC]);
993 /* After possible communicator splitting in make_dd_communicators.
994 * we can set up the intra/inter node communication.
996 gmx_setup_nodecomm(fplog, cr);
999 /* Initialize per-physical-node MPI process/thread ID and counters. */
1000 gmx_init_intranode_counters(cr);
1004 md_print_info(cr, fplog,
1005 "This is simulation %d out of %d running as a composite GROMACS\n"
1006 "multi-simulation job. Setup for this simulation:\n\n",
1007 cr->ms->sim, cr->ms->nsim);
1009 md_print_info(cr, fplog, "Using %d MPI %s\n",
1011 #ifdef GMX_THREAD_MPI
1012 cr->nnodes == 1 ? "thread" : "threads"
1014 cr->nnodes == 1 ? "process" : "processes"
1020 /* Check and update hw_opt for the cut-off scheme */
1021 check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1023 /* Check and update hw_opt for the number of MPI ranks */
1024 check_and_update_hw_opt_3(hw_opt);
1026 gmx_omp_nthreads_init(fplog, cr,
1027 hwinfo->nthreads_hw_avail,
1028 hw_opt->nthreads_omp,
1029 hw_opt->nthreads_omp_pme,
1030 (cr->duty & DUTY_PP) == 0,
1031 inputrec->cutoff_scheme == ecutsVERLET);
1034 if (integrator[inputrec->eI].func != do_tpi &&
1035 inputrec->cutoff_scheme == ecutsVERLET)
1037 gmx_feenableexcept();
1043 /* Select GPU id's to use */
1044 gmx_select_gpu_ids(fplog, cr, &hwinfo->gpu_info, bForceUseGPU,
1049 /* Ignore (potentially) manually selected GPUs */
1050 hw_opt->gpu_opt.n_dev_use = 0;
1053 /* check consistency across ranks of things like SIMD
1054 * support and number of GPUs selected */
1055 gmx_check_hw_runconf_consistency(fplog, hwinfo, cr, hw_opt, bUseGPU);
1057 /* Now that we know the setup is consistent, check for efficiency */
1058 check_resource_division_efficiency(hwinfo, hw_opt, Flags & MD_NTOMPSET,
1061 if (DOMAINDECOMP(cr))
1063 /* When we share GPUs over ranks, we need to know this for the DLB */
1064 dd_setup_dlb_resource_sharing(cr, hwinfo, hw_opt);
1067 /* getting number of PP/PME threads
1068 PME: env variable should be read only on one node to make sure it is
1069 identical everywhere;
1071 /* TODO nthreads_pp is only used for pinning threads.
1072 * This is a temporary solution until we have a hw topology library.
1074 nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
1075 nthreads_pme = gmx_omp_nthreads_get(emntPME);
1077 wcycle = wallcycle_init(fplog, resetstep, cr, nthreads_pp, nthreads_pme);
1081 /* Master synchronizes its value of reset_counters with all nodes
1082 * including PME only nodes */
1083 reset_counters = wcycle_get_reset_counters(wcycle);
1084 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1085 wcycle_set_reset_counters(wcycle, reset_counters);
1089 if (cr->duty & DUTY_PP)
1091 bcast_state(cr, state);
1093 /* Initiate forcerecord */
1095 fr->hwinfo = hwinfo;
1096 fr->gpu_opt = &hw_opt->gpu_opt;
1097 init_forcerec(fplog, oenv, fr, fcd, inputrec, mtop, cr, box,
1098 opt2fn("-table", nfile, fnm),
1099 opt2fn("-tabletf", nfile, fnm),
1100 opt2fn("-tablep", nfile, fnm),
1101 opt2fn("-tableb", nfile, fnm),
1106 /* version for PCA_NOT_READ_NODE (see md.c) */
1107 /*init_forcerec(fplog,fr,fcd,inputrec,mtop,cr,box,FALSE,
1108 "nofile","nofile","nofile","nofile",FALSE,pforce);
1111 /* Initialize QM-MM */
1114 init_QMMMrec(cr, mtop, inputrec, fr);
1117 /* Initialize the mdatoms structure.
1118 * mdatoms is not filled with atom data,
1119 * as this can not be done now with domain decomposition.
1121 mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1123 /* Initialize the virtual site communication */
1124 vsite = init_vsite(mtop, cr, FALSE);
1126 calc_shifts(box, fr->shift_vec);
1128 /* With periodic molecules the charge groups should be whole at start up
1129 * and the virtual sites should not be far from their proper positions.
1131 if (!inputrec->bContinuation && MASTER(cr) &&
1132 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1134 /* Make molecules whole at start of run */
1135 if (fr->ePBC != epbcNONE)
1137 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, state->x);
1141 /* Correct initial vsite positions are required
1142 * for the initial distribution in the domain decomposition
1143 * and for the initial shell prediction.
1145 construct_vsites_mtop(vsite, mtop, state->x);
1149 if (EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype))
1151 ewaldcoeff_q = fr->ewaldcoeff_q;
1152 ewaldcoeff_lj = fr->ewaldcoeff_lj;
1153 pmedata = &fr->pmedata;
1162 /* This is a PME only node */
1164 /* We don't need the state */
1167 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1168 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1172 if (hw_opt->thread_affinity != threadaffOFF)
1174 /* Before setting affinity, check whether the affinity has changed
1175 * - which indicates that probably the OpenMP library has changed it
1176 * since we first checked).
1178 gmx_check_thread_affinity_set(fplog, cr,
1179 hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1181 /* Set the CPU affinity */
1182 gmx_set_thread_affinity(fplog, cr, hw_opt, hwinfo);
1185 /* Initiate PME if necessary,
1186 * either on all nodes or on dedicated PME nodes only. */
1187 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1191 nChargePerturbed = mdatoms->nChargePerturbed;
1192 if (EVDW_PME(inputrec->vdwtype))
1194 nTypePerturbed = mdatoms->nTypePerturbed;
1197 if (cr->npmenodes > 0)
1199 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1200 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1201 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1204 if (cr->duty & DUTY_PME)
1206 status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1207 mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1208 (Flags & MD_REPRODUCIBLE), nthreads_pme);
1211 gmx_fatal(FARGS, "Error %d initializing PME", status);
1217 if (integrator[inputrec->eI].func == do_md)
1219 /* Turn on signal handling on all nodes */
1221 * (A user signal from the PME nodes (if any)
1222 * is communicated to the PP nodes.
1224 signal_handler_install();
1227 if (cr->duty & DUTY_PP)
1229 /* Assumes uniform use of the number of OpenMP threads */
1230 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1232 if (inputrec->bPull)
1234 /* Initialize pull code */
1235 inputrec->pull_work =
1236 init_pull(fplog, inputrec->pull, inputrec, nfile, fnm,
1237 mtop, cr, oenv, inputrec->fepvals->init_lambda,
1238 EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1243 /* Initialize enforced rotation code */
1244 init_rot(fplog, inputrec, nfile, fnm, cr, state->x, box, mtop, oenv,
1248 if (inputrec->eSwapCoords != eswapNO)
1250 /* Initialize ion swapping code */
1251 init_swapcoords(fplog, bVerbose, inputrec, opt2fn_master("-swap", nfile, fnm, cr),
1252 mtop, state->x, state->box, &state->swapstate, cr, oenv, Flags);
1255 constr = init_constraints(fplog, mtop, inputrec, ed, state, cr);
1257 if (DOMAINDECOMP(cr))
1259 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1260 dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1261 Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1263 set_dd_parameters(fplog, cr->dd, dlb_scale, inputrec, &ddbox);
1265 setup_dd_grid(fplog, cr->dd);
1268 /* Now do whatever the user wants us to do (how flexible...) */
1269 integrator[inputrec->eI].func(fplog, cr, nfile, fnm,
1270 oenv, bVerbose, bCompact,
1273 nstepout, inputrec, mtop,
1275 mdatoms, nrnb, wcycle, ed, fr,
1276 repl_ex_nst, repl_ex_nex, repl_ex_seed,
1278 cpt_period, max_hours,
1281 walltime_accounting);
1283 if (inputrec->bPull)
1285 finish_pull(inputrec->pull_work);
1290 finish_rot(inputrec->rot);
1296 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1298 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1299 gmx_pmeonly(*pmedata, cr, nrnb, wcycle, walltime_accounting, ewaldcoeff_q, ewaldcoeff_lj, inputrec);
1302 wallcycle_stop(wcycle, ewcRUN);
1304 /* Finish up, write some stuff
1305 * if rerunMD, don't write last frame again
1307 finish_run(fplog, cr,
1308 inputrec, nrnb, wcycle, walltime_accounting,
1309 fr ? fr->nbv : NULL,
1310 EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1313 /* Free GPU memory and context */
1314 free_gpu_resources(fr, cr, &hwinfo->gpu_info, fr ? fr->gpu_opt : NULL);
1316 if (opt2bSet("-membed", nfile, fnm))
1321 gmx_hardware_info_free(hwinfo);
1323 /* Does what it says */
1324 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1325 walltime_accounting_destroy(walltime_accounting);
1327 /* Close logfile already here if we were appending to it */
1328 if (MASTER(cr) && (Flags & MD_APPENDFILES))
1330 gmx_log_close(fplog);
1333 rc = (int)gmx_get_stop_condition();
1337 #ifdef GMX_THREAD_MPI
1338 /* we need to join all threads. The sub-threads join when they
1339 exit this function, but the master thread needs to be told to
1341 if (PAR(cr) && MASTER(cr))