fixed dispersion correction with parallel EM
authorBerk Hess <hess@kth.se>
Thu, 26 Apr 2012 12:40:13 +0000 (14:40 +0200)
committerBerk Hess <hess@kth.se>
Thu, 26 Apr 2012 12:40:13 +0000 (14:40 +0200)
With EM the energy and pressure dispersion correction terms
were multiplied by the number of nodes. Fixes #901

Change-Id: I289ac0de7d7a9c1f72e939b840c549c1cb49a52a

src/mdlib/minimize.c

index 690f0d98e6f683fea6334b5a9f6a7339451aa2a1..78e922d610d7551687bcd3f61b55c3d4c23190fa 100644 (file)
@@ -697,16 +697,6 @@ static void evaluate_energy(FILE *fplog,gmx_bool bVerbose,t_commrec *cr,
   clear_mat(shake_vir);
   clear_mat(pres);
 
-  /* Calculate long range corrections to pressure and energy */
-  calc_dispcorr(fplog,inputrec,fr,count,top_global->natoms,ems->s.box,ems->s.lambda,
-                pres,force_vir,&prescorr,&enercorr,&dvdlcorr);
-  /* don't think these next 4 lines  can be moved in for now, because we 
-     don't always want to write it -- figure out how to clean this up MRS 8/4/2009 */
-  enerd->term[F_DISPCORR] = enercorr;
-  enerd->term[F_EPOT] += enercorr;
-  enerd->term[F_PRES] += prescorr;
-  enerd->term[F_DVDL] += dvdlcorr;
-
     /* Communicate stuff when parallel */
     if (PAR(cr) && inputrec->eI != eiNM)
     {
@@ -723,6 +713,14 @@ static void evaluate_energy(FILE *fplog,gmx_bool bVerbose,t_commrec *cr,
         wallcycle_stop(wcycle,ewcMoveE);
     }
 
+    /* Calculate long range corrections to pressure and energy */
+    calc_dispcorr(fplog,inputrec,fr,count,top_global->natoms,ems->s.box,ems->s.lambda,
+                  pres,force_vir,&prescorr,&enercorr,&dvdlcorr);
+    enerd->term[F_DISPCORR] = enercorr;
+    enerd->term[F_EPOT] += enercorr;
+    enerd->term[F_PRES] += prescorr;
+    enerd->term[F_DVDL] += dvdlcorr;
+
   ems->epot = enerd->term[F_EPOT];
   
   if (constr) {