Adapt testcases for gapsys softcore
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / tests / nb_free_energy.cpp
index 01af68c47127bb8745a65f8e6e197f0ceda8ea38..e6d42eaab04780cda11853dc4cf8e82aa9d3f66f 100644 (file)
@@ -235,7 +235,7 @@ public:
         fepVals_.sc_sigma                = 0.3;
         fepVals_.sc_sigma_min            = 0.3;
         fepVals_.bScCoul                 = true;
-        fepVals_.scGapsysScaleLinpointLJ = 0.3;
+        fepVals_.scGapsysScaleLinpointLJ = 0.85;
         fepVals_.scGapsysScaleLinpointQ  = 0.3;
         fepVals_.scGapsysSigmaLJ         = 0.3;
         fepVals_.softcoreFunction        = SoftcoreType::Beutler;
@@ -252,8 +252,18 @@ public:
         ljPmeC6Grid_ = makeLJPmeC6GridCorrectionParameters(idef.atnr, idef.iparams, LongRangeVdW::Geom);
     }
 
-    void setSoftcoreAlpha(const real scAlpha) { fepVals_.sc_alpha = scAlpha; }
+    void setSoftcoreAlpha(const real scBeutlerAlphaOrGapsysLinpointScaling)
+    {
+        fepVals_.sc_alpha                = scBeutlerAlphaOrGapsysLinpointScaling;
+        fepVals_.scGapsysScaleLinpointLJ = scBeutlerAlphaOrGapsysLinpointScaling;
+        fepVals_.scGapsysScaleLinpointQ  = scBeutlerAlphaOrGapsysLinpointScaling;
+    }
     void setSoftcoreCoulomb(const bool scCoulomb) { fepVals_.bScCoul = scCoulomb; }
+    void setSoftcoreType(const SoftcoreType softcoreType)
+    {
+        fepVals_.softcoreFunction = softcoreType;
+    }
+
 
     //! get forcerec data as wanted by the nonbonded kernel
     void getForcerec(t_forcerec* fr)
@@ -396,7 +406,7 @@ public:
 };
 
 class NonbondedFepTest :
-    public ::testing::TestWithParam<std::tuple<ListInput, PaddedVector<RVec>, real, real, bool>>
+    public ::testing::TestWithParam<std::tuple<SoftcoreType, ListInput, PaddedVector<RVec>, real, real, bool>>
 {
 protected:
     PaddedVector<RVec>   x_;
@@ -404,16 +414,18 @@ protected:
     real                 lambda_;
     real                 softcoreAlpha_;
     bool                 softcoreCoulomb_;
+    SoftcoreType         softcoreType_;
     TestReferenceData    refData_;
     TestReferenceChecker checker_;
 
     NonbondedFepTest() : checker_(refData_.rootChecker())
     {
-        input_           = std::get<0>(GetParam());
-        x_               = std::get<1>(GetParam());
-        lambda_          = std::get<2>(GetParam());
-        softcoreAlpha_   = std::get<3>(GetParam());
-        softcoreCoulomb_ = std::get<4>(GetParam());
+        softcoreType_    = std::get<0>(GetParam());
+        input_           = std::get<1>(GetParam());
+        x_               = std::get<2>(GetParam());
+        lambda_          = std::get<3>(GetParam());
+        softcoreAlpha_   = std::get<4>(GetParam());
+        softcoreCoulomb_ = std::get<5>(GetParam());
 
         // Note that the reference data for Ewald type interactions has been generated
         // with accurate analytical approximations for the long-range corrections.
@@ -428,6 +440,7 @@ protected:
     {
         input_.frHelper.setSoftcoreAlpha(softcoreAlpha_);
         input_.frHelper.setSoftcoreCoulomb(softcoreCoulomb_);
+        input_.frHelper.setSoftcoreType(softcoreType_);
 
         // get forcerec and interaction_const
         t_forcerec fr;
@@ -496,9 +509,10 @@ std::vector<ListInput> c_interaction = {
 };
 
 //! test parameters
-std::vector<real> c_fepLambdas      = { 0.0, 0.5, 1.0 };
-std::vector<real> c_softcoreAlphas  = { 0.0, 0.3 };
-std::vector<bool> c_softcoreCoulomb = { true, false };
+std::vector<real>         c_fepLambdas                                  = { 0.0, 0.5, 1.0 };
+std::vector<real>         c_softcoreBeutlerAlphaOrGapsysLinpointScaling = { 0.0, 0.3 };
+std::vector<bool>         c_softcoreCoulomb                             = { true, false };
+std::vector<SoftcoreType> c_softcoreType = { SoftcoreType::Beutler, SoftcoreType::Gapsys };
 
 //! Coordinates for testing
 std::vector<PaddedVector<RVec>> c_coordinates = {
@@ -507,10 +521,11 @@ std::vector<PaddedVector<RVec>> c_coordinates = {
 
 INSTANTIATE_TEST_SUITE_P(NBInteraction,
                          NonbondedFepTest,
-                         ::testing::Combine(::testing::ValuesIn(c_interaction),
+                         ::testing::Combine(::testing::ValuesIn(c_softcoreType),
+                                            ::testing::ValuesIn(c_interaction),
                                             ::testing::ValuesIn(c_coordinates),
                                             ::testing::ValuesIn(c_fepLambdas),
-                                            ::testing::ValuesIn(c_softcoreAlphas),
+                                            ::testing::ValuesIn(c_softcoreBeutlerAlphaOrGapsysLinpointScaling),
                                             ::testing::ValuesIn(c_softcoreCoulomb)));
 
 } // namespace