137ebfd4407033245652566b95f8a9e65357571f
[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.3;
170         fepVals_.sc_sigma_min = 0.3;
171         fepVals_.bScCoul      = true;
172
173         fr_.ic = std::make_unique<interaction_const_t>();
174         // set data in ic
175         fr_.ic->softCoreParameters = std::make_unique<interaction_const_t::SoftCoreParameters>(fepVals_);
176
177         // set data in fr
178         real tableRange = 2.9;
179         fr_.pairsTable = make_tables(nullptr, fr_.ic.get(), nullptr, tableRange, GMX_MAKETABLES_14ONLY);
180         fr_.efep       = haveFep_;
181         fr_.fudgeQQ          = 0.5;
182         fr_.use_simd_kernels = useSimd_;
183         fr_.bMolPBC          = haveMolPBC_;
184     }
185
186     t_forcerec* get() { return &fr_; }
187
188 private:
189     FreeEnergyPerturbationType haveFep_    = FreeEnergyPerturbationType::No;
190     bool                       useSimd_    = false;
191     bool                       haveMolPBC_ = false;
192
193     t_lambda   fepVals_;
194     t_forcerec fr_;
195 };
196
197 ForcerecHelper frHelper;
198
199
200 /*! \brief Input structure for listed forces tests
201  */
202 struct ListInput
203 {
204 public:
205     //! Function type
206     int fType = -1;
207     //! do fep
208     bool fep = false;
209     //! Tolerance for float evaluation
210     float floatToler = 1e-6;
211     //! Tolerance for double evaluation
212     double doubleToler = 1e-8;
213     //! Interaction parameters
214     t_iparams iparams = { { 0 } };
215
216     //! Constructor
217     ListInput() {}
218
219     /*! \brief Constructor with tolerance
220      *
221      * \param[in] ftol Single precision tolerance
222      * \param[in] dtol Double precision tolerance
223      */
224     ListInput(float ftol, double dtol)
225     {
226         floatToler  = ftol;
227         doubleToler = dtol;
228     }
229
230     /*! \brief Set parameters for lj14 interaction
231      *
232      * Fep is used if either c6A != c6B or c12A != c12B.
233      *
234      * \param[in] c6A  lj-c6 of state A
235      * \param[in] c12A lj-c12 of state A
236      * \param[in] c6B  lj-c6 of state B
237      * \param[in] c12B lj-c12 of state B
238      */
239     ListInput setLj14Interaction(real c6A, real c12A, real c6B, real c12B)
240     {
241         fType             = F_LJ14;
242         fep               = (c6A != c6B || c12A != c12B);
243         iparams.lj14.c6A  = c6A;
244         iparams.lj14.c12A = c12A;
245         iparams.lj14.c6B  = c6B;
246         iparams.lj14.c12B = c12B;
247
248         return *this;
249     }
250
251     /*! \brief Set parameters for ljc14 interaction
252      *
253      * \param[in] qi     charge i
254      * \param[in] qj     charge j
255      * \param[in] c6     lj-c6
256      * \param[in] c12    lj-c12
257      * \param[in] fudgeQ fudge factor
258      */
259     ListInput setLjc14Interaction(real qi, real qj, real c6, real c12, real fudgeQ)
260     {
261         fType             = F_LJC14_Q;
262         iparams.ljc14.qi  = qi;
263         iparams.ljc14.qj  = qj;
264         iparams.ljc14.c6  = c6;
265         iparams.ljc14.c12 = c12;
266         iparams.ljc14.fqq = fudgeQ;
267
268         return *this;
269     }
270
271     /*! \brief Set parameters for ljcnb interaction
272      *
273      * \param[in] qi  charge i
274      * \param[in] qj  charge j
275      * \param[in] c6  lj-c6
276      * \param[in] c12 lj-c12
277      */
278     ListInput setLjcnbInteraction(real qi, real qj, real c6, real c12)
279     {
280         fType             = F_LJC_PAIRS_NB;
281         iparams.ljcnb.qi  = qi;
282         iparams.ljcnb.qj  = qj;
283         iparams.ljcnb.c6  = c6;
284         iparams.ljcnb.c12 = c12;
285
286         return *this;
287     }
288 };
289
290 class ListedForcesPairsTest :
291     public ::testing::TestWithParam<std::tuple<ListInput, PaddedVector<RVec>, PbcType>>
292 {
293 protected:
294     matrix               box_;
295     t_pbc                pbc_;
296     PaddedVector<RVec>   x_;
297     PbcType              pbcType_;
298     ListInput            input_;
299     TestReferenceData    refData_;
300     TestReferenceChecker checker_;
301
302     ListedForcesPairsTest() : checker_(refData_.rootChecker())
303     {
304         input_   = std::get<0>(GetParam());
305         x_       = std::get<1>(GetParam());
306         pbcType_ = std::get<2>(GetParam());
307         clear_mat(box_);
308         box_[0][0] = box_[1][1] = box_[2][2] = 1.0;
309         set_pbc(&pbc_, pbcType_, box_);
310
311         FloatingPointTolerance tolerance = relativeToleranceAsPrecisionDependentFloatingPoint(
312                 1.0, input_.floatToler, input_.doubleToler);
313
314         checker_.setDefaultTolerance(tolerance);
315     }
316
317     void testOneIfunc(TestReferenceChecker* checker, const real lambda)
318     {
319         SCOPED_TRACE(std::string("Testing PBC type: ") + c_pbcTypeNames[pbcType_]);
320
321         // 'definition of pairs' is a concatenation of #npairs (here 2)
322         // 'nAtomsPerPair+1'-tuples (fType a_0 a_i ... a_nAtomsPerPair)
323         std::vector<t_iatom> iatoms = { 0, 1, 2, 0, 0, 2 };
324
325         std::vector<int>            ddgatindex = { 0, 1, 2 };
326         std::vector<real>           chargeA    = { 1.0, -0.5, -0.5 };
327         std::vector<real>           chargeB    = { 0.0, 0.0, 0.0 };
328         std::vector<unsigned short> egrp       = { 0, 0, 0 };
329         t_mdatoms                   mdatoms    = { 0 };
330
331         mdatoms.chargeA = chargeA.data();
332         mdatoms.chargeB = chargeB.data();
333         mdatoms.cENER   = egrp.data();
334         // nPerturbed is not decisive for fep to be used; it is overruled by
335         // other conditions in do_pairs_general; just here to not segfault
336         // upon query
337         mdatoms.nPerturbed = 0;
338
339         t_forcerec* fr = frHelper.get();
340         fr->efep = input_.fep ? FreeEnergyPerturbationType::Yes : FreeEnergyPerturbationType::No;
341         if (pbcType_ != PbcType::No)
342         {
343             fr->bMolPBC = true;
344         }
345
346         std::vector<BondedKernelFlavor> flavors = { BondedKernelFlavor::ForcesAndVirialAndEnergy };
347
348         if (!input_.fep || lambda == 0)
349         {
350             fr->use_simd_kernels = true;
351             flavors.push_back(BondedKernelFlavor::ForcesSimdWhenAvailable);
352         }
353
354         for (const auto flavor : flavors)
355         {
356             SCOPED_TRACE("Testing bonded kernel flavor: " + c_bondedKernelFlavorStrings[flavor]);
357
358             StepWorkload stepWork;
359             stepWork.computeVirial = computeVirial(flavor);
360             stepWork.computeEnergy = computeEnergy(flavor);
361
362             bool havePerturbedInteractions = input_.fep;
363             if (flavor == BondedKernelFlavor::ForcesSimdWhenAvailable)
364             {
365                 havePerturbedInteractions = false;
366             }
367
368             int numEnergyTerms      = static_cast<int>(NonBondedEnergyTerms::Count);
369             int numFepCouplingTerms = static_cast<int>(FreeEnergyPerturbationCouplingType::Count);
370             OutputQuantities  output(numEnergyTerms);
371             std::vector<real> lambdas(numFepCouplingTerms, lambda);
372
373             do_pairs(input_.fType,
374                      iatoms.size(),
375                      iatoms.data(),
376                      &input_.iparams,
377                      as_rvec_array(x_.data()),
378                      output.f,
379                      output.fShift,
380                      &pbc_,
381                      lambdas.data(),
382                      output.dvdLambda.data(),
383                      gmx::arrayRefFromArray(mdatoms.chargeA, mdatoms.nr),
384                      gmx::arrayRefFromArray(mdatoms.chargeB, mdatoms.nr),
385                      gmx::arrayRefFromArray(mdatoms.bPerturbed, mdatoms.nr),
386                      gmx::arrayRefFromArray(mdatoms.cENER, mdatoms.nr),
387                      mdatoms.nPerturbed,
388                      fr,
389                      havePerturbedInteractions,
390                      stepWork,
391                      &output.energy,
392                      ddgatindex.data());
393
394             checkOutput(checker, output, flavor, input_.fType);
395             auto shiftForcesChecker = checker->checkCompound("Shift-Forces", "Shift-forces");
396
397             if (computeVirial(flavor))
398             {
399                 shiftForcesChecker.checkVector(output.fShift[gmx::c_centralShiftIndex], "Central");
400             }
401             else
402             {
403                 // Permit omitting to compare shift forces with
404                 // reference data when that is useless.
405                 shiftForcesChecker.disableUnusedEntriesCheck();
406             }
407         }
408     }
409
410     void testIfunc()
411     {
412         TestReferenceChecker thisChecker =
413                 checker_.checkCompound("FunctionType", interaction_function[input_.fType].name)
414                         .checkCompound("FEP", (input_.fep ? "Yes" : "No"));
415
416         if (input_.fep)
417         {
418             const int numLambdas = 3;
419             for (int i = 0; i < numLambdas; ++i)
420             {
421                 const real lambda        = i / (numLambdas - 1.0);
422                 auto       lambdaChecker = thisChecker.checkCompound("Lambda", toString(lambda));
423                 testOneIfunc(&lambdaChecker, lambda);
424             }
425         }
426         else
427         {
428             testOneIfunc(&thisChecker, 0.0);
429         }
430     }
431 };
432
433 TEST_P(ListedForcesPairsTest, Ifunc)
434 {
435     testIfunc();
436 }
437
438 //! Function types for testing 1-4 interaction. Add new terms at the end.
439 std::vector<ListInput> c_14Interaction = {
440     { ListInput(1e-5, 1e-7).setLj14Interaction(0.001458, 1.0062882e-6, 0.0, 0.0) },
441     { ListInput(1e-5, 1e-7).setLj14Interaction(0.001458, 1.0062882e-6, 0.001458, 1.0062882e-6) },
442     { ListInput(1e-5, 1e-7).setLjc14Interaction(1.0, -1.0, 0.001458, 1.0062882e-6, 0.5) },
443     { ListInput(1e-5, 1e-7).setLjcnbInteraction(1.0, -1.0, 0.001458, 1.0062882e-6) }
444 };
445
446 //! PBC values for testing
447 std::vector<PbcType> c_pbcForTests = { PbcType::No, PbcType::XY, PbcType::Xyz };
448
449 /*! \brief Coordinates for testing 1-4 interaction
450  *
451  * Define coordinates for 3 atoms here, which will be used in 2 interactions.
452  */
453 std::vector<PaddedVector<RVec>> c_coordinatesFor14Interaction = {
454     { { 0.0, 0.0, 0.0 }, { 1.0, 1.0, 1.0 }, { 1.1, 1.2, 1.3 } }
455 };
456
457 INSTANTIATE_TEST_CASE_P(14Interaction,
458                         ListedForcesPairsTest,
459                         ::testing::Combine(::testing::ValuesIn(c_14Interaction),
460                                            ::testing::ValuesIn(c_coordinatesFor14Interaction),
461                                            ::testing::ValuesIn(c_pbcForTests)));
462
463 } // namespace
464
465 } // namespace test
466
467 } // namespace gmx