1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
53 #include "md_logging.h"
54 #include "md_support.h"
57 #include "pull_rotation.h"
70 #include "checkpoint.h"
71 #include "mtop_util.h"
72 #include "sighandler.h"
75 #include "gmx_detect_hardware.h"
76 #include "gmx_omp_nthreads.h"
77 #include "pull_rotation.h"
78 #include "calc_verletbuf.h"
79 #include "../mdlib/nbnxn_search.h"
80 #include "../mdlib/nbnxn_consts.h"
81 #include "gmx_fatal_collective.h"
85 #include "gmx_thread_affinity.h"
87 #include "gromacs/utility/gmxmpi.h"
93 #include "gpu_utils.h"
94 #include "nbnxn_cuda_data_mgmt.h"
97 gmx_integrator_t *func;
100 /* The array should match the eI array in include/types/enums.h */
101 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}};
103 gmx_large_int_t deform_init_init_step_tpx;
104 matrix deform_init_box_tpx;
105 #ifdef GMX_THREAD_MPI
106 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
110 #ifdef GMX_THREAD_MPI
111 struct mdrunner_arglist
113 gmx_hw_opt_t *hw_opt;
126 const char *dddlb_opt;
131 const char *nbpu_opt;
132 gmx_large_int_t nsteps_cmdline;
142 const char *deviceOptions;
144 int ret; /* return value */
148 /* The function used for spawning threads. Extracts the mdrunner()
149 arguments from its one argument and calls mdrunner(), after making
151 static void mdrunner_start_fn(void *arg)
153 struct mdrunner_arglist *mda = (struct mdrunner_arglist*)arg;
154 struct mdrunner_arglist mc = *mda; /* copy the arg list to make sure
155 that it's thread-local. This doesn't
156 copy pointed-to items, of course,
157 but those are all const. */
158 t_commrec *cr; /* we need a local version of this */
162 fnm = dup_tfn(mc.nfile, mc.fnm);
164 cr = init_par_threads(mc.cr);
171 mda->ret = mdrunner(mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
172 mc.bVerbose, mc.bCompact, mc.nstglobalcomm,
173 mc.ddxyz, mc.dd_node_order, mc.rdd,
174 mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
175 mc.ddcsx, mc.ddcsy, mc.ddcsz,
177 mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
178 mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce,
179 mc.cpt_period, mc.max_hours, mc.deviceOptions, mc.Flags);
182 /* called by mdrunner() to start a specific number of threads (including
183 the main thread) for thread-parallel runs. This in turn calls mdrunner()
185 All options besides nthreads are the same as for mdrunner(). */
186 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
187 FILE *fplog, t_commrec *cr, int nfile,
188 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
189 gmx_bool bCompact, int nstglobalcomm,
190 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
191 const char *dddlb_opt, real dlb_scale,
192 const char *ddcsx, const char *ddcsy, const char *ddcsz,
193 const char *nbpu_opt,
194 gmx_large_int_t nsteps_cmdline,
195 int nstepout, int resetstep,
196 int nmultisim, int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
197 real pforce, real cpt_period, real max_hours,
198 const char *deviceOptions, unsigned long Flags)
201 struct mdrunner_arglist *mda;
202 t_commrec *crn; /* the new commrec */
205 /* first check whether we even need to start tMPI */
206 if (hw_opt->nthreads_tmpi < 2)
211 /* a few small, one-time, almost unavoidable memory leaks: */
213 fnmn = dup_tfn(nfile, fnm);
215 /* fill the data structure to pass as void pointer to thread start fn */
216 mda->hw_opt = hw_opt;
222 mda->bVerbose = bVerbose;
223 mda->bCompact = bCompact;
224 mda->nstglobalcomm = nstglobalcomm;
225 mda->ddxyz[XX] = ddxyz[XX];
226 mda->ddxyz[YY] = ddxyz[YY];
227 mda->ddxyz[ZZ] = ddxyz[ZZ];
228 mda->dd_node_order = dd_node_order;
230 mda->rconstr = rconstr;
231 mda->dddlb_opt = dddlb_opt;
232 mda->dlb_scale = dlb_scale;
236 mda->nbpu_opt = nbpu_opt;
237 mda->nsteps_cmdline = nsteps_cmdline;
238 mda->nstepout = nstepout;
239 mda->resetstep = resetstep;
240 mda->nmultisim = nmultisim;
241 mda->repl_ex_nst = repl_ex_nst;
242 mda->repl_ex_nex = repl_ex_nex;
243 mda->repl_ex_seed = repl_ex_seed;
244 mda->pforce = pforce;
245 mda->cpt_period = cpt_period;
246 mda->max_hours = max_hours;
247 mda->deviceOptions = deviceOptions;
250 /* now spawn new threads that start mdrunner_start_fn(), while
251 the main thread returns, we set thread affinity later */
252 ret = tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi, TMPI_AFFINITY_NONE,
253 mdrunner_start_fn, (void*)(mda) );
254 if (ret != TMPI_SUCCESS)
259 /* make a new comm_rec to reflect the new situation */
260 crn = init_par_threads(cr);
265 static int get_tmpi_omp_thread_division(const gmx_hw_info_t *hwinfo,
266 const gmx_hw_opt_t *hw_opt,
272 /* There are no separate PME nodes here, as we ensured in
273 * check_and_update_hw_opt that nthreads_tmpi>0 with PME nodes
274 * and a conditional ensures we would not have ended up here.
275 * Note that separate PME nodes might be switched on later.
279 nthreads_tmpi = ngpu;
280 if (nthreads_tot > 0 && nthreads_tot < nthreads_tmpi)
282 nthreads_tmpi = nthreads_tot;
285 else if (hw_opt->nthreads_omp > 0)
287 /* Here we could oversubscribe, when we do, we issue a warning later */
288 nthreads_tmpi = max(1, nthreads_tot/hw_opt->nthreads_omp);
292 /* TODO choose nthreads_omp based on hardware topology
293 when we have a hardware topology detection library */
294 /* In general, when running up to 4 threads, OpenMP should be faster.
295 * Note: on AMD Bulldozer we should avoid running OpenMP over two dies.
296 * On Intel>=Nehalem running OpenMP on a single CPU is always faster,
297 * even on two CPUs it's usually faster (but with many OpenMP threads
298 * it could be faster not to use HT, currently we always use HT).
299 * On Nehalem/Westmere we want to avoid running 16 threads over
300 * two CPUs with HT, so we need a limit<16; thus we use 12.
301 * A reasonable limit for Intel Sandy and Ivy bridge,
302 * not knowing the topology, is 16 threads.
304 const int nthreads_omp_always_faster = 4;
305 const int nthreads_omp_always_faster_Nehalem = 12;
306 const int nthreads_omp_always_faster_SandyBridge = 16;
307 const int first_model_Nehalem = 0x1A;
308 const int first_model_SandyBridge = 0x2A;
309 gmx_bool bIntel_Family6;
312 (gmx_cpuid_vendor(hwinfo->cpuid_info) == GMX_CPUID_VENDOR_INTEL &&
313 gmx_cpuid_family(hwinfo->cpuid_info) == 6);
315 if (nthreads_tot <= nthreads_omp_always_faster ||
317 ((gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_Nehalem && nthreads_tot <= nthreads_omp_always_faster_Nehalem) ||
318 (gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_SandyBridge && nthreads_tot <= nthreads_omp_always_faster_SandyBridge))))
320 /* Use pure OpenMP parallelization */
325 /* Don't use OpenMP parallelization */
326 nthreads_tmpi = nthreads_tot;
330 return nthreads_tmpi;
334 /* Get the number of threads to use for thread-MPI based on how many
335 * were requested, which algorithms we're using,
336 * and how many particles there are.
337 * At the point we have already called check_and_update_hw_opt.
338 * Thus all options should be internally consistent and consistent
339 * with the hardware, except that ntmpi could be larger than #GPU.
341 static int get_nthreads_mpi(gmx_hw_info_t *hwinfo,
342 gmx_hw_opt_t *hw_opt,
343 t_inputrec *inputrec, gmx_mtop_t *mtop,
347 int nthreads_hw, nthreads_tot_max, nthreads_tmpi, nthreads_new, ngpu;
348 int min_atoms_per_mpi_thread;
353 if (hw_opt->nthreads_tmpi > 0)
355 /* Trivial, return right away */
356 return hw_opt->nthreads_tmpi;
359 nthreads_hw = hwinfo->nthreads_hw_avail;
361 /* How many total (#tMPI*#OpenMP) threads can we start? */
362 if (hw_opt->nthreads_tot > 0)
364 nthreads_tot_max = hw_opt->nthreads_tot;
368 nthreads_tot_max = nthreads_hw;
371 bCanUseGPU = (inputrec->cutoff_scheme == ecutsVERLET && hwinfo->bCanUseGPU);
374 ngpu = hwinfo->gpu_info.ncuda_dev_use;
382 get_tmpi_omp_thread_division(hwinfo, hw_opt, nthreads_tot_max, ngpu);
384 if (inputrec->eI == eiNM || EI_TPI(inputrec->eI))
386 /* Steps are divided over the nodes iso splitting the atoms */
387 min_atoms_per_mpi_thread = 0;
393 min_atoms_per_mpi_thread = MIN_ATOMS_PER_GPU;
397 min_atoms_per_mpi_thread = MIN_ATOMS_PER_MPI_THREAD;
401 /* Check if an algorithm does not support parallel simulation. */
402 if (nthreads_tmpi != 1 &&
403 ( inputrec->eI == eiLBFGS ||
404 inputrec->coulombtype == eelEWALD ) )
408 md_print_warn(cr, fplog, "The integration or electrostatics algorithm doesn't support parallel runs. Using a single thread-MPI thread.\n");
409 if (hw_opt->nthreads_tmpi > nthreads_tmpi)
411 gmx_fatal(FARGS, "You asked for more than 1 thread-MPI thread, but an algorithm doesn't support that");
414 else if (mtop->natoms/nthreads_tmpi < min_atoms_per_mpi_thread)
416 /* the thread number was chosen automatically, but there are too many
417 threads (too few atoms per thread) */
418 nthreads_new = max(1, mtop->natoms/min_atoms_per_mpi_thread);
420 /* Avoid partial use of Hyper-Threading */
421 if (gmx_cpuid_x86_smt(hwinfo->cpuid_info) == GMX_CPUID_X86_SMT_ENABLED &&
422 nthreads_new > nthreads_hw/2 && nthreads_new < nthreads_hw)
424 nthreads_new = nthreads_hw/2;
427 /* Avoid large prime numbers in the thread count */
428 if (nthreads_new >= 6)
430 /* Use only 6,8,10 with additional factors of 2 */
434 while (3*fac*2 <= nthreads_new)
439 nthreads_new = (nthreads_new/fac)*fac;
444 if (nthreads_new == 5)
450 nthreads_tmpi = nthreads_new;
452 fprintf(stderr, "\n");
453 fprintf(stderr, "NOTE: Parallelization is limited by the small number of atoms,\n");
454 fprintf(stderr, " only starting %d thread-MPI threads.\n", nthreads_tmpi);
455 fprintf(stderr, " You can use the -nt and/or -ntmpi option to optimize the number of threads.\n\n");
458 return nthreads_tmpi;
460 #endif /* GMX_THREAD_MPI */
463 /* Environment variable for setting nstlist */
464 static const char* NSTLIST_ENVVAR = "GMX_NSTLIST";
465 /* Try to increase nstlist when using a GPU with nstlist less than this */
466 static const int NSTLIST_GPU_ENOUGH = 20;
467 /* Increase nstlist until the non-bonded cost increases more than this factor */
468 static const float NBNXN_GPU_LIST_OK_FAC = 1.25;
469 /* Don't increase nstlist beyond a non-bonded cost increases of this factor */
470 static const float NBNXN_GPU_LIST_MAX_FAC = 1.40;
472 /* Try to increase nstlist when running on a GPU */
473 static void increase_nstlist(FILE *fp, t_commrec *cr,
474 t_inputrec *ir, const gmx_mtop_t *mtop, matrix box)
477 int nstlist_orig, nstlist_prev;
478 verletbuf_list_setup_t ls;
479 real rlist_inc, rlist_ok, rlist_max, rlist_new, rlist_prev;
482 gmx_bool bBox, bDD, bCont;
483 const char *nstl_fmt = "\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";
484 const char *vbd_err = "Can not increase nstlist for GPU run because verlet-buffer-drift is not set or used";
485 const char *box_err = "Can not increase nstlist for GPU run because the box is too small";
486 const char *dd_err = "Can not increase nstlist for GPU run because of domain decomposition limitations";
489 /* Number of + nstlist alternative values to try when switching */
490 const int nstl[] = { 20, 25, 40, 50 };
491 #define NNSTL sizeof(nstl)/sizeof(nstl[0])
493 env = getenv(NSTLIST_ENVVAR);
498 fprintf(fp, nstl_fmt, ir->nstlist);
502 if (ir->verletbuf_drift == 0)
504 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");
507 if (ir->verletbuf_drift < 0)
511 fprintf(stderr, "%s\n", vbd_err);
515 fprintf(fp, "%s\n", vbd_err);
521 nstlist_orig = ir->nstlist;
524 sprintf(buf, "Getting nstlist from environment variable GMX_NSTLIST=%s", env);
527 fprintf(stderr, "%s\n", buf);
531 fprintf(fp, "%s\n", buf);
533 sscanf(env, "%d", &ir->nstlist);
536 verletbuf_get_list_setup(TRUE, &ls);
538 /* Allow rlist to make the list double the size of the cut-off sphere */
539 rlist_inc = nbnxn_get_rlist_effective_inc(NBNXN_GPU_CLUSTER_SIZE, mtop->natoms/det(box));
540 rlist_ok = (max(ir->rvdw, ir->rcoulomb) + rlist_inc)*pow(NBNXN_GPU_LIST_OK_FAC, 1.0/3.0) - rlist_inc;
541 rlist_max = (max(ir->rvdw, ir->rcoulomb) + rlist_inc)*pow(NBNXN_GPU_LIST_MAX_FAC, 1.0/3.0) - rlist_inc;
544 fprintf(debug, "GPU nstlist tuning: rlist_inc %.3f rlist_max %.3f\n",
545 rlist_inc, rlist_max);
549 nstlist_prev = nstlist_orig;
550 rlist_prev = ir->rlist;
555 ir->nstlist = nstl[i];
558 /* Set the pair-list buffer size in ir */
559 calc_verlet_buffer_size(mtop, det(box), ir, ir->verletbuf_drift, &ls,
562 /* Does rlist fit in the box? */
563 bBox = (sqr(rlist_new) < max_cutoff2(ir->ePBC, box));
565 if (bBox && DOMAINDECOMP(cr))
567 /* Check if rlist fits in the domain decomposition */
568 if (inputrec2nboundeddim(ir) < DIM)
570 gmx_incons("Changing nstlist with domain decomposition and unbounded dimensions is not implemented yet");
572 copy_mat(box, state_tmp.box);
573 bDD = change_dd_cutoff(cr, &state_tmp, ir, rlist_new);
580 if (bBox && bDD && rlist_new <= rlist_max)
582 /* Increase nstlist */
583 nstlist_prev = ir->nstlist;
584 rlist_prev = rlist_new;
585 bCont = (i+1 < NNSTL && rlist_new < rlist_ok);
589 /* Stick with the previous nstlist */
590 ir->nstlist = nstlist_prev;
591 rlist_new = rlist_prev;
603 gmx_warning(!bBox ? box_err : dd_err);
606 fprintf(fp, "\n%s\n", bBox ? box_err : dd_err);
608 ir->nstlist = nstlist_orig;
610 else if (ir->nstlist != nstlist_orig || rlist_new != ir->rlist)
612 sprintf(buf, "Changing nstlist from %d to %d, rlist from %g to %g",
613 nstlist_orig, ir->nstlist,
614 ir->rlist, rlist_new);
617 fprintf(stderr, "%s\n\n", buf);
621 fprintf(fp, "%s\n\n", buf);
623 ir->rlist = rlist_new;
624 ir->rlistlong = rlist_new;
628 static void prepare_verlet_scheme(FILE *fplog,
629 gmx_hw_info_t *hwinfo,
631 const char *nbpu_opt,
633 const gmx_mtop_t *mtop,
637 /* Here we only check for GPU usage on the MPI master process,
638 * as here we don't know how many GPUs we will use yet.
639 * We check for a GPU on all processes later.
641 *bUseGPU = hwinfo->bCanUseGPU || (getenv("GMX_EMULATE_GPU") != NULL);
643 if (ir->verletbuf_drift > 0)
645 /* Update the Verlet buffer size for the current run setup */
646 verletbuf_list_setup_t ls;
649 /* Here we assume CPU acceleration is on. But as currently
650 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
651 * and 4x2 gives a larger buffer than 4x4, this is ok.
653 verletbuf_get_list_setup(*bUseGPU, &ls);
655 calc_verlet_buffer_size(mtop, det(box), ir,
656 ir->verletbuf_drift, &ls,
658 if (rlist_new != ir->rlist)
662 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
663 ir->rlist, rlist_new,
664 ls.cluster_size_i, ls.cluster_size_j);
666 ir->rlist = rlist_new;
667 ir->rlistlong = rlist_new;
671 /* With GPU or emulation we should check nstlist for performance */
672 if ((EI_DYNAMICS(ir->eI) &&
674 ir->nstlist < NSTLIST_GPU_ENOUGH) ||
675 getenv(NSTLIST_ENVVAR) != NULL)
677 /* Choose a better nstlist */
678 increase_nstlist(fplog, cr, ir, mtop, box);
682 static void convert_to_verlet_scheme(FILE *fplog,
684 gmx_mtop_t *mtop, real box_vol)
686 char *conv_mesg = "Converting input file with group cut-off scheme to the Verlet cut-off scheme";
688 md_print_warn(NULL, fplog, "%s\n", conv_mesg);
690 ir->cutoff_scheme = ecutsVERLET;
691 ir->verletbuf_drift = 0.005;
693 if (ir->rcoulomb != ir->rvdw)
695 gmx_fatal(FARGS, "The VdW and Coulomb cut-offs are different, whereas the Verlet scheme only supports equal cut-offs");
698 if (ir->vdwtype == evdwUSER || EEL_USER(ir->coulombtype))
700 gmx_fatal(FARGS, "User non-bonded potentials are not (yet) supported with the Verlet scheme");
702 else if (EVDW_SWITCHED(ir->vdwtype) || EEL_SWITCHED(ir->coulombtype))
704 md_print_warn(NULL, fplog, "Converting switched or shifted interactions to a shifted potential (without force shift), this will lead to slightly different interaction potentials");
706 if (EVDW_SWITCHED(ir->vdwtype))
708 ir->vdwtype = evdwCUT;
710 if (EEL_SWITCHED(ir->coulombtype))
712 if (EEL_FULL(ir->coulombtype))
714 /* With full electrostatic only PME can be switched */
715 ir->coulombtype = eelPME;
719 md_print_warn(NULL, fplog, "NOTE: Replacing %s electrostatics with reaction-field with epsilon-rf=inf\n", eel_names[ir->coulombtype]);
720 ir->coulombtype = eelRF;
721 ir->epsilon_rf = 0.0;
725 /* We set the target energy drift to a small number.
726 * Note that this is only for testing. For production the user
727 * should think about this and set the mdp options.
729 ir->verletbuf_drift = 1e-4;
732 if (inputrec2nboundeddim(ir) != 3)
734 gmx_fatal(FARGS, "Can only convert old tpr files to the Verlet cut-off scheme with 3D pbc");
737 if (ir->efep != efepNO || ir->implicit_solvent != eisNO)
739 gmx_fatal(FARGS, "Will not convert old tpr files to the Verlet cut-off scheme with free-energy calculations or implicit solvent");
742 if (EI_DYNAMICS(ir->eI) && !(EI_MD(ir->eI) && ir->etc == etcNO))
744 verletbuf_list_setup_t ls;
746 verletbuf_get_list_setup(FALSE, &ls);
747 calc_verlet_buffer_size(mtop, box_vol, ir, ir->verletbuf_drift, &ls,
752 ir->verletbuf_drift = -1;
753 ir->rlist = 1.05*max(ir->rvdw, ir->rcoulomb);
756 gmx_mtop_remove_chargegroups(mtop);
759 static void check_and_update_hw_opt(gmx_hw_opt_t *hw_opt,
761 gmx_bool bIsSimMaster)
763 gmx_omp_nthreads_read_env(&hw_opt->nthreads_omp, bIsSimMaster);
765 #ifndef GMX_THREAD_MPI
766 if (hw_opt->nthreads_tot > 0)
768 gmx_fatal(FARGS, "Setting the total number of threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
770 if (hw_opt->nthreads_tmpi > 0)
772 gmx_fatal(FARGS, "Setting the number of thread-MPI threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
776 if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp_pme <= 0)
778 /* We have the same number of OpenMP threads for PP and PME processes,
779 * thus we can perform several consistency checks.
781 if (hw_opt->nthreads_tmpi > 0 &&
782 hw_opt->nthreads_omp > 0 &&
783 hw_opt->nthreads_tot != hw_opt->nthreads_tmpi*hw_opt->nthreads_omp)
785 gmx_fatal(FARGS, "The total number of threads requested (%d) does not match the thread-MPI threads (%d) times the OpenMP threads (%d) requested",
786 hw_opt->nthreads_tot, hw_opt->nthreads_tmpi, hw_opt->nthreads_omp);
789 if (hw_opt->nthreads_tmpi > 0 &&
790 hw_opt->nthreads_tot % hw_opt->nthreads_tmpi != 0)
792 gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of thread-MPI threads requested (%d)",
793 hw_opt->nthreads_tot, hw_opt->nthreads_tmpi);
796 if (hw_opt->nthreads_omp > 0 &&
797 hw_opt->nthreads_tot % hw_opt->nthreads_omp != 0)
799 gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of OpenMP threads requested (%d)",
800 hw_opt->nthreads_tot, hw_opt->nthreads_omp);
803 if (hw_opt->nthreads_tmpi > 0 &&
804 hw_opt->nthreads_omp <= 0)
806 hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
811 if (hw_opt->nthreads_omp > 1)
813 gmx_fatal(FARGS, "OpenMP threads are requested, but Gromacs was compiled without OpenMP support");
817 if (cutoff_scheme == ecutsGROUP)
819 /* We only have OpenMP support for PME only nodes */
820 if (hw_opt->nthreads_omp > 1)
822 gmx_fatal(FARGS, "OpenMP threads have been requested with cut-off scheme %s, but these are only supported with cut-off scheme %s",
823 ecutscheme_names[cutoff_scheme],
824 ecutscheme_names[ecutsVERLET]);
826 hw_opt->nthreads_omp = 1;
829 if (hw_opt->nthreads_omp_pme > 0 && hw_opt->nthreads_omp <= 0)
831 gmx_fatal(FARGS, "You need to specify -ntomp in addition to -ntomp_pme");
834 if (hw_opt->nthreads_tot == 1)
836 hw_opt->nthreads_tmpi = 1;
838 if (hw_opt->nthreads_omp > 1)
840 gmx_fatal(FARGS, "You requested %d OpenMP threads with %d total threads",
841 hw_opt->nthreads_tmpi, hw_opt->nthreads_tot);
843 hw_opt->nthreads_omp = 1;
846 if (hw_opt->nthreads_omp_pme <= 0 && hw_opt->nthreads_omp > 0)
848 hw_opt->nthreads_omp_pme = hw_opt->nthreads_omp;
853 fprintf(debug, "hw_opt: nt %d ntmpi %d ntomp %d ntomp_pme %d gpu_id '%s'\n",
854 hw_opt->nthreads_tot,
855 hw_opt->nthreads_tmpi,
856 hw_opt->nthreads_omp,
857 hw_opt->nthreads_omp_pme,
858 hw_opt->gpu_id != NULL ? hw_opt->gpu_id : "");
864 /* Override the value in inputrec with value passed on the command line (if any) */
865 static void override_nsteps_cmdline(FILE *fplog,
866 gmx_large_int_t nsteps_cmdline,
870 char sbuf[STEPSTRSIZE];
875 /* override with anything else than the default -2 */
876 if (nsteps_cmdline > -2)
880 ir->nsteps = nsteps_cmdline;
881 if (EI_DYNAMICS(ir->eI))
883 sprintf(stmp, "Overriding nsteps with value passed on the command line: %s steps, %.3f ps",
884 gmx_step_str(nsteps_cmdline, sbuf),
885 nsteps_cmdline*ir->delta_t);
889 sprintf(stmp, "Overriding nsteps with value passed on the command line: %s steps",
890 gmx_step_str(nsteps_cmdline, sbuf));
893 md_print_warn(cr, fplog, "%s\n", stmp);
897 /* Data structure set by SIMMASTER which needs to be passed to all nodes
898 * before the other nodes have read the tpx file and called gmx_detect_hardware.
901 int cutoff_scheme; /* The cutoff scheme from inputrec_t */
902 gmx_bool bUseGPU; /* Use GPU or GPU emulation */
905 int mdrunner(gmx_hw_opt_t *hw_opt,
906 FILE *fplog, t_commrec *cr, int nfile,
907 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
908 gmx_bool bCompact, int nstglobalcomm,
909 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
910 const char *dddlb_opt, real dlb_scale,
911 const char *ddcsx, const char *ddcsy, const char *ddcsz,
912 const char *nbpu_opt,
913 gmx_large_int_t nsteps_cmdline, int nstepout, int resetstep,
914 int nmultisim, int repl_ex_nst, int repl_ex_nex,
915 int repl_ex_seed, real pforce, real cpt_period, real max_hours,
916 const char *deviceOptions, unsigned long Flags)
918 gmx_bool bForceUseGPU, bTryUseGPU;
919 double nodetime = 0, realtime;
920 t_inputrec *inputrec;
921 t_state *state = NULL;
923 gmx_ddbox_t ddbox = {0};
924 int npme_major, npme_minor;
927 gmx_mtop_t *mtop = NULL;
928 t_mdatoms *mdatoms = NULL;
929 t_forcerec *fr = NULL;
930 t_fcdata *fcd = NULL;
932 gmx_pme_t *pmedata = NULL;
933 gmx_vsite_t *vsite = NULL;
935 int i, m, nChargePerturbed = -1, status, nalloc;
937 gmx_wallcycle_t wcycle;
938 gmx_bool bReadRNG, bReadEkin;
940 gmx_runtime_t runtime;
942 gmx_large_int_t reset_counters;
943 gmx_edsam_t ed = NULL;
944 t_commrec *cr_old = cr;
945 int nthreads_pme = 1;
947 gmx_membed_t membed = NULL;
948 gmx_hw_info_t *hwinfo = NULL;
949 master_inf_t minf = {-1, FALSE};
951 /* CAUTION: threads may be started later on in this function, so
952 cr doesn't reflect the final parallel state right now */
956 if (Flags & MD_APPENDFILES)
961 bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
962 bTryUseGPU = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
967 /* Read (nearly) all data required for the simulation */
968 read_tpx_state(ftp2fn(efTPX, nfile, fnm), inputrec, state, NULL, mtop);
970 if (inputrec->cutoff_scheme != ecutsVERLET &&
971 ((Flags & MD_TESTVERLET) || getenv("GMX_VERLET_SCHEME") != NULL))
973 convert_to_verlet_scheme(fplog, inputrec, mtop, det(state->box));
976 /* Detect hardware, gather information. With tMPI only thread 0 does it
977 * and after threads are started broadcasts hwinfo around. */
979 gmx_detect_hardware(fplog, hwinfo, cr,
980 bForceUseGPU, bTryUseGPU, hw_opt->gpu_id);
982 minf.cutoff_scheme = inputrec->cutoff_scheme;
983 minf.bUseGPU = FALSE;
985 if (inputrec->cutoff_scheme == ecutsVERLET)
987 prepare_verlet_scheme(fplog, hwinfo, cr, nbpu_opt,
988 inputrec, mtop, state->box,
991 else if (hwinfo->bCanUseGPU)
993 md_print_warn(cr, fplog,
994 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
995 " To use a GPU, set the mdp option: cutoff-scheme = Verlet\n"
996 " (for quick performance testing you can use the -testverlet option)\n");
1000 gmx_fatal(FARGS, "GPU requested, but can't be used without cutoff-scheme=Verlet");
1004 #ifndef GMX_THREAD_MPI
1007 gmx_bcast_sim(sizeof(minf), &minf, cr);
1010 if (minf.bUseGPU && cr->npmenodes == -1)
1012 /* Don't automatically use PME-only nodes with GPUs */
1016 /* Check for externally set OpenMP affinity and turn off internal
1017 * pinning if any is found. We need to do this check early to tell
1018 * thread-MPI whether it should do pinning when spawning threads.
1019 * TODO: the above no longer holds, we should move these checks down
1021 gmx_omp_check_thread_affinity(fplog, cr, hw_opt);
1023 #ifdef GMX_THREAD_MPI
1024 /* With thread-MPI inputrec is only set here on the master thread */
1028 check_and_update_hw_opt(hw_opt, minf.cutoff_scheme, SIMMASTER(cr));
1030 #ifdef GMX_THREAD_MPI
1031 /* Early check for externally set process affinity. Can't do over all
1032 * MPI processes because hwinfo is not available everywhere, but with
1033 * thread-MPI it's needed as pinning might get turned off which needs
1034 * to be known before starting thread-MPI. */
1035 gmx_check_thread_affinity_set(fplog,
1037 hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1040 #ifdef GMX_THREAD_MPI
1041 if (cr->npmenodes > 0 && hw_opt->nthreads_tmpi <= 0)
1043 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME nodes");
1047 if (hw_opt->nthreads_omp_pme != hw_opt->nthreads_omp &&
1050 gmx_fatal(FARGS, "You need to explicitly specify the number of PME nodes (-npme) when using different number of OpenMP threads for PP and PME nodes");
1054 #ifdef GMX_THREAD_MPI
1057 /* NOW the threads will be started: */
1058 hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
1062 if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp <= 0)
1064 hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
1067 if (hw_opt->nthreads_tmpi > 1)
1069 /* now start the threads. */
1070 cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
1071 oenv, bVerbose, bCompact, nstglobalcomm,
1072 ddxyz, dd_node_order, rdd, rconstr,
1073 dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
1075 nsteps_cmdline, nstepout, resetstep, nmultisim,
1076 repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
1077 cpt_period, max_hours, deviceOptions,
1079 /* the main thread continues here with a new cr. We don't deallocate
1080 the old cr because other threads may still be reading it. */
1083 gmx_comm("Failed to spawn threads");
1088 /* END OF CAUTION: cr is now reliable */
1090 /* g_membed initialisation *
1091 * Because we change the mtop, init_membed is called before the init_parallel *
1092 * (in case we ever want to make it run in parallel) */
1093 if (opt2bSet("-membed", nfile, fnm))
1097 fprintf(stderr, "Initializing membed");
1099 membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
1104 /* now broadcast everything to the non-master nodes/threads: */
1105 init_parallel(fplog, cr, inputrec, mtop);
1107 /* This check needs to happen after get_nthreads_mpi() */
1108 if (inputrec->cutoff_scheme == ecutsVERLET && (Flags & MD_PARTDEC))
1110 gmx_fatal_collective(FARGS, cr, NULL,
1111 "The Verlet cut-off scheme is not supported with particle decomposition.\n"
1112 "You can achieve the same effect as particle decomposition by running in parallel using only OpenMP threads.");
1117 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
1120 #if defined GMX_THREAD_MPI
1121 /* With tMPI we detected on thread 0 and we'll just pass the hwinfo pointer
1122 * to the other threads -- slightly uncool, but works fine, just need to
1123 * make sure that the data doesn't get freed twice. */
1130 gmx_bcast(sizeof(&hwinfo), &hwinfo, cr);
1133 if (PAR(cr) && !SIMMASTER(cr))
1135 /* now we have inputrec on all nodes, can run the detection */
1136 /* TODO: perhaps it's better to propagate within a node instead? */
1138 gmx_detect_hardware(fplog, hwinfo, cr,
1139 bForceUseGPU, bTryUseGPU, hw_opt->gpu_id);
1142 /* Now do the affinity check with MPI/no-MPI (done earlier with thread-MPI). */
1143 gmx_check_thread_affinity_set(fplog, cr,
1144 hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1147 /* now make sure the state is initialized and propagated */
1148 set_state_entries(state, inputrec, cr->nnodes);
1150 /* A parallel command line option consistency check that we can
1151 only do after any threads have started. */
1153 (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
1156 "The -dd or -npme option request a parallel simulation, "
1158 "but %s was compiled without threads or MPI enabled"
1160 #ifdef GMX_THREAD_MPI
1161 "but the number of threads (option -nt) is 1"
1163 "but %s was not started through mpirun/mpiexec or only one process was requested through mpirun/mpiexec"
1170 if ((Flags & MD_RERUN) &&
1171 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1173 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
1176 if (can_use_allvsall(inputrec, mtop, TRUE, cr, fplog) && PAR(cr))
1178 /* Simple neighbour searching and (also?) all-vs-all loops
1179 * do not work with domain decomposition. */
1180 Flags |= MD_PARTDEC;
1183 if (!EEL_PME(inputrec->coulombtype) || (Flags & MD_PARTDEC))
1185 if (cr->npmenodes > 0)
1187 if (!EEL_PME(inputrec->coulombtype))
1189 gmx_fatal_collective(FARGS, cr, NULL,
1190 "PME nodes are requested, but the system does not use PME electrostatics");
1192 if (Flags & MD_PARTDEC)
1194 gmx_fatal_collective(FARGS, cr, NULL,
1195 "PME nodes are requested, but particle decomposition does not support separate PME nodes");
1203 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1206 /* NMR restraints must be initialized before load_checkpoint,
1207 * since with time averaging the history is added to t_state.
1208 * For proper consistency check we therefore need to extend
1210 * So the PME-only nodes (if present) will also initialize
1211 * the distance restraints.
1215 /* This needs to be called before read_checkpoint to extend the state */
1216 init_disres(fplog, mtop, inputrec, cr, Flags & MD_PARTDEC, fcd, state, repl_ex_nst > 0);
1218 if (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0)
1220 if (PAR(cr) && !(Flags & MD_PARTDEC))
1222 gmx_fatal(FARGS, "Orientation restraints do not work (yet) with domain decomposition, use particle decomposition (mdrun option -pd)");
1224 /* Orientation restraints */
1227 init_orires(fplog, mtop, state->x, inputrec, cr->ms, &(fcd->orires),
1232 if (DEFORM(*inputrec))
1234 /* Store the deform reference box before reading the checkpoint */
1237 copy_mat(state->box, box);
1241 gmx_bcast(sizeof(box), box, cr);
1243 /* Because we do not have the update struct available yet
1244 * in which the reference values should be stored,
1245 * we store them temporarily in static variables.
1246 * This should be thread safe, since they are only written once
1247 * and with identical values.
1249 #ifdef GMX_THREAD_MPI
1250 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
1252 deform_init_init_step_tpx = inputrec->init_step;
1253 copy_mat(box, deform_init_box_tpx);
1254 #ifdef GMX_THREAD_MPI
1255 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
1259 if (opt2bSet("-cpi", nfile, fnm))
1261 /* Check if checkpoint file exists before doing continuation.
1262 * This way we can use identical input options for the first and subsequent runs...
1264 if (gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr) )
1266 load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
1267 cr, Flags & MD_PARTDEC, ddxyz,
1268 inputrec, state, &bReadRNG, &bReadEkin,
1269 (Flags & MD_APPENDFILES),
1270 (Flags & MD_APPENDFILESSET));
1274 Flags |= MD_READ_RNG;
1278 Flags |= MD_READ_EKIN;
1283 if (((MASTER(cr) || (Flags & MD_SEPPOT)) && (Flags & MD_APPENDFILES))
1284 #ifdef GMX_THREAD_MPI
1285 /* With thread MPI only the master node/thread exists in mdrun.c,
1286 * therefore non-master nodes need to open the "seppot" log file here.
1288 || (!MASTER(cr) && (Flags & MD_SEPPOT))
1292 gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr, !(Flags & MD_SEPPOT),
1296 /* override nsteps with value from cmdline */
1297 override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
1301 copy_mat(state->box, box);
1306 gmx_bcast(sizeof(box), box, cr);
1309 /* Essential dynamics */
1310 if (opt2bSet("-ei", nfile, fnm))
1312 /* Open input and output files, allocate space for ED data structure */
1313 ed = ed_open(mtop->natoms, &state->edsamstate, nfile, fnm, Flags, oenv, cr);
1316 if (PAR(cr) && !((Flags & MD_PARTDEC) ||
1317 EI_TPI(inputrec->eI) ||
1318 inputrec->eI == eiNM))
1320 cr->dd = init_domain_decomposition(fplog, cr, Flags, ddxyz, rdd, rconstr,
1321 dddlb_opt, dlb_scale,
1322 ddcsx, ddcsy, ddcsz,
1325 &ddbox, &npme_major, &npme_minor);
1327 make_dd_communicators(fplog, cr, dd_node_order);
1329 /* Set overallocation to avoid frequent reallocation of arrays */
1330 set_over_alloc_dd(TRUE);
1334 /* PME, if used, is done on all nodes with 1D decomposition */
1336 cr->duty = (DUTY_PP | DUTY_PME);
1339 if (!EI_TPI(inputrec->eI))
1341 npme_major = cr->nnodes;
1344 if (inputrec->ePBC == epbcSCREW)
1347 "pbc=%s is only implemented with domain decomposition",
1348 epbc_names[inputrec->ePBC]);
1354 /* After possible communicator splitting in make_dd_communicators.
1355 * we can set up the intra/inter node communication.
1357 gmx_setup_nodecomm(fplog, cr);
1360 /* Initialize per-physical-node MPI process/thread ID and counters. */
1361 gmx_init_intranode_counters(cr);
1364 md_print_info(cr, fplog, "Using %d MPI %s\n",
1366 #ifdef GMX_THREAD_MPI
1367 cr->nnodes == 1 ? "thread" : "threads"
1369 cr->nnodes == 1 ? "process" : "processes"
1375 gmx_omp_nthreads_init(fplog, cr,
1376 hwinfo->nthreads_hw_avail,
1377 hw_opt->nthreads_omp,
1378 hw_opt->nthreads_omp_pme,
1379 (cr->duty & DUTY_PP) == 0,
1380 inputrec->cutoff_scheme == ecutsVERLET);
1382 gmx_check_hw_runconf_consistency(fplog, hwinfo, cr, hw_opt->nthreads_tmpi, minf.bUseGPU);
1384 /* getting number of PP/PME threads
1385 PME: env variable should be read only on one node to make sure it is
1386 identical everywhere;
1388 /* TODO nthreads_pp is only used for pinning threads.
1389 * This is a temporary solution until we have a hw topology library.
1391 nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
1392 nthreads_pme = gmx_omp_nthreads_get(emntPME);
1394 wcycle = wallcycle_init(fplog, resetstep, cr, nthreads_pp, nthreads_pme);
1398 /* Master synchronizes its value of reset_counters with all nodes
1399 * including PME only nodes */
1400 reset_counters = wcycle_get_reset_counters(wcycle);
1401 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1402 wcycle_set_reset_counters(wcycle, reset_counters);
1406 if (cr->duty & DUTY_PP)
1408 /* For domain decomposition we allocate dynamically
1409 * in dd_partition_system.
1411 if (DOMAINDECOMP(cr))
1413 bcast_state_setup(cr, state);
1419 bcast_state(cr, state, TRUE);
1423 /* Initiate forcerecord */
1425 fr->hwinfo = hwinfo;
1426 init_forcerec(fplog, oenv, fr, fcd, inputrec, mtop, cr, box, FALSE,
1427 opt2fn("-table", nfile, fnm),
1428 opt2fn("-tabletf", nfile, fnm),
1429 opt2fn("-tablep", nfile, fnm),
1430 opt2fn("-tableb", nfile, fnm),
1434 /* version for PCA_NOT_READ_NODE (see md.c) */
1435 /*init_forcerec(fplog,fr,fcd,inputrec,mtop,cr,box,FALSE,
1436 "nofile","nofile","nofile","nofile",FALSE,pforce);
1438 fr->bSepDVDL = ((Flags & MD_SEPPOT) == MD_SEPPOT);
1440 /* Initialize QM-MM */
1443 init_QMMMrec(cr, box, mtop, inputrec, fr);
1446 /* Initialize the mdatoms structure.
1447 * mdatoms is not filled with atom data,
1448 * as this can not be done now with domain decomposition.
1450 mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1452 if (mdatoms->nPerturbed > 0 && inputrec->cutoff_scheme == ecutsVERLET)
1454 gmx_fatal(FARGS, "The Verlet cut-off scheme does not (yet) support free-energy calculations with perturbed atoms, only perturbed interactions. This will be implemented soon. Use the group scheme for now.");
1457 /* Initialize the virtual site communication */
1458 vsite = init_vsite(mtop, cr, FALSE);
1460 calc_shifts(box, fr->shift_vec);
1462 /* With periodic molecules the charge groups should be whole at start up
1463 * and the virtual sites should not be far from their proper positions.
1465 if (!inputrec->bContinuation && MASTER(cr) &&
1466 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1468 /* Make molecules whole at start of run */
1469 if (fr->ePBC != epbcNONE)
1471 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, state->x);
1475 /* Correct initial vsite positions are required
1476 * for the initial distribution in the domain decomposition
1477 * and for the initial shell prediction.
1479 construct_vsites_mtop(fplog, vsite, mtop, state->x);
1483 if (EEL_PME(fr->eeltype))
1485 ewaldcoeff = fr->ewaldcoeff;
1486 pmedata = &fr->pmedata;
1495 /* This is a PME only node */
1497 /* We don't need the state */
1500 ewaldcoeff = calc_ewaldcoeff(inputrec->rcoulomb, inputrec->ewald_rtol);
1504 if (hw_opt->thread_affinity != threadaffOFF)
1506 /* Before setting affinity, check whether the affinity has changed
1507 * - which indicates that probably the OpenMP library has changed it
1508 * since we first checked).
1510 gmx_check_thread_affinity_set(fplog, cr,
1511 hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1513 /* Set the CPU affinity */
1514 gmx_set_thread_affinity(fplog, cr, hw_opt, nthreads_pme, hwinfo, inputrec);
1517 /* Initiate PME if necessary,
1518 * either on all nodes or on dedicated PME nodes only. */
1519 if (EEL_PME(inputrec->coulombtype))
1523 nChargePerturbed = mdatoms->nChargePerturbed;
1525 if (cr->npmenodes > 0)
1527 /* The PME only nodes need to know nChargePerturbed */
1528 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1531 if (cr->duty & DUTY_PME)
1533 status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1534 mtop ? mtop->natoms : 0, nChargePerturbed,
1535 (Flags & MD_REPRODUCIBLE), nthreads_pme);
1538 gmx_fatal(FARGS, "Error %d initializing PME", status);
1544 if (integrator[inputrec->eI].func == do_md)
1546 /* Turn on signal handling on all nodes */
1548 * (A user signal from the PME nodes (if any)
1549 * is communicated to the PP nodes.
1551 signal_handler_install();
1554 if (cr->duty & DUTY_PP)
1556 if (inputrec->ePull != epullNO)
1558 /* Initialize pull code */
1559 init_pull(fplog, inputrec, nfile, fnm, mtop, cr, oenv, inputrec->fepvals->init_lambda,
1560 EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1565 /* Initialize enforced rotation code */
1566 init_rot(fplog, inputrec, nfile, fnm, cr, state->x, box, mtop, oenv,
1570 constr = init_constraints(fplog, mtop, inputrec, ed, state, cr);
1572 if (DOMAINDECOMP(cr))
1574 dd_init_bondeds(fplog, cr->dd, mtop, vsite, constr, inputrec,
1575 Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1577 set_dd_parameters(fplog, cr->dd, dlb_scale, inputrec, fr, &ddbox);
1579 setup_dd_grid(fplog, cr->dd);
1582 /* Now do whatever the user wants us to do (how flexible...) */
1583 integrator[inputrec->eI].func(fplog, cr, nfile, fnm,
1584 oenv, bVerbose, bCompact,
1587 nstepout, inputrec, mtop,
1589 mdatoms, nrnb, wcycle, ed, fr,
1590 repl_ex_nst, repl_ex_nex, repl_ex_seed,
1592 cpt_period, max_hours,
1597 if (inputrec->ePull != epullNO)
1599 finish_pull(fplog, inputrec->pull);
1604 finish_rot(inputrec->rot);
1611 gmx_pmeonly(*pmedata, cr, nrnb, wcycle, ewaldcoeff, FALSE, inputrec);
1614 if (EI_DYNAMICS(inputrec->eI) || EI_TPI(inputrec->eI))
1616 /* Some timing stats */
1619 if (runtime.proc == 0)
1621 runtime.proc = runtime.real;
1630 wallcycle_stop(wcycle, ewcRUN);
1632 /* Finish up, write some stuff
1633 * if rerunMD, don't write last frame again
1635 finish_run(fplog, cr, ftp2fn(efSTO, nfile, fnm),
1636 inputrec, nrnb, wcycle, &runtime,
1637 fr != NULL && fr->nbv != NULL && fr->nbv->bUseGPU ?
1638 nbnxn_cuda_get_timings(fr->nbv->cu_nbv) : NULL,
1640 EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1642 if ((cr->duty & DUTY_PP) && fr->nbv != NULL && fr->nbv->bUseGPU)
1644 char gpu_err_str[STRLEN];
1646 /* free GPU memory and uninitialize GPU (by destroying the context) */
1647 nbnxn_cuda_free(fplog, fr->nbv->cu_nbv);
1649 if (!free_gpu(gpu_err_str))
1651 gmx_warning("On node %d failed to free GPU #%d: %s",
1652 cr->nodeid, get_current_gpu_device_id(), gpu_err_str);
1656 if (opt2bSet("-membed", nfile, fnm))
1661 #ifdef GMX_THREAD_MPI
1662 if (PAR(cr) && SIMMASTER(cr))
1665 gmx_hardware_info_free(hwinfo);
1668 /* Does what it says */
1669 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", &runtime);
1671 /* Close logfile already here if we were appending to it */
1672 if (MASTER(cr) && (Flags & MD_APPENDFILES))
1674 gmx_log_close(fplog);
1677 rc = (int)gmx_get_stop_condition();
1679 #ifdef GMX_THREAD_MPI
1680 /* we need to join all threads. The sub-threads join when they
1681 exit this function, but the master thread needs to be told to
1683 if (PAR(cr) && MASTER(cr))