OpenMM: added combination rule check, disabled restraint check for now as it's buggy
[gromacs/qmmm-gamess-us.git] / src / kernel / md_openmm.c
blob68162c4de502dbdde4b6de3bf11dd8784e080b33
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
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
32 * And Hey:
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include <signal.h>
41 #include <stdlib.h>
43 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
44 /* _isnan() */
45 #include <float.h>
46 #endif
48 #include "typedefs.h"
49 #include "smalloc.h"
50 #include "sysstuff.h"
51 #include "vec.h"
52 #include "statutil.h"
53 #include "vcm.h"
54 #include "mdebin.h"
55 #include "nrnb.h"
56 #include "calcmu.h"
57 #include "index.h"
58 #include "vsite.h"
59 #include "update.h"
60 #include "ns.h"
61 #include "trnio.h"
62 #include "xtcio.h"
63 #include "mdrun.h"
64 #include "confio.h"
65 #include "network.h"
66 #include "pull.h"
67 #include "xvgr.h"
68 #include "physics.h"
69 #include "names.h"
70 #include "xmdrun.h"
71 #include "ionize.h"
72 #include "disre.h"
73 #include "orires.h"
74 #include "dihre.h"
75 #include "pppm.h"
76 #include "pme.h"
77 #include "mdatoms.h"
78 #include "qmmm.h"
79 #include "mpelogging.h"
80 #include "domdec.h"
81 #include "partdec.h"
82 #include "topsort.h"
83 #include "coulomb.h"
84 #include "constr.h"
85 #include "compute_io.h"
86 #include "mvdata.h"
87 #include "checkpoint.h"
88 #include "mtop_util.h"
89 #include "sighandler.h"
90 #include "genborn.h"
91 #include "string2.h"
93 #ifdef GMX_THREADS
94 #include "tmpi.h"
95 #endif
97 /* include even when OpenMM not used to force compilation of do_md_openmm */
98 #include "openmm_wrapper.h"
100 /* simulation conditions to transmit */
101 enum { eglsNABNSB, eglsCHKPT, eglsSTOPCOND, eglsRESETCOUNTERS, eglsNR };
103 typedef struct
105 int nstms; /* The frequency for intersimulation communication */
106 int sig[eglsNR]; /* The signal set by one process in do_md */
107 int set[eglsNR]; /* The communicated signal, equal for all processes */
108 } globsig_t;
111 static int multisim_min(const gmx_multisim_t *ms,int nmin,int n)
113 int *buf;
114 bool bPos,bEqual;
115 int s,d;
117 snew(buf,ms->nsim);
118 buf[ms->sim] = n;
119 gmx_sumi_sim(ms->nsim,buf,ms);
120 bPos = TRUE;
121 bEqual = TRUE;
122 for (s=0; s<ms->nsim; s++)
124 bPos = bPos && (buf[s] > 0);
125 bEqual = bEqual && (buf[s] == buf[0]);
127 if (bPos)
129 if (bEqual)
131 nmin = min(nmin,buf[0]);
133 else
135 /* Find the least common multiple */
136 for (d=2; d<nmin; d++)
138 s = 0;
139 while (s < ms->nsim && d % buf[s] == 0)
141 s++;
143 if (s == ms->nsim)
145 /* We found the LCM and it is less than nmin */
146 nmin = d;
147 break;
152 sfree(buf);
154 return nmin;
157 static int multisim_nstsimsync(const t_commrec *cr,
158 const t_inputrec *ir,int repl_ex_nst)
160 int nmin;
162 if (MASTER(cr))
164 nmin = INT_MAX;
165 nmin = multisim_min(cr->ms,nmin,ir->nstlist);
166 nmin = multisim_min(cr->ms,nmin,ir->nstcalcenergy);
167 nmin = multisim_min(cr->ms,nmin,repl_ex_nst);
168 if (nmin == INT_MAX)
170 gmx_fatal(FARGS,"Can not find an appropriate interval for inter-simulation communication, since nstlist, nstcalcenergy and -replex are all <= 0");
172 /* Avoid inter-simulation communication at every (second) step */
173 if (nmin <= 2)
175 nmin = 10;
179 gmx_bcast(sizeof(int),&nmin,cr);
181 return nmin;
184 static void init_global_signals(globsig_t *gs,const t_commrec *cr,
185 const t_inputrec *ir,int repl_ex_nst)
187 int i;
189 gs->nstms = 1;
191 for (i=0; i<eglsNR; i++)
193 gs->sig[i] = 0;
194 gs->set[i] = 0;
199 double do_md_openmm(FILE *fplog,t_commrec *cr,int nfile,const t_filenm fnm[],
200 const output_env_t oenv, bool bVerbose,bool bCompact,
201 int nstglobalcomm,
202 gmx_vsite_t *vsite,gmx_constr_t constr,
203 int stepout,t_inputrec *ir,
204 gmx_mtop_t *top_global,
205 t_fcdata *fcd,
206 t_state *state_global,
207 t_mdatoms *mdatoms,
208 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
209 gmx_edsam_t ed,t_forcerec *fr,
210 int repl_ex_nst,int repl_ex_seed,
211 real cpt_period,real max_hours,
212 const char *deviceOptions,
213 unsigned long Flags,
214 gmx_runtime_t *runtime)
216 gmx_mdoutf_t *outf;
217 gmx_large_int_t step,step_rel;
218 double run_time;
219 double t,t0,lam0;
220 bool bSimAnn,
221 bFirstStep,bStateFromTPX,bLastStep,bStartingFromCpt;
222 bool bInitStep=TRUE;
223 bool bNEMD,do_ene,do_log, do_verbose,
224 bX,bV,bF,bCPT;
225 tensor force_vir,shake_vir,total_vir,pres;
226 int i,m;
227 int mdof_flags;
228 rvec mu_tot;
229 t_vcm *vcm;
230 int nchkpt=1;
231 gmx_localtop_t *top;
232 t_mdebin *mdebin=NULL;
233 t_state *state=NULL;
234 rvec *f_global=NULL;
235 int n_xtc=-1;
236 rvec *x_xtc=NULL;
237 gmx_enerdata_t *enerd;
238 rvec *f=NULL;
239 gmx_global_stat_t gstat;
240 gmx_update_t upd=NULL;
241 t_graph *graph=NULL;
242 globsig_t gs;
244 gmx_groups_t *groups;
245 gmx_ekindata_t *ekind, *ekind_save;
246 bool bAppend;
247 int a0,a1;
248 matrix lastbox;
249 real reset_counters=0,reset_counters_now=0;
250 char sbuf[STEPSTRSIZE],sbuf2[STEPSTRSIZE];
251 int handled_stop_condition=gmx_stop_cond_none;
253 const char *ommOptions = NULL;
254 void *openmmData;
256 bAppend = (Flags & MD_APPENDFILES);
257 check_ir_old_tpx_versions(cr,fplog,ir,top_global);
259 groups = &top_global->groups;
261 /* Initial values */
262 init_md(fplog,cr,ir,oenv,&t,&t0,&state_global->lambda,&lam0,
263 nrnb,top_global,&upd,
264 nfile,fnm,&outf,&mdebin,
265 force_vir,shake_vir,mu_tot,&bNEMD,&bSimAnn,&vcm,state_global,Flags);
267 clear_mat(total_vir);
268 clear_mat(pres);
269 /* Energy terms and groups */
270 snew(enerd,1);
271 init_enerdata(top_global->groups.grps[egcENER].nr,ir->n_flambda,enerd);
272 snew(f,top_global->natoms);
274 /* Kinetic energy data */
275 snew(ekind,1);
276 init_ekindata(fplog,top_global,&(ir->opts),ekind);
277 /* needed for iteration of constraints */
278 snew(ekind_save,1);
279 init_ekindata(fplog,top_global,&(ir->opts),ekind_save);
280 /* Copy the cos acceleration to the groups struct */
281 ekind->cosacc.cos_accel = ir->cos_accel;
283 gstat = global_stat_init(ir);
284 debug_gmx();
287 double io = compute_io(ir,top_global->natoms,groups,mdebin->ebin->nener,1);
288 if ((io > 2000) && MASTER(cr))
289 fprintf(stderr,
290 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
291 io);
294 top = gmx_mtop_generate_local_top(top_global,ir);
296 a0 = 0;
297 a1 = top_global->natoms;
299 state = partdec_init_local_state(cr,state_global);
300 f_global = f;
302 atoms2md(top_global,ir,0,NULL,a0,a1-a0,mdatoms);
304 if (vsite)
306 set_vsite_top(vsite,top,mdatoms,cr);
309 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols)
311 graph = mk_graph(fplog,&(top->idef),0,top_global->natoms,FALSE,FALSE);
314 update_mdatoms(mdatoms,state->lambda);
316 if (deviceOptions[0]=='\0')
318 /* empty options, which should default to OpenMM in this build */
319 ommOptions=deviceOptions;
321 else
323 if (gmx_strncasecmp(deviceOptions,"OpenMM",6)!=0)
325 gmx_fatal(FARGS, "This Gromacs version currently only works with OpenMM. Use -device \"OpenMM:<options>\"");
327 else
329 ommOptions=strchr(deviceOptions,':');
330 if (NULL!=ommOptions)
332 /* Increase the pointer to skip the colon */
333 ommOptions++;
338 openmmData = openmm_init(fplog, ommOptions, ir, top_global, top, mdatoms, fr, state);
340 if (MASTER(cr))
342 /* Update mdebin with energy history if appending to output files */
343 if ( Flags & MD_APPENDFILES )
345 restore_energyhistory_from_state(mdebin,&state_global->enerhist);
347 /* Set the initial energy history in state to zero by updating once */
348 update_energyhistory(&state_global->enerhist,mdebin);
351 if (constr)
353 set_constraints(constr,top,ir,mdatoms,cr);
356 if (!ir->bContinuation)
358 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
360 /* Set the velocities of frozen particles to zero */
361 for (i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++)
363 for (m=0; m<DIM; m++)
365 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
367 state->v[i][m] = 0;
373 if (constr)
375 /* Constrain the initial coordinates and velocities */
376 do_constrain_first(fplog,constr,ir,mdatoms,state,f,
377 graph,cr,nrnb,fr,top,shake_vir);
379 if (vsite)
381 /* Construct the virtual sites for the initial configuration */
382 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,NULL,
383 top->idef.iparams,top->idef.il,
384 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
388 debug_gmx();
390 if (MASTER(cr))
392 char tbuf[20];
393 fprintf(fplog,"Initial temperature: %g K\n",enerd->term[F_TEMP]);
394 fprintf(stderr,"starting mdrun '%s'\n",
395 *(top_global->name));
396 if (ir->nsteps >= 0)
398 sprintf(tbuf,"%8.1f",(ir->init_step+ir->nsteps)*ir->delta_t);
400 else
402 sprintf(tbuf,"%s","infinite");
404 if (ir->init_step > 0)
406 fprintf(stderr,"%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
407 gmx_step_str(ir->init_step+ir->nsteps,sbuf),tbuf,
408 gmx_step_str(ir->init_step,sbuf2),
409 ir->init_step*ir->delta_t);
411 else
413 fprintf(stderr,"%s steps, %s ps.\n",
414 gmx_step_str(ir->nsteps,sbuf),tbuf);
418 fprintf(fplog,"\n");
420 /* Set and write start time */
421 runtime_start(runtime);
422 print_date_and_time(fplog,cr->nodeid,"Started mdrun",runtime);
423 wallcycle_start(wcycle,ewcRUN);
424 if (fplog)
425 fprintf(fplog,"\n");
427 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
429 debug_gmx();
430 /***********************************************************
432 * Loop over MD steps
434 ************************************************************/
436 /* loop over MD steps or if rerunMD to end of input trajectory */
437 bFirstStep = TRUE;
438 /* Skip the first Nose-Hoover integration when we get the state from tpx */
439 bStateFromTPX = !opt2bSet("-cpi",nfile,fnm);
440 bInitStep = bFirstStep && bStateFromTPX;
441 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
442 bLastStep = FALSE;
444 init_global_signals(&gs,cr,ir,repl_ex_nst);
446 step = ir->init_step;
447 step_rel = 0;
448 bLastStep = ((ir->nsteps >= 0 && step_rel > ir->nsteps));
450 while (!bLastStep)
452 wallcycle_start(wcycle,ewcSTEP);
454 GMX_MPE_LOG(ev_timestep1);
456 bLastStep = (step_rel == ir->nsteps);
457 t = t0 + step*ir->delta_t;
459 if (gs.set[eglsSTOPCOND] != 0)
461 bLastStep = TRUE;
464 do_log = do_per_step(step,ir->nstlog) || bFirstStep || bLastStep;
465 do_verbose = bVerbose &&
466 (step % stepout == 0 || bFirstStep || bLastStep);
468 if (MASTER(cr) && do_log)
470 print_ebin_header(fplog,step,t,state->lambda);
473 clear_mat(force_vir);
474 GMX_MPE_LOG(ev_timestep2);
476 /* We write a checkpoint at this MD step when:
477 * either when we signalled through gs (in OpenMM NS works different),
478 * or at the last step (but not when we do not want confout),
479 * but never at the first step.
481 bCPT = ((gs.set[eglsCHKPT] ||
482 (bLastStep && (Flags & MD_CONFOUT))) &&
483 step > ir->init_step );
484 if (bCPT)
486 gs.set[eglsCHKPT] = 0;
489 /* Now we have the energies and forces corresponding to the
490 * coordinates at time t. We must output all of this before
491 * the update.
492 * for RerunMD t is read from input trajectory
494 GMX_MPE_LOG(ev_output_start);
496 mdof_flags = 0;
497 if (do_per_step(step,ir->nstxout))
499 mdof_flags |= MDOF_X;
501 if (do_per_step(step,ir->nstvout))
503 mdof_flags |= MDOF_V;
505 if (do_per_step(step,ir->nstfout))
507 mdof_flags |= MDOF_F;
509 if (do_per_step(step,ir->nstxtcout))
511 mdof_flags |= MDOF_XTC;
513 if (bCPT)
515 mdof_flags |= MDOF_CPT;
517 do_ene = do_per_step(step,ir->nstenergy) || (bLastStep && ir->nstenergy);
519 if (mdof_flags != 0)
521 bX = (mdof_flags & (MDOF_X | MDOF_XTC | MDOF_CPT));
522 bV = (mdof_flags & (MDOF_V | MDOF_CPT));
524 wallcycle_start(wcycle,ewcTRAJ);
525 openmm_copy_state(openmmData, state, &t, f, enerd, bX, bV, 0, 0);
526 wallcycle_stop(wcycle,ewcTRAJ);
529 openmm_take_one_step(openmmData);
531 if (mdof_flags != 0 || do_ene)
533 wallcycle_start(wcycle,ewcTRAJ);
534 bF = (mdof_flags & MDOF_F);
535 if (bF || do_ene )
537 openmm_copy_state(openmmData, state, &t, f, enerd, 0, 0, bF, do_ene);
539 upd_mdebin(mdebin, NULL,TRUE,
540 t,mdatoms->tmass,enerd,state,lastbox,
541 shake_vir,force_vir,total_vir,pres,
542 ekind,mu_tot,constr);
543 print_ebin(outf->fp_ene,do_ene,FALSE,FALSE,fplog,step,t,
544 eprAVER,FALSE,mdebin,fcd,groups,&(ir->opts));
546 write_traj(fplog,cr,outf,mdof_flags,top_global,
547 step,t,state,state_global,f,f_global,&n_xtc,&x_xtc);
548 if (bCPT)
550 nchkpt++;
551 bCPT = FALSE;
553 debug_gmx();
554 if (bLastStep && step_rel == ir->nsteps &&
555 (Flags & MD_CONFOUT) && MASTER(cr))
557 /* x and v have been collected in write_traj,
558 * because a checkpoint file will always be written
559 * at the last step.
561 fprintf(stderr,"\nWriting final coordinates.\n");
562 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols)
564 /* Make molecules whole only for confout writing */
565 do_pbc_mtop(fplog,ir->ePBC,state->box,top_global,state_global->x);
567 write_sto_conf_mtop(ftp2fn(efSTO,nfile,fnm),
568 *top_global->name,top_global,
569 state_global->x,state_global->v,
570 ir->ePBC,state->box);
571 debug_gmx();
573 wallcycle_stop(wcycle,ewcTRAJ);
575 GMX_MPE_LOG(ev_output_finish);
578 /* Determine the wallclock run time up till now */
579 run_time = gmx_gettime() - (double)runtime->real;
581 /* Check whether everything is still allright */
582 if (((int)gmx_get_stop_condition() > handled_stop_condition) &&
583 MASTERTHREAD(cr))
585 /* this is just make gs.sig compatible with the hack
586 of sending signals around by MPI_Reduce with together with
587 other floats */
588 /* NOTE: this only works for serial code. For code that allows
589 MPI nodes to propagate their condition, see kernel/md.c*/
590 if ( gmx_get_stop_condition() == gmx_stop_cond_next_ns )
591 gs.set[eglsSTOPCOND]=1;
592 if ( gmx_get_stop_condition() == gmx_stop_cond_next )
593 gs.set[eglsSTOPCOND]=1;
594 /* < 0 means stop at next step, > 0 means stop at next NS step */
595 if (fplog)
597 fprintf(fplog,
598 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
599 gmx_get_signal_name(),
600 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
601 fflush(fplog);
603 fprintf(stderr,
604 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
605 gmx_get_signal_name(),
606 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
607 fflush(stderr);
608 handled_stop_condition=(int)gmx_get_stop_condition();
610 else if (MASTER(cr) &&
611 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
612 gs.set[eglsSTOPCOND] == 0)
614 /* Signal to terminate the run */
615 gs.set[eglsSTOPCOND] = 1;
616 if (fplog)
618 fprintf(fplog,"\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
620 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
623 /* checkpoints */
624 if (MASTER(cr) && (cpt_period >= 0 &&
625 (cpt_period == 0 ||
626 run_time >= nchkpt*cpt_period*60.0)) &&
627 gs.set[eglsCHKPT] == 0)
629 gs.set[eglsCHKPT] = 1;
632 /* Time for performance */
633 if (((step % stepout) == 0) || bLastStep)
635 runtime_upd_proc(runtime);
638 if (do_per_step(step,ir->nstlog))
640 if (fflush(fplog) != 0)
642 gmx_fatal(FARGS,"Cannot flush logfile - maybe you are out of quota?");
646 /* Remaining runtime */
647 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal() ))
649 print_time(stderr,runtime,step,ir,cr);
652 bFirstStep = FALSE;
653 bInitStep = FALSE;
654 bStartingFromCpt = FALSE;
655 step++;
656 step_rel++;
658 /* End of main MD loop */
659 debug_gmx();
661 /* Stop the time */
662 runtime_end(runtime);
664 openmm_cleanup(fplog, openmmData);
666 done_mdoutf(outf);
668 debug_gmx();
670 runtime->nsteps_done = step_rel;
672 return 0;