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