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