30e2e0f3ecf66b268117b8396f1b8e12a73af044
[alexxy/gromacs.git] / src / kernel / md_openmm.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
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.
14
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.
19  * 
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.
26  * 
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.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34  */
35
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <signal.h>
41 #include <stdlib.h>
42
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 "confio.h"
60 #include "network.h"
61 #include "pull.h"
62 #include "xvgr.h"
63 #include "physics.h"
64 #include "names.h"
65 #include "xmdrun.h"
66 #include "ionize.h"
67 #include "disre.h"
68 #include "orires.h"
69 #include "pme.h"
70 #include "mdatoms.h"
71 #include "qmmm.h"
72 #include "mpelogging.h"
73 #include "domdec.h"
74 #include "partdec.h"
75 #include "topsort.h"
76 #include "coulomb.h"
77 #include "constr.h"
78 #include "compute_io.h"
79 #include "mvdata.h"
80 #include "checkpoint.h"
81 #include "mtop_util.h"
82 #include "sighandler.h"
83 #include "genborn.h"
84 #include "string2.h"
85 #include "copyrite.h"
86 #include "membed.h"
87
88 #ifdef GMX_THREAD_MPI
89 #include "tmpi.h"
90 #endif
91
92 /* include even when OpenMM not used to force compilation of do_md_openmm */
93 #include "openmm_wrapper.h"
94
95 double do_md_openmm(FILE *fplog,t_commrec *cr,int nfile,const t_filenm fnm[],
96                     const output_env_t oenv, gmx_bool bVerbose,gmx_bool bCompact,
97                     int nstglobalcomm,
98                     gmx_vsite_t *vsite,gmx_constr_t constr,
99                     int stepout,t_inputrec *ir,
100                     gmx_mtop_t *top_global,
101                     t_fcdata *fcd,
102                     t_state *state_global,
103                     t_mdatoms *mdatoms,
104                     t_nrnb *nrnb,gmx_wallcycle_t wcycle,
105                     gmx_edsam_t ed,t_forcerec *fr,
106                     int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
107                     gmx_membed_t membed,
108                     real cpt_period,real max_hours,
109                     const char *deviceOptions,
110                     unsigned long Flags,
111                     gmx_runtime_t *runtime)
112 {
113     gmx_mdoutf_t *outf;
114     gmx_large_int_t step,step_rel;
115     double     run_time;
116     double     t,t0,lam0;
117     gmx_bool       bSimAnn,
118     bFirstStep,bStateFromTPX,bLastStep,bStartingFromCpt;
119     gmx_bool       bInitStep=TRUE;
120     gmx_bool       do_ene,do_log, do_verbose,
121     bX,bV,bF,bCPT;
122     tensor     force_vir,shake_vir,total_vir,pres;
123     int        i,m;
124     int        mdof_flags;
125     rvec       mu_tot;
126     t_vcm      *vcm;
127     int        nchkpt=1;
128     gmx_localtop_t *top;
129     t_mdebin *mdebin;
130     t_state    *state=NULL;
131     rvec       *f_global=NULL;
132     int        n_xtc=-1;
133     rvec       *x_xtc=NULL;
134     gmx_enerdata_t *enerd;
135     rvec       *f=NULL;
136     gmx_global_stat_t gstat;
137     gmx_update_t upd=NULL;
138     t_graph    *graph=NULL;
139     globsig_t   gs;
140
141     gmx_groups_t *groups;
142     gmx_ekindata_t *ekind, *ekind_save;
143     gmx_bool        bAppend;
144     int         a0,a1;
145     matrix      lastbox;
146     real        reset_counters=0,reset_counters_now=0;
147     char        sbuf[STEPSTRSIZE],sbuf2[STEPSTRSIZE];
148     int         handled_stop_condition=gmx_stop_cond_none; 
149
150     const char *ommOptions = NULL;
151     void   *openmmData;
152
153 #ifdef GMX_DOUBLE
154     /* Checks in cmake should prevent the compilation in double precision
155      * with OpenMM, but just to be sure we check here.
156      */
157     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.");
158 #endif
159
160     bAppend  = (Flags & MD_APPENDFILES);
161     check_ir_old_tpx_versions(cr,fplog,ir,top_global);
162
163     groups = &top_global->groups;
164
165     /* Initial values */
166     init_md(fplog,cr,ir,oenv,&t,&t0,state_global->lambda,
167             &(state_global->fep_state),&lam0,
168             nrnb,top_global,&upd,
169             nfile,fnm,&outf,&mdebin,
170             force_vir,shake_vir,mu_tot,&bSimAnn,&vcm,state_global,Flags);
171
172     clear_mat(total_vir);
173     clear_mat(pres);
174     /* Energy terms and groups */
175     snew(enerd,1);
176     init_enerdata(top_global->groups.grps[egcENER].nr,ir->fepvals->n_lambda,
177                   enerd);
178     snew(f,top_global->natoms);
179
180     /* Kinetic energy data */
181     snew(ekind,1);
182     init_ekindata(fplog,top_global,&(ir->opts),ekind);
183     /* needed for iteration of constraints */
184     snew(ekind_save,1);
185     init_ekindata(fplog,top_global,&(ir->opts),ekind_save);
186     /* Copy the cos acceleration to the groups struct */
187     ekind->cosacc.cos_accel = ir->cos_accel;
188
189     gstat = global_stat_init(ir);
190     debug_gmx();
191
192     {
193         double io = compute_io(ir,top_global->natoms,groups,mdebin->ebin->nener,1);
194         if ((io > 2000) && MASTER(cr))
195             fprintf(stderr,
196                     "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
197                     io);
198     }
199
200     top = gmx_mtop_generate_local_top(top_global,ir);
201
202     a0 = 0;
203     a1 = top_global->natoms;
204
205     state = partdec_init_local_state(cr,state_global);
206     f_global = f;
207
208     atoms2md(top_global,ir,0,NULL,a0,a1-a0,mdatoms);
209
210     if (vsite)
211     {
212         set_vsite_top(vsite,top,mdatoms,cr);
213     }
214
215     if (ir->ePBC != epbcNONE && !ir->bPeriodicMols)
216     {
217         graph = mk_graph(fplog,&(top->idef),0,top_global->natoms,FALSE,FALSE);
218     }
219
220     update_mdatoms(mdatoms,state->lambda[efptMASS]);
221
222     if (deviceOptions[0]=='\0')
223     {
224         /* empty options, which should default to OpenMM in this build */
225         ommOptions=deviceOptions;
226     }
227     else
228     {
229         if (gmx_strncasecmp(deviceOptions,"OpenMM",6)!=0)
230         {
231             gmx_fatal(FARGS, "This Gromacs version currently only works with OpenMM. Use -device \"OpenMM:<options>\"");
232         }
233         else
234         {
235             ommOptions=strchr(deviceOptions,':');
236             if (NULL!=ommOptions)
237             {
238                 /* Increase the pointer to skip the colon */
239                 ommOptions++;
240             }
241         }
242     }
243
244     openmmData = openmm_init(fplog, ommOptions, ir, top_global, top, mdatoms, fr, state);
245     please_cite(fplog,"Friedrichs2009");
246
247     if (MASTER(cr))
248     {
249         /* Update mdebin with energy history if appending to output files */
250         if ( Flags & MD_APPENDFILES )
251         {
252             restore_energyhistory_from_state(mdebin,&state_global->enerhist);
253         }
254         /* Set the initial energy history in state to zero by updating once */
255         update_energyhistory(&state_global->enerhist,mdebin);
256     }
257
258     if (constr)
259     {
260         set_constraints(constr,top,ir,mdatoms,cr);
261     }
262
263     if (!ir->bContinuation)
264     {
265         if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
266         {
267             /* Set the velocities of frozen particles to zero */
268             for (i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++)
269             {
270                 for (m=0; m<DIM; m++)
271                 {
272                     if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
273                     {
274                         state->v[i][m] = 0;
275                     }
276                 }
277             }
278         }
279
280         if (constr)
281         {
282             /* Constrain the initial coordinates and velocities */
283             do_constrain_first(fplog,constr,ir,mdatoms,state,f,
284                                graph,cr,nrnb,fr,top,shake_vir);
285         }
286         if (vsite)
287         {
288             /* Construct the virtual sites for the initial configuration */
289             construct_vsites(fplog,vsite,state->x,nrnb,ir->delta_t,NULL,
290                              top->idef.iparams,top->idef.il,
291                              fr->ePBC,fr->bMolPBC,graph,cr,state->box);
292         }
293     }
294
295     debug_gmx();
296
297     if (MASTER(cr))
298     {
299         char tbuf[20];
300         fprintf(stderr,"starting mdrun '%s'\n",
301                 *(top_global->name));
302         if (ir->nsteps >= 0)
303         {
304             sprintf(tbuf,"%8.1f",(ir->init_step+ir->nsteps)*ir->delta_t);
305         }
306         else
307         {
308             sprintf(tbuf,"%s","infinite");
309         }
310         if (ir->init_step > 0)
311         {
312             fprintf(stderr,"%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
313                     gmx_step_str(ir->init_step+ir->nsteps,sbuf),tbuf,
314                     gmx_step_str(ir->init_step,sbuf2),
315                     ir->init_step*ir->delta_t);
316         }
317         else
318         {
319             fprintf(stderr,"%s steps, %s ps.\n",
320                     gmx_step_str(ir->nsteps,sbuf),tbuf);
321         }
322     }
323
324     fprintf(fplog,"\n");
325
326     /* Set and write start time */
327     runtime_start(runtime);
328     print_date_and_time(fplog,cr->nodeid,"Started mdrun",runtime);
329     wallcycle_start(wcycle,ewcRUN);
330     if (fplog)
331         fprintf(fplog,"\n");
332
333     /* safest point to do file checkpointing is here.  More general point would be immediately before integrator call */
334
335     debug_gmx();
336     /***********************************************************
337      *
338      *             Loop over MD steps
339      *
340      ************************************************************/
341
342     /* loop over MD steps or if rerunMD to end of input trajectory */
343     bFirstStep = TRUE;
344     /* Skip the first Nose-Hoover integration when we get the state from tpx */
345     bStateFromTPX = !opt2bSet("-cpi",nfile,fnm);
346     bInitStep = bFirstStep && bStateFromTPX;
347     bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep;
348     bLastStep = FALSE;
349
350     init_global_signals(&gs,cr,ir,repl_ex_nst);
351
352     step = ir->init_step;
353     step_rel = 0;
354
355     while (!bLastStep)
356     {
357         wallcycle_start(wcycle,ewcSTEP);
358
359         GMX_MPE_LOG(ev_timestep1);
360
361         bLastStep = (step_rel == ir->nsteps);
362         t = t0 + step*ir->delta_t;
363
364         if (gs.set[eglsSTOPCOND] != 0)
365         {
366             bLastStep = TRUE;
367         }
368
369         do_log = do_per_step(step,ir->nstlog) || bFirstStep || bLastStep;
370         do_verbose = bVerbose &&
371                      (step % stepout == 0 || bFirstStep || bLastStep);
372
373         if (MASTER(cr) && do_log)
374         {
375             print_ebin_header(fplog,step,t,state->lambda[efptFEP]);
376         }
377
378         clear_mat(force_vir);
379         GMX_MPE_LOG(ev_timestep2);
380
381         /* We write a checkpoint at this MD step when:
382          * either when we signalled through gs (in OpenMM NS works different),
383          * or at the last step (but not when we do not want confout),
384          * but never at the first step.
385          */
386         bCPT = ((gs.set[eglsCHKPT] ||
387                  (bLastStep && (Flags & MD_CONFOUT))) &&
388                 step > ir->init_step );
389         if (bCPT)
390         {
391             gs.set[eglsCHKPT] = 0;
392         }
393
394         /* Now we have the energies and forces corresponding to the
395          * coordinates at time t. We must output all of this before
396          * the update.
397          * for RerunMD t is read from input trajectory
398          */
399         GMX_MPE_LOG(ev_output_start);
400
401         mdof_flags = 0;
402         if (do_per_step(step,ir->nstxout))
403         {
404             mdof_flags |= MDOF_X;
405         }
406         if (do_per_step(step,ir->nstvout))
407         {
408             mdof_flags |= MDOF_V;
409         }
410         if (do_per_step(step,ir->nstfout))
411         {
412             mdof_flags |= MDOF_F;
413         }
414         if (do_per_step(step,ir->nstxtcout))
415         {
416             mdof_flags |= MDOF_XTC;
417         }
418         if (bCPT)
419         {
420             mdof_flags |= MDOF_CPT;
421         };
422         do_ene = (do_per_step(step,ir->nstenergy) || bLastStep);
423
424         if (mdof_flags != 0 || do_ene || do_log)
425         {
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));
430
431             openmm_copy_state(openmmData, state, &t, f, enerd, bX, bV, bF, do_ene);
432
433             upd_mdebin(mdebin,FALSE,TRUE,
434                        t,mdatoms->tmass,enerd,state,ir->fepvals,ir->expandedvals,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,
438                        step,t,
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);
442             if (bCPT)
443             {
444                 nchkpt++;
445                 bCPT = FALSE;
446             }
447             debug_gmx();
448             if (bLastStep && step_rel == ir->nsteps &&
449                     (Flags & MD_CONFOUT) && MASTER(cr))
450             {
451                 /* x and v have been collected in write_traj,
452                  * because a checkpoint file will always be written
453                  * at the last step.
454                  */
455                 fprintf(stderr,"\nWriting final coordinates.\n");
456                 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols)
457                 {
458                     /* Make molecules whole only for confout writing */
459                     do_pbc_mtop(fplog,ir->ePBC,state->box,top_global,state_global->x);
460                 }
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);
465                 debug_gmx();
466             }
467             wallcycle_stop(wcycle,ewcTRAJ);
468         }
469         GMX_MPE_LOG(ev_output_finish);
470
471
472         /* Determine the wallclock run time up till now */
473         run_time = gmx_gettime() - (double)runtime->real;
474
475         /* Check whether everything is still allright */
476         if (((int)gmx_get_stop_condition() > handled_stop_condition)
477 #ifdef GMX_THREAD_MPI
478             && MASTER(cr)
479 #endif
480             )
481         {
482            /* this is just make gs.sig compatible with the hack 
483                of sending signals around by MPI_Reduce with together with
484                other floats */
485             /* NOTE: this only works for serial code. For code that allows
486                MPI nodes to propagate their condition, see kernel/md.c*/
487             if ( gmx_get_stop_condition() == gmx_stop_cond_next_ns )
488                 gs.set[eglsSTOPCOND]=1;
489             if ( gmx_get_stop_condition() == gmx_stop_cond_next )
490                 gs.set[eglsSTOPCOND]=1;
491             /* < 0 means stop at next step, > 0 means stop at next NS step */
492             if (fplog)
493             {
494                 fprintf(fplog,
495                         "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
496                         gmx_get_signal_name(),
497                         gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
498                 fflush(fplog);
499             }
500             fprintf(stderr,
501                     "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
502                     gmx_get_signal_name(),
503                     gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
504             fflush(stderr);
505             handled_stop_condition=(int)gmx_get_stop_condition();
506         }
507         else if (MASTER(cr) &&
508                  (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
509                  gs.set[eglsSTOPCOND] == 0)
510         {
511             /* Signal to terminate the run */
512             gs.set[eglsSTOPCOND] = 1;
513             if (fplog)
514             {
515                 fprintf(fplog,"\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
516             }
517             fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
518         }
519
520         /* checkpoints */
521         if (MASTER(cr) && (cpt_period >= 0 &&
522                            (cpt_period == 0 ||
523                             run_time >= nchkpt*cpt_period*60.0)) &&
524                 gs.set[eglsCHKPT] == 0)
525         {
526             gs.set[eglsCHKPT] = 1;
527         }
528
529         /* Time for performance */
530         if (((step % stepout) == 0) || bLastStep)
531         {
532             runtime_upd_proc(runtime);
533         }
534
535         if (do_per_step(step,ir->nstlog))
536         {
537             if (fflush(fplog) != 0)
538             {
539                 gmx_fatal(FARGS,"Cannot flush logfile - maybe you are out of disk space?");
540             }
541         }
542
543         /* Remaining runtime */
544         if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal() ))
545         {
546             print_time(stderr,runtime,step,ir,cr);
547         }
548
549         bFirstStep = FALSE;
550         bInitStep = FALSE;
551         bStartingFromCpt = FALSE;
552         step++;
553         step_rel++;
554
555         openmm_take_one_step(openmmData);
556     }
557     /* End of main MD loop */
558     debug_gmx();
559
560     /* Stop the time */
561     runtime_end(runtime);
562
563     if (MASTER(cr))
564     {
565         if (ir->nstcalcenergy > 0) 
566         {
567             print_ebin(outf->fp_ene,FALSE,FALSE,FALSE,fplog,step,t,
568                        eprAVER,FALSE,mdebin,fcd,groups,&(ir->opts));
569         }
570     }
571
572     openmm_cleanup(fplog, openmmData);
573
574     done_mdoutf(outf);
575
576     debug_gmx();
577
578     runtime->nsteps_done = step_rel;
579
580     return 0;
581 }