Fix -Wmaybe-uninitialized warnings in some files
authorAlexey Shvetsov <alexxy@omrb.pnpi.spb.ru>
Fri, 29 Nov 2013 18:54:50 +0000 (22:54 +0400)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Wed, 26 Feb 2014 17:28:23 +0000 (18:28 +0100)
Also clean-up and clarification of clean_vsite_dihs.

Change-Id: Ica4b1d17039e23bc2274a20fe5abca18295a7593
Signed-off-by: Alexey Shvetsov <alexxy@omrb.pnpi.spb.ru>
src/gromacs/gmxpreprocess/vsite_parm.c
src/gromacs/mdlib/forcerec.c

index 3a961a47d9d2a460571fc0285258a627a22bc31c..eaa35b09a16a9f33157100bf109948fbd7b773e5 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2013, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -1404,18 +1404,21 @@ static void clean_vsite_angles(t_params *plist, t_pindex pindex[],
 static void clean_vsite_dihs(t_params *plist, t_pindex pindex[],
                              int cftype, int vsite_type[])
 {
-    int          ftype, i, parnr, k, l, m, n, nvsite, kept_i, vsnral;
-    atom_id      atom, constr;
-    atom_id      vsiteatoms[4];
-    gmx_bool     bKeep, bUsed, bPresent;
-    t_params    *ps;
+    int       i, kept_i;
+    t_params *ps;
 
     ps = &(plist[cftype]);
 
-    vsnral = 0;
     kept_i = 0;
     for (i = 0; (i < ps->nr); i++) /* for all dihedrals in the plist */
     {
+        int      ftype, parnr, k, l, m, n, nvsite;
+        int      vsnral = 0;            /* keep the compiler happy */
+        atom_id  atom, constr;
+        atom_id  vsiteatoms[4] = { 0 }; /* init to zero to make gcc4.8 happy */
+        gmx_bool bKeep, bUsed, bPresent;
+
+
         bKeep = FALSE;
         /* check if all virtual sites are constructed from the same atoms */
         nvsite = 0;
@@ -1424,8 +1427,7 @@ static void clean_vsite_dihs(t_params *plist, t_pindex pindex[],
             atom = ps->param[i].a[k];
             if (vsite_type[atom] != NOTSET)
             {
-                nvsite++;
-                if (nvsite == 1)
+                if (nvsite == 0)
                 {
                     /* store construction atoms of first vsite */
                     vsnral = NRAL(pindex[atom].ftype)-1;
@@ -1445,27 +1447,30 @@ static void clean_vsite_dihs(t_params *plist, t_pindex pindex[],
                     }
                 }
                 else
-                /* check if this vsite is constructed from the same atoms */
-                if (vsnral == NRAL(pindex[atom].ftype)-1)
                 {
-                    for (m = 0; (m < vsnral) && !bKeep; m++)
+                    /* check if this vsite is constructed from the same atoms */
+                    if (vsnral == NRAL(pindex[atom].ftype)-1)
                     {
-                        bPresent = FALSE;
-                        constr   =
-                            plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
-                        for (n = 0; (n < vsnral) && !bPresent; n++)
+                        for (m = 0; (m < vsnral) && !bKeep; m++)
                         {
-                            if (constr == vsiteatoms[n])
+                            bPresent = FALSE;
+                            constr   =
+                                plist[pindex[atom].ftype].param[pindex[atom].parnr].a[m+1];
+                            for (n = 0; (n < vsnral) && !bPresent; n++)
                             {
-                                bPresent = TRUE;
+                                if (constr == vsiteatoms[n])
+                                {
+                                    bPresent = TRUE;
+                                }
+                            }
+                            if (!bPresent)
+                            {
+                                bKeep = TRUE;
                             }
-                        }
-                        if (!bPresent)
-                        {
-                            bKeep = TRUE;
                         }
                     }
                 }
+                nvsite++;
             }
         }
 
@@ -1482,6 +1487,7 @@ static void clean_vsite_dihs(t_params *plist, t_pindex pindex[],
             atom = ps->param[i].a[k];
             if (vsite_type[atom] == NOTSET)
             {
+                /* vsnral will be set here, we don't get here with nvsite==0 */
                 bUsed = FALSE;
                 for (m = 0; (m < vsnral) && !bUsed; m++)
                 {
index 0c7db6ed40fbdadb4eb34f350c6366292927c389..8d18cb7a249fe54facacb4a37c29d9d14995f49f 100644 (file)
@@ -228,15 +228,15 @@ check_solvent_cg(const gmx_moltype_t    *molt,
                  int                     cginfo,
                  int                    *cg_sp)
 {
-    const t_blocka     *  excl;
+    const t_blocka       *excl;
     t_atom               *atom;
     int                   j, k;
     int                   j0, j1, nj;
     gmx_bool              perturbed;
     gmx_bool              has_vdw[4];
     gmx_bool              match;
-    real                  tmp_charge[4];
-    int                   tmp_vdwtype[4];
+    real                  tmp_charge[4]  = { 0.0 }; /* init to zero to make gcc4.8 happy */
+    int                   tmp_vdwtype[4] = { 0 };   /* init to zero to make gcc4.8 happy */
     int                   tjA;
     gmx_bool              qm;
     solvent_parameters_t *solvent_parameters;