Fixing a problem when dhdl and replex are not multiples
authorMichael Shirts <michael.shirts@virginia.edu>
Wed, 19 Jun 2013 02:31:32 +0000 (22:31 -0400)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Thu, 20 Jun 2013 14:21:32 +0000 (16:21 +0200)
Problems is fixed by guaranteeing the free energy information
is now updated at a frequency that is the greatest common
denominator of the relevant internal variables.

An existing gcd function was taken from gmx_genion and put
in the include/maths.h file.

Change-Id: I7473185ab02c7c9a7557fbaa165bb20b63bafb7a

include/maths.h
src/kernel/md.c
src/tools/gmx_genion.c

index bf95c7ca6c8f47df4a31e9132d61959a1bf54b2d..74e937a9d011cc0b603832b339ea7490675cb463 100644 (file)
@@ -210,6 +210,18 @@ check_int_multiply_for_overflow(gmx_large_int_t  a,
                                 gmx_large_int_t  b,
                                 gmx_large_int_t *result);
 
+static int gmx_greatest_common_divisor(int p, int q)
+{
+    int tmp;
+    while (q != 0)
+    {
+        tmp = q;
+        q = p % q;
+        p = tmp;
+    }
+    return p;
+}
+
 #ifdef __cplusplus
 }
 #endif
index 4bf12a67d321324ae7fb22a25d2fb39f1aa35d53..f82e6974ca5095ff0eb4b7268ddf9b63798ca64d 100644 (file)
@@ -592,15 +592,16 @@ double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
 
     debug_gmx();
 
-    /* set free energy calculation frequency as the minimum of nstdhdl, nstexpanded, and nstrepl_ex_nst*/
+    /* set free energy calculation frequency as the minimum
+       greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
     nstfep = ir->fepvals->nstdhdl;
-    if (ir->bExpanded && (nstfep > ir->expandedvals->nstexpanded))
+    if (ir->bExpanded)
     {
-        nstfep = ir->expandedvals->nstexpanded;
+        nstfep = gmx_greatest_common_divisor(ir->fepvals->nstdhdl,nstfep);
     }
-    if (repl_ex_nst > 0 && nstfep > repl_ex_nst)
+    if (repl_ex_nst > 0)
     {
-        nstfep = repl_ex_nst;
+        nstfep = gmx_greatest_common_divisor(repl_ex_nst,nstfep);
     }
 
     /* I'm assuming we need global communication the first time! MRS */
index 227e8c99329dbdb180384c583cec12d12a94e184..f691be21b8c64afd8d876d245205b39991d77748 100644 (file)
 #include "mtop_util.h"
 #include "gmx_ana.h"
 
-static int greatest_common_divisor(int p, int q)
-{
-    int tmp;
-    while (q != 0)
-    {
-        tmp = q;
-        q = p % q;
-        p = tmp;
-    }
-    return p;
-}
-
 static void insert_ion(int nsa, int *nwater,
                        gmx_bool bSet[], int repl[], atom_id index[],
                        rvec x[], t_pbc *pbc,
@@ -459,7 +447,7 @@ int gmx_genion(int argc, char *argv[])
 
         /* Check if the system is neutralizable
          * is (qdelta == p_q*p_num + n_q*n_num) solvable for p_num and n_num? */
-        int gcd = greatest_common_divisor(n_q, p_q);
+        int gcd = gmx_greatest_common_divisor(n_q, p_q);
         if ((qdelta % gcd) != 0)
         {
             gmx_fatal(FARGS, "Can't neutralize this system using -nq %d and"