SYCL: Use acc.bind(cgh) instead of cgh.require(acc)
[alexxy/gromacs.git] / src / gromacs / mdlib / wholemoleculetransform.h
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 /*! \libinternal \file
36  *
37  * \brief Declares the WholeMolecules class for generating whole molecules
38  *
39  * \author Berk Hess <hess@kth.se>
40  * \ingroup module_mdlib
41  * \inlibraryapi
42  */
43 #ifndef GMX_MDLIB_WHOLEMOLECULETRANSFORM_H
44 #define GMX_MDLIB_WHOLEMOLECULETRANSFORM_H
45
46 #include <vector>
47
48 #include "gromacs/math/vectypes.h"
49 #include "gromacs/pbcutil/mshift.h"
50 #include "gromacs/utility/arrayref.h"
51
52 class gmx_ga2la_t;
53 struct gmx_mtop_t;
54 enum class PbcType : int;
55
56 namespace gmx
57 {
58
59 /*! \libinternal
60  * \brief This class manages a coordinate buffer with molecules not split
61  * over periodic boundary conditions for use in force calculations
62  * which require whole molecules.
63  *
64  * Note: This class should not be used for computation of forces which
65  *       have virial contributions through shift forces.
66  */
67 class WholeMoleculeTransform
68 {
69 public:
70     /*! \brief Constructor
71      *
72      * \param[in] mtop               The global topology use for getting the connections between atoms
73      * \param[in] pbcType            The type of PBC
74      * \param[in] useAtomReordering  Whether we will use atom reordering
75      */
76     WholeMoleculeTransform(const gmx_mtop_t& mtop, PbcType pbcType, bool useAtomReordering);
77
78     /*! \brief Changes the atom order to the one provided
79      *
80      * This method is called after domain repartitioning.
81      * The object should have been constructed with \p useAtomReordering set to \p true.
82      *
83      * \param[in] globalAtomIndices  The global atom indices for the local atoms, size should be the system size
84      * \param[in] ga2la              Global to local atom index lookup (the inverse of \p globalAtomIndices)
85      */
86     void updateAtomOrder(ArrayRef<const int> globalAtomIndices, const gmx_ga2la_t& ga2la);
87
88     /*! \brief Updates the graph when atoms have been shifted by periodic vectors */
89     void updateForAtomPbcJumps(ArrayRef<const RVec> x, const matrix box);
90
91     /*! \brief Create and return coordinates with whole molecules for input coordinates \p x
92      *
93      * \param[in] x  Input coordinates, should not have periodic displacement compared
94      *               with the coordinates passed in the last call to \p updateForAtomPbcJumps().
95      * \param[in] box  The current periodic image vectors
96      *
97      * Note: this operation is not free. If you need whole molecules coordinates
98      * more than once during the force calculation, store the result and reuse it.
99      */
100     ArrayRef<const RVec> wholeMoleculeCoordinates(ArrayRef<const RVec> x, const matrix box);
101
102 private:
103     //! The type of PBC
104     PbcType pbcType_;
105     //! The graph
106     t_graph graph_;
107     //! The atom index at which graphGlobalAtomOrderEdges_ starts
108     int globalEdgeAtomBegin_;
109     //! The edges for the global atom order
110     ListOfLists<int> graphGlobalAtomOrderEdges_;
111     //! Buffer for storing coordinates for whole molecules
112     std::vector<RVec> wholeMoleculeCoordinates_;
113 };
114
115 } // namespace gmx
116
117 #endif