Fix clang 3.5 warnings regarding *abs*
authorAlexey Shvetsov <alexxy@omrb.pnpi.spb.ru>
Sun, 20 Apr 2014 13:01:54 +0000 (17:01 +0400)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Fri, 9 May 2014 19:06:10 +0000 (21:06 +0200)
Fix abs() family functions type usage mismatch.

Change-Id: I85ed2931d681aa1ad024678b4209a524abc2cc61
Signed-off-by: Alexey Shvetsov <alexxy@omrb.pnpi.spb.ru>
src/gmxlib/libxdrf.c
src/tools/geminate.c
src/tools/gmx_anaeig.c
src/tools/gmx_current.c
src/tools/gmx_mindist.c
src/tools/gmx_rms.c
src/tools/gmx_xpm2ps.c

index 3194f94528524b1d9d8a147a5da74f20db5b1872..d3209280241d9c746c01700b42f4c1aca5de0869 100644 (file)
@@ -1361,6 +1361,12 @@ xdr_xtc_estimate_dt(FILE *fp, XDR *xdrs, int natoms, gmx_bool * bOK)
     return res;
 }
 
+/* can be replaced by llabs in 5.0 (allows C99)*/
+static gmx_off_t
+gmx_llabs(gmx_off_t x)
+{
+    return x>0?x:-x;
+}
 
 int
 xdr_xtc_seek_frame(int frame, FILE *fp, XDR *xdrs, int natoms)
@@ -1399,7 +1405,7 @@ xdr_xtc_seek_frame(int frame, FILE *fp, XDR *xdrs, int natoms)
         {
             return -1;
         }
-        if (fr != frame && abs(low-high) > header_size)
+        if (fr != frame && gmx_llabs(low-high) > header_size)
         {
             if (fr < frame)
             {
@@ -1533,10 +1539,9 @@ int xdr_xtc_seek_time(real time, FILE *fp, XDR *xdrs, int natoms, gmx_bool bSeek
            the current time and the target time is bigger than dt and above all the distance between high
            and low is bigger than 1 frame, then do another step of binary search. Otherwise stop and check
            if we reached the solution */
-        if ((((t < time && dt_sign >= 0) || (t > time && dt_sign == -1)) || ((t
-                                                                              - time) >= dt && dt_sign >= 0)
-             || ((time - t) >= -dt && dt_sign < 0)) && (abs(low - high)
-                                                        > header_size))
+        if ((((t < time && dt_sign >= 0) || (t > time && dt_sign == -1)) ||
+             ((t - time) >= dt && dt_sign >= 0) || ((time - t) >= -dt && dt_sign < 0)) &&
+            (gmx_llabs(low - high) > header_size))
         {
             if (dt >= 0 && dt_sign != -1)
             {
@@ -1574,7 +1579,7 @@ int xdr_xtc_seek_time(real time, FILE *fp, XDR *xdrs, int natoms, gmx_bool bSeek
         }
         else
         {
-            if (abs(low - high) <= header_size)
+            if (gmx_llabs(low - high) <= header_size)
             {
                 break;
             }
index d9e69984502dc1f6ad5855c63de2ac8aa62625fe..ee008df124b9effc8475c9ea32e4ba132242d966 100644 (file)
@@ -1195,10 +1195,10 @@ extern void fixGemACF(double *ct, int len)
     bBad = FALSE;
 
     /* An acf of binary data must be one at t=0. */
-    if (abs(ct[0]-1.0) > 1e-6)
+    if (fabs(ct[0]-1.0) > 1e-6)
     {
         ct[0] = 1.0;
-        fprintf(stderr, "|ct[0]-1.0| = %1.6d. Setting ct[0] to 1.0.\n", abs(ct[0]-1.0));
+        fprintf(stderr, "|ct[0]-1.0| = %1.6f. Setting ct[0] to 1.0.\n", fabs(ct[0]-1.0));
     }
 
     for (i = 0; i < len; i++)
index bda9de7c879867ac77c2a0a473af8430980d744b..b01c5754baad4f2cde943cbaa65f874eba034257 100644 (file)
@@ -258,7 +258,7 @@ static void write_xvgr_graphs(const char *file, int ngraphs, int nsetspergraph,
         {
             for (i = 0; i < n; i++)
             {
-                if (bSplit && i > 0 && abs(x[i]) < 1e-5)
+                if (bSplit && i > 0 && fabs(x[i]) < 1e-5)
                 {
                     if (output_env_get_print_xvgr_codes(oenv))
                     {
@@ -671,7 +671,7 @@ static void project(const char *trajfile, t_topology *top, int ePBC, matrix topb
                            oenv);
         for (i = 0; i < nframes; i++)
         {
-            if (bSplit && i > 0 && abs(inprod[noutvec][i]) < 1e-5)
+            if (bSplit && i > 0 && fabs(inprod[noutvec][i]) < 1e-5)
             {
                 fprintf(xvgrout, "&\n");
             }
@@ -759,7 +759,7 @@ static void project(const char *trajfile, t_topology *top, int ePBC, matrix topb
             j = 0;
             for (i = 0; i < atoms.nr; i++)
             {
-                if (j > 0 && bSplit && abs(inprod[noutvec][i]) < 1e-5)
+                if (j > 0 && bSplit && fabs(inprod[noutvec][i]) < 1e-5)
                 {
                     fprintf(out, "TER\n");
                     j = 0;
index ad2ab1903568b30b19f94e043cfc0328cb93d1c2..efa25a3d4a22e7d0df3ad3802e58d6625c3beeb7 100644 (file)
@@ -141,7 +141,7 @@ static gmx_bool precalc(t_topology top, real mass2[], real qmol[])
         }
     }
 
-    if (abs(qall) > 0.01)
+    if (fabs(qall) > 0.01)
     {
         printf("\n\nSystem not neutral (q=%f) will not calculate translational part of the dipole moment.\n", qall);
         bNEU = FALSE;
index c44c181412ac7763938d82d6a2ba0ae59263ca68..30d5abbff38811151395b16a2313c6d92cb5b248 100644 (file)
@@ -172,7 +172,7 @@ static void periodic_mindist_plot(const char *trxfn, const char *outfn,
             ind_mini = ind_min[0];
             ind_minj = ind_min[1];
         }
-        if (bSplit && !bFirst && abs(t/output_env_get_time_factor(oenv)) < 1e-5)
+        if (bSplit && !bFirst && fabs(t/output_env_get_time_factor(oenv)) < 1e-5)
         {
             fprintf(out, "&\n");
         }
@@ -423,7 +423,7 @@ void dist_plot(const char *fn, const char *afile, const char *dfile,
     bFirst = TRUE;
     do
     {
-        if (bSplit && !bFirst && abs(t/output_env_get_time_factor(oenv)) < 1e-5)
+        if (bSplit && !bFirst && fabs(t/output_env_get_time_factor(oenv)) < 1e-5)
         {
             fprintf(dist, "&\n");
             if (num)
index f4109f13fdfa5f4977af5fb20738803e72c21e50..628ed7f0e49fac3ddf0387bc671a9421f8eb3db5 100644 (file)
@@ -1135,7 +1135,7 @@ int gmx_rms(int argc, char *argv[])
     for (i = 0; (i < teller); i++)
     {
         if (bSplit && i > 0 &&
-            abs(time[bPrev ? freq*i : i]/output_env_get_time_factor(oenv)) < 1e-5)
+            fabs(time[bPrev ? freq*i : i]/output_env_get_time_factor(oenv)) < 1e-5)
         {
             fprintf(fp, "&\n");
         }
@@ -1177,7 +1177,7 @@ int gmx_rms(int argc, char *argv[])
         }
         for (i = 0; (i < teller); i++)
         {
-            if (bSplit && i > 0 && abs(time[i]) < 1e-5)
+            if (bSplit && i > 0 && fabs(time[i]) < 1e-5)
             {
                 fprintf(fp, "&\n");
             }
index efb73b420839911376252cf64f255084972cb3cf..f61e6a65f49877770a4b77ebcf9321276fdfbc3b 100644 (file)
@@ -525,8 +525,8 @@ static void draw_zerolines(t_psdata out, real x0, real y0, real w,
                 xx = xx00+(x+0.7)*psr->xboxsize;
                 /* draw lines whenever tick label almost zero (e.g. next trajectory) */
                 if (x != 0 && x < mat[i].nx-1 &&
-                    abs(mat[i].axis_x[x]) <
-                    0.1*abs(mat[i].axis_x[x+1]-mat[i].axis_x[x]) )
+                    fabs(mat[i].axis_x[x]) <
+                    0.1*fabs(mat[i].axis_x[x+1]-mat[i].axis_x[x]) )
                 {
                     ps_line (out, xx, yy00, xx, yy00+dy+2);
                 }
@@ -540,8 +540,8 @@ static void draw_zerolines(t_psdata out, real x0, real y0, real w,
                 yy = yy00+(y+0.7)*psr->yboxsize;
                 /* draw lines whenever tick label almost zero (e.g. next trajectory) */
                 if (y != 0 && y < mat[i].ny-1 &&
-                    abs(mat[i].axis_y[y]) <
-                    0.1*abs(mat[i].axis_y[y+1]-mat[i].axis_y[y]) )
+                    fabs(mat[i].axis_y[y]) <
+                    0.1*fabs(mat[i].axis_y[y+1]-mat[i].axis_y[y]) )
                 {
                     ps_line (out, xx00, yy, xx00+w+2, yy);
                 }
@@ -1137,7 +1137,7 @@ void zero_lines(int nmat, t_matrix *mat, t_matrix *mat2)
             }
             for (x = 0; x < mats[i].nx-1; x++)
             {
-                if (abs(mats[i].axis_x[x+1]) < 1e-5)
+                if (fabs(mats[i].axis_x[x+1]) < 1e-5)
                 {
                     for (y = 0; y < mats[i].ny; y++)
                     {
@@ -1147,7 +1147,7 @@ void zero_lines(int nmat, t_matrix *mat, t_matrix *mat2)
             }
             for (y = 0; y < mats[i].ny-1; y++)
             {
-                if (abs(mats[i].axis_y[y+1]) < 1e-5)
+                if (fabs(mats[i].axis_y[y+1]) < 1e-5)
                 {
                     for (x = 0; x < mats[i].nx; x++)
                     {