2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 * This tests that sample code can run
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>
47 #include "gromacs/utility/arrayref.h"
49 // The entire nblib public API can be included with a single header or individual components
50 // can be included via their respective headers.
51 #include "nblib/nblib.h"
55 // Create an argon particle with a name and a mass.
56 nblib::ParticleType argonAtom(nblib::ParticleTypeName("Ar"), nblib::Mass(39.94800));
57 // Create an argon molecule.
58 nblib::Molecule argonMolecule(nblib::MoleculeName("AR"));
59 // Add the argon particle to a molecule. The names are for bookkeeping and need not match.
60 argonMolecule.addParticle(nblib::ParticleName("Argon"), argonAtom);
61 // Define Lennard-Jones params for argon (parameters from gromos43A1).
62 nblib::C6 ArC6{ 0.0062647225 }; // C6 parameter
63 nblib::C12 ArC12{ 9.847044e-06 }; // C12 parameter
64 // Holder for non-bonded interactions.
65 nblib::ParticleTypesInteractions interactions;
66 // Add non-bonded interactions for argon.
67 interactions.add(argonAtom.name(), ArC6, ArC12);
68 // The TopologyBuilder builds the Topology!
69 nblib::TopologyBuilder topologyBuilder;
70 // Number of Argon particles (molecules) in the system.
71 int numParticles = 12;
72 // Add the requested number of argon molecules to a topology.
73 topologyBuilder.addMolecule(argonMolecule, numParticles);
74 // Add the argon interactions to the topology.
75 topologyBuilder.addParticleTypesInteractions(interactions);
76 // Build the topology.
77 nblib::Topology topology = topologyBuilder.buildTopology();
78 // The system needs a bounding box. Only cubic and rectangular boxes are supported.
79 nblib::Box box(6.05449);
80 // User defined coordinates.
81 std::vector<nblib::Vec3> coordinates = {
82 { 0.794, 1.439, 0.610 }, { 1.397, 0.673, 1.916 }, { 0.659, 1.080, 0.573 },
83 { 1.105, 0.090, 3.431 }, { 1.741, 1.291, 3.432 }, { 1.936, 1.441, 5.873 },
84 { 0.960, 2.246, 1.659 }, { 0.382, 3.023, 2.793 }, { 0.053, 4.857, 4.242 },
85 { 2.655, 5.057, 2.211 }, { 4.114, 0.737, 0.614 }, { 5.977, 5.104, 5.217 },
87 // User defined velocities.
88 std::vector<nblib::Vec3> velocities = {
89 { 0.0055, -0.1400, 0.2127 }, { 0.0930, -0.0160, -0.0086 }, { 0.1678, 0.2476, -0.0660 },
90 { 0.1591, -0.0934, -0.0835 }, { -0.0317, 0.0573, 0.1453 }, { 0.0597, 0.0013, -0.0462 },
91 { 0.0484, -0.0357, 0.0168 }, { 0.0530, 0.0295, -0.2694 }, { -0.0550, -0.0896, 0.0494 },
92 { -0.0799, -0.2534, -0.0079 }, { 0.0436, -0.1557, 0.1849 }, { -0.0214, 0.0446, 0.0758 },
94 // Force buffer initialization for each particle.
95 std::vector<nblib::Vec3> forces = {
96 { 0.0000, 0.0000, 0.0000 }, { 0.0000, 0.0000, 0.0000 }, { 0.0000, 0.0000, 0.0000 },
97 { 0.0000, 0.0000, 0.0000 }, { 0.0000, 0.0000, 0.0000 }, { 0.0000, 0.0000, 0.0000 },
98 { 0.0000, 0.0000, 0.0000 }, { 0.0000, 0.0000, 0.0000 }, { 0.0000, 0.0000, 0.0000 },
99 { 0.0000, 0.0000, 0.0000 }, { 0.0000, 0.0000, 0.0000 }, { 0.0000, 0.0000, 0.0000 },
101 // A simulation state contains all the molecular information about the system.
102 nblib::SimulationState simState(coordinates, velocities, forces, box, topology);
103 // Kernel options are flags needed for force calculation.
104 nblib::NBKernelOptions options = nblib::NBKernelOptions();
105 // Use a simple cutoff rule for Coulomb
106 options.coulombType = nblib::CoulombType::Cutoff;
107 // Some performance flags can be set a run time
108 options.nbnxmSimd = nblib::SimdKernels::SimdNo;
109 // The force calculator contains all the data needed to compute forces.
110 nblib::ForceCalculator forceCalculator(simState, options);
111 // Integration requires masses, positions, and forces
112 nblib::LeapFrog integrator(simState.topology(), simState.box());
113 // Print some diagnostic info
114 printf("initial forces on particle 0: x %4f y %4f z %4f\n", forces[0][0], forces[0][1], forces[0][2]);
115 // The forces are computed for the user
116 gmx::ArrayRef<nblib::Vec3> userForces(simState.forces());
117 forceCalculator.compute(simState.coordinates(), userForces);
118 // Print some diagnostic info
119 printf(" final forces on particle 0: x %4f y %4f z %4f\n",
123 // User may modify forces stored in simState.forces() if needed
124 // Print some diagnostic info
125 printf("initial position of particle 0: x %4f y %4f z %4f\n",
126 simState.coordinates()[0][0],
127 simState.coordinates()[0][1],
128 simState.coordinates()[0][2]);
129 // Integrate with a time step of 1 fs
130 integrator.integrate(1.0, simState.coordinates(), simState.velocities(), simState.forces());
131 // Print some diagnostic info
133 printf(" final position of particle 0: x %4f y %4f z %4f\n",
134 simState.coordinates()[0][0],
135 simState.coordinates()[0][1],
136 simState.coordinates()[0][2]);