bdbf8bee30aea76b107e8788b314d10357a00ce9
[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, 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 struct gmx_ga2la_t;
56 struct t_commrec;
57
58 #ifdef __cplusplus
59 extern "C" {
60 #endif
61
62 /*! \brief Select local atoms of a group.
63  *
64  * Selects the indices of local atoms of a group and stores them in anrs_loc[0..nr_loc].
65  * If you need the positions of the group's atoms on all nodes, provide a coll_ind[0..nr]
66  * array and pass it on to communicate_group_positions. Thus the collective array
67  * will always have the same atom order (ascending indices).
68  *
69  *  \param[in]     ga2la      Global to local atom index conversion data.
70  *  \param[in]     nr         The total number of atoms that the group contains.
71  *  \param[in]     anrs       The global atom number of the group's atoms.
72  *  \param[out]    nr_loc     The number of group atoms present on the local node.
73  *  \param[out]    anrs_loc   The local atom numbers of the group.
74  *  \param[in,out] nalloc_loc Local allocation size of anrs_loc array.
75  *  \param[out]    coll_ind   If not NULL this array must be of size nr. It stores
76  *                            for each local atom where it belongs in the global
77  *                            (collective) array such that it can be gmx_summed
78  *                            in the communicate_group_positions routine.
79  */
80
81 extern void dd_make_local_group_indices(gmx_ga2la_t *ga2la,
82                                         const int nr, int anrs[], int *nr_loc,
83                                         int *anrs_loc[], int *nalloc_loc,
84                                         int coll_ind[]);
85
86
87 /*! \brief Assemble local positions into a collective array present on all nodes.
88  *
89  * Communicate the positions of the group's atoms such that every node has all of
90  * them. Unless running on huge number of cores, this is not a big performance impact
91  * as long as the collective subset [0..nr] is kept small. The atom indices are
92  * retrieved from anrs_loc[0..nr_loc]. If you call the routine for the serial case,
93  * provide an array coll_ind[i] = i for i in 1..nr.
94  *
95  * If shifts != NULL, the PBC representation of each atom is chosen such that a
96  * continuous trajectory results. Therefore, if the group is whole at the start
97  * of the simulation, it will always stay whole.
98  * If shifts = NULL, the group positions are not made whole again, but assembled
99  * and distributed to all nodes. The variables marked "optional" are not used in
100  * that case.
101  *
102  * \param[in]     cr           Pointer to MPI communication data.
103  * \param[out]    xcoll        Collective array of positions, identical on all nodes
104  *                             after this routine has been called.
105  * \param[in,out] shifts       Collective array of shifts for xcoll, needed to make
106  *                             the group whole. This array remembers the shifts
107  *                             since the start of the simulation (where the group
108  *                             is whole) and must therefore not be changed outside
109  *                             of this routine! If NULL, the group will not be made
110  *                             whole and the optional variables are ignored.
111  * \param[out]    extra_shifts Extra shifts since last time step, only needed as
112  *                             buffer variable [0..nr] (optional).
113  * \param[in]     bNS          Neighbor searching / domain re-decomposition has been
114  *                             performed at the begin of this time step such that
115  *                             the shifts have changed and need to be updated
116  *                             (optional).
117  * \param[in]     x_loc        Pointer to the local atom positions this node has.
118  * \param[in]     nr           Total number of atoms in the group.
119  * \param[in]     nr_loc       Number of group atoms on the local node.
120  * \param[in]     anrs_loc     Array of the local atom indices.
121  * \param[in]     coll_ind     This array of size nr stores for each local atom where
122  *                             it belongs in the collective array so that the local
123  *                             contributions can be gmx_summed. It is provided by
124  *                             dd_make_local_group_indices.
125  * \param[in,out] xcoll_old    Positions from the last time step, used to make the
126  *                             group whole (optional).
127  * \param[in]     box          Simulation box matrix, needed to shift xcoll such that
128  *                             the group becomes whole (optional).
129  */
130 extern void communicate_group_positions(const t_commrec *cr, rvec *xcoll, ivec *shifts,
131                                         ivec *extra_shifts, const gmx_bool bNS,
132                                         const rvec *x_loc, const int nr, const int nr_loc,
133                                         const int *anrs_loc, const int *coll_ind, rvec *xcoll_old,
134                                         matrix box);
135
136 /*! \brief Calculates the center of the positions x locally.
137  *
138  * Calculates the center of mass (if masses are given in the weight array) or
139  * the geometrical center (if NULL is passed as weight).
140  *
141  * \param[in]   x            Positions.
142  * \param[in]   weight       Can be NULL or an array of weights. If masses are
143  *                           given as weights, the COM is calculated.
144  * \param[in]   nr           Number of positions and weights if present.
145  * \param[out]  center       The (weighted) center of the positions.
146  *
147  */
148 extern void get_center(rvec x[], real weight[], const int nr, rvec center);
149
150
151 /*! \brief Calculates the sum of the positions x locally.
152  *
153  * Calculates the (weighted) sum of position vectors and returns the sum of
154  * weights, which is needed when local contributions shall be summed to a
155  * global weighted center.
156  *
157  * \param[in]   x            Array of positions.
158  * \param[in]   weight       Can be NULL or an array of weights.
159  * \param[in]   nr           Number of positions and weights if present.
160  * \param[out]  dsumvec      The (weighted) sum of the positions.
161  * \return Sum of weights.
162  *
163  */
164 extern double get_sum_of_positions(rvec x[], real weight[], const int nr, dvec dsumvec);
165
166
167 /*! \brief Calculates the global center of all local arrays x.
168  *
169  * Get the center from local positions [0..nr_loc], this involves communication.
170  * Not that the positions must already have the correct PBC representation. Use
171  * this routine if no collective coordinates are assembled from which the center
172  * could be calculated without communication.
173  *
174  * \param[in]   cr           Pointer to MPI communication data.
175  * \param[in]   x_loc        Array of local positions [0..nr_loc].
176  * \param[in]   weight_loc   Array of local weights, these are the masses if the
177  *                           center of mass is to be calculated.
178  * \param[in]   nr_loc       The number of positions on the local node.
179  * \param[in]   nr_group     The number of positions in the whole group. Since
180  *                           this is known anyway, we do not need to communicate
181  *                           and sum nr_loc if we pass it over.
182  * \param[out]  center       The (weighted) center of all x_loc from all the
183  *                           nodes.
184  */
185 extern void get_center_comm(const t_commrec *cr, rvec x_loc[], real weight_loc[],
186                             int nr_loc, int nr_group, rvec center);
187
188
189 /*! \brief Translate positions.
190  *
191  * Add a translation vector to the positions x.
192  *
193  * \param[in,out] x          Array of positions.
194  * \param[in]     nr         Number of entries in the position array.
195  * \param[in]     transvec   Translation vector to be added to all positions.
196  *
197  */
198 extern void translate_x(rvec x[], const int nr, const rvec transvec);
199
200
201 /*! \brief Rotate positions.
202  *
203  * Rotate the positions with the rotation matrix.
204  *
205  * \param[in,out] x          Array of positions.
206  * \param[in]     nr         Number of entries in the position array.
207  * \param[in]     rmat       Rotation matrix to operate on all positions.
208  *
209  */
210 extern void rotate_x(rvec x[], const int nr, matrix rmat);
211
212 #ifdef __cplusplus
213 }
214 #endif