Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / 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 "gmx_random.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,dvdlambda[efptNR];
433     real   dvdl_dum,lambda_dum;
434     t_pbc  pbc;
435     float  cycles_ppdpme,cycles_pme,cycles_seppme,cycles_force;
436
437     start  = mdatoms->start;
438     homenr = mdatoms->homenr;
439
440     bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
441
442     clear_mat(vir_force);
443
444     if (PARTDECOMP(cr))
445     {
446         pd_cg_range(cr,&cg0,&cg1);
447     }
448     else
449     {
450         cg0 = 0;
451         if (DOMAINDECOMP(cr))
452         {
453             cg1 = cr->dd->ncg_tot;
454         }
455         else
456         {
457             cg1 = top->cgs.nr;
458         }
459         if (fr->n_tpi > 0)
460         {
461             cg1--;
462         }
463     }
464
465     bStateChanged = (flags & GMX_FORCE_STATECHANGED);
466     bNS           = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE);
467     bFillGrid     = (bNS && bStateChanged);
468     bCalcCGCM     = (bFillGrid && !DOMAINDECOMP(cr));
469     bDoLongRange  = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DOLR));
470     bDoForces     = (flags & GMX_FORCE_FORCES);
471     bSepLRF       = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
472     /* should probably move this to the forcerec since it doesn't change */
473     bDoAdressWF   = ((fr->adress_type!=eAdressOff));
474
475     if (bStateChanged)
476     {
477         update_forcerec(fplog,fr,box);
478
479         /* Calculate total (local) dipole moment in a temporary common array.
480          * This makes it possible to sum them over nodes faster.
481          */
482         calc_mu(start,homenr,
483                 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
484                 mu,mu+DIM);
485     }
486
487   if (fr->ePBC != epbcNONE) {
488     /* Compute shift vectors every step,
489      * because of pressure coupling or box deformation!
490      */
491     if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
492       calc_shifts(box,fr->shift_vec);
493
494     if (bCalcCGCM) {
495       put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,box,
496                                &(top->cgs),x,fr->cg_cm);
497       inc_nrnb(nrnb,eNR_CGCM,homenr);
498       inc_nrnb(nrnb,eNR_RESETX,cg1-cg0);
499     }
500     else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
501       unshift_self(graph,box,x);
502     }
503   }
504   else if (bCalcCGCM) {
505     calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm);
506     inc_nrnb(nrnb,eNR_CGCM,homenr);
507   }
508
509   if (bCalcCGCM) {
510     if (PAR(cr)) {
511       move_cgcm(fplog,cr,fr->cg_cm);
512     }
513     if (gmx_debug_at)
514       pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr);
515   }
516
517 #ifdef GMX_MPI
518   if (!(cr->duty & DUTY_PME)) {
519     /* Send particle coordinates to the pme nodes.
520      * Since this is only implemented for domain decomposition
521      * and domain decomposition does not use the graph,
522      * we do not need to worry about shifting.
523      */
524
525     wallcycle_start(wcycle,ewcPP_PMESENDX);
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[efptCOUL],
535                    ( flags & GMX_FORCE_VIRIAL),step);
536
537     wallcycle_stop(wcycle,ewcPP_PMESENDX);
538   }
539 #endif /* GMX_MPI */
540
541     /* Communicate coordinates and sum dipole if necessary */
542     if (PAR(cr))
543     {
544         wallcycle_start(wcycle,ewcMOVEX);
545         if (DOMAINDECOMP(cr))
546         {
547             dd_move_x(cr->dd,box,x);
548         }
549         else
550         {
551             move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb);
552         }
553         /* When we don't need the total dipole we sum it in global_stat */
554         if (bStateChanged && NEED_MUTOT(*inputrec))
555         {
556             gmx_sumd(2*DIM,mu,cr);
557         }
558         wallcycle_stop(wcycle,ewcMOVEX);
559     }
560     if (bStateChanged)
561     {
562
563         /* update adress weight beforehand */
564         if(bDoAdressWF)
565         {
566             /* need pbc for adress weight calculation with pbc_dx */
567             set_pbc(&pbc,inputrec->ePBC,box);
568             if(fr->adress_site == eAdressSITEcog)
569             {
570                 update_adress_weights_cog(top->idef.iparams,top->idef.il,x,fr,mdatoms,
571                                           inputrec->ePBC==epbcNONE ? NULL : &pbc);
572             }
573             else if (fr->adress_site == eAdressSITEcom)
574             {
575                 update_adress_weights_com(fplog,cg0,cg1,&(top->cgs),x,fr,mdatoms,
576                                           inputrec->ePBC==epbcNONE ? NULL : &pbc);
577             }
578             else if (fr->adress_site == eAdressSITEatomatom){
579                 update_adress_weights_atom_per_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
580                                           inputrec->ePBC==epbcNONE ? NULL : &pbc);
581             }
582             else
583             {
584                 update_adress_weights_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
585                                            inputrec->ePBC==epbcNONE ? NULL : &pbc);
586             }
587         }
588
589         for(i=0; i<2; i++)
590         {
591             for(j=0;j<DIM;j++)
592             {
593                 fr->mu_tot[i][j] = mu[i*DIM + j];
594             }
595         }
596     }
597     if (fr->efep == efepNO)
598     {
599         copy_rvec(fr->mu_tot[0],mu_tot);
600     }
601     else
602     {
603         for(j=0; j<DIM; j++)
604         {
605             mu_tot[j] =
606                 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
607         }
608     }
609
610     /* Reset energies */
611     reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
612     clear_rvecs(SHIFTS,fr->fshift);
613
614     if (bNS)
615     {
616         wallcycle_start(wcycle,ewcNS);
617
618         if (graph && bStateChanged)
619         {
620             /* Calculate intramolecular shift vectors to make molecules whole */
621             mk_mshift(fplog,graph,fr->ePBC,box,x);
622         }
623
624         /* Reset long range forces if necessary */
625         if (fr->bTwinRange)
626         {
627             /* Reset the (long-range) forces if necessary */
628             clear_rvecs(fr->natoms_force_constr,bSepLRF ? fr->f_twin : f);
629         }
630
631         /* Do the actual neighbour searching and if twin range electrostatics
632          * also do the calculation of long range forces and energies.
633          */
634         for (i=0;i<efptNR;i++) {dvdlambda[i] = 0;}
635         ns(fplog,fr,x,box,
636            groups,&(inputrec->opts),top,mdatoms,
637            cr,nrnb,lambda,dvdlambda,&enerd->grpp,bFillGrid,
638            bDoLongRange,bDoForces,bSepLRF ? fr->f_twin : f);
639         if (bSepDVDL)
640         {
641             fprintf(fplog,sepdvdlformat,"LR non-bonded",0.0,dvdlambda);
642         }
643         enerd->dvdl_lin[efptVDW] += dvdlambda[efptVDW];
644         enerd->dvdl_lin[efptCOUL] += dvdlambda[efptCOUL];
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                 if (fr->bDomDec)
691                 {
692                     clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
693                 }
694                 else
695                 {
696                     clear_rvecs(homenr,fr->f_novirsum+start);
697                 }
698             }
699             else
700             {
701                 /* We are not calculating the pressure so we do not need
702                  * a separate array for forces that do not contribute
703                  * to the pressure.
704                  */
705                 fr->f_novirsum = f;
706             }
707         }
708
709         if (bSepLRF)
710         {
711             /* Add the long range forces to the short range forces */
712             for(i=0; i<fr->natoms_force_constr; i++)
713             {
714                 copy_rvec(fr->f_twin[i],f[i]);
715             }
716         }
717         else if (!(fr->bTwinRange && bNS))
718         {
719             /* Clear the short-range forces */
720             clear_rvecs(fr->natoms_force_constr,f);
721         }
722
723         clear_rvec(fr->vir_diag_posres);
724     }
725     if (inputrec->ePull == epullCONSTRAINT)
726     {
727         clear_pull_forces(inputrec->pull);
728     }
729
730     /* update QMMMrec, if necessary */
731     if(fr->bQMMM)
732     {
733         update_QMMMrec(cr,fr,x,mdatoms,box,top);
734     }
735
736     if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
737     {
738         /* Position restraints always require full pbc. Check if we already did it for Adress */
739         if(!(bStateChanged && bDoAdressWF))
740         {
741             set_pbc(&pbc,inputrec->ePBC,box);
742         }
743         v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
744                    top->idef.iparams_posres,
745                    (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres,
746                    inputrec->ePBC==epbcNONE ? NULL : &pbc,lambda[efptRESTRAINT],&(dvdlambda[efptRESTRAINT]),
747                    fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
748         if (bSepDVDL)
749         {
750             fprintf(fplog,sepdvdlformat,
751                     interaction_function[F_POSRES].longname,v,dvdlambda);
752         }
753         enerd->term[F_POSRES] += v;
754         /* This linear lambda dependence assumption is only correct
755          * when only k depends on lambda,
756          * not when the reference position depends on lambda.
757          * grompp checks for this.  (verify this is still the case?)
758          */
759         enerd->dvdl_nonlin[efptRESTRAINT] += dvdlambda[efptRESTRAINT]; /* if just the force constant changes, this is linear,
760                                                                           but we can't be sure w/o additional checking that is
761                                                                           hard to do at this level of code. Otherwise,
762                                                                           the dvdl is not differentiable */
763         inc_nrnb(nrnb,eNR_POSRES,top->idef.il[F_POSRES].nr/2);
764         if ((inputrec->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
765         {
766             for(i=0; i<enerd->n_lambda; i++)
767             {
768                 lambda_dum = (i==0 ? lambda[efptRESTRAINT] : inputrec->fepvals->all_lambda[efptRESTRAINT][i-1]);
769                 v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
770                            top->idef.iparams_posres,
771                            (const rvec*)x,NULL,NULL,
772                            inputrec->ePBC==epbcNONE ? NULL : &pbc,lambda_dum,&dvdl_dum,
773                            fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
774                 enerd->enerpart_lambda[i] += v;
775             }
776         }
777    }
778
779     if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
780     {
781         /* Flat-bottomed position restraints always require full pbc */
782         if(!(bStateChanged && bDoAdressWF))
783         {
784             set_pbc(&pbc,inputrec->ePBC,box);
785         }
786         v = fbposres(top->idef.il[F_FBPOSRES].nr,top->idef.il[F_FBPOSRES].iatoms,
787                      top->idef.iparams_fbposres,
788                      (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres,
789                      inputrec->ePBC==epbcNONE ? NULL : &pbc,
790                      fr->rc_scaling,fr->ePBC,fr->posres_com);
791         enerd->term[F_FBPOSRES] += v;
792         inc_nrnb(nrnb,eNR_FBPOSRES,top->idef.il[F_FBPOSRES].nr/2);
793     }
794
795     /* Compute the bonded and non-bonded energies and optionally forces */
796     do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
797                       cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
798                       x,hist,f,enerd,fcd,mtop,top,fr->born,
799                       &(top->atomtypes),bBornRadii,box,
800                       inputrec->fepvals,lambda,graph,&(top->excls),fr->mu_tot,
801                       flags,&cycles_pme);
802
803     cycles_force = wallcycle_stop(wcycle,ewcFORCE);
804
805     if (ed)
806     {
807         do_flood(fplog,cr,x,f,ed,box,step,bNS);
808     }
809
810     if (DOMAINDECOMP(cr))
811     {
812         dd_force_flop_stop(cr->dd,nrnb);
813         if (wcycle)
814         {
815             dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
816         }
817     }
818
819     if (bDoForces)
820     {
821         if (IR_ELEC_FIELD(*inputrec))
822         {
823             /* Compute forces due to electric field */
824             calc_f_el(MASTER(cr) ? field : NULL,
825                       start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
826                       inputrec->ex,inputrec->et,t);
827         }
828
829         if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
830         {
831             /* Compute thermodynamic force in hybrid AdResS region */
832             adress_thermo_force(start,homenr,&(top->cgs),x,fr->f_novirsum,fr,mdatoms,
833                                 inputrec->ePBC==epbcNONE ? NULL : &pbc);
834         }
835
836         /* Communicate the forces */
837         if (PAR(cr))
838         {
839             wallcycle_start(wcycle,ewcMOVEF);
840             if (DOMAINDECOMP(cr))
841             {
842                 dd_move_f(cr->dd,f,fr->fshift);
843                 /* Do we need to communicate the separate force array
844                  * for terms that do not contribute to the single sum virial?
845                  * Position restraints and electric fields do not introduce
846                  * inter-cg forces, only full electrostatics methods do.
847                  * When we do not calculate the virial, fr->f_novirsum = f,
848                  * so we have already communicated these forces.
849                  */
850                 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
851                     (flags & GMX_FORCE_VIRIAL))
852                 {
853                     dd_move_f(cr->dd,fr->f_novirsum,NULL);
854                 }
855                 if (bSepLRF)
856                 {
857                     /* We should not update the shift forces here,
858                      * since f_twin is already included in f.
859                      */
860                     dd_move_f(cr->dd,fr->f_twin,NULL);
861                 }
862             }
863             else
864             {
865                 pd_move_f(cr,f,nrnb);
866                 if (bSepLRF)
867                 {
868                     pd_move_f(cr,fr->f_twin,nrnb);
869                 }
870             }
871             wallcycle_stop(wcycle,ewcMOVEF);
872         }
873
874         /* If we have NoVirSum forces, but we do not calculate the virial,
875          * we sum fr->f_novirum=f later.
876          */
877         if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
878         {
879             wallcycle_start(wcycle,ewcVSITESPREAD);
880             spread_vsite_f(fplog,vsite,x,f,fr->fshift,FALSE,NULL,nrnb,
881                            &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
882             wallcycle_stop(wcycle,ewcVSITESPREAD);
883
884             if (bSepLRF)
885             {
886                 wallcycle_start(wcycle,ewcVSITESPREAD);
887                 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,FALSE,NULL,
888                                nrnb,
889                                &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
890                 wallcycle_stop(wcycle,ewcVSITESPREAD);
891             }
892         }
893
894         if (flags & GMX_FORCE_VIRIAL)
895         {
896             /* Calculation of the virial must be done after vsites! */
897             calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
898                         vir_force,graph,box,nrnb,fr,inputrec->ePBC);
899         }
900     }
901
902     enerd->term[F_COM_PULL] = 0;
903     if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
904     {
905         /* Calculate the center of mass forces, this requires communication,
906          * which is why pull_potential is called close to other communication.
907          * The virial contribution is calculated directly,
908          * which is why we call pull_potential after calc_virial.
909          */
910         set_pbc(&pbc,inputrec->ePBC,box);
911         dvdlambda[efptRESTRAINT] = 0;
912         enerd->term[F_COM_PULL] +=
913             pull_potential(inputrec->ePull,inputrec->pull,mdatoms,&pbc,
914                            cr,t,lambda[efptRESTRAINT],x,f,vir_force,&(dvdlambda[efptRESTRAINT]));
915         if (bSepDVDL)
916         {
917             fprintf(fplog,sepdvdlformat,"Com pull",enerd->term[F_COM_PULL],dvdlambda[efptRESTRAINT]);
918         }
919         enerd->dvdl_lin[efptRESTRAINT] += dvdlambda[efptRESTRAINT];
920     }
921
922     /* Add the forces from enforced rotation potentials (if any) */
923     if (inputrec->bRot)
924     {
925         wallcycle_start(wcycle,ewcROTadd);
926         enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr,step,t);
927         wallcycle_stop(wcycle,ewcROTadd);
928     }
929
930     if (PAR(cr) && !(cr->duty & DUTY_PME))
931     {
932         cycles_ppdpme = wallcycle_stop(wcycle,ewcPPDURINGPME);
933         dd_cycles_add(cr->dd,cycles_ppdpme,ddCyclPPduringPME);
934
935         /* In case of node-splitting, the PP nodes receive the long-range
936          * forces, virial and energy from the PME nodes here.
937          */
938         wallcycle_start(wcycle,ewcPP_PMEWAITRECVF);
939         dvdlambda[efptCOUL] = 0;
940         gmx_pme_receive_f(cr,fr->f_novirsum,fr->vir_el_recip,&e,&dvdlambda[efptCOUL],
941                           &cycles_seppme);
942         if (bSepDVDL)
943         {
944             fprintf(fplog,sepdvdlformat,"PME mesh",e,dvdlambda[efptCOUL]);
945         }
946         enerd->term[F_COUL_RECIP] += e;
947         enerd->dvdl_lin[efptCOUL] += dvdlambda[efptCOUL];
948         if (wcycle)
949         {
950             dd_cycles_add(cr->dd,cycles_seppme,ddCyclPME);
951         }
952         wallcycle_stop(wcycle,ewcPP_PMEWAITRECVF);
953     }
954
955     if (bDoForces && fr->bF_NoVirSum)
956     {
957         if (vsite)
958         {
959             /* Spread the mesh force on virtual sites to the other particles...
960              * This is parallellized. MPI communication is performed
961              * if the constructing atoms aren't local.
962              */
963             wallcycle_start(wcycle,ewcVSITESPREAD);
964             spread_vsite_f(fplog,vsite,x,fr->f_novirsum,NULL,
965                            (flags & GMX_FORCE_VIRIAL),fr->vir_el_recip,
966                            nrnb,
967                            &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
968             wallcycle_stop(wcycle,ewcVSITESPREAD);
969         }
970         if (flags & GMX_FORCE_VIRIAL)
971         {
972             /* Now add the forces, this is local */
973             if (fr->bDomDec)
974             {
975                 sum_forces(0,fr->f_novirsum_n,f,fr->f_novirsum);
976             }
977             else
978             {
979                 sum_forces(start,start+homenr,f,fr->f_novirsum);
980             }
981             if (EEL_FULL(fr->eeltype))
982             {
983                 /* Add the mesh contribution to the virial */
984                 m_add(vir_force,fr->vir_el_recip,vir_force);
985             }
986             if (debug)
987             {
988                 pr_rvecs(debug,0,"vir_force",vir_force,DIM);
989             }
990         }
991     }
992
993     /* Sum the potential energy terms from group contributions */
994     sum_epot(&(inputrec->opts),enerd);
995
996     if (fr->print_force >= 0 && bDoForces)
997     {
998         print_large_forces(stderr,mdatoms,cr,step,fr->print_force,x,f);
999     }
1000 }
1001
1002 void do_constrain_first(FILE *fplog,gmx_constr_t constr,
1003                         t_inputrec *ir,t_mdatoms *md,
1004                         t_state *state,rvec *f,
1005                         t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
1006                         t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
1007 {
1008     int    i,m,start,end;
1009     gmx_large_int_t step;
1010     real   dt=ir->delta_t;
1011     real   dvdl_dum;
1012     rvec   *savex;
1013
1014     snew(savex,state->natoms);
1015
1016     start = md->start;
1017     end   = md->homenr + start;
1018
1019     if (debug)
1020         fprintf(debug,"vcm: start=%d, homenr=%d, end=%d\n",
1021                 start,md->homenr,end);
1022     /* Do a first constrain to reset particles... */
1023     step = ir->init_step;
1024     if (fplog)
1025     {
1026         char buf[STEPSTRSIZE];
1027         fprintf(fplog,"\nConstraining the starting coordinates (step %s)\n",
1028                 gmx_step_str(step,buf));
1029     }
1030     dvdl_dum = 0;
1031
1032     /* constrain the current position */
1033     constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1034               ir,NULL,cr,step,0,md,
1035               state->x,state->x,NULL,
1036               state->box,state->lambda[efptBONDED],&dvdl_dum,
1037               NULL,NULL,nrnb,econqCoord,ir->epc==epcMTTK,state->veta,state->veta);
1038     if (EI_VV(ir->eI))
1039     {
1040         /* constrain the inital velocity, and save it */
1041         /* also may be useful if we need the ekin from the halfstep for velocity verlet */
1042         /* might not yet treat veta correctly */
1043         constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1044                   ir,NULL,cr,step,0,md,
1045                   state->x,state->v,state->v,
1046                   state->box,state->lambda[efptBONDED],&dvdl_dum,
1047                   NULL,NULL,nrnb,econqVeloc,ir->epc==epcMTTK,state->veta,state->veta);
1048     }
1049     /* constrain the inital velocities at t-dt/2 */
1050     if (EI_STATE_VELOCITY(ir->eI) && ir->eI!=eiVV)
1051     {
1052         for(i=start; (i<end); i++)
1053         {
1054             for(m=0; (m<DIM); m++)
1055             {
1056                 /* Reverse the velocity */
1057                 state->v[i][m] = -state->v[i][m];
1058                 /* Store the position at t-dt in buf */
1059                 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
1060             }
1061         }
1062     /* Shake the positions at t=-dt with the positions at t=0
1063      * as reference coordinates.
1064          */
1065         if (fplog)
1066         {
1067             char buf[STEPSTRSIZE];
1068             fprintf(fplog,"\nConstraining the coordinates at t0-dt (step %s)\n",
1069                     gmx_step_str(step,buf));
1070         }
1071         dvdl_dum = 0;
1072         constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1073                   ir,NULL,cr,step,-1,md,
1074                   state->x,savex,NULL,
1075                   state->box,state->lambda[efptBONDED],&dvdl_dum,
1076                   state->v,NULL,nrnb,econqCoord,ir->epc==epcMTTK,state->veta,state->veta);
1077
1078         for(i=start; i<end; i++) {
1079             for(m=0; m<DIM; m++) {
1080                 /* Re-reverse the velocities */
1081                 state->v[i][m] = -state->v[i][m];
1082             }
1083         }
1084     }
1085     sfree(savex);
1086 }
1087
1088 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr)
1089 {
1090   double eners[2],virs[2],enersum,virsum,y0,f,g,h;
1091   double r0,r1,r,rc3,rc9,ea,eb,ec,pa,pb,pc,pd;
1092   double invscale,invscale2,invscale3;
1093   int    ri0,ri1,ri,i,offstart,offset;
1094   real   scale,*vdwtab;
1095
1096   fr->enershiftsix = 0;
1097   fr->enershifttwelve = 0;
1098   fr->enerdiffsix = 0;
1099   fr->enerdifftwelve = 0;
1100   fr->virdiffsix = 0;
1101   fr->virdifftwelve = 0;
1102
1103   if (eDispCorr != edispcNO) {
1104     for(i=0; i<2; i++) {
1105       eners[i] = 0;
1106       virs[i]  = 0;
1107     }
1108     if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT)) {
1109       if (fr->rvdw_switch == 0)
1110         gmx_fatal(FARGS,
1111                   "With dispersion correction rvdw-switch can not be zero "
1112                   "for vdw-type = %s",evdw_names[fr->vdwtype]);
1113
1114       scale  = fr->nblists[0].tab.scale;
1115       vdwtab = fr->nblists[0].vdwtab;
1116
1117       /* Round the cut-offs to exact table values for precision */
1118       ri0 = floor(fr->rvdw_switch*scale);
1119       ri1 = ceil(fr->rvdw*scale);
1120       r0  = ri0/scale;
1121       r1  = ri1/scale;
1122       rc3 = r0*r0*r0;
1123       rc9  = rc3*rc3*rc3;
1124
1125       if (fr->vdwtype == evdwSHIFT) {
1126         /* Determine the constant energy shift below rvdw_switch */
1127         fr->enershiftsix    = (real)(-1.0/(rc3*rc3)) - vdwtab[8*ri0];
1128         fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - vdwtab[8*ri0 + 4];
1129       }
1130       /* Add the constant part from 0 to rvdw_switch.
1131        * This integration from 0 to rvdw_switch overcounts the number
1132        * of interactions by 1, as it also counts the self interaction.
1133        * We will correct for this later.
1134        */
1135       eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
1136       eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
1137
1138       invscale = 1.0/(scale);
1139       invscale2 = invscale*invscale;
1140       invscale3 = invscale*invscale2;
1141
1142       /* following summation derived from cubic spline definition,
1143         Numerical Recipies in C, second edition, p. 113-116.  Exact
1144         for the cubic spline.  We first calculate the negative of
1145         the energy from rvdw to rvdw_switch, assuming that g(r)=1,
1146         and then add the more standard, abrupt cutoff correction to
1147         that result, yielding the long-range correction for a
1148         switched function.  We perform both the pressure and energy
1149         loops at the same time for simplicity, as the computational
1150         cost is low. */
1151
1152       for (i=0;i<2;i++) {
1153         enersum = 0.0; virsum = 0.0;
1154         if (i==0)
1155           offstart = 0;
1156         else
1157           offstart = 4;
1158         for (ri=ri0; ri<ri1; ri++) {
1159           r = ri*invscale;
1160           ea = invscale3;
1161           eb = 2.0*invscale2*r;
1162           ec = invscale*r*r;
1163
1164           pa = invscale3;
1165           pb = 3.0*invscale2*r;
1166           pc = 3.0*invscale*r*r;
1167           pd = r*r*r;
1168
1169           /* this "8" is from the packing in the vdwtab array - perhaps
1170             should be #define'ed? */
1171           offset = 8*ri + offstart;
1172           y0 = vdwtab[offset];
1173           f = vdwtab[offset+1];
1174           g = vdwtab[offset+2];
1175           h = vdwtab[offset+3];
1176
1177           enersum += y0*(ea/3 + eb/2 + ec) + f*(ea/4 + eb/3 + ec/2)+
1178             g*(ea/5 + eb/4 + ec/3) + h*(ea/6 + eb/5 + ec/4);
1179           virsum  +=  f*(pa/4 + pb/3 + pc/2 + pd) +
1180             2*g*(pa/5 + pb/4 + pc/3 + pd/2) + 3*h*(pa/6 + pb/5 + pc/4 + pd/3);
1181
1182         }
1183         enersum *= 4.0*M_PI;
1184         virsum  *= 4.0*M_PI;
1185         eners[i] -= enersum;
1186         virs[i]  -= virsum;
1187       }
1188
1189       /* now add the correction for rvdw_switch to infinity */
1190       eners[0] += -4.0*M_PI/(3.0*rc3);
1191       eners[1] +=  4.0*M_PI/(9.0*rc9);
1192       virs[0]  +=  8.0*M_PI/rc3;
1193       virs[1]  += -16.0*M_PI/(3.0*rc9);
1194     }
1195     else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER)) {
1196       if (fr->vdwtype == evdwUSER && fplog)
1197         fprintf(fplog,
1198                 "WARNING: using dispersion correction with user tables\n");
1199       rc3  = fr->rvdw*fr->rvdw*fr->rvdw;
1200       rc9  = rc3*rc3*rc3;
1201       eners[0] += -4.0*M_PI/(3.0*rc3);
1202       eners[1] +=  4.0*M_PI/(9.0*rc9);
1203       virs[0]  +=  8.0*M_PI/rc3;
1204       virs[1]  += -16.0*M_PI/(3.0*rc9);
1205     } else {
1206       gmx_fatal(FARGS,
1207                 "Dispersion correction is not implemented for vdw-type = %s",
1208                 evdw_names[fr->vdwtype]);
1209     }
1210     fr->enerdiffsix    = eners[0];
1211     fr->enerdifftwelve = eners[1];
1212     /* The 0.5 is due to the Gromacs definition of the virial */
1213     fr->virdiffsix     = 0.5*virs[0];
1214     fr->virdifftwelve  = 0.5*virs[1];
1215   }
1216 }
1217
1218 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
1219                    gmx_large_int_t step,int natoms,
1220                    matrix box,real lambda,tensor pres,tensor virial,
1221                    real *prescorr, real *enercorr, real *dvdlcorr)
1222 {
1223     gmx_bool bCorrAll,bCorrPres;
1224     real dvdlambda,invvol,dens,ninter,avcsix,avctwelve,enerdiff,svir=0,spres=0;
1225     int  m;
1226
1227     *prescorr = 0;
1228     *enercorr = 0;
1229     *dvdlcorr = 0;
1230
1231     clear_mat(virial);
1232     clear_mat(pres);
1233
1234     if (ir->eDispCorr != edispcNO) {
1235         bCorrAll  = (ir->eDispCorr == edispcAllEner ||
1236                      ir->eDispCorr == edispcAllEnerPres);
1237         bCorrPres = (ir->eDispCorr == edispcEnerPres ||
1238                      ir->eDispCorr == edispcAllEnerPres);
1239
1240         invvol = 1/det(box);
1241         if (fr->n_tpi)
1242         {
1243             /* Only correct for the interactions with the inserted molecule */
1244             dens = (natoms - fr->n_tpi)*invvol;
1245             ninter = fr->n_tpi;
1246         }
1247         else
1248         {
1249             dens = natoms*invvol;
1250             ninter = 0.5*natoms;
1251         }
1252
1253         if (ir->efep == efepNO)
1254         {
1255             avcsix    = fr->avcsix[0];
1256             avctwelve = fr->avctwelve[0];
1257         }
1258         else
1259         {
1260             avcsix    = (1 - lambda)*fr->avcsix[0]    + lambda*fr->avcsix[1];
1261             avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
1262         }
1263
1264         enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
1265         *enercorr += avcsix*enerdiff;
1266         dvdlambda = 0.0;
1267         if (ir->efep != efepNO)
1268         {
1269             dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
1270         }
1271         if (bCorrAll)
1272         {
1273             enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
1274             *enercorr += avctwelve*enerdiff;
1275             if (fr->efep != efepNO)
1276             {
1277                 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
1278             }
1279         }
1280
1281         if (bCorrPres)
1282         {
1283             svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
1284             if (ir->eDispCorr == edispcAllEnerPres)
1285             {
1286                 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
1287             }
1288             /* The factor 2 is because of the Gromacs virial definition */
1289             spres = -2.0*invvol*svir*PRESFAC;
1290
1291             for(m=0; m<DIM; m++) {
1292                 virial[m][m] += svir;
1293                 pres[m][m] += spres;
1294             }
1295             *prescorr += spres;
1296         }
1297
1298         /* Can't currently control when it prints, for now, just print when degugging */
1299         if (debug)
1300         {
1301             if (bCorrAll) {
1302                 fprintf(debug,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1303                         avcsix,avctwelve);
1304             }
1305             if (bCorrPres)
1306             {
1307                 fprintf(debug,
1308                         "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
1309                         *enercorr,spres,svir);
1310             }
1311             else
1312             {
1313                 fprintf(debug,"Long Range LJ corr.: Epot %10g\n",*enercorr);
1314             }
1315         }
1316
1317         if (fr->bSepDVDL && do_per_step(step,ir->nstlog))
1318         {
1319             fprintf(fplog,sepdvdlformat,"Dispersion correction",
1320                     *enercorr,dvdlambda);
1321         }
1322         if (fr->efep != efepNO)
1323         {
1324             *dvdlcorr += dvdlambda;
1325         }
1326     }
1327 }
1328
1329 void do_pbc_first(FILE *fplog,matrix box,t_forcerec *fr,
1330                   t_graph *graph,rvec x[])
1331 {
1332   if (fplog)
1333     fprintf(fplog,"Removing pbc first time\n");
1334   calc_shifts(box,fr->shift_vec);
1335   if (graph) {
1336     mk_mshift(fplog,graph,fr->ePBC,box,x);
1337     if (gmx_debug_at)
1338       p_graph(debug,"do_pbc_first 1",graph);
1339     shift_self(graph,box,x);
1340     /* By doing an extra mk_mshift the molecules that are broken
1341      * because they were e.g. imported from another software
1342      * will be made whole again. Such are the healing powers
1343      * of GROMACS.
1344      */
1345     mk_mshift(fplog,graph,fr->ePBC,box,x);
1346     if (gmx_debug_at)
1347       p_graph(debug,"do_pbc_first 2",graph);
1348   }
1349   if (fplog)
1350     fprintf(fplog,"Done rmpbc\n");
1351 }
1352
1353 static void low_do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1354                             gmx_mtop_t *mtop,rvec x[],
1355                             gmx_bool bFirst)
1356 {
1357   t_graph *graph;
1358   int mb,as,mol;
1359   gmx_molblock_t *molb;
1360
1361   if (bFirst && fplog)
1362     fprintf(fplog,"Removing pbc first time\n");
1363
1364   snew(graph,1);
1365   as = 0;
1366   for(mb=0; mb<mtop->nmolblock; mb++) {
1367     molb = &mtop->molblock[mb];
1368     if (molb->natoms_mol == 1 ||
1369         (!bFirst && mtop->moltype[molb->type].cgs.nr == 1)) {
1370       /* Just one atom or charge group in the molecule, no PBC required */
1371       as += molb->nmol*molb->natoms_mol;
1372     } else {
1373       /* Pass NULL iso fplog to avoid graph prints for each molecule type */
1374       mk_graph_ilist(NULL,mtop->moltype[molb->type].ilist,
1375                      0,molb->natoms_mol,FALSE,FALSE,graph);
1376
1377       for(mol=0; mol<molb->nmol; mol++) {
1378         mk_mshift(fplog,graph,ePBC,box,x+as);
1379
1380         shift_self(graph,box,x+as);
1381         /* The molecule is whole now.
1382          * We don't need the second mk_mshift call as in do_pbc_first,
1383          * since we no longer need this graph.
1384          */
1385
1386         as += molb->natoms_mol;
1387       }
1388       done_graph(graph);
1389     }
1390   }
1391   sfree(graph);
1392 }
1393
1394 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
1395                        gmx_mtop_t *mtop,rvec x[])
1396 {
1397   low_do_pbc_mtop(fplog,ePBC,box,mtop,x,TRUE);
1398 }
1399
1400 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1401                  gmx_mtop_t *mtop,rvec x[])
1402 {
1403   low_do_pbc_mtop(fplog,ePBC,box,mtop,x,FALSE);
1404 }
1405
1406 void finish_run(FILE *fplog,t_commrec *cr,const char *confout,
1407                 t_inputrec *inputrec,
1408                 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
1409                 gmx_runtime_t *runtime,
1410                 gmx_bool bWriteStat)
1411 {
1412   int    i,j;
1413   t_nrnb *nrnb_tot=NULL;
1414   real   delta_t;
1415   double nbfs,mflop;
1416   double cycles[ewcNR];
1417
1418   wallcycle_sum(cr,wcycle,cycles);
1419
1420   if (cr->nnodes > 1) {
1421     if (SIMMASTER(cr))
1422       snew(nrnb_tot,1);
1423 #ifdef GMX_MPI
1424     MPI_Reduce(nrnb->n,nrnb_tot->n,eNRNB,MPI_DOUBLE,MPI_SUM,
1425                MASTERRANK(cr),cr->mpi_comm_mysim);
1426 #endif
1427   } else {
1428     nrnb_tot = nrnb;
1429   }
1430
1431   if (SIMMASTER(cr)) {
1432     print_flop(fplog,nrnb_tot,&nbfs,&mflop);
1433     if (cr->nnodes > 1) {
1434       sfree(nrnb_tot);
1435     }
1436   }
1437
1438   if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr)) {
1439     print_dd_statistics(cr,inputrec,fplog);
1440   }
1441
1442 #ifdef GMX_MPI
1443     if (PARTDECOMP(cr))
1444     {
1445         if (MASTER(cr))
1446         {
1447             t_nrnb     *nrnb_all;
1448             int        s;
1449             MPI_Status stat;
1450
1451             snew(nrnb_all,cr->nnodes);
1452             nrnb_all[0] = *nrnb;
1453             for(s=1; s<cr->nnodes; s++)
1454             {
1455                 MPI_Recv(nrnb_all[s].n,eNRNB,MPI_DOUBLE,s,0,
1456                          cr->mpi_comm_mysim,&stat);
1457             }
1458             pr_load(fplog,cr,nrnb_all);
1459             sfree(nrnb_all);
1460         }
1461         else
1462         {
1463             MPI_Send(nrnb->n,eNRNB,MPI_DOUBLE,MASTERRANK(cr),0,
1464                      cr->mpi_comm_mysim);
1465         }
1466     }
1467 #endif
1468
1469   if (SIMMASTER(cr)) {
1470     wallcycle_print(fplog,cr->nnodes,cr->npmenodes,runtime->realtime,
1471                     wcycle,cycles);
1472
1473     if (EI_DYNAMICS(inputrec->eI)) {
1474       delta_t = inputrec->delta_t;
1475     } else {
1476       delta_t = 0;
1477     }
1478
1479     if (fplog) {
1480         print_perf(fplog,runtime->proctime,runtime->realtime,
1481                    cr->nnodes-cr->npmenodes,
1482                    runtime->nsteps_done,delta_t,nbfs,mflop);
1483     }
1484     if (bWriteStat) {
1485         print_perf(stderr,runtime->proctime,runtime->realtime,
1486                    cr->nnodes-cr->npmenodes,
1487                    runtime->nsteps_done,delta_t,nbfs,mflop);
1488     }
1489
1490     /*
1491     runtime=inputrec->nsteps*inputrec->delta_t;
1492     if (bWriteStat) {
1493       if (cr->nnodes == 1)
1494         fprintf(stderr,"\n\n");
1495       print_perf(stderr,nodetime,realtime,runtime,&ntot,
1496                  cr->nnodes-cr->npmenodes,FALSE);
1497     }
1498     wallcycle_print(fplog,cr->nnodes,cr->npmenodes,realtime,wcycle,cycles);
1499     print_perf(fplog,nodetime,realtime,runtime,&ntot,cr->nnodes-cr->npmenodes,
1500                TRUE);
1501     if (PARTDECOMP(cr))
1502       pr_load(fplog,cr,nrnb_all);
1503     if (cr->nnodes > 1)
1504       sfree(nrnb_all);
1505     */
1506   }
1507 }
1508
1509 extern void initialize_lambdas(FILE *fplog,t_inputrec *ir,int *fep_state,real *lambda,double *lam0)
1510 {
1511     /* this function works, but could probably use a logic rewrite to keep all the different
1512        types of efep straight. */
1513
1514     int i;
1515     t_lambda *fep = ir->fepvals;
1516
1517     if ((ir->efep==efepNO) && (ir->bSimTemp == FALSE)) {
1518         for (i=0;i<efptNR;i++)  {
1519             lambda[i] = 0.0;
1520             if (lam0)
1521             {
1522                 lam0[i] = 0.0;
1523             }
1524         }
1525         return;
1526     } else {
1527         *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
1528                                              if checkpoint is set -- a kludge is in for now
1529                                              to prevent this.*/
1530         for (i=0;i<efptNR;i++)
1531         {
1532             /* overwrite lambda state with init_lambda for now for backwards compatibility */
1533             if (fep->init_lambda>=0) /* if it's -1, it was never initializd */
1534             {
1535                 lambda[i] = fep->init_lambda;
1536                 if (lam0) {
1537                     lam0[i] = lambda[i];
1538                 }
1539             }
1540             else
1541             {
1542                 lambda[i] = fep->all_lambda[i][*fep_state];
1543                 if (lam0) {
1544                     lam0[i] = lambda[i];
1545                 }
1546             }
1547         }
1548         if (ir->bSimTemp) {
1549             /* need to rescale control temperatures to match current state */
1550             for (i=0;i<ir->opts.ngtc;i++) {
1551                 if (ir->opts.ref_t[i] > 0) {
1552                     ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
1553                 }
1554             }
1555         }
1556     }
1557
1558     /* Send to the log the information on the current lambdas */
1559     if (fplog != NULL)
1560     {
1561         fprintf(fplog,"Initial vector of lambda components:[ ");
1562         for (i=0;i<efptNR;i++)
1563         {
1564             fprintf(fplog,"%10.4f ",lambda[i]);
1565         }
1566         fprintf(fplog,"]\n");
1567     }
1568     return;
1569 }
1570
1571
1572 void init_md(FILE *fplog,
1573              t_commrec *cr,t_inputrec *ir,const output_env_t oenv,
1574              double *t,double *t0,
1575              real *lambda, int *fep_state, double *lam0,
1576              t_nrnb *nrnb,gmx_mtop_t *mtop,
1577              gmx_update_t *upd,
1578              int nfile,const t_filenm fnm[],
1579              gmx_mdoutf_t **outf,t_mdebin **mdebin,
1580              tensor force_vir,tensor shake_vir,rvec mu_tot,
1581              gmx_bool *bSimAnn,t_vcm **vcm, t_state *state, unsigned long Flags)
1582 {
1583     int  i,j,n;
1584     real tmpt,mod;
1585
1586     /* Initial values */
1587     *t = *t0       = ir->init_t;
1588
1589     *bSimAnn=FALSE;
1590     for(i=0;i<ir->opts.ngtc;i++)
1591     {
1592         /* set bSimAnn if any group is being annealed */
1593         if(ir->opts.annealing[i]!=eannNO)
1594         {
1595             *bSimAnn = TRUE;
1596         }
1597     }
1598     if (*bSimAnn)
1599     {
1600         update_annealing_target_temp(&(ir->opts),ir->init_t);
1601     }
1602
1603     /* Initialize lambda variables */
1604     initialize_lambdas(fplog,ir,fep_state,lambda,lam0);
1605
1606     if (upd)
1607     {
1608         *upd = init_update(fplog,ir);
1609     }
1610
1611
1612     if (vcm != NULL)
1613     {
1614         *vcm = init_vcm(fplog,&mtop->groups,ir);
1615     }
1616
1617     if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
1618     {
1619         if (ir->etc == etcBERENDSEN)
1620         {
1621             please_cite(fplog,"Berendsen84a");
1622         }
1623         if (ir->etc == etcVRESCALE)
1624         {
1625             please_cite(fplog,"Bussi2007a");
1626         }
1627     }
1628
1629     init_nrnb(nrnb);
1630
1631     if (nfile != -1)
1632     {
1633         *outf = init_mdoutf(nfile,fnm,Flags,cr,ir,oenv);
1634
1635         *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
1636                               mtop,ir, (*outf)->fp_dhdl);
1637     }
1638
1639     if (ir->bAdress)
1640     {
1641       please_cite(fplog,"Fritsch12");
1642       please_cite(fplog,"Junghans10");
1643     }
1644     /* Initiate variables */
1645     clear_mat(force_vir);
1646     clear_mat(shake_vir);
1647     clear_rvec(mu_tot);
1648
1649     debug_gmx();
1650 }
1651
1652
1653
1654