Rework -Weverything
[alexxy/gromacs.git] / api / nblib / tests / topology.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2020,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
37  * This implements topology setup tests
38  *
39  * \author Victor Holanda <victor.holanda@cscs.ch>
40  * \author Joe Jordan <ejjordan@kth.se>
41  * \author Prashanth Kanduri <kanduri@cscs.ch>
42  * \author Sebastian Keller <keller@cscs.ch>
43  */
44 #include <gmock/gmock.h>
45 #include <gtest/gtest.h>
46
47 #include "gromacs/topology/exclusionblocks.h"
48 #include "gromacs/utility/listoflists.h"
49 #include "nblib/exception.h"
50 #include "nblib/particletype.h"
51 #include "nblib/sequencing.hpp"
52 #include "nblib/tests/testsystems.h"
53 #include "nblib/topology.h"
54 #include "nblib/topologyhelpers.h"
55 #include "nblib/particlesequencer.h"
56
57 namespace nblib
58 {
59 namespace test
60 {
61 namespace
62 {
63
64 using ::testing::Eq;
65 using ::testing::Pointwise;
66
67 //! Compares all element between two lists of lists
68 //! Todo: unify this with the identical function in nbkernelsystem test make this a method
69 //!       of ListOfLists<>
70 template<typename T>
71 void compareLists(const gmx::ListOfLists<T>& list, const std::vector<std::vector<T>>& v)
72 {
73     ASSERT_EQ(list.size(), v.size());
74     for (std::size_t i = 0; i < list.size(); i++)
75     {
76         ASSERT_EQ(list[i].size(), v[i].size());
77         EXPECT_THAT(list[i], Pointwise(Eq(), v[i]));
78     }
79 }
80
81 TEST(NBlibTest, TopologyHasNumParticles)
82 {
83     WaterTopologyBuilder waters;
84     Topology             watersTopology = waters.buildTopology(2);
85     const int            test           = watersTopology.numParticles();
86     const int            ref            = 6;
87     EXPECT_EQ(ref, test);
88 }
89
90 TEST(NBlibTest, TopologyHasCharges)
91 {
92     WaterTopologyBuilder     waters;
93     Topology                 watersTopology = waters.buildTopology(2);
94     const std::vector<real>& test           = watersTopology.getCharges();
95     const std::vector<real>& ref = { Charges.at("Ow"), Charges.at("Hw"), Charges.at("Hw"),
96                                      Charges.at("Ow"), Charges.at("Hw"), Charges.at("Hw") };
97     EXPECT_EQ(ref, test);
98 }
99
100 TEST(NBlibTest, TopologyHasMasses)
101 {
102     WaterTopologyBuilder waters;
103     Topology             watersTopology = waters.buildTopology(2);
104
105     const Mass              refOwMass = waters.water().at("Ow").mass();
106     const Mass              refHwMass = waters.water().at("H").mass();
107     const std::vector<Mass> ref = { refOwMass, refHwMass, refHwMass, refOwMass, refHwMass, refHwMass };
108     const std::vector<Mass> test = expandQuantity(watersTopology, &ParticleType::mass);
109     EXPECT_EQ(ref, test);
110 }
111
112 TEST(NBlibTest, TopologyHasParticleTypes)
113 {
114     WaterTopologyBuilder             waters;
115     Topology                         watersTopology = waters.buildTopology(2);
116     const std::vector<ParticleType>& test           = watersTopology.getParticleTypes();
117     const ParticleType               refOw          = waters.water().at("Ow");
118     const ParticleType               refHw          = waters.water().at("H");
119     const std::vector<ParticleType>& ref            = { refOw, refHw };
120     const std::vector<ParticleType>& ref2           = { refHw, refOw };
121     EXPECT_TRUE(ref == test || ref2 == test);
122 }
123
124 TEST(NBlibTest, TopologyHasParticleTypeIds)
125 {
126     WaterTopologyBuilder waters;
127     Topology             watersTopology = waters.buildTopology(2);
128
129     const std::vector<int>&          testIds   = watersTopology.getParticleTypeIdOfAllParticles();
130     const std::vector<ParticleType>& testTypes = watersTopology.getParticleTypes();
131
132     std::vector<ParticleType> testTypesExpanded;
133     testTypesExpanded.reserve(testTypes.size());
134     for (int i : testIds)
135     {
136         testTypesExpanded.push_back(testTypes[i]);
137     }
138
139     const ParticleType              refOw = waters.water().at("Ow");
140     const ParticleType              refHw = waters.water().at("H");
141     const std::vector<ParticleType> ref   = { refOw, refHw, refHw, refOw, refHw, refHw };
142
143     EXPECT_TRUE(ref == testTypesExpanded);
144 }
145
146 TEST(NBlibTest, TopologyThrowsIdenticalParticleType)
147 {
148     //! User error: Two different ParticleTypes with the same name
149     ParticleType U235(ParticleTypeName("Uranium"), Mass(235));
150     ParticleType U238(ParticleTypeName("Uranium"), Mass(238));
151
152     Molecule ud235(MoleculeName("UraniumDimer235"));
153     ud235.addParticle(ParticleName("U1"), U235);
154     ud235.addParticle(ParticleName("U2"), U235);
155
156     Molecule ud238(MoleculeName("UraniumDimer238"));
157     ud238.addParticle(ParticleName("U1"), U238);
158     ud238.addParticle(ParticleName("U2"), U238);
159
160     TopologyBuilder topologyBuilder;
161     topologyBuilder.addMolecule(ud235, 1);
162     EXPECT_THROW(topologyBuilder.addMolecule(ud238, 1), InputException);
163 }
164
165 TEST(NBlibTest, TopologyHasExclusions)
166 {
167     WaterTopologyBuilder        waters;
168     Topology                    watersTopology = waters.buildTopology(2);
169     ExclusionLists<int>         exclusionLists = watersTopology.exclusionLists();
170     const gmx::ListOfLists<int> testExclusions(std::move(exclusionLists.ListRanges),
171                                                std::move(exclusionLists.ListElements));
172
173     const std::vector<std::vector<int>>& refExclusions = { { 0, 1, 2 }, { 0, 1, 2 }, { 0, 1, 2 },
174                                                            { 3, 4, 5 }, { 3, 4, 5 }, { 3, 4, 5 } };
175
176     compareLists(testExclusions, refExclusions);
177 }
178
179 TEST(NBlibTest, TopologyHasSequencing)
180 {
181     WaterTopologyBuilder waters;
182     Topology             watersTopology = waters.buildTopology(2);
183
184     EXPECT_EQ(0, watersTopology.sequenceID(MoleculeName("SOL"), 0, ResidueName("SOL"), ParticleName("Oxygen")));
185     EXPECT_EQ(1, watersTopology.sequenceID(MoleculeName("SOL"), 0, ResidueName("SOL"), ParticleName("H1")));
186     EXPECT_EQ(2, watersTopology.sequenceID(MoleculeName("SOL"), 0, ResidueName("SOL"), ParticleName("H2")));
187     EXPECT_EQ(3, watersTopology.sequenceID(MoleculeName("SOL"), 1, ResidueName("SOL"), ParticleName("Oxygen")));
188     EXPECT_EQ(4, watersTopology.sequenceID(MoleculeName("SOL"), 1, ResidueName("SOL"), ParticleName("H1")));
189     EXPECT_EQ(5, watersTopology.sequenceID(MoleculeName("SOL"), 1, ResidueName("SOL"), ParticleName("H2")));
190 }
191
192 TEST(NBlibTest, TopologyCanAggregateBonds)
193 {
194     Molecule water    = WaterMoleculeBuilder{}.waterMolecule();
195     Molecule methanol = MethanolMoleculeBuilder{}.methanolMolecule();
196
197     std::vector<std::tuple<Molecule, int>> molecules{ std::make_tuple(water, 2),
198                                                       std::make_tuple(methanol, 1) };
199     std::vector<HarmonicBondType>          bonds;
200     std::vector<size_t>                    bondsExpansion;
201     std::tie(bondsExpansion, bonds) = detail::collectInteractions<HarmonicBondType>(molecules);
202
203     std::vector<HarmonicBondType> bondsTest;
204     // use the expansionArray (bondsExpansion) to expand to the full list if bonds
205     std::transform(begin(bondsExpansion),
206                    end(bondsExpansion),
207                    std::back_inserter(bondsTest),
208                    [&bonds](size_t i) { return bonds[i]; });
209
210     std::vector<HarmonicBondType> waterBonds =
211             pickType<HarmonicBondType>(water.interactionData()).interactionTypes_;
212     std::vector<HarmonicBondType> methanolBonds =
213             pickType<HarmonicBondType>(methanol.interactionData()).interactionTypes_;
214
215     std::vector<HarmonicBondType> bondsReference;
216     std::copy(begin(waterBonds), end(waterBonds), std::back_inserter(bondsReference));
217     std::copy(begin(waterBonds), end(waterBonds), std::back_inserter(bondsReference));
218     std::copy(begin(methanolBonds), end(methanolBonds), std::back_inserter(bondsReference));
219
220     EXPECT_EQ(bondsTest, bondsReference);
221 }
222
223 TEST(NBlibTest, TopologyCanSequencePairIDs)
224 {
225     Molecule water    = WaterMoleculeBuilder{}.waterMolecule();
226     Molecule methanol = MethanolMoleculeBuilder{}.methanolMolecule();
227
228     std::vector<std::tuple<Molecule, int>> molecules{ std::make_tuple(water, 2),
229                                                       std::make_tuple(methanol, 1) };
230     ParticleSequencer                      particleSequencer;
231     particleSequencer.build(molecules);
232     auto pairs = detail::sequenceIDs<HarmonicBondType>(molecules, particleSequencer);
233
234     int Ow1 = particleSequencer(MoleculeName("SOL"), 0, ResidueName("SOL"), ParticleName("Oxygen"));
235     int H11 = particleSequencer(MoleculeName("SOL"), 0, ResidueName("SOL"), ParticleName("H1"));
236     int H12 = particleSequencer(MoleculeName("SOL"), 0, ResidueName("SOL"), ParticleName("H2"));
237     int Ow2 = particleSequencer(MoleculeName("SOL"), 1, ResidueName("SOL"), ParticleName("Oxygen"));
238     int H21 = particleSequencer(MoleculeName("SOL"), 1, ResidueName("SOL"), ParticleName("H1"));
239     int H22 = particleSequencer(MoleculeName("SOL"), 1, ResidueName("SOL"), ParticleName("H2"));
240
241     int Me  = particleSequencer(MoleculeName("MeOH"), 0, ResidueName("MeOH"), ParticleName("Me1"));
242     int MeO = particleSequencer(MoleculeName("MeOH"), 0, ResidueName("MeOH"), ParticleName("O2"));
243     int MeH = particleSequencer(MoleculeName("MeOH"), 0, ResidueName("MeOH"), ParticleName("H3"));
244
245     /// \cond DO_NOT_DOCUMENT
246 #define SORT(i, j) (i < j) ? i : j, (i < j) ? j : i
247
248     std::vector<CoordinateIndex<HarmonicBondType>> pairs_reference{
249         { SORT(Ow1, H11) }, { SORT(Ow1, H12) }, { SORT(Ow2, H21) },
250         { SORT(Ow2, H22) }, { SORT(MeO, MeH) }, { SORT(MeO, Me) }
251     };
252 #undef SORT
253     /// \endcond
254
255     EXPECT_EQ(pairs, pairs_reference);
256 }
257
258 TEST(NBlibTest, TopologySequenceIdThrows)
259 {
260     Molecule water    = WaterMoleculeBuilder{}.waterMolecule();
261     Molecule methanol = MethanolMoleculeBuilder{}.methanolMolecule();
262
263     std::vector<std::tuple<Molecule, int>> molecules{ std::make_tuple(water, 2),
264                                                       std::make_tuple(methanol, 1) };
265     ParticleSequencer                      particleSequencer;
266     particleSequencer.build(molecules);
267     auto pairs = detail::sequenceIDs<HarmonicBondType>(molecules, particleSequencer);
268
269     // Input error: no particle called O-Atom in molecule "water"
270     EXPECT_THROW(particleSequencer(MoleculeName("SOL"), 0, ResidueName("SOL"), ParticleName("O-Atom")),
271                  InputException);
272 }
273
274 TEST(NBlibTest, TopologyCanEliminateDuplicateBonds)
275 {
276     HarmonicBondType b1(1.0, 2.0);
277     HarmonicBondType b2(1.1, 2.1);
278     HarmonicBondType b3(1.2, 2.2);
279
280     // can be compressed to {b1,b2,b3} + {1,1,2,0,1,0,2,2}
281     std::vector<HarmonicBondType> bonds{ b2, b2, b3, b1, b2, b1, b3, b3 };
282
283     // expected output
284     std::vector<HarmonicBondType> uniqueBondsReference{ b1, b2, b3 };
285     std::vector<size_t>           indicesReference{ 1, 1, 2, 0, 1, 0, 2, 2 };
286
287     std::tuple<std::vector<size_t>, std::vector<HarmonicBondType>> bondData =
288             detail::eliminateDuplicateInteractions(bonds);
289
290     auto indices     = std::get<0>(bondData);
291     auto uniqueBonds = std::get<1>(bondData);
292
293     EXPECT_EQ(uniqueBondsReference, uniqueBonds);
294     EXPECT_EQ(indicesReference, indices);
295 }
296
297 TEST(NBlibTest, TopologyListedInteractions)
298 {
299     // Todo: add angles to SpcMethanolTopologyBuilder and extend interactions_reference below to
300     // Todo: include angles
301
302     Topology spcTopology = SpcMethanolTopologyBuilder{}.buildTopology(1, 2);
303
304     auto  interactionData = spcTopology.getInteractionData();
305     auto& harmonicBonds   = pickType<HarmonicBondType>(interactionData);
306
307     auto& indices = harmonicBonds.indices;
308     auto& bonds   = harmonicBonds.parameters;
309
310     std::map<std::tuple<int, int>, HarmonicBondType> interactions_test;
311     for (auto& ituple : indices)
312     {
313         interactions_test[std::make_tuple(std::get<0>(ituple), std::get<1>(ituple))] =
314                 bonds[std::get<2>(ituple)];
315     }
316
317     // there should be 3 unique HarmonicBondType instances
318     EXPECT_EQ(bonds.size(), 3);
319     // and 6 interaction pairs (bonds)
320     EXPECT_EQ(indices.size(), 6);
321
322     HarmonicBondType ohBond(1., 1.);
323     HarmonicBondType ohBondMethanol(1.01, 1.02);
324     HarmonicBondType ometBond(1.1, 1.2);
325
326     std::map<std::tuple<int, int>, HarmonicBondType> interactions_reference;
327
328     int Ow = spcTopology.sequenceID(MoleculeName("SOL"), 0, ResidueName("SOL"), ParticleName("Oxygen"));
329     int H1 = spcTopology.sequenceID(MoleculeName("SOL"), 0, ResidueName("SOL"), ParticleName("H1"));
330     int H2 = spcTopology.sequenceID(MoleculeName("SOL"), 0, ResidueName("SOL"), ParticleName("H2"));
331
332     int Me1 = spcTopology.sequenceID(MoleculeName("MeOH"), 0, ResidueName("MeOH"), ParticleName("Me1"));
333     int MeO1 = spcTopology.sequenceID(MoleculeName("MeOH"), 0, ResidueName("MeOH"), ParticleName("O2"));
334     int MeH1 = spcTopology.sequenceID(MoleculeName("MeOH"), 0, ResidueName("MeOH"), ParticleName("H3"));
335
336     int Me2 = spcTopology.sequenceID(MoleculeName("MeOH"), 1, ResidueName("MeOH"), ParticleName("Me1"));
337     int MeO2 = spcTopology.sequenceID(MoleculeName("MeOH"), 1, ResidueName("MeOH"), ParticleName("O2"));
338     int MeH2 = spcTopology.sequenceID(MoleculeName("MeOH"), 1, ResidueName("MeOH"), ParticleName("H3"));
339
340     /// \cond DO_NOT_DOCUMENT
341 #define SORT(i, j) (i < j) ? i : j, (i < j) ? j : i
342     interactions_reference[std::make_tuple(SORT(Ow, H1))]     = ohBond;
343     interactions_reference[std::make_tuple(SORT(Ow, H2))]     = ohBond;
344     interactions_reference[std::make_tuple(SORT(MeO1, MeH1))] = ohBondMethanol;
345     interactions_reference[std::make_tuple(SORT(MeO1, Me1))]  = ometBond;
346     interactions_reference[std::make_tuple(SORT(MeO2, MeH2))] = ohBondMethanol;
347     interactions_reference[std::make_tuple(SORT(MeO2, Me2))]  = ometBond;
348 #undef SORT
349     /// \endcond
350
351     EXPECT_TRUE(std::equal(
352             begin(interactions_reference), end(interactions_reference), begin(interactions_test)));
353 }
354
355 TEST(NBlibTest, TopologyListedInteractionsMultipleTypes)
356 {
357     // Todo: add an angle type here
358
359     Molecule water    = WaterMoleculeBuilder{}.waterMolecule();
360     Molecule methanol = MethanolMoleculeBuilder{}.methanolMolecule();
361
362     CubicBondType testBond(1., 1., 1.);
363     HarmonicAngle testAngle(Degrees(1), 1);
364
365     water.addInteraction(ParticleName("H1"), ParticleName("H2"), testBond);
366     water.addInteraction(ParticleName("H1"), ParticleName("Oxygen"), ParticleName("H2"), testAngle);
367
368     ParticleTypesInteractions nbInteractions;
369     std::vector<std::string>  particleTypeNames = { "Ow", "H", "OMet", "CMet" };
370     for (const auto& name : particleTypeNames)
371     {
372         nbInteractions.add(ParticleTypeName(name), C6(0), C12(0));
373     }
374
375     TopologyBuilder topologyBuilder;
376     topologyBuilder.addMolecule(water, 1);
377     topologyBuilder.addMolecule(methanol, 1);
378     topologyBuilder.addParticleTypesInteractions(nbInteractions);
379
380     Topology topology = topologyBuilder.buildTopology();
381
382     auto  interactionData   = topology.getInteractionData();
383     auto& harmonicBonds     = pickType<HarmonicBondType>(interactionData);
384     auto& cubicBonds        = pickType<CubicBondType>(interactionData);
385     auto& angleInteractions = pickType<HarmonicAngle>(interactionData);
386
387     HarmonicBondType              ohBond(1., 1.);
388     HarmonicBondType              ohBondMethanol(1.01, 1.02);
389     HarmonicBondType              ometBond(1.1, 1.2);
390     std::vector<HarmonicBondType> harmonicBondsReference{ ohBond, ohBondMethanol, ometBond };
391
392     EXPECT_EQ(harmonicBonds.parameters, harmonicBondsReference);
393
394     int H1 = topology.sequenceID(MoleculeName("SOL"), 0, ResidueName("SOL"), ParticleName("H1"));
395     int H2 = topology.sequenceID(MoleculeName("SOL"), 0, ResidueName("SOL"), ParticleName("H2"));
396     int Ow = topology.sequenceID(MoleculeName("SOL"), 0, ResidueName("SOL"), ParticleName("Oxygen"));
397
398     int Me1 = topology.sequenceID(MoleculeName("MeOH"), 0, ResidueName("MeOH"), ParticleName("Me1"));
399     int MeO1 = topology.sequenceID(MoleculeName("MeOH"), 0, ResidueName("MeOH"), ParticleName("O2"));
400     int MeH1 = topology.sequenceID(MoleculeName("MeOH"), 0, ResidueName("MeOH"), ParticleName("H3"));
401
402     std::vector<CubicBondType>                   cubicBondsReference{ testBond };
403     std::vector<InteractionIndex<CubicBondType>> cubicIndicesReference{
404         { std::min(H1, H2), std::max(H1, H2), 0 }
405     };
406     EXPECT_EQ(cubicBondsReference, cubicBonds.parameters);
407     EXPECT_EQ(cubicIndicesReference, cubicBonds.indices);
408
409     HarmonicAngle                                methanolAngle(Degrees(108.52), 397.5);
410     std::vector<HarmonicAngle>                   angleReference{ testAngle, methanolAngle };
411     std::vector<InteractionIndex<HarmonicAngle>> angleIndicesReference{
412         { std::min(H1, H2), Ow, std::max(H1, H2), 0 }, { std::min(MeH1, MeO1), Me1, std::max(MeO1, MeH1), 1 }
413     };
414     EXPECT_EQ(angleReference, angleInteractions.parameters);
415     EXPECT_EQ(angleIndicesReference, angleInteractions.indices);
416 }
417
418 TEST(NBlibTest, TopologyInvalidParticleInInteractionThrows)
419 {
420     Molecule water    = WaterMoleculeBuilder{}.waterMolecule();
421     Molecule methanol = MethanolMoleculeBuilder{}.methanolMolecule();
422
423     HarmonicBondType testBond(1., 1.);
424
425     // Invalid input: no particle named "Iron" in molecule water
426     water.addInteraction(ParticleName("H1"), ParticleName("Iron"), testBond);
427
428     ParticleTypesInteractions nbInteractions;
429     std::vector<std::string>  particleTypeNames = { "Ow", "H", "OMet", "CMet" };
430     for (const auto& name : particleTypeNames)
431     {
432         nbInteractions.add(ParticleTypeName(name), C6(0), C12(0));
433     }
434
435     TopologyBuilder topologyBuilder;
436     topologyBuilder.addMolecule(water, 1);
437     topologyBuilder.addMolecule(methanol, 1);
438     topologyBuilder.addParticleTypesInteractions(nbInteractions);
439
440     EXPECT_THROW(topologyBuilder.buildTopology(), InputException);
441 }
442
443
444 TEST(NBlibTest, toGmxExclusionBlockWorks)
445 {
446     std::vector<std::tuple<int, int>> testInput{ { 0, 0 }, { 0, 1 }, { 0, 2 }, { 1, 0 }, { 1, 1 },
447                                                  { 1, 2 }, { 2, 0 }, { 2, 1 }, { 2, 2 } };
448
449     std::vector<gmx::ExclusionBlock> reference;
450
451     gmx::ExclusionBlock localBlock;
452     localBlock.atomNumber.push_back(0);
453     localBlock.atomNumber.push_back(1);
454     localBlock.atomNumber.push_back(2);
455
456     reference.push_back(localBlock);
457     reference.push_back(localBlock);
458     reference.push_back(localBlock);
459
460     std::vector<gmx::ExclusionBlock> probe = toGmxExclusionBlock(testInput);
461
462     ASSERT_EQ(reference.size(), probe.size());
463     for (size_t i = 0; i < reference.size(); ++i)
464     {
465         ASSERT_EQ(reference[i].atomNumber.size(), probe[i].atomNumber.size());
466         for (size_t j = 0; j < reference[i].atomNumber.size(); ++j)
467         {
468             EXPECT_EQ(reference[i].atomNumber[j], probe[i].atomNumber[j]);
469         }
470     }
471 }
472
473 } // namespace
474 } // namespace test
475 } // namespace nblib