Clean up CMAP placement in parameter list
authorErik Lindahl <erik@kth.se>
Mon, 23 Jun 2014 12:03:17 +0000 (14:03 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Tue, 24 Jun 2014 12:11:32 +0000 (14:11 +0200)
Some of the CMAP variables were always placed into and
read from the first entry of the parameter list. While this
did not result in any errors, this patch now places them
correctly in the F_CMAP position.

Fixes #1345.

Change-Id: Ic6e09e46c352976cfc1f57ff903c082b4ec43df8

src/gromacs/gmxpreprocess/grompp.c
src/gromacs/gmxpreprocess/topio.c
src/gromacs/gmxpreprocess/toppush.c

index b562b5188fd5bc4dbdcb96297521fd9558395e56..fc7b2c5555164fff59a550c796094b1cc27723c4 100644 (file)
@@ -1724,10 +1724,10 @@ int gmx_grompp(int argc, char *argv[])
     }
 
     /* If we are using CMAP, setup the pre-interpolation grid */
-    if (plist->ncmap > 0)
+    if (plist[F_CMAP].ncmap > 0)
     {
-        init_cmap_grid(&sys->ffparams.cmap_grid, plist->nc, plist->grid_spacing);
-        setup_cmap(plist->grid_spacing, plist->nc, plist->cmap, &sys->ffparams.cmap_grid);
+        init_cmap_grid(&sys->ffparams.cmap_grid, plist[F_CMAP].nc, plist[F_CMAP].grid_spacing);
+        setup_cmap(plist[F_CMAP].grid_spacing, plist[F_CMAP].nc, plist[F_CMAP].cmap, &sys->ffparams.cmap_grid);
     }
 
     set_wall_atomtype(atype, opts, ir, wi);
index e5b97a4d4320dd813de6c441f60b6f57a50867b7..00addd94532138dce23061c1826a464a478e8eca 100644 (file)
@@ -610,8 +610,8 @@ static char **read_topol(const char *infile, const char *outfile,
     comb     = 0;
 
     /* Init the number of CMAP torsion angles  and grid spacing */
-    plist->grid_spacing = 0;
-    plist->nc           = 0;
+    plist[F_CMAP].grid_spacing = 0;
+    plist[F_CMAP].nc           = 0;
 
     bWarn_copy_A_B = bFEP;
 
index 8f25139fd2cb8cb44a320eefad4886d1c69b254c..bc7f299e551dd888f80cf00673cd510384c3abbb 100644 (file)
@@ -1113,8 +1113,8 @@ push_cmaptype(directive d, t_params bt[], int nral, gpp_atomtype_t at,
     nrfp   = nrfpA+nrfpB;
 
     /* Allocate memory for the CMAP grid */
-    bt->ncmap += nrfp;
-    srenew(bt->cmap, bt->ncmap);
+    bt[F_CMAP].ncmap += nrfp;
+    srenew(bt[F_CMAP].cmap, bt[F_CMAP].ncmap);
 
     /* Read in CMAP parameters */
     sl = 0;
@@ -1126,7 +1126,7 @@ push_cmaptype(directive d, t_params bt[], int nral, gpp_atomtype_t at,
         }
         nn  = sscanf(line+start+sl, " %s ", s);
         sl += strlen(s);
-        bt->cmap[i+(bt->ncmap)-nrfp] = strtod(s, NULL);
+        bt[F_CMAP].cmap[i+(bt[F_CMAP].ncmap)-nrfp] = strtod(s, NULL);
 
         if (nn == 1)
         {
@@ -1144,7 +1144,7 @@ push_cmaptype(directive d, t_params bt[], int nral, gpp_atomtype_t at,
     {
         for (i = 0; i < ncmap; i++)
         {
-            bt->cmap[i+ncmap] = bt->cmap[i];
+            bt[F_CMAP].cmap[i+ncmap] = bt[F_CMAP].cmap[i];
         }
     }
     else
@@ -1167,12 +1167,12 @@ push_cmaptype(directive d, t_params bt[], int nral, gpp_atomtype_t at,
     /* Set grid spacing and the number of grids (we assume these numbers to be the same for all grids
      * so we can safely assign them each time
      */
-    bt->grid_spacing = nxcmap;     /* Or nycmap, they need to be equal */
-    bt->nc           = bt->nc + 1; /* Since we are incrementing here, we need to subtract later, see (*****) */
-    nct              = (nral+1) * bt->nc;
+    bt[F_CMAP].grid_spacing = nxcmap;            /* Or nycmap, they need to be equal */
+    bt[F_CMAP].nc           = bt[F_CMAP].nc + 1; /* Since we are incrementing here, we need to subtract later, see (*****) */
+    nct                     = (nral+1) * bt[F_CMAP].nc;
 
     /* Allocate memory for the cmap_types information */
-    srenew(bt->cmap_types, nct);
+    srenew(bt[F_CMAP].cmap_types, nct);
 
     for (i = 0; (i < nral); i++)
     {
@@ -1186,16 +1186,16 @@ push_cmaptype(directive d, t_params bt[], int nral, gpp_atomtype_t at,
         }
 
         /* Assign a grid number to each cmap_type */
-        bt->cmap_types[bt->nct++] = get_bond_atomtype_type(alc[i], bat);
+        bt[F_CMAP].cmap_types[bt[F_CMAP].nct++] = get_bond_atomtype_type(alc[i], bat);
     }
 
     /* Assign a type number to this cmap */
-    bt->cmap_types[bt->nct++] = bt->nc-1; /* Since we inremented earlier, we need to subtrac here, to get the types right (****) */
+    bt[F_CMAP].cmap_types[bt[F_CMAP].nct++] = bt[F_CMAP].nc-1; /* Since we inremented earlier, we need to subtrac here, to get the types right (****) */
 
     /* Check for the correct number of atoms (again) */
-    if (bt->nct != nct)
+    if (bt[F_CMAP].nct != nct)
     {
-        gmx_fatal(FARGS, "Incorrect number of atom types (%d) in cmap type %d\n", nct, bt->nc);
+        gmx_fatal(FARGS, "Incorrect number of atom types (%d) in cmap type %d\n", nct, bt[F_CMAP].nc);
     }
 
     /* Is this correct?? */
@@ -1533,7 +1533,7 @@ static gmx_bool default_cmap_params(t_params bondtype[],
     ct           = 0;
 
     /* Match the current cmap angle against the list of cmap_types */
-    for (i = 0; i < bondtype->nct && !bFound; i += 6)
+    for (i = 0; i < bondtype[F_CMAP].nct && !bFound; i += 6)
     {
         if (bB)
         {
@@ -1542,15 +1542,15 @@ static gmx_bool default_cmap_params(t_params bondtype[],
         else
         {
             if (
-                (get_atomtype_batype(at->atom[p->a[0]].type, atype) == bondtype->cmap_types[i])   &&
-                (get_atomtype_batype(at->atom[p->a[1]].type, atype) == bondtype->cmap_types[i+1]) &&
-                (get_atomtype_batype(at->atom[p->a[2]].type, atype) == bondtype->cmap_types[i+2]) &&
-                (get_atomtype_batype(at->atom[p->a[3]].type, atype) == bondtype->cmap_types[i+3]) &&
-                (get_atomtype_batype(at->atom[p->a[4]].type, atype) == bondtype->cmap_types[i+4]))
+                (get_atomtype_batype(at->atom[p->a[0]].type, atype) == bondtype[F_CMAP].cmap_types[i])   &&
+                (get_atomtype_batype(at->atom[p->a[1]].type, atype) == bondtype[F_CMAP].cmap_types[i+1]) &&
+                (get_atomtype_batype(at->atom[p->a[2]].type, atype) == bondtype[F_CMAP].cmap_types[i+2]) &&
+                (get_atomtype_batype(at->atom[p->a[3]].type, atype) == bondtype[F_CMAP].cmap_types[i+3]) &&
+                (get_atomtype_batype(at->atom[p->a[4]].type, atype) == bondtype[F_CMAP].cmap_types[i+4]))
             {
                 /* Found cmap torsion */
                 bFound       = TRUE;
-                ct           = bondtype->cmap_types[i+5];
+                ct           = bondtype[F_CMAP].cmap_types[i+5];
                 nparam_found = 1;
             }
         }