Merge release-4-5-patches into release-4-6
[alexxy/gromacs.git] / src / mdlib / sim_util.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  *                        VERSION 3.2.0
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2004, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  * 
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  * 
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  * 
31  * For more info, check our website at http://www.gromacs.org
32  * 
33  * And Hey:
34  * GROwing Monsters And Cloning Shrimps
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #ifdef GMX_CRAY_XT3
41 #include<catamount/dclock.h>
42 #endif
43
44
45 #include <stdio.h>
46 #include <time.h>
47 #ifdef HAVE_SYS_TIME_H
48 #include <sys/time.h>
49 #endif
50 #include <math.h>
51 #include "typedefs.h"
52 #include "string2.h"
53 #include "gmxfio.h"
54 #include "smalloc.h"
55 #include "names.h"
56 #include "confio.h"
57 #include "mvdata.h"
58 #include "txtdump.h"
59 #include "pbc.h"
60 #include "chargegroup.h"
61 #include "vec.h"
62 #include <time.h>
63 #include "nrnb.h"
64 #include "mshift.h"
65 #include "mdrun.h"
66 #include "update.h"
67 #include "physics.h"
68 #include "main.h"
69 #include "mdatoms.h"
70 #include "force.h"
71 #include "bondf.h"
72 #include "pme.h"
73 #include "disre.h"
74 #include "orires.h"
75 #include "network.h"
76 #include "calcmu.h"
77 #include "constr.h"
78 #include "xvgr.h"
79 #include "trnio.h"
80 #include "xtcio.h"
81 #include "copyrite.h"
82 #include "pull_rotation.h"
83 #include "mpelogging.h"
84 #include "domdec.h"
85 #include "partdec.h"
86 #include "gmx_wallcycle.h"
87 #include "genborn.h"
88
89 #ifdef GMX_LIB_MPI
90 #include <mpi.h>
91 #endif
92 #ifdef GMX_THREAD_MPI
93 #include "tmpi.h"
94 #endif
95
96 #include "adress.h"
97 #include "qmmm.h"
98
99 #if 0
100 typedef struct gmx_timeprint {
101     
102 } t_gmx_timeprint;
103 #endif
104
105 /* Portable version of ctime_r implemented in src/gmxlib/string2.c, but we do not want it declared in public installed headers */
106 char *
107 gmx_ctime_r(const time_t *clock,char *buf, int n);
108
109
110 double
111 gmx_gettime()
112 {
113 #ifdef HAVE_GETTIMEOFDAY
114         struct timeval t;
115         double seconds;
116         
117         gettimeofday(&t,NULL);
118         
119         seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
120         
121         return seconds;
122 #else
123         double  seconds;
124         
125         seconds = time(NULL);
126         
127         return seconds;
128 #endif
129 }
130
131
132 #define difftime(end,start) ((double)(end)-(double)(start))
133
134 void print_time(FILE *out,gmx_runtime_t *runtime,gmx_large_int_t step,   
135                 t_inputrec *ir, t_commrec *cr)
136 {
137     time_t finish;
138     char   timebuf[STRLEN];
139     double dt;
140     char buf[48];
141     
142 #ifndef GMX_THREAD_MPI
143     if (!PAR(cr))
144 #endif
145     {
146         fprintf(out,"\r");
147     }
148     fprintf(out,"step %s",gmx_step_str(step,buf));
149     if ((step >= ir->nstlist))
150     {
151         if ((ir->nstlist == 0) || ((step % ir->nstlist) == 0))
152         {
153             /* We have done a full cycle let's update time_per_step */
154             runtime->last = gmx_gettime();
155             dt = difftime(runtime->last,runtime->real);
156             runtime->time_per_step = dt/(step - ir->init_step + 1);
157         }
158         dt = (ir->nsteps + ir->init_step - step)*runtime->time_per_step;
159         
160         if (ir->nsteps >= 0)
161         {
162             if (dt >= 300)
163             {    
164                 finish = (time_t) (runtime->last + dt);
165                 gmx_ctime_r(&finish,timebuf,STRLEN);
166                 sprintf(buf,"%s",timebuf);
167                 buf[strlen(buf)-1]='\0';
168                 fprintf(out,", will finish %s",buf);
169             }
170             else
171                 fprintf(out,", remaining runtime: %5d s          ",(int)dt);
172         }
173         else
174         {
175             fprintf(out," performance: %.1f ns/day    ",
176                     ir->delta_t/1000*24*60*60/runtime->time_per_step);
177         }
178     }
179 #ifndef GMX_THREAD_MPI
180     if (PAR(cr))
181     {
182         fprintf(out,"\n");
183     }
184 #endif
185
186     fflush(out);
187 }
188
189 #ifdef NO_CLOCK 
190 #define clock() -1
191 #endif
192
193 static double set_proctime(gmx_runtime_t *runtime)
194 {
195     double diff;
196 #ifdef GMX_CRAY_XT3
197     double prev;
198
199     prev = runtime->proc;
200     runtime->proc = dclock();
201     
202     diff = runtime->proc - prev;
203 #else
204     clock_t prev;
205
206     prev = runtime->proc;
207     runtime->proc = clock();
208
209     diff = (double)(runtime->proc - prev)/(double)CLOCKS_PER_SEC;
210 #endif
211     if (diff < 0)
212     {
213         /* The counter has probably looped, ignore this data */
214         diff = 0;
215     }
216
217     return diff;
218 }
219
220 void runtime_start(gmx_runtime_t *runtime)
221 {
222     runtime->real = gmx_gettime();
223     runtime->proc          = 0;
224     set_proctime(runtime);
225     runtime->realtime      = 0;
226     runtime->proctime      = 0;
227     runtime->last          = 0;
228     runtime->time_per_step = 0;
229 }
230
231 void runtime_end(gmx_runtime_t *runtime)
232 {
233     double now;
234     
235     now = gmx_gettime();
236     
237     runtime->proctime += set_proctime(runtime);
238     runtime->realtime  = now - runtime->real;
239     runtime->real      = now;
240 }
241
242 void runtime_upd_proc(gmx_runtime_t *runtime)
243 {
244     runtime->proctime += set_proctime(runtime);
245 }
246
247 void print_date_and_time(FILE *fplog,int nodeid,const char *title,
248                          const gmx_runtime_t *runtime)
249 {
250     int i;
251     char timebuf[STRLEN];
252     char time_string[STRLEN];
253     time_t tmptime;
254
255     if (fplog)
256     {
257         if (runtime != NULL)
258         {
259             tmptime = (time_t) runtime->real;
260             gmx_ctime_r(&tmptime,timebuf,STRLEN);
261         }
262         else
263         {
264             tmptime = (time_t) gmx_gettime();
265             gmx_ctime_r(&tmptime,timebuf,STRLEN);
266         }
267         for(i=0; timebuf[i]>=' '; i++)
268         {
269             time_string[i]=timebuf[i];
270         }
271         time_string[i]='\0';
272
273         fprintf(fplog,"%s on node %d %s\n",title,nodeid,time_string);
274     }
275 }
276
277 static void sum_forces(int start,int end,rvec f[],rvec flr[])
278 {
279   int i;
280   
281   if (gmx_debug_at) {
282     pr_rvecs(debug,0,"fsr",f+start,end-start);
283     pr_rvecs(debug,0,"flr",flr+start,end-start);
284   }
285   for(i=start; (i<end); i++)
286     rvec_inc(f[i],flr[i]);
287 }
288
289 /* 
290  * calc_f_el calculates forces due to an electric field.
291  *
292  * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e 
293  *
294  * Et[] contains the parameters for the time dependent 
295  * part of the field (not yet used). 
296  * Ex[] contains the parameters for
297  * the spatial dependent part of the field. You can have cool periodic
298  * fields in principle, but only a constant field is supported
299  * now. 
300  * The function should return the energy due to the electric field
301  * (if any) but for now returns 0.
302  *
303  * WARNING:
304  * There can be problems with the virial.
305  * Since the field is not self-consistent this is unavoidable.
306  * For neutral molecules the virial is correct within this approximation.
307  * For neutral systems with many charged molecules the error is small.
308  * But for systems with a net charge or a few charged molecules
309  * the error can be significant when the field is high.
310  * Solution: implement a self-consitent electric field into PME.
311  */
312 static void calc_f_el(FILE *fp,int  start,int homenr,
313                       real charge[],rvec x[],rvec f[],
314                       t_cosines Ex[],t_cosines Et[],double t)
315 {
316     rvec Ext;
317     real t0;
318     int  i,m;
319     
320     for(m=0; (m<DIM); m++)
321     {
322         if (Et[m].n > 0)
323         {
324             if (Et[m].n == 3)
325             {
326                 t0 = Et[m].a[1];
327                 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
328             }
329             else
330             {
331                 Ext[m] = cos(Et[m].a[0]*t);
332             }
333         }
334         else
335         {
336             Ext[m] = 1.0;
337         }
338         if (Ex[m].n > 0)
339         {
340             /* Convert the field strength from V/nm to MD-units */
341             Ext[m] *= Ex[m].a[0]*FIELDFAC;
342             for(i=start; (i<start+homenr); i++)
343                 f[i][m] += charge[i]*Ext[m];
344         }
345         else
346         {
347             Ext[m] = 0;
348         }
349     }
350     if (fp != NULL)
351     {
352         fprintf(fp,"%10g  %10g  %10g  %10g #FIELD\n",t,
353                 Ext[XX]/FIELDFAC,Ext[YY]/FIELDFAC,Ext[ZZ]/FIELDFAC);
354     }
355 }
356
357 static void calc_virial(FILE *fplog,int start,int homenr,rvec x[],rvec f[],
358                         tensor vir_part,t_graph *graph,matrix box,
359                         t_nrnb *nrnb,const t_forcerec *fr,int ePBC)
360 {
361   int i,j;
362   tensor virtest;
363
364   /* The short-range virial from surrounding boxes */
365   clear_mat(vir_part);
366   calc_vir(fplog,SHIFTS,fr->shift_vec,fr->fshift,vir_part,ePBC==epbcSCREW,box);
367   inc_nrnb(nrnb,eNR_VIRIAL,SHIFTS);
368   
369   /* Calculate partial virial, for local atoms only, based on short range. 
370    * Total virial is computed in global_stat, called from do_md 
371    */
372   f_calc_vir(fplog,start,start+homenr,x,f,vir_part,graph,box);
373   inc_nrnb(nrnb,eNR_VIRIAL,homenr);
374
375   /* Add position restraint contribution */
376   for(i=0; i<DIM; i++) {
377     vir_part[i][i] += fr->vir_diag_posres[i];
378   }
379
380   /* Add wall contribution */
381   for(i=0; i<DIM; i++) {
382     vir_part[i][ZZ] += fr->vir_wall_z[i];
383   }
384
385   if (debug)
386     pr_rvecs(debug,0,"vir_part",vir_part,DIM);
387 }
388
389 static void print_large_forces(FILE *fp,t_mdatoms *md,t_commrec *cr,
390                                gmx_large_int_t step,real pforce,rvec *x,rvec *f)
391 {
392   int  i;
393   real pf2,fn2;
394   char buf[STEPSTRSIZE];
395
396   pf2 = sqr(pforce);
397   for(i=md->start; i<md->start+md->homenr; i++) {
398     fn2 = norm2(f[i]);
399     /* We also catch NAN, if the compiler does not optimize this away. */
400     if (fn2 >= pf2 || fn2 != fn2) {
401       fprintf(fp,"step %s  atom %6d  x %8.3f %8.3f %8.3f  force %12.5e\n",
402               gmx_step_str(step,buf),
403               ddglatnr(cr->dd,i),x[i][XX],x[i][YY],x[i][ZZ],sqrt(fn2));
404     }
405   }
406 }
407
408 void do_force(FILE *fplog,t_commrec *cr,
409               t_inputrec *inputrec,
410               gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
411               gmx_localtop_t *top,
412               gmx_mtop_t *mtop,
413               gmx_groups_t *groups,
414               matrix box,rvec x[],history_t *hist,
415               rvec f[],
416               tensor vir_force,
417               t_mdatoms *mdatoms,
418               gmx_enerdata_t *enerd,t_fcdata *fcd,
419               real lambda,t_graph *graph,
420               t_forcerec *fr,gmx_vsite_t *vsite,rvec mu_tot,
421               double t,FILE *field,gmx_edsam_t ed,
422               gmx_bool bBornRadii,
423               int flags)
424 {
425     int    cg0,cg1,i,j;
426     int    start,homenr;
427     double mu[2*DIM]; 
428     gmx_bool   bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
429     gmx_bool   bDoLongRange,bDoForces,bSepLRF;
430     gmx_bool   bDoAdressWF;
431     matrix boxs;
432     real   e,v,dvdl;
433     t_pbc  pbc;
434     float  cycles_ppdpme,cycles_pme,cycles_seppme,cycles_force;
435   
436     start  = mdatoms->start;
437     homenr = mdatoms->homenr;
438
439     bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
440
441     clear_mat(vir_force);
442
443     if (PARTDECOMP(cr))
444     {
445         pd_cg_range(cr,&cg0,&cg1);
446     }
447     else
448     {
449         cg0 = 0;
450         if (DOMAINDECOMP(cr))
451         {
452             cg1 = cr->dd->ncg_tot;
453         }
454         else
455         {
456             cg1 = top->cgs.nr;
457         }
458         if (fr->n_tpi > 0)
459         {
460             cg1--;
461         }
462     }
463
464     bStateChanged = (flags & GMX_FORCE_STATECHANGED);
465     bNS           = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE); 
466     bFillGrid     = (bNS && bStateChanged);
467     bCalcCGCM     = (bFillGrid && !DOMAINDECOMP(cr));
468     bDoLongRange  = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DOLR));
469     bDoForces     = (flags & GMX_FORCE_FORCES);
470     bSepLRF       = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
471     /* should probably move this to the forcerec since it doesn't change */
472     bDoAdressWF   = ((fr->adress_type!=eAdressOff));
473
474     if (bStateChanged)
475     {
476         update_forcerec(fplog,fr,box);
477         
478         /* Calculate total (local) dipole moment in a temporary common array. 
479          * This makes it possible to sum them over nodes faster.
480          */
481         calc_mu(start,homenr,
482                 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
483                 mu,mu+DIM);
484     }
485   
486   if (fr->ePBC != epbcNONE) { 
487     /* Compute shift vectors every step,
488      * because of pressure coupling or box deformation!
489      */
490     if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
491       calc_shifts(box,fr->shift_vec);
492     
493     if (bCalcCGCM) { 
494       put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,box,
495                                &(top->cgs),x,fr->cg_cm);
496       inc_nrnb(nrnb,eNR_CGCM,homenr);
497       inc_nrnb(nrnb,eNR_RESETX,cg1-cg0);
498     } 
499     else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
500       unshift_self(graph,box,x);
501     }
502   } 
503   else if (bCalcCGCM) {
504     calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm);
505     inc_nrnb(nrnb,eNR_CGCM,homenr);
506   }
507   
508   if (bCalcCGCM) {
509     if (PAR(cr)) {
510       move_cgcm(fplog,cr,fr->cg_cm);
511     }
512     if (gmx_debug_at)
513       pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr);
514   }
515
516 #ifdef GMX_MPI
517   if (!(cr->duty & DUTY_PME)) {
518     /* Send particle coordinates to the pme nodes.
519      * Since this is only implemented for domain decomposition
520      * and domain decomposition does not use the graph,
521      * we do not need to worry about shifting.
522      */    
523
524     wallcycle_start(wcycle,ewcPP_PMESENDX);
525     GMX_MPE_LOG(ev_send_coordinates_start);
526
527     bBS = (inputrec->nwall == 2);
528     if (bBS) {
529       copy_mat(box,boxs);
530       svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
531     }
532
533     gmx_pme_send_x(cr,bBS ? boxs : box,x,
534                    mdatoms->nChargePerturbed,lambda,
535                    ( flags & GMX_FORCE_VIRIAL),step);
536
537     GMX_MPE_LOG(ev_send_coordinates_finish);
538     wallcycle_stop(wcycle,ewcPP_PMESENDX);
539   }
540 #endif /* GMX_MPI */
541
542     /* Communicate coordinates and sum dipole if necessary */
543     if (PAR(cr))
544     {
545         wallcycle_start(wcycle,ewcMOVEX);
546         if (DOMAINDECOMP(cr))
547         {
548             dd_move_x(cr->dd,box,x);
549         }
550         else
551         {
552             move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb);
553         }
554         /* When we don't need the total dipole we sum it in global_stat */
555         if (bStateChanged && NEED_MUTOT(*inputrec))
556         {
557             gmx_sumd(2*DIM,mu,cr);
558         }
559         wallcycle_stop(wcycle,ewcMOVEX);
560     }
561     if (bStateChanged)
562     {
563
564         /* update adress weight beforehand */
565         if(bDoAdressWF)
566         {
567             /* need pbc for adress weight calculation with pbc_dx */
568             set_pbc(&pbc,inputrec->ePBC,box);
569             if(fr->adress_site == eAdressSITEcog)
570             {
571                 update_adress_weights_cog(top->idef.iparams,top->idef.il,x,fr,mdatoms,
572                                           inputrec->ePBC==epbcNONE ? NULL : &pbc);
573             }
574             else if (fr->adress_site == eAdressSITEcom)
575             {
576                 update_adress_weights_com(fplog,cg0,cg1,&(top->cgs),x,fr,mdatoms,
577                                           inputrec->ePBC==epbcNONE ? NULL : &pbc);
578             }
579             else if (fr->adress_site == eAdressSITEatomatom){
580                 update_adress_weights_atom_per_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
581                                           inputrec->ePBC==epbcNONE ? NULL : &pbc);
582             }
583             else
584             {
585                 update_adress_weights_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
586                                            inputrec->ePBC==epbcNONE ? NULL : &pbc);
587             }
588         }
589
590         for(i=0; i<2; i++)
591         {
592             for(j=0;j<DIM;j++)
593             {
594                 fr->mu_tot[i][j] = mu[i*DIM + j];
595             }
596         }
597     }
598     if (fr->efep == efepNO)
599     {
600         copy_rvec(fr->mu_tot[0],mu_tot);
601     }
602     else
603     {
604         for(j=0; j<DIM; j++)
605         {
606             mu_tot[j] =
607                 (1.0 - lambda)*fr->mu_tot[0][j] + lambda*fr->mu_tot[1][j];
608         }
609     }
610
611     /* Reset energies */
612     reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
613     clear_rvecs(SHIFTS,fr->fshift);
614
615     if (bNS)
616     {
617         wallcycle_start(wcycle,ewcNS);
618         
619         if (graph && bStateChanged)
620         {
621             /* Calculate intramolecular shift vectors to make molecules whole */
622             mk_mshift(fplog,graph,fr->ePBC,box,x);
623         }
624
625         /* Reset long range forces if necessary */
626         if (fr->bTwinRange)
627         {
628             /* Reset the (long-range) forces if necessary */
629             clear_rvecs(fr->natoms_force_constr,bSepLRF ? fr->f_twin : f);
630         }
631
632         /* Do the actual neighbour searching and if twin range electrostatics
633          * also do the calculation of long range forces and energies.
634          */
635         dvdl = 0; 
636         ns(fplog,fr,x,box,
637            groups,&(inputrec->opts),top,mdatoms,
638            cr,nrnb,lambda,&dvdl,&enerd->grpp,bFillGrid,
639            bDoLongRange,bDoForces,bSepLRF ? fr->f_twin : f);
640         if (bSepDVDL)
641         {
642             fprintf(fplog,sepdvdlformat,"LR non-bonded",0.0,dvdl);
643         }
644         enerd->dvdl_lin += dvdl;
645         
646         wallcycle_stop(wcycle,ewcNS);
647     }
648         
649     if (inputrec->implicit_solvent && bNS) 
650     {
651         make_gb_nblist(cr,inputrec->gb_algorithm,inputrec->rlist,
652                        x,box,fr,&top->idef,graph,fr->born);
653     }
654         
655     if (DOMAINDECOMP(cr))
656     {
657         if (!(cr->duty & DUTY_PME))
658         {
659             wallcycle_start(wcycle,ewcPPDURINGPME);
660             dd_force_flop_start(cr->dd,nrnb);
661         }
662     }
663     
664     if (inputrec->bRot)
665     {
666         /* Enforced rotation has its own cycle counter that starts after the collective
667          * coordinates have been communicated. It is added to ddCyclF to allow
668          * for proper load-balancing */
669         wallcycle_start(wcycle,ewcROT);
670         do_rotation(cr,inputrec,box,x,t,step,wcycle,bNS);
671         wallcycle_stop(wcycle,ewcROT);
672     }
673
674     /* Start the force cycle counter.
675      * This counter is stopped in do_forcelow_level.
676      * No parallel communication should occur while this counter is running,
677      * since that will interfere with the dynamic load balancing.
678      */
679     wallcycle_start(wcycle,ewcFORCE);
680
681     if (bDoForces)
682     {
683         /* Reset forces for which the virial is calculated separately:
684          * PME/Ewald forces if necessary */
685         if (fr->bF_NoVirSum) 
686         {
687             if (flags & GMX_FORCE_VIRIAL)
688             {
689                 fr->f_novirsum = fr->f_novirsum_alloc;
690                 GMX_BARRIER(cr->mpi_comm_mygroup);
691                 if (fr->bDomDec)
692                 {
693                     clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
694                 }
695                 else
696                 {
697                     clear_rvecs(homenr,fr->f_novirsum+start);
698                 }
699                 GMX_BARRIER(cr->mpi_comm_mygroup);
700             }
701             else
702             {
703                 /* We are not calculating the pressure so we do not need
704                  * a separate array for forces that do not contribute
705                  * to the pressure.
706                  */
707                 fr->f_novirsum = f;
708             }
709         }
710
711         if (bSepLRF)
712         {
713             /* Add the long range forces to the short range forces */
714             for(i=0; i<fr->natoms_force_constr; i++)
715             {
716                 copy_rvec(fr->f_twin[i],f[i]);
717             }
718         }
719         else if (!(fr->bTwinRange && bNS))
720         {
721             /* Clear the short-range forces */
722             clear_rvecs(fr->natoms_force_constr,f);
723         }
724
725         clear_rvec(fr->vir_diag_posres);
726
727         GMX_BARRIER(cr->mpi_comm_mygroup);
728     }
729     if (inputrec->ePull == epullCONSTRAINT)
730     {
731         clear_pull_forces(inputrec->pull);
732     }
733
734     /* update QMMMrec, if necessary */
735     if(fr->bQMMM)
736     {
737         update_QMMMrec(cr,fr,x,mdatoms,box,top);
738     }
739
740     if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
741     {
742         /* Position restraints always require full pbc. Check if we already did it for Adress */
743         if(!(bStateChanged && bDoAdressWF))
744         {
745             set_pbc(&pbc,inputrec->ePBC,box);
746         }
747         v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
748                    top->idef.iparams_posres,
749                    (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres,
750                    inputrec->ePBC==epbcNONE ? NULL : &pbc,lambda,&dvdl,
751                    fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
752         if (bSepDVDL)
753         {
754             fprintf(fplog,sepdvdlformat,
755                     interaction_function[F_POSRES].longname,v,dvdl);
756         }
757         enerd->term[F_POSRES] += v;
758         /* This linear lambda dependence assumption is only correct
759          * when only k depends on lambda,
760          * not when the reference position depends on lambda.
761          * grompp checks for this.
762          */
763         enerd->dvdl_lin += dvdl;
764         inc_nrnb(nrnb,eNR_POSRES,top->idef.il[F_POSRES].nr/2);
765     }
766
767     /* Compute the bonded and non-bonded energies and optionally forces */    
768     do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
769                       cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
770                       x,hist,f,enerd,fcd,mtop,top,fr->born,
771                       &(top->atomtypes),bBornRadii,box,
772                       lambda,graph,&(top->excls),fr->mu_tot,
773                       flags,&cycles_pme);
774     
775     cycles_force = wallcycle_stop(wcycle,ewcFORCE);
776     GMX_BARRIER(cr->mpi_comm_mygroup);
777     
778     if (ed)
779     {
780         do_flood(fplog,cr,x,f,ed,box,step);
781     }
782         
783     if (DOMAINDECOMP(cr))
784     {
785         dd_force_flop_stop(cr->dd,nrnb);
786         if (wcycle)
787         {
788             dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
789         }
790     }
791     
792     if (bDoForces)
793     {
794         if (IR_ELEC_FIELD(*inputrec))
795         {
796             /* Compute forces due to electric field */
797             calc_f_el(MASTER(cr) ? field : NULL,
798                       start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
799                       inputrec->ex,inputrec->et,t);
800         }
801
802         if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
803         {
804             /* Compute thermodynamic force in hybrid AdResS region */
805             adress_thermo_force(start,homenr,&(top->cgs),x,fr->f_novirsum,fr,mdatoms,
806                                 inputrec->ePBC==epbcNONE ? NULL : &pbc);
807         }
808         
809         /* Communicate the forces */
810         if (PAR(cr))
811         {
812             wallcycle_start(wcycle,ewcMOVEF);
813             if (DOMAINDECOMP(cr))
814             {
815                 dd_move_f(cr->dd,f,fr->fshift);
816                 /* Do we need to communicate the separate force array
817                  * for terms that do not contribute to the single sum virial?
818                  * Position restraints and electric fields do not introduce
819                  * inter-cg forces, only full electrostatics methods do.
820                  * When we do not calculate the virial, fr->f_novirsum = f,
821                  * so we have already communicated these forces.
822                  */
823                 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
824                     (flags & GMX_FORCE_VIRIAL))
825                 {
826                     dd_move_f(cr->dd,fr->f_novirsum,NULL);
827                 }
828                 if (bSepLRF)
829                 {
830                     /* We should not update the shift forces here,
831                      * since f_twin is already included in f.
832                      */
833                     dd_move_f(cr->dd,fr->f_twin,NULL);
834                 }
835             }
836             else
837             {
838                 pd_move_f(cr,f,nrnb);
839                 if (bSepLRF)
840                 {
841                     pd_move_f(cr,fr->f_twin,nrnb);
842                 }
843             }
844             wallcycle_stop(wcycle,ewcMOVEF);
845         }
846
847         /* If we have NoVirSum forces, but we do not calculate the virial,
848          * we sum fr->f_novirum=f later.
849          */
850         if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
851         {
852             wallcycle_start(wcycle,ewcVSITESPREAD);
853             spread_vsite_f(fplog,vsite,x,f,fr->fshift,FALSE,NULL,nrnb,
854                            &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
855             wallcycle_stop(wcycle,ewcVSITESPREAD);
856
857             if (bSepLRF)
858             {
859                 wallcycle_start(wcycle,ewcVSITESPREAD);
860                 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,FALSE,NULL,
861                                nrnb,
862                                &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
863                 wallcycle_stop(wcycle,ewcVSITESPREAD);
864             }
865         }
866         
867         if (flags & GMX_FORCE_VIRIAL)
868         {
869             /* Calculation of the virial must be done after vsites! */
870             calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
871                         vir_force,graph,box,nrnb,fr,inputrec->ePBC);
872         }
873     }
874
875     enerd->term[F_COM_PULL] = 0;
876     if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
877     {
878         /* Calculate the center of mass forces, this requires communication,
879          * which is why pull_potential is called close to other communication.
880          * The virial contribution is calculated directly,
881          * which is why we call pull_potential after calc_virial.
882          */
883         set_pbc(&pbc,inputrec->ePBC,box);
884         dvdl = 0; 
885         enerd->term[F_COM_PULL] +=
886             pull_potential(inputrec->ePull,inputrec->pull,mdatoms,&pbc,
887                            cr,t,lambda,x,f,vir_force,&dvdl);
888         if (bSepDVDL)
889         {
890             fprintf(fplog,sepdvdlformat,"Com pull",enerd->term[F_COM_PULL],dvdl);
891         }
892         enerd->dvdl_lin += dvdl;
893     }
894     
895     /* Add the forces from enforced rotation potentials (if any) */
896     if (inputrec->bRot)
897     {
898         wallcycle_start(wcycle,ewcROTadd);
899         enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr,step,t);
900         wallcycle_stop(wcycle,ewcROTadd);
901     }
902
903     if (PAR(cr) && !(cr->duty & DUTY_PME))
904     {
905         cycles_ppdpme = wallcycle_stop(wcycle,ewcPPDURINGPME);
906         dd_cycles_add(cr->dd,cycles_ppdpme,ddCyclPPduringPME);
907
908         /* In case of node-splitting, the PP nodes receive the long-range 
909          * forces, virial and energy from the PME nodes here.
910          */    
911         wallcycle_start(wcycle,ewcPP_PMEWAITRECVF);
912         dvdl = 0;
913         gmx_pme_receive_f(cr,fr->f_novirsum,fr->vir_el_recip,&e,&dvdl,
914                           &cycles_seppme);
915         if (bSepDVDL)
916         {
917             fprintf(fplog,sepdvdlformat,"PME mesh",e,dvdl);
918         }
919         enerd->term[F_COUL_RECIP] += e;
920         enerd->dvdl_lin += dvdl;
921         if (wcycle)
922         {
923             dd_cycles_add(cr->dd,cycles_seppme,ddCyclPME);
924         }
925         wallcycle_stop(wcycle,ewcPP_PMEWAITRECVF);
926     }
927
928     if (bDoForces && fr->bF_NoVirSum)
929     {
930         if (vsite)
931         {
932             /* Spread the mesh force on virtual sites to the other particles... 
933              * This is parallellized. MPI communication is performed
934              * if the constructing atoms aren't local.
935              */
936             wallcycle_start(wcycle,ewcVSITESPREAD);
937             spread_vsite_f(fplog,vsite,x,fr->f_novirsum,NULL,
938                            (flags & GMX_FORCE_VIRIAL),fr->vir_el_recip,
939                            nrnb,
940                            &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
941             wallcycle_stop(wcycle,ewcVSITESPREAD);
942         }
943         if (flags & GMX_FORCE_VIRIAL)
944         {
945             /* Now add the forces, this is local */
946             if (fr->bDomDec)
947             {
948                 sum_forces(0,fr->f_novirsum_n,f,fr->f_novirsum);
949             }
950             else
951             {
952                 sum_forces(start,start+homenr,f,fr->f_novirsum);
953             }
954             if (EEL_FULL(fr->eeltype))
955             {
956                 /* Add the mesh contribution to the virial */
957                 m_add(vir_force,fr->vir_el_recip,vir_force);
958             }
959             if (debug)
960             {
961                 pr_rvecs(debug,0,"vir_force",vir_force,DIM);
962             }
963         }
964     }
965     
966     /* Sum the potential energy terms from group contributions */
967     sum_epot(&(inputrec->opts),enerd);
968     
969     if (fr->print_force >= 0 && bDoForces)
970     {
971         print_large_forces(stderr,mdatoms,cr,step,fr->print_force,x,f);
972     }
973 }
974
975 void do_constrain_first(FILE *fplog,gmx_constr_t constr,
976                         t_inputrec *ir,t_mdatoms *md,
977                         t_state *state,rvec *f,
978                         t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
979                         t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
980 {
981     int    i,m,start,end;
982     gmx_large_int_t step;
983     real   dt=ir->delta_t;
984     real   dvdlambda;
985     rvec   *savex;
986     
987     snew(savex,state->natoms);
988
989     start = md->start;
990     end   = md->homenr + start;
991     
992     if (debug)
993         fprintf(debug,"vcm: start=%d, homenr=%d, end=%d\n",
994                 start,md->homenr,end);
995     /* Do a first constrain to reset particles... */
996     step = ir->init_step;
997     if (fplog)
998     {
999         char buf[STEPSTRSIZE];
1000         fprintf(fplog,"\nConstraining the starting coordinates (step %s)\n",
1001                 gmx_step_str(step,buf));
1002     }
1003     dvdlambda = 0;
1004     
1005     /* constrain the current position */
1006     constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1007               ir,NULL,cr,step,0,md,
1008               state->x,state->x,NULL,
1009               state->box,state->lambda,&dvdlambda,
1010               NULL,NULL,nrnb,econqCoord,ir->epc==epcMTTK,state->veta,state->veta);
1011     if (EI_VV(ir->eI)) 
1012     {
1013         /* constrain the inital velocity, and save it */
1014         /* also may be useful if we need the ekin from the halfstep for velocity verlet */
1015         /* might not yet treat veta correctly */
1016         constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1017                   ir,NULL,cr,step,0,md,
1018                   state->x,state->v,state->v,
1019                   state->box,state->lambda,&dvdlambda,
1020                   NULL,NULL,nrnb,econqVeloc,ir->epc==epcMTTK,state->veta,state->veta);
1021     }
1022     /* constrain the inital velocities at t-dt/2 */
1023     if (EI_STATE_VELOCITY(ir->eI) && ir->eI!=eiVV)
1024     {
1025         for(i=start; (i<end); i++) 
1026         {
1027             for(m=0; (m<DIM); m++) 
1028             {
1029                 /* Reverse the velocity */
1030                 state->v[i][m] = -state->v[i][m];
1031                 /* Store the position at t-dt in buf */
1032                 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
1033             }
1034         }
1035     /* Shake the positions at t=-dt with the positions at t=0                        
1036      * as reference coordinates.                                                     
1037          */
1038         if (fplog)
1039         {
1040             char buf[STEPSTRSIZE];
1041             fprintf(fplog,"\nConstraining the coordinates at t0-dt (step %s)\n",
1042                     gmx_step_str(step,buf));
1043         }
1044         dvdlambda = 0;
1045         constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1046                   ir,NULL,cr,step,-1,md,
1047                   state->x,savex,NULL,
1048                   state->box,state->lambda,&dvdlambda,
1049                   state->v,NULL,nrnb,econqCoord,ir->epc==epcMTTK,state->veta,state->veta);
1050         
1051         for(i=start; i<end; i++) {
1052             for(m=0; m<DIM; m++) {
1053                 /* Re-reverse the velocities */
1054                 state->v[i][m] = -state->v[i][m];
1055             }
1056         }
1057     }
1058     
1059     sfree(savex);
1060 }
1061
1062 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr)
1063 {
1064   double eners[2],virs[2],enersum,virsum,y0,f,g,h;
1065   double r0,r1,r,rc3,rc9,ea,eb,ec,pa,pb,pc,pd;
1066   double invscale,invscale2,invscale3;
1067   int    ri0,ri1,ri,i,offstart,offset;
1068   real   scale,*vdwtab; 
1069
1070   fr->enershiftsix = 0;
1071   fr->enershifttwelve = 0;
1072   fr->enerdiffsix = 0;
1073   fr->enerdifftwelve = 0;
1074   fr->virdiffsix = 0;
1075   fr->virdifftwelve = 0;
1076
1077   if (eDispCorr != edispcNO) {
1078     for(i=0; i<2; i++) {
1079       eners[i] = 0;
1080       virs[i]  = 0;
1081     }
1082     if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT)) {
1083       if (fr->rvdw_switch == 0)
1084         gmx_fatal(FARGS,
1085                   "With dispersion correction rvdw-switch can not be zero "
1086                   "for vdw-type = %s",evdw_names[fr->vdwtype]);
1087
1088       scale  = fr->nblists[0].tab.scale;
1089       vdwtab = fr->nblists[0].vdwtab;
1090
1091       /* Round the cut-offs to exact table values for precision */
1092       ri0 = floor(fr->rvdw_switch*scale);
1093       ri1 = ceil(fr->rvdw*scale);
1094       r0  = ri0/scale;
1095       r1  = ri1/scale;
1096       rc3 = r0*r0*r0;
1097       rc9  = rc3*rc3*rc3;
1098
1099       if (fr->vdwtype == evdwSHIFT) {
1100         /* Determine the constant energy shift below rvdw_switch */
1101         fr->enershiftsix    = (real)(-1.0/(rc3*rc3)) - vdwtab[8*ri0];
1102         fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - vdwtab[8*ri0 + 4];
1103       }
1104       /* Add the constant part from 0 to rvdw_switch.
1105        * This integration from 0 to rvdw_switch overcounts the number
1106        * of interactions by 1, as it also counts the self interaction.
1107        * We will correct for this later.
1108        */
1109       eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
1110       eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
1111       
1112       invscale = 1.0/(scale);  
1113       invscale2 = invscale*invscale;
1114       invscale3 = invscale*invscale2;
1115
1116       /* following summation derived from cubic spline definition,
1117         Numerical Recipies in C, second edition, p. 113-116.  Exact
1118         for the cubic spline.  We first calculate the negative of
1119         the energy from rvdw to rvdw_switch, assuming that g(r)=1,
1120         and then add the more standard, abrupt cutoff correction to
1121         that result, yielding the long-range correction for a
1122         switched function.  We perform both the pressure and energy
1123         loops at the same time for simplicity, as the computational
1124         cost is low. */
1125       
1126       for (i=0;i<2;i++) {
1127         enersum = 0.0; virsum = 0.0;
1128         if (i==0)
1129           offstart = 0;
1130         else
1131           offstart = 4;
1132         for (ri=ri0; ri<ri1; ri++) {
1133           r = ri*invscale;
1134           ea = invscale3;
1135           eb = 2.0*invscale2*r;
1136           ec = invscale*r*r;
1137           
1138           pa = invscale3;
1139           pb = 3.0*invscale2*r;
1140           pc = 3.0*invscale*r*r;
1141           pd = r*r*r;
1142           
1143           /* this "8" is from the packing in the vdwtab array - perhaps
1144             should be #define'ed? */
1145           offset = 8*ri + offstart;
1146           y0 = vdwtab[offset];
1147           f = vdwtab[offset+1];
1148           g = vdwtab[offset+2];
1149           h = vdwtab[offset+3];
1150           
1151           enersum += y0*(ea/3 + eb/2 + ec) + f*(ea/4 + eb/3 + ec/2)+
1152             g*(ea/5 + eb/4 + ec/3) + h*(ea/6 + eb/5 + ec/4);  
1153           virsum  +=  f*(pa/4 + pb/3 + pc/2 + pd) + 
1154             2*g*(pa/5 + pb/4 + pc/3 + pd/2) + 3*h*(pa/6 + pb/5 + pc/4 + pd/3);
1155           
1156         }
1157         enersum *= 4.0*M_PI;
1158         virsum  *= 4.0*M_PI; 
1159         eners[i] -= enersum;
1160         virs[i]  -= virsum;
1161       }
1162
1163       /* now add the correction for rvdw_switch to infinity */
1164       eners[0] += -4.0*M_PI/(3.0*rc3);
1165       eners[1] +=  4.0*M_PI/(9.0*rc9);
1166       virs[0]  +=  8.0*M_PI/rc3;
1167       virs[1]  += -16.0*M_PI/(3.0*rc9);
1168     } 
1169     else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER)) {
1170       if (fr->vdwtype == evdwUSER && fplog)
1171         fprintf(fplog,
1172                 "WARNING: using dispersion correction with user tables\n");
1173       rc3  = fr->rvdw*fr->rvdw*fr->rvdw;
1174       rc9  = rc3*rc3*rc3;
1175       eners[0] += -4.0*M_PI/(3.0*rc3);
1176       eners[1] +=  4.0*M_PI/(9.0*rc9);
1177       virs[0]  +=  8.0*M_PI/rc3;
1178       virs[1]  += -16.0*M_PI/(3.0*rc9);
1179     } else {
1180       gmx_fatal(FARGS,
1181                 "Dispersion correction is not implemented for vdw-type = %s",
1182                 evdw_names[fr->vdwtype]);
1183     }
1184     fr->enerdiffsix    = eners[0];
1185     fr->enerdifftwelve = eners[1];
1186     /* The 0.5 is due to the Gromacs definition of the virial */
1187     fr->virdiffsix     = 0.5*virs[0];
1188     fr->virdifftwelve  = 0.5*virs[1];
1189   }
1190 }
1191
1192 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
1193                    gmx_large_int_t step,int natoms,
1194                    matrix box,real lambda,tensor pres,tensor virial,
1195                    real *prescorr, real *enercorr, real *dvdlcorr)
1196 {
1197     gmx_bool bCorrAll,bCorrPres;
1198     real dvdlambda,invvol,dens,ninter,avcsix,avctwelve,enerdiff,svir=0,spres=0;
1199     int  m;
1200     
1201     *prescorr = 0;
1202     *enercorr = 0;
1203     *dvdlcorr = 0;
1204     
1205     clear_mat(virial);
1206     clear_mat(pres);
1207     
1208     if (ir->eDispCorr != edispcNO) {
1209         bCorrAll  = (ir->eDispCorr == edispcAllEner ||
1210                      ir->eDispCorr == edispcAllEnerPres);
1211         bCorrPres = (ir->eDispCorr == edispcEnerPres ||
1212                      ir->eDispCorr == edispcAllEnerPres);
1213         
1214         invvol = 1/det(box);
1215         if (fr->n_tpi) 
1216         {
1217             /* Only correct for the interactions with the inserted molecule */
1218             dens = (natoms - fr->n_tpi)*invvol;
1219             ninter = fr->n_tpi;
1220         } 
1221         else 
1222         {
1223             dens = natoms*invvol;
1224             ninter = 0.5*natoms;
1225         }
1226         
1227         if (ir->efep == efepNO) 
1228         {
1229             avcsix    = fr->avcsix[0];
1230             avctwelve = fr->avctwelve[0];
1231         } 
1232         else 
1233         {
1234             avcsix    = (1 - lambda)*fr->avcsix[0]    + lambda*fr->avcsix[1];
1235             avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
1236         }
1237         
1238         enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
1239         *enercorr += avcsix*enerdiff;
1240         dvdlambda = 0.0;
1241         if (ir->efep != efepNO) 
1242         {
1243             dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
1244         }
1245         if (bCorrAll) 
1246         {
1247             enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
1248             *enercorr += avctwelve*enerdiff;
1249             if (fr->efep != efepNO) 
1250             {
1251                 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
1252             }
1253         }
1254         
1255         if (bCorrPres) 
1256         {
1257             svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
1258             if (ir->eDispCorr == edispcAllEnerPres)
1259             {
1260                 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
1261             }
1262             /* The factor 2 is because of the Gromacs virial definition */
1263             spres = -2.0*invvol*svir*PRESFAC;
1264             
1265             for(m=0; m<DIM; m++) {
1266                 virial[m][m] += svir;
1267                 pres[m][m] += spres;
1268             }
1269             *prescorr += spres;
1270         }
1271         
1272         /* Can't currently control when it prints, for now, just print when degugging */
1273         if (debug)
1274         {
1275             if (bCorrAll) {
1276                 fprintf(debug,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1277                         avcsix,avctwelve);
1278             }
1279             if (bCorrPres) 
1280             {
1281                 fprintf(debug,
1282                         "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
1283                         *enercorr,spres,svir);
1284             }
1285             else
1286             {
1287                 fprintf(debug,"Long Range LJ corr.: Epot %10g\n",*enercorr);
1288             }
1289         }
1290         
1291         if (fr->bSepDVDL && do_per_step(step,ir->nstlog))
1292         {
1293             fprintf(fplog,sepdvdlformat,"Dispersion correction",
1294                     *enercorr,dvdlambda);
1295         }
1296         if (fr->efep != efepNO) 
1297         {
1298             *dvdlcorr += dvdlambda;
1299         }
1300     }
1301 }
1302
1303 void do_pbc_first(FILE *fplog,matrix box,t_forcerec *fr,
1304                   t_graph *graph,rvec x[])
1305 {
1306   if (fplog)
1307     fprintf(fplog,"Removing pbc first time\n");
1308   calc_shifts(box,fr->shift_vec);
1309   if (graph) {
1310     mk_mshift(fplog,graph,fr->ePBC,box,x);
1311     if (gmx_debug_at)
1312       p_graph(debug,"do_pbc_first 1",graph);
1313     shift_self(graph,box,x);
1314     /* By doing an extra mk_mshift the molecules that are broken
1315      * because they were e.g. imported from another software
1316      * will be made whole again. Such are the healing powers
1317      * of GROMACS.
1318      */
1319     mk_mshift(fplog,graph,fr->ePBC,box,x);
1320     if (gmx_debug_at)
1321       p_graph(debug,"do_pbc_first 2",graph);
1322   }
1323   if (fplog)
1324     fprintf(fplog,"Done rmpbc\n");
1325 }
1326
1327 static void low_do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1328                             gmx_mtop_t *mtop,rvec x[],
1329                             gmx_bool bFirst)
1330 {
1331   t_graph *graph;
1332   int mb,as,mol;
1333   gmx_molblock_t *molb;
1334
1335   if (bFirst && fplog)
1336     fprintf(fplog,"Removing pbc first time\n");
1337
1338   snew(graph,1);
1339   as = 0;
1340   for(mb=0; mb<mtop->nmolblock; mb++) {
1341     molb = &mtop->molblock[mb];
1342     if (molb->natoms_mol == 1 || 
1343         (!bFirst && mtop->moltype[molb->type].cgs.nr == 1)) {
1344       /* Just one atom or charge group in the molecule, no PBC required */
1345       as += molb->nmol*molb->natoms_mol;
1346     } else {
1347       /* Pass NULL iso fplog to avoid graph prints for each molecule type */
1348       mk_graph_ilist(NULL,mtop->moltype[molb->type].ilist,
1349                      0,molb->natoms_mol,FALSE,FALSE,graph);
1350       
1351       for(mol=0; mol<molb->nmol; mol++) {
1352         mk_mshift(fplog,graph,ePBC,box,x+as);
1353         
1354         shift_self(graph,box,x+as);
1355         /* The molecule is whole now.
1356          * We don't need the second mk_mshift call as in do_pbc_first,
1357          * since we no longer need this graph.
1358          */
1359         
1360         as += molb->natoms_mol;
1361       }
1362       done_graph(graph);
1363     }
1364   }
1365   sfree(graph);
1366 }
1367
1368 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
1369                        gmx_mtop_t *mtop,rvec x[])
1370 {
1371   low_do_pbc_mtop(fplog,ePBC,box,mtop,x,TRUE);
1372 }
1373
1374 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1375                  gmx_mtop_t *mtop,rvec x[])
1376 {
1377   low_do_pbc_mtop(fplog,ePBC,box,mtop,x,FALSE);
1378 }
1379
1380 void finish_run(FILE *fplog,t_commrec *cr,const char *confout,
1381                 t_inputrec *inputrec,
1382                 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
1383                 gmx_runtime_t *runtime,
1384                 gmx_bool bWriteStat)
1385 {
1386   int    i,j;
1387   t_nrnb *nrnb_tot=NULL;
1388   real   delta_t;
1389   double nbfs,mflop;
1390   double cycles[ewcNR];
1391
1392   wallcycle_sum(cr,wcycle,cycles);
1393
1394   if (cr->nnodes > 1) {
1395     if (SIMMASTER(cr))
1396       snew(nrnb_tot,1);
1397 #ifdef GMX_MPI
1398     MPI_Reduce(nrnb->n,nrnb_tot->n,eNRNB,MPI_DOUBLE,MPI_SUM,
1399                MASTERRANK(cr),cr->mpi_comm_mysim);
1400 #endif  
1401   } else {
1402     nrnb_tot = nrnb;
1403   }
1404     
1405   if (SIMMASTER(cr)) {
1406     print_flop(fplog,nrnb_tot,&nbfs,&mflop);
1407     if (cr->nnodes > 1) {
1408       sfree(nrnb_tot);
1409     }
1410   }
1411
1412   if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr)) {
1413     print_dd_statistics(cr,inputrec,fplog);
1414   }
1415
1416 #ifdef GMX_MPI
1417     if (PARTDECOMP(cr))
1418     {
1419         if (MASTER(cr))
1420         {
1421             t_nrnb     *nrnb_all;
1422             int        s;
1423             MPI_Status stat;
1424
1425             snew(nrnb_all,cr->nnodes);
1426             nrnb_all[0] = *nrnb;
1427             for(s=1; s<cr->nnodes; s++)
1428             {
1429                 MPI_Recv(nrnb_all[s].n,eNRNB,MPI_DOUBLE,s,0,
1430                          cr->mpi_comm_mysim,&stat);
1431             }
1432             pr_load(fplog,cr,nrnb_all);
1433             sfree(nrnb_all);
1434         }
1435         else
1436         {
1437             MPI_Send(nrnb->n,eNRNB,MPI_DOUBLE,MASTERRANK(cr),0,
1438                      cr->mpi_comm_mysim);
1439         }
1440     }
1441 #endif  
1442
1443   if (SIMMASTER(cr)) {
1444     wallcycle_print(fplog,cr->nnodes,cr->npmenodes,runtime->realtime,
1445                     wcycle,cycles);
1446
1447     if (EI_DYNAMICS(inputrec->eI)) {
1448       delta_t = inputrec->delta_t;
1449     } else {
1450       delta_t = 0;
1451     }
1452     
1453     if (fplog) {
1454         print_perf(fplog,runtime->proctime,runtime->realtime,
1455                    cr->nnodes-cr->npmenodes,
1456                    runtime->nsteps_done,delta_t,nbfs,mflop);
1457     }
1458     if (bWriteStat) {
1459         print_perf(stderr,runtime->proctime,runtime->realtime,
1460                    cr->nnodes-cr->npmenodes,
1461                    runtime->nsteps_done,delta_t,nbfs,mflop);
1462     }
1463
1464     /*
1465     runtime=inputrec->nsteps*inputrec->delta_t;
1466     if (bWriteStat) {
1467       if (cr->nnodes == 1)
1468         fprintf(stderr,"\n\n");
1469       print_perf(stderr,nodetime,realtime,runtime,&ntot,
1470                  cr->nnodes-cr->npmenodes,FALSE);
1471     }
1472     wallcycle_print(fplog,cr->nnodes,cr->npmenodes,realtime,wcycle,cycles);
1473     print_perf(fplog,nodetime,realtime,runtime,&ntot,cr->nnodes-cr->npmenodes,
1474                TRUE);
1475     if (PARTDECOMP(cr))
1476       pr_load(fplog,cr,nrnb_all);
1477     if (cr->nnodes > 1)
1478       sfree(nrnb_all);
1479     */
1480   }
1481 }
1482
1483 void init_md(FILE *fplog,
1484              t_commrec *cr,t_inputrec *ir,const output_env_t oenv,
1485              double *t,double *t0,
1486              real *lambda,double *lam0,
1487              t_nrnb *nrnb,gmx_mtop_t *mtop,
1488              gmx_update_t *upd,
1489              int nfile,const t_filenm fnm[],
1490              gmx_mdoutf_t **outf,t_mdebin **mdebin,
1491              tensor force_vir,tensor shake_vir,rvec mu_tot,
1492              gmx_bool *bSimAnn,t_vcm **vcm, t_state *state, unsigned long Flags)
1493 {
1494     int  i,j,n;
1495     real tmpt,mod;
1496         
1497     /* Initial values */
1498     *t = *t0       = ir->init_t;
1499     if (ir->efep != efepNO)
1500     {
1501         *lam0 = ir->init_lambda;
1502         *lambda = *lam0 + ir->init_step*ir->delta_lambda;
1503     }
1504     else
1505     {
1506         *lambda = *lam0   = 0.0;
1507     } 
1508
1509     *bSimAnn=FALSE;
1510     for(i=0;i<ir->opts.ngtc;i++)
1511     {
1512         /* set bSimAnn if any group is being annealed */
1513         if(ir->opts.annealing[i]!=eannNO)
1514         {
1515             *bSimAnn = TRUE;
1516         }
1517     }
1518     if (*bSimAnn)
1519     {
1520         update_annealing_target_temp(&(ir->opts),ir->init_t);
1521     }
1522     
1523     if (upd)
1524     {
1525         *upd = init_update(fplog,ir);
1526     }
1527     
1528     if (vcm != NULL)
1529     {
1530         *vcm = init_vcm(fplog,&mtop->groups,ir);
1531     }
1532     
1533     if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
1534     {
1535         if (ir->etc == etcBERENDSEN)
1536         {
1537             please_cite(fplog,"Berendsen84a");
1538         }
1539         if (ir->etc == etcVRESCALE)
1540         {
1541             please_cite(fplog,"Bussi2007a");
1542         }
1543     }
1544     
1545     init_nrnb(nrnb);
1546     
1547     if (nfile != -1)
1548     {
1549         *outf = init_mdoutf(nfile,fnm,Flags,cr,ir,oenv);
1550
1551         *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
1552                               mtop,ir, (*outf)->fp_dhdl);
1553     }
1554     
1555     if (ir->bAdress)
1556     {
1557       please_cite(fplog,"Fritsch12");
1558       please_cite(fplog,"Junghans10");
1559     }
1560     /* Initiate variables */  
1561     clear_mat(force_vir);
1562     clear_mat(shake_vir);
1563     clear_rvec(mu_tot);
1564     
1565     debug_gmx();
1566 }
1567