Removed old-fashioned MPE logging code
[alexxy/gromacs.git] / src / programs / mdrun / 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 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
44 /* _isnan() */
45 #include <float.h>
46 #endif
47
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 "domdec.h"
80 #include "partdec.h"
81 #include "topsort.h"
82 #include "coulomb.h"
83 #include "constr.h"
84 #include "compute_io.h"
85 #include "mvdata.h"
86 #include "checkpoint.h"
87 #include "mtop_util.h"
88 #include "sighandler.h"
89 #include "genborn.h"
90 #include "string2.h"
91 #include "copyrite.h"
92 #include "membed.h"
93
94 #ifdef GMX_THREADS
95 #include "tmpi.h"
96 #endif
97
98 /* include even when OpenMM not used to force compilation of do_md_openmm */
99 #include "openmm_wrapper.h"
100
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,
103                     int nstglobalcomm,
104                     gmx_vsite_t *vsite,gmx_constr_t constr,
105                     int stepout,t_inputrec *ir,
106                     gmx_mtop_t *top_global,
107                     t_fcdata *fcd,
108                     t_state *state_global,
109                     t_mdatoms *mdatoms,
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,
116                     unsigned long Flags,
117                     gmx_runtime_t *runtime)
118 {
119     gmx_mdoutf_t *outf;
120     gmx_large_int_t step,step_rel;
121     double     run_time;
122     double     t,t0,lam0;
123     gmx_bool       bSimAnn,
124     bFirstStep,bStateFromTPX,bLastStep,bStartingFromCpt;
125     gmx_bool       bInitStep=TRUE;
126     gmx_bool       do_ene,do_log, do_verbose,
127     bX,bV,bF,bCPT;
128     tensor     force_vir,shake_vir,total_vir,pres;
129     int        i,m;
130     int        mdof_flags;
131     rvec       mu_tot;
132     t_vcm      *vcm;
133     int        nchkpt=1;
134     gmx_localtop_t *top;
135     t_mdebin *mdebin=NULL;
136     t_state    *state=NULL;
137     rvec       *f_global=NULL;
138     int        n_xtc=-1;
139     rvec       *x_xtc=NULL;
140     gmx_enerdata_t *enerd;
141     rvec       *f=NULL;
142     gmx_global_stat_t gstat;
143     gmx_update_t upd=NULL;
144     t_graph    *graph=NULL;
145     globsig_t   gs;
146
147     gmx_groups_t *groups;
148     gmx_ekindata_t *ekind, *ekind_save;
149     gmx_bool        bAppend;
150     int         a0,a1;
151     matrix      lastbox;
152     real        reset_counters=0,reset_counters_now=0;
153     char        sbuf[STEPSTRSIZE],sbuf2[STEPSTRSIZE];
154     int         handled_stop_condition=gmx_stop_cond_none; 
155
156     const char *ommOptions = NULL;
157     void   *openmmData;
158
159 #ifdef GMX_DOUBLE
160     /* Checks in cmake should prevent the compilation in double precision
161      * with OpenMM, but just to be sure we check here.
162      */
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.");
164 #endif
165
166     bAppend  = (Flags & MD_APPENDFILES);
167     check_ir_old_tpx_versions(cr,fplog,ir,top_global);
168
169     groups = &top_global->groups;
170
171     /* Initial values */
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);
176
177     clear_mat(total_vir);
178     clear_mat(pres);
179     /* Energy terms and groups */
180     snew(enerd,1);
181     init_enerdata(top_global->groups.grps[egcENER].nr,ir->n_flambda,enerd);
182     snew(f,top_global->natoms);
183
184     /* Kinetic energy data */
185     snew(ekind,1);
186     init_ekindata(fplog,top_global,&(ir->opts),ekind);
187     /* needed for iteration of constraints */
188     snew(ekind_save,1);
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;
192
193     gstat = global_stat_init(ir);
194     debug_gmx();
195
196     {
197         double io = compute_io(ir,top_global->natoms,groups,mdebin->ebin->nener,1);
198         if ((io > 2000) && MASTER(cr))
199             fprintf(stderr,
200                     "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
201                     io);
202     }
203
204     top = gmx_mtop_generate_local_top(top_global,ir);
205
206     a0 = 0;
207     a1 = top_global->natoms;
208
209     state = partdec_init_local_state(cr,state_global);
210     f_global = f;
211
212     atoms2md(top_global,ir,0,NULL,a0,a1-a0,mdatoms);
213
214     if (vsite)
215     {
216         set_vsite_top(vsite,top,mdatoms,cr);
217     }
218
219     if (ir->ePBC != epbcNONE && !ir->bPeriodicMols)
220     {
221         graph = mk_graph(fplog,&(top->idef),0,top_global->natoms,FALSE,FALSE);
222     }
223
224     update_mdatoms(mdatoms,state->lambda);
225
226     if (deviceOptions[0]=='\0')
227     {
228         /* empty options, which should default to OpenMM in this build */
229         ommOptions=deviceOptions;
230     }
231     else
232     {
233         if (gmx_strncasecmp(deviceOptions,"OpenMM",6)!=0)
234         {
235             gmx_fatal(FARGS, "This Gromacs version currently only works with OpenMM. Use -device \"OpenMM:<options>\"");
236         }
237         else
238         {
239             ommOptions=strchr(deviceOptions,':');
240             if (NULL!=ommOptions)
241             {
242                 /* Increase the pointer to skip the colon */
243                 ommOptions++;
244             }
245         }
246     }
247
248     openmmData = openmm_init(fplog, ommOptions, ir, top_global, top, mdatoms, fr, state);
249     please_cite(fplog,"Friedrichs2009");
250
251     if (MASTER(cr))
252     {
253         /* Update mdebin with energy history if appending to output files */
254         if ( Flags & MD_APPENDFILES )
255         {
256             restore_energyhistory_from_state(mdebin,&state_global->enerhist);
257         }
258         /* Set the initial energy history in state to zero by updating once */
259         update_energyhistory(&state_global->enerhist,mdebin);
260     }
261
262     if (constr)
263     {
264         set_constraints(constr,top,ir,mdatoms,cr);
265     }
266
267     if (!ir->bContinuation)
268     {
269         if (mdatoms->cFREEZE && (state->flags & (1<<estV)))
270         {
271             /* Set the velocities of frozen particles to zero */
272             for (i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++)
273             {
274                 for (m=0; m<DIM; m++)
275                 {
276                     if (ir->opts.nFreeze[mdatoms->cFREEZE[i]][m])
277                     {
278                         state->v[i][m] = 0;
279                     }
280                 }
281             }
282         }
283
284         if (constr)
285         {
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);
289         }
290         if (vsite)
291         {
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);
296         }
297     }
298
299     debug_gmx();
300
301     if (MASTER(cr))
302     {
303         char tbuf[20];
304         fprintf(fplog,"Initial temperature: %g K\n",enerd->term[F_TEMP]);
305         fprintf(stderr,"starting mdrun '%s'\n",
306                 *(top_global->name));
307         if (ir->nsteps >= 0)
308         {
309             sprintf(tbuf,"%8.1f",(ir->init_step+ir->nsteps)*ir->delta_t);
310         }
311         else
312         {
313             sprintf(tbuf,"%s","infinite");
314         }
315         if (ir->init_step > 0)
316         {
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);
321         }
322         else
323         {
324             fprintf(stderr,"%s steps, %s ps.\n",
325                     gmx_step_str(ir->nsteps,sbuf),tbuf);
326         }
327     }
328
329     fprintf(fplog,"\n");
330
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);
335     if (fplog)
336         fprintf(fplog,"\n");
337
338     /* safest point to do file checkpointing is here.  More general point would be immediately before integrator call */
339
340     debug_gmx();
341     /***********************************************************
342      *
343      *             Loop over MD steps
344      *
345      ************************************************************/
346
347     /* loop over MD steps or if rerunMD to end of input trajectory */
348     bFirstStep = TRUE;
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;
353     bLastStep = FALSE;
354
355     init_global_signals(&gs,cr,ir,repl_ex_nst);
356
357     step = ir->init_step;
358     step_rel = 0;
359
360     while (!bLastStep)
361     {
362         wallcycle_start(wcycle,ewcSTEP);
363
364         bLastStep = (step_rel == ir->nsteps);
365         t = t0 + step*ir->delta_t;
366
367         if (gs.set[eglsSTOPCOND] != 0)
368         {
369             bLastStep = TRUE;
370         }
371
372         do_log = do_per_step(step,ir->nstlog) || bFirstStep || bLastStep;
373         do_verbose = bVerbose &&
374                      (step % stepout == 0 || bFirstStep || bLastStep);
375
376         if (MASTER(cr) && do_log)
377         {
378             print_ebin_header(fplog,step,t,state->lambda);
379         }
380
381         clear_mat(force_vir);
382
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.
387          */
388         bCPT = ((gs.set[eglsCHKPT] ||
389                  (bLastStep && (Flags & MD_CONFOUT))) &&
390                 step > ir->init_step );
391         if (bCPT)
392         {
393             gs.set[eglsCHKPT] = 0;
394         }
395
396         /* Now we have the energies and forces corresponding to the
397          * coordinates at time t. We must output all of this before
398          * the update.
399          * for RerunMD t is read from input trajectory
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,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
470         /* Determine the wallclock run time up till now */
471         run_time = gmx_gettime() - (double)runtime->real;
472
473         /* Check whether everything is still allright */
474         if (((int)gmx_get_stop_condition() > handled_stop_condition)
475 #ifdef GMX_THREADS
476             && MASTER(cr)
477 #endif
478             )
479         {
480            /* this is just make gs.sig compatible with the hack 
481                of sending signals around by MPI_Reduce with together with
482                other floats */
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 */
490             if (fplog)
491             {
492                 fprintf(fplog,
493                         "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
494                         gmx_get_signal_name(),
495                         gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
496                 fflush(fplog);
497             }
498             fprintf(stderr,
499                     "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
500                     gmx_get_signal_name(),
501                     gs.sig[eglsSTOPCOND]==1 ? "NS " : "");
502             fflush(stderr);
503             handled_stop_condition=(int)gmx_get_stop_condition();
504         }
505         else if (MASTER(cr) &&
506                  (max_hours > 0 && run_time > max_hours*60.0*60.0*0.99) &&
507                  gs.set[eglsSTOPCOND] == 0)
508         {
509             /* Signal to terminate the run */
510             gs.set[eglsSTOPCOND] = 1;
511             if (fplog)
512             {
513                 fprintf(fplog,"\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
514             }
515             fprintf(stderr, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n",gmx_step_str(step,sbuf),max_hours*0.99);
516         }
517
518         /* checkpoints */
519         if (MASTER(cr) && (cpt_period >= 0 &&
520                            (cpt_period == 0 ||
521                             run_time >= nchkpt*cpt_period*60.0)) &&
522                 gs.set[eglsCHKPT] == 0)
523         {
524             gs.set[eglsCHKPT] = 1;
525         }
526
527         /* Time for performance */
528         if (((step % stepout) == 0) || bLastStep)
529         {
530             runtime_upd_proc(runtime);
531         }
532
533         if (do_per_step(step,ir->nstlog))
534         {
535             if (fflush(fplog) != 0)
536             {
537                 gmx_fatal(FARGS,"Cannot flush logfile - maybe you are out of quota?");
538             }
539         }
540
541         /* Remaining runtime */
542         if (MULTIMASTER(cr) && (do_verbose || gmx_got_usr_signal() ))
543         {
544             print_time(stderr,runtime,step,ir,cr);
545         }
546
547         bFirstStep = FALSE;
548         bInitStep = FALSE;
549         bStartingFromCpt = FALSE;
550         step++;
551         step_rel++;
552
553         openmm_take_one_step(openmmData);
554     }
555     /* End of main MD loop */
556     debug_gmx();
557
558     /* Stop the time */
559     runtime_end(runtime);
560
561     if (MASTER(cr))
562     {
563         if (ir->nstcalcenergy > 0) 
564         {
565             print_ebin(outf->fp_ene,FALSE,FALSE,FALSE,fplog,step,t,
566                        eprAVER,FALSE,mdebin,fcd,groups,&(ir->opts));
567         }
568     }
569
570     openmm_cleanup(fplog, openmmData);
571
572     done_mdoutf(outf);
573
574     debug_gmx();
575
576     runtime->nsteps_done = step_rel;
577
578     return 0;
579 }