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