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, 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.
54 #include "md_logging.h"
55 #include "md_support.h"
69 #include "checkpoint.h"
70 #include "mtop_util.h"
71 #include "sighandler.h"
73 #include "gmx_detect_hardware.h"
74 #include "gmx_omp_nthreads.h"
75 #include "calc_verletbuf.h"
76 #include "gmx_fatal_collective.h"
79 #include "gmx_thread_affinity.h"
81 #include "gromacs/fileio/tpxio.h"
82 #include "gromacs/mdlib/nbnxn_search.h"
83 #include "gromacs/mdlib/nbnxn_consts.h"
84 #include "gromacs/timing/wallcycle.h"
85 #include "gromacs/utility/gmxmpi.h"
86 #include "gromacs/utility/gmxomp.h"
87 #include "gromacs/swap/swapcoords.h"
88 #include "gromacs/essentialdynamics/edsam.h"
89 #include "gromacs/pulling/pull.h"
90 #include "gromacs/pulling/pull_rotation.h"
96 #include "gpu_utils.h"
97 #include "nbnxn_cuda_data_mgmt.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 struct mdrunner_arglist
127 const char *dddlb_opt;
132 const char *nbpu_opt;
134 gmx_int64_t nsteps_cmdline;
144 const char *deviceOptions;
146 int ret; /* return value */
150 /* The function used for spawning threads. Extracts the mdrunner()
151 arguments from its one argument and calls mdrunner(), after making
153 static void mdrunner_start_fn(void *arg)
155 struct mdrunner_arglist *mda = (struct mdrunner_arglist*)arg;
156 struct mdrunner_arglist mc = *mda; /* copy the arg list to make sure
157 that it's thread-local. This doesn't
158 copy pointed-to items, of course,
159 but those are all const. */
160 t_commrec *cr; /* we need a local version of this */
164 fnm = dup_tfn(mc.nfile, mc.fnm);
166 cr = reinitialize_commrec_for_this_thread(mc.cr);
173 mda->ret = mdrunner(&mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
174 mc.bVerbose, mc.bCompact, mc.nstglobalcomm,
175 mc.ddxyz, mc.dd_node_order, mc.rdd,
176 mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
177 mc.ddcsx, mc.ddcsy, mc.ddcsz,
178 mc.nbpu_opt, mc.nstlist_cmdline,
179 mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
180 mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce,
181 mc.cpt_period, mc.max_hours, mc.deviceOptions, mc.Flags);
184 /* called by mdrunner() to start a specific number of threads (including
185 the main thread) for thread-parallel runs. This in turn calls mdrunner()
187 All options besides nthreads are the same as for mdrunner(). */
188 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
189 FILE *fplog, t_commrec *cr, int nfile,
190 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
191 gmx_bool bCompact, int nstglobalcomm,
192 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
193 const char *dddlb_opt, real dlb_scale,
194 const char *ddcsx, const char *ddcsy, const char *ddcsz,
195 const char *nbpu_opt, int nstlist_cmdline,
196 gmx_int64_t nsteps_cmdline,
197 int nstepout, int resetstep,
198 int nmultisim, int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
199 real pforce, real cpt_period, real max_hours,
200 const char *deviceOptions, unsigned long Flags)
203 struct mdrunner_arglist *mda;
204 t_commrec *crn; /* the new commrec */
207 /* first check whether we even need to start tMPI */
208 if (hw_opt->nthreads_tmpi < 2)
213 /* a few small, one-time, almost unavoidable memory leaks: */
215 fnmn = dup_tfn(nfile, fnm);
217 /* fill the data structure to pass as void pointer to thread start fn */
218 /* hw_opt contains pointers, which should all be NULL at this stage */
219 mda->hw_opt = *hw_opt;
225 mda->bVerbose = bVerbose;
226 mda->bCompact = bCompact;
227 mda->nstglobalcomm = nstglobalcomm;
228 mda->ddxyz[XX] = ddxyz[XX];
229 mda->ddxyz[YY] = ddxyz[YY];
230 mda->ddxyz[ZZ] = ddxyz[ZZ];
231 mda->dd_node_order = dd_node_order;
233 mda->rconstr = rconstr;
234 mda->dddlb_opt = dddlb_opt;
235 mda->dlb_scale = dlb_scale;
239 mda->nbpu_opt = nbpu_opt;
240 mda->nstlist_cmdline = nstlist_cmdline;
241 mda->nsteps_cmdline = nsteps_cmdline;
242 mda->nstepout = nstepout;
243 mda->resetstep = resetstep;
244 mda->nmultisim = nmultisim;
245 mda->repl_ex_nst = repl_ex_nst;
246 mda->repl_ex_nex = repl_ex_nex;
247 mda->repl_ex_seed = repl_ex_seed;
248 mda->pforce = pforce;
249 mda->cpt_period = cpt_period;
250 mda->max_hours = max_hours;
251 mda->deviceOptions = deviceOptions;
254 /* now spawn new threads that start mdrunner_start_fn(), while
255 the main thread returns, we set thread affinity later */
256 ret = tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi, TMPI_AFFINITY_NONE,
257 mdrunner_start_fn, (void*)(mda) );
258 if (ret != TMPI_SUCCESS)
263 crn = reinitialize_commrec_for_this_thread(cr);
268 static int get_tmpi_omp_thread_division(const gmx_hw_info_t *hwinfo,
269 const gmx_hw_opt_t *hw_opt,
275 /* There are no separate PME nodes here, as we ensured in
276 * check_and_update_hw_opt that nthreads_tmpi>0 with PME nodes
277 * and a conditional ensures we would not have ended up here.
278 * Note that separate PME nodes might be switched on later.
282 nthreads_tmpi = ngpu;
283 if (nthreads_tot > 0 && nthreads_tot < nthreads_tmpi)
285 nthreads_tmpi = nthreads_tot;
288 else if (hw_opt->nthreads_omp > 0)
290 /* Here we could oversubscribe, when we do, we issue a warning later */
291 nthreads_tmpi = max(1, nthreads_tot/hw_opt->nthreads_omp);
295 /* TODO choose nthreads_omp based on hardware topology
296 when we have a hardware topology detection library */
297 /* In general, when running up to 4 threads, OpenMP should be faster.
298 * Note: on AMD Bulldozer we should avoid running OpenMP over two dies.
299 * On Intel>=Nehalem running OpenMP on a single CPU is always faster,
300 * even on two CPUs it's usually faster (but with many OpenMP threads
301 * it could be faster not to use HT, currently we always use HT).
302 * On Nehalem/Westmere we want to avoid running 16 threads over
303 * two CPUs with HT, so we need a limit<16; thus we use 12.
304 * A reasonable limit for Intel Sandy and Ivy bridge,
305 * not knowing the topology, is 16 threads.
307 const int nthreads_omp_always_faster = 4;
308 const int nthreads_omp_always_faster_Nehalem = 12;
309 const int nthreads_omp_always_faster_SandyBridge = 16;
310 const int first_model_Nehalem = 0x1A;
311 const int first_model_SandyBridge = 0x2A;
312 gmx_bool bIntel_Family6;
315 (gmx_cpuid_vendor(hwinfo->cpuid_info) == GMX_CPUID_VENDOR_INTEL &&
316 gmx_cpuid_family(hwinfo->cpuid_info) == 6);
318 if (nthreads_tot <= nthreads_omp_always_faster ||
320 ((gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_Nehalem && nthreads_tot <= nthreads_omp_always_faster_Nehalem) ||
321 (gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_SandyBridge && nthreads_tot <= nthreads_omp_always_faster_SandyBridge))))
323 /* Use pure OpenMP parallelization */
328 /* Don't use OpenMP parallelization */
329 nthreads_tmpi = nthreads_tot;
333 return nthreads_tmpi;
337 /* Get the number of threads to use for thread-MPI based on how many
338 * were requested, which algorithms we're using,
339 * and how many particles there are.
340 * At the point we have already called check_and_update_hw_opt.
341 * Thus all options should be internally consistent and consistent
342 * with the hardware, except that ntmpi could be larger than #GPU.
344 static int get_nthreads_mpi(const gmx_hw_info_t *hwinfo,
345 gmx_hw_opt_t *hw_opt,
346 t_inputrec *inputrec, gmx_mtop_t *mtop,
350 int nthreads_hw, nthreads_tot_max, nthreads_tmpi, nthreads_new, ngpu;
351 int min_atoms_per_mpi_thread;
356 if (hw_opt->nthreads_tmpi > 0)
358 /* Trivial, return right away */
359 return hw_opt->nthreads_tmpi;
362 nthreads_hw = hwinfo->nthreads_hw_avail;
364 /* How many total (#tMPI*#OpenMP) threads can we start? */
365 if (hw_opt->nthreads_tot > 0)
367 nthreads_tot_max = hw_opt->nthreads_tot;
371 nthreads_tot_max = nthreads_hw;
374 bCanUseGPU = (inputrec->cutoff_scheme == ecutsVERLET &&
375 hwinfo->gpu_info.ncuda_dev_compatible > 0);
378 ngpu = hwinfo->gpu_info.ncuda_dev_compatible;
385 if (inputrec->cutoff_scheme == ecutsGROUP)
387 /* We checked this before, but it doesn't hurt to do it once more */
388 assert(hw_opt->nthreads_omp == 1);
392 get_tmpi_omp_thread_division(hwinfo, hw_opt, nthreads_tot_max, ngpu);
394 if (inputrec->eI == eiNM || EI_TPI(inputrec->eI))
396 /* Dims/steps are divided over the nodes iso splitting the atoms */
397 min_atoms_per_mpi_thread = 0;
403 min_atoms_per_mpi_thread = MIN_ATOMS_PER_GPU;
407 min_atoms_per_mpi_thread = MIN_ATOMS_PER_MPI_THREAD;
411 /* Check if an algorithm does not support parallel simulation. */
412 if (nthreads_tmpi != 1 &&
413 ( inputrec->eI == eiLBFGS ||
414 inputrec->coulombtype == eelEWALD ) )
418 md_print_warn(cr, fplog, "The integration or electrostatics algorithm doesn't support parallel runs. Using a single thread-MPI thread.\n");
419 if (hw_opt->nthreads_tmpi > nthreads_tmpi)
421 gmx_fatal(FARGS, "You asked for more than 1 thread-MPI thread, but an algorithm doesn't support that");
424 else if (mtop->natoms/nthreads_tmpi < min_atoms_per_mpi_thread)
426 /* the thread number was chosen automatically, but there are too many
427 threads (too few atoms per thread) */
428 nthreads_new = max(1, mtop->natoms/min_atoms_per_mpi_thread);
430 /* Avoid partial use of Hyper-Threading */
431 if (gmx_cpuid_x86_smt(hwinfo->cpuid_info) == GMX_CPUID_X86_SMT_ENABLED &&
432 nthreads_new > nthreads_hw/2 && nthreads_new < nthreads_hw)
434 nthreads_new = nthreads_hw/2;
437 /* Avoid large prime numbers in the thread count */
438 if (nthreads_new >= 6)
440 /* Use only 6,8,10 with additional factors of 2 */
444 while (3*fac*2 <= nthreads_new)
449 nthreads_new = (nthreads_new/fac)*fac;
454 if (nthreads_new == 5)
460 nthreads_tmpi = nthreads_new;
462 fprintf(stderr, "\n");
463 fprintf(stderr, "NOTE: Parallelization is limited by the small number of atoms,\n");
464 fprintf(stderr, " only starting %d thread-MPI threads.\n", nthreads_tmpi);
465 fprintf(stderr, " You can use the -nt and/or -ntmpi option to optimize the number of threads.\n\n");
468 return nthreads_tmpi;
470 #endif /* GMX_THREAD_MPI */
473 /* We determine the extra cost of the non-bonded kernels compared to
474 * a reference nstlist value of 10 (which is the default in grompp).
476 static const int nbnxn_reference_nstlist = 10;
477 /* The values to try when switching */
478 const int nstlist_try[] = { 20, 25, 40 };
479 #define NNSTL sizeof(nstlist_try)/sizeof(nstlist_try[0])
480 /* Increase nstlist until the non-bonded cost increases more than listfac_ok,
481 * but never more than listfac_max.
482 * A standard (protein+)water system at 300K with PME ewald_rtol=1e-5
483 * needs 1.28 at rcoulomb=0.9 and 1.24 at rcoulomb=1.0 to get to nstlist=40.
484 * Note that both CPU and GPU factors are conservative. Performance should
485 * not go down due to this tuning, except with a relatively slow GPU.
486 * On the other hand, at medium/high parallelization or with fast GPUs
487 * nstlist will not be increased enough to reach optimal performance.
489 /* CPU: pair-search is about a factor 1.5 slower than the non-bonded kernel */
490 static const float nbnxn_cpu_listfac_ok = 1.05;
491 static const float nbnxn_cpu_listfac_max = 1.09;
492 /* GPU: pair-search is a factor 1.5-3 slower than the non-bonded kernel */
493 static const float nbnxn_gpu_listfac_ok = 1.20;
494 static const float nbnxn_gpu_listfac_max = 1.30;
496 /* Try to increase nstlist when using the Verlet cut-off scheme */
497 static void increase_nstlist(FILE *fp, t_commrec *cr,
498 t_inputrec *ir, int nstlist_cmdline,
499 const gmx_mtop_t *mtop, matrix box,
502 float listfac_ok, listfac_max;
503 int nstlist_orig, nstlist_prev;
504 verletbuf_list_setup_t ls;
505 real rlist_nstlist10, rlist_inc, rlist_ok, rlist_max;
506 real rlist_new, rlist_prev;
509 gmx_bool bBox, bDD, bCont;
510 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";
511 const char *vbd_err = "Can not increase nstlist because verlet-buffer-tolerance is not set or used";
512 const char *box_err = "Can not increase nstlist because the box is too small";
513 const char *dd_err = "Can not increase nstlist because of domain decomposition limitations";
516 if (nstlist_cmdline <= 0)
518 if (fp != NULL && bGPU && ir->nstlist < nstlist_try[0])
520 fprintf(fp, nstl_gpu, ir->nstlist);
523 while (nstlist_ind < NNSTL && ir->nstlist >= nstlist_try[nstlist_ind])
527 if (nstlist_ind == NNSTL)
529 /* There are no larger nstlist value to try */
534 if (ir->verletbuf_tol == 0 && bGPU)
536 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");
539 if (ir->verletbuf_tol < 0)
543 fprintf(stderr, "%s\n", vbd_err);
547 fprintf(fp, "%s\n", vbd_err);
555 listfac_ok = nbnxn_gpu_listfac_ok;
556 listfac_max = nbnxn_gpu_listfac_max;
560 listfac_ok = nbnxn_cpu_listfac_ok;
561 listfac_max = nbnxn_cpu_listfac_max;
564 nstlist_orig = ir->nstlist;
565 if (nstlist_cmdline > 0)
569 sprintf(buf, "Getting nstlist=%d from command line option",
572 ir->nstlist = nstlist_cmdline;
575 verletbuf_get_list_setup(bGPU, &ls);
577 /* Allow rlist to make the list a given factor larger than the list
578 * would be with nstlist=10.
580 nstlist_prev = ir->nstlist;
582 calc_verlet_buffer_size(mtop, det(box), ir, &ls, NULL, &rlist_nstlist10);
583 ir->nstlist = nstlist_prev;
585 /* Determine the pair list size increase due to zero interactions */
586 rlist_inc = nbnxn_get_rlist_effective_inc(ls.cluster_size_j,
587 mtop->natoms/det(box));
588 rlist_ok = (rlist_nstlist10 + rlist_inc)*pow(listfac_ok, 1.0/3.0) - rlist_inc;
589 rlist_max = (rlist_nstlist10 + rlist_inc)*pow(listfac_max, 1.0/3.0) - rlist_inc;
592 fprintf(debug, "nstlist tuning: rlist_inc %.3f rlist_ok %.3f rlist_max %.3f\n",
593 rlist_inc, rlist_ok, rlist_max);
596 nstlist_prev = nstlist_orig;
597 rlist_prev = ir->rlist;
600 if (nstlist_cmdline <= 0)
602 ir->nstlist = nstlist_try[nstlist_ind];
605 /* Set the pair-list buffer size in ir */
606 calc_verlet_buffer_size(mtop, det(box), ir, &ls, NULL, &rlist_new);
608 /* Does rlist fit in the box? */
609 bBox = (sqr(rlist_new) < max_cutoff2(ir->ePBC, box));
611 if (bBox && DOMAINDECOMP(cr))
613 /* Check if rlist fits in the domain decomposition */
614 if (inputrec2nboundeddim(ir) < DIM)
616 gmx_incons("Changing nstlist with domain decomposition and unbounded dimensions is not implemented yet");
618 copy_mat(box, state_tmp.box);
619 bDD = change_dd_cutoff(cr, &state_tmp, ir, rlist_new);
624 fprintf(debug, "nstlist %d rlist %.3f bBox %d bDD %d\n",
625 ir->nstlist, rlist_new, bBox, bDD);
630 if (nstlist_cmdline <= 0)
632 if (bBox && bDD && rlist_new <= rlist_max)
634 /* Increase nstlist */
635 nstlist_prev = ir->nstlist;
636 rlist_prev = rlist_new;
637 bCont = (nstlist_ind+1 < NNSTL && rlist_new < rlist_ok);
641 /* Stick with the previous nstlist */
642 ir->nstlist = nstlist_prev;
643 rlist_new = rlist_prev;
655 gmx_warning(!bBox ? box_err : dd_err);
658 fprintf(fp, "\n%s\n", bBox ? box_err : dd_err);
660 ir->nstlist = nstlist_orig;
662 else if (ir->nstlist != nstlist_orig || rlist_new != ir->rlist)
664 sprintf(buf, "Changing nstlist from %d to %d, rlist from %g to %g",
665 nstlist_orig, ir->nstlist,
666 ir->rlist, rlist_new);
669 fprintf(stderr, "%s\n\n", buf);
673 fprintf(fp, "%s\n\n", buf);
675 ir->rlist = rlist_new;
676 ir->rlistlong = rlist_new;
680 static void prepare_verlet_scheme(FILE *fplog,
684 const gmx_mtop_t *mtop,
688 if (ir->verletbuf_tol > 0)
690 /* Update the Verlet buffer size for the current run setup */
691 verletbuf_list_setup_t ls;
694 /* Here we assume CPU acceleration is on. But as currently
695 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
696 * and 4x2 gives a larger buffer than 4x4, this is ok.
698 verletbuf_get_list_setup(bUseGPU, &ls);
700 calc_verlet_buffer_size(mtop, det(box), ir, &ls, NULL, &rlist_new);
702 if (rlist_new != ir->rlist)
706 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
707 ir->rlist, rlist_new,
708 ls.cluster_size_i, ls.cluster_size_j);
710 ir->rlist = rlist_new;
711 ir->rlistlong = rlist_new;
715 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
717 gmx_fatal(FARGS, "Can not set nstlist without %s",
718 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
721 if (EI_DYNAMICS(ir->eI))
723 /* Set or try nstlist values */
724 increase_nstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, bUseGPU);
728 static void convert_to_verlet_scheme(FILE *fplog,
730 gmx_mtop_t *mtop, real box_vol)
732 char *conv_mesg = "Converting input file with group cut-off scheme to the Verlet cut-off scheme";
734 md_print_warn(NULL, fplog, "%s\n", conv_mesg);
736 ir->cutoff_scheme = ecutsVERLET;
737 ir->verletbuf_tol = 0.005;
739 if (ir->rcoulomb != ir->rvdw)
741 gmx_fatal(FARGS, "The VdW and Coulomb cut-offs are different, whereas the Verlet scheme only supports equal cut-offs");
744 if (ir->vdwtype == evdwUSER || EEL_USER(ir->coulombtype))
746 gmx_fatal(FARGS, "User non-bonded potentials are not (yet) supported with the Verlet scheme");
748 else if (EVDW_SWITCHED(ir->vdwtype) || EEL_SWITCHED(ir->coulombtype))
750 md_print_warn(NULL, fplog, "Converting switched or shifted interactions to a shifted potential (without force shift), this will lead to slightly different interaction potentials");
752 if (EVDW_SWITCHED(ir->vdwtype))
754 ir->vdwtype = evdwCUT;
756 if (EEL_SWITCHED(ir->coulombtype))
758 if (EEL_FULL(ir->coulombtype))
760 /* With full electrostatic only PME can be switched */
761 ir->coulombtype = eelPME;
765 md_print_warn(NULL, fplog, "NOTE: Replacing %s electrostatics with reaction-field with epsilon-rf=inf\n", eel_names[ir->coulombtype]);
766 ir->coulombtype = eelRF;
767 ir->epsilon_rf = 0.0;
771 /* We set the pair energy error tolerance to a small number.
772 * Note that this is only for testing. For production the user
773 * should think about this and set the mdp options.
775 ir->verletbuf_tol = 1e-4;
778 if (inputrec2nboundeddim(ir) != 3)
780 gmx_fatal(FARGS, "Can only convert old tpr files to the Verlet cut-off scheme with 3D pbc");
783 if (ir->efep != efepNO || ir->implicit_solvent != eisNO)
785 gmx_fatal(FARGS, "Will not convert old tpr files to the Verlet cut-off scheme with free-energy calculations or implicit solvent");
788 if (EI_DYNAMICS(ir->eI) && !(EI_MD(ir->eI) && ir->etc == etcNO))
790 verletbuf_list_setup_t ls;
792 verletbuf_get_list_setup(FALSE, &ls);
793 calc_verlet_buffer_size(mtop, box_vol, ir, &ls, NULL, &ir->rlist);
797 ir->verletbuf_tol = -1;
798 ir->rlist = 1.05*max(ir->rvdw, ir->rcoulomb);
801 gmx_mtop_remove_chargegroups(mtop);
804 static void print_hw_opt(FILE *fp, const gmx_hw_opt_t *hw_opt)
806 fprintf(fp, "hw_opt: nt %d ntmpi %d ntomp %d ntomp_pme %d gpu_id '%s'\n",
807 hw_opt->nthreads_tot,
808 hw_opt->nthreads_tmpi,
809 hw_opt->nthreads_omp,
810 hw_opt->nthreads_omp_pme,
811 hw_opt->gpu_opt.gpu_id != NULL ? hw_opt->gpu_opt.gpu_id : "");
814 /* Checks we can do when we don't (yet) know the cut-off scheme */
815 static void check_and_update_hw_opt_1(gmx_hw_opt_t *hw_opt,
816 gmx_bool bIsSimMaster)
818 gmx_omp_nthreads_read_env(&hw_opt->nthreads_omp, bIsSimMaster);
820 #ifndef GMX_THREAD_MPI
821 if (hw_opt->nthreads_tot > 0)
823 gmx_fatal(FARGS, "Setting the total number of threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
825 if (hw_opt->nthreads_tmpi > 0)
827 gmx_fatal(FARGS, "Setting the number of thread-MPI threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
832 if (hw_opt->nthreads_omp > 1)
834 gmx_fatal(FARGS, "More than 1 OpenMP thread requested, but Gromacs was compiled without OpenMP support");
836 hw_opt->nthreads_omp = 1;
839 if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp_pme <= 0)
841 /* We have the same number of OpenMP threads for PP and PME processes,
842 * thus we can perform several consistency checks.
844 if (hw_opt->nthreads_tmpi > 0 &&
845 hw_opt->nthreads_omp > 0 &&
846 hw_opt->nthreads_tot != hw_opt->nthreads_tmpi*hw_opt->nthreads_omp)
848 gmx_fatal(FARGS, "The total number of threads requested (%d) does not match the thread-MPI threads (%d) times the OpenMP threads (%d) requested",
849 hw_opt->nthreads_tot, hw_opt->nthreads_tmpi, hw_opt->nthreads_omp);
852 if (hw_opt->nthreads_tmpi > 0 &&
853 hw_opt->nthreads_tot % hw_opt->nthreads_tmpi != 0)
855 gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of thread-MPI threads requested (%d)",
856 hw_opt->nthreads_tot, hw_opt->nthreads_tmpi);
859 if (hw_opt->nthreads_omp > 0 &&
860 hw_opt->nthreads_tot % hw_opt->nthreads_omp != 0)
862 gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of OpenMP threads requested (%d)",
863 hw_opt->nthreads_tot, hw_opt->nthreads_omp);
866 if (hw_opt->nthreads_tmpi > 0 &&
867 hw_opt->nthreads_omp <= 0)
869 hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
874 if (hw_opt->nthreads_omp > 1)
876 gmx_fatal(FARGS, "OpenMP threads are requested, but Gromacs was compiled without OpenMP support");
880 if (hw_opt->nthreads_omp_pme > 0 && hw_opt->nthreads_omp <= 0)
882 gmx_fatal(FARGS, "You need to specify -ntomp in addition to -ntomp_pme");
885 if (hw_opt->nthreads_tot == 1)
887 hw_opt->nthreads_tmpi = 1;
889 if (hw_opt->nthreads_omp > 1)
891 gmx_fatal(FARGS, "You requested %d OpenMP threads with %d total threads",
892 hw_opt->nthreads_tmpi, hw_opt->nthreads_tot);
894 hw_opt->nthreads_omp = 1;
897 if (hw_opt->nthreads_omp_pme <= 0 && hw_opt->nthreads_omp > 0)
899 hw_opt->nthreads_omp_pme = hw_opt->nthreads_omp;
902 /* Parse GPU IDs, if provided.
903 * We check consistency with the tMPI thread count later.
905 gmx_parse_gpu_ids(&hw_opt->gpu_opt);
907 #ifdef GMX_THREAD_MPI
908 if (hw_opt->gpu_opt.ncuda_dev_use > 0 && hw_opt->nthreads_tmpi == 0)
910 /* Set the number of MPI threads equal to the number of GPUs */
911 hw_opt->nthreads_tmpi = hw_opt->gpu_opt.ncuda_dev_use;
913 if (hw_opt->nthreads_tot > 0 &&
914 hw_opt->nthreads_tmpi > hw_opt->nthreads_tot)
916 /* We have more GPUs than total threads requested.
917 * We choose to (later) generate a mismatch error,
918 * instead of launching more threads than requested.
920 hw_opt->nthreads_tmpi = hw_opt->nthreads_tot;
927 print_hw_opt(debug, hw_opt);
931 /* Checks we can do when we know the cut-off scheme */
932 static void check_and_update_hw_opt_2(gmx_hw_opt_t *hw_opt,
935 if (cutoff_scheme == ecutsGROUP)
937 /* We only have OpenMP support for PME only nodes */
938 if (hw_opt->nthreads_omp > 1)
940 gmx_fatal(FARGS, "OpenMP threads have been requested with cut-off scheme %s, but these are only supported with cut-off scheme %s",
941 ecutscheme_names[cutoff_scheme],
942 ecutscheme_names[ecutsVERLET]);
944 hw_opt->nthreads_omp = 1;
947 if (hw_opt->nthreads_omp_pme <= 0 && hw_opt->nthreads_omp > 0)
949 hw_opt->nthreads_omp_pme = hw_opt->nthreads_omp;
954 print_hw_opt(debug, hw_opt);
959 /* Override the value in inputrec with value passed on the command line (if any) */
960 static void override_nsteps_cmdline(FILE *fplog,
961 gmx_int64_t nsteps_cmdline,
965 char sbuf[STEPSTRSIZE];
970 /* override with anything else than the default -2 */
971 if (nsteps_cmdline > -2)
975 ir->nsteps = nsteps_cmdline;
976 if (EI_DYNAMICS(ir->eI))
978 sprintf(stmp, "Overriding nsteps with value passed on the command line: %s steps, %.3f ps",
979 gmx_step_str(nsteps_cmdline, sbuf),
980 nsteps_cmdline*ir->delta_t);
984 sprintf(stmp, "Overriding nsteps with value passed on the command line: %s steps",
985 gmx_step_str(nsteps_cmdline, sbuf));
988 md_print_warn(cr, fplog, "%s\n", stmp);
992 /* Frees GPU memory and destroys the CUDA context.
994 * Note that this function needs to be called even if GPUs are not used
995 * in this run because the PME ranks have no knowledge of whether GPUs
996 * are used or not, but all ranks need to enter the barrier below.
998 static void free_gpu_resources(const t_forcerec *fr,
1001 gmx_bool bIsPPrankUsingGPU;
1002 char gpu_err_str[STRLEN];
1004 bIsPPrankUsingGPU = (cr->duty & DUTY_PP) && fr->nbv != NULL && fr->nbv->bUseGPU;
1006 if (bIsPPrankUsingGPU)
1008 /* free nbnxn data in GPU memory */
1009 nbnxn_cuda_free(fr->nbv->cu_nbv);
1011 /* With tMPI we need to wait for all ranks to finish deallocation before
1012 * destroying the context in free_gpu() as some ranks may be sharing
1014 * Note: as only PP ranks need to free GPU resources, so it is safe to
1015 * not call the barrier on PME ranks.
1017 #ifdef GMX_THREAD_MPI
1022 #endif /* GMX_THREAD_MPI */
1024 /* uninitialize GPU (by destroying the context) */
1025 if (!free_gpu(gpu_err_str))
1027 gmx_warning("On node %d failed to free GPU #%d: %s",
1028 cr->nodeid, get_current_gpu_device_id(), gpu_err_str);
1033 int mdrunner(gmx_hw_opt_t *hw_opt,
1034 FILE *fplog, t_commrec *cr, int nfile,
1035 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
1036 gmx_bool bCompact, int nstglobalcomm,
1037 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
1038 const char *dddlb_opt, real dlb_scale,
1039 const char *ddcsx, const char *ddcsy, const char *ddcsz,
1040 const char *nbpu_opt, int nstlist_cmdline,
1041 gmx_int64_t nsteps_cmdline, int nstepout, int resetstep,
1042 int gmx_unused nmultisim, int repl_ex_nst, int repl_ex_nex,
1043 int repl_ex_seed, real pforce, real cpt_period, real max_hours,
1044 const char *deviceOptions, unsigned long Flags)
1046 gmx_bool bForceUseGPU, bTryUseGPU;
1047 double nodetime = 0, realtime;
1048 t_inputrec *inputrec;
1049 t_state *state = NULL;
1051 gmx_ddbox_t ddbox = {0};
1052 int npme_major, npme_minor;
1055 gmx_mtop_t *mtop = NULL;
1056 t_mdatoms *mdatoms = NULL;
1057 t_forcerec *fr = NULL;
1058 t_fcdata *fcd = NULL;
1059 real ewaldcoeff_q = 0;
1060 real ewaldcoeff_lj = 0;
1061 gmx_pme_t *pmedata = NULL;
1062 gmx_vsite_t *vsite = NULL;
1063 gmx_constr_t constr;
1064 int i, m, nChargePerturbed = -1, nTypePerturbed = 0, status, nalloc;
1066 gmx_wallcycle_t wcycle;
1067 gmx_bool bReadRNG, bReadEkin;
1069 gmx_walltime_accounting_t walltime_accounting = NULL;
1071 gmx_int64_t reset_counters;
1072 gmx_edsam_t ed = NULL;
1073 t_commrec *cr_old = cr;
1074 int nthreads_pme = 1;
1075 int nthreads_pp = 1;
1076 gmx_membed_t membed = NULL;
1077 gmx_hw_info_t *hwinfo = NULL;
1078 /* The master rank decides early on bUseGPU and broadcasts this later */
1079 gmx_bool bUseGPU = FALSE;
1081 /* CAUTION: threads may be started later on in this function, so
1082 cr doesn't reflect the final parallel state right now */
1086 if (Flags & MD_APPENDFILES)
1091 bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
1092 bTryUseGPU = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
1094 /* Detect hardware, gather information. This is an operation that is
1095 * global for this process (MPI rank). */
1096 hwinfo = gmx_detect_hardware(fplog, cr, bTryUseGPU);
1102 /* Read (nearly) all data required for the simulation */
1103 read_tpx_state(ftp2fn(efTPX, nfile, fnm), inputrec, state, NULL, mtop);
1105 if (inputrec->cutoff_scheme != ecutsVERLET &&
1106 ((Flags & MD_TESTVERLET) || getenv("GMX_VERLET_SCHEME") != NULL))
1108 convert_to_verlet_scheme(fplog, inputrec, mtop, det(state->box));
1111 if (inputrec->cutoff_scheme == ecutsVERLET)
1113 /* Here the master rank decides if all ranks will use GPUs */
1114 bUseGPU = (hwinfo->gpu_info.ncuda_dev_compatible > 0 ||
1115 getenv("GMX_EMULATE_GPU") != NULL);
1117 prepare_verlet_scheme(fplog, cr,
1118 inputrec, nstlist_cmdline, mtop, state->box,
1123 if (nstlist_cmdline > 0)
1125 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
1128 if (hwinfo->gpu_info.ncuda_dev_compatible > 0)
1130 md_print_warn(cr, fplog,
1131 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
1132 " To use a GPU, set the mdp option: cutoff-scheme = Verlet\n"
1133 " (for quick performance testing you can use the -testverlet option)\n");
1138 gmx_fatal(FARGS, "GPU requested, but can't be used without cutoff-scheme=Verlet");
1141 #ifdef GMX_TARGET_BGQ
1142 md_print_warn(cr, fplog,
1143 "NOTE: There is no SIMD implementation of the group scheme kernels on\n"
1144 " BlueGene/Q. You will observe better performance from using the\n"
1145 " Verlet cut-off scheme.\n");
1150 /* Check for externally set OpenMP affinity and turn off internal
1151 * pinning if any is found. We need to do this check early to tell
1152 * thread-MPI whether it should do pinning when spawning threads.
1153 * TODO: the above no longer holds, we should move these checks down
1155 gmx_omp_check_thread_affinity(fplog, cr, hw_opt);
1157 /* Check and update the hardware options for internal consistency */
1158 check_and_update_hw_opt_1(hw_opt, SIMMASTER(cr));
1162 #ifdef GMX_THREAD_MPI
1163 /* Early check for externally set process affinity.
1164 * With thread-MPI this is needed as pinning might get turned off,
1165 * which needs to be known before starting thread-MPI.
1166 * With thread-MPI hw_opt is processed here on the master rank
1167 * and passed to the other ranks later, so we only do this on master.
1169 gmx_check_thread_affinity_set(fplog,
1171 hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1174 #ifdef GMX_THREAD_MPI
1175 if (cr->npmenodes > 0 && hw_opt->nthreads_tmpi <= 0)
1177 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME nodes");
1181 if (hw_opt->nthreads_omp_pme != hw_opt->nthreads_omp &&
1184 gmx_fatal(FARGS, "You need to explicitly specify the number of PME nodes (-npme) when using different number of OpenMP threads for PP and PME nodes");
1188 #ifdef GMX_THREAD_MPI
1191 /* Since the master knows the cut-off scheme, update hw_opt for this.
1192 * This is done later for normal MPI and also once more with tMPI
1193 * for all tMPI ranks.
1195 check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1197 /* NOW the threads will be started: */
1198 hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
1202 if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp <= 0)
1204 hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
1207 if (hw_opt->nthreads_tmpi > 1)
1209 /* now start the threads. */
1210 cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
1211 oenv, bVerbose, bCompact, nstglobalcomm,
1212 ddxyz, dd_node_order, rdd, rconstr,
1213 dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
1214 nbpu_opt, nstlist_cmdline,
1215 nsteps_cmdline, nstepout, resetstep, nmultisim,
1216 repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
1217 cpt_period, max_hours, deviceOptions,
1219 /* the main thread continues here with a new cr. We don't deallocate
1220 the old cr because other threads may still be reading it. */
1223 gmx_comm("Failed to spawn threads");
1228 /* END OF CAUTION: cr is now reliable */
1230 /* g_membed initialisation *
1231 * Because we change the mtop, init_membed is called before the init_parallel *
1232 * (in case we ever want to make it run in parallel) */
1233 if (opt2bSet("-membed", nfile, fnm))
1237 fprintf(stderr, "Initializing membed");
1239 membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
1244 /* now broadcast everything to the non-master nodes/threads: */
1245 init_parallel(cr, inputrec, mtop);
1247 /* This check needs to happen after get_nthreads_mpi() */
1248 if (inputrec->cutoff_scheme == ecutsVERLET && (Flags & MD_PARTDEC))
1250 gmx_fatal_collective(FARGS, cr, NULL,
1251 "The Verlet cut-off scheme is not supported with particle decomposition.\n"
1252 "You can achieve the same effect as particle decomposition by running in parallel using only OpenMP threads.");
1257 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
1260 /* now make sure the state is initialized and propagated */
1261 set_state_entries(state, inputrec, cr->nnodes);
1263 /* A parallel command line option consistency check that we can
1264 only do after any threads have started. */
1266 (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
1269 "The -dd or -npme option request a parallel simulation, "
1271 "but %s was compiled without threads or MPI enabled"
1273 #ifdef GMX_THREAD_MPI
1274 "but the number of threads (option -nt) is 1"
1276 "but %s was not started through mpirun/mpiexec or only one process was requested through mpirun/mpiexec"
1283 if ((Flags & MD_RERUN) &&
1284 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1286 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
1289 if (can_use_allvsall(inputrec, TRUE, cr, fplog) && PAR(cr))
1291 /* Simple neighbour searching and (also?) all-vs-all loops
1292 * do not work with domain decomposition. */
1293 Flags |= MD_PARTDEC;
1296 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)) || (Flags & MD_PARTDEC))
1298 if (cr->npmenodes > 0)
1300 if (!EEL_PME(inputrec->coulombtype) && !EVDW_PME(inputrec->vdwtype))
1302 gmx_fatal_collective(FARGS, cr, NULL,
1303 "PME nodes are requested, but the system does not use PME electrostatics or LJ-PME");
1305 if (Flags & MD_PARTDEC)
1307 gmx_fatal_collective(FARGS, cr, NULL,
1308 "PME nodes are requested, but particle decomposition does not support separate PME nodes");
1318 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1322 /* NMR restraints must be initialized before load_checkpoint,
1323 * since with time averaging the history is added to t_state.
1324 * For proper consistency check we therefore need to extend
1326 * So the PME-only nodes (if present) will also initialize
1327 * the distance restraints.
1331 /* This needs to be called before read_checkpoint to extend the state */
1332 init_disres(fplog, mtop, inputrec, cr, Flags & MD_PARTDEC, fcd, state, repl_ex_nst > 0);
1334 if (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0)
1336 if (PAR(cr) && !(Flags & MD_PARTDEC))
1338 gmx_fatal(FARGS, "Orientation restraints do not work (yet) with domain decomposition, use particle decomposition (mdrun option -pd)");
1340 /* Orientation restraints */
1343 init_orires(fplog, mtop, state->x, inputrec, cr->ms, &(fcd->orires),
1348 if (DEFORM(*inputrec))
1350 /* Store the deform reference box before reading the checkpoint */
1353 copy_mat(state->box, box);
1357 gmx_bcast(sizeof(box), box, cr);
1359 /* Because we do not have the update struct available yet
1360 * in which the reference values should be stored,
1361 * we store them temporarily in static variables.
1362 * This should be thread safe, since they are only written once
1363 * and with identical values.
1365 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
1366 deform_init_init_step_tpx = inputrec->init_step;
1367 copy_mat(box, deform_init_box_tpx);
1368 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
1371 if (opt2bSet("-cpi", nfile, fnm))
1373 /* Check if checkpoint file exists before doing continuation.
1374 * This way we can use identical input options for the first and subsequent runs...
1376 if (gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr) )
1378 load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
1379 cr, Flags & MD_PARTDEC, ddxyz,
1380 inputrec, state, &bReadRNG, &bReadEkin,
1381 (Flags & MD_APPENDFILES),
1382 (Flags & MD_APPENDFILESSET));
1386 Flags |= MD_READ_RNG;
1390 Flags |= MD_READ_EKIN;
1395 if (((MASTER(cr) || (Flags & MD_SEPPOT)) && (Flags & MD_APPENDFILES))
1396 #ifdef GMX_THREAD_MPI
1397 /* With thread MPI only the master node/thread exists in mdrun.c,
1398 * therefore non-master nodes need to open the "seppot" log file here.
1400 || (!MASTER(cr) && (Flags & MD_SEPPOT))
1404 gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr, !(Flags & MD_SEPPOT),
1408 /* override nsteps with value from cmdline */
1409 override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
1413 copy_mat(state->box, box);
1418 gmx_bcast(sizeof(box), box, cr);
1421 /* Essential dynamics */
1422 if (opt2bSet("-ei", nfile, fnm))
1424 /* Open input and output files, allocate space for ED data structure */
1425 ed = ed_open(mtop->natoms, &state->edsamstate, nfile, fnm, Flags, oenv, cr);
1428 if (PAR(cr) && !((Flags & MD_PARTDEC) ||
1429 EI_TPI(inputrec->eI) ||
1430 inputrec->eI == eiNM))
1432 cr->dd = init_domain_decomposition(fplog, cr, Flags, ddxyz, rdd, rconstr,
1433 dddlb_opt, dlb_scale,
1434 ddcsx, ddcsy, ddcsz,
1437 &ddbox, &npme_major, &npme_minor);
1439 make_dd_communicators(fplog, cr, dd_node_order);
1441 /* Set overallocation to avoid frequent reallocation of arrays */
1442 set_over_alloc_dd(TRUE);
1446 /* PME, if used, is done on all nodes with 1D decomposition */
1448 cr->duty = (DUTY_PP | DUTY_PME);
1451 /* NM and TPI perform single node energy calculations in parallel */
1452 if (!(inputrec->eI == eiNM || EI_TPI(inputrec->eI)))
1454 npme_major = cr->nnodes;
1457 if (inputrec->ePBC == epbcSCREW)
1460 "pbc=%s is only implemented with domain decomposition",
1461 epbc_names[inputrec->ePBC]);
1467 /* After possible communicator splitting in make_dd_communicators.
1468 * we can set up the intra/inter node communication.
1470 gmx_setup_nodecomm(fplog, cr);
1473 /* Initialize per-physical-node MPI process/thread ID and counters. */
1474 gmx_init_intranode_counters(cr);
1477 md_print_info(cr, fplog, "Using %d MPI %s\n",
1479 #ifdef GMX_THREAD_MPI
1480 cr->nnodes == 1 ? "thread" : "threads"
1482 cr->nnodes == 1 ? "process" : "processes"
1488 /* Check and update hw_opt for the cut-off scheme */
1489 check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1491 gmx_omp_nthreads_init(fplog, cr,
1492 hwinfo->nthreads_hw_avail,
1493 hw_opt->nthreads_omp,
1494 hw_opt->nthreads_omp_pme,
1495 (cr->duty & DUTY_PP) == 0,
1496 inputrec->cutoff_scheme == ecutsVERLET);
1500 /* The master rank decided on the use of GPUs,
1501 * broadcast this information to all ranks.
1503 gmx_bcast_sim(sizeof(bUseGPU), &bUseGPU, cr);
1508 if (cr->npmenodes == -1)
1510 /* Don't automatically use PME-only nodes with GPUs */
1514 /* Select GPU id's to use */
1515 gmx_select_gpu_ids(fplog, cr, &hwinfo->gpu_info, bForceUseGPU,
1520 /* Ignore (potentially) manually selected GPUs */
1521 hw_opt->gpu_opt.ncuda_dev_use = 0;
1524 /* check consistency of CPU acceleration and number of GPUs selected */
1525 gmx_check_hw_runconf_consistency(fplog, hwinfo, cr, hw_opt, bUseGPU);
1527 if (DOMAINDECOMP(cr))
1529 /* When we share GPUs over ranks, we need to know this for the DLB */
1530 dd_setup_dlb_resource_sharing(cr, hwinfo, hw_opt);
1533 /* getting number of PP/PME threads
1534 PME: env variable should be read only on one node to make sure it is
1535 identical everywhere;
1537 /* TODO nthreads_pp is only used for pinning threads.
1538 * This is a temporary solution until we have a hw topology library.
1540 nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
1541 nthreads_pme = gmx_omp_nthreads_get(emntPME);
1543 wcycle = wallcycle_init(fplog, resetstep, cr, nthreads_pp, nthreads_pme);
1547 /* Master synchronizes its value of reset_counters with all nodes
1548 * including PME only nodes */
1549 reset_counters = wcycle_get_reset_counters(wcycle);
1550 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1551 wcycle_set_reset_counters(wcycle, reset_counters);
1555 if (cr->duty & DUTY_PP)
1557 /* For domain decomposition we allocate dynamically
1558 * in dd_partition_system.
1560 if (DOMAINDECOMP(cr))
1562 bcast_state_setup(cr, state);
1568 bcast_state(cr, state, TRUE);
1572 /* Initiate forcerecord */
1574 fr->hwinfo = hwinfo;
1575 fr->gpu_opt = &hw_opt->gpu_opt;
1576 init_forcerec(fplog, oenv, fr, fcd, inputrec, mtop, cr, box,
1577 opt2fn("-table", nfile, fnm),
1578 opt2fn("-tabletf", nfile, fnm),
1579 opt2fn("-tablep", nfile, fnm),
1580 opt2fn("-tableb", nfile, fnm),
1585 /* version for PCA_NOT_READ_NODE (see md.c) */
1586 /*init_forcerec(fplog,fr,fcd,inputrec,mtop,cr,box,FALSE,
1587 "nofile","nofile","nofile","nofile",FALSE,pforce);
1589 fr->bSepDVDL = ((Flags & MD_SEPPOT) == MD_SEPPOT);
1591 /* Initialize QM-MM */
1594 init_QMMMrec(cr, mtop, inputrec, fr);
1597 /* Initialize the mdatoms structure.
1598 * mdatoms is not filled with atom data,
1599 * as this can not be done now with domain decomposition.
1601 mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1603 if (mdatoms->nPerturbed > 0 && inputrec->cutoff_scheme == ecutsVERLET)
1605 gmx_fatal(FARGS, "The Verlet cut-off scheme does not (yet) support free-energy calculations with perturbed atoms, only perturbed interactions. This will be implemented soon. Use the group scheme for now.");
1608 /* Initialize the virtual site communication */
1609 vsite = init_vsite(mtop, cr, FALSE);
1611 calc_shifts(box, fr->shift_vec);
1613 /* With periodic molecules the charge groups should be whole at start up
1614 * and the virtual sites should not be far from their proper positions.
1616 if (!inputrec->bContinuation && MASTER(cr) &&
1617 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1619 /* Make molecules whole at start of run */
1620 if (fr->ePBC != epbcNONE)
1622 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, state->x);
1626 /* Correct initial vsite positions are required
1627 * for the initial distribution in the domain decomposition
1628 * and for the initial shell prediction.
1630 construct_vsites_mtop(vsite, mtop, state->x);
1634 if (EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype))
1636 ewaldcoeff_q = fr->ewaldcoeff_q;
1637 ewaldcoeff_lj = fr->ewaldcoeff_lj;
1638 pmedata = &fr->pmedata;
1647 /* This is a PME only node */
1649 /* We don't need the state */
1652 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1653 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1657 if (hw_opt->thread_affinity != threadaffOFF)
1659 /* Before setting affinity, check whether the affinity has changed
1660 * - which indicates that probably the OpenMP library has changed it
1661 * since we first checked).
1663 gmx_check_thread_affinity_set(fplog, cr,
1664 hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1666 /* Set the CPU affinity */
1667 gmx_set_thread_affinity(fplog, cr, hw_opt, hwinfo);
1670 /* Initiate PME if necessary,
1671 * either on all nodes or on dedicated PME nodes only. */
1672 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1676 nChargePerturbed = mdatoms->nChargePerturbed;
1677 if (EVDW_PME(inputrec->vdwtype))
1679 nTypePerturbed = mdatoms->nTypePerturbed;
1682 if (cr->npmenodes > 0)
1684 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1685 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1686 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1689 if (cr->duty & DUTY_PME)
1691 status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1692 mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1693 (Flags & MD_REPRODUCIBLE), nthreads_pme);
1696 gmx_fatal(FARGS, "Error %d initializing PME", status);
1702 if (integrator[inputrec->eI].func == do_md)
1704 /* Turn on signal handling on all nodes */
1706 * (A user signal from the PME nodes (if any)
1707 * is communicated to the PP nodes.
1709 signal_handler_install();
1712 if (cr->duty & DUTY_PP)
1714 /* Assumes uniform use of the number of OpenMP threads */
1715 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1717 if (inputrec->ePull != epullNO)
1719 /* Initialize pull code */
1720 init_pull(fplog, inputrec, nfile, fnm, mtop, cr, oenv, inputrec->fepvals->init_lambda,
1721 EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1726 /* Initialize enforced rotation code */
1727 init_rot(fplog, inputrec, nfile, fnm, cr, state->x, box, mtop, oenv,
1731 if (inputrec->eSwapCoords != eswapNO)
1733 /* Initialize ion swapping code */
1734 init_swapcoords(fplog, bVerbose, inputrec, opt2fn_master("-swap", nfile, fnm, cr),
1735 mtop, state->x, state->box, &state->swapstate, cr, oenv, Flags);
1738 constr = init_constraints(fplog, mtop, inputrec, ed, state, cr);
1740 if (DOMAINDECOMP(cr))
1742 dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1743 Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1745 set_dd_parameters(fplog, cr->dd, dlb_scale, inputrec, &ddbox);
1747 setup_dd_grid(fplog, cr->dd);
1750 /* Now do whatever the user wants us to do (how flexible...) */
1751 integrator[inputrec->eI].func(fplog, cr, nfile, fnm,
1752 oenv, bVerbose, bCompact,
1755 nstepout, inputrec, mtop,
1757 mdatoms, nrnb, wcycle, ed, fr,
1758 repl_ex_nst, repl_ex_nex, repl_ex_seed,
1760 cpt_period, max_hours,
1763 walltime_accounting);
1765 if (inputrec->ePull != epullNO)
1767 finish_pull(inputrec->pull);
1772 finish_rot(inputrec->rot);
1779 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1780 gmx_pmeonly(*pmedata, cr, nrnb, wcycle, walltime_accounting, ewaldcoeff_q, ewaldcoeff_lj, inputrec);
1783 wallcycle_stop(wcycle, ewcRUN);
1785 /* Finish up, write some stuff
1786 * if rerunMD, don't write last frame again
1788 finish_run(fplog, cr,
1789 inputrec, nrnb, wcycle, walltime_accounting,
1790 fr != NULL && fr->nbv != NULL && fr->nbv->bUseGPU ?
1791 nbnxn_cuda_get_timings(fr->nbv->cu_nbv) : NULL,
1792 EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1795 /* Free GPU memory and context */
1796 free_gpu_resources(fr, cr);
1798 if (opt2bSet("-membed", nfile, fnm))
1803 gmx_hardware_info_free(hwinfo);
1805 /* Does what it says */
1806 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", walltime_accounting);
1807 walltime_accounting_destroy(walltime_accounting);
1809 /* Close logfile already here if we were appending to it */
1810 if (MASTER(cr) && (Flags & MD_APPENDFILES))
1812 gmx_log_close(fplog);
1815 rc = (int)gmx_get_stop_condition();
1817 #ifdef GMX_THREAD_MPI
1818 /* we need to join all threads. The sub-threads join when they
1819 exit this function, but the master thread needs to be told to
1821 if (PAR(cr) && MASTER(cr))