Remove unused sign parameter from dih_angle()
authorBerk Hess <hess@kth.se>
Sun, 29 Oct 2017 21:20:54 +0000 (22:20 +0100)
committerMark Abraham <mark.j.abraham@gmail.com>
Mon, 30 Oct 2017 15:25:21 +0000 (16:25 +0100)
Change-Id: I88a73ca49b6acfc59b4baf0d847aa81542a870ca

src/gromacs/gmxana/anadih.cpp
src/gromacs/gmxana/gmx_dipoles.cpp
src/gromacs/gmxana/hxprops.cpp
src/gromacs/gmxana/nrama.cpp
src/gromacs/gmxpreprocess/x2top.cpp
src/gromacs/listed-forces/bonded.cpp
src/gromacs/listed-forces/bonded.h
src/gromacs/listed-forces/tests/bonded.cpp
src/gromacs/listed-forces/tests/refdata/BondedTest_DihedralAnglePbcNone.xml
src/gromacs/listed-forces/tests/refdata/BondedTest_DihedralAnglePbcXy.xml
src/gromacs/listed-forces/tests/refdata/BondedTest_DihedralAnglePbcXyz.xml

index 2393f556f012acf64784ff3e536f7cbd27c3ee85..142d1917aefd7ad83d01a1ed7f7bab85b09b0beb 100644 (file)
@@ -734,14 +734,14 @@ static void calc_dihs(struct t_pbc *pbc,
 {
     int  i, ix, t1, t2, t3;
     rvec r_ij, r_kj, r_kl, m, n;
-    real sign, aaa;
+    real aaa;
 
     for (i = ix = 0; (ix < n4); i++, ix += 4)
     {
         aaa = dih_angle(x_s[index[ix]], x_s[index[ix+1]], x_s[index[ix+2]],
                         x_s[index[ix+3]], pbc,
                         r_ij, r_kj, r_kl, m, n,
-                        &sign, &t1, &t2, &t3);
+                        &t1, &t2, &t3);
 
         ang[i] = aaa; /* not taking into account ryckaert bellemans yet */
     }
index 32dd9125ccf5aad2f612cf74ee615fe4a8b04a66..169a028d2f7aaf2c108a894f140479d7a203fd3a 100644 (file)
@@ -238,7 +238,6 @@ static void do_gkr(t_gkrbin *gb, int ncos, int *ngrp, int *molindex[],
                 {
                     rvec xi, xj, xk, xl;
                     rvec r_ij, r_kj, r_kl, mm, nn;
-                    real sign;
                     int  t1, t2, t3;
 
                     copy_rvec(xcm[grp0][i], xj);
@@ -247,7 +246,7 @@ static void do_gkr(t_gkrbin *gb, int ncos, int *ngrp, int *molindex[],
                     rvec_add(xk, mu[gj], xl);
                     phi = dih_angle(xi, xj, xk, xl, &pbc,
                                     r_ij, r_kj, r_kl, mm, nn, /* out */
-                                    &sign, &t1, &t2, &t3);
+                                    &t1, &t2, &t3);
                     cosa = std::cos(phi);
                 }
                 else
index 497281d13b6f21030fd99832d22f60cfe1d917b5..9c6d7aea0a09b7951e7febf03d70b8ab29a10d41 100644 (file)
@@ -185,7 +185,6 @@ real ca_phi(int gnx, int index[], rvec x[])
     real phi, phitot;
     int  i, ai, aj, ak, al, t1, t2, t3;
     rvec r_ij, r_kj, r_kl, m, n;
-    real sign;
 
     if (gnx <= 4)
     {
@@ -202,7 +201,7 @@ real ca_phi(int gnx, int index[], rvec x[])
         phi = RAD2DEG*
             dih_angle(x[ai], x[aj], x[ak], x[al], nullptr,
                       r_ij, r_kj, r_kl, m, n,
-                      &sign, &t1, &t2, &t3);
+                      &t1, &t2, &t3);
         phitot += phi;
     }
 
@@ -487,7 +486,6 @@ void calc_hxprops(int nres, t_bb bb[], const rvec x[])
 {
     int  i, ao, an, t1, t2, t3;
     rvec dx, r_ij, r_kj, r_kl, m, n;
-    real sign;
 
     for (i = 0; (i < nres); i++)
     {
@@ -515,11 +513,11 @@ void calc_hxprops(int nres, t_bb bb[], const rvec x[])
         bb[i].phi = RAD2DEG*
             dih_angle(x[bb[i].Cprev], x[bb[i].N], x[bb[i].CA], x[bb[i].C], nullptr,
                       r_ij, r_kj, r_kl, m, n,
-                      &sign, &t1, &t2, &t3);
+                      &t1, &t2, &t3);
         bb[i].psi = RAD2DEG*
             dih_angle(x[bb[i].N], x[bb[i].CA], x[bb[i].C], x[bb[i].Nnext], nullptr,
                       r_ij, r_kj, r_kl, m, n,
-                      &sign, &t1, &t2, &t3);
+                      &t1, &t2, &t3);
         bb[i].pprms2 = gmx::square(bb[i].phi-PHI_AHX)+gmx::square(bb[i].psi-PSI_AHX);
 
         bb[i].jcaha +=
index 7a95fe66a471ced5f9ddd6b15edb4a61e3d79743..36dcc8a9b1b933f6d13f29d1c46dbe536ad4a7ad 100644 (file)
@@ -74,7 +74,6 @@ static void calc_dihs(t_xrama *xr)
 {
     int          i, t1, t2, t3;
     rvec         r_ij, r_kj, r_kl, m, n;
-    real         sign;
     t_dih       *dd;
     gmx_rmpbc_t  gpbc = nullptr;
 
@@ -88,7 +87,7 @@ static void calc_dihs(t_xrama *xr)
         dd->ang = dih_angle(xr->x[dd->ai[0]], xr->x[dd->ai[1]],
                             xr->x[dd->ai[2]], xr->x[dd->ai[3]],
                             nullptr,
-                            r_ij, r_kj, r_kl, m, n, &sign, &t1, &t2, &t3);
+                            r_ij, r_kj, r_kl, m, n, &t1, &t2, &t3);
     }
 }
 
index e9e85fb4a0dfdb6a3d602e5c1317fdc6d901d2ff..7ca13e3733d5c7730b1bd5c5fa6ff88cac493fed 100644 (file)
@@ -282,7 +282,7 @@ static void calc_angles_dihs(t_params *ang, t_params *dih, const rvec x[], gmx_b
 {
     int    i, ai, aj, ak, al, t1, t2, t3;
     rvec   r_ij, r_kj, r_kl, m, n;
-    real   sign, th, costh, ph;
+    real   th, costh, ph;
     t_pbc  pbc;
 
     if (bPBC)
@@ -314,7 +314,7 @@ static void calc_angles_dihs(t_params *ang, t_params *dih, const rvec x[], gmx_b
         ak = dih->param[i].ak();
         al = dih->param[i].al();
         ph = RAD2DEG*dih_angle(x[ai], x[aj], x[ak], x[al], bPBC ? &pbc : nullptr,
-                               r_ij, r_kj, r_kl, m, n, &sign, &t1, &t2, &t3);
+                               r_ij, r_kj, r_kl, m, n, &t1, &t2, &t3);
         if (debug)
         {
             fprintf(debug, "X2TOP: ai=%3d aj=%3d ak=%3d al=%3d r_ij=%8.3f r_kj=%8.3f r_kl=%8.3f ph=%8.3f\n",
index 90512069e9bede4a7b9c5681af366580db7e72dd..f752d8b2c2f8e591215791b98952e4f79b5bf1f3 100644 (file)
@@ -1389,20 +1389,18 @@ real quartic_angles(int nbonds,
 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
                const t_pbc *pbc,
                rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n,
-               real *sign, int *t1, int *t2, int *t3)
+               int *t1, int *t2, int *t3)
 {
-    real ipr, phi;
-
     *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /*  3        */
     *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /*  3               */
     *t3 = pbc_rvec_sub(pbc, xk, xl, r_kl); /*  3               */
 
     cprod(r_ij, r_kj, m);                  /*  9        */
     cprod(r_kj, r_kl, n);                  /*  9               */
-    phi     = gmx_angle(m, n);             /* 49 (assuming 25 for atan2) */
-    ipr     = iprod(r_ij, n);              /*  5        */
-    (*sign) = (ipr < 0.0) ? -1.0 : 1.0;
-    phi     = (*sign)*phi;                 /*  1               */
+    real phi  = gmx_angle(m, n);           /* 49 (assuming 25 for atan2) */
+    real ipr  = iprod(r_ij, n);            /*  5        */
+    real sign = (ipr < 0.0) ? -1.0 : 1.0;
+    phi       = sign*phi;                  /*  1               */
     /* 82 TOTAL        */
     return phi;
 }
@@ -1740,7 +1738,7 @@ real pdihs(int nbonds,
     int  i, type, ai, aj, ak, al;
     int  t1, t2, t3;
     rvec r_ij, r_kj, r_kl, m, n;
-    real phi, sign, ddphi, vpd, vtot;
+    real phi, ddphi, vpd, vtot;
 
     vtot = 0.0;
 
@@ -1753,7 +1751,7 @@ real pdihs(int nbonds,
         al   = forceatoms[i++];
 
         phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
-                        &sign, &t1, &t2, &t3);  /*  84      */
+                        &t1, &t2, &t3);  /*  84      */
         *dvdlambda += dopdihs(forceparams[type].pdihs.cpA,
                               forceparams[type].pdihs.cpB,
                               forceparams[type].pdihs.phiA,
@@ -1801,7 +1799,7 @@ pdihs_noener(int nbonds,
     int  i, type, ai, aj, ak, al;
     int  t1, t2, t3;
     rvec r_ij, r_kj, r_kl, m, n;
-    real phi, sign, ddphi_tot, ddphi;
+    real phi, ddphi_tot, ddphi;
 
     for (i = 0; (i < nbonds); )
     {
@@ -1811,7 +1809,7 @@ pdihs_noener(int nbonds,
         al   = forceatoms[i+4];
 
         phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
-                        &sign, &t1, &t2, &t3);
+                        &t1, &t2, &t3);
 
         ddphi_tot = 0;
 
@@ -2099,7 +2097,7 @@ real idihs(int nbonds,
 {
     int  i, type, ai, aj, ak, al;
     int  t1, t2, t3;
-    real phi, phi0, dphi0, ddphi, sign, vtot;
+    real phi, phi0, dphi0, ddphi, vtot;
     rvec r_ij, r_kj, r_kl, m, n;
     real L1, kk, dp, dp2, kA, kB, pA, pB, dvdl_term;
 
@@ -2115,7 +2113,7 @@ real idihs(int nbonds,
         al   = forceatoms[i++];
 
         phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
-                        &sign, &t1, &t2, &t3);  /*  84         */
+                        &t1, &t2, &t3);  /*  84                */
 
         /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
          * force changes if we just apply a normal harmonic.
@@ -2292,7 +2290,7 @@ real dihres(int nbonds,
     real vtot = 0;
     int  ai, aj, ak, al, i, k, type, t1, t2, t3;
     real phi0A, phi0B, dphiA, dphiB, kfacA, kfacB, phi0, dphi, kfac;
-    real phi, ddphi, ddp, ddp2, dp, sign, d2r, L1;
+    real phi, ddphi, ddp, ddp2, dp, d2r, L1;
     rvec r_ij, r_kj, r_kl, m, n;
 
     L1 = 1.0-lambda;
@@ -2321,7 +2319,7 @@ real dihres(int nbonds,
         kfac  = L1*kfacA + lambda*kfacB;
 
         phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
-                        &sign, &t1, &t2, &t3);
+                        &t1, &t2, &t3);
         /* 84 flops */
 
         if (debug)
@@ -2724,7 +2722,7 @@ real rbdihs(int nbonds,
     real       parmB[NR_RBDIHS];
     real       parm[NR_RBDIHS];
     real       cos_phi, phi, rbp, rbpBA;
-    real       v, sign, ddphi, sin_phi;
+    real       v, ddphi, sin_phi;
     real       cosfac, vtot;
     real       L1        = 1.0-lambda;
     real       dvdl_term = 0;
@@ -2739,7 +2737,7 @@ real rbdihs(int nbonds,
         al   = forceatoms[i++];
 
         phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
-                        &sign, &t1, &t2, &t3);  /*  84         */
+                        &t1, &t2, &t3);  /*  84                */
 
         /* Change to polymer convention */
         if (phi < c0)
@@ -2876,8 +2874,8 @@ cmap_dihs(int nbonds,
     int         pos1, pos2, pos3, pos4;
 
     real        ty[4], ty1[4], ty2[4], ty12[4], tc[16], tx[16];
-    real        phi1, cos_phi1, sin_phi1, sign1, xphi1;
-    real        phi2, cos_phi2, sin_phi2, sign2, xphi2;
+    real        phi1, cos_phi1, sin_phi1, xphi1;
+    real        phi2, cos_phi2, sin_phi2, xphi2;
     real        dx, xx, tt, tu, e, df1, df2, vtot;
     real        ra21, rb21, rg21, rg1, rgr1, ra2r1, rb2r1, rabr1;
     real        ra22, rb22, rg22, rg2, rgr2, ra2r2, rb2r2, rabr2;
@@ -2928,7 +2926,7 @@ cmap_dihs(int nbonds,
         a1l   = al;
 
         phi1  = dih_angle(x[a1i], x[a1j], x[a1k], x[a1l], pbc, r1_ij, r1_kj, r1_kl, m1, n1,
-                          &sign1, &t11, &t21, &t31);  /* 84 */
+                          &t11, &t21, &t31);  /* 84 */
 
         cos_phi1 = std::cos(phi1);
 
@@ -2989,7 +2987,7 @@ cmap_dihs(int nbonds,
         a2l   = am;
 
         phi2  = dih_angle(x[a2i], x[a2j], x[a2k], x[a2l], pbc, r2_ij, r2_kj, r2_kl, m2, n2,
-                          &sign2, &t12, &t22, &t32); /* 84 */
+                          &t12, &t22, &t32); /* 84 */
 
         cos_phi2 = std::cos(phi2);
 
@@ -3784,7 +3782,7 @@ real tab_dihs(int nbonds,
     int  i, type, ai, aj, ak, al, table;
     int  t1, t2, t3;
     rvec r_ij, r_kj, r_kl, m, n;
-    real phi, sign, ddphi, vpd, vtot;
+    real phi, ddphi, vpd, vtot;
 
     vtot = 0.0;
     for (i = 0; (i < nbonds); )
@@ -3796,7 +3794,7 @@ real tab_dihs(int nbonds,
         al   = forceatoms[i++];
 
         phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
-                        &sign, &t1, &t2, &t3);  /*  84  */
+                        &t1, &t2, &t3);  /*  84  */
 
         table = forceparams[type].tab.table;
 
index 0c86e8a127f55a1ec02e8d9409409fa54e2b37bc..2b5d737296b9e0a28d1745cf8b5d67a8101e4f00 100644 (file)
@@ -73,7 +73,6 @@ real bond_angle(const rvec xi, const rvec xj, const rvec xk,
 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
                const struct t_pbc *pbc,
                rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n, /* out */
-               real *sign,
                int *t1, int *t2, int *t3);
 
 /*! \brief Do an update of the forces for dihedral potentials */
index e7b72be85ffadc8f44472926e23c359025141356..d3d80bd1585a94d7e0043f2c1fd601dd0d36c19d 100644 (file)
@@ -106,17 +106,16 @@ class BondedTest : public ::testing::Test
         void testDihedralAngle(int epbc)
         {
             rvec  r_ij, r_kj, r_kl, m, n;
-            real  cosine_angle, angle;
+            real  angle;
             int   t1, t2, t3;
             t_pbc pbc;
 
             set_pbc(&pbc, epbc, box);
             angle = dih_angle(x[0], x[1], x[2], x[3], &pbc,
-                              r_ij, r_kj, r_kl, m, n, &cosine_angle,
+                              r_ij, r_kj, r_kl, m, n,
                               &t1, &t2, &t3);
 
             checker_.checkReal(angle, "angle");
-            checker_.checkReal(cosine_angle, "cosine_angle");
             checker_.checkInteger(t1, "t1");
             checker_.checkInteger(t2, "t2");
             checker_.checkInteger(t3, "t3");
index 5458d88aaf04da95415c24decea7e734c90f3caf..4cbb82d0f65099f1d9ec3ffa8d6ac3c4f9251810 100644 (file)
@@ -2,7 +2,6 @@
 <?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
 <ReferenceData>
   <Real Name="angle">-1.5707964</Real>
-  <Real Name="cosine_angle">-1</Real>
   <Int Name="t1">22</Int>
   <Int Name="t2">22</Int>
   <Int Name="t3">22</Int>
index 3bf5d9b43b39d2d9a5490daf41651f6d09d44ef6..36301d1e061418a794b900a4045852d37bebeb10 100644 (file)
@@ -2,7 +2,6 @@
 <?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
 <ReferenceData>
   <Real Name="angle">-1.5707964</Real>
-  <Real Name="cosine_angle">-1</Real>
   <Int Name="t1">22</Int>
   <Int Name="t2">17</Int>
   <Int Name="t3">23</Int>
index 8941df4caf5947ec0e9adc2e5c0aa37d79d92c63..2852780783ad416d58cde8ac6139d458b97ed661 100644 (file)
@@ -2,7 +2,6 @@
 <?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
 <ReferenceData>
   <Real Name="angle">1.5707964</Real>
-  <Real Name="cosine_angle">1</Real>
   <Int Name="t1">37</Int>
   <Int Name="t2">17</Int>
   <Int Name="t3">23</Int>