f72d6b46d1fecc2059e7cba2479afcace1ca0a5d
[alexxy/gromacs.git] / api / nblib / integrator.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2020, 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 /*! \inpublicapi \file
36  * \brief
37  * Implements nblib integrator
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 #ifndef NBLIB_INTEGRATOR_H
46 #define NBLIB_INTEGRATOR_H
47
48 #include <vector>
49
50 #include "nblib/box.h"
51 #include "nblib/vector.h"
52
53 namespace gmx
54 {
55 template<typename T>
56 class ArrayRef;
57 } // namespace gmx
58
59 namespace nblib
60 {
61
62 class Topology;
63
64 /*! \brief Simple integrator
65  *
66  */
67 class LeapFrog final
68 {
69 public:
70     /*! \brief Constructor.
71      *
72      * \param[in] topology  Topology object to build list of inverse masses.
73      * \param[in] box       Box object for ensuring that coordinates remain within bounds
74      */
75     LeapFrog(const Topology& topology, const Box& box);
76
77     /*! \brief Integrate
78      *
79      * Integrates the equation of motion using Leap-Frog algorithm.
80      * Updates coordinates and velocities.
81      *
82      * \param[in]  dt          Timestep.
83      * \param[out] coordinates Coordinate array that would be modified in-place.
84      * \param[out] velocities  Velocity array that would be modified in-place.
85      * \param[in]  forces      Force array to be read.
86      *
87      */
88     void integrate(real                      dt,
89                    gmx::ArrayRef<Vec3>       coordinates,
90                    gmx::ArrayRef<Vec3>       velocities,
91                    gmx::ArrayRef<const Vec3> forces);
92
93 private:
94     //! 1/mass for all atoms
95     std::vector<real> inverseMasses_;
96     //! Box for PBC conformity
97     Box box_;
98 };
99
100 } // namespace nblib
101
102 #endif // NBLIB_INTEGRATOR_H