removed group non-boneded call with verlet scheme
[gromacs/AngularHB.git] / src / kernel / md_openmm.c
blob80df1d1827d6a91d8f3e8d2cd5b9be912dbbfdcf
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 #include "typedefs.h"
44 #include "smalloc.h"
45 #include "sysstuff.h"
46 #include "vec.h"
47 #include "statutil.h"
48 #include "vcm.h"
49 #include "mdebin.h"
50 #include "nrnb.h"
51 #include "calcmu.h"
52 #include "index.h"
53 #include "vsite.h"
54 #include "update.h"
55 #include "ns.h"
56 #include "trnio.h"
57 #include "xtcio.h"
58 #include "mdrun.h"
59 #include "md_support.h"
60 #include "confio.h"
61 #include "network.h"
62 #include "pull.h"
63 #include "xvgr.h"
64 #include "physics.h"
65 #include "names.h"
66 #include "xmdrun.h"
67 #include "ionize.h"
68 #include "disre.h"
69 #include "orires.h"
70 #include "pme.h"
71 #include "mdatoms.h"
72 #include "qmmm.h"
73 #include "mpelogging.h"
74 #include "domdec.h"
75 #include "partdec.h"
76 #include "topsort.h"
77 #include "coulomb.h"
78 #include "constr.h"
79 #include "compute_io.h"
80 #include "mvdata.h"
81 #include "checkpoint.h"
82 #include "mtop_util.h"
83 #include "sighandler.h"
84 #include "genborn.h"
85 #include "string2.h"
86 #include "copyrite.h"
87 #include "membed.h"
89 #ifdef GMX_THREAD_MPI
90 #include "tmpi.h"
91 #endif
93 /* include even when OpenMM not used to force compilation of do_md_openmm */
94 #include "openmm_wrapper.h"
96 double do_md_openmm(FILE *fplog,t_commrec *cr,int nfile,const t_filenm fnm[],
97 const output_env_t oenv, gmx_bool bVerbose,gmx_bool bCompact,
98 int nstglobalcomm,
99 gmx_vsite_t *vsite,gmx_constr_t constr,
100 int stepout,t_inputrec *ir,
101 gmx_mtop_t *top_global,
102 t_fcdata *fcd,
103 t_state *state_global,
104 t_mdatoms *mdatoms,
105 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
106 gmx_edsam_t ed,t_forcerec *fr,
107 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
108 gmx_membed_t membed,
109 real cpt_period,real max_hours,
110 const char *deviceOptions,
111 unsigned long Flags,
112 gmx_runtime_t *runtime)
114 gmx_mdoutf_t *outf;
115 gmx_large_int_t step,step_rel;
116 double run_time;
117 double t,t0,lam0;
118 gmx_bool bSimAnn,
119 bFirstStep,bStateFromTPX,bLastStep,bStartingFromCpt;
120 gmx_bool bInitStep=TRUE;
121 gmx_bool do_ene,do_log, do_verbose,
122 bX,bV,bF,bCPT;
123 tensor force_vir,shake_vir,total_vir,pres;
124 int i,m;
125 int mdof_flags;
126 rvec mu_tot;
127 t_vcm *vcm;
128 int nchkpt=1;
129 gmx_localtop_t *top;
130 t_mdebin *mdebin;
131 t_state *state=NULL;
132 rvec *f_global=NULL;
133 int n_xtc=-1;
134 rvec *x_xtc=NULL;
135 gmx_enerdata_t *enerd;
136 rvec *f=NULL;
137 gmx_global_stat_t gstat;
138 gmx_update_t upd=NULL;
139 t_graph *graph=NULL;
140 globsig_t gs;
142 gmx_groups_t *groups;
143 gmx_ekindata_t *ekind, *ekind_save;
144 gmx_bool bAppend;
145 int a0,a1;
146 matrix lastbox;
147 real reset_counters=0,reset_counters_now=0;
148 char sbuf[STEPSTRSIZE],sbuf2[STEPSTRSIZE];
149 int handled_stop_condition=gmx_stop_cond_none;
151 const char *ommOptions = NULL;
152 void *openmmData;
154 #ifdef GMX_DOUBLE
155 /* Checks in cmake should prevent the compilation in double precision
156 * with OpenMM, but just to be sure we check here.
158 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.");
159 #endif
161 bAppend = (Flags & MD_APPENDFILES);
162 check_ir_old_tpx_versions(cr,fplog,ir,top_global);
164 groups = &top_global->groups;
166 /* Initial values */
167 init_md(fplog,cr,ir,oenv,&t,&t0,state_global->lambda,
168 &(state_global->fep_state),&lam0,
169 nrnb,top_global,&upd,
170 nfile,fnm,&outf,&mdebin,
171 force_vir,shake_vir,mu_tot,&bSimAnn,&vcm,state_global,Flags);
173 clear_mat(total_vir);
174 clear_mat(pres);
175 /* Energy terms and groups */
176 snew(enerd,1);
177 init_enerdata(top_global->groups.grps[egcENER].nr,ir->fepvals->n_lambda,
178 enerd);
179 snew(f,top_global->natoms);
181 /* Kinetic energy data */
182 snew(ekind,1);
183 init_ekindata(fplog,top_global,&(ir->opts),ekind);
184 /* needed for iteration of constraints */
185 snew(ekind_save,1);
186 init_ekindata(fplog,top_global,&(ir->opts),ekind_save);
187 /* Copy the cos acceleration to the groups struct */
188 ekind->cosacc.cos_accel = ir->cos_accel;
190 gstat = global_stat_init(ir);
191 debug_gmx();
194 double io = compute_io(ir,top_global->natoms,groups,mdebin->ebin->nener,1);
195 if ((io > 2000) && MASTER(cr))
196 fprintf(stderr,
197 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
198 io);
201 top = gmx_mtop_generate_local_top(top_global,ir);
203 a0 = 0;
204 a1 = top_global->natoms;
206 state = partdec_init_local_state(cr,state_global);
207 f_global = f;
209 atoms2md(top_global,ir,0,NULL,a0,a1-a0,mdatoms);
211 if (vsite)
213 set_vsite_top(vsite,top,mdatoms,cr);
216 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols)
218 graph = mk_graph(fplog,&(top->idef),0,top_global->natoms,FALSE,FALSE);
221 update_mdatoms(mdatoms,state->lambda[efptMASS]);
223 if (deviceOptions[0]=='\0')
225 /* empty options, which should default to OpenMM in this build */
226 ommOptions=deviceOptions;
228 else
230 if (gmx_strncasecmp(deviceOptions,"OpenMM",6)!=0)
232 gmx_fatal(FARGS, "This Gromacs version currently only works with OpenMM. Use -device \"OpenMM:<options>\"");
234 else
236 ommOptions=strchr(deviceOptions,':');
237 if (NULL!=ommOptions)
239 /* Increase the pointer to skip the colon */
240 ommOptions++;
245 openmmData = openmm_init(fplog, ommOptions, ir, top_global, top, mdatoms, fr, state);
246 please_cite(fplog,"Friedrichs2009");
248 if (MASTER(cr))
250 /* Update mdebin with energy history if appending to output files */
251 if ( Flags & MD_APPENDFILES )
253 restore_energyhistory_from_state(mdebin,&state_global->enerhist);
255 /* Set the initial energy history in state to zero by updating once */
256 update_energyhistory(&state_global->enerhist,mdebin);
259 if (constr)
261 set_constraints(constr,top,ir,mdatoms,cr);
264 if (!ir->bContinuation)
266 if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
268 /* Set the velocities of frozen particles to zero */
269 for (i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++)
271 for (m=0; m<DIM; m++)
273 if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
275 state->v[i][m] = 0;
281 if (constr)
283 /* Constrain the initial coordinates and velocities */
284 do_constrain_first(fplog,constr,ir,mdatoms,state,f,
285 graph,cr,nrnb,fr,top,shake_vir);
287 if (vsite)
289 /* Construct the virtual sites for the initial configuration */
290 construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,NULL,
291 top->idef.iparams,top->idef.il,
292 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
296 debug_gmx();
298 if (MASTER(cr))
300 char tbuf[20];
301 fprintf(stderr,"starting mdrun '%s'\n",
302 *(top_global->name));
303 if (ir->nsteps >= 0)
305 sprintf(tbuf,"%8.1f",(ir->init_step+ir->nsteps)*ir->delta_t);
307 else
309 sprintf(tbuf,"%s","infinite");
311 if (ir->init_step > 0)
313 fprintf(stderr,"%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
314 gmx_step_str(ir->init_step+ir->nsteps,sbuf),tbuf,
315 gmx_step_str(ir->init_step,sbuf2),
316 ir->init_step*ir->delta_t);
318 else
320 fprintf(stderr,"%s steps, %s ps.\n",
321 gmx_step_str(ir->nsteps,sbuf),tbuf);
325 fprintf(fplog,"\n");
327 /* Set and write start time */
328 runtime_start(runtime);
329 print_date_and_time(fplog,cr->nodeid,"Started mdrun",runtime);
330 wallcycle_start(wcycle,ewcRUN);
331 if (fplog)
332 fprintf(fplog,"\n");
334 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
336 debug_gmx();
337 /***********************************************************
339 * Loop over MD steps
341 ************************************************************/
343 /* loop over MD steps or if rerunMD to end of input trajectory */
344 bFirstStep = TRUE;
345 /* Skip the first Nose-Hoover integration when we get the state from tpx */
346 bStateFromTPX = !opt2bSet("-cpi",nfile,fnm);
347 bInitStep = bFirstStep && bStateFromTPX;
348 bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
349 bLastStep = FALSE;
351 init_global_signals(&gs,cr,ir,repl_ex_nst);
353 step = ir->init_step;
354 step_rel = 0;
356 while (!bLastStep)
358 wallcycle_start(wcycle,ewcSTEP);
360 GMX_MPE_LOG(ev_timestep1);
362 bLastStep = (step_rel == ir->nsteps);
363 t = t0 + step*ir->delta_t;
365 if (gs.set[eglsSTOPCOND] != 0)
367 bLastStep = TRUE;
370 do_log = do_per_step(step,ir->nstlog) || bFirstStep || bLastStep;
371 do_verbose = bVerbose &&
372 (step % stepout == 0 || bFirstStep || bLastStep);
374 if (MASTER(cr) && do_log)
376 print_ebin_header(fplog,step,t,state->lambda[efptFEP]);
379 clear_mat(force_vir);
380 GMX_MPE_LOG(ev_timestep2);
382 /* We write a checkpoint at this MD step when:
383 * either when we signalled through gs (in OpenMM NS works different),
384 * or at the last step (but not when we do not want confout),
385 * but never at the first step.
387 bCPT = ((gs.set[eglsCHKPT] ||
388 (bLastStep && (Flags & MD_CONFOUT))) &&
389 step > ir->init_step );
390 if (bCPT)
392 gs.set[eglsCHKPT] = 0;
395 /* Now we have the energies and forces corresponding to the
396 * coordinates at time t. We must output all of this before
397 * the update.
398 * for RerunMD t is read from input trajectory
400 GMX_MPE_LOG(ev_output_start);
402 mdof_flags = 0;
403 if (do_per_step(step,ir->nstxout))
405 mdof_flags |= MDOF_X;
407 if (do_per_step(step,ir->nstvout))
409 mdof_flags |= MDOF_V;
411 if (do_per_step(step,ir->nstfout))
413 mdof_flags |= MDOF_F;
415 if (do_per_step(step,ir->nstxtcout))
417 mdof_flags |= MDOF_XTC;
419 if (bCPT)
421 mdof_flags |= MDOF_CPT;
423 do_ene = (do_per_step(step,ir->nstenergy) || bLastStep);
425 if (mdof_flags != 0 || do_ene || do_log)
427 wallcycle_start(wcycle,ewcTRAJ);
428 bF = (mdof_flags & MDOF_F);
429 bX = (mdof_flags & (MDOF_X | MDOF_XTC | MDOF_CPT));
430 bV = (mdof_flags & (MDOF_V | MDOF_CPT));
432 openmm_copy_state(openmmData, state, &t, f, enerd, bX, bV, bF, do_ene);
434 upd_mdebin(mdebin,FALSE,TRUE,
435 t,mdatoms->tmass,enerd,state,ir->fepvals,ir->expandedvals,lastbox,
436 shake_vir,force_vir,total_vir,pres,
437 ekind,mu_tot,constr);
438 print_ebin(outf->fp_ene,do_ene,FALSE,FALSE,do_log?fplog:NULL,
439 step,t,
440 eprNORMAL,bCompact,mdebin,fcd,groups,&(ir->opts));
441 write_traj(fplog,cr,outf,mdof_flags,top_global,
442 step,t,state,state_global,f,f_global,&n_xtc,&x_xtc);
443 if (bCPT)
445 nchkpt++;
446 bCPT = FALSE;
448 debug_gmx();
449 if (bLastStep && step_rel == ir->nsteps &&
450 (Flags & MD_CONFOUT) && MASTER(cr))
452 /* x and v have been collected in write_traj,
453 * because a checkpoint file will always be written
454 * at the last step.
456 fprintf(stderr,"\nWriting final coordinates.\n");
457 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols)
459 /* Make molecules whole only for confout writing */
460 do_pbc_mtop(fplog,ir->ePBC,state->box,top_global,state_global->x);
462 write_sto_conf_mtop(ftp2fn(efSTO,nfile,fnm),
463 *top_global->name,top_global,
464 state_global->x,state_global->v,
465 ir->ePBC,state->box);
466 debug_gmx();
468 wallcycle_stop(wcycle,ewcTRAJ);
470 GMX_MPE_LOG(ev_output_finish);
473 /* Determine the wallclock run time up till now */
474 run_time = gmx_gettime() - (double)runtime->real;
476 /* Check whether everything is still allright */
477 if (((int)gmx_get_stop_condition() > handled_stop_condition)
478 #ifdef GMX_THREAD_MPI
479 && MASTER(cr)
480 #endif
483 /* this is just make gs.sig compatible with the hack
484 of sending signals around by MPI_Reduce with together with
485 other floats */
486 /* NOTE: this only works for serial code. For code that allows
487 MPI nodes to propagate their condition, see kernel/md.c*/
488 if ( gmx_get_stop_condition() == gmx_stop_cond_next_ns )
489 gs.set[eglsSTOPCOND]=1;
490 if ( gmx_get_stop_condition() == gmx_stop_cond_next )
491 gs.set[eglsSTOPCOND]=1;
492 /* < 0 means stop at next step, > 0 means stop at next NS step */
493 if (fplog)
495 fprintf(fplog,
496 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
497 gmx_get_signal_name(),
498 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
499 fflush(fplog);
501 fprintf(stderr,
502 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
503 gmx_get_signal_name(),
504 gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
505 fflush(stderr);
506 handled_stop_condition=(int)gmx_get_stop_condition();
508 else if (MASTER(cr) &&
509 (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
510 gs.set[eglsSTOPCOND] == 0)
512 /* Signal to terminate the run */
513 gs.set[eglsSTOPCOND] = 1;
514 if (fplog)
516 fprintf(fplog,"\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
518 fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
521 /* checkpoints */
522 if (MASTER(cr) && (cpt_period >= 0 &&
523 (cpt_period == 0 ||
524 run_time >= nchkpt*cpt_period*60.0)) &&
525 gs.set[eglsCHKPT] == 0)
527 gs.set[eglsCHKPT] = 1;
530 /* Time for performance */
531 if (((step % stepout) == 0) || bLastStep)
533 runtime_upd_proc(runtime);
536 if (do_per_step(step,ir->nstlog))
538 if (fflush(fplog) != 0)
540 gmx_fatal(FARGS,"Cannot flush logfile - maybe you are out of disk space?");
544 /* Remaining runtime */
545 if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal() ))
547 print_time(stderr,runtime,step,ir,cr);
550 bFirstStep = FALSE;
551 bInitStep = FALSE;
552 bStartingFromCpt = FALSE;
553 step++;
554 step_rel++;
556 openmm_take_one_step(openmmData);
558 /* End of main MD loop */
559 debug_gmx();
561 /* Stop the time */
562 runtime_end(runtime);
564 if (MASTER(cr))
566 if (ir->nstcalcenergy > 0)
568 print_ebin(outf->fp_ene,FALSE,FALSE,FALSE,fplog,step,t,
569 eprAVER,FALSE,mdebin,fcd,groups,&(ir->opts));
573 openmm_cleanup(fplog, openmmData);
575 done_mdoutf(outf);
577 debug_gmx();
579 runtime->nsteps_done = step_rel;
581 return 0;