SYCL: Avoid using no_init read accessor in rocFFT
[alexxy/gromacs.git] / src / gromacs / mdlib / groupcoord.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2008, The GROMACS development team.
6  * Copyright (c) 2012,2014,2015,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37
38 /*! \libinternal \file
39  * \brief Assemble atomic positions of a (small) subset of atoms and distribute to all nodes.
40  *
41  *  This file contains functions to assemble the positions of a subset of the
42  *  atoms and to do operations on it like determining the center of mass, or
43  *  doing translations and rotations. These functions are useful when
44  *  a subset of the positions needs to be compared to some set of reference
45  *  positions, as e.g. done for essential dynamics.
46  *
47  * \inlibraryapi
48  */
49
50 #include <stdio.h>
51
52 #include "gromacs/math/vectypes.h"
53 #include "gromacs/utility/basedefinitions.h"
54
55 class gmx_ga2la_t;
56 struct t_commrec;
57
58 /*! \brief Select local atoms of a group.
59  *
60  * Selects the indices of local atoms of a group and stores them in anrs_loc[0..nr_loc].
61  * If you need the positions of the group's atoms on all nodes, provide a coll_ind[0..nr]
62  * array and pass it on to communicate_group_positions. Thus the collective array
63  * will always have the same atom order (ascending indices).
64  *
65  *  \param[in]     ga2la      Global to local atom index conversion data.
66  *  \param[in]     nr         The total number of atoms that the group contains.
67  *  \param[in]     anrs       The global atom number of the group's atoms.
68  *  \param[out]    nr_loc     The number of group atoms present on the local node.
69  *  \param[out]    anrs_loc   The local atom numbers of the group.
70  *  \param[in,out] nalloc_loc Local allocation size of anrs_loc array.
71  *  \param[out]    coll_ind   If not NULL this array must be of size nr. It stores
72  *                            for each local atom where it belongs in the global
73  *                            (collective) array such that it can be gmx_summed
74  *                            in the communicate_group_positions routine.
75  */
76
77 extern void dd_make_local_group_indices(const gmx_ga2la_t* ga2la,
78                                         int                nr,
79                                         int                anrs[],
80                                         int*               nr_loc,
81                                         int*               anrs_loc[],
82                                         int*               nalloc_loc,
83                                         int                coll_ind[]);
84
85
86 /*! \brief Assemble local positions into a collective array present on all nodes.
87  *
88  * Communicate the positions of the group's atoms such that every node has all of
89  * them. Unless running on huge number of cores, this is not a big performance impact
90  * as long as the collective subset [0..nr] is kept small. The atom indices are
91  * retrieved from anrs_loc[0..nr_loc]. If you call the routine for the serial case,
92  * provide an array coll_ind[i] = i for i in 1..nr.
93  *
94  * If shifts != NULL, the PBC representation of each atom is chosen such that a
95  * continuous trajectory results. Therefore, if the group is whole at the start
96  * of the simulation, it will always stay whole.
97  * If shifts = NULL, the group positions are not made whole again, but assembled
98  * and distributed to all nodes. The variables marked "optional" are not used in
99  * that case.
100  *
101  * \param[in]     cr           Pointer to MPI communication data.
102  * \param[out]    xcoll        Collective array of positions, identical on all nodes
103  *                             after this routine has been called.
104  * \param[in,out] shifts       Collective array of shifts for xcoll, needed to make
105  *                             the group whole. This array remembers the shifts
106  *                             since the start of the simulation (where the group
107  *                             is whole) and must therefore not be changed outside
108  *                             of this routine! If NULL, the group will not be made
109  *                             whole and the optional variables are ignored.
110  * \param[out]    extra_shifts Extra shifts since last time step, only needed as
111  *                             buffer variable [0..nr] (optional).
112  * \param[in]     bNS          Neighbor searching / domain re-decomposition has been
113  *                             performed at the begin of this time step such that
114  *                             the shifts have changed and need to be updated
115  *                             (optional).
116  * \param[in]     x_loc        Pointer to the local atom positions this node has.
117  * \param[in]     nr           Total number of atoms in the group.
118  * \param[in]     nr_loc       Number of group atoms on the local node.
119  * \param[in]     anrs_loc     Array of the local atom indices.
120  * \param[in]     coll_ind     This array of size nr stores for each local atom where
121  *                             it belongs in the collective array so that the local
122  *                             contributions can be gmx_summed. It is provided by
123  *                             dd_make_local_group_indices.
124  * \param[in,out] xcoll_old    Positions from the last time step, used to make the
125  *                             group whole (optional).
126  * \param[in]     box          Simulation box matrix, needed to shift xcoll such that
127  *                             the group becomes whole (optional).
128  */
129 extern void communicate_group_positions(const t_commrec* cr,
130                                         rvec*            xcoll,
131                                         ivec*            shifts,
132                                         ivec*            extra_shifts,
133                                         gmx_bool         bNS,
134                                         const rvec*      x_loc,
135                                         int              nr,
136                                         int              nr_loc,
137                                         const int*       anrs_loc,
138                                         const int*       coll_ind,
139                                         rvec*            xcoll_old,
140                                         const matrix     box);
141
142 /*! \brief Calculates the center of the positions x locally.
143  *
144  * Calculates the center of mass (if masses are given in the weight array) or
145  * the geometrical center (if NULL is passed as weight).
146  *
147  * \param[in]   x            Positions.
148  * \param[in]   weight       Can be NULL or an array of weights. If masses are
149  *                           given as weights, the COM is calculated.
150  * \param[in]   nr           Number of positions and weights if present.
151  * \param[out]  center       The (weighted) center of the positions.
152  *
153  */
154 extern void get_center(rvec x[], real weight[], int nr, rvec center);
155
156
157 /*! \brief Calculates the sum of the positions x locally.
158  *
159  * Calculates the (weighted) sum of position vectors and returns the sum of
160  * weights, which is needed when local contributions shall be summed to a
161  * global weighted center.
162  *
163  * \param[in]   x            Array of positions.
164  * \param[in]   weight       Can be NULL or an array of weights.
165  * \param[in]   nr           Number of positions and weights if present.
166  * \param[out]  dsumvec      The (weighted) sum of the positions.
167  * \return Sum of weights.
168  *
169  */
170 extern double get_sum_of_positions(rvec x[], real weight[], int nr, dvec dsumvec);
171
172
173 /*! \brief Calculates the global center of all local arrays x.
174  *
175  * Get the center from local positions [0..nr_loc], this involves communication.
176  * Not that the positions must already have the correct PBC representation. Use
177  * this routine if no collective coordinates are assembled from which the center
178  * could be calculated without communication.
179  *
180  * \param[in]   cr           Pointer to MPI communication data.
181  * \param[in]   x_loc        Array of local positions [0..nr_loc].
182  * \param[in]   weight_loc   Array of local weights, these are the masses if the
183  *                           center of mass is to be calculated.
184  * \param[in]   nr_loc       The number of positions on the local node.
185  * \param[in]   nr_group     The number of positions in the whole group. Since
186  *                           this is known anyway, we do not need to communicate
187  *                           and sum nr_loc if we pass it over.
188  * \param[out]  center       The (weighted) center of all x_loc from all the
189  *                           nodes.
190  */
191 extern void get_center_comm(const t_commrec* cr, rvec x_loc[], real weight_loc[], int nr_loc, int nr_group, rvec center);
192
193
194 /*! \brief Translate positions.
195  *
196  * Add a translation vector to the positions x.
197  *
198  * \param[in,out] x          Array of positions.
199  * \param[in]     nr         Number of entries in the position array.
200  * \param[in]     transvec   Translation vector to be added to all positions.
201  *
202  */
203 extern void translate_x(rvec x[], int nr, const rvec transvec);
204
205
206 /*! \brief Rotate positions.
207  *
208  * Rotate the positions with the rotation matrix.
209  *
210  * \param[in,out] x          Array of positions.
211  * \param[in]     nr         Number of entries in the position array.
212  * \param[in]     rmat       Rotation matrix to operate on all positions.
213  *
214  */
215 extern void rotate_x(rvec x[], int nr, matrix rmat);