Replace qsort/bsearch with std::sort/lower_bound
authorRoland Schulz <roland.schulz@intel.com>
Fri, 21 Oct 2016 10:17:38 +0000 (03:17 -0700)
committerErik Lindahl <erik.lindahl@gmail.com>
Sun, 30 Oct 2016 09:35:04 +0000 (10:35 +0100)
Change-Id: Ic8dc3d86f9738751e463065a012ff1adc02bdf7c

src/gromacs/gmxana/gmx_cluster.cpp
src/gromacs/gmxana/nrama.cpp
src/gromacs/gmxpreprocess/resall.cpp

index dd621875495fab073a1f59de75415015ece86f76..f8a90473237217a3d612ec00d44377a889b5ab10 100644 (file)
@@ -330,43 +330,20 @@ static real rms_dist(int isize, real **d, real **d_r)
     return std::sqrt(r2);
 }
 
-static int rms_dist_comp(const void *a, const void *b)
+static bool rms_dist_comp(const t_dist &a, const t_dist &b)
 {
-    const t_dist *da, *db;
-
-    da = static_cast<const t_dist *>(a);
-    db = static_cast<const t_dist *>(b);
-
-    if (da->dist - db->dist < 0)
-    {
-        return -1;
-    }
-    else if (da->dist - db->dist > 0)
-    {
-        return 1;
-    }
-    return 0;
+    return a.dist < b.dist;
 }
 
-static int clust_id_comp(const void *a, const void *b)
+static bool clust_id_comp(const t_clustid &a, const t_clustid &b)
 {
-    const t_clustid *da, *db;
-
-    da = static_cast<const t_clustid *>(a);
-    db = static_cast<const t_clustid *>(b);
-
-    return da->clust - db->clust;
+    return a.clust < b.clust;
 }
 
-static int nrnb_comp(const void *a, const void *b)
+static bool nrnb_comp(const t_nnb &a, const t_nnb &b)
 {
-    const t_nnb *da, *db;
-
-    da = static_cast<const t_nnb *>(a);
-    db = static_cast<const t_nnb *>(b);
-
-    /* return the b-a, we want highest first */
-    return db->nr - da->nr;
+    /* return b<a, we want highest first */
+    return b.nr < a.nr;
 }
 
 void gather(t_mat *m, real cutoff, t_clusters *clust)
@@ -393,7 +370,7 @@ void gather(t_mat *m, real cutoff, t_clusters *clust)
     {
         gmx_incons("gather algortihm");
     }
-    qsort(d, nn, sizeof(d[0]), rms_dist_comp);
+    std::sort(d, d+nn, rms_dist_comp);
 
     /* Now we make a cluster index for all of the conformations */
     c = new_clustid(n1);
@@ -424,7 +401,7 @@ void gather(t_mat *m, real cutoff, t_clusters *clust)
     while (bChange);
     fprintf(stderr, "\nSorting and renumbering clusters\n");
     /* Sort on cluster number */
-    qsort(c, n1, sizeof(c[0]), clust_id_comp);
+    std::sort(c, c+n1, clust_id_comp);
 
     /* Renumber clusters */
     cid = 1;
@@ -526,7 +503,7 @@ static void jarvis_patrick(int n1, real **mat, int M, int P,
             row[j].j    = j;
             row[j].dist = mat[i][j];
         }
-        qsort(row, n1, sizeof(row[0]), rms_dist_comp);
+        std::sort(row, row+n1, rms_dist_comp);
         if (M > 0)
         {
             /* Put the M nearest neighbors in the list */
@@ -623,7 +600,7 @@ static void jarvis_patrick(int n1, real **mat, int M, int P,
 
     fprintf(stderr, "\nSorting and renumbering clusters\n");
     /* Sort on cluster number */
-    qsort(c, n1, sizeof(c[0]), clust_id_comp);
+    std::sort(c, c+n1, clust_id_comp);
 
     /* Renumber clusters */
     cid = 1;
@@ -720,7 +697,7 @@ static void gromos(int n1, real **mat, real rmsdcut, t_clusters *clust)
     sfree(row);
 
     /* sort neighbor list on number of neighbors, largest first */
-    qsort(nnb, n1, sizeof(nnb[0]), nrnb_comp);
+    std::sort(nnb, nnb+n1, nrnb_comp);
 
     if (debug)
     {
@@ -763,7 +740,7 @@ static void gromos(int n1, real **mat, real rmsdcut, t_clusters *clust)
         }
         /* sort again on nnb[].nr, because we have new # neighbors: */
         /* but we only need to sort upto i, i.e. when nnb[].nr>0 */
-        qsort(nnb, i, sizeof(nnb[0]), nrnb_comp);
+        std::sort(nnb, nnb+i, nrnb_comp);
 
         fprintf(stderr, "\b\b\b\b%4d", k);
         /* new cluster id */
index eb202583e6466272054afd05d7bf2ce10ffe6e72..5759488a5472586857e67802534b50b4da614bbc 100644 (file)
@@ -41,6 +41,8 @@
 #include <cstdlib>
 #include <cstring>
 
+#include <algorithm>
+
 #include "gromacs/listed-forces/bonded.h"
 #include "gromacs/pbcutil/rmpbc.h"
 #include "gromacs/utility/cstringutil.h"
 static const char *pp_pat[] = { "C", "N", "CA", "C", "N" };
 #define NPP (sizeof(pp_pat)/sizeof(pp_pat[0]))
 
-static int d_comp(const void *a, const void *b)
+static bool d_comp(const t_dih &a, const t_dih &b)
 {
-    const t_dih *da, *db;
-
-    da = static_cast<const t_dih *>(a);
-    db = static_cast<const t_dih *>(b);
-
-    if (da->ai[1] < db->ai[1])
+    if (a.ai[1] < b.ai[1])
     {
-        return -1;
+        return true;
     }
-    else if (da->ai[1] == db->ai[1])
+    else if (a.ai[1] == b.ai[1])
     {
-        return (da->ai[2] - db->ai[2]);
+        return a.ai[2] < b.ai[2];
     }
     else
     {
-        return 1;
+        return false;
     }
 }
 
@@ -217,8 +214,8 @@ static void get_dih_props(t_xrama *xr, const t_idef *idef, int mult)
 
         key.ai[1] = ia[2];
         key.ai[2] = ia[3];
-        if ((dd = static_cast<t_dih *>(bsearch(&key, xr->dih, xr->ndih, (size_t)sizeof(key), d_comp)))
-            != NULL)
+        dd        = std::lower_bound(xr->dih, xr->dih+xr->ndih, key, d_comp);
+        if (dd < xr->dih+xr->ndih && !d_comp(key, *dd))
         {
             dd->mult = idef->iparams[ft].pdihs.mult;
             dd->phi0 = idef->iparams[ft].pdihs.phiA;
index e262876d1abf65ea9c624ed8ba20ef7c599c24ba..44dcc5f4b61e5fa61f190334af2c103aa5d1e791 100644 (file)
@@ -260,16 +260,6 @@ static void check_rtp(int nrtp, t_restp rtp[], char *libfn)
     }
 }
 
-static int comprtp(const void *a, const void *b)
-{
-    const t_restp *ra, *rb;
-
-    ra = static_cast<const t_restp *>(a);
-    rb = static_cast<const t_restp *>(b);
-
-    return gmx_strcasecmp(ra->resname, rb->resname);
-}
-
 int get_bt(char* header)
 {
     int i;
@@ -551,7 +541,7 @@ void read_resall(char *rrdb, int *nrtpptr, t_restp **rtp,
     srenew(rrtp, nrtp);
 
     fprintf(stderr, "\nSorting it all out...\n");
-    qsort(rrtp, nrtp, (size_t)sizeof(rrtp[0]), comprtp);
+    std::sort(rrtp, rrtp+nrtp, [](const t_restp &a, const t_restp &b) {return gmx_strcasecmp(a.resname, b.resname) < 0; });
 
     check_rtp(nrtp, rrtp, rrdb);