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