Add gapsys softcore function.
[alexxy/gromacs.git] / src / gromacs / listed_forces / tests / pairs.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2021, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \file
36  * \brief Implements test of 1-4 interactions
37  *
38  * This test is copied from the bonded interactions test and slightly
39  * modified since 'do_pairs' takes a different set of arguments than
40  * 'calculateSimpleBond'. To keep the test setup uncluttered this test is
41  * therefore not merged into the bonded test but implemented standalone.
42  *
43  * The test setup consists of 2 atom pairs that are tested in an fep setting
44  * (vanishing charge and lennard-jones parameters of one atom) and without
45  * fep. Placement of the atoms in the box is such that shift-forces and pbc
46  * paths in do_pairs are covered.
47  *
48  * \author Sebastian Kehl <sebastian.kehl@mpcdf.mpg.de>
49  * \ingroup module_listed_forces
50  */
51 #include "gmxpre.h"
52
53 #include "gromacs/listed_forces/bonded.h"
54
55 #include <cmath>
56
57 #include <memory>
58 #include <unordered_map>
59
60 #include <gtest/gtest.h>
61
62 #include "gromacs/listed_forces/listed_forces.h"
63 #include "gromacs/listed_forces/pairs.h"
64 #include "gromacs/math/paddedvector.h"
65 #include "gromacs/math/units.h"
66 #include "gromacs/math/vec.h"
67 #include "gromacs/math/vectypes.h"
68 #include "gromacs/mdtypes/mdatom.h"
69 #include "gromacs/mdtypes/simulation_workload.h"
70 #include "gromacs/mdtypes/enerdata.h"
71 #include "gromacs/mdtypes/forcerec.h"
72 #include "gromacs/mdtypes/inputrec.h"
73 #include "gromacs/mdtypes/interaction_const.h"
74 #include "gromacs/mdtypes/nblist.h"
75 #include "gromacs/tables/forcetable.h"
76 #include "gromacs/pbcutil/ishift.h"
77 #include "gromacs/pbcutil/pbc.h"
78 #include "gromacs/topology/idef.h"
79 #include "gromacs/utility/enumerationhelpers.h"
80 #include "gromacs/utility/strconvert.h"
81 #include "gromacs/utility/stringstream.h"
82 #include "gromacs/utility/textwriter.h"
83 #include "gromacs/utility/fatalerror.h"
84
85 #include "testutils/refdata.h"
86 #include "testutils/testasserts.h"
87
88 namespace gmx
89 {
90 namespace test
91 {
92 namespace
93 {
94
95 //! Number of atoms used in these tests.
96 constexpr int c_numAtoms = 3;
97
98 /*! \brief Output from pairs kernels
99  *
100  */
101 struct OutputQuantities
102 {
103     OutputQuantities(int energyGroup) :
104         energy(energyGroup), dvdLambda(static_cast<int>(FreeEnergyPerturbationCouplingType::Count), 0.0)
105     {
106     }
107
108     //! Energy of this interaction
109     gmx_grppairener_t energy;
110     //! Derivative with respect to lambda
111     std::vector<real> dvdLambda;
112     //! Shift vectors
113     rvec fShift[gmx::detail::c_numIvecs] = { { 0 } };
114     //! Forces
115     alignas(GMX_REAL_MAX_SIMD_WIDTH * sizeof(real)) rvec4 f[c_numAtoms] = { { 0 } };
116 };
117
118 /*! \brief Utility to check the output from pairs tests
119  *
120  * \param[in] checker Reference checker
121  * \param[in] output  The output from the test to check
122  * \param[in] bondedKernelFlavor  Flavor for determining what output to check
123  * \param[in] functionType type of the interaction
124  */
125 void checkOutput(TestReferenceChecker*    checker,
126                  const OutputQuantities&  output,
127                  const BondedKernelFlavor bondedKernelFlavor,
128                  const int                functionType)
129 {
130     if (computeEnergy(bondedKernelFlavor))
131     {
132         switch (functionType)
133         {
134             case F_LJ14:
135             case F_LJC14_Q:
136                 checker->checkReal(output.energy.energyGroupPairTerms[NonBondedEnergyTerms::Coulomb14][0],
137                                    "Epot Coulomb14");
138                 checker->checkReal(output.energy.energyGroupPairTerms[NonBondedEnergyTerms::LJ14][0],
139                                    "Epot LJ14");
140                 break;
141             case F_LJC_PAIRS_NB:
142                 checker->checkReal(output.energy.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR][0],
143                                    "Epot Coulomb14");
144                 checker->checkReal(output.energy.energyGroupPairTerms[NonBondedEnergyTerms::LJSR][0],
145                                    "Epot LJ14");
146                 break;
147             default: gmx_fatal(FARGS, "Unknown function type %d in do_nonbonded14", functionType);
148         }
149         checker->checkReal(output.dvdLambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
150                            "dVdlCoul ");
151         checker->checkReal(output.dvdLambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
152                            "dVdlVdw ");
153     }
154     checker->checkSequence(std::begin(output.f), std::end(output.f), "Forces");
155 }
156
157 /* \brief Utility class to setup forcerec and interaction parameters
158  *
159  * Data is only initialized as necessary for the 1-4 interactions to work!
160  */
161 class ForcerecHelper
162 {
163 public:
164     ForcerecHelper()
165     {
166         fepVals_.sc_alpha         = 0.3;
167         fepVals_.sc_power         = 1;
168         fepVals_.sc_r_power       = 6.0;
169         fepVals_.sc_sigma         = 0.85;
170         fepVals_.sc_sigma_min     = 0.3;
171         fepVals_.bScCoul          = true;
172         fepVals_.softcoreFunction = SoftcoreType::Beutler;
173
174         fr_.ic = std::make_unique<interaction_const_t>();
175         // set data in ic
176         fr_.ic->softCoreParameters = std::make_unique<interaction_const_t::SoftCoreParameters>(fepVals_);
177
178         // set data in fr
179         real tableRange = 2.9;
180         fr_.pairsTable = make_tables(nullptr, fr_.ic.get(), nullptr, tableRange, GMX_MAKETABLES_14ONLY);
181         fr_.efep       = haveFep_;
182         fr_.fudgeQQ          = 0.5;
183         fr_.use_simd_kernels = useSimd_;
184         fr_.bMolPBC          = haveMolPBC_;
185     }
186
187     t_forcerec* get() { return &fr_; }
188
189 private:
190     FreeEnergyPerturbationType haveFep_    = FreeEnergyPerturbationType::No;
191     bool                       useSimd_    = false;
192     bool                       haveMolPBC_ = false;
193
194     t_lambda   fepVals_;
195     t_forcerec fr_;
196 };
197
198 ForcerecHelper frHelper;
199
200
201 /*! \brief Input structure for listed forces tests
202  */
203 struct ListInput
204 {
205 public:
206     //! Function type
207     int fType = -1;
208     //! do fep
209     bool fep = false;
210     //! Tolerance for float evaluation
211     float floatToler = 1e-6;
212     //! Tolerance for double evaluation
213     double doubleToler = 1e-8;
214     //! Interaction parameters
215     t_iparams iparams = { { 0 } };
216
217     //! Constructor
218     ListInput() {}
219
220     /*! \brief Constructor with tolerance
221      *
222      * \param[in] ftol Single precision tolerance
223      * \param[in] dtol Double precision tolerance
224      */
225     ListInput(float ftol, double dtol)
226     {
227         floatToler  = ftol;
228         doubleToler = dtol;
229     }
230
231     /*! \brief Set parameters for lj14 interaction
232      *
233      * Fep is used if either c6A != c6B or c12A != c12B.
234      *
235      * \param[in] c6A  lj-c6 of state A
236      * \param[in] c12A lj-c12 of state A
237      * \param[in] c6B  lj-c6 of state B
238      * \param[in] c12B lj-c12 of state B
239      */
240     ListInput setLj14Interaction(real c6A, real c12A, real c6B, real c12B)
241     {
242         fType             = F_LJ14;
243         fep               = (c6A != c6B || c12A != c12B);
244         iparams.lj14.c6A  = c6A;
245         iparams.lj14.c12A = c12A;
246         iparams.lj14.c6B  = c6B;
247         iparams.lj14.c12B = c12B;
248
249         return *this;
250     }
251
252     /*! \brief Set parameters for ljc14 interaction
253      *
254      * \param[in] qi     charge i
255      * \param[in] qj     charge j
256      * \param[in] c6     lj-c6
257      * \param[in] c12    lj-c12
258      * \param[in] fudgeQ fudge factor
259      */
260     ListInput setLjc14Interaction(real qi, real qj, real c6, real c12, real fudgeQ)
261     {
262         fType             = F_LJC14_Q;
263         iparams.ljc14.qi  = qi;
264         iparams.ljc14.qj  = qj;
265         iparams.ljc14.c6  = c6;
266         iparams.ljc14.c12 = c12;
267         iparams.ljc14.fqq = fudgeQ;
268
269         return *this;
270     }
271
272     /*! \brief Set parameters for ljcnb interaction
273      *
274      * \param[in] qi  charge i
275      * \param[in] qj  charge j
276      * \param[in] c6  lj-c6
277      * \param[in] c12 lj-c12
278      */
279     ListInput setLjcnbInteraction(real qi, real qj, real c6, real c12)
280     {
281         fType             = F_LJC_PAIRS_NB;
282         iparams.ljcnb.qi  = qi;
283         iparams.ljcnb.qj  = qj;
284         iparams.ljcnb.c6  = c6;
285         iparams.ljcnb.c12 = c12;
286
287         return *this;
288     }
289 };
290
291 class ListedForcesPairsTest :
292     public ::testing::TestWithParam<std::tuple<ListInput, PaddedVector<RVec>, PbcType>>
293 {
294 protected:
295     matrix               box_;
296     t_pbc                pbc_;
297     PaddedVector<RVec>   x_;
298     PbcType              pbcType_;
299     ListInput            input_;
300     TestReferenceData    refData_;
301     TestReferenceChecker checker_;
302
303     ListedForcesPairsTest() : checker_(refData_.rootChecker())
304     {
305         input_   = std::get<0>(GetParam());
306         x_       = std::get<1>(GetParam());
307         pbcType_ = std::get<2>(GetParam());
308         clear_mat(box_);
309         box_[0][0] = box_[1][1] = box_[2][2] = 1.0;
310         set_pbc(&pbc_, pbcType_, box_);
311
312         FloatingPointTolerance tolerance = relativeToleranceAsPrecisionDependentFloatingPoint(
313                 1.0, input_.floatToler, input_.doubleToler);
314
315         checker_.setDefaultTolerance(tolerance);
316     }
317
318     void testOneIfunc(TestReferenceChecker* checker, const real lambda, const SoftcoreType softcoreType)
319     {
320         SCOPED_TRACE(std::string("Testing PBC type: ") + c_pbcTypeNames[pbcType_]);
321
322         // 'definition of pairs' is a concatenation of #npairs (here 2)
323         // 'nAtomsPerPair+1'-tuples (fType a_0 a_i ... a_nAtomsPerPair)
324         std::vector<t_iatom> iatoms = { 0, 1, 2, 0, 0, 2 };
325
326         std::vector<int>            ddgatindex = { 0, 1, 2 };
327         std::vector<real>           chargeA    = { 1.0, -0.5, -0.5 };
328         std::vector<real>           chargeB    = { 0.0, 0.0, 0.0 };
329         std::vector<unsigned short> egrp       = { 0, 0, 0 };
330         t_mdatoms                   mdatoms    = { 0 };
331
332         mdatoms.chargeA = chargeA.data();
333         mdatoms.chargeB = chargeB.data();
334         mdatoms.cENER   = egrp.data();
335         // nPerturbed is not decisive for fep to be used; it is overruled by
336         // other conditions in do_pairs_general; just here to not segfault
337         // upon query
338         mdatoms.nPerturbed = 0;
339
340         t_forcerec* fr = frHelper.get();
341         fr->efep = input_.fep ? FreeEnergyPerturbationType::Yes : FreeEnergyPerturbationType::No;
342         fr->ic->softCoreParameters->softcoreType = softcoreType;
343         if (pbcType_ != PbcType::No)
344         {
345             fr->bMolPBC = true;
346         }
347
348         std::vector<BondedKernelFlavor> flavors = { BondedKernelFlavor::ForcesAndVirialAndEnergy };
349
350         if (!input_.fep || lambda == 0)
351         {
352             fr->use_simd_kernels = true;
353             flavors.push_back(BondedKernelFlavor::ForcesSimdWhenAvailable);
354         }
355
356         for (const auto flavor : flavors)
357         {
358             SCOPED_TRACE("Testing bonded kernel flavor: " + c_bondedKernelFlavorStrings[flavor]);
359
360             StepWorkload stepWork;
361             stepWork.computeVirial = computeVirial(flavor);
362             stepWork.computeEnergy = computeEnergy(flavor);
363
364             bool havePerturbedInteractions = input_.fep;
365             if (flavor == BondedKernelFlavor::ForcesSimdWhenAvailable)
366             {
367                 havePerturbedInteractions = false;
368             }
369
370             int numEnergyTerms      = static_cast<int>(NonBondedEnergyTerms::Count);
371             int numFepCouplingTerms = static_cast<int>(FreeEnergyPerturbationCouplingType::Count);
372             OutputQuantities  output(numEnergyTerms);
373             std::vector<real> lambdas(numFepCouplingTerms, lambda);
374
375             do_pairs(input_.fType,
376                      iatoms.size(),
377                      iatoms.data(),
378                      &input_.iparams,
379                      as_rvec_array(x_.data()),
380                      output.f,
381                      output.fShift,
382                      &pbc_,
383                      lambdas.data(),
384                      output.dvdLambda.data(),
385                      gmx::arrayRefFromArray(mdatoms.chargeA, mdatoms.nr),
386                      gmx::arrayRefFromArray(mdatoms.chargeB, mdatoms.nr),
387                      gmx::arrayRefFromArray(mdatoms.bPerturbed, mdatoms.nr),
388                      gmx::arrayRefFromArray(mdatoms.cENER, mdatoms.nr),
389                      mdatoms.nPerturbed,
390                      fr,
391                      havePerturbedInteractions,
392                      stepWork,
393                      &output.energy,
394                      ddgatindex.data());
395
396             checkOutput(checker, output, flavor, input_.fType);
397             auto shiftForcesChecker = checker->checkCompound("Shift-Forces", "Shift-forces");
398
399             if (computeVirial(flavor))
400             {
401                 shiftForcesChecker.checkVector(output.fShift[gmx::c_centralShiftIndex], "Central");
402             }
403             else
404             {
405                 // Permit omitting to compare shift forces with
406                 // reference data when that is useless.
407                 shiftForcesChecker.disableUnusedEntriesCheck();
408             }
409         }
410     }
411
412     void testIfunc()
413     {
414         TestReferenceChecker thisChecker =
415                 checker_.checkCompound("FunctionType", interaction_function[input_.fType].name)
416                         .checkCompound("FEP", (input_.fep ? "Yes" : "No"));
417
418         if (input_.fep)
419         {
420             const int numLambdas = 3;
421             for (int i = 0; i < numLambdas; ++i)
422             {
423                 const real lambda = i / (numLambdas - 1.0);
424                 for (SoftcoreType c : EnumerationWrapper<SoftcoreType>{})
425                 {
426                     auto lambdaChecker = thisChecker.checkCompound("Lambda", toString(lambda));
427                     auto softcoreChecker = lambdaChecker.checkCompound("Sofcore", enumValueToString(c));
428                     testOneIfunc(&softcoreChecker, lambda, c);
429                 }
430             }
431         }
432         else
433         {
434             testOneIfunc(&thisChecker, 0.0, SoftcoreType::None);
435         }
436     }
437 };
438
439 TEST_P(ListedForcesPairsTest, Ifunc)
440 {
441     testIfunc();
442 }
443
444 //! Function types for testing 1-4 interaction. Add new terms at the end.
445 std::vector<ListInput> c_14Interaction = {
446     { ListInput(1e-5, 1e-7).setLj14Interaction(0.001458, 1.0062882e-6, 0.0, 0.0) },
447     { ListInput(1e-5, 1e-7).setLj14Interaction(0.001458, 1.0062882e-6, 0.001458, 1.0062882e-6) },
448     { ListInput(1e-5, 1e-7).setLjc14Interaction(1.0, -1.0, 0.001458, 1.0062882e-6, 0.5) },
449     { ListInput(1e-5, 1e-7).setLjcnbInteraction(1.0, -1.0, 0.001458, 1.0062882e-6) }
450 };
451
452 //! PBC values for testing
453 std::vector<PbcType> c_pbcForTests = { PbcType::No, PbcType::XY, PbcType::Xyz };
454
455 /*! \brief Coordinates for testing 1-4 interaction
456  *
457  * Define coordinates for 3 atoms here, which will be used in 2 interactions.
458  */
459 std::vector<PaddedVector<RVec>> c_coordinatesFor14Interaction = {
460     { { 0.0, 0.0, 0.0 }, { 1.0, 1.0, 1.0 }, { 1.1, 1.2, 1.3 } }
461 };
462
463 INSTANTIATE_TEST_SUITE_P(14Interaction,
464                          ListedForcesPairsTest,
465                          ::testing::Combine(::testing::ValuesIn(c_14Interaction),
466                                             ::testing::ValuesIn(c_coordinatesFor14Interaction),
467                                             ::testing::ValuesIn(c_pbcForTests)));
468
469 } // namespace
470
471 } // namespace test
472
473 } // namespace gmx