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 "gromacs/topology/mtop_util.h"
71 #include "sighandler.h"
73 #include "gmx_detect_hardware.h"
74 #include "gmx_omp_nthreads.h"
75 #include "gromacs/gmxpreprocess/calc_verletbuf.h"
77 #include "gmx_thread_affinity.h"
81 #include "gromacs/essentialdynamics/edsam.h"
82 #include "gromacs/fileio/tpxio.h"
83 #include "gromacs/math/vec.h"
84 #include "gromacs/mdlib/nbnxn_search.h"
85 #include "gromacs/mdlib/nbnxn_consts.h"
86 #include "gromacs/pbcutil/pbc.h"
87 #include "gromacs/pulling/pull.h"
88 #include "gromacs/pulling/pull_rotation.h"
89 #include "gromacs/swap/swapcoords.h"
90 #include "gromacs/timing/wallcycle.h"
91 #include "gromacs/utility/gmxassert.h"
92 #include "gromacs/utility/gmxmpi.h"
93 #include "gromacs/utility/smalloc.h"
99 #include "gpu_utils.h"
100 #include "nbnxn_cuda_data_mgmt.h"
103 gmx_integrator_t *func;
106 /* The array should match the eI array in include/types/enums.h */
107 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}};
109 gmx_int64_t deform_init_init_step_tpx;
110 matrix deform_init_box_tpx;
111 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
114 #ifdef GMX_THREAD_MPI
115 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
116 * the number of threads will get lowered.
118 #define MIN_ATOMS_PER_MPI_THREAD 90
119 #define MIN_ATOMS_PER_GPU 900
121 struct mdrunner_arglist
136 const char *dddlb_opt;
141 const char *nbpu_opt;
143 gmx_int64_t nsteps_cmdline;
153 const char *deviceOptions;
159 /* The function used for spawning threads. Extracts the mdrunner()
160 arguments from its one argument and calls mdrunner(), after making
162 static void mdrunner_start_fn(void *arg)
164 struct mdrunner_arglist *mda = (struct mdrunner_arglist*)arg;
165 struct mdrunner_arglist mc = *mda; /* copy the arg list to make sure
166 that it's thread-local. This doesn't
167 copy pointed-to items, of course,
168 but those are all const. */
169 t_commrec *cr; /* we need a local version of this */
173 fnm = dup_tfn(mc.nfile, mc.fnm);
175 cr = reinitialize_commrec_for_this_thread(mc.cr);
182 mdrunner(&mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
183 mc.bVerbose, mc.bCompact, mc.nstglobalcomm,
184 mc.ddxyz, mc.dd_node_order, mc.rdd,
185 mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
186 mc.ddcsx, mc.ddcsy, mc.ddcsz,
187 mc.nbpu_opt, mc.nstlist_cmdline,
188 mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
189 mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce,
190 mc.cpt_period, mc.max_hours, mc.deviceOptions, mc.imdport, mc.Flags);
193 /* called by mdrunner() to start a specific number of threads (including
194 the main thread) for thread-parallel runs. This in turn calls mdrunner()
196 All options besides nthreads are the same as for mdrunner(). */
197 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
198 FILE *fplog, t_commrec *cr, int nfile,
199 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
200 gmx_bool bCompact, int nstglobalcomm,
201 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
202 const char *dddlb_opt, real dlb_scale,
203 const char *ddcsx, const char *ddcsy, const char *ddcsz,
204 const char *nbpu_opt, int nstlist_cmdline,
205 gmx_int64_t nsteps_cmdline,
206 int nstepout, int resetstep,
207 int nmultisim, int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
208 real pforce, real cpt_period, real max_hours,
209 const char *deviceOptions, unsigned long Flags)
212 struct mdrunner_arglist *mda;
213 t_commrec *crn; /* the new commrec */
216 /* first check whether we even need to start tMPI */
217 if (hw_opt->nthreads_tmpi < 2)
222 /* a few small, one-time, almost unavoidable memory leaks: */
224 fnmn = dup_tfn(nfile, fnm);
226 /* fill the data structure to pass as void pointer to thread start fn */
227 /* hw_opt contains pointers, which should all be NULL at this stage */
228 mda->hw_opt = *hw_opt;
234 mda->bVerbose = bVerbose;
235 mda->bCompact = bCompact;
236 mda->nstglobalcomm = nstglobalcomm;
237 mda->ddxyz[XX] = ddxyz[XX];
238 mda->ddxyz[YY] = ddxyz[YY];
239 mda->ddxyz[ZZ] = ddxyz[ZZ];
240 mda->dd_node_order = dd_node_order;
242 mda->rconstr = rconstr;
243 mda->dddlb_opt = dddlb_opt;
244 mda->dlb_scale = dlb_scale;
248 mda->nbpu_opt = nbpu_opt;
249 mda->nstlist_cmdline = nstlist_cmdline;
250 mda->nsteps_cmdline = nsteps_cmdline;
251 mda->nstepout = nstepout;
252 mda->resetstep = resetstep;
253 mda->nmultisim = nmultisim;
254 mda->repl_ex_nst = repl_ex_nst;
255 mda->repl_ex_nex = repl_ex_nex;
256 mda->repl_ex_seed = repl_ex_seed;
257 mda->pforce = pforce;
258 mda->cpt_period = cpt_period;
259 mda->max_hours = max_hours;
260 mda->deviceOptions = deviceOptions;
263 /* now spawn new threads that start mdrunner_start_fn(), while
264 the main thread returns, we set thread affinity later */
265 ret = tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi, TMPI_AFFINITY_NONE,
266 mdrunner_start_fn, (void*)(mda) );
267 if (ret != TMPI_SUCCESS)
272 crn = reinitialize_commrec_for_this_thread(cr);
277 static int get_tmpi_omp_thread_division(const gmx_hw_info_t *hwinfo,
278 const gmx_hw_opt_t *hw_opt,
284 /* There are no separate PME nodes here, as we ensured in
285 * check_and_update_hw_opt that nthreads_tmpi>0 with PME nodes
286 * and a conditional ensures we would not have ended up here.
287 * Note that separate PME nodes might be switched on later.
291 nthreads_tmpi = ngpu;
292 if (nthreads_tot > 0 && nthreads_tot < nthreads_tmpi)
294 nthreads_tmpi = nthreads_tot;
297 else if (hw_opt->nthreads_omp > 0)
299 /* Here we could oversubscribe, when we do, we issue a warning later */
300 nthreads_tmpi = std::max(1, nthreads_tot/hw_opt->nthreads_omp);
304 /* TODO choose nthreads_omp based on hardware topology
305 when we have a hardware topology detection library */
306 /* In general, when running up to 4 threads, OpenMP should be faster.
307 * Note: on AMD Bulldozer we should avoid running OpenMP over two dies.
308 * On Intel>=Nehalem running OpenMP on a single CPU is always faster,
309 * even on two CPUs it's usually faster (but with many OpenMP threads
310 * it could be faster not to use HT, currently we always use HT).
311 * On Nehalem/Westmere we want to avoid running 16 threads over
312 * two CPUs with HT, so we need a limit<16; thus we use 12.
313 * A reasonable limit for Intel Sandy and Ivy bridge,
314 * not knowing the topology, is 16 threads.
316 const int nthreads_omp_always_faster = 4;
317 const int nthreads_omp_always_faster_Nehalem = 12;
318 const int nthreads_omp_always_faster_SandyBridge = 16;
319 gmx_bool bIntel_Family6;
322 (gmx_cpuid_vendor(hwinfo->cpuid_info) == GMX_CPUID_VENDOR_INTEL &&
323 gmx_cpuid_family(hwinfo->cpuid_info) == 6);
325 if (nthreads_tot <= nthreads_omp_always_faster ||
327 ((gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_Nehalem && nthreads_tot <= nthreads_omp_always_faster_Nehalem) ||
328 (gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_SandyBridge && nthreads_tot <= nthreads_omp_always_faster_SandyBridge))))
330 /* Use pure OpenMP parallelization */
335 /* Don't use OpenMP parallelization */
336 nthreads_tmpi = nthreads_tot;
340 return nthreads_tmpi;
344 /* Get the number of threads to use for thread-MPI based on how many
345 * were requested, which algorithms we're using,
346 * and how many particles there are.
347 * At the point we have already called check_and_update_hw_opt.
348 * Thus all options should be internally consistent and consistent
349 * with the hardware, except that ntmpi could be larger than #GPU.
351 static int get_nthreads_mpi(const gmx_hw_info_t *hwinfo,
352 gmx_hw_opt_t *hw_opt,
353 t_inputrec *inputrec, gmx_mtop_t *mtop,
357 int nthreads_hw, nthreads_tot_max, nthreads_tmpi, nthreads_new, ngpu;
358 int min_atoms_per_mpi_thread;
361 if (hw_opt->nthreads_tmpi > 0)
363 /* Trivial, return right away */
364 return hw_opt->nthreads_tmpi;
367 nthreads_hw = hwinfo->nthreads_hw_avail;
369 /* How many total (#tMPI*#OpenMP) threads can we start? */
370 if (hw_opt->nthreads_tot > 0)
372 nthreads_tot_max = hw_opt->nthreads_tot;
376 nthreads_tot_max = nthreads_hw;
379 bCanUseGPU = (inputrec->cutoff_scheme == ecutsVERLET &&
380 hwinfo->gpu_info.ncuda_dev_compatible > 0);
383 ngpu = hwinfo->gpu_info.ncuda_dev_compatible;
390 if (inputrec->cutoff_scheme == ecutsGROUP)
392 /* We checked this before, but it doesn't hurt to do it once more */
393 assert(hw_opt->nthreads_omp == 1);
397 get_tmpi_omp_thread_division(hwinfo, hw_opt, nthreads_tot_max, ngpu);
399 if (inputrec->eI == eiNM || EI_TPI(inputrec->eI))
401 /* Dims/steps are divided over the nodes iso splitting the atoms */
402 min_atoms_per_mpi_thread = 0;
408 min_atoms_per_mpi_thread = MIN_ATOMS_PER_GPU;
412 min_atoms_per_mpi_thread = MIN_ATOMS_PER_MPI_THREAD;
416 /* Check if an algorithm does not support parallel simulation. */
417 if (nthreads_tmpi != 1 &&
418 ( inputrec->eI == eiLBFGS ||
419 inputrec->coulombtype == eelEWALD ) )
423 md_print_warn(cr, fplog, "The integration or electrostatics algorithm doesn't support parallel runs. Using a single thread-MPI thread.\n");
424 if (hw_opt->nthreads_tmpi > nthreads_tmpi)
426 gmx_fatal(FARGS, "You asked for more than 1 thread-MPI thread, but an algorithm doesn't support that");
429 else if (mtop->natoms/nthreads_tmpi < min_atoms_per_mpi_thread)
431 /* the thread number was chosen automatically, but there are too many
432 threads (too few atoms per thread) */
433 nthreads_new = std::max(1, mtop->natoms/min_atoms_per_mpi_thread);
435 /* Avoid partial use of Hyper-Threading */
436 if (gmx_cpuid_x86_smt(hwinfo->cpuid_info) == GMX_CPUID_X86_SMT_ENABLED &&
437 nthreads_new > nthreads_hw/2 && nthreads_new < nthreads_hw)
439 nthreads_new = nthreads_hw/2;
442 /* Avoid large prime numbers in the thread count */
443 if (nthreads_new >= 6)
445 /* Use only 6,8,10 with additional factors of 2 */
449 while (3*fac*2 <= nthreads_new)
454 nthreads_new = (nthreads_new/fac)*fac;
459 if (nthreads_new == 5)
465 nthreads_tmpi = nthreads_new;
467 fprintf(stderr, "\n");
468 fprintf(stderr, "NOTE: Parallelization is limited by the small number of atoms,\n");
469 fprintf(stderr, " only starting %d thread-MPI threads.\n", nthreads_tmpi);
470 fprintf(stderr, " You can use the -nt and/or -ntmpi option to optimize the number of threads.\n\n");
473 return nthreads_tmpi;
475 #endif /* GMX_THREAD_MPI */
478 /* We determine the extra cost of the non-bonded kernels compared to
479 * a reference nstlist value of 10 (which is the default in grompp).
481 static const int nbnxnReferenceNstlist = 10;
482 /* The values to try when switching */
483 const int nstlist_try[] = { 20, 25, 40 };
484 #define NNSTL sizeof(nstlist_try)/sizeof(nstlist_try[0])
485 /* Increase nstlist until the non-bonded cost increases more than listfac_ok,
486 * but never more than listfac_max.
487 * A standard (protein+)water system at 300K with PME ewald_rtol=1e-5
488 * needs 1.28 at rcoulomb=0.9 and 1.24 at rcoulomb=1.0 to get to nstlist=40.
489 * Note that both CPU and GPU factors are conservative. Performance should
490 * not go down due to this tuning, except with a relatively slow GPU.
491 * On the other hand, at medium/high parallelization or with fast GPUs
492 * nstlist will not be increased enough to reach optimal performance.
494 /* CPU: pair-search is about a factor 1.5 slower than the non-bonded kernel */
495 static const float nbnxn_cpu_listfac_ok = 1.05;
496 static const float nbnxn_cpu_listfac_max = 1.09;
497 /* GPU: pair-search is a factor 1.5-3 slower than the non-bonded kernel */
498 static const float nbnxn_gpu_listfac_ok = 1.20;
499 static const float nbnxn_gpu_listfac_max = 1.30;
501 /* Try to increase nstlist when using the Verlet cut-off scheme */
502 static void increase_nstlist(FILE *fp, t_commrec *cr,
503 t_inputrec *ir, int nstlist_cmdline,
504 const gmx_mtop_t *mtop, matrix box,
507 float listfac_ok, listfac_max;
508 int nstlist_orig, nstlist_prev;
509 verletbuf_list_setup_t ls;
510 real rlistWithReferenceNstlist, rlist_inc, rlist_ok, rlist_max;
511 real rlist_new, rlist_prev;
512 size_t nstlist_ind = 0;
514 gmx_bool bBox, bDD, bCont;
515 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";
516 const char *nve_err = "Can not increase nstlist because an NVE ensemble is used";
517 const char *vbd_err = "Can not increase nstlist because verlet-buffer-tolerance is not set or used";
518 const char *box_err = "Can not increase nstlist because the box is too small";
519 const char *dd_err = "Can not increase nstlist because of domain decomposition limitations";
521 const float oneThird = 1.0f / 3.0f;
523 if (nstlist_cmdline <= 0)
525 if (ir->nstlist == 1)
527 /* The user probably set nstlist=1 for a reason,
528 * don't mess with the settings.
533 if (fp != NULL && bGPU && ir->nstlist < nstlist_try[0])
535 fprintf(fp, nstl_gpu, ir->nstlist);
538 while (nstlist_ind < NNSTL && ir->nstlist >= nstlist_try[nstlist_ind])
542 if (nstlist_ind == NNSTL)
544 /* There are no larger nstlist value to try */
549 if (EI_MD(ir->eI) && ir->etc == etcNO)
553 fprintf(stderr, "%s\n", nve_err);
557 fprintf(fp, "%s\n", nve_err);
563 if (ir->verletbuf_tol == 0 && bGPU)
565 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");
568 if (ir->verletbuf_tol < 0)
572 fprintf(stderr, "%s\n", vbd_err);
576 fprintf(fp, "%s\n", vbd_err);
584 listfac_ok = nbnxn_gpu_listfac_ok;
585 listfac_max = nbnxn_gpu_listfac_max;
589 listfac_ok = nbnxn_cpu_listfac_ok;
590 listfac_max = nbnxn_cpu_listfac_max;
593 nstlist_orig = ir->nstlist;
594 if (nstlist_cmdline > 0)
598 sprintf(buf, "Getting nstlist=%d from command line option",
601 ir->nstlist = nstlist_cmdline;
604 verletbuf_get_list_setup(bGPU, &ls);
606 /* Allow rlist to make the list a given factor larger than the list
607 * would be with the reference value for nstlist (10).
609 nstlist_prev = ir->nstlist;
610 ir->nstlist = nbnxnReferenceNstlist;
611 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL,
612 &rlistWithReferenceNstlist);
613 ir->nstlist = nstlist_prev;
615 /* Determine the pair list size increase due to zero interactions */
616 rlist_inc = nbnxn_get_rlist_effective_inc(ls.cluster_size_j,
617 mtop->natoms/det(box));
618 rlist_ok = (rlistWithReferenceNstlist + rlist_inc)*pow(listfac_ok, oneThird) - rlist_inc;
619 rlist_max = (rlistWithReferenceNstlist + rlist_inc)*pow(listfac_max, oneThird) - rlist_inc;
622 fprintf(debug, "nstlist tuning: rlist_inc %.3f rlist_ok %.3f rlist_max %.3f\n",
623 rlist_inc, rlist_ok, rlist_max);
626 nstlist_prev = nstlist_orig;
627 rlist_prev = ir->rlist;
630 if (nstlist_cmdline <= 0)
632 ir->nstlist = nstlist_try[nstlist_ind];
635 /* Set the pair-list buffer size in ir */
636 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
638 /* Does rlist fit in the box? */
639 bBox = (sqr(rlist_new) < max_cutoff2(ir->ePBC, box));
641 if (bBox && DOMAINDECOMP(cr))
643 /* Check if rlist fits in the domain decomposition */
644 if (inputrec2nboundeddim(ir) < DIM)
646 gmx_incons("Changing nstlist with domain decomposition and unbounded dimensions is not implemented yet");
648 copy_mat(box, state_tmp.box);
649 bDD = change_dd_cutoff(cr, &state_tmp, ir, rlist_new);
654 fprintf(debug, "nstlist %d rlist %.3f bBox %d bDD %d\n",
655 ir->nstlist, rlist_new, bBox, bDD);
660 if (nstlist_cmdline <= 0)
662 if (bBox && bDD && rlist_new <= rlist_max)
664 /* Increase nstlist */
665 nstlist_prev = ir->nstlist;
666 rlist_prev = rlist_new;
667 bCont = (nstlist_ind+1 < NNSTL && rlist_new < rlist_ok);
671 /* Stick with the previous nstlist */
672 ir->nstlist = nstlist_prev;
673 rlist_new = rlist_prev;
685 gmx_warning(!bBox ? box_err : dd_err);
688 fprintf(fp, "\n%s\n", bBox ? box_err : dd_err);
690 ir->nstlist = nstlist_orig;
692 else if (ir->nstlist != nstlist_orig || rlist_new != ir->rlist)
694 sprintf(buf, "Changing nstlist from %d to %d, rlist from %g to %g",
695 nstlist_orig, ir->nstlist,
696 ir->rlist, rlist_new);
699 fprintf(stderr, "%s\n\n", buf);
703 fprintf(fp, "%s\n\n", buf);
705 ir->rlist = rlist_new;
706 ir->rlistlong = rlist_new;
710 static void prepare_verlet_scheme(FILE *fplog,
714 const gmx_mtop_t *mtop,
718 /* For NVE simulations, we will retain the initial list buffer */
719 if (ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
721 /* Update the Verlet buffer size for the current run setup */
722 verletbuf_list_setup_t ls;
725 /* Here we assume SIMD-enabled kernels are being used. But as currently
726 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
727 * and 4x2 gives a larger buffer than 4x4, this is ok.
729 verletbuf_get_list_setup(bUseGPU, &ls);
731 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
733 if (rlist_new != ir->rlist)
737 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
738 ir->rlist, rlist_new,
739 ls.cluster_size_i, ls.cluster_size_j);
741 ir->rlist = rlist_new;
742 ir->rlistlong = rlist_new;
746 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
748 gmx_fatal(FARGS, "Can not set nstlist without %s",
749 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
752 if (EI_DYNAMICS(ir->eI))
754 /* Set or try nstlist values */
755 increase_nstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, bUseGPU);
759 static void convert_to_verlet_scheme(FILE *fplog,
761 gmx_mtop_t *mtop, real box_vol)
763 const char *conv_mesg = "Converting input file with group cut-off scheme to the Verlet cut-off scheme";
765 md_print_warn(NULL, fplog, "%s\n", conv_mesg);
767 ir->cutoff_scheme = ecutsVERLET;
768 ir->verletbuf_tol = 0.005;
770 if (ir->rcoulomb != ir->rvdw)
772 gmx_fatal(FARGS, "The VdW and Coulomb cut-offs are different, whereas the Verlet scheme only supports equal cut-offs");
775 if (ir->vdwtype == evdwUSER || EEL_USER(ir->coulombtype))
777 gmx_fatal(FARGS, "User non-bonded potentials are not (yet) supported with the Verlet scheme");
779 else if (ir_vdw_switched(ir) || ir_coulomb_switched(ir))
781 if (ir_vdw_switched(ir) && ir->vdw_modifier == eintmodNONE)
783 ir->vdwtype = evdwCUT;
787 case evdwSHIFT: ir->vdw_modifier = eintmodFORCESWITCH; break;
788 case evdwSWITCH: ir->vdw_modifier = eintmodPOTSWITCH; break;
789 default: gmx_fatal(FARGS, "The Verlet scheme does not support Van der Waals interactions of type '%s'", evdw_names[ir->vdwtype]);
792 if (ir_coulomb_switched(ir) && ir->coulomb_modifier == eintmodNONE)
794 if (EEL_FULL(ir->coulombtype))
796 /* With full electrostatic only PME can be switched */
797 ir->coulombtype = eelPME;
798 ir->coulomb_modifier = eintmodPOTSHIFT;
802 md_print_warn(NULL, fplog, "NOTE: Replacing %s electrostatics with reaction-field with epsilon-rf=inf\n", eel_names[ir->coulombtype]);
803 ir->coulombtype = eelRF;
804 ir->epsilon_rf = 0.0;
805 ir->coulomb_modifier = eintmodPOTSHIFT;
809 /* We set the pair energy error tolerance to a small number.
810 * Note that this is only for testing. For production the user
811 * should think about this and set the mdp options.
813 ir->verletbuf_tol = 1e-4;
816 if (inputrec2nboundeddim(ir) != 3)
818 gmx_fatal(FARGS, "Can only convert old tpr files to the Verlet cut-off scheme with 3D pbc");
821 if (ir->efep != efepNO || ir->implicit_solvent != eisNO)
823 gmx_fatal(FARGS, "Will not convert old tpr files to the Verlet cut-off scheme with free-energy calculations or implicit solvent");
826 if (EI_DYNAMICS(ir->eI) && !(EI_MD(ir->eI) && ir->etc == etcNO))
828 verletbuf_list_setup_t ls;
830 verletbuf_get_list_setup(FALSE, &ls);
831 calc_verlet_buffer_size(mtop, box_vol, ir, -1, &ls, NULL, &ir->rlist);
839 rlist_fac = 1 + verlet_buffer_ratio_NVE_T0;
843 rlist_fac = 1 + verlet_buffer_ratio_nodynamics;
845 ir->verletbuf_tol = -1;
846 ir->rlist = rlist_fac*std::max(ir->rvdw, ir->rcoulomb);
849 gmx_mtop_remove_chargegroups(mtop);
852 static void print_hw_opt(FILE *fp, const gmx_hw_opt_t *hw_opt)
854 fprintf(fp, "hw_opt: nt %d ntmpi %d ntomp %d ntomp_pme %d gpu_id '%s'\n",
855 hw_opt->nthreads_tot,
856 hw_opt->nthreads_tmpi,
857 hw_opt->nthreads_omp,
858 hw_opt->nthreads_omp_pme,
859 hw_opt->gpu_opt.gpu_id != NULL ? hw_opt->gpu_opt.gpu_id : "");
862 /* Checks we can do when we don't (yet) know the cut-off scheme */
863 static void check_and_update_hw_opt_1(gmx_hw_opt_t *hw_opt,
864 gmx_bool bIsSimMaster)
866 gmx_omp_nthreads_read_env(&hw_opt->nthreads_omp, bIsSimMaster);
868 #ifndef GMX_THREAD_MPI
869 if (hw_opt->nthreads_tot > 0)
871 gmx_fatal(FARGS, "Setting the total number of threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
873 if (hw_opt->nthreads_tmpi > 0)
875 gmx_fatal(FARGS, "Setting the number of thread-MPI threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
880 if (hw_opt->nthreads_omp > 1)
882 gmx_fatal(FARGS, "More than 1 OpenMP thread requested, but Gromacs was compiled without OpenMP support");
884 hw_opt->nthreads_omp = 1;
887 if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp_pme <= 0)
889 /* We have the same number of OpenMP threads for PP and PME processes,
890 * thus we can perform several consistency checks.
892 if (hw_opt->nthreads_tmpi > 0 &&
893 hw_opt->nthreads_omp > 0 &&
894 hw_opt->nthreads_tot != hw_opt->nthreads_tmpi*hw_opt->nthreads_omp)
896 gmx_fatal(FARGS, "The total number of threads requested (%d) does not match the thread-MPI threads (%d) times the OpenMP threads (%d) requested",
897 hw_opt->nthreads_tot, hw_opt->nthreads_tmpi, hw_opt->nthreads_omp);
900 if (hw_opt->nthreads_tmpi > 0 &&
901 hw_opt->nthreads_tot % hw_opt->nthreads_tmpi != 0)
903 gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of thread-MPI threads requested (%d)",
904 hw_opt->nthreads_tot, hw_opt->nthreads_tmpi);
907 if (hw_opt->nthreads_omp > 0 &&
908 hw_opt->nthreads_tot % hw_opt->nthreads_omp != 0)
910 gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of OpenMP threads requested (%d)",
911 hw_opt->nthreads_tot, hw_opt->nthreads_omp);
914 if (hw_opt->nthreads_tmpi > 0 &&
915 hw_opt->nthreads_omp <= 0)
917 hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
922 if (hw_opt->nthreads_omp > 1)
924 gmx_fatal(FARGS, "OpenMP threads are requested, but Gromacs was compiled without OpenMP support");
928 if (hw_opt->nthreads_omp_pme > 0 && hw_opt->nthreads_omp <= 0)
930 gmx_fatal(FARGS, "You need to specify -ntomp in addition to -ntomp_pme");
933 if (hw_opt->nthreads_tot == 1)
935 hw_opt->nthreads_tmpi = 1;
937 if (hw_opt->nthreads_omp > 1)
939 gmx_fatal(FARGS, "You requested %d OpenMP threads with %d total threads",
940 hw_opt->nthreads_tmpi, hw_opt->nthreads_tot);
942 hw_opt->nthreads_omp = 1;
945 if (hw_opt->nthreads_omp_pme <= 0 && hw_opt->nthreads_omp > 0)
947 hw_opt->nthreads_omp_pme = hw_opt->nthreads_omp;
950 /* Parse GPU IDs, if provided.
951 * We check consistency with the tMPI thread count later.
953 gmx_parse_gpu_ids(&hw_opt->gpu_opt);
955 #ifdef GMX_THREAD_MPI
956 if (hw_opt->gpu_opt.ncuda_dev_use > 0 && hw_opt->nthreads_tmpi == 0)
958 /* Set the number of MPI threads equal to the number of GPUs */
959 hw_opt->nthreads_tmpi = hw_opt->gpu_opt.ncuda_dev_use;
961 if (hw_opt->nthreads_tot > 0 &&
962 hw_opt->nthreads_tmpi > hw_opt->nthreads_tot)
964 /* We have more GPUs than total threads requested.
965 * We choose to (later) generate a mismatch error,
966 * instead of launching more threads than requested.
968 hw_opt->nthreads_tmpi = hw_opt->nthreads_tot;
975 print_hw_opt(debug, hw_opt);
979 /* Checks we can do when we know the cut-off scheme */
980 static void check_and_update_hw_opt_2(gmx_hw_opt_t *hw_opt,
983 if (cutoff_scheme == ecutsGROUP)
985 /* We only have OpenMP support for PME only nodes */
986 if (hw_opt->nthreads_omp > 1)
988 gmx_fatal(FARGS, "OpenMP threads have been requested with cut-off scheme %s, but these are only supported with cut-off scheme %s",
989 ecutscheme_names[cutoff_scheme],
990 ecutscheme_names[ecutsVERLET]);
992 hw_opt->nthreads_omp = 1;
995 if (hw_opt->nthreads_omp_pme <= 0 && hw_opt->nthreads_omp > 0)
997 hw_opt->nthreads_omp_pme = hw_opt->nthreads_omp;
1002 print_hw_opt(debug, hw_opt);
1007 /* Override the value in inputrec with value passed on the command line (if any) */
1008 static void override_nsteps_cmdline(FILE *fplog,
1009 gmx_int64_t nsteps_cmdline,
1011 const t_commrec *cr)
1013 char sbuf[STEPSTRSIZE];
1018 /* override with anything else than the default -2 */
1019 if (nsteps_cmdline > -2)
1023 ir->nsteps = nsteps_cmdline;
1024 if (EI_DYNAMICS(ir->eI))
1026 sprintf(stmp, "Overriding nsteps with value passed on the command line: %s steps, %.3f ps",
1027 gmx_step_str(nsteps_cmdline, sbuf),
1028 nsteps_cmdline*ir->delta_t);
1032 sprintf(stmp, "Overriding nsteps with value passed on the command line: %s steps",
1033 gmx_step_str(nsteps_cmdline, sbuf));
1036 md_print_warn(cr, fplog, "%s\n", stmp);
1040 /* Frees GPU memory and destroys the CUDA context.
1042 * Note that this function needs to be called even if GPUs are not used
1043 * in this run because the PME ranks have no knowledge of whether GPUs
1044 * are used or not, but all ranks need to enter the barrier below.
1046 static void free_gpu_resources(const t_forcerec *fr,
1047 const t_commrec *cr)
1049 gmx_bool bIsPPrankUsingGPU;
1050 char gpu_err_str[STRLEN];
1052 bIsPPrankUsingGPU = (cr->duty & DUTY_PP) && fr != NULL && fr->nbv != NULL && fr->nbv->bUseGPU;
1054 if (bIsPPrankUsingGPU)
1056 /* free nbnxn data in GPU memory */
1057 nbnxn_cuda_free(fr->nbv->cu_nbv);
1059 /* With tMPI we need to wait for all ranks to finish deallocation before
1060 * destroying the context in free_gpu() as some ranks may be sharing
1062 * Note: as only PP ranks need to free GPU resources, so it is safe to
1063 * not call the barrier on PME ranks.
1065 #ifdef GMX_THREAD_MPI
1070 #endif /* GMX_THREAD_MPI */
1072 /* uninitialize GPU (by destroying the context) */
1073 if (!free_gpu(gpu_err_str))
1075 gmx_warning("On rank %d failed to free GPU #%d: %s",
1076 cr->nodeid, get_current_gpu_device_id(), gpu_err_str);
1081 int mdrunner(gmx_hw_opt_t *hw_opt,
1082 FILE *fplog, t_commrec *cr, int nfile,
1083 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
1084 gmx_bool bCompact, int nstglobalcomm,
1085 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
1086 const char *dddlb_opt, real dlb_scale,
1087 const char *ddcsx, const char *ddcsy, const char *ddcsz,
1088 const char *nbpu_opt, int nstlist_cmdline,
1089 gmx_int64_t nsteps_cmdline, int nstepout, int resetstep,
1090 int gmx_unused nmultisim, int repl_ex_nst, int repl_ex_nex,
1091 int repl_ex_seed, real pforce, real cpt_period, real max_hours,
1092 const char *deviceOptions, int imdport, unsigned long Flags)
1094 gmx_bool bForceUseGPU, bTryUseGPU;
1095 t_inputrec *inputrec;
1096 t_state *state = NULL;
1098 gmx_ddbox_t ddbox = {0};
1099 int npme_major, npme_minor;
1101 gmx_mtop_t *mtop = NULL;
1102 t_mdatoms *mdatoms = NULL;
1103 t_forcerec *fr = NULL;
1104 t_fcdata *fcd = NULL;
1105 real ewaldcoeff_q = 0;
1106 real ewaldcoeff_lj = 0;
1107 gmx_pme_t *pmedata = NULL;
1108 gmx_vsite_t *vsite = NULL;
1109 gmx_constr_t constr;
1110 int nChargePerturbed = -1, nTypePerturbed = 0, status;
1111 gmx_wallcycle_t wcycle;
1113 gmx_walltime_accounting_t walltime_accounting = NULL;
1115 gmx_int64_t reset_counters;
1116 gmx_edsam_t ed = NULL;
1117 int nthreads_pme = 1;
1118 int nthreads_pp = 1;
1119 gmx_membed_t membed = NULL;
1120 gmx_hw_info_t *hwinfo = NULL;
1121 /* The master rank decides early on bUseGPU and broadcasts this later */
1122 gmx_bool bUseGPU = FALSE;
1124 /* CAUTION: threads may be started later on in this function, so
1125 cr doesn't reflect the final parallel state right now */
1129 if (Flags & MD_APPENDFILES)
1134 bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
1135 bTryUseGPU = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
1137 /* Detect hardware, gather information. This is an operation that is
1138 * global for this process (MPI rank). */
1139 hwinfo = gmx_detect_hardware(fplog, cr, bTryUseGPU);
1145 /* Read (nearly) all data required for the simulation */
1146 read_tpx_state(ftp2fn(efTPX, nfile, fnm), inputrec, state, NULL, mtop);
1148 if (inputrec->cutoff_scheme != ecutsVERLET &&
1149 ((Flags & MD_TESTVERLET) || getenv("GMX_VERLET_SCHEME") != NULL))
1151 convert_to_verlet_scheme(fplog, inputrec, mtop, det(state->box));
1154 if (inputrec->cutoff_scheme == ecutsVERLET)
1156 /* Here the master rank decides if all ranks will use GPUs */
1157 bUseGPU = (hwinfo->gpu_info.ncuda_dev_compatible > 0 ||
1158 getenv("GMX_EMULATE_GPU") != NULL);
1160 /* TODO add GPU kernels for this and replace this check by:
1161 * (bUseGPU && (ir->vdwtype == evdwPME &&
1162 * ir->ljpme_combination_rule == eljpmeLB))
1163 * update the message text and the content of nbnxn_acceleration_supported.
1166 !nbnxn_acceleration_supported(fplog, cr, inputrec, bUseGPU))
1168 /* Fallback message printed by nbnxn_acceleration_supported */
1171 gmx_fatal(FARGS, "GPU acceleration requested, but not supported with the given input settings");
1176 prepare_verlet_scheme(fplog, cr,
1177 inputrec, nstlist_cmdline, mtop, state->box,
1182 if (nstlist_cmdline > 0)
1184 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
1187 if (hwinfo->gpu_info.ncuda_dev_compatible > 0)
1189 md_print_warn(cr, fplog,
1190 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
1191 " To use a GPU, set the mdp option: cutoff-scheme = Verlet\n"
1192 " (for quick performance testing you can use the -testverlet option)\n");
1197 gmx_fatal(FARGS, "GPU requested, but can't be used without cutoff-scheme=Verlet");
1200 #ifdef GMX_TARGET_BGQ
1201 md_print_warn(cr, fplog,
1202 "NOTE: There is no SIMD implementation of the group scheme kernels on\n"
1203 " BlueGene/Q. You will observe better performance from using the\n"
1204 " Verlet cut-off scheme.\n");
1208 if (inputrec->eI == eiSD2)
1210 md_print_warn(cr, fplog, "The stochastic dynamics integrator %s is deprecated, since\n"
1211 "it is slower than integrator %s and is slightly less accurate\n"
1212 "with constraints. Use the %s integrator.",
1213 ei_names[inputrec->eI], ei_names[eiSD1], ei_names[eiSD1]);
1217 /* Check and update the hardware options for internal consistency */
1218 check_and_update_hw_opt_1(hw_opt, SIMMASTER(cr));
1220 /* Early check for externally set process affinity. */
1221 gmx_check_thread_affinity_set(fplog, cr,
1222 hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1226 #ifdef GMX_THREAD_MPI
1227 if (cr->npmenodes > 0 && hw_opt->nthreads_tmpi <= 0)
1229 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
1233 if (hw_opt->nthreads_omp_pme != hw_opt->nthreads_omp &&
1236 gmx_fatal(FARGS, "You need to explicitly specify the number of PME ranks (-npme) when using different number of OpenMP threads for PP and PME ranks");
1240 #ifdef GMX_THREAD_MPI
1243 /* Since the master knows the cut-off scheme, update hw_opt for this.
1244 * This is done later for normal MPI and also once more with tMPI
1245 * for all tMPI ranks.
1247 check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1249 /* NOW the threads will be started: */
1250 hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
1254 if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp <= 0)
1256 hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
1259 if (hw_opt->nthreads_tmpi > 1)
1261 t_commrec *cr_old = cr;
1262 /* now start the threads. */
1263 cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
1264 oenv, bVerbose, bCompact, nstglobalcomm,
1265 ddxyz, dd_node_order, rdd, rconstr,
1266 dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
1267 nbpu_opt, nstlist_cmdline,
1268 nsteps_cmdline, nstepout, resetstep, nmultisim,
1269 repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
1270 cpt_period, max_hours, deviceOptions,
1272 /* the main thread continues here with a new cr. We don't deallocate
1273 the old cr because other threads may still be reading it. */
1276 gmx_comm("Failed to spawn threads");
1281 /* END OF CAUTION: cr is now reliable */
1283 /* g_membed initialisation *
1284 * Because we change the mtop, init_membed is called before the init_parallel *
1285 * (in case we ever want to make it run in parallel) */
1286 if (opt2bSet("-membed", nfile, fnm))
1290 fprintf(stderr, "Initializing membed");
1292 membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
1297 /* now broadcast everything to the non-master nodes/threads: */
1298 init_parallel(cr, inputrec, mtop);
1302 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
1305 /* now make sure the state is initialized and propagated */
1306 set_state_entries(state, inputrec);
1308 /* A parallel command line option consistency check that we can
1309 only do after any threads have started. */
1311 (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
1314 "The -dd or -npme option request a parallel simulation, "
1316 "but %s was compiled without threads or MPI enabled"
1318 #ifdef GMX_THREAD_MPI
1319 "but the number of threads (option -nt) is 1"
1321 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
1324 , output_env_get_program_display_name(oenv)
1328 if ((Flags & MD_RERUN) &&
1329 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1331 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
1334 if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
1336 gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
1339 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
1341 if (cr->npmenodes > 0)
1343 gmx_fatal_collective(FARGS, cr, NULL,
1344 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
1353 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1357 /* NMR restraints must be initialized before load_checkpoint,
1358 * since with time averaging the history is added to t_state.
1359 * For proper consistency check we therefore need to extend
1361 * So the PME-only nodes (if present) will also initialize
1362 * the distance restraints.
1366 /* This needs to be called before read_checkpoint to extend the state */
1367 init_disres(fplog, mtop, inputrec, cr, fcd, state, repl_ex_nst > 0);
1369 init_orires(fplog, mtop, state->x, inputrec, cr, &(fcd->orires),
1372 if (DEFORM(*inputrec))
1374 /* Store the deform reference box before reading the checkpoint */
1377 copy_mat(state->box, box);
1381 gmx_bcast(sizeof(box), box, cr);
1383 /* Because we do not have the update struct available yet
1384 * in which the reference values should be stored,
1385 * we store them temporarily in static variables.
1386 * This should be thread safe, since they are only written once
1387 * and with identical values.
1389 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
1390 deform_init_init_step_tpx = inputrec->init_step;
1391 copy_mat(box, deform_init_box_tpx);
1392 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
1395 if (opt2bSet("-cpi", nfile, fnm))
1397 /* Check if checkpoint file exists before doing continuation.
1398 * This way we can use identical input options for the first and subsequent runs...
1400 if (gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr) )
1402 load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
1404 inputrec, state, &bReadEkin,
1405 (Flags & MD_APPENDFILES),
1406 (Flags & MD_APPENDFILESSET));
1410 Flags |= MD_READ_EKIN;
1415 if (((MASTER(cr) || (Flags & MD_SEPPOT)) && (Flags & MD_APPENDFILES))
1416 #ifdef GMX_THREAD_MPI
1417 /* With thread MPI only the master node/thread exists in mdrun.c,
1418 * therefore non-master nodes need to open the "seppot" log file here.
1420 || (!MASTER(cr) && (Flags & MD_SEPPOT))
1424 gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr, !(Flags & MD_SEPPOT),
1428 /* override nsteps with value from cmdline */
1429 override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
1433 copy_mat(state->box, box);
1438 gmx_bcast(sizeof(box), box, cr);
1441 /* Essential dynamics */
1442 if (opt2bSet("-ei", nfile, fnm))
1444 /* Open input and output files, allocate space for ED data structure */
1445 ed = ed_open(mtop->natoms, &state->edsamstate, nfile, fnm, Flags, oenv, cr);
1448 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
1449 inputrec->eI == eiNM))
1451 cr->dd = init_domain_decomposition(fplog, cr, Flags, ddxyz, rdd, rconstr,
1452 dddlb_opt, dlb_scale,
1453 ddcsx, ddcsy, ddcsz,
1456 &ddbox, &npme_major, &npme_minor);
1458 make_dd_communicators(fplog, cr, dd_node_order);
1460 /* Set overallocation to avoid frequent reallocation of arrays */
1461 set_over_alloc_dd(TRUE);
1465 /* PME, if used, is done on all nodes with 1D decomposition */
1467 cr->duty = (DUTY_PP | DUTY_PME);
1471 if (inputrec->ePBC == epbcSCREW)
1474 "pbc=%s is only implemented with domain decomposition",
1475 epbc_names[inputrec->ePBC]);
1481 /* After possible communicator splitting in make_dd_communicators.
1482 * we can set up the intra/inter node communication.
1484 gmx_setup_nodecomm(fplog, cr);
1487 /* Initialize per-physical-node MPI process/thread ID and counters. */
1488 gmx_init_intranode_counters(cr);
1492 md_print_info(cr, fplog,
1493 "This is simulation %d out of %d running as a composite Gromacs\n"
1494 "multi-simulation job. Setup for this simulation:\n\n",
1495 cr->ms->sim, cr->ms->nsim);
1497 md_print_info(cr, fplog, "Using %d MPI %s\n",
1499 #ifdef GMX_THREAD_MPI
1500 cr->nnodes == 1 ? "thread" : "threads"
1502 cr->nnodes == 1 ? "process" : "processes"
1508 /* Check and update hw_opt for the cut-off scheme */
1509 check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1511 gmx_omp_nthreads_init(fplog, cr,
1512 hwinfo->nthreads_hw_avail,
1513 hw_opt->nthreads_omp,
1514 hw_opt->nthreads_omp_pme,
1515 (cr->duty & DUTY_PP) == 0,
1516 inputrec->cutoff_scheme == ecutsVERLET);
1520 /* The master rank decided on the use of GPUs,
1521 * broadcast this information to all ranks.
1523 gmx_bcast_sim(sizeof(bUseGPU), &bUseGPU, cr);
1528 if (cr->npmenodes == -1)
1530 /* Don't automatically use PME-only nodes with GPUs */
1534 /* Select GPU id's to use */
1535 gmx_select_gpu_ids(fplog, cr, &hwinfo->gpu_info, bForceUseGPU,
1540 /* Ignore (potentially) manually selected GPUs */
1541 hw_opt->gpu_opt.ncuda_dev_use = 0;
1544 /* check consistency across ranks of things like SIMD
1545 * support and number of GPUs selected */
1546 gmx_check_hw_runconf_consistency(fplog, hwinfo, cr, hw_opt, bUseGPU);
1548 if (DOMAINDECOMP(cr))
1550 /* When we share GPUs over ranks, we need to know this for the DLB */
1551 dd_setup_dlb_resource_sharing(cr, hwinfo, hw_opt);
1554 /* getting number of PP/PME threads
1555 PME: env variable should be read only on one node to make sure it is
1556 identical everywhere;
1558 /* TODO nthreads_pp is only used for pinning threads.
1559 * This is a temporary solution until we have a hw topology library.
1561 nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
1562 nthreads_pme = gmx_omp_nthreads_get(emntPME);
1564 wcycle = wallcycle_init(fplog, resetstep, cr, nthreads_pp, nthreads_pme);
1568 /* Master synchronizes its value of reset_counters with all nodes
1569 * including PME only nodes */
1570 reset_counters = wcycle_get_reset_counters(wcycle);
1571 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1572 wcycle_set_reset_counters(wcycle, reset_counters);
1576 if (cr->duty & DUTY_PP)
1578 bcast_state(cr, state);
1580 /* Initiate forcerecord */
1582 fr->hwinfo = hwinfo;
1583 fr->gpu_opt = &hw_opt->gpu_opt;
1584 init_forcerec(fplog, oenv, fr, fcd, inputrec, mtop, cr, box,
1585 opt2fn("-table", nfile, fnm),
1586 opt2fn("-tabletf", nfile, fnm),
1587 opt2fn("-tablep", nfile, fnm),
1588 opt2fn("-tableb", nfile, fnm),
1593 /* version for PCA_NOT_READ_NODE (see md.c) */
1594 /*init_forcerec(fplog,fr,fcd,inputrec,mtop,cr,box,FALSE,
1595 "nofile","nofile","nofile","nofile",FALSE,pforce);
1597 fr->bSepDVDL = ((Flags & MD_SEPPOT) == MD_SEPPOT);
1599 /* Initialize QM-MM */
1602 init_QMMMrec(cr, mtop, inputrec, fr);
1605 /* Initialize the mdatoms structure.
1606 * mdatoms is not filled with atom data,
1607 * as this can not be done now with domain decomposition.
1609 mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1611 /* Initialize the virtual site communication */
1612 vsite = init_vsite(mtop, cr, FALSE);
1614 calc_shifts(box, fr->shift_vec);
1616 /* With periodic molecules the charge groups should be whole at start up
1617 * and the virtual sites should not be far from their proper positions.
1619 if (!inputrec->bContinuation && MASTER(cr) &&
1620 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1622 /* Make molecules whole at start of run */
1623 if (fr->ePBC != epbcNONE)
1625 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, state->x);
1629 /* Correct initial vsite positions are required
1630 * for the initial distribution in the domain decomposition
1631 * and for the initial shell prediction.
1633 construct_vsites_mtop(vsite, mtop, state->x);
1637 if (EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype))
1639 ewaldcoeff_q = fr->ewaldcoeff_q;
1640 ewaldcoeff_lj = fr->ewaldcoeff_lj;
1641 pmedata = &fr->pmedata;
1650 /* This is a PME only node */
1652 /* We don't need the state */
1655 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1656 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1660 if (hw_opt->thread_affinity != threadaffOFF)
1662 /* Before setting affinity, check whether the affinity has changed
1663 * - which indicates that probably the OpenMP library has changed it
1664 * since we first checked).
1666 gmx_check_thread_affinity_set(fplog, cr,
1667 hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1669 /* Set the CPU affinity */
1670 gmx_set_thread_affinity(fplog, cr, hw_opt, hwinfo);
1673 /* Initiate PME if necessary,
1674 * either on all nodes or on dedicated PME nodes only. */
1675 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1679 nChargePerturbed = mdatoms->nChargePerturbed;
1680 if (EVDW_PME(inputrec->vdwtype))
1682 nTypePerturbed = mdatoms->nTypePerturbed;
1685 if (cr->npmenodes > 0)
1687 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1688 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1689 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1692 if (cr->duty & DUTY_PME)
1694 status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1695 mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1696 (Flags & MD_REPRODUCIBLE), nthreads_pme);
1699 gmx_fatal(FARGS, "Error %d initializing PME", status);
1705 if (integrator[inputrec->eI].func == do_md)
1707 /* Turn on signal handling on all nodes */
1709 * (A user signal from the PME nodes (if any)
1710 * is communicated to the PP nodes.
1712 signal_handler_install();
1715 if (cr->duty & DUTY_PP)
1717 /* Assumes uniform use of the number of OpenMP threads */
1718 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1720 if (inputrec->ePull != epullNO)
1722 /* Initialize pull code */
1723 init_pull(fplog, inputrec, nfile, fnm, mtop, cr, oenv, inputrec->fepvals->init_lambda,
1724 EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1729 /* Initialize enforced rotation code */
1730 init_rot(fplog, inputrec, nfile, fnm, cr, state->x, box, mtop, oenv,
1734 if (inputrec->eSwapCoords != eswapNO)
1736 /* Initialize ion swapping code */
1737 init_swapcoords(fplog, bVerbose, inputrec, opt2fn_master("-swap", nfile, fnm, cr),
1738 mtop, state->x, state->box, &state->swapstate, cr, oenv, Flags);
1741 constr = init_constraints(fplog, mtop, inputrec, ed, state, cr);
1743 if (DOMAINDECOMP(cr))
1745 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1746 dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1747 Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1749 set_dd_parameters(fplog, cr->dd, dlb_scale, inputrec, &ddbox);
1751 setup_dd_grid(fplog, cr->dd);
1754 /* Now do whatever the user wants us to do (how flexible...) */
1755 integrator[inputrec->eI].func(fplog, cr, nfile, fnm,
1756 oenv, bVerbose, bCompact,
1759 nstepout, inputrec, mtop,
1761 mdatoms, nrnb, wcycle, ed, fr,
1762 repl_ex_nst, repl_ex_nex, repl_ex_seed,
1764 cpt_period, max_hours,
1768 walltime_accounting);
1770 if (inputrec->ePull != epullNO)
1772 finish_pull(inputrec->pull);
1777 finish_rot(inputrec->rot);
1783 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1785 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1786 gmx_pmeonly(*pmedata, cr, nrnb, wcycle, walltime_accounting, ewaldcoeff_q, ewaldcoeff_lj, inputrec);
1789 wallcycle_stop(wcycle, ewcRUN);
1791 /* Finish up, write some stuff
1792 * if rerunMD, don't write last frame again
1794 finish_run(fplog, cr,
1795 inputrec, nrnb, wcycle, walltime_accounting,
1796 fr != NULL && fr->nbv != NULL && fr->nbv->bUseGPU ?
1797 nbnxn_cuda_get_timings(fr->nbv->cu_nbv) : NULL,
1798 EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1801 /* Free GPU memory and context */
1802 free_gpu_resources(fr, cr);
1804 if (opt2bSet("-membed", nfile, fnm))
1809 gmx_hardware_info_free(hwinfo);
1811 /* Does what it says */
1812 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1813 walltime_accounting_destroy(walltime_accounting);
1815 /* Close logfile already here if we were appending to it */
1816 if (MASTER(cr) && (Flags & MD_APPENDFILES))
1818 gmx_log_close(fplog);
1821 rc = (int)gmx_get_stop_condition();
1823 #ifdef GMX_THREAD_MPI
1824 /* we need to join all threads. The sub-threads join when they
1825 exit this function, but the master thread needs to be told to
1827 if (PAR(cr) && MASTER(cr))