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
52 #include "md_logging.h"
53 #include "md_support.h"
56 #include "pull_rotation.h"
69 #include "checkpoint.h"
70 #include "mtop_util.h"
71 #include "sighandler.h"
74 #include "gmx_detect_hardware.h"
75 #include "gmx_omp_nthreads.h"
76 #include "pull_rotation.h"
77 #include "calc_verletbuf.h"
78 #include "../mdlib/nbnxn_search.h"
79 #include "../mdlib/nbnxn_consts.h"
80 #include "gmx_fatal_collective.h"
84 #include "gmx_thread_affinity.h"
86 #include "gromacs/utility/gmxmpi.h"
92 #include "gpu_utils.h"
93 #include "nbnxn_cuda_data_mgmt.h"
96 gmx_integrator_t *func;
99 /* The array should match the eI array in include/types/enums.h */
100 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}};
102 gmx_large_int_t deform_init_init_step_tpx;
103 matrix deform_init_box_tpx;
104 #ifdef GMX_THREAD_MPI
105 tMPI_Thread_mutex_t deform_init_box_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
109 #ifdef GMX_THREAD_MPI
110 struct mdrunner_arglist
112 gmx_hw_opt_t *hw_opt;
125 const char *dddlb_opt;
130 const char *nbpu_opt;
131 gmx_large_int_t nsteps_cmdline;
141 const char *deviceOptions;
143 int ret; /* return value */
147 /* The function used for spawning threads. Extracts the mdrunner()
148 arguments from its one argument and calls mdrunner(), after making
150 static void mdrunner_start_fn(void *arg)
152 struct mdrunner_arglist *mda = (struct mdrunner_arglist*)arg;
153 struct mdrunner_arglist mc = *mda; /* copy the arg list to make sure
154 that it's thread-local. This doesn't
155 copy pointed-to items, of course,
156 but those are all const. */
157 t_commrec *cr; /* we need a local version of this */
161 fnm = dup_tfn(mc.nfile, mc.fnm);
163 cr = init_par_threads(mc.cr);
170 mda->ret = mdrunner(mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
171 mc.bVerbose, mc.bCompact, mc.nstglobalcomm,
172 mc.ddxyz, mc.dd_node_order, mc.rdd,
173 mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
174 mc.ddcsx, mc.ddcsy, mc.ddcsz,
176 mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
177 mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce,
178 mc.cpt_period, mc.max_hours, mc.deviceOptions, mc.Flags);
181 /* called by mdrunner() to start a specific number of threads (including
182 the main thread) for thread-parallel runs. This in turn calls mdrunner()
184 All options besides nthreads are the same as for mdrunner(). */
185 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
186 FILE *fplog, t_commrec *cr, int nfile,
187 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
188 gmx_bool bCompact, int nstglobalcomm,
189 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
190 const char *dddlb_opt, real dlb_scale,
191 const char *ddcsx, const char *ddcsy, const char *ddcsz,
192 const char *nbpu_opt,
193 gmx_large_int_t nsteps_cmdline,
194 int nstepout, int resetstep,
195 int nmultisim, int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
196 real pforce, real cpt_period, real max_hours,
197 const char *deviceOptions, unsigned long Flags)
200 struct mdrunner_arglist *mda;
201 t_commrec *crn; /* the new commrec */
204 /* first check whether we even need to start tMPI */
205 if (hw_opt->nthreads_tmpi < 2)
210 /* a few small, one-time, almost unavoidable memory leaks: */
212 fnmn = dup_tfn(nfile, fnm);
214 /* fill the data structure to pass as void pointer to thread start fn */
215 mda->hw_opt = hw_opt;
221 mda->bVerbose = bVerbose;
222 mda->bCompact = bCompact;
223 mda->nstglobalcomm = nstglobalcomm;
224 mda->ddxyz[XX] = ddxyz[XX];
225 mda->ddxyz[YY] = ddxyz[YY];
226 mda->ddxyz[ZZ] = ddxyz[ZZ];
227 mda->dd_node_order = dd_node_order;
229 mda->rconstr = rconstr;
230 mda->dddlb_opt = dddlb_opt;
231 mda->dlb_scale = dlb_scale;
235 mda->nbpu_opt = nbpu_opt;
236 mda->nsteps_cmdline = nsteps_cmdline;
237 mda->nstepout = nstepout;
238 mda->resetstep = resetstep;
239 mda->nmultisim = nmultisim;
240 mda->repl_ex_nst = repl_ex_nst;
241 mda->repl_ex_nex = repl_ex_nex;
242 mda->repl_ex_seed = repl_ex_seed;
243 mda->pforce = pforce;
244 mda->cpt_period = cpt_period;
245 mda->max_hours = max_hours;
246 mda->deviceOptions = deviceOptions;
249 /* now spawn new threads that start mdrunner_start_fn(), while
250 the main thread returns, we set thread affinity later */
251 ret = tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi, TMPI_AFFINITY_NONE,
252 mdrunner_start_fn, (void*)(mda) );
253 if (ret != TMPI_SUCCESS)
258 /* make a new comm_rec to reflect the new situation */
259 crn = init_par_threads(cr);
264 static int get_tmpi_omp_thread_division(const gmx_hw_info_t *hwinfo,
265 const gmx_hw_opt_t *hw_opt,
271 /* There are no separate PME nodes here, as we ensured in
272 * check_and_update_hw_opt that nthreads_tmpi>0 with PME nodes
273 * and a conditional ensures we would not have ended up here.
274 * Note that separate PME nodes might be switched on later.
278 nthreads_tmpi = ngpu;
279 if (nthreads_tot > 0 && nthreads_tot < nthreads_tmpi)
281 nthreads_tmpi = nthreads_tot;
284 else if (hw_opt->nthreads_omp > 0)
286 /* Here we could oversubscribe, when we do, we issue a warning later */
287 nthreads_tmpi = max(1, nthreads_tot/hw_opt->nthreads_omp);
291 /* TODO choose nthreads_omp based on hardware topology
292 when we have a hardware topology detection library */
293 /* In general, when running up to 4 threads, OpenMP should be faster.
294 * Note: on AMD Bulldozer we should avoid running OpenMP over two dies.
295 * On Intel>=Nehalem running OpenMP on a single CPU is always faster,
296 * even on two CPUs it's usually faster (but with many OpenMP threads
297 * it could be faster not to use HT, currently we always use HT).
298 * On Nehalem/Westmere we want to avoid running 16 threads over
299 * two CPUs with HT, so we need a limit<16; thus we use 12.
300 * A reasonable limit for Intel Sandy and Ivy bridge,
301 * not knowing the topology, is 16 threads.
303 const int nthreads_omp_always_faster = 4;
304 const int nthreads_omp_always_faster_Nehalem = 12;
305 const int nthreads_omp_always_faster_SandyBridge = 16;
306 const int first_model_Nehalem = 0x1A;
307 const int first_model_SandyBridge = 0x2A;
308 gmx_bool bIntel_Family6;
311 (gmx_cpuid_vendor(hwinfo->cpuid_info) == GMX_CPUID_VENDOR_INTEL &&
312 gmx_cpuid_family(hwinfo->cpuid_info) == 6);
314 if (nthreads_tot <= nthreads_omp_always_faster ||
316 ((gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_Nehalem && nthreads_tot <= nthreads_omp_always_faster_Nehalem) ||
317 (gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_SandyBridge && nthreads_tot <= nthreads_omp_always_faster_SandyBridge))))
319 /* Use pure OpenMP parallelization */
324 /* Don't use OpenMP parallelization */
325 nthreads_tmpi = nthreads_tot;
329 return nthreads_tmpi;
333 /* Get the number of threads to use for thread-MPI based on how many
334 * were requested, which algorithms we're using,
335 * and how many particles there are.
336 * At the point we have already called check_and_update_hw_opt.
337 * Thus all options should be internally consistent and consistent
338 * with the hardware, except that ntmpi could be larger than #GPU.
340 static int get_nthreads_mpi(gmx_hw_info_t *hwinfo,
341 gmx_hw_opt_t *hw_opt,
342 t_inputrec *inputrec, gmx_mtop_t *mtop,
346 int nthreads_hw, nthreads_tot_max, nthreads_tmpi, nthreads_new, ngpu;
347 int min_atoms_per_mpi_thread;
352 if (hw_opt->nthreads_tmpi > 0)
354 /* Trivial, return right away */
355 return hw_opt->nthreads_tmpi;
358 nthreads_hw = hwinfo->nthreads_hw_avail;
360 /* How many total (#tMPI*#OpenMP) threads can we start? */
361 if (hw_opt->nthreads_tot > 0)
363 nthreads_tot_max = hw_opt->nthreads_tot;
367 nthreads_tot_max = nthreads_hw;
370 bCanUseGPU = (inputrec->cutoff_scheme == ecutsVERLET && hwinfo->bCanUseGPU);
373 ngpu = hwinfo->gpu_info.ncuda_dev_use;
381 get_tmpi_omp_thread_division(hwinfo, hw_opt, nthreads_tot_max, ngpu);
383 if (inputrec->eI == eiNM || EI_TPI(inputrec->eI))
385 /* Steps are divided over the nodes iso splitting the atoms */
386 min_atoms_per_mpi_thread = 0;
392 min_atoms_per_mpi_thread = MIN_ATOMS_PER_GPU;
396 min_atoms_per_mpi_thread = MIN_ATOMS_PER_MPI_THREAD;
400 /* Check if an algorithm does not support parallel simulation. */
401 if (nthreads_tmpi != 1 &&
402 ( inputrec->eI == eiLBFGS ||
403 inputrec->coulombtype == eelEWALD ) )
407 md_print_warn(cr, fplog, "The integration or electrostatics algorithm doesn't support parallel runs. Using a single thread-MPI thread.\n");
408 if (hw_opt->nthreads_tmpi > nthreads_tmpi)
410 gmx_fatal(FARGS, "You asked for more than 1 thread-MPI thread, but an algorithm doesn't support that");
413 else if (mtop->natoms/nthreads_tmpi < min_atoms_per_mpi_thread)
415 /* the thread number was chosen automatically, but there are too many
416 threads (too few atoms per thread) */
417 nthreads_new = max(1, mtop->natoms/min_atoms_per_mpi_thread);
419 /* Avoid partial use of Hyper-Threading */
420 if (gmx_cpuid_x86_smt(hwinfo->cpuid_info) == GMX_CPUID_X86_SMT_ENABLED &&
421 nthreads_new > nthreads_hw/2 && nthreads_new < nthreads_hw)
423 nthreads_new = nthreads_hw/2;
426 /* Avoid large prime numbers in the thread count */
427 if (nthreads_new >= 6)
429 /* Use only 6,8,10 with additional factors of 2 */
433 while (3*fac*2 <= nthreads_new)
438 nthreads_new = (nthreads_new/fac)*fac;
443 if (nthreads_new == 5)
449 nthreads_tmpi = nthreads_new;
451 fprintf(stderr, "\n");
452 fprintf(stderr, "NOTE: Parallelization is limited by the small number of atoms,\n");
453 fprintf(stderr, " only starting %d thread-MPI threads.\n", nthreads_tmpi);
454 fprintf(stderr, " You can use the -nt and/or -ntmpi option to optimize the number of threads.\n\n");
457 return nthreads_tmpi;
459 #endif /* GMX_THREAD_MPI */
462 /* Environment variable for setting nstlist */
463 static const char* NSTLIST_ENVVAR = "GMX_NSTLIST";
464 /* Try to increase nstlist when using a GPU with nstlist less than this */
465 static const int NSTLIST_GPU_ENOUGH = 20;
466 /* Increase nstlist until the non-bonded cost increases more than this factor */
467 static const float NBNXN_GPU_LIST_OK_FAC = 1.25;
468 /* Don't increase nstlist beyond a non-bonded cost increases of this factor */
469 static const float NBNXN_GPU_LIST_MAX_FAC = 1.40;
471 /* Try to increase nstlist when running on a GPU */
472 static void increase_nstlist(FILE *fp, t_commrec *cr,
473 t_inputrec *ir, const gmx_mtop_t *mtop, matrix box)
476 int nstlist_orig, nstlist_prev;
477 verletbuf_list_setup_t ls;
478 real rlist_inc, rlist_ok, rlist_max, rlist_new, rlist_prev;
481 gmx_bool bBox, bDD, bCont;
482 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";
483 const char *vbd_err = "Can not increase nstlist for GPU run because verlet-buffer-drift is not set or used";
484 const char *box_err = "Can not increase nstlist for GPU run because the box is too small";
485 const char *dd_err = "Can not increase nstlist for GPU run because of domain decomposition limitations";
488 /* Number of + nstlist alternative values to try when switching */
489 const int nstl[] = { 20, 25, 40, 50 };
490 #define NNSTL sizeof(nstl)/sizeof(nstl[0])
492 env = getenv(NSTLIST_ENVVAR);
497 fprintf(fp, nstl_fmt, ir->nstlist);
501 if (ir->verletbuf_drift == 0)
503 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");
506 if (ir->verletbuf_drift < 0)
510 fprintf(stderr, "%s\n", vbd_err);
514 fprintf(fp, "%s\n", vbd_err);
520 nstlist_orig = ir->nstlist;
523 sprintf(buf, "Getting nstlist from environment variable GMX_NSTLIST=%s", env);
526 fprintf(stderr, "%s\n", buf);
530 fprintf(fp, "%s\n", buf);
532 sscanf(env, "%d", &ir->nstlist);
535 verletbuf_get_list_setup(TRUE, &ls);
537 /* Allow rlist to make the list double the size of the cut-off sphere */
538 rlist_inc = nbnxn_get_rlist_effective_inc(NBNXN_GPU_CLUSTER_SIZE, mtop->natoms/det(box));
539 rlist_ok = (max(ir->rvdw, ir->rcoulomb) + rlist_inc)*pow(NBNXN_GPU_LIST_OK_FAC, 1.0/3.0) - rlist_inc;
540 rlist_max = (max(ir->rvdw, ir->rcoulomb) + rlist_inc)*pow(NBNXN_GPU_LIST_MAX_FAC, 1.0/3.0) - rlist_inc;
543 fprintf(debug, "GPU nstlist tuning: rlist_inc %.3f rlist_max %.3f\n",
544 rlist_inc, rlist_max);
548 nstlist_prev = nstlist_orig;
549 rlist_prev = ir->rlist;
554 ir->nstlist = nstl[i];
557 /* Set the pair-list buffer size in ir */
558 calc_verlet_buffer_size(mtop, det(box), ir, ir->verletbuf_drift, &ls,
561 /* Does rlist fit in the box? */
562 bBox = (sqr(rlist_new) < max_cutoff2(ir->ePBC, box));
564 if (bBox && DOMAINDECOMP(cr))
566 /* Check if rlist fits in the domain decomposition */
567 if (inputrec2nboundeddim(ir) < DIM)
569 gmx_incons("Changing nstlist with domain decomposition and unbounded dimensions is not implemented yet");
571 copy_mat(box, state_tmp.box);
572 bDD = change_dd_cutoff(cr, &state_tmp, ir, rlist_new);
579 if (bBox && bDD && rlist_new <= rlist_max)
581 /* Increase nstlist */
582 nstlist_prev = ir->nstlist;
583 rlist_prev = rlist_new;
584 bCont = (i+1 < NNSTL && rlist_new < rlist_ok);
588 /* Stick with the previous nstlist */
589 ir->nstlist = nstlist_prev;
590 rlist_new = rlist_prev;
602 gmx_warning(!bBox ? box_err : dd_err);
605 fprintf(fp, "\n%s\n", bBox ? box_err : dd_err);
607 ir->nstlist = nstlist_orig;
609 else if (ir->nstlist != nstlist_orig || rlist_new != ir->rlist)
611 sprintf(buf, "Changing nstlist from %d to %d, rlist from %g to %g",
612 nstlist_orig, ir->nstlist,
613 ir->rlist, rlist_new);
616 fprintf(stderr, "%s\n\n", buf);
620 fprintf(fp, "%s\n\n", buf);
622 ir->rlist = rlist_new;
623 ir->rlistlong = rlist_new;
627 static void prepare_verlet_scheme(FILE *fplog,
628 gmx_hw_info_t *hwinfo,
630 const char *nbpu_opt,
632 const gmx_mtop_t *mtop,
636 /* Here we only check for GPU usage on the MPI master process,
637 * as here we don't know how many GPUs we will use yet.
638 * We check for a GPU on all processes later.
640 *bUseGPU = hwinfo->bCanUseGPU || (getenv("GMX_EMULATE_GPU") != NULL);
642 if (ir->verletbuf_drift > 0)
644 /* Update the Verlet buffer size for the current run setup */
645 verletbuf_list_setup_t ls;
648 /* Here we assume CPU acceleration is on. But as currently
649 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
650 * and 4x2 gives a larger buffer than 4x4, this is ok.
652 verletbuf_get_list_setup(*bUseGPU, &ls);
654 calc_verlet_buffer_size(mtop, det(box), ir,
655 ir->verletbuf_drift, &ls,
657 if (rlist_new != ir->rlist)
661 fprintf(fplog, "\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
662 ir->rlist, rlist_new,
663 ls.cluster_size_i, ls.cluster_size_j);
665 ir->rlist = rlist_new;
666 ir->rlistlong = rlist_new;
670 /* With GPU or emulation we should check nstlist for performance */
671 if ((EI_DYNAMICS(ir->eI) &&
673 ir->nstlist < NSTLIST_GPU_ENOUGH) ||
674 getenv(NSTLIST_ENVVAR) != NULL)
676 /* Choose a better nstlist */
677 increase_nstlist(fplog, cr, ir, mtop, box);
681 static void convert_to_verlet_scheme(FILE *fplog,
683 gmx_mtop_t *mtop, real box_vol)
685 char *conv_mesg = "Converting input file with group cut-off scheme to the Verlet cut-off scheme";
687 md_print_warn(NULL, fplog, "%s\n", conv_mesg);
689 ir->cutoff_scheme = ecutsVERLET;
690 ir->verletbuf_drift = 0.005;
692 if (ir->rcoulomb != ir->rvdw)
694 gmx_fatal(FARGS, "The VdW and Coulomb cut-offs are different, whereas the Verlet scheme only supports equal cut-offs");
697 if (ir->vdwtype == evdwUSER || EEL_USER(ir->coulombtype))
699 gmx_fatal(FARGS, "User non-bonded potentials are not (yet) supported with the Verlet scheme");
701 else if (EVDW_SWITCHED(ir->vdwtype) || EEL_SWITCHED(ir->coulombtype))
703 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");
705 if (EVDW_SWITCHED(ir->vdwtype))
707 ir->vdwtype = evdwCUT;
709 if (EEL_SWITCHED(ir->coulombtype))
711 if (EEL_FULL(ir->coulombtype))
713 /* With full electrostatic only PME can be switched */
714 ir->coulombtype = eelPME;
718 md_print_warn(NULL, fplog, "NOTE: Replacing %s electrostatics with reaction-field with epsilon-rf=inf\n", eel_names[ir->coulombtype]);
719 ir->coulombtype = eelRF;
720 ir->epsilon_rf = 0.0;
724 /* We set the target energy drift to a small number.
725 * Note that this is only for testing. For production the user
726 * should think about this and set the mdp options.
728 ir->verletbuf_drift = 1e-4;
731 if (inputrec2nboundeddim(ir) != 3)
733 gmx_fatal(FARGS, "Can only convert old tpr files to the Verlet cut-off scheme with 3D pbc");
736 if (ir->efep != efepNO || ir->implicit_solvent != eisNO)
738 gmx_fatal(FARGS, "Will not convert old tpr files to the Verlet cut-off scheme with free-energy calculations or implicit solvent");
741 if (EI_DYNAMICS(ir->eI) && !(EI_MD(ir->eI) && ir->etc == etcNO))
743 verletbuf_list_setup_t ls;
745 verletbuf_get_list_setup(FALSE, &ls);
746 calc_verlet_buffer_size(mtop, box_vol, ir, ir->verletbuf_drift, &ls,
751 ir->verletbuf_drift = -1;
752 ir->rlist = 1.05*max(ir->rvdw, ir->rcoulomb);
755 gmx_mtop_remove_chargegroups(mtop);
758 static void check_and_update_hw_opt(gmx_hw_opt_t *hw_opt,
760 gmx_bool bIsSimMaster)
762 gmx_omp_nthreads_read_env(&hw_opt->nthreads_omp, bIsSimMaster);
764 #ifndef GMX_THREAD_MPI
765 if (hw_opt->nthreads_tot > 0)
767 gmx_fatal(FARGS, "Setting the total number of threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
769 if (hw_opt->nthreads_tmpi > 0)
771 gmx_fatal(FARGS, "Setting the number of thread-MPI threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
775 if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp_pme <= 0)
777 /* We have the same number of OpenMP threads for PP and PME processes,
778 * thus we can perform several consistency checks.
780 if (hw_opt->nthreads_tmpi > 0 &&
781 hw_opt->nthreads_omp > 0 &&
782 hw_opt->nthreads_tot != hw_opt->nthreads_tmpi*hw_opt->nthreads_omp)
784 gmx_fatal(FARGS, "The total number of threads requested (%d) does not match the thread-MPI threads (%d) times the OpenMP threads (%d) requested",
785 hw_opt->nthreads_tot, hw_opt->nthreads_tmpi, hw_opt->nthreads_omp);
788 if (hw_opt->nthreads_tmpi > 0 &&
789 hw_opt->nthreads_tot % hw_opt->nthreads_tmpi != 0)
791 gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of thread-MPI threads requested (%d)",
792 hw_opt->nthreads_tot, hw_opt->nthreads_tmpi);
795 if (hw_opt->nthreads_omp > 0 &&
796 hw_opt->nthreads_tot % hw_opt->nthreads_omp != 0)
798 gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of OpenMP threads requested (%d)",
799 hw_opt->nthreads_tot, hw_opt->nthreads_omp);
802 if (hw_opt->nthreads_tmpi > 0 &&
803 hw_opt->nthreads_omp <= 0)
805 hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
810 if (hw_opt->nthreads_omp > 1)
812 gmx_fatal(FARGS, "OpenMP threads are requested, but Gromacs was compiled without OpenMP support");
816 if (cutoff_scheme == ecutsGROUP)
818 /* We only have OpenMP support for PME only nodes */
819 if (hw_opt->nthreads_omp > 1)
821 gmx_fatal(FARGS, "OpenMP threads have been requested with cut-off scheme %s, but these are only supported with cut-off scheme %s",
822 ecutscheme_names[cutoff_scheme],
823 ecutscheme_names[ecutsVERLET]);
825 hw_opt->nthreads_omp = 1;
828 if (hw_opt->nthreads_omp_pme > 0 && hw_opt->nthreads_omp <= 0)
830 gmx_fatal(FARGS, "You need to specify -ntomp in addition to -ntomp_pme");
833 if (hw_opt->nthreads_tot == 1)
835 hw_opt->nthreads_tmpi = 1;
837 if (hw_opt->nthreads_omp > 1)
839 gmx_fatal(FARGS, "You requested %d OpenMP threads with %d total threads",
840 hw_opt->nthreads_tmpi, hw_opt->nthreads_tot);
842 hw_opt->nthreads_omp = 1;
845 if (hw_opt->nthreads_omp_pme <= 0 && hw_opt->nthreads_omp > 0)
847 hw_opt->nthreads_omp_pme = hw_opt->nthreads_omp;
852 fprintf(debug, "hw_opt: nt %d ntmpi %d ntomp %d ntomp_pme %d gpu_id '%s'\n",
853 hw_opt->nthreads_tot,
854 hw_opt->nthreads_tmpi,
855 hw_opt->nthreads_omp,
856 hw_opt->nthreads_omp_pme,
857 hw_opt->gpu_id != NULL ? hw_opt->gpu_id : "");
863 /* Override the value in inputrec with value passed on the command line (if any) */
864 static void override_nsteps_cmdline(FILE *fplog,
865 gmx_large_int_t nsteps_cmdline,
869 char sbuf[STEPSTRSIZE];
874 /* override with anything else than the default -2 */
875 if (nsteps_cmdline > -2)
879 ir->nsteps = nsteps_cmdline;
880 if (EI_DYNAMICS(ir->eI))
882 sprintf(stmp, "Overriding nsteps with value passed on the command line: %s steps, %.3f ps",
883 gmx_step_str(nsteps_cmdline, sbuf),
884 nsteps_cmdline*ir->delta_t);
888 sprintf(stmp, "Overriding nsteps with value passed on the command line: %s steps",
889 gmx_step_str(nsteps_cmdline, sbuf));
892 md_print_warn(cr, fplog, "%s\n", stmp);
896 /* Data structure set by SIMMASTER which needs to be passed to all nodes
897 * before the other nodes have read the tpx file and called gmx_detect_hardware.
900 int cutoff_scheme; /* The cutoff scheme from inputrec_t */
901 gmx_bool bUseGPU; /* Use GPU or GPU emulation */
904 int mdrunner(gmx_hw_opt_t *hw_opt,
905 FILE *fplog, t_commrec *cr, int nfile,
906 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
907 gmx_bool bCompact, int nstglobalcomm,
908 ivec ddxyz, int dd_node_order, real rdd, real rconstr,
909 const char *dddlb_opt, real dlb_scale,
910 const char *ddcsx, const char *ddcsy, const char *ddcsz,
911 const char *nbpu_opt,
912 gmx_large_int_t nsteps_cmdline, int nstepout, int resetstep,
913 int nmultisim, int repl_ex_nst, int repl_ex_nex,
914 int repl_ex_seed, real pforce, real cpt_period, real max_hours,
915 const char *deviceOptions, unsigned long Flags)
917 gmx_bool bForceUseGPU, bTryUseGPU;
918 double nodetime = 0, realtime;
919 t_inputrec *inputrec;
920 t_state *state = NULL;
922 gmx_ddbox_t ddbox = {0};
923 int npme_major, npme_minor;
926 gmx_mtop_t *mtop = NULL;
927 t_mdatoms *mdatoms = NULL;
928 t_forcerec *fr = NULL;
929 t_fcdata *fcd = NULL;
931 gmx_pme_t *pmedata = NULL;
932 gmx_vsite_t *vsite = NULL;
934 int i, m, nChargePerturbed = -1, status, nalloc;
936 gmx_wallcycle_t wcycle;
937 gmx_bool bReadRNG, bReadEkin;
939 gmx_runtime_t runtime;
941 gmx_large_int_t reset_counters;
942 gmx_edsam_t ed = NULL;
943 t_commrec *cr_old = cr;
944 int nthreads_pme = 1;
946 gmx_membed_t membed = NULL;
947 gmx_hw_info_t *hwinfo = NULL;
948 master_inf_t minf = {-1, FALSE};
950 /* CAUTION: threads may be started later on in this function, so
951 cr doesn't reflect the final parallel state right now */
955 if (Flags & MD_APPENDFILES)
960 bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
961 bTryUseGPU = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
966 /* Read (nearly) all data required for the simulation */
967 read_tpx_state(ftp2fn(efTPX, nfile, fnm), inputrec, state, NULL, mtop);
969 if (inputrec->cutoff_scheme != ecutsVERLET &&
970 ((Flags & MD_TESTVERLET) || getenv("GMX_VERLET_SCHEME") != NULL))
972 convert_to_verlet_scheme(fplog, inputrec, mtop, det(state->box));
975 /* Detect hardware, gather information. With tMPI only thread 0 does it
976 * and after threads are started broadcasts hwinfo around. */
978 gmx_detect_hardware(fplog, hwinfo, cr,
979 bForceUseGPU, bTryUseGPU, hw_opt->gpu_id);
981 minf.cutoff_scheme = inputrec->cutoff_scheme;
982 minf.bUseGPU = FALSE;
984 if (inputrec->cutoff_scheme == ecutsVERLET)
986 prepare_verlet_scheme(fplog, hwinfo, cr, nbpu_opt,
987 inputrec, mtop, state->box,
990 else if (hwinfo->bCanUseGPU)
992 md_print_warn(cr, fplog,
993 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
994 " To use a GPU, set the mdp option: cutoff-scheme = Verlet\n"
995 " (for quick performance testing you can use the -testverlet option)\n");
999 gmx_fatal(FARGS, "GPU requested, but can't be used without cutoff-scheme=Verlet");
1003 #ifndef GMX_THREAD_MPI
1006 gmx_bcast_sim(sizeof(minf), &minf, cr);
1009 if (minf.bUseGPU && cr->npmenodes == -1)
1011 /* Don't automatically use PME-only nodes with GPUs */
1015 /* Check for externally set OpenMP affinity and turn off internal
1016 * pinning if any is found. We need to do this check early to tell
1017 * thread-MPI whether it should do pinning when spawning threads.
1018 * TODO: the above no longer holds, we should move these checks down
1020 gmx_omp_check_thread_affinity(fplog, cr, hw_opt);
1022 #ifdef GMX_THREAD_MPI
1023 /* With thread-MPI inputrec is only set here on the master thread */
1027 check_and_update_hw_opt(hw_opt, minf.cutoff_scheme, SIMMASTER(cr));
1029 #ifdef GMX_THREAD_MPI
1030 /* Early check for externally set process affinity. Can't do over all
1031 * MPI processes because hwinfo is not available everywhere, but with
1032 * thread-MPI it's needed as pinning might get turned off which needs
1033 * to be known before starting thread-MPI. */
1034 gmx_check_thread_affinity_set(fplog,
1036 hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1039 #ifdef GMX_THREAD_MPI
1040 if (cr->npmenodes > 0 && hw_opt->nthreads_tmpi <= 0)
1042 gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME nodes");
1046 if (hw_opt->nthreads_omp_pme != hw_opt->nthreads_omp &&
1049 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");
1053 #ifdef GMX_THREAD_MPI
1056 /* NOW the threads will be started: */
1057 hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
1061 if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp <= 0)
1063 hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
1066 if (hw_opt->nthreads_tmpi > 1)
1068 /* now start the threads. */
1069 cr = mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
1070 oenv, bVerbose, bCompact, nstglobalcomm,
1071 ddxyz, dd_node_order, rdd, rconstr,
1072 dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
1074 nsteps_cmdline, nstepout, resetstep, nmultisim,
1075 repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
1076 cpt_period, max_hours, deviceOptions,
1078 /* the main thread continues here with a new cr. We don't deallocate
1079 the old cr because other threads may still be reading it. */
1082 gmx_comm("Failed to spawn threads");
1087 /* END OF CAUTION: cr is now reliable */
1089 /* g_membed initialisation *
1090 * Because we change the mtop, init_membed is called before the init_parallel *
1091 * (in case we ever want to make it run in parallel) */
1092 if (opt2bSet("-membed", nfile, fnm))
1096 fprintf(stderr, "Initializing membed");
1098 membed = init_membed(fplog, nfile, fnm, mtop, inputrec, state, cr, &cpt_period);
1103 /* now broadcast everything to the non-master nodes/threads: */
1104 init_parallel(fplog, cr, inputrec, mtop);
1106 /* This check needs to happen after get_nthreads_mpi() */
1107 if (inputrec->cutoff_scheme == ecutsVERLET && (Flags & MD_PARTDEC))
1109 gmx_fatal_collective(FARGS, cr, NULL,
1110 "The Verlet cut-off scheme is not supported with particle decomposition.\n"
1111 "You can achieve the same effect as particle decomposition by running in parallel using only OpenMP threads.");
1116 pr_inputrec(fplog, 0, "Input Parameters", inputrec, FALSE);
1119 #if defined GMX_THREAD_MPI
1120 /* With tMPI we detected on thread 0 and we'll just pass the hwinfo pointer
1121 * to the other threads -- slightly uncool, but works fine, just need to
1122 * make sure that the data doesn't get freed twice. */
1129 gmx_bcast(sizeof(&hwinfo), &hwinfo, cr);
1132 if (PAR(cr) && !SIMMASTER(cr))
1134 /* now we have inputrec on all nodes, can run the detection */
1135 /* TODO: perhaps it's better to propagate within a node instead? */
1137 gmx_detect_hardware(fplog, hwinfo, cr,
1138 bForceUseGPU, bTryUseGPU, hw_opt->gpu_id);
1141 /* Now do the affinity check with MPI/no-MPI (done earlier with thread-MPI). */
1142 gmx_check_thread_affinity_set(fplog, cr,
1143 hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1146 /* now make sure the state is initialized and propagated */
1147 set_state_entries(state, inputrec, cr->nnodes);
1149 /* A parallel command line option consistency check that we can
1150 only do after any threads have started. */
1152 (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
1155 "The -dd or -npme option request a parallel simulation, "
1157 "but %s was compiled without threads or MPI enabled"
1159 #ifdef GMX_THREAD_MPI
1160 "but the number of threads (option -nt) is 1"
1162 "but %s was not started through mpirun/mpiexec or only one process was requested through mpirun/mpiexec"
1169 if ((Flags & MD_RERUN) &&
1170 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1172 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
1175 if (can_use_allvsall(inputrec, mtop, TRUE, cr, fplog) && PAR(cr))
1177 /* Simple neighbour searching and (also?) all-vs-all loops
1178 * do not work with domain decomposition. */
1179 Flags |= MD_PARTDEC;
1182 if (!EEL_PME(inputrec->coulombtype) || (Flags & MD_PARTDEC))
1184 if (cr->npmenodes > 0)
1186 if (!EEL_PME(inputrec->coulombtype))
1188 gmx_fatal_collective(FARGS, cr, NULL,
1189 "PME nodes are requested, but the system does not use PME electrostatics");
1191 if (Flags & MD_PARTDEC)
1193 gmx_fatal_collective(FARGS, cr, NULL,
1194 "PME nodes are requested, but particle decomposition does not support separate PME nodes");
1202 fcRegisterSteps(inputrec->nsteps, inputrec->init_step);
1205 /* NMR restraints must be initialized before load_checkpoint,
1206 * since with time averaging the history is added to t_state.
1207 * For proper consistency check we therefore need to extend
1209 * So the PME-only nodes (if present) will also initialize
1210 * the distance restraints.
1214 /* This needs to be called before read_checkpoint to extend the state */
1215 init_disres(fplog, mtop, inputrec, cr, Flags & MD_PARTDEC, fcd, state, repl_ex_nst > 0);
1217 if (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0)
1219 if (PAR(cr) && !(Flags & MD_PARTDEC))
1221 gmx_fatal(FARGS, "Orientation restraints do not work (yet) with domain decomposition, use particle decomposition (mdrun option -pd)");
1223 /* Orientation restraints */
1226 init_orires(fplog, mtop, state->x, inputrec, cr->ms, &(fcd->orires),
1231 if (DEFORM(*inputrec))
1233 /* Store the deform reference box before reading the checkpoint */
1236 copy_mat(state->box, box);
1240 gmx_bcast(sizeof(box), box, cr);
1242 /* Because we do not have the update struct available yet
1243 * in which the reference values should be stored,
1244 * we store them temporarily in static variables.
1245 * This should be thread safe, since they are only written once
1246 * and with identical values.
1248 #ifdef GMX_THREAD_MPI
1249 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
1251 deform_init_init_step_tpx = inputrec->init_step;
1252 copy_mat(box, deform_init_box_tpx);
1253 #ifdef GMX_THREAD_MPI
1254 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
1258 if (opt2bSet("-cpi", nfile, fnm))
1260 /* Check if checkpoint file exists before doing continuation.
1261 * This way we can use identical input options for the first and subsequent runs...
1263 if (gmx_fexist_master(opt2fn_master("-cpi", nfile, fnm, cr), cr) )
1265 load_checkpoint(opt2fn_master("-cpi", nfile, fnm, cr), &fplog,
1266 cr, Flags & MD_PARTDEC, ddxyz,
1267 inputrec, state, &bReadRNG, &bReadEkin,
1268 (Flags & MD_APPENDFILES),
1269 (Flags & MD_APPENDFILESSET));
1273 Flags |= MD_READ_RNG;
1277 Flags |= MD_READ_EKIN;
1282 if (((MASTER(cr) || (Flags & MD_SEPPOT)) && (Flags & MD_APPENDFILES))
1283 #ifdef GMX_THREAD_MPI
1284 /* With thread MPI only the master node/thread exists in mdrun.c,
1285 * therefore non-master nodes need to open the "seppot" log file here.
1287 || (!MASTER(cr) && (Flags & MD_SEPPOT))
1291 gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr, !(Flags & MD_SEPPOT),
1295 /* override nsteps with value from cmdline */
1296 override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
1300 copy_mat(state->box, box);
1305 gmx_bcast(sizeof(box), box, cr);
1308 /* Essential dynamics */
1309 if (opt2bSet("-ei", nfile, fnm))
1311 /* Open input and output files, allocate space for ED data structure */
1312 ed = ed_open(mtop->natoms, &state->edsamstate, nfile, fnm, Flags, oenv, cr);
1315 if (PAR(cr) && !((Flags & MD_PARTDEC) ||
1316 EI_TPI(inputrec->eI) ||
1317 inputrec->eI == eiNM))
1319 cr->dd = init_domain_decomposition(fplog, cr, Flags, ddxyz, rdd, rconstr,
1320 dddlb_opt, dlb_scale,
1321 ddcsx, ddcsy, ddcsz,
1324 &ddbox, &npme_major, &npme_minor);
1326 make_dd_communicators(fplog, cr, dd_node_order);
1328 /* Set overallocation to avoid frequent reallocation of arrays */
1329 set_over_alloc_dd(TRUE);
1333 /* PME, if used, is done on all nodes with 1D decomposition */
1335 cr->duty = (DUTY_PP | DUTY_PME);
1338 if (!EI_TPI(inputrec->eI))
1340 npme_major = cr->nnodes;
1343 if (inputrec->ePBC == epbcSCREW)
1346 "pbc=%s is only implemented with domain decomposition",
1347 epbc_names[inputrec->ePBC]);
1353 /* After possible communicator splitting in make_dd_communicators.
1354 * we can set up the intra/inter node communication.
1356 gmx_setup_nodecomm(fplog, cr);
1359 /* Initialize per-physical-node MPI process/thread ID and counters. */
1360 gmx_init_intranode_counters(cr);
1363 md_print_info(cr, fplog, "Using %d MPI %s\n",
1365 #ifdef GMX_THREAD_MPI
1366 cr->nnodes == 1 ? "thread" : "threads"
1368 cr->nnodes == 1 ? "process" : "processes"
1374 gmx_omp_nthreads_init(fplog, cr,
1375 hwinfo->nthreads_hw_avail,
1376 hw_opt->nthreads_omp,
1377 hw_opt->nthreads_omp_pme,
1378 (cr->duty & DUTY_PP) == 0,
1379 inputrec->cutoff_scheme == ecutsVERLET);
1381 gmx_check_hw_runconf_consistency(fplog, hwinfo, cr, hw_opt->nthreads_tmpi, minf.bUseGPU);
1383 /* getting number of PP/PME threads
1384 PME: env variable should be read only on one node to make sure it is
1385 identical everywhere;
1387 /* TODO nthreads_pp is only used for pinning threads.
1388 * This is a temporary solution until we have a hw topology library.
1390 nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
1391 nthreads_pme = gmx_omp_nthreads_get(emntPME);
1393 wcycle = wallcycle_init(fplog, resetstep, cr, nthreads_pp, nthreads_pme);
1397 /* Master synchronizes its value of reset_counters with all nodes
1398 * including PME only nodes */
1399 reset_counters = wcycle_get_reset_counters(wcycle);
1400 gmx_bcast_sim(sizeof(reset_counters), &reset_counters, cr);
1401 wcycle_set_reset_counters(wcycle, reset_counters);
1405 if (cr->duty & DUTY_PP)
1407 /* For domain decomposition we allocate dynamically
1408 * in dd_partition_system.
1410 if (DOMAINDECOMP(cr))
1412 bcast_state_setup(cr, state);
1418 bcast_state(cr, state, TRUE);
1422 /* Initiate forcerecord */
1424 fr->hwinfo = hwinfo;
1425 init_forcerec(fplog, oenv, fr, fcd, inputrec, mtop, cr, box, FALSE,
1426 opt2fn("-table", nfile, fnm),
1427 opt2fn("-tabletf", nfile, fnm),
1428 opt2fn("-tablep", nfile, fnm),
1429 opt2fn("-tableb", nfile, fnm),
1433 /* version for PCA_NOT_READ_NODE (see md.c) */
1434 /*init_forcerec(fplog,fr,fcd,inputrec,mtop,cr,box,FALSE,
1435 "nofile","nofile","nofile","nofile",FALSE,pforce);
1437 fr->bSepDVDL = ((Flags & MD_SEPPOT) == MD_SEPPOT);
1439 /* Initialize QM-MM */
1442 init_QMMMrec(cr, box, mtop, inputrec, fr);
1445 /* Initialize the mdatoms structure.
1446 * mdatoms is not filled with atom data,
1447 * as this can not be done now with domain decomposition.
1449 mdatoms = init_mdatoms(fplog, mtop, inputrec->efep != efepNO);
1451 if (mdatoms->nPerturbed > 0 && inputrec->cutoff_scheme == ecutsVERLET)
1453 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.");
1456 /* Initialize the virtual site communication */
1457 vsite = init_vsite(mtop, cr, FALSE);
1459 calc_shifts(box, fr->shift_vec);
1461 /* With periodic molecules the charge groups should be whole at start up
1462 * and the virtual sites should not be far from their proper positions.
1464 if (!inputrec->bContinuation && MASTER(cr) &&
1465 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1467 /* Make molecules whole at start of run */
1468 if (fr->ePBC != epbcNONE)
1470 do_pbc_first_mtop(fplog, inputrec->ePBC, box, mtop, state->x);
1474 /* Correct initial vsite positions are required
1475 * for the initial distribution in the domain decomposition
1476 * and for the initial shell prediction.
1478 construct_vsites_mtop(fplog, vsite, mtop, state->x);
1482 if (EEL_PME(fr->eeltype))
1484 ewaldcoeff = fr->ewaldcoeff;
1485 pmedata = &fr->pmedata;
1494 /* This is a PME only node */
1496 /* We don't need the state */
1499 ewaldcoeff = calc_ewaldcoeff(inputrec->rcoulomb, inputrec->ewald_rtol);
1503 if (hw_opt->thread_affinity != threadaffOFF)
1505 /* Before setting affinity, check whether the affinity has changed
1506 * - which indicates that probably the OpenMP library has changed it
1507 * since we first checked).
1509 gmx_check_thread_affinity_set(fplog, cr,
1510 hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1512 /* Set the CPU affinity */
1513 gmx_set_thread_affinity(fplog, cr, hw_opt, nthreads_pme, hwinfo, inputrec);
1516 /* Initiate PME if necessary,
1517 * either on all nodes or on dedicated PME nodes only. */
1518 if (EEL_PME(inputrec->coulombtype))
1522 nChargePerturbed = mdatoms->nChargePerturbed;
1524 if (cr->npmenodes > 0)
1526 /* The PME only nodes need to know nChargePerturbed */
1527 gmx_bcast_sim(sizeof(nChargePerturbed), &nChargePerturbed, cr);
1530 if (cr->duty & DUTY_PME)
1532 status = gmx_pme_init(pmedata, cr, npme_major, npme_minor, inputrec,
1533 mtop ? mtop->natoms : 0, nChargePerturbed,
1534 (Flags & MD_REPRODUCIBLE), nthreads_pme);
1537 gmx_fatal(FARGS, "Error %d initializing PME", status);
1543 if (integrator[inputrec->eI].func == do_md)
1545 /* Turn on signal handling on all nodes */
1547 * (A user signal from the PME nodes (if any)
1548 * is communicated to the PP nodes.
1550 signal_handler_install();
1553 if (cr->duty & DUTY_PP)
1555 if (inputrec->ePull != epullNO)
1557 /* Initialize pull code */
1558 init_pull(fplog, inputrec, nfile, fnm, mtop, cr, oenv, inputrec->fepvals->init_lambda,
1559 EI_DYNAMICS(inputrec->eI) && MASTER(cr), Flags);
1564 /* Initialize enforced rotation code */
1565 init_rot(fplog, inputrec, nfile, fnm, cr, state->x, box, mtop, oenv,
1569 constr = init_constraints(fplog, mtop, inputrec, ed, state, cr);
1571 if (DOMAINDECOMP(cr))
1573 dd_init_bondeds(fplog, cr->dd, mtop, vsite, constr, inputrec,
1574 Flags & MD_DDBONDCHECK, fr->cginfo_mb);
1576 set_dd_parameters(fplog, cr->dd, dlb_scale, inputrec, fr, &ddbox);
1578 setup_dd_grid(fplog, cr->dd);
1581 /* Now do whatever the user wants us to do (how flexible...) */
1582 integrator[inputrec->eI].func(fplog, cr, nfile, fnm,
1583 oenv, bVerbose, bCompact,
1586 nstepout, inputrec, mtop,
1588 mdatoms, nrnb, wcycle, ed, fr,
1589 repl_ex_nst, repl_ex_nex, repl_ex_seed,
1591 cpt_period, max_hours,
1596 if (inputrec->ePull != epullNO)
1598 finish_pull(fplog, inputrec->pull);
1603 finish_rot(inputrec->rot);
1610 gmx_pmeonly(*pmedata, cr, nrnb, wcycle, ewaldcoeff, FALSE, inputrec);
1613 if (EI_DYNAMICS(inputrec->eI) || EI_TPI(inputrec->eI))
1615 /* Some timing stats */
1618 if (runtime.proc == 0)
1620 runtime.proc = runtime.real;
1629 wallcycle_stop(wcycle, ewcRUN);
1631 /* Finish up, write some stuff
1632 * if rerunMD, don't write last frame again
1634 finish_run(fplog, cr, ftp2fn(efSTO, nfile, fnm),
1635 inputrec, nrnb, wcycle, &runtime,
1636 fr != NULL && fr->nbv != NULL && fr->nbv->bUseGPU ?
1637 nbnxn_cuda_get_timings(fr->nbv->cu_nbv) : NULL,
1639 EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1641 if ((cr->duty & DUTY_PP) && fr->nbv != NULL && fr->nbv->bUseGPU)
1643 char gpu_err_str[STRLEN];
1645 /* free GPU memory and uninitialize GPU (by destroying the context) */
1646 nbnxn_cuda_free(fplog, fr->nbv->cu_nbv);
1648 if (!free_gpu(gpu_err_str))
1650 gmx_warning("On node %d failed to free GPU #%d: %s",
1651 cr->nodeid, get_current_gpu_device_id(), gpu_err_str);
1655 if (opt2bSet("-membed", nfile, fnm))
1660 #ifdef GMX_THREAD_MPI
1661 if (PAR(cr) && SIMMASTER(cr))
1664 gmx_hardware_info_free(hwinfo);
1667 /* Does what it says */
1668 print_date_and_time(fplog, cr->nodeid, "Finished mdrun", &runtime);
1670 /* Close logfile already here if we were appending to it */
1671 if (MASTER(cr) && (Flags & MD_APPENDFILES))
1673 gmx_log_close(fplog);
1676 rc = (int)gmx_get_stop_condition();
1678 #ifdef GMX_THREAD_MPI
1679 /* we need to join all threads. The sub-threads join when they
1680 exit this function, but the master thread needs to be told to
1682 if (PAR(cr) && MASTER(cr))