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.
56 #include "md_logging.h"
57 #include "md_support.h"
71 #include "checkpoint.h"
72 #include "gromacs/topology/mtop_util.h"
73 #include "sighandler.h"
75 #include "gmx_detect_hardware.h"
76 #include "gmx_omp_nthreads.h"
77 #include "gromacs/gmxpreprocess/calc_verletbuf.h"
79 #include "gmx_thread_affinity.h"
83 #include "gromacs/essentialdynamics/edsam.h"
84 #include "gromacs/fileio/tpxio.h"
85 #include "gromacs/math/vec.h"
86 #include "gromacs/mdlib/nbnxn_search.h"
87 #include "gromacs/mdlib/nbnxn_consts.h"
88 #include "gromacs/pbcutil/pbc.h"
89 #include "gromacs/pulling/pull.h"
90 #include "gromacs/pulling/pull_rotation.h"
91 #include "gromacs/swap/swapcoords.h"
92 #include "gromacs/timing/wallcycle.h"
93 #include "gromacs/utility/gmxassert.h"
94 #include "gromacs/utility/gmxmpi.h"
95 #include "gromacs/utility/smalloc.h"
101 #include "gpu_utils.h"
102 #include "nbnxn_cuda_data_mgmt.h"
105 gmx_integrator_t *func;
108 /* The array should match the eI array in include/types/enums.h */
109 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}};
111 gmx_int64_t deform_init_init_step_tpx;
112 matrix deform_init_box_tpx;
113 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
116 #ifdef GMX_THREAD_MPI
117 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
118 * the number of threads will get lowered.
120 #define MIN_ATOMS_PER_MPI_THREAD 90
121 #define MIN_ATOMS_PER_GPU 900
123 struct mdrunner_arglist
138 const char *dddlb_opt;
143 const char *nbpu_opt;
145 gmx_int64_t nsteps_cmdline;
155 const char *deviceOptions;
161 /* The function used for spawning threads. Extracts the mdrunner()
162 arguments from its one argument and calls mdrunner(), after making
164 static void mdrunner_start_fn(void *arg)
166 struct mdrunner_arglist *mda = (struct mdrunner_arglist*)arg;
167 struct mdrunner_arglist mc = *mda; /* copy the arg list to make sure
168 that it's thread-local. This doesn't
169 copy pointed-to items, of course,
170 but those are all const. */
171 t_commrec *cr; /* we need a local version of this */
175 fnm = dup_tfn(mc.nfile, mc.fnm);
177 cr = reinitialize_commrec_for_this_thread(mc.cr);
184 mdrunner(&mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
185 mc.bVerbose, mc.bCompact, mc.nstglobalcomm,
186 mc.ddxyz, mc.dd_node_order, mc.rdd,
187 mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
188 mc.ddcsx, mc.ddcsy, mc.ddcsz,
189 mc.nbpu_opt, mc.nstlist_cmdline,
190 mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
191 mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce,
192 mc.cpt_period, mc.max_hours, mc.deviceOptions, mc.imdport, mc.Flags);
195 /* called by mdrunner() to start a specific number of threads (including
196 the main thread) for thread-parallel runs. This in turn calls mdrunner()
198 All options besides nthreads are the same as for mdrunner(). */
199 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
200 FILE *fplog, t_commrec *cr, int nfile,
201 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
202 gmx_bool bCompact, int nstglobalcomm,
203 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
204 const char *dddlb_opt, real dlb_scale,
205 const char *ddcsx, const char *ddcsy, const char *ddcsz,
206 const char *nbpu_opt, int nstlist_cmdline,
207 gmx_int64_t nsteps_cmdline,
208 int nstepout, int resetstep,
209 int nmultisim, int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
210 real pforce, real cpt_period, real max_hours,
211 const char *deviceOptions, unsigned long Flags)
214 struct mdrunner_arglist *mda;
215 t_commrec *crn; /* the new commrec */
218 /* first check whether we even need to start tMPI */
219 if (hw_opt->nthreads_tmpi < 2)
224 /* a few small, one-time, almost unavoidable memory leaks: */
226 fnmn = dup_tfn(nfile, fnm);
228 /* fill the data structure to pass as void pointer to thread start fn */
229 /* hw_opt contains pointers, which should all be NULL at this stage */
230 mda->hw_opt = *hw_opt;
236 mda->bVerbose = bVerbose;
237 mda->bCompact = bCompact;
238 mda->nstglobalcomm = nstglobalcomm;
239 mda->ddxyz[XX] = ddxyz[XX];
240 mda->ddxyz[YY] = ddxyz[YY];
241 mda->ddxyz[ZZ] = ddxyz[ZZ];
242 mda->dd_node_order = dd_node_order;
244 mda->rconstr = rconstr;
245 mda->dddlb_opt = dddlb_opt;
246 mda->dlb_scale = dlb_scale;
250 mda->nbpu_opt = nbpu_opt;
251 mda->nstlist_cmdline = nstlist_cmdline;
252 mda->nsteps_cmdline = nsteps_cmdline;
253 mda->nstepout = nstepout;
254 mda->resetstep = resetstep;
255 mda->nmultisim = nmultisim;
256 mda->repl_ex_nst = repl_ex_nst;
257 mda->repl_ex_nex = repl_ex_nex;
258 mda->repl_ex_seed = repl_ex_seed;
259 mda->pforce = pforce;
260 mda->cpt_period = cpt_period;
261 mda->max_hours = max_hours;
262 mda->deviceOptions = deviceOptions;
265 /* now spawn new threads that start mdrunner_start_fn(), while
266 the main thread returns, we set thread affinity later */
267 ret = tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi, TMPI_AFFINITY_NONE,
268 mdrunner_start_fn, (void*)(mda) );
269 if (ret != TMPI_SUCCESS)
274 crn = reinitialize_commrec_for_this_thread(cr);
279 static int get_tmpi_omp_thread_division(const gmx_hw_info_t *hwinfo,
280 const gmx_hw_opt_t *hw_opt,
286 /* There are no separate PME nodes here, as we ensured in
287 * check_and_update_hw_opt that nthreads_tmpi>0 with PME nodes
288 * and a conditional ensures we would not have ended up here.
289 * Note that separate PME nodes might be switched on later.
293 nthreads_tmpi = ngpu;
294 if (nthreads_tot > 0 && nthreads_tot < nthreads_tmpi)
296 nthreads_tmpi = nthreads_tot;
299 else if (hw_opt->nthreads_omp > 0)
301 /* Here we could oversubscribe, when we do, we issue a warning later */
302 nthreads_tmpi = std::max(1, nthreads_tot/hw_opt->nthreads_omp);
306 /* TODO choose nthreads_omp based on hardware topology
307 when we have a hardware topology detection library */
308 /* In general, when running up to 4 threads, OpenMP should be faster.
309 * Note: on AMD Bulldozer we should avoid running OpenMP over two dies.
310 * On Intel>=Nehalem running OpenMP on a single CPU is always faster,
311 * even on two CPUs it's usually faster (but with many OpenMP threads
312 * it could be faster not to use HT, currently we always use HT).
313 * On Nehalem/Westmere we want to avoid running 16 threads over
314 * two CPUs with HT, so we need a limit<16; thus we use 12.
315 * A reasonable limit for Intel Sandy and Ivy bridge,
316 * not knowing the topology, is 16 threads.
318 const int nthreads_omp_always_faster = 4;
319 const int nthreads_omp_always_faster_Nehalem = 12;
320 const int nthreads_omp_always_faster_SandyBridge = 16;
321 gmx_bool bIntel_Family6;
324 (gmx_cpuid_vendor(hwinfo->cpuid_info) == GMX_CPUID_VENDOR_INTEL &&
325 gmx_cpuid_family(hwinfo->cpuid_info) == 6);
327 if (nthreads_tot <= nthreads_omp_always_faster ||
329 ((gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_Nehalem && nthreads_tot <= nthreads_omp_always_faster_Nehalem) ||
330 (gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_SandyBridge && nthreads_tot <= nthreads_omp_always_faster_SandyBridge))))
332 /* Use pure OpenMP parallelization */
337 /* Don't use OpenMP parallelization */
338 nthreads_tmpi = nthreads_tot;
342 return nthreads_tmpi;
346 /* Get the number of threads to use for thread-MPI based on how many
347 * were requested, which algorithms we're using,
348 * and how many particles there are.
349 * At the point we have already called check_and_update_hw_opt.
350 * Thus all options should be internally consistent and consistent
351 * with the hardware, except that ntmpi could be larger than #GPU.
353 static int get_nthreads_mpi(const gmx_hw_info_t *hwinfo,
354 gmx_hw_opt_t *hw_opt,
355 t_inputrec *inputrec, gmx_mtop_t *mtop,
359 int nthreads_hw, nthreads_tot_max, nthreads_tmpi, nthreads_new, ngpu;
360 int min_atoms_per_mpi_thread;
363 if (hw_opt->nthreads_tmpi > 0)
365 /* Trivial, return right away */
366 return hw_opt->nthreads_tmpi;
369 nthreads_hw = hwinfo->nthreads_hw_avail;
371 /* How many total (#tMPI*#OpenMP) threads can we start? */
372 if (hw_opt->nthreads_tot > 0)
374 nthreads_tot_max = hw_opt->nthreads_tot;
378 nthreads_tot_max = nthreads_hw;
381 bCanUseGPU = (inputrec->cutoff_scheme == ecutsVERLET &&
382 hwinfo->gpu_info.ncuda_dev_compatible > 0);
385 ngpu = hwinfo->gpu_info.ncuda_dev_compatible;
392 if (inputrec->cutoff_scheme == ecutsGROUP)
394 /* We checked this before, but it doesn't hurt to do it once more */
395 assert(hw_opt->nthreads_omp == 1);
399 get_tmpi_omp_thread_division(hwinfo, hw_opt, nthreads_tot_max, ngpu);
401 if (inputrec->eI == eiNM || EI_TPI(inputrec->eI))
403 /* Dims/steps are divided over the nodes iso splitting the atoms */
404 min_atoms_per_mpi_thread = 0;
410 min_atoms_per_mpi_thread = MIN_ATOMS_PER_GPU;
414 min_atoms_per_mpi_thread = MIN_ATOMS_PER_MPI_THREAD;
418 /* Check if an algorithm does not support parallel simulation. */
419 if (nthreads_tmpi != 1 &&
420 ( inputrec->eI == eiLBFGS ||
421 inputrec->coulombtype == eelEWALD ) )
425 md_print_warn(cr, fplog, "The integration or electrostatics algorithm doesn't support parallel runs. Using a single thread-MPI thread.\n");
426 if (hw_opt->nthreads_tmpi > nthreads_tmpi)
428 gmx_fatal(FARGS, "You asked for more than 1 thread-MPI thread, but an algorithm doesn't support that");
431 else if (mtop->natoms/nthreads_tmpi < min_atoms_per_mpi_thread)
433 /* the thread number was chosen automatically, but there are too many
434 threads (too few atoms per thread) */
435 nthreads_new = std::max(1, mtop->natoms/min_atoms_per_mpi_thread);
437 /* Avoid partial use of Hyper-Threading */
438 if (gmx_cpuid_x86_smt(hwinfo->cpuid_info) == GMX_CPUID_X86_SMT_ENABLED &&
439 nthreads_new > nthreads_hw/2 && nthreads_new < nthreads_hw)
441 nthreads_new = nthreads_hw/2;
444 /* Avoid large prime numbers in the thread count */
445 if (nthreads_new >= 6)
447 /* Use only 6,8,10 with additional factors of 2 */
451 while (3*fac*2 <= nthreads_new)
456 nthreads_new = (nthreads_new/fac)*fac;
461 if (nthreads_new == 5)
467 nthreads_tmpi = nthreads_new;
469 fprintf(stderr, "\n");
470 fprintf(stderr, "NOTE: Parallelization is limited by the small number of atoms,\n");
471 fprintf(stderr, " only starting %d thread-MPI threads.\n", nthreads_tmpi);
472 fprintf(stderr, " You can use the -nt and/or -ntmpi option to optimize the number of threads.\n\n");
475 return nthreads_tmpi;
477 #endif /* GMX_THREAD_MPI */
480 /* We determine the extra cost of the non-bonded kernels compared to
481 * a reference nstlist value of 10 (which is the default in grompp).
483 static const int nbnxnReferenceNstlist = 10;
484 /* The values to try when switching */
485 const int nstlist_try[] = { 20, 25, 40 };
486 #define NNSTL sizeof(nstlist_try)/sizeof(nstlist_try[0])
487 /* Increase nstlist until the non-bonded cost increases more than listfac_ok,
488 * but never more than listfac_max.
489 * A standard (protein+)water system at 300K with PME ewald_rtol=1e-5
490 * needs 1.28 at rcoulomb=0.9 and 1.24 at rcoulomb=1.0 to get to nstlist=40.
491 * Note that both CPU and GPU factors are conservative. Performance should
492 * not go down due to this tuning, except with a relatively slow GPU.
493 * On the other hand, at medium/high parallelization or with fast GPUs
494 * nstlist will not be increased enough to reach optimal performance.
496 /* CPU: pair-search is about a factor 1.5 slower than the non-bonded kernel */
497 static const float nbnxn_cpu_listfac_ok = 1.05;
498 static const float nbnxn_cpu_listfac_max = 1.09;
499 /* GPU: pair-search is a factor 1.5-3 slower than the non-bonded kernel */
500 static const float nbnxn_gpu_listfac_ok = 1.20;
501 static const float nbnxn_gpu_listfac_max = 1.30;
503 /* Try to increase nstlist when using the Verlet cut-off scheme */
504 static void increase_nstlist(FILE *fp, t_commrec *cr,
505 t_inputrec *ir, int nstlist_cmdline,
506 const gmx_mtop_t *mtop, matrix box,
509 float listfac_ok, listfac_max;
510 int nstlist_orig, nstlist_prev;
511 verletbuf_list_setup_t ls;
512 real rlistWithReferenceNstlist, rlist_inc, rlist_ok, rlist_max;
513 real rlist_new, rlist_prev;
514 size_t nstlist_ind = 0;
516 gmx_bool bBox, bDD, bCont;
517 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";
518 const char *nve_err = "Can not increase nstlist because an NVE ensemble is used";
519 const char *vbd_err = "Can not increase nstlist because verlet-buffer-tolerance is not set or used";
520 const char *box_err = "Can not increase nstlist because the box is too small";
521 const char *dd_err = "Can not increase nstlist because of domain decomposition limitations";
523 const float oneThird = 1.0f / 3.0f;
525 if (nstlist_cmdline <= 0)
527 if (ir->nstlist == 1)
529 /* The user probably set nstlist=1 for a reason,
530 * don't mess with the settings.
535 if (fp != NULL && bGPU && ir->nstlist < nstlist_try[0])
537 fprintf(fp, nstl_gpu, ir->nstlist);
540 while (nstlist_ind < NNSTL && ir->nstlist >= nstlist_try[nstlist_ind])
544 if (nstlist_ind == NNSTL)
546 /* There are no larger nstlist value to try */
551 if (EI_MD(ir->eI) && ir->etc == etcNO)
555 fprintf(stderr, "%s\n", nve_err);
559 fprintf(fp, "%s\n", nve_err);
565 if (ir->verletbuf_tol == 0 && bGPU)
567 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");
570 if (ir->verletbuf_tol < 0)
574 fprintf(stderr, "%s\n", vbd_err);
578 fprintf(fp, "%s\n", vbd_err);
586 listfac_ok = nbnxn_gpu_listfac_ok;
587 listfac_max = nbnxn_gpu_listfac_max;
591 listfac_ok = nbnxn_cpu_listfac_ok;
592 listfac_max = nbnxn_cpu_listfac_max;
595 nstlist_orig = ir->nstlist;
596 if (nstlist_cmdline > 0)
600 sprintf(buf, "Getting nstlist=%d from command line option",
603 ir->nstlist = nstlist_cmdline;
606 verletbuf_get_list_setup(bGPU, &ls);
608 /* Allow rlist to make the list a given factor larger than the list
609 * would be with the reference value for nstlist (10).
611 nstlist_prev = ir->nstlist;
612 ir->nstlist = nbnxnReferenceNstlist;
613 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL,
614 &rlistWithReferenceNstlist);
615 ir->nstlist = nstlist_prev;
617 /* Determine the pair list size increase due to zero interactions */
618 rlist_inc = nbnxn_get_rlist_effective_inc(ls.cluster_size_j,
619 mtop->natoms/det(box));
620 rlist_ok = (rlistWithReferenceNstlist + rlist_inc)*pow(listfac_ok, oneThird) - rlist_inc;
621 rlist_max = (rlistWithReferenceNstlist + rlist_inc)*pow(listfac_max, oneThird) - rlist_inc;
624 fprintf(debug, "nstlist tuning: rlist_inc %.3f rlist_ok %.3f rlist_max %.3f\n",
625 rlist_inc, rlist_ok, rlist_max);
628 nstlist_prev = nstlist_orig;
629 rlist_prev = ir->rlist;
632 if (nstlist_cmdline <= 0)
634 ir->nstlist = nstlist_try[nstlist_ind];
637 /* Set the pair-list buffer size in ir */
638 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
640 /* Does rlist fit in the box? */
641 bBox = (sqr(rlist_new) < max_cutoff2(ir->ePBC, box));
643 if (bBox && DOMAINDECOMP(cr))
645 /* Check if rlist fits in the domain decomposition */
646 if (inputrec2nboundeddim(ir) < DIM)
648 gmx_incons("Changing nstlist with domain decomposition and unbounded dimensions is not implemented yet");
650 copy_mat(box, state_tmp.box);
651 bDD = change_dd_cutoff(cr, &state_tmp, ir, rlist_new);
656 fprintf(debug, "nstlist %d rlist %.3f bBox %d bDD %d\n",
657 ir->nstlist, rlist_new, bBox, bDD);
662 if (nstlist_cmdline <= 0)
664 if (bBox && bDD && rlist_new <= rlist_max)
666 /* Increase nstlist */
667 nstlist_prev = ir->nstlist;
668 rlist_prev = rlist_new;
669 bCont = (nstlist_ind+1 < NNSTL && rlist_new < rlist_ok);
673 /* Stick with the previous nstlist */
674 ir->nstlist = nstlist_prev;
675 rlist_new = rlist_prev;
687 gmx_warning(!bBox ? box_err : dd_err);
690 fprintf(fp, "\n%s\n", bBox ? box_err : dd_err);
692 ir->nstlist = nstlist_orig;
694 else if (ir->nstlist != nstlist_orig || rlist_new != ir->rlist)
696 sprintf(buf, "Changing nstlist from %d to %d, rlist from %g to %g",
697 nstlist_orig, ir->nstlist,
698 ir->rlist, rlist_new);
701 fprintf(stderr, "%s\n\n", buf);
705 fprintf(fp, "%s\n\n", buf);
707 ir->rlist = rlist_new;
708 ir->rlistlong = rlist_new;
712 static void prepare_verlet_scheme(FILE *fplog,
716 const gmx_mtop_t *mtop,
720 /* For NVE simulations, we will retain the initial list buffer */
721 if (ir->verletbuf_tol > 0 && !(EI_MD(ir->eI) && ir->etc == etcNO))
723 /* Update the Verlet buffer size for the current run setup */
724 verletbuf_list_setup_t ls;
727 /* Here we assume SIMD-enabled kernels are being used. But as currently
728 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
729 * and 4x2 gives a larger buffer than 4x4, this is ok.
731 verletbuf_get_list_setup(bUseGPU, &ls);
733 calc_verlet_buffer_size(mtop, det(box), ir, -1, &ls, NULL, &rlist_new);
735 if (rlist_new != ir->rlist)
739 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
740 ir->rlist, rlist_new,
741 ls.cluster_size_i, ls.cluster_size_j);
743 ir->rlist = rlist_new;
744 ir->rlistlong = rlist_new;
748 if (nstlist_cmdline > 0 && (!EI_DYNAMICS(ir->eI) || ir->verletbuf_tol <= 0))
750 gmx_fatal(FARGS, "Can not set nstlist without %s",
751 !EI_DYNAMICS(ir->eI) ? "dynamics" : "verlet-buffer-tolerance");
754 if (EI_DYNAMICS(ir->eI))
756 /* Set or try nstlist values */
757 increase_nstlist(fplog, cr, ir, nstlist_cmdline, mtop, box, bUseGPU);
761 static void convert_to_verlet_scheme(FILE *fplog,
763 gmx_mtop_t *mtop, real box_vol)
765 const char *conv_mesg = "Converting input file with group cut-off scheme to the Verlet cut-off scheme";
767 md_print_warn(NULL, fplog, "%s\n", conv_mesg);
769 ir->cutoff_scheme = ecutsVERLET;
770 ir->verletbuf_tol = 0.005;
772 if (ir->rcoulomb != ir->rvdw)
774 gmx_fatal(FARGS, "The VdW and Coulomb cut-offs are different, whereas the Verlet scheme only supports equal cut-offs");
777 if (ir->vdwtype == evdwUSER || EEL_USER(ir->coulombtype))
779 gmx_fatal(FARGS, "User non-bonded potentials are not (yet) supported with the Verlet scheme");
781 else if (ir_vdw_switched(ir) || ir_coulomb_switched(ir))
783 if (ir_vdw_switched(ir) && ir->vdw_modifier == eintmodNONE)
785 ir->vdwtype = evdwCUT;
789 case evdwSHIFT: ir->vdw_modifier = eintmodFORCESWITCH; break;
790 case evdwSWITCH: ir->vdw_modifier = eintmodPOTSWITCH; break;
791 default: gmx_fatal(FARGS, "The Verlet scheme does not support Van der Waals interactions of type '%s'", evdw_names[ir->vdwtype]);
794 if (ir_coulomb_switched(ir) && ir->coulomb_modifier == eintmodNONE)
796 if (EEL_FULL(ir->coulombtype))
798 /* With full electrostatic only PME can be switched */
799 ir->coulombtype = eelPME;
800 ir->coulomb_modifier = eintmodPOTSHIFT;
804 md_print_warn(NULL, fplog, "NOTE: Replacing %s electrostatics with reaction-field with epsilon-rf=inf\n", eel_names[ir->coulombtype]);
805 ir->coulombtype = eelRF;
806 ir->epsilon_rf = 0.0;
807 ir->coulomb_modifier = eintmodPOTSHIFT;
811 /* We set the pair energy error tolerance to a small number.
812 * Note that this is only for testing. For production the user
813 * should think about this and set the mdp options.
815 ir->verletbuf_tol = 1e-4;
818 if (inputrec2nboundeddim(ir) != 3)
820 gmx_fatal(FARGS, "Can only convert old tpr files to the Verlet cut-off scheme with 3D pbc");
823 if (ir->efep != efepNO || ir->implicit_solvent != eisNO)
825 gmx_fatal(FARGS, "Will not convert old tpr files to the Verlet cut-off scheme with free-energy calculations or implicit solvent");
828 if (EI_DYNAMICS(ir->eI) && !(EI_MD(ir->eI) && ir->etc == etcNO))
830 verletbuf_list_setup_t ls;
832 verletbuf_get_list_setup(FALSE, &ls);
833 calc_verlet_buffer_size(mtop, box_vol, ir, -1, &ls, NULL, &ir->rlist);
841 rlist_fac = 1 + verlet_buffer_ratio_NVE_T0;
845 rlist_fac = 1 + verlet_buffer_ratio_nodynamics;
847 ir->verletbuf_tol = -1;
848 ir->rlist = rlist_fac*std::max(ir->rvdw, ir->rcoulomb);
851 gmx_mtop_remove_chargegroups(mtop);
854 static void print_hw_opt(FILE *fp, const gmx_hw_opt_t *hw_opt)
856 fprintf(fp, "hw_opt: nt %d ntmpi %d ntomp %d ntomp_pme %d gpu_id '%s'\n",
857 hw_opt->nthreads_tot,
858 hw_opt->nthreads_tmpi,
859 hw_opt->nthreads_omp,
860 hw_opt->nthreads_omp_pme,
861 hw_opt->gpu_opt.gpu_id != NULL ? hw_opt->gpu_opt.gpu_id : "");
864 /* Checks we can do when we don't (yet) know the cut-off scheme */
865 static void check_and_update_hw_opt_1(gmx_hw_opt_t *hw_opt,
866 gmx_bool bIsSimMaster)
868 gmx_omp_nthreads_read_env(&hw_opt->nthreads_omp, bIsSimMaster);
870 #ifndef GMX_THREAD_MPI
871 if (hw_opt->nthreads_tot > 0)
873 gmx_fatal(FARGS, "Setting the total number of threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
875 if (hw_opt->nthreads_tmpi > 0)
877 gmx_fatal(FARGS, "Setting the number of thread-MPI threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
882 if (hw_opt->nthreads_omp > 1)
884 gmx_fatal(FARGS, "More than 1 OpenMP thread requested, but Gromacs was compiled without OpenMP support");
886 hw_opt->nthreads_omp = 1;
889 if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp_pme <= 0)
891 /* We have the same number of OpenMP threads for PP and PME processes,
892 * thus we can perform several consistency checks.
894 if (hw_opt->nthreads_tmpi > 0 &&
895 hw_opt->nthreads_omp > 0 &&
896 hw_opt->nthreads_tot != hw_opt->nthreads_tmpi*hw_opt->nthreads_omp)
898 gmx_fatal(FARGS, "The total number of threads requested (%d) does not match the thread-MPI threads (%d) times the OpenMP threads (%d) requested",
899 hw_opt->nthreads_tot, hw_opt->nthreads_tmpi, hw_opt->nthreads_omp);
902 if (hw_opt->nthreads_tmpi > 0 &&
903 hw_opt->nthreads_tot % hw_opt->nthreads_tmpi != 0)
905 gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of thread-MPI threads requested (%d)",
906 hw_opt->nthreads_tot, hw_opt->nthreads_tmpi);
909 if (hw_opt->nthreads_omp > 0 &&
910 hw_opt->nthreads_tot % hw_opt->nthreads_omp != 0)
912 gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of OpenMP threads requested (%d)",
913 hw_opt->nthreads_tot, hw_opt->nthreads_omp);
916 if (hw_opt->nthreads_tmpi > 0 &&
917 hw_opt->nthreads_omp <= 0)
919 hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
924 if (hw_opt->nthreads_omp > 1)
926 gmx_fatal(FARGS, "OpenMP threads are requested, but Gromacs was compiled without OpenMP support");
930 if (hw_opt->nthreads_omp_pme > 0 && hw_opt->nthreads_omp <= 0)
932 gmx_fatal(FARGS, "You need to specify -ntomp in addition to -ntomp_pme");
935 if (hw_opt->nthreads_tot == 1)
937 hw_opt->nthreads_tmpi = 1;
939 if (hw_opt->nthreads_omp > 1)
941 gmx_fatal(FARGS, "You requested %d OpenMP threads with %d total threads",
942 hw_opt->nthreads_tmpi, hw_opt->nthreads_tot);
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;
952 /* Parse GPU IDs, if provided.
953 * We check consistency with the tMPI thread count later.
955 gmx_parse_gpu_ids(&hw_opt->gpu_opt);
957 #ifdef GMX_THREAD_MPI
958 if (hw_opt->gpu_opt.ncuda_dev_use > 0 && hw_opt->nthreads_tmpi == 0)
960 /* Set the number of MPI threads equal to the number of GPUs */
961 hw_opt->nthreads_tmpi = hw_opt->gpu_opt.ncuda_dev_use;
963 if (hw_opt->nthreads_tot > 0 &&
964 hw_opt->nthreads_tmpi > hw_opt->nthreads_tot)
966 /* We have more GPUs than total threads requested.
967 * We choose to (later) generate a mismatch error,
968 * instead of launching more threads than requested.
970 hw_opt->nthreads_tmpi = hw_opt->nthreads_tot;
977 print_hw_opt(debug, hw_opt);
981 /* Checks we can do when we know the cut-off scheme */
982 static void check_and_update_hw_opt_2(gmx_hw_opt_t *hw_opt,
985 if (cutoff_scheme == ecutsGROUP)
987 /* We only have OpenMP support for PME only nodes */
988 if (hw_opt->nthreads_omp > 1)
990 gmx_fatal(FARGS, "OpenMP threads have been requested with cut-off scheme %s, but these are only supported with cut-off scheme %s",
991 ecutscheme_names[cutoff_scheme],
992 ecutscheme_names[ecutsVERLET]);
994 hw_opt->nthreads_omp = 1;
997 if (hw_opt->nthreads_omp_pme <= 0 && hw_opt->nthreads_omp > 0)
999 hw_opt->nthreads_omp_pme = hw_opt->nthreads_omp;
1004 print_hw_opt(debug, hw_opt);
1009 /* Override the value in inputrec with value passed on the command line (if any) */
1010 static void override_nsteps_cmdline(FILE *fplog,
1011 gmx_int64_t nsteps_cmdline,
1013 const t_commrec *cr)
1015 char sbuf[STEPSTRSIZE];
1020 /* override with anything else than the default -2 */
1021 if (nsteps_cmdline > -2)
1025 ir->nsteps = nsteps_cmdline;
1026 if (EI_DYNAMICS(ir->eI))
1028 sprintf(stmp, "Overriding nsteps with value passed on the command line: %s steps, %.3f ps",
1029 gmx_step_str(nsteps_cmdline, sbuf),
1030 nsteps_cmdline*ir->delta_t);
1034 sprintf(stmp, "Overriding nsteps with value passed on the command line: %s steps",
1035 gmx_step_str(nsteps_cmdline, sbuf));
1038 md_print_warn(cr, fplog, "%s\n", stmp);
1042 /* Frees GPU memory and destroys the CUDA context.
1044 * Note that this function needs to be called even if GPUs are not used
1045 * in this run because the PME ranks have no knowledge of whether GPUs
1046 * are used or not, but all ranks need to enter the barrier below.
1048 static void free_gpu_resources(const t_forcerec *fr,
1049 const t_commrec *cr)
1051 gmx_bool bIsPPrankUsingGPU;
1052 char gpu_err_str[STRLEN];
1054 bIsPPrankUsingGPU = (cr->duty & DUTY_PP) && fr != NULL && fr->nbv != NULL && fr->nbv->bUseGPU;
1056 if (bIsPPrankUsingGPU)
1058 /* free nbnxn data in GPU memory */
1059 nbnxn_cuda_free(fr->nbv->cu_nbv);
1061 /* With tMPI we need to wait for all ranks to finish deallocation before
1062 * destroying the context in free_gpu() as some ranks may be sharing
1064 * Note: as only PP ranks need to free GPU resources, so it is safe to
1065 * not call the barrier on PME ranks.
1067 #ifdef GMX_THREAD_MPI
1072 #endif /* GMX_THREAD_MPI */
1074 /* uninitialize GPU (by destroying the context) */
1075 if (!free_gpu(gpu_err_str))
1077 gmx_warning("On rank %d failed to free GPU #%d: %s",
1078 cr->nodeid, get_current_gpu_device_id(), gpu_err_str);
1083 int mdrunner(gmx_hw_opt_t *hw_opt,
1084 FILE *fplog, t_commrec *cr, int nfile,
1085 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
1086 gmx_bool bCompact, int nstglobalcomm,
1087 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
1088 const char *dddlb_opt, real dlb_scale,
1089 const char *ddcsx, const char *ddcsy, const char *ddcsz,
1090 const char *nbpu_opt, int nstlist_cmdline,
1091 gmx_int64_t nsteps_cmdline, int nstepout, int resetstep,
1092 int gmx_unused nmultisim, int repl_ex_nst, int repl_ex_nex,
1093 int repl_ex_seed, real pforce, real cpt_period, real max_hours,
1094 const char *deviceOptions, int imdport, unsigned long Flags)
1096 gmx_bool bForceUseGPU, bTryUseGPU;
1097 t_inputrec *inputrec;
1098 t_state *state = NULL;
1100 gmx_ddbox_t ddbox = {0};
1101 int npme_major, npme_minor;
1103 gmx_mtop_t *mtop = NULL;
1104 t_mdatoms *mdatoms = NULL;
1105 t_forcerec *fr = NULL;
1106 t_fcdata *fcd = NULL;
1107 real ewaldcoeff_q = 0;
1108 real ewaldcoeff_lj = 0;
1109 gmx_pme_t *pmedata = NULL;
1110 gmx_vsite_t *vsite = NULL;
1111 gmx_constr_t constr;
1112 int nChargePerturbed = -1, nTypePerturbed = 0, status;
1113 gmx_wallcycle_t wcycle;
1115 gmx_walltime_accounting_t walltime_accounting = NULL;
1117 gmx_int64_t reset_counters;
1118 gmx_edsam_t ed = NULL;
1119 int nthreads_pme = 1;
1120 int nthreads_pp = 1;
1121 gmx_membed_t membed = NULL;
1122 gmx_hw_info_t *hwinfo = NULL;
1123 /* The master rank decides early on bUseGPU and broadcasts this later */
1124 gmx_bool bUseGPU = FALSE;
1126 /* CAUTION: threads may be started later on in this function, so
1127 cr doesn't reflect the final parallel state right now */
1131 if (Flags & MD_APPENDFILES)
1136 bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
1137 bTryUseGPU = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
1139 /* Detect hardware, gather information. This is an operation that is
1140 * global for this process (MPI rank). */
1141 hwinfo = gmx_detect_hardware(fplog, cr, bTryUseGPU);
1147 /* Read (nearly) all data required for the simulation */
1148 read_tpx_state(ftp2fn(efTPX, nfile, fnm), inputrec, state, NULL, mtop);
1150 if (inputrec->cutoff_scheme != ecutsVERLET &&
1151 ((Flags & MD_TESTVERLET) || getenv("GMX_VERLET_SCHEME") != NULL))
1153 convert_to_verlet_scheme(fplog, inputrec, mtop, det(state->box));
1156 if (inputrec->cutoff_scheme == ecutsVERLET)
1158 /* Here the master rank decides if all ranks will use GPUs */
1159 bUseGPU = (hwinfo->gpu_info.ncuda_dev_compatible > 0 ||
1160 getenv("GMX_EMULATE_GPU") != NULL);
1162 /* TODO add GPU kernels for this and replace this check by:
1163 * (bUseGPU && (ir->vdwtype == evdwPME &&
1164 * ir->ljpme_combination_rule == eljpmeLB))
1165 * update the message text and the content of nbnxn_acceleration_supported.
1168 !nbnxn_acceleration_supported(fplog, cr, inputrec, bUseGPU))
1170 /* Fallback message printed by nbnxn_acceleration_supported */
1173 gmx_fatal(FARGS, "GPU acceleration requested, but not supported with the given input settings");
1178 prepare_verlet_scheme(fplog, cr,
1179 inputrec, nstlist_cmdline, mtop, state->box,
1184 if (nstlist_cmdline > 0)
1186 gmx_fatal(FARGS, "Can not set nstlist with the group cut-off scheme");
1189 if (hwinfo->gpu_info.ncuda_dev_compatible > 0)
1191 md_print_warn(cr, fplog,
1192 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
1193 " To use a GPU, set the mdp option: cutoff-scheme = Verlet\n"
1194 " (for quick performance testing you can use the -testverlet option)\n");
1199 gmx_fatal(FARGS, "GPU requested, but can't be used without cutoff-scheme=Verlet");
1202 #ifdef GMX_TARGET_BGQ
1203 md_print_warn(cr, fplog,
1204 "NOTE: There is no SIMD implementation of the group scheme kernels on\n"
1205 " BlueGene/Q. You will observe better performance from using the\n"
1206 " Verlet cut-off scheme.\n");
1210 if (inputrec->eI == eiSD2)
1212 md_print_warn(cr, fplog, "The stochastic dynamics integrator %s is deprecated, since\n"
1213 "it is slower than integrator %s and is slightly less accurate\n"
1214 "with constraints. Use the %s integrator.",
1215 ei_names[inputrec->eI], ei_names[eiSD1], ei_names[eiSD1]);
1219 /* Check and update the hardware options for internal consistency */
1220 check_and_update_hw_opt_1(hw_opt, SIMMASTER(cr));
1222 /* Early check for externally set process affinity. */
1223 gmx_check_thread_affinity_set(fplog, cr,
1224 hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1228 #ifdef GMX_THREAD_MPI
1229 if (cr->npmenodes > 0 && hw_opt->nthreads_tmpi <= 0)
1231 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
1235 if (hw_opt->nthreads_omp_pme != hw_opt->nthreads_omp &&
1238 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");
1242 #ifdef GMX_THREAD_MPI
1245 /* Since the master knows the cut-off scheme, update hw_opt for this.
1246 * This is done later for normal MPI and also once more with tMPI
1247 * for all tMPI ranks.
1249 check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1251 /* NOW the threads will be started: */
1252 hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
1256 if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp <= 0)
1258 hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
1261 if (hw_opt->nthreads_tmpi > 1)
1263 t_commrec *cr_old = cr;
1264 /* now start the threads. */
1265 cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
1266 oenv, bVerbose, bCompact, nstglobalcomm,
1267 ddxyz, dd_node_order, rdd, rconstr,
1268 dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
1269 nbpu_opt, nstlist_cmdline,
1270 nsteps_cmdline, nstepout, resetstep, nmultisim,
1271 repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
1272 cpt_period, max_hours, deviceOptions,
1274 /* the main thread continues here with a new cr. We don't deallocate
1275 the old cr because other threads may still be reading it. */
1278 gmx_comm("Failed to spawn threads");
1283 /* END OF CAUTION: cr is now reliable */
1285 /* g_membed initialisation *
1286 * Because we change the mtop, init_membed is called before the init_parallel *
1287 * (in case we ever want to make it run in parallel) */
1288 if (opt2bSet("-membed", nfile, fnm))
1292 fprintf(stderr, "Initializing membed");
1294 membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
1299 /* now broadcast everything to the non-master nodes/threads: */
1300 init_parallel(cr, inputrec, mtop);
1304 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
1307 /* now make sure the state is initialized and propagated */
1308 set_state_entries(state, inputrec);
1310 /* A parallel command line option consistency check that we can
1311 only do after any threads have started. */
1313 (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
1316 "The -dd or -npme option request a parallel simulation, "
1318 "but %s was compiled without threads or MPI enabled"
1320 #ifdef GMX_THREAD_MPI
1321 "but the number of threads (option -nt) is 1"
1323 "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
1330 if ((Flags & MD_RERUN) &&
1331 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1333 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
1336 if (can_use_allvsall(inputrec, TRUE, cr, fplog) && DOMAINDECOMP(cr))
1338 gmx_fatal(FARGS, "All-vs-all loops do not work with domain decomposition, use a single MPI rank");
1341 if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
1343 if (cr->npmenodes > 0)
1345 gmx_fatal_collective(FARGS, cr, NULL,
1346 "PME-only ranks are requested, but the system does not use PME for electrostatics or LJ");
1355 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1359 /* NMR restraints must be initialized before load_checkpoint,
1360 * since with time averaging the history is added to t_state.
1361 * For proper consistency check we therefore need to extend
1363 * So the PME-only nodes (if present) will also initialize
1364 * the distance restraints.
1368 /* This needs to be called before read_checkpoint to extend the state */
1369 init_disres(fplog, mtop, inputrec, cr, fcd, state, repl_ex_nst > 0);
1371 init_orires(fplog, mtop, state->x, inputrec, cr, &(fcd->orires),
1374 if (DEFORM(*inputrec))
1376 /* Store the deform reference box before reading the checkpoint */
1379 copy_mat(state->box, box);
1383 gmx_bcast(sizeof(box), box, cr);
1385 /* Because we do not have the update struct available yet
1386 * in which the reference values should be stored,
1387 * we store them temporarily in static variables.
1388 * This should be thread safe, since they are only written once
1389 * and with identical values.
1391 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
1392 deform_init_init_step_tpx = inputrec->init_step;
1393 copy_mat(box, deform_init_box_tpx);
1394 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
1397 if (opt2bSet("-cpi", nfile, fnm))
1399 /* Check if checkpoint file exists before doing continuation.
1400 * This way we can use identical input options for the first and subsequent runs...
1402 if (gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr) )
1404 load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
1406 inputrec, state, &bReadEkin,
1407 (Flags & MD_APPENDFILES),
1408 (Flags & MD_APPENDFILESSET));
1412 Flags |= MD_READ_EKIN;
1417 if (((MASTER(cr) || (Flags & MD_SEPPOT)) && (Flags & MD_APPENDFILES))
1418 #ifdef GMX_THREAD_MPI
1419 /* With thread MPI only the master node/thread exists in mdrun.c,
1420 * therefore non-master nodes need to open the "seppot" log file here.
1422 || (!MASTER(cr) && (Flags & MD_SEPPOT))
1426 gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr, !(Flags & MD_SEPPOT),
1430 /* override nsteps with value from cmdline */
1431 override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
1435 copy_mat(state->box, box);
1440 gmx_bcast(sizeof(box), box, cr);
1443 /* Essential dynamics */
1444 if (opt2bSet("-ei", nfile, fnm))
1446 /* Open input and output files, allocate space for ED data structure */
1447 ed = ed_open(mtop->natoms, &state->edsamstate, nfile, fnm, Flags, oenv, cr);
1450 if (PAR(cr) && !(EI_TPI(inputrec->eI) ||
1451 inputrec->eI == eiNM))
1453 cr->dd = init_domain_decomposition(fplog, cr, Flags, ddxyz, rdd, rconstr,
1454 dddlb_opt, dlb_scale,
1455 ddcsx, ddcsy, ddcsz,
1458 &ddbox, &npme_major, &npme_minor);
1460 make_dd_communicators(fplog, cr, dd_node_order);
1462 /* Set overallocation to avoid frequent reallocation of arrays */
1463 set_over_alloc_dd(TRUE);
1467 /* PME, if used, is done on all nodes with 1D decomposition */
1469 cr->duty = (DUTY_PP | DUTY_PME);
1473 if (inputrec->ePBC == epbcSCREW)
1476 "pbc=%s is only implemented with domain decomposition",
1477 epbc_names[inputrec->ePBC]);
1483 /* After possible communicator splitting in make_dd_communicators.
1484 * we can set up the intra/inter node communication.
1486 gmx_setup_nodecomm(fplog, cr);
1489 /* Initialize per-physical-node MPI process/thread ID and counters. */
1490 gmx_init_intranode_counters(cr);
1494 md_print_info(cr, fplog,
1495 "This is simulation %d out of %d running as a composite Gromacs\n"
1496 "multi-simulation job. Setup for this simulation:\n\n",
1497 cr->ms->sim, cr->ms->nsim);
1499 md_print_info(cr, fplog, "Using %d MPI %s\n",
1501 #ifdef GMX_THREAD_MPI
1502 cr->nnodes == 1 ? "thread" : "threads"
1504 cr->nnodes == 1 ? "process" : "processes"
1510 /* Check and update hw_opt for the cut-off scheme */
1511 check_and_update_hw_opt_2(hw_opt, inputrec->cutoff_scheme);
1513 gmx_omp_nthreads_init(fplog, cr,
1514 hwinfo->nthreads_hw_avail,
1515 hw_opt->nthreads_omp,
1516 hw_opt->nthreads_omp_pme,
1517 (cr->duty & DUTY_PP) == 0,
1518 inputrec->cutoff_scheme == ecutsVERLET);
1522 /* The master rank decided on the use of GPUs,
1523 * broadcast this information to all ranks.
1525 gmx_bcast_sim(sizeof(bUseGPU), &bUseGPU, cr);
1530 if (cr->npmenodes == -1)
1532 /* Don't automatically use PME-only nodes with GPUs */
1536 /* Select GPU id's to use */
1537 gmx_select_gpu_ids(fplog, cr, &hwinfo->gpu_info, bForceUseGPU,
1542 /* Ignore (potentially) manually selected GPUs */
1543 hw_opt->gpu_opt.ncuda_dev_use = 0;
1546 /* check consistency across ranks of things like SIMD
1547 * support and number of GPUs selected */
1548 gmx_check_hw_runconf_consistency(fplog, hwinfo, cr, hw_opt, bUseGPU);
1550 if (DOMAINDECOMP(cr))
1552 /* When we share GPUs over ranks, we need to know this for the DLB */
1553 dd_setup_dlb_resource_sharing(cr, hwinfo, hw_opt);
1556 /* getting number of PP/PME threads
1557 PME: env variable should be read only on one node to make sure it is
1558 identical everywhere;
1560 /* TODO nthreads_pp is only used for pinning threads.
1561 * This is a temporary solution until we have a hw topology library.
1563 nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
1564 nthreads_pme = gmx_omp_nthreads_get(emntPME);
1566 wcycle = wallcycle_init(fplog, resetstep, cr, nthreads_pp, nthreads_pme);
1570 /* Master synchronizes its value of reset_counters with all nodes
1571 * including PME only nodes */
1572 reset_counters = wcycle_get_reset_counters(wcycle);
1573 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1574 wcycle_set_reset_counters(wcycle, reset_counters);
1578 if (cr->duty & DUTY_PP)
1580 bcast_state(cr, state);
1582 /* Initiate forcerecord */
1584 fr->hwinfo = hwinfo;
1585 fr->gpu_opt = &hw_opt->gpu_opt;
1586 init_forcerec(fplog, oenv, fr, fcd, inputrec, mtop, cr, box,
1587 opt2fn("-table", nfile, fnm),
1588 opt2fn("-tabletf", nfile, fnm),
1589 opt2fn("-tablep", nfile, fnm),
1590 opt2fn("-tableb", nfile, fnm),
1595 /* version for PCA_NOT_READ_NODE (see md.c) */
1596 /*init_forcerec(fplog,fr,fcd,inputrec,mtop,cr,box,FALSE,
1597 "nofile","nofile","nofile","nofile",FALSE,pforce);
1599 fr->bSepDVDL = ((Flags & MD_SEPPOT) == MD_SEPPOT);
1601 /* Initialize QM-MM */
1604 init_QMMMrec(cr, mtop, inputrec, fr);
1607 /* Initialize the mdatoms structure.
1608 * mdatoms is not filled with atom data,
1609 * as this can not be done now with domain decomposition.
1611 mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1613 /* Initialize the virtual site communication */
1614 vsite = init_vsite(mtop, cr, FALSE);
1616 calc_shifts(box, fr->shift_vec);
1618 /* With periodic molecules the charge groups should be whole at start up
1619 * and the virtual sites should not be far from their proper positions.
1621 if (!inputrec->bContinuation && MASTER(cr) &&
1622 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1624 /* Make molecules whole at start of run */
1625 if (fr->ePBC != epbcNONE)
1627 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, state->x);
1631 /* Correct initial vsite positions are required
1632 * for the initial distribution in the domain decomposition
1633 * and for the initial shell prediction.
1635 construct_vsites_mtop(vsite, mtop, state->x);
1639 if (EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype))
1641 ewaldcoeff_q = fr->ewaldcoeff_q;
1642 ewaldcoeff_lj = fr->ewaldcoeff_lj;
1643 pmedata = &fr->pmedata;
1652 /* This is a PME only node */
1654 /* We don't need the state */
1657 ewaldcoeff_q = calc_ewaldcoeff_q(inputrec->rcoulomb, inputrec->ewald_rtol);
1658 ewaldcoeff_lj = calc_ewaldcoeff_lj(inputrec->rvdw, inputrec->ewald_rtol_lj);
1662 if (hw_opt->thread_affinity != threadaffOFF)
1664 /* Before setting affinity, check whether the affinity has changed
1665 * - which indicates that probably the OpenMP library has changed it
1666 * since we first checked).
1668 gmx_check_thread_affinity_set(fplog, cr,
1669 hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1671 /* Set the CPU affinity */
1672 gmx_set_thread_affinity(fplog, cr, hw_opt, hwinfo);
1675 /* Initiate PME if necessary,
1676 * either on all nodes or on dedicated PME nodes only. */
1677 if (EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype))
1681 nChargePerturbed = mdatoms->nChargePerturbed;
1682 if (EVDW_PME(inputrec->vdwtype))
1684 nTypePerturbed = mdatoms->nTypePerturbed;
1687 if (cr->npmenodes > 0)
1689 /* The PME only nodes need to know nChargePerturbed(FEP on Q) and nTypePerturbed(FEP on LJ)*/
1690 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1691 gmx_bcast_sim(sizeof(nTypePerturbed), &nTypePerturbed, cr);
1694 if (cr->duty & DUTY_PME)
1696 status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1697 mtop ? mtop->natoms : 0, nChargePerturbed, nTypePerturbed,
1698 (Flags & MD_REPRODUCIBLE), nthreads_pme);
1701 gmx_fatal(FARGS, "Error %d initializing PME", status);
1707 if (integrator[inputrec->eI].func == do_md)
1709 /* Turn on signal handling on all nodes */
1711 * (A user signal from the PME nodes (if any)
1712 * is communicated to the PP nodes.
1714 signal_handler_install();
1717 if (cr->duty & DUTY_PP)
1719 /* Assumes uniform use of the number of OpenMP threads */
1720 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntDefault));
1722 if (inputrec->ePull != epullNO)
1724 /* Initialize pull code */
1725 init_pull(fplog, inputrec, nfile, fnm, mtop, cr, oenv, inputrec->fepvals->init_lambda,
1726 EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1731 /* Initialize enforced rotation code */
1732 init_rot(fplog, inputrec, nfile, fnm, cr, state->x, box, mtop, oenv,
1736 if (inputrec->eSwapCoords != eswapNO)
1738 /* Initialize ion swapping code */
1739 init_swapcoords(fplog, bVerbose, inputrec, opt2fn_master("-swap", nfile, fnm, cr),
1740 mtop, state->x, state->box, &state->swapstate, cr, oenv, Flags);
1743 constr = init_constraints(fplog, mtop, inputrec, ed, state, cr);
1745 if (DOMAINDECOMP(cr))
1747 GMX_RELEASE_ASSERT(fr, "fr was NULL while cr->duty was DUTY_PP");
1748 dd_init_bondeds(fplog, cr->dd, mtop, vsite, inputrec,
1749 Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1751 set_dd_parameters(fplog, cr->dd, dlb_scale, inputrec, &ddbox);
1753 setup_dd_grid(fplog, cr->dd);
1756 /* Now do whatever the user wants us to do (how flexible...) */
1757 integrator[inputrec->eI].func(fplog, cr, nfile, fnm,
1758 oenv, bVerbose, bCompact,
1761 nstepout, inputrec, mtop,
1763 mdatoms, nrnb, wcycle, ed, fr,
1764 repl_ex_nst, repl_ex_nex, repl_ex_seed,
1766 cpt_period, max_hours,
1770 walltime_accounting);
1772 if (inputrec->ePull != epullNO)
1774 finish_pull(inputrec->pull);
1779 finish_rot(inputrec->rot);
1785 GMX_RELEASE_ASSERT(pmedata, "pmedata was NULL while cr->duty was not DUTY_PP");
1787 walltime_accounting = walltime_accounting_init(gmx_omp_nthreads_get(emntPME));
1788 gmx_pmeonly(*pmedata, cr, nrnb, wcycle, walltime_accounting, ewaldcoeff_q, ewaldcoeff_lj, inputrec);
1791 wallcycle_stop(wcycle, ewcRUN);
1793 /* Finish up, write some stuff
1794 * if rerunMD, don't write last frame again
1796 finish_run(fplog, cr,
1797 inputrec, nrnb, wcycle, walltime_accounting,
1798 fr != NULL && fr->nbv != NULL && fr->nbv->bUseGPU ?
1799 nbnxn_cuda_get_timings(fr->nbv->cu_nbv) : NULL,
1800 EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1803 /* Free GPU memory and context */
1804 free_gpu_resources(fr, cr);
1806 if (opt2bSet("-membed", nfile, fnm))
1811 gmx_hardware_info_free(hwinfo);
1813 /* Does what it says */
1814 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", gmx_gettime());
1815 walltime_accounting_destroy(walltime_accounting);
1817 /* Close logfile already here if we were appending to it */
1818 if (MASTER(cr) && (Flags & MD_APPENDFILES))
1820 gmx_log_close(fplog);
1823 rc = (int)gmx_get_stop_condition();
1825 #ifdef GMX_THREAD_MPI
1826 /* we need to join all threads. The sub-threads join when they
1827 exit this function, but the master thread needs to be told to
1829 if (PAR(cr) && MASTER(cr))