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
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2010, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
43 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
84 #include "compute_io.h"
86 #include "checkpoint.h"
87 #include "mtop_util.h"
88 #include "sighandler.h"
98 /* include even when OpenMM not used to force compilation of do_md_openmm */
99 #include "openmm_wrapper.h"
101 double do_md_openmm(FILE *fplog,t_commrec *cr,int nfile,const t_filenm fnm[],
102 const output_env_t oenv, gmx_bool bVerbose,gmx_bool bCompact,
104 gmx_vsite_t *vsite,gmx_constr_t constr,
105 int stepout,t_inputrec *ir,
106 gmx_mtop_t *top_global,
108 t_state *state_global,
110 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
111 gmx_edsam_t ed,t_forcerec *fr,
112 int repl_ex_nst,int repl_ex_seed,
113 gmx_membed_t *membed,
114 real cpt_period,real max_hours,
115 const char *deviceOptions,
117 gmx_runtime_t *runtime)
120 gmx_large_int_t step,step_rel;
124 bFirstStep,bStateFromTPX,bLastStep,bStartingFromCpt;
125 gmx_bool bInitStep=TRUE;
126 gmx_bool do_ene,do_log, do_verbose,
128 tensor force_vir,shake_vir,total_vir,pres;
135 t_mdebin *mdebin=NULL;
140 gmx_enerdata_t *enerd;
142 gmx_global_stat_t gstat;
143 gmx_update_t upd=NULL;
147 gmx_groups_t *groups;
148 gmx_ekindata_t *ekind, *ekind_save;
152 real reset_counters=0,reset_counters_now=0;
153 char sbuf[STEPSTRSIZE],sbuf2[STEPSTRSIZE];
154 int handled_stop_condition=gmx_stop_cond_none;
156 const char *ommOptions = NULL;
160 /* Checks in cmake should prevent the compilation in double precision
161 * with OpenMM, but just to be sure we check here.
163 gmx_fatal(FARGS,"Compilation was performed in double precision, but OpenMM only supports single precision. If you want to use to OpenMM, compile in single precision.");
166 bAppend = (Flags & MD_APPENDFILES);
167 check_ir_old_tpx_versions(cr,fplog,ir,top_global);
169 groups = &top_global->groups;
172 init_md(fplog,cr,ir,oenv,&t,&t0,&state_global->lambda,&lam0,
173 nrnb,top_global,&upd,
174 nfile,fnm,&outf,&mdebin,
175 force_vir,shake_vir,mu_tot,&bSimAnn,&vcm,state_global,Flags);
177 clear_mat(total_vir);
179 /* Energy terms and groups */
181 init_enerdata(top_global->groups.grps[egcENER].nr,ir->n_flambda,enerd);
182 snew(f,top_global->natoms);
184 /* Kinetic energy data */
186 init_ekindata(fplog,top_global,&(ir->opts),ekind);
187 /* needed for iteration of constraints */
189 init_ekindata(fplog,top_global,&(ir->opts),ekind_save);
190 /* Copy the cos acceleration to the groups struct */
191 ekind->cosacc.cos_accel = ir->cos_accel;
193 gstat = global_stat_init(ir);
197 double io = compute_io(ir,top_global->natoms,groups,mdebin->ebin->nener,1);
198 if ((io > 2000) && MASTER(cr))
200 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
204 top = gmx_mtop_generate_local_top(top_global,ir);
207 a1 = top_global->natoms;
209 state = partdec_init_local_state(cr,state_global);
212 atoms2md(top_global,ir,0,NULL,a0,a1-a0,mdatoms);
216 set_vsite_top(vsite,top,mdatoms,cr);
219 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols)
221 graph = mk_graph(fplog,&(top->idef),0,top_global->natoms,FALSE,FALSE);
224 update_mdatoms(mdatoms,state->lambda);
226 if (deviceOptions[0]=='\0')
228 /* empty options, which should default to OpenMM in this build */
229 ommOptions=deviceOptions;
233 if (gmx_strncasecmp(deviceOptions,"OpenMM",6)!=0)
235 gmx_fatal(FARGS, "This Gromacs version currently only works with OpenMM. Use -device \"OpenMM:<options>\"");
239 ommOptions=strchr(deviceOptions,':');
240 if (NULL!=ommOptions)
242 /* Increase the pointer to skip the colon */
248 openmmData = openmm_init(fplog, ommOptions, ir, top_global, top, mdatoms, fr, state);
249 please_cite(fplog,"Friedrichs2009");
253 /* Update mdebin with energy history if appending to output files */
254 if ( Flags & MD_APPENDFILES )
256 restore_energyhistory_from_state(mdebin,&state_global->enerhist);
258 /* Set the initial energy history in state to zero by updating once */
259 update_energyhistory(&state_global->enerhist,mdebin);
264 set_constraints(constr,top,ir,mdatoms,cr);
267 if (!ir->bContinuation)
269 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
271 /* Set the velocities of frozen particles to zero */
272 for (i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++)
274 for (m=0; m<DIM; m++)
276 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
286 /* Constrain the initial coordinates and velocities */
287 do_constrain_first(fplog,constr,ir,mdatoms,state,f,
288 graph,cr,nrnb,fr,top,shake_vir);
292 /* Construct the virtual sites for the initial configuration */
293 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,NULL,
294 top->idef.iparams,top->idef.il,
295 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
304 fprintf(fplog,"Initial temperature: %g K\n",enerd->term[F_TEMP]);
305 fprintf(stderr,"starting mdrun '%s'\n",
306 *(top_global->name));
309 sprintf(tbuf,"%8.1f",(ir->init_step+ir->nsteps)*ir->delta_t);
313 sprintf(tbuf,"%s","infinite");
315 if (ir->init_step > 0)
317 fprintf(stderr,"%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
318 gmx_step_str(ir->init_step+ir->nsteps,sbuf),tbuf,
319 gmx_step_str(ir->init_step,sbuf2),
320 ir->init_step*ir->delta_t);
324 fprintf(stderr,"%s steps, %s ps.\n",
325 gmx_step_str(ir->nsteps,sbuf),tbuf);
331 /* Set and write start time */
332 runtime_start(runtime);
333 print_date_and_time(fplog,cr->nodeid,"Started mdrun",runtime);
334 wallcycle_start(wcycle,ewcRUN);
338 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
341 /***********************************************************
345 ************************************************************/
347 /* loop over MD steps or if rerunMD to end of input trajectory */
349 /* Skip the first Nose-Hoover integration when we get the state from tpx */
350 bStateFromTPX = !opt2bSet("-cpi",nfile,fnm);
351 bInitStep = bFirstStep && bStateFromTPX;
352 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
355 init_global_signals(&gs,cr,ir,repl_ex_nst);
357 step = ir->init_step;
362 wallcycle_start(wcycle,ewcSTEP);
364 bLastStep = (step_rel == ir->nsteps);
365 t = t0 + step*ir->delta_t;
367 if (gs.set[eglsSTOPCOND] != 0)
372 do_log = do_per_step(step,ir->nstlog) || bFirstStep || bLastStep;
373 do_verbose = bVerbose &&
374 (step % stepout == 0 || bFirstStep || bLastStep);
376 if (MASTER(cr) && do_log)
378 print_ebin_header(fplog,step,t,state->lambda);
381 clear_mat(force_vir);
383 /* We write a checkpoint at this MD step when:
384 * either when we signalled through gs (in OpenMM NS works different),
385 * or at the last step (but not when we do not want confout),
386 * but never at the first step.
388 bCPT = ((gs.set[eglsCHKPT] ||
389 (bLastStep && (Flags & MD_CONFOUT))) &&
390 step > ir->init_step );
393 gs.set[eglsCHKPT] = 0;
396 /* Now we have the energies and forces corresponding to the
397 * coordinates at time t. We must output all of this before
399 * for RerunMD t is read from input trajectory
402 if (do_per_step(step,ir->nstxout))
404 mdof_flags |= MDOF_X;
406 if (do_per_step(step,ir->nstvout))
408 mdof_flags |= MDOF_V;
410 if (do_per_step(step,ir->nstfout))
412 mdof_flags |= MDOF_F;
414 if (do_per_step(step,ir->nstxtcout))
416 mdof_flags |= MDOF_XTC;
420 mdof_flags |= MDOF_CPT;
422 do_ene = (do_per_step(step,ir->nstenergy) || bLastStep);
424 if (mdof_flags != 0 || do_ene || do_log)
426 wallcycle_start(wcycle,ewcTRAJ);
427 bF = (mdof_flags & MDOF_F);
428 bX = (mdof_flags & (MDOF_X | MDOF_XTC | MDOF_CPT));
429 bV = (mdof_flags & (MDOF_V | MDOF_CPT));
431 openmm_copy_state(openmmData, state, &t, f, enerd, bX, bV, bF, do_ene);
433 upd_mdebin(mdebin,FALSE,TRUE,
434 t,mdatoms->tmass,enerd,state,lastbox,
435 shake_vir,force_vir,total_vir,pres,
436 ekind,mu_tot,constr);
437 print_ebin(outf->fp_ene,do_ene,FALSE,FALSE,do_log?fplog:NULL,
439 eprNORMAL,bCompact,mdebin,fcd,groups,&(ir->opts));
440 write_traj(fplog,cr,outf,mdof_flags,top_global,
441 step,t,state,state_global,f,f_global,&n_xtc,&x_xtc);
448 if (bLastStep && step_rel == ir->nsteps &&
449 (Flags & MD_CONFOUT) && MASTER(cr))
451 /* x and v have been collected in write_traj,
452 * because a checkpoint file will always be written
455 fprintf(stderr,"\nWriting final coordinates.\n");
456 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols)
458 /* Make molecules whole only for confout writing */
459 do_pbc_mtop(fplog,ir->ePBC,state->box,top_global,state_global->x);
461 write_sto_conf_mtop(ftp2fn(efSTO,nfile,fnm),
462 *top_global->name,top_global,
463 state_global->x,state_global->v,
464 ir->ePBC,state->box);
467 wallcycle_stop(wcycle,ewcTRAJ);
470 /* Determine the wallclock run time up till now */
471 run_time = gmx_gettime() - (double)runtime->real;
473 /* Check whether everything is still allright */
474 if (((int)gmx_get_stop_condition() > handled_stop_condition)
480 /* this is just make gs.sig compatible with the hack
481 of sending signals around by MPI_Reduce with together with
483 /* NOTE: this only works for serial code. For code that allows
484 MPI nodes to propagate their condition, see kernel/md.c*/
485 if ( gmx_get_stop_condition() == gmx_stop_cond_next_ns )
486 gs.set[eglsSTOPCOND]=1;
487 if ( gmx_get_stop_condition() == gmx_stop_cond_next )
488 gs.set[eglsSTOPCOND]=1;
489 /* < 0 means stop at next step, > 0 means stop at next NS step */
493 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
494 gmx_get_signal_name(),
495 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
499 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
500 gmx_get_signal_name(),
501 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
503 handled_stop_condition=(int)gmx_get_stop_condition();
505 else if (MASTER(cr) &&
506 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
507 gs.set[eglsSTOPCOND] == 0)
509 /* Signal to terminate the run */
510 gs.set[eglsSTOPCOND] = 1;
513 fprintf(fplog,"\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
515 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
519 if (MASTER(cr) && (cpt_period >= 0 &&
521 run_time >= nchkpt*cpt_period*60.0)) &&
522 gs.set[eglsCHKPT] == 0)
524 gs.set[eglsCHKPT] = 1;
527 /* Time for performance */
528 if (((step % stepout) == 0) || bLastStep)
530 runtime_upd_proc(runtime);
533 if (do_per_step(step,ir->nstlog))
535 if (fflush(fplog) != 0)
537 gmx_fatal(FARGS,"Cannot flush logfile - maybe you are out of quota?");
541 /* Remaining runtime */
542 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal() ))
544 print_time(stderr,runtime,step,ir,cr);
549 bStartingFromCpt = FALSE;
553 openmm_take_one_step(openmmData);
555 /* End of main MD loop */
559 runtime_end(runtime);
563 if (ir->nstcalcenergy > 0)
565 print_ebin(outf->fp_ene,FALSE,FALSE,FALSE,fplog,step,t,
566 eprAVER,FALSE,mdebin,fcd,groups,&(ir->opts));
570 openmm_cleanup(fplog, openmmData);
576 runtime->nsteps_done = step_rel;