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
43 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
64 #include "mpelogging.h"
70 #include "checkpoint.h"
71 #include "mtop_util.h"
72 #include "sighandler.h"
77 #include "md_openmm.h"
91 #include "md_openmm.h"
96 gmx_integrator_t *func;
99 /* The array should match the eI array in include/types/enums.h */
100 #ifdef GMX_OPENMM /* FIXME do_md_openmm needs fixing */
101 const gmx_intp_t integrator[eiNR] = { {do_md_openmm}, {do_md_openmm}, {do_md_openmm}, {do_md_openmm}, {do_md_openmm}, {do_md_openmm}, {do_md_openmm}, {do_md_openmm}, {do_md_openmm}, {do_md_openmm}, {do_md_openmm},{do_md_openmm}};
103 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}};
106 gmx_large_int_t deform_init_init_step_tpx;
107 matrix deform_init_box_tpx;
109 tMPI_Thread_mutex_t deform_init_box_mutex=TMPI_THREAD_MUTEX_INITIALIZER;
114 struct mdrunner_arglist
128 const char *dddlb_opt;
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(cr->nnodes, 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, mc.nstepout, mc.resetstep,
175 mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_seed, mc.pforce,
176 mc.cpt_period, mc.max_hours, mc.deviceOptions, mc.Flags);
179 /* called by mdrunner() to start a specific number of threads (including
180 the main thread) for thread-parallel runs. This in turn calls mdrunner()
182 All options besides nthreads are the same as for mdrunner(). */
183 static t_commrec *mdrunner_start_threads(int nthreads,
184 FILE *fplog,t_commrec *cr,int nfile,
185 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
186 gmx_bool bCompact, int nstglobalcomm,
187 ivec ddxyz,int dd_node_order,real rdd,real rconstr,
188 const char *dddlb_opt,real dlb_scale,
189 const char *ddcsx,const char *ddcsy,const char *ddcsz,
190 int nstepout,int resetstep,int nmultisim,int repl_ex_nst,
191 int repl_ex_seed, real pforce,real cpt_period, real max_hours,
192 const char *deviceOptions, unsigned long Flags)
195 struct mdrunner_arglist *mda;
196 t_commrec *crn; /* the new commrec */
199 /* first check whether we even need to start tMPI */
203 /* a few small, one-time, almost unavoidable memory leaks: */
205 fnmn=dup_tfn(nfile, fnm);
207 /* fill the data structure to pass as void pointer to thread start fn */
213 mda->bVerbose=bVerbose;
214 mda->bCompact=bCompact;
215 mda->nstglobalcomm=nstglobalcomm;
216 mda->ddxyz[XX]=ddxyz[XX];
217 mda->ddxyz[YY]=ddxyz[YY];
218 mda->ddxyz[ZZ]=ddxyz[ZZ];
219 mda->dd_node_order=dd_node_order;
221 mda->rconstr=rconstr;
222 mda->dddlb_opt=dddlb_opt;
223 mda->dlb_scale=dlb_scale;
227 mda->nstepout=nstepout;
228 mda->resetstep=resetstep;
229 mda->nmultisim=nmultisim;
230 mda->repl_ex_nst=repl_ex_nst;
231 mda->repl_ex_seed=repl_ex_seed;
233 mda->cpt_period=cpt_period;
234 mda->max_hours=max_hours;
235 mda->deviceOptions=deviceOptions;
238 fprintf(stderr, "Starting %d threads\n",nthreads);
240 /* now spawn new threads that start mdrunner_start_fn(), while
241 the main thread returns */
242 ret=tMPI_Init_fn(TRUE, nthreads, mdrunner_start_fn, (void*)(mda) );
243 if (ret!=TMPI_SUCCESS)
246 /* make a new comm_rec to reflect the new situation */
247 crn=init_par_threads(cr);
252 /* get the number of threads based on how many there were requested,
253 which algorithms we're using, and how many particles there are. */
254 static int get_nthreads(int nthreads_requested, t_inputrec *inputrec,
257 int nthreads,nthreads_new;
258 int min_atoms_per_thread;
261 nthreads = nthreads_requested;
263 /* determine # of hardware threads. */
264 if (nthreads_requested < 1)
266 if ((env = getenv("GMX_MAX_THREADS")) != NULL)
269 sscanf(env,"%d",&nthreads);
272 gmx_fatal(FARGS,"GMX_MAX_THREADS (%d) should be larger than 0",
278 nthreads = tMPI_Get_recommended_nthreads();
282 if (inputrec->eI == eiNM || EI_TPI(inputrec->eI))
284 /* Steps are divided over the nodes iso splitting the atoms */
285 min_atoms_per_thread = 0;
289 min_atoms_per_thread = MIN_ATOMS_PER_THREAD;
292 /* Check if an algorithm does not support parallel simulation. */
294 ( inputrec->eI == eiLBFGS ||
295 inputrec->coulombtype == eelEWALD ) )
297 fprintf(stderr,"\nThe integration or electrostatics algorithm doesn't support parallel runs. Not starting any threads.\n");
300 else if (nthreads_requested < 1 &&
301 mtop->natoms/nthreads < min_atoms_per_thread)
303 /* the thread number was chosen automatically, but there are too many
304 threads (too few atoms per thread) */
305 nthreads_new = max(1,mtop->natoms/min_atoms_per_thread);
307 if (nthreads_new > 8 || (nthreads == 8 && nthreads_new > 4))
309 /* Use only multiples of 4 above 8 threads
310 * or with an 8-core processor
311 * (to avoid 6 threads on 8 core processors with 4 real cores).
313 nthreads_new = (nthreads_new/4)*4;
315 else if (nthreads_new > 4)
317 /* Avoid 5 or 7 threads */
318 nthreads_new = (nthreads_new/2)*2;
321 nthreads = nthreads_new;
323 fprintf(stderr,"\n");
324 fprintf(stderr,"NOTE: Parallelization is limited by the small number of atoms,\n");
325 fprintf(stderr," only starting %d threads.\n",nthreads);
326 fprintf(stderr," You can use the -nt option to optimize the number of threads.\n\n");
333 int mdrunner(int nthreads_requested, FILE *fplog,t_commrec *cr,int nfile,
334 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
335 gmx_bool bCompact, int nstglobalcomm,
336 ivec ddxyz,int dd_node_order,real rdd,real rconstr,
337 const char *dddlb_opt,real dlb_scale,
338 const char *ddcsx,const char *ddcsy,const char *ddcsz,
339 int nstepout,int resetstep,int nmultisim,int repl_ex_nst,
340 int repl_ex_seed, real pforce,real cpt_period,real max_hours,
341 const char *deviceOptions, unsigned long Flags)
343 double nodetime=0,realtime;
344 t_inputrec *inputrec;
347 gmx_ddbox_t ddbox={0};
348 int npme_major,npme_minor;
351 gmx_mtop_t *mtop=NULL;
352 t_mdatoms *mdatoms=NULL;
356 gmx_pme_t *pmedata=NULL;
357 gmx_vsite_t *vsite=NULL;
359 int i,m,nChargePerturbed=-1,status,nalloc;
361 gmx_wallcycle_t wcycle;
362 gmx_bool bReadRNG,bReadEkin;
364 gmx_runtime_t runtime;
366 gmx_large_int_t reset_counters;
368 t_commrec *cr_old=cr;
370 gmx_membed_t *membed=NULL;
372 /* CAUTION: threads may be started later on in this function, so
373 cr doesn't reflect the final parallel state right now */
377 if (bVerbose && SIMMASTER(cr))
379 fprintf(stderr,"Getting Loaded...\n");
382 if (Flags & MD_APPENDFILES)
390 /* Read (nearly) all data required for the simulation */
391 read_tpx_state(ftp2fn(efTPX,nfile,fnm),inputrec,state,NULL,mtop);
393 /* NOW the threads will be started: */
395 nthreads = get_nthreads(nthreads_requested, inputrec, mtop);
399 /* now start the threads. */
400 cr=mdrunner_start_threads(nthreads, fplog, cr_old, nfile, fnm,
401 oenv, bVerbose, bCompact, nstglobalcomm,
402 ddxyz, dd_node_order, rdd, rconstr,
403 dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
404 nstepout, resetstep, nmultisim,
405 repl_ex_nst, repl_ex_seed, pforce,
406 cpt_period, max_hours, deviceOptions,
408 /* the main thread continues here with a new cr. We don't deallocate
409 the old cr because other threads may still be reading it. */
412 gmx_comm("Failed to spawn threads");
417 /* END OF CAUTION: cr is now reliable */
419 /* g_membed initialisation *
420 * Because we change the mtop, init_membed is called before the init_parallel *
421 * (in case we ever want to make it run in parallel) */
422 if (opt2bSet("-membed",nfile,fnm))
424 fprintf(stderr,"Entering membed code");
426 init_membed(fplog,membed,nfile,fnm,mtop,inputrec,state,cr,&cpt_period);
431 /* now broadcast everything to the non-master nodes/threads: */
432 init_parallel(fplog, cr, inputrec, mtop);
436 pr_inputrec(fplog,0,"Input Parameters",inputrec,FALSE);
439 /* now make sure the state is initialized and propagated */
440 set_state_entries(state,inputrec,cr->nnodes);
442 /* A parallel command line option consistency check that we can
443 only do after any threads have started. */
445 (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
448 "The -dd or -npme option request a parallel simulation, "
450 "but mdrun was compiled without threads or MPI enabled"
453 "but the number of threads (option -nt) is 1"
455 "but mdrun was not started through mpirun/mpiexec or only one process was requested through mpirun/mpiexec"
461 if ((Flags & MD_RERUN) &&
462 (EI_ENERGY_MINIMIZATION(inputrec->eI) || eiNM == inputrec->eI))
464 gmx_fatal(FARGS, "The .mdp file specified an energy mininization or normal mode algorithm, and these are not compatible with mdrun -rerun");
467 if (can_use_allvsall(inputrec,mtop,TRUE,cr,fplog))
469 /* All-vs-all loops do not work with domain decomposition */
473 if (!EEL_PME(inputrec->coulombtype) || (Flags & MD_PARTDEC))
479 fcRegisterSteps(inputrec->nsteps,inputrec->init_step);
482 /* NMR restraints must be initialized before load_checkpoint,
483 * since with time averaging the history is added to t_state.
484 * For proper consistency check we therefore need to extend
486 * So the PME-only nodes (if present) will also initialize
487 * the distance restraints.
491 /* This needs to be called before read_checkpoint to extend the state */
492 init_disres(fplog,mtop,inputrec,cr,Flags & MD_PARTDEC,fcd,state);
494 if (gmx_mtop_ftype_count(mtop,F_ORIRES) > 0)
496 if (PAR(cr) && !(Flags & MD_PARTDEC))
498 gmx_fatal(FARGS,"Orientation restraints do not work (yet) with domain decomposition, use particle decomposition (mdrun option -pd)");
500 /* Orientation restraints */
503 init_orires(fplog,mtop,state->x,inputrec,cr->ms,&(fcd->orires),
508 if (DEFORM(*inputrec))
510 /* Store the deform reference box before reading the checkpoint */
513 copy_mat(state->box,box);
517 gmx_bcast(sizeof(box),box,cr);
519 /* Because we do not have the update struct available yet
520 * in which the reference values should be stored,
521 * we store them temporarily in static variables.
522 * This should be thread safe, since they are only written once
523 * and with identical values.
526 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
528 deform_init_init_step_tpx = inputrec->init_step;
529 copy_mat(box,deform_init_box_tpx);
531 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
535 if (opt2bSet("-cpi",nfile,fnm))
537 /* Check if checkpoint file exists before doing continuation.
538 * This way we can use identical input options for the first and subsequent runs...
540 if( gmx_fexist_master(opt2fn_master("-cpi",nfile,fnm,cr),cr) )
542 load_checkpoint(opt2fn_master("-cpi",nfile,fnm,cr),&fplog,
543 cr,Flags & MD_PARTDEC,ddxyz,
544 inputrec,state,&bReadRNG,&bReadEkin,
545 (Flags & MD_APPENDFILES));
549 Flags |= MD_READ_RNG;
553 Flags |= MD_READ_EKIN;
558 if (((MASTER(cr) || (Flags & MD_SEPPOT)) && (Flags & MD_APPENDFILES))
560 /* With thread MPI only the master node/thread exists in mdrun.c,
561 * therefore non-master nodes need to open the "seppot" log file here.
563 || (!MASTER(cr) && (Flags & MD_SEPPOT))
567 gmx_log_open(ftp2fn(efLOG,nfile,fnm),cr,!(Flags & MD_SEPPOT),
573 copy_mat(state->box,box);
578 gmx_bcast(sizeof(box),box,cr);
581 /* Essential dynamics */
582 if (opt2bSet("-ei",nfile,fnm))
584 /* Open input and output files, allocate space for ED data structure */
585 ed = ed_open(nfile,fnm,Flags,cr);
588 if (bVerbose && SIMMASTER(cr))
590 fprintf(stderr,"Loaded with Money\n\n");
593 if (PAR(cr) && !((Flags & MD_PARTDEC) ||
594 EI_TPI(inputrec->eI) ||
595 inputrec->eI == eiNM))
597 cr->dd = init_domain_decomposition(fplog,cr,Flags,ddxyz,rdd,rconstr,
602 &ddbox,&npme_major,&npme_minor);
604 make_dd_communicators(fplog,cr,dd_node_order);
606 /* Set overallocation to avoid frequent reallocation of arrays */
607 set_over_alloc_dd(TRUE);
611 /* PME, if used, is done on all nodes with 1D decomposition */
613 cr->duty = (DUTY_PP | DUTY_PME);
616 if (!EI_TPI(inputrec->eI))
618 npme_major = cr->nnodes;
621 if (inputrec->ePBC == epbcSCREW)
624 "pbc=%s is only implemented with domain decomposition",
625 epbc_names[inputrec->ePBC]);
631 /* After possible communicator splitting in make_dd_communicators.
632 * we can set up the intra/inter node communication.
634 gmx_setup_nodecomm(fplog,cr);
637 wcycle = wallcycle_init(fplog,resetstep,cr);
640 /* Master synchronizes its value of reset_counters with all nodes
641 * including PME only nodes */
642 reset_counters = wcycle_get_reset_counters(wcycle);
643 gmx_bcast_sim(sizeof(reset_counters),&reset_counters,cr);
644 wcycle_set_reset_counters(wcycle, reset_counters);
649 if (cr->duty & DUTY_PP)
651 /* For domain decomposition we allocate dynamically
652 * in dd_partition_system.
654 if (DOMAINDECOMP(cr))
656 bcast_state_setup(cr,state);
666 bcast_state(cr,state,TRUE);
670 /* Dihedral Restraints */
671 if (gmx_mtop_ftype_count(mtop,F_DIHRES) > 0)
673 init_dihres(fplog,mtop,inputrec,fcd);
676 /* Initiate forcerecord */
678 init_forcerec(fplog,oenv,fr,fcd,inputrec,mtop,cr,box,FALSE,
679 opt2fn("-table",nfile,fnm),
680 opt2fn("-tablep",nfile,fnm),
681 opt2fn("-tableb",nfile,fnm),FALSE,pforce);
683 /* version for PCA_NOT_READ_NODE (see md.c) */
684 /*init_forcerec(fplog,fr,fcd,inputrec,mtop,cr,box,FALSE,
685 "nofile","nofile","nofile",FALSE,pforce);
687 fr->bSepDVDL = ((Flags & MD_SEPPOT) == MD_SEPPOT);
689 /* Initialize QM-MM */
692 init_QMMMrec(cr,box,mtop,inputrec,fr);
695 /* Initialize the mdatoms structure.
696 * mdatoms is not filled with atom data,
697 * as this can not be done now with domain decomposition.
699 mdatoms = init_mdatoms(fplog,mtop,inputrec->efep!=efepNO);
701 /* Initialize the virtual site communication */
702 vsite = init_vsite(mtop,cr);
704 calc_shifts(box,fr->shift_vec);
706 /* With periodic molecules the charge groups should be whole at start up
707 * and the virtual sites should not be far from their proper positions.
709 if (!inputrec->bContinuation && MASTER(cr) &&
710 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
712 /* Make molecules whole at start of run */
713 if (fr->ePBC != epbcNONE)
715 do_pbc_first_mtop(fplog,inputrec->ePBC,box,mtop,state->x);
719 /* Correct initial vsite positions are required
720 * for the initial distribution in the domain decomposition
721 * and for the initial shell prediction.
723 construct_vsites_mtop(fplog,vsite,mtop,state->x);
727 /* Initiate PPPM if necessary */
728 if (fr->eeltype == eelPPPM)
730 if (mdatoms->nChargePerturbed)
732 gmx_fatal(FARGS,"Free energy with %s is not implemented",
733 eel_names[fr->eeltype]);
735 status = gmx_pppm_init(fplog,cr,oenv,FALSE,TRUE,box,
736 getenv("GMXGHAT"),inputrec, (Flags & MD_REPRODUCIBLE));
739 gmx_fatal(FARGS,"Error %d initializing PPPM",status);
743 if (EEL_PME(fr->eeltype))
745 ewaldcoeff = fr->ewaldcoeff;
746 pmedata = &fr->pmedata;
755 /* This is a PME only node */
757 /* We don't need the state */
760 ewaldcoeff = calc_ewaldcoeff(inputrec->rcoulomb, inputrec->ewald_rtol);
764 /* Initiate PME if necessary,
765 * either on all nodes or on dedicated PME nodes only. */
766 if (EEL_PME(inputrec->coulombtype))
770 nChargePerturbed = mdatoms->nChargePerturbed;
772 if (cr->npmenodes > 0)
774 /* The PME only nodes need to know nChargePerturbed */
775 gmx_bcast_sim(sizeof(nChargePerturbed),&nChargePerturbed,cr);
777 if (cr->duty & DUTY_PME)
779 status = gmx_pme_init(pmedata,cr,npme_major,npme_minor,inputrec,
780 mtop ? mtop->natoms : 0,nChargePerturbed,
781 (Flags & MD_REPRODUCIBLE));
784 gmx_fatal(FARGS,"Error %d initializing PME",status);
790 if (integrator[inputrec->eI].func == do_md
793 integrator[inputrec->eI].func == do_md_openmm
797 /* Turn on signal handling on all nodes */
799 * (A user signal from the PME nodes (if any)
800 * is communicated to the PP nodes.
802 signal_handler_install();
805 if (cr->duty & DUTY_PP)
807 if (inputrec->ePull != epullNO)
809 /* Initialize pull code */
810 init_pull(fplog,inputrec,nfile,fnm,mtop,cr,oenv,
811 EI_DYNAMICS(inputrec->eI) && MASTER(cr),Flags);
814 constr = init_constraints(fplog,mtop,inputrec,ed,state,cr);
816 if (DOMAINDECOMP(cr))
818 dd_init_bondeds(fplog,cr->dd,mtop,vsite,constr,inputrec,
819 Flags & MD_DDBONDCHECK,fr->cginfo_mb);
821 set_dd_parameters(fplog,cr->dd,dlb_scale,inputrec,fr,&ddbox);
823 setup_dd_grid(fplog,cr->dd);
826 /* Now do whatever the user wants us to do (how flexible...) */
827 integrator[inputrec->eI].func(fplog,cr,nfile,fnm,
828 oenv,bVerbose,bCompact,
831 nstepout,inputrec,mtop,
833 mdatoms,nrnb,wcycle,ed,fr,
834 repl_ex_nst,repl_ex_seed,
836 cpt_period,max_hours,
841 if (inputrec->ePull != epullNO)
843 finish_pull(fplog,inputrec->pull);
849 gmx_pmeonly(*pmedata,cr,nrnb,wcycle,ewaldcoeff,FALSE,inputrec);
852 if (EI_DYNAMICS(inputrec->eI) || EI_TPI(inputrec->eI))
854 /* Some timing stats */
857 if (runtime.proc == 0)
859 runtime.proc = runtime.real;
868 wallcycle_stop(wcycle,ewcRUN);
870 /* Finish up, write some stuff
871 * if rerunMD, don't write last frame again
873 finish_run(fplog,cr,ftp2fn(efSTO,nfile,fnm),
874 inputrec,nrnb,wcycle,&runtime,
875 EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
877 if (opt2bSet("-membed",nfile,fnm))
882 /* Does what it says */
883 print_date_and_time(fplog,cr->nodeid,"Finished mdrun",&runtime);
885 /* Close logfile already here if we were appending to it */
886 if (MASTER(cr) && (Flags & MD_APPENDFILES))
888 gmx_log_close(fplog);
891 rc=(int)gmx_get_stop_condition();
894 /* we need to join all threads. The sub-threads join when they
895 exit this function, but the master thread needs to be told to
897 if (PAR(cr) && MASTER(cr))
906 void md_print_warning(const t_commrec *cr,FILE *fplog,const char *buf)
910 fprintf(stderr,"\n%s\n",buf);
914 fprintf(fplog,"\n%s\n",buf);