Introduce AngleInteractionType
[alexxy/gromacs.git] / api / nblib / tests / testsystems.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 nblib test systems
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  * \author Artem Zhmurov <zhmurov@gmail.com>
44  */
45 #include "nblib/tests/testsystems.h"
46 #include "nblib/exception.h"
47
48 namespace nblib
49 {
50
51 //! User class to hold ParticleTypes
52 //! Note: this is not part of NBLIB, users should write their own
53 class ParticleLibrary
54 {
55 public:
56     ParticleLibrary()
57     {
58         ParticleType Ow(ParticleTypeName("Ow"), Mass(15.99940));
59         ParticleType H(ParticleTypeName("H"), Mass(1.008));
60         ParticleType OMet(ParticleTypeName("OMet"), Mass(15.999));
61         ParticleType CMet(ParticleTypeName("CMet"), Mass(15.035));
62         ParticleType Ar_gromos(ParticleTypeName("Ar_gromos"), Mass(39.94800));
63         ParticleType Ar_opls(ParticleTypeName("Ar_opls"), Mass(39.94800));
64
65         particles_.insert(std::make_pair(Ow.name(), Ow));
66         particles_.insert(std::make_pair(H.name(), H));
67         particles_.insert(std::make_pair(OMet.name(), OMet));
68         particles_.insert(std::make_pair(CMet.name(), CMet));
69         particles_.insert(std::make_pair(Ar_gromos.name(), Ar_gromos));
70         particles_.insert(std::make_pair(Ar_opls.name(), Ar_opls));
71
72         c6_[Ow.name()]        = 0.0026173456;
73         c6_[H.name()]         = 0;
74         c6_[OMet.name()]      = 0.0022619536;
75         c6_[CMet.name()]      = 0.0088755241;
76         c6_[Ar_gromos.name()] = 0.0062647225;
77         c6_[Ar_opls.name()]   = 0.0058560692;
78
79         c12_[Ow.name()]        = 2.634129e-06;
80         c12_[H.name()]         = 0;
81         c12_[OMet.name()]      = 1.505529e-06;
82         c12_[CMet.name()]      = 2.0852922e-05;
83         c12_[Ar_gromos.name()] = 9.847044e-06;
84         c12_[Ar_opls.name()]   = 8.203193e-06;
85     }
86
87     //! Get particle type using the string identifier
88     [[nodiscard]] ParticleType type(const std::string& particleTypeName) const
89     {
90         return particles_.at(ParticleTypeName(particleTypeName));
91     }
92
93     //! Get C6 parameter of a given particle
94     [[nodiscard]] C6 c6(const ParticleName& particleName) const
95     {
96         return c6_.at(ParticleTypeName(particleName.value()));
97     }
98
99     //! Get C12 parameter of a given particle
100     [[nodiscard]] C12 c12(const ParticleName& particleName) const
101     {
102         return c12_.at(ParticleTypeName(particleName.value()));
103     }
104
105 private:
106     std::map<ParticleTypeName, ParticleType> particles_;
107     std::map<ParticleTypeName, C6>           c6_;
108     std::map<ParticleTypeName, C12>          c12_;
109 };
110
111 std::unordered_map<std::string, Charge> Charges{ { "Ow", Charge(-0.82) },
112                                                  { "Hw", Charge(+0.41) },
113                                                  { "OMet", Charge(-0.574) },
114                                                  { "CMet", Charge(+0.176) },
115                                                  { "HMet", Charge(+0.398) } };
116
117 WaterMoleculeBuilder::WaterMoleculeBuilder() : water_(MoleculeName("SOL"))
118 {
119     ParticleLibrary plib;
120
121     //! Add the particles
122     water_.addParticle(ParticleName("Oxygen"), Charges.at("Ow"), plib.type("Ow"));
123     water_.addParticle(ParticleName("H1"), Charges.at("Hw"), plib.type("H"));
124     water_.addParticle(ParticleName("H2"), Charges.at("Hw"), plib.type("H"));
125
126     HarmonicBondType ohBond(1., 1.);
127     water_.addInteraction(ParticleName("Oxygen"), ParticleName("H1"), ohBond);
128     water_.addInteraction(ParticleName("Oxygen"), ParticleName("H2"), ohBond);
129 }
130
131 Molecule WaterMoleculeBuilder::waterMolecule()
132 {
133     addExclusionsFromNames();
134     return water_;
135 }
136
137 Molecule WaterMoleculeBuilder::waterMoleculeWithoutExclusions()
138 {
139     return water_;
140 }
141
142 void WaterMoleculeBuilder::addExclusionsFromNames()
143 {
144     water_.addExclusion(ParticleName("H1"), ParticleName("Oxygen"));
145     water_.addExclusion(ParticleName("H2"), ParticleName("Oxygen"));
146     water_.addExclusion(ParticleName("H1"), ParticleName("H2"));
147 }
148
149 MethanolMoleculeBuilder::MethanolMoleculeBuilder() : methanol_(MoleculeName("MeOH"))
150 {
151     ParticleLibrary library;
152
153     //! Add the particles
154     methanol_.addParticle(ParticleName("Me1"), Charges.at("CMet"), library.type("CMet"));
155     methanol_.addParticle(ParticleName("O2"), Charges.at("OMet"), library.type("OMet"));
156     methanol_.addParticle(ParticleName("H3"), Charges.at("HMet"), library.type("H"));
157
158     // Add the exclusions
159     methanol_.addExclusion(ParticleName("Me1"), ParticleName("O2"));
160     methanol_.addExclusion(ParticleName("Me1"), ParticleName("H3"));
161     methanol_.addExclusion(ParticleName("H3"), ParticleName("O2"));
162
163     HarmonicBondType ohBond(1.01, 1.02);
164     methanol_.addInteraction(ParticleName("O2"), ParticleName("H3"), ohBond);
165
166     HarmonicBondType ometBond(1.1, 1.2);
167     methanol_.addInteraction(ParticleName("O2"), ParticleName("Me1"), ometBond);
168
169     HarmonicAngle ochAngle(397.5, Degrees(108.52));
170     methanol_.addInteraction(ParticleName("O2"), ParticleName("Me1"), ParticleName("H3"), ochAngle);
171 }
172
173 Molecule MethanolMoleculeBuilder::methanolMolecule()
174 {
175     return methanol_;
176 }
177
178
179 Topology WaterTopologyBuilder::buildTopology(int numMolecules)
180 {
181     ParticleLibrary library;
182
183     ParticleTypesInteractions interactions;
184     std::vector<std::string>  typeNames = { "Ow", "H" };
185     for (const auto& name : typeNames)
186     {
187         interactions.add(
188                 ParticleTypeName(name), library.c6(ParticleName(name)), library.c12(ParticleName(name)));
189     }
190
191     // Add some molecules to the topology
192     TopologyBuilder topologyBuilder;
193     topologyBuilder.addMolecule(water(), numMolecules);
194
195     // Add non-bonded interaction information
196     topologyBuilder.addParticleTypesInteractions(interactions);
197
198     Topology topology = topologyBuilder.buildTopology();
199     return topology;
200 }
201
202 Molecule WaterTopologyBuilder::water()
203 {
204     return waterMolecule_.waterMolecule();
205 }
206
207 Topology SpcMethanolTopologyBuilder::buildTopology(int numWater, int numMethanol)
208 {
209     ParticleLibrary library;
210
211     ParticleTypesInteractions interactions;
212     std::vector<std::string>  typeNames = { "Ow", "H", "OMet", "CMet" };
213     for (const auto& name : typeNames)
214     {
215         interactions.add(
216                 ParticleTypeName(name), library.c6(ParticleName(name)), library.c12(ParticleName(name)));
217     }
218
219     // Add some molecules to the topology
220     TopologyBuilder topologyBuilder;
221     topologyBuilder.addMolecule(methanol(), numMethanol);
222     topologyBuilder.addMolecule(water(), numWater);
223
224     // Add non-bonded interaction information
225     topologyBuilder.addParticleTypesInteractions(interactions);
226
227     Topology topology = topologyBuilder.buildTopology();
228     return topology;
229 }
230
231 Molecule SpcMethanolTopologyBuilder::methanol()
232 {
233     return methanolMolecule_.methanolMolecule();
234 }
235
236 Molecule SpcMethanolTopologyBuilder::water()
237 {
238     return waterMolecule_.waterMolecule();
239 }
240
241 ArgonTopologyBuilder::ArgonTopologyBuilder(const int& numParticles, const fftypes forceField)
242 {
243     ParticleLibrary library;
244
245     ParticleTypesInteractions nbinteractions;
246
247     Molecule argonMolecule(MoleculeName("AR"));
248
249     std::string particleName;
250
251     if (forceField == fftypes::GROMOS43A1)
252     {
253         particleName = "Ar_gromos";
254     }
255     else if (forceField == fftypes::OPLSA)
256     {
257         particleName = "Ar_opls";
258     }
259
260     nbinteractions.add(ParticleTypeName(particleName),
261                        library.c6(ParticleName(particleName)),
262                        library.c12(ParticleName(particleName)));
263     argonMolecule.addParticle(ParticleName("AR"), library.type(particleName));
264
265     topologyBuilder_.addMolecule(argonMolecule, numParticles);
266     topologyBuilder_.addParticleTypesInteractions((nbinteractions));
267 }
268
269 Topology ArgonTopologyBuilder::argonTopology()
270 {
271     return topologyBuilder_.buildTopology();
272 }
273
274 ArgonSimulationStateBuilder::ArgonSimulationStateBuilder(const fftypes forceField) :
275     box_(6.05449), topology_(ArgonTopologyBuilder(12, forceField).argonTopology())
276 {
277
278     coordinates_ = {
279         { 0.794, 1.439, 0.610 }, { 1.397, 0.673, 1.916 }, { 0.659, 1.080, 0.573 },
280         { 1.105, 0.090, 3.431 }, { 1.741, 1.291, 3.432 }, { 1.936, 1.441, 5.873 },
281         { 0.960, 2.246, 1.659 }, { 0.382, 3.023, 2.793 }, { 0.053, 4.857, 4.242 },
282         { 2.655, 5.057, 2.211 }, { 4.114, 0.737, 0.614 }, { 5.977, 5.104, 5.217 },
283     };
284
285     velocities_ = {
286         { 0.0055, -0.1400, 0.2127 },   { 0.0930, -0.0160, -0.0086 }, { 0.1678, 0.2476, -0.0660 },
287         { 0.1591, -0.0934, -0.0835 },  { -0.0317, 0.0573, 0.1453 },  { 0.0597, 0.0013, -0.0462 },
288         { 0.0484, -0.0357, 0.0168 },   { 0.0530, 0.0295, -0.2694 },  { -0.0550, -0.0896, 0.0494 },
289         { -0.0799, -0.2534, -0.0079 }, { 0.0436, -0.1557, 0.1849 },  { -0.0214, 0.0446, 0.0758 },
290     };
291     forces_ = {
292         { 0.0000, 0.0000, 0.0000 }, { 0.0000, 0.0000, 0.0000 }, { 0.0000, 0.0000, 0.0000 },
293         { 0.0000, 0.0000, 0.0000 }, { 0.0000, 0.0000, 0.0000 }, { 0.0000, 0.0000, 0.0000 },
294         { 0.0000, 0.0000, 0.0000 }, { 0.0000, 0.0000, 0.0000 }, { 0.0000, 0.0000, 0.0000 },
295         { 0.0000, 0.0000, 0.0000 }, { 0.0000, 0.0000, 0.0000 }, { 0.0000, 0.0000, 0.0000 },
296     };
297 }
298
299 void ArgonSimulationStateBuilder::setCoordinate(int particleNum, int dimension, real value)
300 {
301     if (dimension < 0 or dimension > 2)
302     {
303         throw InputException("Must provide a valid dimension\n");
304     }
305     coordinates_.at(particleNum)[dimension] = value;
306 }
307
308 void ArgonSimulationStateBuilder::setVelocity(int particleNum, int dimension, real value)
309 {
310     if (dimension < 0 or dimension > 2)
311     {
312         throw InputException("Must provide a valid dimension\n");
313     }
314     velocities_.at(particleNum)[dimension] = value;
315 }
316
317 SimulationState ArgonSimulationStateBuilder::setupSimulationState()
318 {
319     return SimulationState(coordinates_, velocities_, forces_, box_, topology_);
320 }
321
322 const Topology& ArgonSimulationStateBuilder::topology() const
323 {
324     return topology_;
325 }
326
327 Box& ArgonSimulationStateBuilder::box()
328 {
329     return box_;
330 }
331
332 std::vector<Vec3>& ArgonSimulationStateBuilder::coordinates()
333 {
334     return coordinates_;
335 }
336
337 std::vector<Vec3>& ArgonSimulationStateBuilder::velocities()
338 {
339     return velocities_;
340 }
341
342 SpcMethanolSimulationStateBuilder::SpcMethanolSimulationStateBuilder() :
343     box_(3.01000), topology_(SpcMethanolTopologyBuilder().buildTopology(1, 1))
344 {
345     coordinates_ = {
346         { 1.970, 1.460, 1.209 }, // Me1
347         { 1.978, 1.415, 1.082 }, // O2
348         { 1.905, 1.460, 1.030 }, // H3
349         { 1.555, 1.511, 0.703 }, // Ow
350         { 1.498, 1.495, 0.784 }, // Hw1
351         { 1.496, 1.521, 0.623 }, // Hw2
352     };
353
354     velocities_ = {
355         { -0.8587, -0.1344, -0.0643 }, { 0.0623, -0.1787, 0.0036 }, { -0.5020, -0.9564, 0.0997 },
356         { 0.869, 1.245, 1.665 },       { 0.169, 0.275, 1.565 },     { 0.269, 2.275, 1.465 },
357     };
358
359     forces_ = {
360         { 0.000, 0.000, 0.000 }, { 0.000, 0.000, 0.000 }, { 0.000, 0.000, 0.000 },
361         { 0.000, 0.000, 0.000 }, { 0.000, 0.000, 0.000 }, { 0.000, 0.000, 0.000 },
362     };
363 }
364
365 SimulationState SpcMethanolSimulationStateBuilder::setupSimulationState()
366 {
367     return SimulationState(coordinates_, velocities_, forces_, box_, topology_);
368 }
369
370 std::vector<Vec3>& SpcMethanolSimulationStateBuilder::coordinates()
371 {
372     return coordinates_;
373 }
374
375 std::vector<Vec3>& SpcMethanolSimulationStateBuilder::velocities()
376 {
377     return velocities_;
378 }
379
380 } // namespace nblib