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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
41 #if defined(HAVE_SCHED_H) && defined(HAVE_SCHED_GETAFFINITY)
44 #include <sys/syscall.h>
56 #include "md_logging.h"
57 #include "md_support.h"
67 #include "mpelogging.h"
73 #include "checkpoint.h"
74 #include "mtop_util.h"
75 #include "sighandler.h"
78 #include "gmx_detect_hardware.h"
79 #include "gmx_omp_nthreads.h"
80 #include "pull_rotation.h"
81 #include "calc_verletbuf.h"
82 #include "../mdlib/nbnxn_search.h"
83 #include "../mdlib/nbnxn_consts.h"
84 #include "gmx_fatal_collective.h"
88 #include "thread_mpi/threads.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_large_int_t deform_init_init_step_tpx;
112 matrix deform_init_box_tpx;
113 #ifdef GMX_THREAD_MPI
114 tMPI_Thread_mutex_t deform_init_box_mutex=TMPI_THREAD_MUTEX_INITIALIZER;
118 #ifdef GMX_THREAD_MPI
119 struct mdrunner_arglist
121 gmx_hw_opt_t *hw_opt;
134 const char *dddlb_opt;
139 const char *nbpu_opt;
150 const char *deviceOptions;
152 int ret; /* return value */
156 /* The function used for spawning threads. Extracts the mdrunner()
157 arguments from its one argument and calls mdrunner(), after making
159 static void mdrunner_start_fn(void *arg)
161 struct mdrunner_arglist *mda=(struct mdrunner_arglist*)arg;
162 struct mdrunner_arglist mc=*mda; /* copy the arg list to make sure
163 that it's thread-local. This doesn't
164 copy pointed-to items, of course,
165 but those are all const. */
166 t_commrec *cr; /* we need a local version of this */
170 fnm = dup_tfn(mc.nfile, mc.fnm);
172 cr = init_par_threads(mc.cr);
179 mda->ret=mdrunner(mc.hw_opt, fplog, cr, mc.nfile, fnm, mc.oenv,
180 mc.bVerbose, mc.bCompact, mc.nstglobalcomm,
181 mc.ddxyz, mc.dd_node_order, mc.rdd,
182 mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
183 mc.ddcsx, mc.ddcsy, mc.ddcsz,
185 mc.nsteps_cmdline, mc.nstepout, mc.resetstep,
186 mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_nex, mc.repl_ex_seed, mc.pforce,
187 mc.cpt_period, mc.max_hours, mc.deviceOptions, mc.Flags);
190 /* called by mdrunner() to start a specific number of threads (including
191 the main thread) for thread-parallel runs. This in turn calls mdrunner()
193 All options besides nthreads are the same as for mdrunner(). */
194 static t_commrec *mdrunner_start_threads(gmx_hw_opt_t *hw_opt,
195 FILE *fplog,t_commrec *cr,int nfile,
196 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
197 gmx_bool bCompact, int nstglobalcomm,
198 ivec ddxyz,int dd_node_order,real rdd,real rconstr,
199 const char *dddlb_opt,real dlb_scale,
200 const char *ddcsx,const char *ddcsy,const char *ddcsz,
201 const char *nbpu_opt,
202 int nsteps_cmdline, int nstepout,int resetstep,
203 int nmultisim,int repl_ex_nst,int repl_ex_nex, int repl_ex_seed,
204 real pforce,real cpt_period, real max_hours,
205 const char *deviceOptions, unsigned long Flags)
208 struct mdrunner_arglist *mda;
209 t_commrec *crn; /* the new commrec */
212 /* first check whether we even need to start tMPI */
213 if (hw_opt->nthreads_tmpi < 2)
218 /* a few small, one-time, almost unavoidable memory leaks: */
220 fnmn=dup_tfn(nfile, fnm);
222 /* fill the data structure to pass as void pointer to thread start fn */
229 mda->bVerbose=bVerbose;
230 mda->bCompact=bCompact;
231 mda->nstglobalcomm=nstglobalcomm;
232 mda->ddxyz[XX]=ddxyz[XX];
233 mda->ddxyz[YY]=ddxyz[YY];
234 mda->ddxyz[ZZ]=ddxyz[ZZ];
235 mda->dd_node_order=dd_node_order;
237 mda->rconstr=rconstr;
238 mda->dddlb_opt=dddlb_opt;
239 mda->dlb_scale=dlb_scale;
243 mda->nbpu_opt=nbpu_opt;
244 mda->nsteps_cmdline=nsteps_cmdline;
245 mda->nstepout=nstepout;
246 mda->resetstep=resetstep;
247 mda->nmultisim=nmultisim;
248 mda->repl_ex_nst=repl_ex_nst;
249 mda->repl_ex_nex=repl_ex_nex;
250 mda->repl_ex_seed=repl_ex_seed;
252 mda->cpt_period=cpt_period;
253 mda->max_hours=max_hours;
254 mda->deviceOptions=deviceOptions;
257 /* now spawn new threads that start mdrunner_start_fn(), while
258 the main thread returns */
259 ret=tMPI_Init_fn(TRUE, hw_opt->nthreads_tmpi,
260 (hw_opt->bThreadPinning ? TMPI_AFFINITY_ALL_CORES : TMPI_AFFINITY_NONE),
261 mdrunner_start_fn, (void*)(mda) );
262 if (ret!=TMPI_SUCCESS)
265 /* make a new comm_rec to reflect the new situation */
266 crn=init_par_threads(cr);
271 static int get_tmpi_omp_thread_division(const gmx_hw_info_t *hwinfo,
272 const gmx_hw_opt_t *hw_opt,
278 /* There are no separate PME nodes here, as we ensured in
279 * check_and_update_hw_opt that nthreads_tmpi>0 with PME nodes
280 * and a conditional ensures we would not have ended up here.
281 * Note that separate PME nodes might be switched on later.
285 nthreads_tmpi = ngpu;
286 if (nthreads_tot > 0 && nthreads_tot < nthreads_tmpi)
288 nthreads_tmpi = nthreads_tot;
291 else if (hw_opt->nthreads_omp > 0)
293 /* Here we could oversubscribe, when we do, we issue a warning later */
294 nthreads_tmpi = max(1,nthreads_tot/hw_opt->nthreads_omp);
298 /* TODO choose nthreads_omp based on hardware topology
299 when we have a hardware topology detection library */
300 /* In general, when running up to 4 threads, OpenMP should be faster.
301 * Note: on AMD Bulldozer we should avoid running OpenMP over two dies.
302 * On Intel>=Nehalem running OpenMP on a single CPU is always faster,
303 * even on two CPUs it's usually faster (but with many OpenMP threads
304 * it could be faster not to use HT, currently we always use HT).
305 * On Nehalem/Westmere we want to avoid running 16 threads over
306 * two CPUs with HT, so we need a limit<16; thus we use 12.
307 * A reasonable limit for Intel Sandy and Ivy bridge,
308 * not knowing the topology, is 16 threads.
310 const int nthreads_omp_always_faster = 4;
311 const int nthreads_omp_always_faster_Nehalem = 12;
312 const int nthreads_omp_always_faster_SandyBridge = 16;
313 const int first_model_Nehalem = 0x1A;
314 const int first_model_SandyBridge = 0x2A;
315 gmx_bool bIntel_Family6;
318 (gmx_cpuid_vendor(hwinfo->cpuid_info) == GMX_CPUID_VENDOR_INTEL &&
319 gmx_cpuid_family(hwinfo->cpuid_info) == 6);
321 if (nthreads_tot <= nthreads_omp_always_faster ||
323 ((gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_Nehalem && nthreads_tot <= nthreads_omp_always_faster_Nehalem) ||
324 (gmx_cpuid_model(hwinfo->cpuid_info) >= nthreads_omp_always_faster_SandyBridge && nthreads_tot <= nthreads_omp_always_faster_SandyBridge))))
326 /* Use pure OpenMP parallelization */
331 /* Don't use OpenMP parallelization */
332 nthreads_tmpi = nthreads_tot;
336 return nthreads_tmpi;
340 /* Get the number of threads to use for thread-MPI based on how many
341 * were requested, which algorithms we're using,
342 * and how many particles there are.
343 * At the point we have already called check_and_update_hw_opt.
344 * Thus all options should be internally consistent and consistent
345 * with the hardware, except that ntmpi could be larger than #GPU.
347 static int get_nthreads_mpi(gmx_hw_info_t *hwinfo,
348 gmx_hw_opt_t *hw_opt,
349 t_inputrec *inputrec, gmx_mtop_t *mtop,
353 int nthreads_hw,nthreads_tot_max,nthreads_tmpi,nthreads_new,ngpu;
354 int min_atoms_per_mpi_thread;
359 if (hw_opt->nthreads_tmpi > 0)
361 /* Trivial, return right away */
362 return hw_opt->nthreads_tmpi;
365 nthreads_hw = hwinfo->nthreads_hw_avail;
367 /* How many total (#tMPI*#OpenMP) threads can we start? */
368 if (hw_opt->nthreads_tot > 0)
370 nthreads_tot_max = hw_opt->nthreads_tot;
374 nthreads_tot_max = nthreads_hw;
377 bCanUseGPU = (inputrec->cutoff_scheme == ecutsVERLET && hwinfo->bCanUseGPU);
380 ngpu = hwinfo->gpu_info.ncuda_dev_use;
388 get_tmpi_omp_thread_division(hwinfo,hw_opt,nthreads_tot_max,ngpu);
390 if (inputrec->eI == eiNM || EI_TPI(inputrec->eI))
392 /* Steps are divided over the nodes iso splitting the atoms */
393 min_atoms_per_mpi_thread = 0;
399 min_atoms_per_mpi_thread = MIN_ATOMS_PER_GPU;
403 min_atoms_per_mpi_thread = MIN_ATOMS_PER_MPI_THREAD;
407 /* Check if an algorithm does not support parallel simulation. */
408 if (nthreads_tmpi != 1 &&
409 ( inputrec->eI == eiLBFGS ||
410 inputrec->coulombtype == eelEWALD ) )
414 md_print_warn(cr,fplog,"The integration or electrostatics algorithm doesn't support parallel runs. Using a single thread-MPI thread.\n");
415 if (hw_opt->nthreads_tmpi > nthreads_tmpi)
417 gmx_fatal(FARGS,"You asked for more than 1 thread-MPI thread, but an algorithm doesn't support that");
420 else if (mtop->natoms/nthreads_tmpi < min_atoms_per_mpi_thread)
422 /* the thread number was chosen automatically, but there are too many
423 threads (too few atoms per thread) */
424 nthreads_new = max(1,mtop->natoms/min_atoms_per_mpi_thread);
426 /* Avoid partial use of Hyper-Threading */
427 if (gmx_cpuid_x86_smt(hwinfo->cpuid_info) == GMX_CPUID_X86_SMT_ENABLED &&
428 nthreads_new > nthreads_hw/2 && nthreads_new < nthreads_hw)
430 nthreads_new = nthreads_hw/2;
433 /* Avoid large prime numbers in the thread count */
434 if (nthreads_new >= 6)
436 /* Use only 6,8,10 with additional factors of 2 */
440 while (3*fac*2 <= nthreads_new)
445 nthreads_new = (nthreads_new/fac)*fac;
450 if (nthreads_new == 5)
456 nthreads_tmpi = nthreads_new;
458 fprintf(stderr,"\n");
459 fprintf(stderr,"NOTE: Parallelization is limited by the small number of atoms,\n");
460 fprintf(stderr," only starting %d thread-MPI threads.\n",nthreads_tmpi);
461 fprintf(stderr," You can use the -nt and/or -ntmpi option to optimize the number of threads.\n\n");
464 return nthreads_tmpi;
466 #endif /* GMX_THREAD_MPI */
469 /* Environment variable for setting nstlist */
470 static const char* NSTLIST_ENVVAR = "GMX_NSTLIST";
471 /* Try to increase nstlist when using a GPU with nstlist less than this */
472 static const int NSTLIST_GPU_ENOUGH = 20;
473 /* Increase nstlist until the non-bonded cost increases more than this factor */
474 static const float NBNXN_GPU_LIST_OK_FAC = 1.25;
475 /* Don't increase nstlist beyond a non-bonded cost increases of this factor */
476 static const float NBNXN_GPU_LIST_MAX_FAC = 1.40;
478 /* Try to increase nstlist when running on a GPU */
479 static void increase_nstlist(FILE *fp,t_commrec *cr,
480 t_inputrec *ir,const gmx_mtop_t *mtop,matrix box)
483 int nstlist_orig,nstlist_prev;
484 verletbuf_list_setup_t ls;
485 real rlist_inc,rlist_ok,rlist_max,rlist_new,rlist_prev;
488 gmx_bool bBox,bDD,bCont;
489 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";
490 const char *vbd_err="Can not increase nstlist for GPU run because verlet-buffer-drift is not set or used";
491 const char *box_err="Can not increase nstlist for GPU run because the box is too small";
492 const char *dd_err ="Can not increase nstlist for GPU run because of domain decomposition limitations";
495 /* Number of + nstlist alternative values to try when switching */
496 const int nstl[]={ 20, 25, 40, 50 };
497 #define NNSTL sizeof(nstl)/sizeof(nstl[0])
499 env = getenv(NSTLIST_ENVVAR);
504 fprintf(fp,nstl_fmt,ir->nstlist);
508 if (ir->verletbuf_drift == 0)
510 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");
513 if (ir->verletbuf_drift < 0)
517 fprintf(stderr,"%s\n",vbd_err);
521 fprintf(fp,"%s\n",vbd_err);
527 nstlist_orig = ir->nstlist;
530 sprintf(buf,"Getting nstlist from environment variable GMX_NSTLIST=%s",env);
533 fprintf(stderr,"%s\n",buf);
537 fprintf(fp,"%s\n",buf);
539 sscanf(env,"%d",&ir->nstlist);
542 verletbuf_get_list_setup(TRUE,&ls);
544 /* Allow rlist to make the list double the size of the cut-off sphere */
545 rlist_inc = nbnxn_get_rlist_effective_inc(NBNXN_GPU_CLUSTER_SIZE,mtop->natoms/det(box));
546 rlist_ok = (max(ir->rvdw,ir->rcoulomb) + rlist_inc)*pow(NBNXN_GPU_LIST_OK_FAC,1.0/3.0) - rlist_inc;
547 rlist_max = (max(ir->rvdw,ir->rcoulomb) + rlist_inc)*pow(NBNXN_GPU_LIST_MAX_FAC,1.0/3.0) - rlist_inc;
550 fprintf(debug,"GPU nstlist tuning: rlist_inc %.3f rlist_max %.3f\n",
551 rlist_inc,rlist_max);
555 nstlist_prev = nstlist_orig;
556 rlist_prev = ir->rlist;
561 ir->nstlist = nstl[i];
564 /* Set the pair-list buffer size in ir */
565 calc_verlet_buffer_size(mtop,det(box),ir,ir->verletbuf_drift,&ls,
568 /* Does rlist fit in the box? */
569 bBox = (sqr(rlist_new) < max_cutoff2(ir->ePBC,box));
571 if (bBox && DOMAINDECOMP(cr))
573 /* Check if rlist fits in the domain decomposition */
574 if (inputrec2nboundeddim(ir) < DIM)
576 gmx_incons("Changing nstlist with domain decomposition and unbounded dimensions is not implemented yet");
578 copy_mat(box,state_tmp.box);
579 bDD = change_dd_cutoff(cr,&state_tmp,ir,rlist_new);
586 if (bBox && bDD && rlist_new <= rlist_max)
588 /* Increase nstlist */
589 nstlist_prev = ir->nstlist;
590 rlist_prev = rlist_new;
591 bCont = (i+1 < NNSTL && rlist_new < rlist_ok);
595 /* Stick with the previous nstlist */
596 ir->nstlist = nstlist_prev;
597 rlist_new = rlist_prev;
609 gmx_warning(!bBox ? box_err : dd_err);
612 fprintf(fp,"\n%s\n",bBox ? box_err : dd_err);
614 ir->nstlist = nstlist_orig;
616 else if (ir->nstlist != nstlist_orig || rlist_new != ir->rlist)
618 sprintf(buf,"Changing nstlist from %d to %d, rlist from %g to %g",
619 nstlist_orig,ir->nstlist,
620 ir->rlist,rlist_new);
623 fprintf(stderr,"%s\n\n",buf);
627 fprintf(fp,"%s\n\n",buf);
629 ir->rlist = rlist_new;
630 ir->rlistlong = rlist_new;
634 static void prepare_verlet_scheme(FILE *fplog,
635 gmx_hw_info_t *hwinfo,
637 gmx_hw_opt_t *hw_opt,
638 const char *nbpu_opt,
640 const gmx_mtop_t *mtop,
644 /* Here we only check for GPU usage on the MPI master process,
645 * as here we don't know how many GPUs we will use yet.
646 * We check for a GPU on all processes later.
648 *bUseGPU = hwinfo->bCanUseGPU || (getenv("GMX_EMULATE_GPU") != NULL);
650 if (ir->verletbuf_drift > 0)
652 /* Update the Verlet buffer size for the current run setup */
653 verletbuf_list_setup_t ls;
656 /* Here we assume CPU acceleration is on. But as currently
657 * calc_verlet_buffer_size gives the same results for 4x8 and 4x4
658 * and 4x2 gives a larger buffer than 4x4, this is ok.
660 verletbuf_get_list_setup(*bUseGPU,&ls);
662 calc_verlet_buffer_size(mtop,det(box),ir,
663 ir->verletbuf_drift,&ls,
665 if (rlist_new != ir->rlist)
669 fprintf(fplog,"\nChanging rlist from %g to %g for non-bonded %dx%d atom kernels\n\n",
671 ls.cluster_size_i,ls.cluster_size_j);
673 ir->rlist = rlist_new;
674 ir->rlistlong = rlist_new;
678 /* With GPU or emulation we should check nstlist for performance */
679 if ((EI_DYNAMICS(ir->eI) &&
681 ir->nstlist < NSTLIST_GPU_ENOUGH) ||
682 getenv(NSTLIST_ENVVAR) != NULL)
684 /* Choose a better nstlist */
685 increase_nstlist(fplog,cr,ir,mtop,box);
689 static void convert_to_verlet_scheme(FILE *fplog,
691 gmx_mtop_t *mtop,real box_vol)
693 char *conv_mesg="Converting input file with group cut-off scheme to the Verlet cut-off scheme";
695 md_print_warn(NULL,fplog,"%s\n",conv_mesg);
697 ir->cutoff_scheme = ecutsVERLET;
698 ir->verletbuf_drift = 0.005;
700 if (ir->rcoulomb != ir->rvdw)
702 gmx_fatal(FARGS,"The VdW and Coulomb cut-offs are different, whereas the Verlet scheme only supports equal cut-offs");
705 if (ir->vdwtype == evdwUSER || EEL_USER(ir->coulombtype))
707 gmx_fatal(FARGS,"User non-bonded potentials are not (yet) supported with the Verlet scheme");
709 else if (EVDW_SWITCHED(ir->vdwtype) || EEL_SWITCHED(ir->coulombtype))
711 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");
713 if (EVDW_SWITCHED(ir->vdwtype))
715 ir->vdwtype = evdwCUT;
717 if (EEL_SWITCHED(ir->coulombtype))
719 if (EEL_FULL(ir->coulombtype))
721 /* With full electrostatic only PME can be switched */
722 ir->coulombtype = eelPME;
726 md_print_warn(NULL,fplog,"NOTE: Replacing %s electrostatics with reaction-field with epsilon-rf=inf\n",eel_names[ir->coulombtype]);
727 ir->coulombtype = eelRF;
728 ir->epsilon_rf = 0.0;
732 /* We set the target energy drift to a small number.
733 * Note that this is only for testing. For production the user
734 * should think about this and set the mdp options.
736 ir->verletbuf_drift = 1e-4;
739 if (inputrec2nboundeddim(ir) != 3)
741 gmx_fatal(FARGS,"Can only convert old tpr files to the Verlet cut-off scheme with 3D pbc");
744 if (ir->efep != efepNO || ir->implicit_solvent != eisNO)
746 gmx_fatal(FARGS,"Will not convert old tpr files to the Verlet cut-off scheme with free-energy calculations or implicit solvent");
749 if (EI_DYNAMICS(ir->eI) && !(EI_MD(ir->eI) && ir->etc == etcNO))
751 verletbuf_list_setup_t ls;
753 verletbuf_get_list_setup(FALSE,&ls);
754 calc_verlet_buffer_size(mtop,box_vol,ir,ir->verletbuf_drift,&ls,
759 ir->verletbuf_drift = -1;
760 ir->rlist = 1.05*max(ir->rvdw,ir->rcoulomb);
763 gmx_mtop_remove_chargegroups(mtop);
766 /* Check the process affinity mask. If it is non-zero, something
767 * else has set the affinity, and mdrun should honor that and
768 * not attempt to do its own thread pinning.
770 * This function should be called twice. Once before the OpenMP
771 * library gets initialized with bAfterOpenMPInit=FALSE (which will
772 * detect affinity set by external tools like taskset), and again
773 * later, after the OpenMP initialization, with bAfterOpenMPInit=TRUE
774 * (which will detect affinity changes made by the OpenMP library).
776 * Note that this will only work on Linux, because we use a GNU
778 static void check_cpu_affinity_set(FILE *fplog, const t_commrec *cr,
779 gmx_hw_opt_t *hw_opt, int ncpus,
780 gmx_bool bAfterOpenmpInit)
782 #ifdef HAVE_SCHED_GETAFFINITY
783 cpu_set_t mask_current;
784 int i, ret, cpu_count, cpu_set;
788 if (!hw_opt->bThreadPinning)
790 /* internal affinity setting is off, don't bother checking process affinity */
794 CPU_ZERO(&mask_current);
795 if ((ret = sched_getaffinity(0, sizeof(cpu_set_t), &mask_current)) != 0)
797 /* failed to query affinity mask, will just return */
800 fprintf(debug, "Failed to query affinity mask (error %d)", ret);
805 /* Before proceeding with the actual check, make sure that the number of
806 * detected CPUs is >= the CPUs in the current set.
807 * We need to check for CPU_COUNT as it was added only in glibc 2.6. */
809 if (ncpus < CPU_COUNT(&mask_current))
813 fprintf(debug, "%d CPUs detected, but %d was returned by CPU_COUNT",
814 ncpus, CPU_COUNT(&mask_current));
818 #endif /* CPU_COUNT */
821 for (i = 0; (i < ncpus && i < CPU_SETSIZE); i++)
823 bAllSet = bAllSet && (CPU_ISSET(i, &mask_current) != 0);
828 if (!bAfterOpenmpInit)
830 md_print_warn(cr, fplog,
831 "%s detected a non-default process affinity, "
832 "so it will not attempt to pin its threads", ShortProgram());
836 md_print_warn(cr, fplog,
837 "%s detected a non-default process affinity, "
838 "probably set by the OpenMP library, "
839 "so it will not attempt to pin its threads", ShortProgram());
841 hw_opt->bThreadPinning = FALSE;
845 fprintf(debug, "Non-default affinity mask found, mdrun will not pin threads\n");
852 fprintf(debug, "Default affinity mask found\n");
855 #endif /* HAVE_SCHED_GETAFFINITY */
858 /* Set CPU affinity. Can be important for performance.
859 On some systems (e.g. Cray) CPU Affinity is set by default.
860 But default assigning doesn't work (well) with only some ranks
861 having threads. This causes very low performance.
862 External tools have cumbersome syntax for setting affinity
863 in the case that only some ranks have threads.
864 Thus it is important that GROMACS sets the affinity internally
865 if only PME is using threads.
867 static void set_cpu_affinity(FILE *fplog,
869 gmx_hw_opt_t *hw_opt,
871 const gmx_hw_info_t *hwinfo,
872 const t_inputrec *inputrec)
874 #if defined GMX_THREAD_MPI
875 /* With the number of TMPI threads equal to the number of cores
876 * we already pinned in thread-MPI, so don't pin again here.
878 if (hw_opt->nthreads_tmpi == tMPI_Thread_get_hw_number())
885 /* If the tMPI thread affinity setting is not supported encourage the user
886 * to report it as it's either a bug or an exotic platform which we might
887 * want to support. */
888 if (tMPI_Thread_setaffinity_support() != TMPI_SETAFFINITY_SUPPORT_YES)
890 md_print_warn(NULL, fplog,
891 "Can not set thread affinities on the current plarform. On NUMA systems this\n"
892 "can cause performance degradation. If you think your platform should support\n"
893 "setting affinities, contact the GROMACS developers.");
896 #endif /* __APPLE__ */
898 if (hw_opt->bThreadPinning)
900 int nth_affinity_set, thread_id_node, thread_id,
901 nthread_local, nthread_node, nthread_hw_max, nphyscore;
905 /* threads on this MPI process or TMPI thread */
906 if (cr->duty & DUTY_PP)
908 nthread_local = gmx_omp_nthreads_get(emntNonbonded);
912 nthread_local = gmx_omp_nthreads_get(emntPME);
915 /* map the current process to cores */
917 nthread_node = nthread_local;
919 if (PAR(cr) || MULTISIM(cr))
921 /* We need to determine a scan of the thread counts in this
926 MPI_Comm_split(MPI_COMM_WORLD,gmx_hostname_num(),cr->rank_intranode,
928 MPI_Scan(&nthread_local,&thread_id_node,1,MPI_INT,MPI_SUM,comm_intra);
929 /* MPI_Scan is inclusive, but here we need exclusive */
930 thread_id_node -= nthread_local;
931 /* Get the total number of threads on this physical node */
932 MPI_Allreduce(&nthread_local,&nthread_node,1,MPI_INT,MPI_SUM,comm_intra);
933 MPI_Comm_free(&comm_intra);
938 if (hw_opt->core_pinning_offset > 0)
940 offset = hw_opt->core_pinning_offset;
943 fprintf(stderr, "Applying core pinning offset %d\n", offset);
947 fprintf(fplog, "Applying core pinning offset %d\n", offset);
951 /* With Intel Hyper-Threading enabled, we want to pin consecutive
952 * threads to physical cores when using more threads than physical
953 * cores or when the user requests so.
955 nthread_hw_max = hwinfo->nthreads_hw_avail;
957 if (hw_opt->bPinHyperthreading ||
958 (gmx_cpuid_x86_smt(hwinfo->cpuid_info) == GMX_CPUID_X86_SMT_ENABLED &&
959 nthread_node > nthread_hw_max/2 && getenv("GMX_DISABLE_PINHT") == NULL))
961 if (gmx_cpuid_x86_smt(hwinfo->cpuid_info) != GMX_CPUID_X86_SMT_ENABLED)
963 /* We print to stderr on all processes, as we might have
964 * different settings on different physical nodes.
966 if (gmx_cpuid_vendor(hwinfo->cpuid_info) != GMX_CPUID_VENDOR_INTEL)
968 md_print_warn(NULL, fplog, "Pinning for Hyper-Threading layout requested, "
969 "but non-Intel CPU detected (vendor: %s)\n",
970 gmx_cpuid_vendor_string[gmx_cpuid_vendor(hwinfo->cpuid_info)]);
974 md_print_warn(NULL, fplog, "Pinning for Hyper-Threading layout requested, "
975 "but the CPU detected does not have Intel Hyper-Threading support "
976 "(or it is turned off)\n");
979 nphyscore = nthread_hw_max/2;
983 fprintf(stderr, "Pinning to Hyper-Threading cores with %d physical cores in a compute node\n",
988 fprintf(fplog, "Pinning to Hyper-Threading cores with %d physical cores in a compute node\n",
993 /* Set the per-thread affinity. In order to be able to check the success
994 * of affinity settings, we will set nth_affinity_set to 1 on threads
995 * where the affinity setting succeded and to 0 where it failed.
996 * Reducing these 0/1 values over the threads will give the total number
997 * of threads on which we succeeded.
999 nth_affinity_set = 0;
1000 #pragma omp parallel firstprivate(thread_id_node) num_threads(nthread_local) \
1001 reduction(+:nth_affinity_set)
1004 gmx_bool setaffinity_ret;
1006 thread_id = gmx_omp_get_thread_num();
1007 thread_id_node += thread_id;
1010 core = offset + thread_id_node;
1014 /* Lock pairs of threads to the same hyperthreaded core */
1015 core = offset + thread_id_node/2 + (thread_id_node % 2)*nphyscore;
1018 setaffinity_ret = tMPI_Thread_setaffinity_single(tMPI_Thread_self(), core);
1020 /* store the per-thread success-values of the setaffinity */
1021 nth_affinity_set = (setaffinity_ret == 0);
1025 fprintf(debug, "On rank %2d, thread %2d, core %2d the affinity setting returned %d\n",
1026 cr->nodeid, gmx_omp_get_thread_num(), core, setaffinity_ret);
1030 if (nth_affinity_set > nthread_local)
1034 sprintf(msg, "Looks like we have set affinity for more threads than "
1035 "we have (%d > %d)!\n", nth_affinity_set, nthread_local);
1040 /* check & warn if some threads failed to set their affinities */
1041 if (nth_affinity_set != nthread_local)
1043 char sbuf1[STRLEN], sbuf2[STRLEN];
1045 /* sbuf1 contains rank info, while sbuf2 OpenMP thread info */
1046 sbuf1[0] = sbuf2[0] = '\0';
1048 #ifdef GMX_THREAD_MPI
1049 sprintf(sbuf1, "In thread-MPI thread #%d: ", cr->nodeid);
1050 #else /* GMX_LIB_MPI */
1051 sprintf(sbuf1, "In MPI process #%d: ", cr->nodeid);
1053 #endif /* GMX_MPI */
1055 if (nthread_local > 1)
1057 sprintf(sbuf2, "of %d/%d thread%s ",
1058 nthread_local - nth_affinity_set, nthread_local,
1059 (nthread_local - nth_affinity_set) > 1 ? "s" : "");
1062 md_print_warn(NULL, fplog,
1063 "NOTE: %sAffinity setting %sfailed.\n"
1064 " This can cause performance degradation!",
1072 static void check_and_update_hw_opt(gmx_hw_opt_t *hw_opt,
1074 gmx_bool bIsSimMaster)
1076 gmx_omp_nthreads_read_env(&hw_opt->nthreads_omp, bIsSimMaster);
1078 #ifndef GMX_THREAD_MPI
1079 if (hw_opt->nthreads_tot > 0)
1081 gmx_fatal(FARGS,"Setting the total number of threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
1083 if (hw_opt->nthreads_tmpi > 0)
1085 gmx_fatal(FARGS,"Setting the number of thread-MPI threads is only supported with thread-MPI and Gromacs was compiled without thread-MPI");
1089 if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp_pme <= 0)
1091 /* We have the same number of OpenMP threads for PP and PME processes,
1092 * thus we can perform several consistency checks.
1094 if (hw_opt->nthreads_tmpi > 0 &&
1095 hw_opt->nthreads_omp > 0 &&
1096 hw_opt->nthreads_tot != hw_opt->nthreads_tmpi*hw_opt->nthreads_omp)
1098 gmx_fatal(FARGS,"The total number of threads requested (%d) does not match the thread-MPI threads (%d) times the OpenMP threads (%d) requested",
1099 hw_opt->nthreads_tot,hw_opt->nthreads_tmpi,hw_opt->nthreads_omp);
1102 if (hw_opt->nthreads_tmpi > 0 &&
1103 hw_opt->nthreads_tot % hw_opt->nthreads_tmpi != 0)
1105 gmx_fatal(FARGS,"The total number of threads requested (%d) is not divisible by the number of thread-MPI threads requested (%d)",
1106 hw_opt->nthreads_tot,hw_opt->nthreads_tmpi);
1109 if (hw_opt->nthreads_omp > 0 &&
1110 hw_opt->nthreads_tot % hw_opt->nthreads_omp != 0)
1112 gmx_fatal(FARGS,"The total number of threads requested (%d) is not divisible by the number of OpenMP threads requested (%d)",
1113 hw_opt->nthreads_tot,hw_opt->nthreads_omp);
1116 if (hw_opt->nthreads_tmpi > 0 &&
1117 hw_opt->nthreads_omp <= 0)
1119 hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
1124 if (hw_opt->nthreads_omp > 1)
1126 gmx_fatal(FARGS,"OpenMP threads are requested, but Gromacs was compiled without OpenMP support");
1130 if (cutoff_scheme == ecutsGROUP)
1132 /* We only have OpenMP support for PME only nodes */
1133 if (hw_opt->nthreads_omp > 1)
1135 gmx_fatal(FARGS,"OpenMP threads have been requested with cut-off scheme %s, but these are only supported with cut-off scheme %s",
1136 ecutscheme_names[cutoff_scheme],
1137 ecutscheme_names[ecutsVERLET]);
1139 hw_opt->nthreads_omp = 1;
1142 if (hw_opt->nthreads_omp_pme > 0 && hw_opt->nthreads_omp <= 0)
1144 gmx_fatal(FARGS,"You need to specify -ntomp in addition to -ntomp_pme");
1147 if (hw_opt->nthreads_tot == 1)
1149 hw_opt->nthreads_tmpi = 1;
1151 if (hw_opt->nthreads_omp > 1)
1153 gmx_fatal(FARGS,"You requested %d OpenMP threads with %d total threads",
1154 hw_opt->nthreads_tmpi,hw_opt->nthreads_tot);
1156 hw_opt->nthreads_omp = 1;
1159 if (hw_opt->nthreads_omp_pme <= 0 && hw_opt->nthreads_omp > 0)
1161 hw_opt->nthreads_omp_pme = hw_opt->nthreads_omp;
1166 fprintf(debug,"hw_opt: nt %d ntmpi %d ntomp %d ntomp_pme %d gpu_id '%s'\n",
1167 hw_opt->nthreads_tot,
1168 hw_opt->nthreads_tmpi,
1169 hw_opt->nthreads_omp,
1170 hw_opt->nthreads_omp_pme,
1171 hw_opt->gpu_id!=NULL ? hw_opt->gpu_id : "");
1177 /* Override the value in inputrec with value passed on the command line (if any) */
1178 static void override_nsteps_cmdline(FILE *fplog,
1181 const t_commrec *cr)
1186 /* override with anything else than the default -2 */
1187 if (nsteps_cmdline > -2)
1191 ir->nsteps = nsteps_cmdline;
1192 if (EI_DYNAMICS(ir->eI))
1194 sprintf(stmp, "Overriding nsteps with value passed on the command line: %d steps, %.3f ps",
1195 nsteps_cmdline, nsteps_cmdline*ir->delta_t);
1199 sprintf(stmp, "Overriding nsteps with value passed on the command line: %d steps",
1203 md_print_warn(cr, fplog, "%s\n", stmp);
1207 /* Data structure set by SIMMASTER which needs to be passed to all nodes
1208 * before the other nodes have read the tpx file and called gmx_detect_hardware.
1211 int cutoff_scheme; /* The cutoff scheme from inputrec_t */
1212 gmx_bool bUseGPU; /* Use GPU or GPU emulation */
1215 int mdrunner(gmx_hw_opt_t *hw_opt,
1216 FILE *fplog,t_commrec *cr,int nfile,
1217 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
1218 gmx_bool bCompact, int nstglobalcomm,
1219 ivec ddxyz,int dd_node_order,real rdd,real rconstr,
1220 const char *dddlb_opt,real dlb_scale,
1221 const char *ddcsx,const char *ddcsy,const char *ddcsz,
1222 const char *nbpu_opt,
1223 int nsteps_cmdline, int nstepout,int resetstep,
1224 int nmultisim,int repl_ex_nst,int repl_ex_nex,
1225 int repl_ex_seed, real pforce,real cpt_period,real max_hours,
1226 const char *deviceOptions, unsigned long Flags)
1228 gmx_bool bForceUseGPU,bTryUseGPU;
1229 double nodetime=0,realtime;
1230 t_inputrec *inputrec;
1231 t_state *state=NULL;
1233 gmx_ddbox_t ddbox={0};
1234 int npme_major,npme_minor;
1237 gmx_mtop_t *mtop=NULL;
1238 t_mdatoms *mdatoms=NULL;
1239 t_forcerec *fr=NULL;
1242 gmx_pme_t *pmedata=NULL;
1243 gmx_vsite_t *vsite=NULL;
1244 gmx_constr_t constr;
1245 int i,m,nChargePerturbed=-1,status,nalloc;
1247 gmx_wallcycle_t wcycle;
1248 gmx_bool bReadRNG,bReadEkin;
1250 gmx_runtime_t runtime;
1252 gmx_large_int_t reset_counters;
1253 gmx_edsam_t ed=NULL;
1254 t_commrec *cr_old=cr;
1257 gmx_membed_t membed=NULL;
1258 gmx_hw_info_t *hwinfo=NULL;
1259 master_inf_t minf={-1,FALSE};
1261 /* CAUTION: threads may be started later on in this function, so
1262 cr doesn't reflect the final parallel state right now */
1266 if (Flags & MD_APPENDFILES)
1271 bForceUseGPU = (strncmp(nbpu_opt, "gpu", 3) == 0);
1272 bTryUseGPU = (strncmp(nbpu_opt, "auto", 4) == 0) || bForceUseGPU;
1277 /* Read (nearly) all data required for the simulation */
1278 read_tpx_state(ftp2fn(efTPX,nfile,fnm),inputrec,state,NULL,mtop);
1280 if (inputrec->cutoff_scheme != ecutsVERLET &&
1281 ((Flags & MD_TESTVERLET) || getenv("GMX_VERLET_SCHEME") != NULL))
1283 convert_to_verlet_scheme(fplog,inputrec,mtop,det(state->box));
1286 /* Detect hardware, gather information. With tMPI only thread 0 does it
1287 * and after threads are started broadcasts hwinfo around. */
1289 gmx_detect_hardware(fplog, hwinfo, cr,
1290 bForceUseGPU, bTryUseGPU, hw_opt->gpu_id);
1292 minf.cutoff_scheme = inputrec->cutoff_scheme;
1293 minf.bUseGPU = FALSE;
1295 if (inputrec->cutoff_scheme == ecutsVERLET)
1297 prepare_verlet_scheme(fplog,hwinfo,cr,hw_opt,nbpu_opt,
1298 inputrec,mtop,state->box,
1301 else if (hwinfo->bCanUseGPU)
1303 md_print_warn(cr,fplog,
1304 "NOTE: GPU(s) found, but the current simulation can not use GPUs\n"
1305 " To use a GPU, set the mdp option: cutoff-scheme = Verlet\n"
1306 " (for quick performance testing you can use the -testverlet option)\n");
1310 gmx_fatal(FARGS,"GPU requested, but can't be used without cutoff-scheme=Verlet");
1314 #ifndef GMX_THREAD_MPI
1317 gmx_bcast_sim(sizeof(minf),&minf,cr);
1320 if (minf.bUseGPU && cr->npmenodes == -1)
1322 /* Don't automatically use PME-only nodes with GPUs */
1326 /* Check for externally set OpenMP affinity and turn off internal
1327 * pinning if any is found. We need to do this check early to tell
1328 * thread-MPI whether it should do pinning when spawning threads.
1330 gmx_omp_check_thread_affinity(fplog, cr, hw_opt);
1332 #ifdef GMX_THREAD_MPI
1333 /* With thread-MPI inputrec is only set here on the master thread */
1337 check_and_update_hw_opt(hw_opt,minf.cutoff_scheme,SIMMASTER(cr));
1339 #ifdef GMX_THREAD_MPI
1340 /* Early check for externally set process affinity. Can't do over all
1341 * MPI processes because hwinfo is not available everywhere, but with
1342 * thread-MPI it's needed as pinning might get turned off which needs
1343 * to be known before starting thread-MPI. */
1344 check_cpu_affinity_set(fplog,
1346 hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1349 #ifdef GMX_THREAD_MPI
1350 if (cr->npmenodes > 0 && hw_opt->nthreads_tmpi <= 0)
1352 gmx_fatal(FARGS,"You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME nodes");
1356 if (hw_opt->nthreads_omp_pme != hw_opt->nthreads_omp &&
1359 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");
1363 #ifdef GMX_THREAD_MPI
1366 /* NOW the threads will be started: */
1367 hw_opt->nthreads_tmpi = get_nthreads_mpi(hwinfo,
1371 if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp <= 0)
1373 hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
1376 if (hw_opt->nthreads_tmpi > 1)
1378 /* now start the threads. */
1379 cr=mdrunner_start_threads(hw_opt, fplog, cr_old, nfile, fnm,
1380 oenv, bVerbose, bCompact, nstglobalcomm,
1381 ddxyz, dd_node_order, rdd, rconstr,
1382 dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
1384 nsteps_cmdline, nstepout, resetstep, nmultisim,
1385 repl_ex_nst, repl_ex_nex, repl_ex_seed, pforce,
1386 cpt_period, max_hours, deviceOptions,
1388 /* the main thread continues here with a new cr. We don't deallocate
1389 the old cr because other threads may still be reading it. */
1392 gmx_comm("Failed to spawn threads");
1397 /* END OF CAUTION: cr is now reliable */
1399 /* g_membed initialisation *
1400 * Because we change the mtop, init_membed is called before the init_parallel *
1401 * (in case we ever want to make it run in parallel) */
1402 if (opt2bSet("-membed",nfile,fnm))
1406 fprintf(stderr,"Initializing membed");
1408 membed = init_membed(fplog,nfile,fnm,mtop,inputrec,state,cr,&cpt_period);
1413 /* now broadcast everything to the non-master nodes/threads: */
1414 init_parallel(fplog, cr, inputrec, mtop);
1416 /* This check needs to happen after get_nthreads_mpi() */
1417 if (inputrec->cutoff_scheme == ecutsVERLET && (Flags & MD_PARTDEC))
1419 gmx_fatal_collective(FARGS,cr,NULL,
1420 "The Verlet cut-off scheme is not supported with particle decomposition.\n"
1421 "You can achieve the same effect as particle decomposition by running in parallel using only OpenMP threads.");
1426 pr_inputrec(fplog,0,"Input Parameters",inputrec,FALSE);
1429 #if defined GMX_THREAD_MPI
1430 /* With tMPI we detected on thread 0 and we'll just pass the hwinfo pointer
1431 * to the other threads -- slightly uncool, but works fine, just need to
1432 * make sure that the data doesn't get freed twice. */
1439 gmx_bcast(sizeof(&hwinfo), &hwinfo, cr);
1442 if (PAR(cr) && !SIMMASTER(cr))
1444 /* now we have inputrec on all nodes, can run the detection */
1445 /* TODO: perhaps it's better to propagate within a node instead? */
1447 gmx_detect_hardware(fplog, hwinfo, cr,
1448 bForceUseGPU, bTryUseGPU, hw_opt->gpu_id);
1451 /* Now do the affinity check with MPI/no-MPI (done earlier with thread-MPI). */
1452 check_cpu_affinity_set(fplog, cr,
1453 hw_opt, hwinfo->nthreads_hw_avail, FALSE);
1456 /* now make sure the state is initialized and propagated */
1457 set_state_entries(state,inputrec,cr->nnodes);
1459 /* A parallel command line option consistency check that we can
1460 only do after any threads have started. */
1462 (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
1465 "The -dd or -npme option request a parallel simulation, "
1467 "but %s was compiled without threads or MPI enabled"
1469 #ifdef GMX_THREAD_MPI
1470 "but the number of threads (option -nt) is 1"
1472 "but %s was not started through mpirun/mpiexec or only one process was requested through mpirun/mpiexec"
1479 if ((Flags & MD_RERUN) &&
1480 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
1482 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
1485 if (can_use_allvsall(inputrec,mtop,TRUE,cr,fplog) && PAR(cr))
1487 /* All-vs-all loops do not work with domain decomposition */
1488 Flags |= MD_PARTDEC;
1491 if (!EEL_PME(inputrec->coulombtype) || (Flags & MD_PARTDEC))
1493 if (cr->npmenodes > 0)
1495 if (!EEL_PME(inputrec->coulombtype))
1497 gmx_fatal_collective(FARGS,cr,NULL,
1498 "PME nodes are requested, but the system does not use PME electrostatics");
1500 if (Flags & MD_PARTDEC)
1502 gmx_fatal_collective(FARGS,cr,NULL,
1503 "PME nodes are requested, but particle decomposition does not support separate PME nodes");
1511 fcRegisterSteps(inputrec->nsteps,inputrec->init_step);
1514 /* NMR restraints must be initialized before load_checkpoint,
1515 * since with time averaging the history is added to t_state.
1516 * For proper consistency check we therefore need to extend
1518 * So the PME-only nodes (if present) will also initialize
1519 * the distance restraints.
1523 /* This needs to be called before read_checkpoint to extend the state */
1524 init_disres(fplog,mtop,inputrec,cr,Flags & MD_PARTDEC,fcd,state, repl_ex_nst > 0);
1526 if (gmx_mtop_ftype_count(mtop,F_ORIRES) > 0)
1528 if (PAR(cr) && !(Flags & MD_PARTDEC))
1530 gmx_fatal(FARGS,"Orientation restraints do not work (yet) with domain decomposition, use particle decomposition (mdrun option -pd)");
1532 /* Orientation restraints */
1535 init_orires(fplog,mtop,state->x,inputrec,cr->ms,&(fcd->orires),
1540 if (DEFORM(*inputrec))
1542 /* Store the deform reference box before reading the checkpoint */
1545 copy_mat(state->box,box);
1549 gmx_bcast(sizeof(box),box,cr);
1551 /* Because we do not have the update struct available yet
1552 * in which the reference values should be stored,
1553 * we store them temporarily in static variables.
1554 * This should be thread safe, since they are only written once
1555 * and with identical values.
1557 #ifdef GMX_THREAD_MPI
1558 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
1560 deform_init_init_step_tpx = inputrec->init_step;
1561 copy_mat(box,deform_init_box_tpx);
1562 #ifdef GMX_THREAD_MPI
1563 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
1567 if (opt2bSet("-cpi",nfile,fnm))
1569 /* Check if checkpoint file exists before doing continuation.
1570 * This way we can use identical input options for the first and subsequent runs...
1572 if( gmx_fexist_master(opt2fn_master("-cpi",nfile,fnm,cr),cr) )
1574 load_checkpoint(opt2fn_master("-cpi",nfile,fnm,cr),&fplog,
1575 cr,Flags & MD_PARTDEC,ddxyz,
1576 inputrec,state,&bReadRNG,&bReadEkin,
1577 (Flags & MD_APPENDFILES),
1578 (Flags & MD_APPENDFILESSET));
1582 Flags |= MD_READ_RNG;
1586 Flags |= MD_READ_EKIN;
1591 if (((MASTER(cr) || (Flags & MD_SEPPOT)) && (Flags & MD_APPENDFILES))
1592 #ifdef GMX_THREAD_MPI
1593 /* With thread MPI only the master node/thread exists in mdrun.c,
1594 * therefore non-master nodes need to open the "seppot" log file here.
1596 || (!MASTER(cr) && (Flags & MD_SEPPOT))
1600 gmx_log_open(ftp2fn(efLOG,nfile,fnm),cr,!(Flags & MD_SEPPOT),
1604 /* override nsteps with value from cmdline */
1605 override_nsteps_cmdline(fplog, nsteps_cmdline, inputrec, cr);
1609 copy_mat(state->box,box);
1614 gmx_bcast(sizeof(box),box,cr);
1617 /* Essential dynamics */
1618 if (opt2bSet("-ei",nfile,fnm))
1620 /* Open input and output files, allocate space for ED data structure */
1621 ed = ed_open(mtop->natoms,&state->edsamstate,nfile,fnm,Flags,oenv,cr);
1624 if (PAR(cr) && !((Flags & MD_PARTDEC) ||
1625 EI_TPI(inputrec->eI) ||
1626 inputrec->eI == eiNM))
1628 cr->dd = init_domain_decomposition(fplog,cr,Flags,ddxyz,rdd,rconstr,
1629 dddlb_opt,dlb_scale,
1633 &ddbox,&npme_major,&npme_minor);
1635 make_dd_communicators(fplog,cr,dd_node_order);
1637 /* Set overallocation to avoid frequent reallocation of arrays */
1638 set_over_alloc_dd(TRUE);
1642 /* PME, if used, is done on all nodes with 1D decomposition */
1644 cr->duty = (DUTY_PP | DUTY_PME);
1647 if (!EI_TPI(inputrec->eI))
1649 npme_major = cr->nnodes;
1652 if (inputrec->ePBC == epbcSCREW)
1655 "pbc=%s is only implemented with domain decomposition",
1656 epbc_names[inputrec->ePBC]);
1662 /* After possible communicator splitting in make_dd_communicators.
1663 * we can set up the intra/inter node communication.
1665 gmx_setup_nodecomm(fplog,cr);
1668 /* Initialize per-physical-node MPI process/thread ID and counters. */
1669 gmx_init_intranode_counters(cr);
1672 md_print_info(cr,fplog,"Using %d MPI %s\n",
1674 #ifdef GMX_THREAD_MPI
1675 cr->nnodes==1 ? "thread" : "threads"
1677 cr->nnodes==1 ? "process" : "processes"
1683 gmx_omp_nthreads_init(fplog, cr,
1684 hwinfo->nthreads_hw_avail,
1685 hw_opt->nthreads_omp,
1686 hw_opt->nthreads_omp_pme,
1687 (cr->duty & DUTY_PP) == 0,
1688 inputrec->cutoff_scheme == ecutsVERLET);
1690 gmx_check_hw_runconf_consistency(fplog, hwinfo, cr, hw_opt->nthreads_tmpi, minf.bUseGPU);
1692 /* getting number of PP/PME threads
1693 PME: env variable should be read only on one node to make sure it is
1694 identical everywhere;
1696 /* TODO nthreads_pp is only used for pinning threads.
1697 * This is a temporary solution until we have a hw topology library.
1699 nthreads_pp = gmx_omp_nthreads_get(emntNonbonded);
1700 nthreads_pme = gmx_omp_nthreads_get(emntPME);
1702 wcycle = wallcycle_init(fplog,resetstep,cr,nthreads_pp,nthreads_pme);
1706 /* Master synchronizes its value of reset_counters with all nodes
1707 * including PME only nodes */
1708 reset_counters = wcycle_get_reset_counters(wcycle);
1709 gmx_bcast_sim(sizeof(reset_counters),&reset_counters,cr);
1710 wcycle_set_reset_counters(wcycle, reset_counters);
1714 if (cr->duty & DUTY_PP)
1716 /* For domain decomposition we allocate dynamically
1717 * in dd_partition_system.
1719 if (DOMAINDECOMP(cr))
1721 bcast_state_setup(cr,state);
1727 bcast_state(cr,state,TRUE);
1731 /* Initiate forcerecord */
1733 fr->hwinfo = hwinfo;
1734 init_forcerec(fplog,oenv,fr,fcd,inputrec,mtop,cr,box,FALSE,
1735 opt2fn("-table",nfile,fnm),
1736 opt2fn("-tabletf",nfile,fnm),
1737 opt2fn("-tablep",nfile,fnm),
1738 opt2fn("-tableb",nfile,fnm),
1742 /* version for PCA_NOT_READ_NODE (see md.c) */
1743 /*init_forcerec(fplog,fr,fcd,inputrec,mtop,cr,box,FALSE,
1744 "nofile","nofile","nofile","nofile",FALSE,pforce);
1746 fr->bSepDVDL = ((Flags & MD_SEPPOT) == MD_SEPPOT);
1748 /* Initialize QM-MM */
1751 init_QMMMrec(cr,box,mtop,inputrec,fr);
1754 /* Initialize the mdatoms structure.
1755 * mdatoms is not filled with atom data,
1756 * as this can not be done now with domain decomposition.
1758 mdatoms = init_mdatoms(fplog,mtop,inputrec->efep!=efepNO);
1760 /* Initialize the virtual site communication */
1761 vsite = init_vsite(mtop,cr,FALSE);
1763 calc_shifts(box,fr->shift_vec);
1765 /* With periodic molecules the charge groups should be whole at start up
1766 * and the virtual sites should not be far from their proper positions.
1768 if (!inputrec->bContinuation && MASTER(cr) &&
1769 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
1771 /* Make molecules whole at start of run */
1772 if (fr->ePBC != epbcNONE)
1774 do_pbc_first_mtop(fplog,inputrec->ePBC,box,mtop,state->x);
1778 /* Correct initial vsite positions are required
1779 * for the initial distribution in the domain decomposition
1780 * and for the initial shell prediction.
1782 construct_vsites_mtop(fplog,vsite,mtop,state->x);
1786 if (EEL_PME(fr->eeltype))
1788 ewaldcoeff = fr->ewaldcoeff;
1789 pmedata = &fr->pmedata;
1798 /* This is a PME only node */
1800 /* We don't need the state */
1803 ewaldcoeff = calc_ewaldcoeff(inputrec->rcoulomb, inputrec->ewald_rtol);
1807 /* Before setting affinity, check whether the affinity has changed
1808 * - which indicates that probably the OpenMP library has changed it since
1809 * we first checked). */
1810 check_cpu_affinity_set(fplog, cr, hw_opt, hwinfo->nthreads_hw_avail, TRUE);
1812 /* Set the CPU affinity */
1813 set_cpu_affinity(fplog,cr,hw_opt,nthreads_pme,hwinfo,inputrec);
1815 /* Initiate PME if necessary,
1816 * either on all nodes or on dedicated PME nodes only. */
1817 if (EEL_PME(inputrec->coulombtype))
1821 nChargePerturbed = mdatoms->nChargePerturbed;
1823 if (cr->npmenodes > 0)
1825 /* The PME only nodes need to know nChargePerturbed */
1826 gmx_bcast_sim(sizeof(nChargePerturbed),&nChargePerturbed,cr);
1829 if (cr->duty & DUTY_PME)
1831 status = gmx_pme_init(pmedata,cr,npme_major,npme_minor,inputrec,
1832 mtop ? mtop->natoms : 0,nChargePerturbed,
1833 (Flags & MD_REPRODUCIBLE),nthreads_pme);
1836 gmx_fatal(FARGS,"Error %d initializing PME",status);
1842 if (integrator[inputrec->eI].func == do_md)
1844 /* Turn on signal handling on all nodes */
1846 * (A user signal from the PME nodes (if any)
1847 * is communicated to the PP nodes.
1849 signal_handler_install();
1852 if (cr->duty & DUTY_PP)
1854 if (inputrec->ePull != epullNO)
1856 /* Initialize pull code */
1857 init_pull(fplog,inputrec,nfile,fnm,mtop,cr,oenv, inputrec->fepvals->init_lambda,
1858 EI_DYNAMICS(inputrec->eI) && MASTER(cr),Flags);
1863 /* Initialize enforced rotation code */
1864 init_rot(fplog,inputrec,nfile,fnm,cr,state->x,box,mtop,oenv,
1868 constr = init_constraints(fplog,mtop,inputrec,ed,state,cr);
1870 if (DOMAINDECOMP(cr))
1872 dd_init_bondeds(fplog,cr->dd,mtop,vsite,constr,inputrec,
1873 Flags & MD_DDBONDCHECK,fr->cginfo_mb);
1875 set_dd_parameters(fplog,cr->dd,dlb_scale,inputrec,fr,&ddbox);
1877 setup_dd_grid(fplog,cr->dd);
1880 /* Now do whatever the user wants us to do (how flexible...) */
1881 integrator[inputrec->eI].func(fplog,cr,nfile,fnm,
1882 oenv,bVerbose,bCompact,
1885 nstepout,inputrec,mtop,
1887 mdatoms,nrnb,wcycle,ed,fr,
1888 repl_ex_nst,repl_ex_nex,repl_ex_seed,
1890 cpt_period,max_hours,
1895 if (inputrec->ePull != epullNO)
1897 finish_pull(fplog,inputrec->pull);
1902 finish_rot(fplog,inputrec->rot);
1909 gmx_pmeonly(*pmedata,cr,nrnb,wcycle,ewaldcoeff,FALSE,inputrec);
1912 if (EI_DYNAMICS(inputrec->eI) || EI_TPI(inputrec->eI))
1914 /* Some timing stats */
1917 if (runtime.proc == 0)
1919 runtime.proc = runtime.real;
1928 wallcycle_stop(wcycle,ewcRUN);
1930 /* Finish up, write some stuff
1931 * if rerunMD, don't write last frame again
1933 finish_run(fplog,cr,ftp2fn(efSTO,nfile,fnm),
1934 inputrec,nrnb,wcycle,&runtime,
1935 fr != NULL && fr->nbv != NULL && fr->nbv->bUseGPU ?
1936 nbnxn_cuda_get_timings(fr->nbv->cu_nbv) : NULL,
1938 EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
1940 if ((cr->duty & DUTY_PP) && fr->nbv != NULL && fr->nbv->bUseGPU)
1942 char gpu_err_str[STRLEN];
1944 /* free GPU memory and uninitialize GPU (by destroying the context) */
1945 nbnxn_cuda_free(fplog, fr->nbv->cu_nbv);
1947 if (!free_gpu(gpu_err_str))
1949 gmx_warning("On node %d failed to free GPU #%d: %s",
1950 cr->nodeid, get_current_gpu_device_id(), gpu_err_str);
1954 if (opt2bSet("-membed",nfile,fnm))
1959 #ifdef GMX_THREAD_MPI
1960 if (PAR(cr) && SIMMASTER(cr))
1963 gmx_hardware_info_free(hwinfo);
1966 /* Does what it says */
1967 print_date_and_time(fplog,cr->nodeid,"Finished mdrun",&runtime);
1969 /* Close logfile already here if we were appending to it */
1970 if (MASTER(cr) && (Flags & MD_APPENDFILES))
1972 gmx_log_close(fplog);
1975 rc=(int)gmx_get_stop_condition();
1977 #ifdef GMX_THREAD_MPI
1978 /* we need to join all threads. The sub-threads join when they
1979 exit this function, but the master thread needs to be told to
1981 if (PAR(cr) && MASTER(cr))