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"
76 #include "md_openmm.h"
90 #include "md_openmm.h"
95 gmx_integrator_t *func;
98 /* The array should match the eI array in include/types/enums.h */
99 #ifdef GMX_OPENMM /* FIXME do_md_openmm needs fixing */
100 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}};
102 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}};
105 gmx_large_int_t deform_init_init_step_tpx;
106 matrix deform_init_box_tpx;
108 tMPI_Thread_mutex_t deform_init_box_mutex=TMPI_THREAD_MUTEX_INITIALIZER;
113 struct mdrunner_arglist
127 const char *dddlb_opt;
140 const char *deviceOptions;
142 int ret; /* return value */
146 /* The function used for spawning threads. Extracts the mdrunner()
147 arguments from its one argument and calls mdrunner(), after making
149 static void mdrunner_start_fn(void *arg)
151 struct mdrunner_arglist *mda=(struct mdrunner_arglist*)arg;
152 struct mdrunner_arglist mc=*mda; /* copy the arg list to make sure
153 that it's thread-local. This doesn't
154 copy pointed-to items, of course,
155 but those are all const. */
156 t_commrec *cr; /* we need a local version of this */
160 fnm = dup_tfn(mc.nfile, mc.fnm);
162 cr = init_par_threads(mc.cr);
169 mda->ret=mdrunner(cr->nnodes, fplog, cr, mc.nfile, fnm, mc.oenv,
170 mc.bVerbose, mc.bCompact, mc.nstglobalcomm,
171 mc.ddxyz, mc.dd_node_order, mc.rdd,
172 mc.rconstr, mc.dddlb_opt, mc.dlb_scale,
173 mc.ddcsx, mc.ddcsy, mc.ddcsz, mc.nstepout, mc.resetstep,
174 mc.nmultisim, mc.repl_ex_nst, mc.repl_ex_seed, mc.pforce,
175 mc.cpt_period, mc.max_hours, mc.deviceOptions, mc.Flags);
178 /* called by mdrunner() to start a specific number of threads (including
179 the main thread) for thread-parallel runs. This in turn calls mdrunner()
181 All options besides nthreads are the same as for mdrunner(). */
182 static t_commrec *mdrunner_start_threads(int nthreads,
183 FILE *fplog,t_commrec *cr,int nfile,
184 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
185 gmx_bool bCompact, int nstglobalcomm,
186 ivec ddxyz,int dd_node_order,real rdd,real rconstr,
187 const char *dddlb_opt,real dlb_scale,
188 const char *ddcsx,const char *ddcsy,const char *ddcsz,
189 int nstepout,int resetstep,int nmultisim,int repl_ex_nst,
190 int repl_ex_seed, real pforce,real cpt_period, real max_hours,
191 const char *deviceOptions, unsigned long Flags)
194 struct mdrunner_arglist *mda;
195 t_commrec *crn; /* the new commrec */
198 /* first check whether we even need to start tMPI */
202 /* a few small, one-time, almost unavoidable memory leaks: */
204 fnmn=dup_tfn(nfile, fnm);
206 /* fill the data structure to pass as void pointer to thread start fn */
212 mda->bVerbose=bVerbose;
213 mda->bCompact=bCompact;
214 mda->nstglobalcomm=nstglobalcomm;
215 mda->ddxyz[XX]=ddxyz[XX];
216 mda->ddxyz[YY]=ddxyz[YY];
217 mda->ddxyz[ZZ]=ddxyz[ZZ];
218 mda->dd_node_order=dd_node_order;
220 mda->rconstr=rconstr;
221 mda->dddlb_opt=dddlb_opt;
222 mda->dlb_scale=dlb_scale;
226 mda->nstepout=nstepout;
227 mda->resetstep=resetstep;
228 mda->nmultisim=nmultisim;
229 mda->repl_ex_nst=repl_ex_nst;
230 mda->repl_ex_seed=repl_ex_seed;
232 mda->cpt_period=cpt_period;
233 mda->max_hours=max_hours;
234 mda->deviceOptions=deviceOptions;
237 fprintf(stderr, "Starting %d threads\n",nthreads);
239 /* now spawn new threads that start mdrunner_start_fn(), while
240 the main thread returns */
241 ret=tMPI_Init_fn(TRUE, nthreads, mdrunner_start_fn, (void*)(mda) );
242 if (ret!=TMPI_SUCCESS)
245 /* make a new comm_rec to reflect the new situation */
246 crn=init_par_threads(cr);
251 /* get the number of threads based on how many there were requested,
252 which algorithms we're using, and how many particles there are. */
253 static int get_nthreads(int nthreads_requested, t_inputrec *inputrec,
256 int nthreads,nthreads_new;
257 int min_atoms_per_thread;
259 nthreads = nthreads_requested;
261 /* determine # of hardware threads. */
262 if (nthreads_requested < 1)
264 nthreads = tMPI_Get_recommended_nthreads();
267 if (inputrec->eI == eiNM || EI_TPI(inputrec->eI))
269 /* Steps are divided over the nodes iso splitting the atoms */
270 min_atoms_per_thread = 0;
274 min_atoms_per_thread = MIN_ATOMS_PER_THREAD;
277 /* Check if an algorithm does not support parallel simulation. */
279 ( inputrec->eI == eiLBFGS ||
280 inputrec->coulombtype == eelEWALD ) )
282 fprintf(stderr,"\nThe integration or electrostatics algorithm doesn't support parallel runs. Not starting any threads.\n");
285 else if (nthreads_requested < 1 &&
286 mtop->natoms/nthreads < min_atoms_per_thread)
288 /* the thread number was chosen automatically, but there are too many
289 threads (too few atoms per thread) */
290 nthreads_new = max(1,mtop->natoms/min_atoms_per_thread);
292 if (nthreads_new > 8 || (nthreads == 8 && nthreads_new > 4))
294 /* Use only multiples of 4 above 8 threads
295 * or with an 8-core processor
296 * (to avoid 6 threads on 8 core processors with 4 real cores).
298 nthreads_new = (nthreads_new/4)*4;
300 else if (nthreads_new > 4)
302 /* Avoid 5 or 7 threads */
303 nthreads_new = (nthreads_new/2)*2;
306 nthreads = nthreads_new;
308 fprintf(stderr,"\n");
309 fprintf(stderr,"NOTE: Parallelization is limited by the small number of atoms,\n");
310 fprintf(stderr," only starting %d threads.\n",nthreads);
311 fprintf(stderr," You can use the -nt option to optimize the number of threads.\n\n");
318 int mdrunner(int nthreads_requested, FILE *fplog,t_commrec *cr,int nfile,
319 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
320 gmx_bool bCompact, int nstglobalcomm,
321 ivec ddxyz,int dd_node_order,real rdd,real rconstr,
322 const char *dddlb_opt,real dlb_scale,
323 const char *ddcsx,const char *ddcsy,const char *ddcsz,
324 int nstepout,int resetstep,int nmultisim,int repl_ex_nst,
325 int repl_ex_seed, real pforce,real cpt_period,real max_hours,
326 const char *deviceOptions, unsigned long Flags)
328 double nodetime=0,realtime;
329 t_inputrec *inputrec;
332 gmx_ddbox_t ddbox={0};
333 int npme_major,npme_minor;
336 gmx_mtop_t *mtop=NULL;
337 t_mdatoms *mdatoms=NULL;
341 gmx_pme_t *pmedata=NULL;
342 gmx_vsite_t *vsite=NULL;
344 int i,m,nChargePerturbed=-1,status,nalloc;
346 gmx_wallcycle_t wcycle;
347 gmx_bool bReadRNG,bReadEkin;
349 gmx_runtime_t runtime;
351 gmx_large_int_t reset_counters;
353 t_commrec *cr_old=cr;
356 /* CAUTION: threads may be started later on in this function, so
357 cr doesn't reflect the final parallel state right now */
361 if (bVerbose && SIMMASTER(cr))
363 fprintf(stderr,"Getting Loaded...\n");
366 if (Flags & MD_APPENDFILES)
374 /* Read (nearly) all data required for the simulation */
375 read_tpx_state(ftp2fn(efTPX,nfile,fnm),inputrec,state,NULL,mtop);
377 /* NOW the threads will be started: */
379 nthreads = get_nthreads(nthreads_requested, inputrec, mtop);
383 /* now start the threads. */
384 cr=mdrunner_start_threads(nthreads, fplog, cr_old, nfile, fnm,
385 oenv, bVerbose, bCompact, nstglobalcomm,
386 ddxyz, dd_node_order, rdd, rconstr,
387 dddlb_opt, dlb_scale, ddcsx, ddcsy, ddcsz,
388 nstepout, resetstep, nmultisim,
389 repl_ex_nst, repl_ex_seed, pforce,
390 cpt_period, max_hours, deviceOptions,
392 /* the main thread continues here with a new cr. We don't deallocate
393 the old cr because other threads may still be reading it. */
396 gmx_comm("Failed to spawn threads");
401 /* END OF CAUTION: cr is now reliable */
405 /* now broadcast everything to the non-master nodes/threads: */
406 init_parallel(fplog, cr, inputrec, mtop);
410 pr_inputrec(fplog,0,"Input Parameters",inputrec,FALSE);
413 /* now make sure the state is initialized and propagated */
414 set_state_entries(state,inputrec,cr->nnodes);
416 /* A parallel command line option consistency check that we can
417 only do after any threads have started. */
419 (ddxyz[XX] > 1 || ddxyz[YY] > 1 || ddxyz[ZZ] > 1 || cr->npmenodes > 0))
422 "The -dd or -npme option request a parallel simulation, "
424 "but mdrun was compiled without threads or MPI enabled"
427 "but the number of threads (option -nt) is 1"
429 "but mdrun was not started through mpirun/mpiexec or only one process was requested through mpirun/mpiexec"
435 if (can_use_allvsall(inputrec,mtop,TRUE,cr,fplog))
437 /* All-vs-all loops do not work with domain decomposition */
441 if (!EEL_PME(inputrec->coulombtype) || (Flags & MD_PARTDEC))
447 fcRegisterSteps(inputrec->nsteps,inputrec->init_step);
450 /* NMR restraints must be initialized before load_checkpoint,
451 * since with time averaging the history is added to t_state.
452 * For proper consistency check we therefore need to extend
454 * So the PME-only nodes (if present) will also initialize
455 * the distance restraints.
459 /* This needs to be called before read_checkpoint to extend the state */
460 init_disres(fplog,mtop,inputrec,cr,Flags & MD_PARTDEC,fcd,state);
462 if (gmx_mtop_ftype_count(mtop,F_ORIRES) > 0)
464 if (PAR(cr) && !(Flags & MD_PARTDEC))
466 gmx_fatal(FARGS,"Orientation restraints do not work (yet) with domain decomposition, use particle decomposition (mdrun option -pd)");
468 /* Orientation restraints */
471 init_orires(fplog,mtop,state->x,inputrec,cr->ms,&(fcd->orires),
476 if (DEFORM(*inputrec))
478 /* Store the deform reference box before reading the checkpoint */
481 copy_mat(state->box,box);
485 gmx_bcast(sizeof(box),box,cr);
487 /* Because we do not have the update struct available yet
488 * in which the reference values should be stored,
489 * we store them temporarily in static variables.
490 * This should be thread safe, since they are only written once
491 * and with identical values.
494 tMPI_Thread_mutex_lock(&deform_init_box_mutex);
496 deform_init_init_step_tpx = inputrec->init_step;
497 copy_mat(box,deform_init_box_tpx);
499 tMPI_Thread_mutex_unlock(&deform_init_box_mutex);
503 if (opt2bSet("-cpi",nfile,fnm))
505 /* Check if checkpoint file exists before doing continuation.
506 * This way we can use identical input options for the first and subsequent runs...
508 if( gmx_fexist_master(opt2fn_master("-cpi",nfile,fnm,cr),cr) )
510 load_checkpoint(opt2fn_master("-cpi",nfile,fnm,cr),&fplog,
511 cr,Flags & MD_PARTDEC,ddxyz,
512 inputrec,state,&bReadRNG,&bReadEkin,
513 (Flags & MD_APPENDFILES));
517 Flags |= MD_READ_RNG;
521 Flags |= MD_READ_EKIN;
526 if (((MASTER(cr) || (Flags & MD_SEPPOT)) && (Flags & MD_APPENDFILES))
528 /* With thread MPI only the master node/thread exists in mdrun.c,
529 * therefore non-master nodes need to open the "seppot" log file here.
531 || (!MASTER(cr) && (Flags & MD_SEPPOT))
535 gmx_log_open(ftp2fn(efLOG,nfile,fnm),cr,!(Flags & MD_SEPPOT),
541 copy_mat(state->box,box);
546 gmx_bcast(sizeof(box),box,cr);
549 /* Essential dynamics */
550 if (opt2bSet("-ei",nfile,fnm))
552 /* Open input and output files, allocate space for ED data structure */
553 ed = ed_open(nfile,fnm,Flags,cr);
556 if (bVerbose && SIMMASTER(cr))
558 fprintf(stderr,"Loaded with Money\n\n");
561 if (PAR(cr) && !((Flags & MD_PARTDEC) || EI_TPI(inputrec->eI)))
563 cr->dd = init_domain_decomposition(fplog,cr,Flags,ddxyz,rdd,rconstr,
568 &ddbox,&npme_major,&npme_minor);
570 make_dd_communicators(fplog,cr,dd_node_order);
572 /* Set overallocation to avoid frequent reallocation of arrays */
573 set_over_alloc_dd(TRUE);
577 /* PME, if used, is done on all nodes with 1D decomposition */
579 cr->duty = (DUTY_PP | DUTY_PME);
580 npme_major = cr->nnodes;
583 if (inputrec->ePBC == epbcSCREW)
586 "pbc=%s is only implemented with domain decomposition",
587 epbc_names[inputrec->ePBC]);
593 /* After possible communicator splitting in make_dd_communicators.
594 * we can set up the intra/inter node communication.
596 gmx_setup_nodecomm(fplog,cr);
599 wcycle = wallcycle_init(fplog,resetstep,cr);
602 /* Master synchronizes its value of reset_counters with all nodes
603 * including PME only nodes */
604 reset_counters = wcycle_get_reset_counters(wcycle);
605 gmx_bcast_sim(sizeof(reset_counters),&reset_counters,cr);
606 wcycle_set_reset_counters(wcycle, reset_counters);
611 if (cr->duty & DUTY_PP)
613 /* For domain decomposition we allocate dynamically
614 * in dd_partition_system.
616 if (DOMAINDECOMP(cr))
618 bcast_state_setup(cr,state);
628 bcast_state(cr,state,TRUE);
632 /* Dihedral Restraints */
633 if (gmx_mtop_ftype_count(mtop,F_DIHRES) > 0)
635 init_dihres(fplog,mtop,inputrec,fcd);
638 /* Initiate forcerecord */
640 init_forcerec(fplog,oenv,fr,fcd,inputrec,mtop,cr,box,FALSE,
641 opt2fn("-table",nfile,fnm),
642 opt2fn("-tablep",nfile,fnm),
643 opt2fn("-tableb",nfile,fnm),FALSE,pforce);
645 /* version for PCA_NOT_READ_NODE (see md.c) */
646 /*init_forcerec(fplog,fr,fcd,inputrec,mtop,cr,box,FALSE,
647 "nofile","nofile","nofile",FALSE,pforce);
649 fr->bSepDVDL = ((Flags & MD_SEPPOT) == MD_SEPPOT);
651 /* Initialize QM-MM */
654 init_QMMMrec(cr,box,mtop,inputrec,fr);
657 /* Initialize the mdatoms structure.
658 * mdatoms is not filled with atom data,
659 * as this can not be done now with domain decomposition.
661 mdatoms = init_mdatoms(fplog,mtop,inputrec->efep!=efepNO);
663 /* Initialize the virtual site communication */
664 vsite = init_vsite(mtop,cr);
666 calc_shifts(box,fr->shift_vec);
668 /* With periodic molecules the charge groups should be whole at start up
669 * and the virtual sites should not be far from their proper positions.
671 if (!inputrec->bContinuation && MASTER(cr) &&
672 !(inputrec->ePBC != epbcNONE && inputrec->bPeriodicMols))
674 /* Make molecules whole at start of run */
675 if (fr->ePBC != epbcNONE)
677 do_pbc_first_mtop(fplog,inputrec->ePBC,box,mtop,state->x);
681 /* Correct initial vsite positions are required
682 * for the initial distribution in the domain decomposition
683 * and for the initial shell prediction.
685 construct_vsites_mtop(fplog,vsite,mtop,state->x);
689 /* Initiate PPPM if necessary */
690 if (fr->eeltype == eelPPPM)
692 if (mdatoms->nChargePerturbed)
694 gmx_fatal(FARGS,"Free energy with %s is not implemented",
695 eel_names[fr->eeltype]);
697 status = gmx_pppm_init(fplog,cr,oenv,FALSE,TRUE,box,
698 getenv("GMXGHAT"),inputrec, (Flags & MD_REPRODUCIBLE));
701 gmx_fatal(FARGS,"Error %d initializing PPPM",status);
705 if (EEL_PME(fr->eeltype))
707 ewaldcoeff = fr->ewaldcoeff;
708 pmedata = &fr->pmedata;
717 /* This is a PME only node */
719 /* We don't need the state */
722 ewaldcoeff = calc_ewaldcoeff(inputrec->rcoulomb, inputrec->ewald_rtol);
726 /* Initiate PME if necessary,
727 * either on all nodes or on dedicated PME nodes only. */
728 if (EEL_PME(inputrec->coulombtype))
732 nChargePerturbed = mdatoms->nChargePerturbed;
734 if (cr->npmenodes > 0)
736 /* The PME only nodes need to know nChargePerturbed */
737 gmx_bcast_sim(sizeof(nChargePerturbed),&nChargePerturbed,cr);
739 if (cr->duty & DUTY_PME)
741 status = gmx_pme_init(pmedata,cr,npme_major,npme_minor,inputrec,
742 mtop ? mtop->natoms : 0,nChargePerturbed,
743 (Flags & MD_REPRODUCIBLE));
746 gmx_fatal(FARGS,"Error %d initializing PME",status);
752 if (integrator[inputrec->eI].func == do_md
755 integrator[inputrec->eI].func == do_md_openmm
759 /* Turn on signal handling on all nodes */
761 * (A user signal from the PME nodes (if any)
762 * is communicated to the PP nodes.
764 signal_handler_install();
767 if (cr->duty & DUTY_PP)
769 if (inputrec->ePull != epullNO)
771 /* Initialize pull code */
772 init_pull(fplog,inputrec,nfile,fnm,mtop,cr,oenv,
773 EI_DYNAMICS(inputrec->eI) && MASTER(cr),Flags);
776 constr = init_constraints(fplog,mtop,inputrec,ed,state,cr);
778 if (DOMAINDECOMP(cr))
780 dd_init_bondeds(fplog,cr->dd,mtop,vsite,constr,inputrec,
781 Flags & MD_DDBONDCHECK,fr->cginfo_mb);
783 set_dd_parameters(fplog,cr->dd,dlb_scale,inputrec,fr,&ddbox);
785 setup_dd_grid(fplog,cr->dd);
788 /* Now do whatever the user wants us to do (how flexible...) */
789 integrator[inputrec->eI].func(fplog,cr,nfile,fnm,
790 oenv,bVerbose,bCompact,
793 nstepout,inputrec,mtop,
795 mdatoms,nrnb,wcycle,ed,fr,
796 repl_ex_nst,repl_ex_seed,
797 cpt_period,max_hours,
802 if (inputrec->ePull != epullNO)
804 finish_pull(fplog,inputrec->pull);
810 gmx_pmeonly(*pmedata,cr,nrnb,wcycle,ewaldcoeff,FALSE,inputrec);
813 if (EI_DYNAMICS(inputrec->eI) || EI_TPI(inputrec->eI))
815 /* Some timing stats */
818 if (runtime.proc == 0)
820 runtime.proc = runtime.real;
829 wallcycle_stop(wcycle,ewcRUN);
831 /* Finish up, write some stuff
832 * if rerunMD, don't write last frame again
834 finish_run(fplog,cr,ftp2fn(efSTO,nfile,fnm),
835 inputrec,nrnb,wcycle,&runtime,
836 EI_DYNAMICS(inputrec->eI) && !MULTISIM(cr));
838 /* Does what it says */
839 print_date_and_time(fplog,cr->nodeid,"Finished mdrun",&runtime);
841 /* Close logfile already here if we were appending to it */
842 if (MASTER(cr) && (Flags & MD_APPENDFILES))
844 gmx_log_close(fplog);
847 rc=(int)gmx_get_stop_condition();
850 /* we need to join all threads. The sub-threads join when they
851 exit this function, but the master thread needs to be told to
853 if (PAR(cr) && MASTER(cr))
862 void md_print_warning(const t_commrec *cr,FILE *fplog,const char *buf)
866 fprintf(stderr,"\n%s\n",buf);
870 fprintf(fplog,"\n%s\n",buf);