SIMD support for nonbonded free-energy kernels
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / tests / nb_free_energy.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
36 /*! \internal \file
37  * \brief Implements test of nonbonded fep kernel
38  *
39  * Implements the test logic from the bonded interactions also for the
40  * nonbonded fep kernel. This requires setting up some more input
41  * structures that in the bonded case.
42  *
43  * The test setup consists of an atom pair that is evaluated in an fep setting
44  * (vanishing charge and lennard-jones parameters of atom #2) with and without
45  * softcore Potentials.
46  *
47  * \author Sebastian Kehl <sebastian.kehl@mpcdf.mpg.de>
48  * \ingroup module_gmxlib_nonbonded
49  */
50 #include "gmxpre.h"
51
52 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
53 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
54
55 #include <cmath>
56
57 #include <gtest/gtest.h>
58
59 #include "gromacs/math/paddedvector.h"
60 #include "gromacs/math/units.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/math/vectypes.h"
63 #include "gromacs/math/arrayrefwithpadding.h"
64 #include "gromacs/mdtypes/mdatom.h"
65 #include "gromacs/mdtypes/enerdata.h"
66 #include "gromacs/mdtypes/forcerec.h"
67 #include "gromacs/mdtypes/inputrec.h"
68 #include "gromacs/mdtypes/interaction_const.h"
69 #include "gromacs/mdtypes/nblist.h"
70 #include "gromacs/mdtypes/forceoutput.h"
71 #include "gromacs/mdlib/forcerec.h"
72 #include "gromacs/tables/forcetable.h"
73 #include "gromacs/pbcutil/ishift.h"
74 #include "gromacs/pbcutil/pbc.h"
75 #include "gromacs/topology/idef.h"
76 #include "gromacs/topology/forcefieldparameters.h"
77 #include "gromacs/utility/enumerationhelpers.h"
78 #include "gromacs/utility/strconvert.h"
79 #include "gromacs/utility/stringstream.h"
80 #include "gromacs/utility/textwriter.h"
81 #include "gromacs/utility/arrayref.h"
82 #include "gromacs/utility/smalloc.h"
83 #include "gromacs/ewald/ewald_utils.h"
84 #include "gromacs/gmxlib/nrnb.h"
85
86 #include "testutils/refdata.h"
87 #include "testutils/testasserts.h"
88
89 namespace gmx
90 {
91 namespace test
92 {
93 namespace
94 {
95
96 //! Number of atoms used in these tests.
97 constexpr int c_numAtoms     = 4;
98 constexpr int c_numAtomTypes = 3;
99
100 /*! \brief Output from nonbonded fep kernel
101  *
102  */
103 struct OutputQuantities
104 {
105     OutputQuantities() :
106         energy(static_cast<int>(NonBondedEnergyTerms::Count)),
107         dvdLambda(static_cast<int>(FreeEnergyPerturbationCouplingType::Count), 0.0),
108         fShift(1, { 0.0, 0.0, 0.0 }),
109         f(c_numAtoms, { 0.0, 0.0, 0.0 })
110     {
111     }
112
113     //! Energies of this interaction (size EgNR)
114     gmx_grppairener_t energy;
115     //! Derivative with respect to lambda (size efptNR)
116     std::vector<real> dvdLambda;
117     //! Shift force vectors (size N_IVEC but in this test only 1)
118     std::vector<RVec> fShift;
119     //! Forces (size c_numAtoms)
120     PaddedVector<RVec> f;
121 };
122
123 /*! \brief Utility to check the output from nonbonded test
124  *
125  * \param[in] checker Reference checker
126  * \param[in] output  The output from the test to check
127  */
128 void checkOutput(TestReferenceChecker* checker, const OutputQuantities& output)
129 {
130     checker->checkReal(output.energy.energyGroupPairTerms[NonBondedEnergyTerms::LJSR][0], "EVdw ");
131     checker->checkReal(output.energy.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR][0],
132                        "ECoul ");
133     checker->checkReal(output.dvdLambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
134                        "dVdlCoul ");
135     checker->checkReal(output.dvdLambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
136                        "dVdlVdw ");
137
138     checker->checkSequence(std::begin(output.f), std::end(output.f), "Forces");
139
140     auto shiftForcesChecker = checker->checkCompound("Shift-Forces", "Shift-forces");
141     shiftForcesChecker.checkVector(output.fShift[0], "Central");
142 }
143
144 class InteractionConstHelper
145 {
146 public:
147     InteractionConstHelper() {}
148
149     //! init data to construct interaction_const
150     void initInteractionConst(CoulombInteractionType coulType, VanDerWaalsType vdwType, InteractionModifiers vdwMod)
151     {
152         coulType_ = coulType;
153         vdwType_  = vdwType;
154         vdwMod_   = vdwMod;
155
156         // initialize correction tables
157         interaction_const_t tmp;
158         tmp.ewaldcoeff_q       = calc_ewaldcoeff_q(1.0, 1.0e-5);
159         coulEwaldCoeff_        = tmp.ewaldcoeff_q;
160         tmp.ewaldcoeff_lj      = calc_ewaldcoeff_lj(1.0, 1.0e-5);
161         vdwEwaldCoeff_         = tmp.ewaldcoeff_lj;
162         tmp.eeltype            = coulType;
163         tmp.vdwtype            = vdwType;
164         tmp.coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
165         tmp.vdwEwaldTables     = std::make_unique<EwaldCorrectionTables>();
166
167         init_interaction_const_tables(nullptr, &tmp, 1.0, 0.0);
168         coulombTables_ = *tmp.coulombEwaldTables;
169         vdwTables_     = *tmp.vdwEwaldTables;
170     }
171
172     /*! \brief Setup interaction_const_t
173      *
174      * \param[in]  fepVals t_lambda struct of fep values
175      * \param[out] ic      interaction_const_t pointer with data
176      */
177     void getInteractionConst(const t_lambda& fepVals, interaction_const_t* ic)
178     {
179         ic->softCoreParameters = std::unique_ptr<interaction_const_t::SoftCoreParameters>(
180                 new interaction_const_t::SoftCoreParameters(fepVals));
181
182         ic->coulombEwaldTables  = std::unique_ptr<EwaldCorrectionTables>(new EwaldCorrectionTables);
183         *ic->coulombEwaldTables = coulombTables_;
184
185         ic->vdwEwaldTables  = std::unique_ptr<EwaldCorrectionTables>(new EwaldCorrectionTables);
186         *ic->vdwEwaldTables = vdwTables_;
187
188         // set coulomb and vdw types
189         ic->eeltype      = coulType_;
190         ic->vdwtype      = vdwType_;
191         ic->vdw_modifier = vdwMod_;
192
193         // some non default parameters used in this testcase
194         ic->epsfac                   = gmx::c_one4PiEps0 * 0.25;
195         ic->reactionFieldCoefficient = 0.0; // former k_rf
196         ic->reactionFieldShift       = 1.0; // former c_rf
197         ic->ewaldcoeff_q             = coulEwaldCoeff_;
198         ic->ewaldcoeff_lj            = vdwEwaldCoeff_;
199         ic->sh_ewald                 = 1.0e-5;
200         ic->sh_lj_ewald              = -1.0;
201         ic->dispersion_shift.cpot    = -1.0;
202         ic->repulsion_shift.cpot     = -1.0;
203     }
204
205 private:
206     //! correction tables
207     EwaldCorrectionTables coulombTables_;
208     EwaldCorrectionTables vdwTables_;
209
210     //! coulomb and vdw type specifiers
211     CoulombInteractionType coulType_;
212     real                   coulEwaldCoeff_;
213     VanDerWaalsType        vdwType_;
214     InteractionModifiers   vdwMod_;
215     real                   vdwEwaldCoeff_;
216 };
217
218
219 /* \brief Utility class to setup forcerec
220  *
221  * This helper takes care of handling the neccessary data that are kept
222  * at various places and in various forms in the forcerec hierarchy such
223  * that this class can safely be used.
224  *
225  * Data is only initialized as necessary for the nonbonded kernel to work!
226  */
227 class ForcerecHelper
228 {
229 public:
230     ForcerecHelper()
231     {
232         fepVals_.sc_alpha     = 0.3;
233         fepVals_.sc_power     = 1;
234         fepVals_.sc_r_power   = 6.0;
235         fepVals_.sc_sigma     = 0.3;
236         fepVals_.sc_sigma_min = 0.3;
237         fepVals_.bScCoul      = true;
238     }
239
240     //! initialize data structure to construct forcerec
241     void initForcerec(const gmx_ffparams_t&  idef,
242                       CoulombInteractionType coulType,
243                       VanDerWaalsType        vdwType,
244                       InteractionModifiers   vdwMod)
245     {
246         icHelper_.initInteractionConst(coulType, vdwType, vdwMod);
247         nbfp_        = makeNonBondedParameterLists(idef.atnr, idef.iparams, false);
248         ljPmeC6Grid_ = makeLJPmeC6GridCorrectionParameters(idef.atnr, idef.iparams, LongRangeVdW::Geom);
249     }
250
251     void setSoftcoreAlpha(const real scAlpha) { fepVals_.sc_alpha = scAlpha; }
252     void setSoftcoreCoulomb(const bool scCoulomb) { fepVals_.bScCoul = scCoulomb; }
253
254     //! get forcerec data as wanted by the nonbonded kernel
255     void getForcerec(t_forcerec* fr)
256     {
257         fr->ic = std::make_unique<interaction_const_t>();
258
259         // set data in ic
260         icHelper_.getInteractionConst(fepVals_, fr->ic.get());
261
262         // set data in fr
263         fr->ljpme_c6grid = ljPmeC6Grid_;
264         fr->nbfp         = nbfp_;
265         fr->shift_vec    = { { 0.0, 0.0, 0.0 } };
266         fr->rlist        = rlist_;
267         fr->ntype        = c_numAtomTypes;
268
269         // simd
270         fr->use_simd_kernels = GMX_USE_SIMD_KERNELS;
271     }
272
273 private:
274     InteractionConstHelper icHelper_;
275     std::vector<real>      ljPmeC6Grid_;
276     std::vector<real>      nbfp_;
277     t_lambda               fepVals_;
278     real                   rlist_ = 1.0;
279 };
280
281 /*! \brief Utility structure to hold atoms data
282  *
283  * A system having 4 interactions with the following perturbation pattern:
284  *  - no perturbation
285  *  - vdw- and coulomb-perturbation
286  *  - coulomb-perturbation only
287  *  - vdw-perturbation only
288  *
289  * This is realized by defining 3 different atom types that control
290  * the vdw-perturbation. The coulomb-perturbation is controlled by directly
291  * setting the charge of the atoms at the lambda states A/B.
292  */
293 struct AtomData
294 {
295     AtomData()
296     {
297         idef.atnr = c_numAtomTypes;
298         idef.iparams.resize(2 * c_numAtomTypes * c_numAtomTypes);
299
300         // set interaction parameters for different combinations of types
301         idef.iparams[0].lj = { 0.001458, 1.0062882e-6 }; // 0-0
302         idef.iparams[1].lj = { 0.0, 0.0 };               // 0-1
303         idef.iparams[2].lj = { 0.001458, 1.0062882e-6 }; // 0-2
304         idef.iparams[3].lj = { 0.0, 0.0 };               // 1-0
305         idef.iparams[4].lj = { 0.0, 0.0 };               // 1-1
306         idef.iparams[5].lj = { 0.0, 0.0 };               // 1-2
307         idef.iparams[6].lj = { 0.001458, 1.0062882e-6 }; // 2-0
308         idef.iparams[7].lj = { 0.0, 0.0 };               // 2-1
309         idef.iparams[8].lj = { 0.001458, 1.0062882e-6 }; // 2-2
310
311         GMX_ASSERT(chargeA.size() == c_numAtoms, "This test wants 4 atoms");
312         GMX_ASSERT(chargeB.size() == c_numAtoms, "This test wants 4 atoms");
313         GMX_ASSERT(typeA.size() == c_numAtoms, "This test wants 4 atoms");
314         GMX_ASSERT(typeB.size() == c_numAtoms, "This test wants 4 atoms");
315     }
316
317     // forcefield parameters
318     gmx_ffparams_t idef;
319
320     // atom data
321     std::vector<real> chargeA = { 1.0, -1.0, -1.0, 1.0 };
322     std::vector<real> chargeB = { 1.0, 0.0, 0.0, 1.0 };
323     std::vector<int>  typeA   = { 0, 0, 0, 0 };
324     std::vector<int>  typeB   = { 0, 1, 2, 1 };
325     // perturbation pattern: {no-pert, vdw- and coul-pert, coul-pert, vdw-pert}
326
327     // neighbourhood information
328     std::vector<int> iAtoms  = { 0 };
329     std::vector<int> jAtoms  = { 0, 1, 2, 3 };
330     std::vector<int> jIndex  = { 0, 4 };
331     std::vector<int> shift   = { 0 };
332     std::vector<int> gid     = { 0 };
333     std::vector<int> exclFep = { 0, 1, 1, 1 };
334
335     // construct t_nblist
336     t_nblist getNbList()
337     {
338         t_nblist nbl;
339         nbl.nri      = 1;
340         nbl.nrj      = 4;
341         nbl.iinr     = iAtoms;
342         nbl.jindex   = jIndex;
343         nbl.jjnr     = jAtoms;
344         nbl.shift    = shift;
345         nbl.gid      = gid;
346         nbl.excl_fep = exclFep;
347         return nbl;
348     }
349 };
350
351 /*! \brief Input structure for nonbonded fep kernel
352  */
353 struct ListInput
354 {
355 public:
356     //! Function type
357     int fType = F_LJ;
358     //! Tolerance for float evaluation
359     float floatToler = 1e-6;
360     //! Tolerance for double evaluation
361     double doubleToler = 1e-8;
362     //! atom parameters
363     AtomData atoms;
364     //! forcerec helper
365     ForcerecHelper frHelper;
366
367     //! Constructor
368     ListInput() {}
369
370     /*! \brief Constructor with tolerance
371      *
372      * \param[in] ftol Single precision tolerance
373      * \param[in] dtol Double precision tolerance
374      */
375     ListInput(float ftol, double dtol)
376     {
377         floatToler  = ftol;
378         doubleToler = dtol;
379     }
380
381     /*! \brief Set parameters for nonbonded interaction
382      *
383      * \param[in] coulType coulomb type
384      * \param[in] vdwType  vdw type
385      * \param[in] vdwMod   vdw potential modifier
386      */
387     ListInput setInteraction(CoulombInteractionType coulType, VanDerWaalsType vdwType, InteractionModifiers vdwMod)
388     {
389         frHelper.initForcerec(atoms.idef, coulType, vdwType, vdwMod);
390         return *this;
391     }
392 };
393
394 class NonbondedFepTest :
395     public ::testing::TestWithParam<std::tuple<ListInput, PaddedVector<RVec>, real, real, bool>>
396 {
397 protected:
398     PaddedVector<RVec>   x_;
399     ListInput            input_;
400     real                 lambda_;
401     real                 softcoreAlpha_;
402     bool                 softcoreCoulomb_;
403     TestReferenceData    refData_;
404     TestReferenceChecker checker_;
405
406     NonbondedFepTest() : checker_(refData_.rootChecker())
407     {
408         input_           = std::get<0>(GetParam());
409         x_               = std::get<1>(GetParam());
410         lambda_          = std::get<2>(GetParam());
411         softcoreAlpha_   = std::get<3>(GetParam());
412         softcoreCoulomb_ = std::get<4>(GetParam());
413
414         // Note that the reference data for Ewald type interactions has been generated
415         // with accurate analytical approximations for the long-range corrections.
416         // When the free-energy kernel switches from tabulated to analytical corrections,
417         // the double precision tolerance can be tightend to 1e-11.
418         test::FloatingPointTolerance tolerance(
419                 input_.floatToler, input_.doubleToler, 1.0e-6, 1.0e-11, 10000, 100, false);
420         checker_.setDefaultTolerance(tolerance);
421     }
422
423     void testKernel()
424     {
425         input_.frHelper.setSoftcoreAlpha(softcoreAlpha_);
426         input_.frHelper.setSoftcoreCoulomb(softcoreCoulomb_);
427
428         // get forcerec and interaction_const
429         t_forcerec fr;
430         input_.frHelper.getForcerec(&fr);
431
432         // t_nblist
433         t_nblist nbl = input_.atoms.getNbList();
434
435         // output buffers
436         OutputQuantities output;
437
438         // lambda vector
439         int numFepCouplingTerms = static_cast<int>(FreeEnergyPerturbationCouplingType::Count);
440         std::vector<real> lambdas(numFepCouplingTerms, lambda_);
441
442         // fep kernel data
443         int doNBFlags = 0;
444         doNBFlags |= GMX_NONBONDED_DO_FORCE;
445         doNBFlags |= GMX_NONBONDED_DO_SHIFTFORCE;
446         doNBFlags |= GMX_NONBONDED_DO_POTENTIAL;
447
448         // force buffers
449         bool                      unusedBool = true; // this bool has no effect in the kernel
450         gmx::ForceWithShiftForces forces(output.f.arrayRefWithPadding(), unusedBool, output.fShift);
451
452         // dummy counter
453         t_nrnb nrnb;
454
455         // run fep kernel
456         gmx_nb_free_energy_kernel(nbl,
457                                   x_.arrayRefWithPadding(),
458                                   fr.use_simd_kernels,
459                                   fr.ntype,
460                                   fr.rlist,
461                                   *fr.ic,
462                                   fr.shift_vec,
463                                   fr.nbfp,
464                                   fr.ljpme_c6grid,
465                                   input_.atoms.chargeA,
466                                   input_.atoms.chargeB,
467                                   input_.atoms.typeA,
468                                   input_.atoms.typeB,
469                                   doNBFlags,
470                                   lambdas,
471                                   &nrnb,
472                                   output.f.arrayRefWithPadding().paddedArrayRef().data(),
473                                   as_rvec_array(output.fShift.data()),
474                                   output.energy.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR],
475                                   output.energy.energyGroupPairTerms[NonBondedEnergyTerms::LJSR],
476                                   output.dvdLambda);
477
478         checkOutput(&checker_, output);
479     }
480 };
481
482 TEST_P(NonbondedFepTest, testKernel)
483 {
484     testKernel();
485 }
486
487 //! configurations to test
488 std::vector<ListInput> c_interaction = {
489     { ListInput(1e-6, 1e-8).setInteraction(CoulombInteractionType::Cut, VanDerWaalsType::Cut, InteractionModifiers::None) },
490     { ListInput(1e-6, 1e-8).setInteraction(CoulombInteractionType::Cut, VanDerWaalsType::Cut, InteractionModifiers::PotSwitch) },
491     { ListInput(1e-6, 1e-8).setInteraction(CoulombInteractionType::Pme, VanDerWaalsType::Pme, InteractionModifiers::None) }
492 };
493
494 //! test parameters
495 std::vector<real> c_fepLambdas      = { 0.0, 0.5, 1.0 };
496 std::vector<real> c_softcoreAlphas  = { 0.0, 0.3 };
497 std::vector<bool> c_softcoreCoulomb = { true, false };
498
499 //! Coordinates for testing
500 std::vector<PaddedVector<RVec>> c_coordinates = {
501     { { 1.0, 1.0, 1.0 }, { 1.1, 1.15, 1.2 }, { 0.9, 0.85, 0.8 }, { 1.1, 1.15, 0.8 } }
502 };
503
504 INSTANTIATE_TEST_SUITE_P(NBInteraction,
505                          NonbondedFepTest,
506                          ::testing::Combine(::testing::ValuesIn(c_interaction),
507                                             ::testing::ValuesIn(c_coordinates),
508                                             ::testing::ValuesIn(c_fepLambdas),
509                                             ::testing::ValuesIn(c_softcoreAlphas),
510                                             ::testing::ValuesIn(c_softcoreCoulomb)));
511
512 } // namespace
513
514 } // namespace test
515
516 } // namespace gmx