SYCL: Use acc.bind(cgh) instead of cgh.require(acc)
[alexxy/gromacs.git] / src / gromacs / mdlib / wholemoleculetransform.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
36 #include "gmxpre.h"
37
38 #include "wholemoleculetransform.h"
39
40 #include "gromacs/domdec/ga2la.h"
41 #include "gromacs/topology/mtop_util.h"
42
43 namespace gmx
44 {
45
46 WholeMoleculeTransform::WholeMoleculeTransform(const gmx_mtop_t& mtop,
47                                                const PbcType     pbcType,
48                                                const bool        useAtomReordering) :
49     pbcType_(pbcType)
50 {
51     gmx_localtop_t localTop(mtop.ffparams);
52
53     gmx_mtop_generate_local_top(mtop, &localTop, false);
54
55     graph_ = mk_graph(localTop.idef, mtop.natoms);
56
57     wholeMoleculeCoordinates_.resize(mtop.natoms);
58
59     if (useAtomReordering)
60     {
61         // Store a copy of the global edges
62         globalEdgeAtomBegin_       = graph_.edgeAtomBegin;
63         graphGlobalAtomOrderEdges_ = graph_.edges;
64         // Reset the graph range to cover the whole system
65         graph_.edgeAtomBegin = 0;
66         graph_.edgeAtomEnd   = graph_.shiftAtomEnd;
67
68         // Resize the edge color list for potential addition of non-connected atoms
69         graph_.edgeColor.resize(graph_.edgeAtomEnd);
70     }
71 }
72
73 void WholeMoleculeTransform::updateAtomOrder(ArrayRef<const int> globalAtomIndices, const gmx_ga2la_t& ga2la)
74 {
75     GMX_RELEASE_ASSERT(!graphGlobalAtomOrderEdges_.empty(),
76                        "We need the edges in the global atom order");
77
78     GMX_RELEASE_ASSERT(globalAtomIndices.ssize() == graph_.shiftAtomEnd,
79                        "The number of indices should match the number of atoms we shift");
80
81     const int globalEdgeAtomEnd = globalEdgeAtomBegin_ + graphGlobalAtomOrderEdges_.size();
82     graph_.edges.clear();
83     for (const int globalAtomIndex : globalAtomIndices)
84     {
85         if (globalAtomIndex >= globalEdgeAtomBegin_ && globalAtomIndex < globalEdgeAtomEnd)
86         {
87             // Get the list of global atoms linked to this atom and add their local indices
88             const ArrayRef<const int> globalList =
89                     graphGlobalAtomOrderEdges_[globalAtomIndex - globalEdgeAtomBegin_];
90             graph_.edges.pushBackListOfSize(globalList.size());
91             ArrayRef<int> lastList = graph_.edges.back();
92             std::transform(globalList.begin(), globalList.end(), lastList.begin(), [&ga2la](int i) -> int {
93                 return ga2la.find(i)->la;
94             });
95         }
96         else
97         {
98             // The atom has no edges, push back an empty list
99             graph_.edges.pushBackListOfSize(0);
100         }
101     }
102
103     GMX_RELEASE_ASSERT(int(graph_.edges.size()) == graph_.shiftAtomEnd,
104                        "We should have as many lists of edges as the system (shift) size");
105 }
106
107 void WholeMoleculeTransform::updateForAtomPbcJumps(ArrayRef<const RVec> x, const matrix box)
108 {
109     mk_mshift(nullptr, &graph_, pbcType_, box, as_rvec_array(x.data()));
110 }
111
112 ArrayRef<const RVec> WholeMoleculeTransform::wholeMoleculeCoordinates(ArrayRef<const RVec> x, const matrix box)
113 {
114     shift_x(&graph_, box, as_rvec_array(x.data()), as_rvec_array(wholeMoleculeCoordinates_.data()));
115
116     return wholeMoleculeCoordinates_;
117 }
118
119 } // namespace gmx