Remove unused thole polarization rfac parameter
authorJoe Jordan <ejjordan12@gmail.com>
Wed, 3 Nov 2021 06:24:23 +0000 (06:24 +0000)
committerMark Abraham <mark.j.abraham@gmail.com>
Wed, 3 Nov 2021 06:24:23 +0000 (06:24 +0000)
src/gromacs/fileio/tpxio.cpp
src/gromacs/gmxpreprocess/convparm.cpp
src/gromacs/listed_forces/tests/bonded.cpp
src/gromacs/topology/idef.cpp
src/gromacs/topology/idef.h

index d190d663ed65439cde1a942b82d33e4199ba6572..f8b0ced58defbeb726d921d29b5b9480923a4e01 100644 (file)
@@ -140,6 +140,7 @@ enum tpxv
     tpxv_TransformationPullCoord,     /**< Support for transformation pull coordinates */
     tpxv_SoftcoreGapsys,              /**< Added gapsys softcore function */
     tpxv_ReaddedConstantAcceleration, /**< Re-added support for constant acceleration NEMD. */
+    tpxv_RemoveTholeRfac,             /**< Remove unused rfac parameter from thole listed force */
     tpxv_Count                        /**< the total number of tpxv versions */
 };
 
@@ -1899,7 +1900,12 @@ static void do_iparams(gmx::ISerializer* serializer, t_functype ftype, t_iparams
             serializer->doReal(&iparams->thole.a);
             serializer->doReal(&iparams->thole.alpha1);
             serializer->doReal(&iparams->thole.alpha2);
-            serializer->doReal(&iparams->thole.rfac);
+            if (file_version < tpxv_RemoveTholeRfac)
+            {
+                real noRfac = 0;
+                serializer->doReal(&noRfac);
+            }
+
             break;
         case F_LJ:
             serializer->doReal(&iparams->lj.c6);
index 37b59f1eb5bce38d43afe26566dab8f77d221dcf..e8e5fc897f5aa85f7126ffa5261c383158b6327b 100644 (file)
@@ -258,14 +258,6 @@ static int assign_param(t_functype                ftype,
             newparam->thole.a      = old[0];
             newparam->thole.alpha1 = old[1];
             newparam->thole.alpha2 = old[2];
-            if ((old[1] > 0) && (old[2] > 0))
-            {
-                newparam->thole.rfac = old[0] * gmx::invsixthroot(old[1] * old[2]);
-            }
-            else
-            {
-                newparam->thole.rfac = 1;
-            }
             break;
         case F_BHAM:
             newparam->bham.a = old[0];
index 6abade7c01b027d9b54f6076989a4909e69215ee..04bf2e1eb7a9adf4951edb18783b09fd838871c4 100644 (file)
@@ -451,17 +451,15 @@ public:
      * \param[in] a      Thole factor
      * \param[in] alpha1 Polarizability 1 (nm^3)
      * \param[in] alpha2 Polarizability 2 (nm^3)
-     * \param[in] rfac   Distance factor
      * \return The structure itself.
      */
-    iListInput setTholePolarization(real a, real alpha1, real alpha2, real rfac)
+    iListInput setTholePolarization(real a, real alpha1, real alpha2)
     {
         ftype                = F_THOLE_POL;
         fep                  = false;
         iparams.thole.a      = a;
         iparams.thole.alpha1 = alpha1;
         iparams.thole.alpha2 = alpha2;
-        iparams.thole.rfac   = rfac;
         return *this;
     }
     /*! \brief Set parameters for Water Polarization
@@ -708,7 +706,7 @@ std::vector<iListInput> c_InputDihs = {
 std::vector<iListInput> c_InputPols = {
     { iListInput(2e-5, 1e-8).setPolarization(0.12) },
     { iListInput(2e-3, 1e-8).setAnharmPolarization(0.0013, 0.02, 1235.6) },
-    { iListInput(1.4e-3, 1e-8).setTholePolarization(0.26, 0.07, 0.09, 1.6) },
+    { iListInput(1.4e-3, 1e-8).setTholePolarization(0.26, 0.07, 0.09) },
     { iListInput(2e-3, 1e-8).setWaterPolarization(0.001, 0.0012, 0.0016, 0.095, 0.15, 0.02) },
 };
 
index b54d73c4f5cb11c8ffa1e124e4ca04f30b2fc902..543a0f41135fc173604f8325950a2025100a5ed0 100644 (file)
@@ -180,11 +180,10 @@ void printInteractionParameters(gmx::TextWriter* writer, t_functype ftype, const
                                        iparams.anharm_polarize.khyp);
             break;
         case F_THOLE_POL:
-            writer->writeLineFormatted("a=%15.8e, alpha1=%15.8e, alpha2=%15.8e, rfac=%15.8e",
+            writer->writeLineFormatted("a=%15.8e, alpha1=%15.8e, alpha2=%15.8e",
                                        iparams.thole.a,
                                        iparams.thole.alpha1,
-                                       iparams.thole.alpha2,
-                                       iparams.thole.rfac);
+                                       iparams.thole.alpha2);
             break;
         case F_WATER_POL:
             writer->writeLineFormatted(
index 7c68e2c606b629d4b5d540efde9f89bd1a8d71f4..60cf6bee0e122ce330e3cd199655bf9be1552800 100644 (file)
@@ -112,7 +112,7 @@ typedef union t_iparams
     } wpol;
     struct
     {
-        real a, alpha1, alpha2, rfac;
+        real a, alpha1, alpha2;
     } thole;
     struct
     {