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