2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2008, The GROMACS development team.
6 * Copyright (c) 2012,2014, 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.
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.
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.
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.
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.
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.
37 /*! \libinternal \file groupcoord.h
39 * Assemble atom positions for comparison with a reference set.
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.
56 #include "types/commrec.h"
59 /*! \brief Select local atoms of a group.
61 * Selects the indices of local atoms of a group and stores them in anrs_loc[0..nr_loc].
62 * If you need the positions of the group's atoms on all nodes, provide a coll_ind[0..nr]
63 * array and pass it on to communicate_group_positions. Thus the collective array
64 * will always have the same atom order (ascending indices).
66 * \param[in] ga2la Global to local atom index conversion data.
67 * \param[in] nr The total number of atoms that the group contains.
68 * \param[in] anrs The global atom number of the group's atoms.
69 * \param[out] nr_loc The number of group atoms present on the local node.
70 * \param[out] anrs_loc The local atom numbers of the group.
71 * \param[in,out] nalloc_loc Local allocation size of anrs_loc array.
72 * \param[out] coll_ind If not NULL this array must be of size nr. It stores
73 * for each local atom where it belongs in the global
74 * (collective) array such that it can be gmx_summed
75 * in the communicate_group_positions routine.
77 extern void dd_make_local_group_indices(gmx_ga2la_t ga2la,
78 const int nr, int anrs[], int *nr_loc,
79 int *anrs_loc[], int *nalloc_loc,
83 /*! \brief Assemble local positions into a collective array present on all nodes.
85 * Communicate the positions of the group's atoms such that every node has all of
86 * them. Unless running on huge number of cores, this is not a big performance impact
87 * as long as the collective subset [0..nr] is kept small. The atom indices are
88 * retrieved from anrs_loc[0..nr_loc]. If you call the routine for the serial case,
89 * provide an array coll_ind[i] = i for i in 1..nr.
91 * \param[in] cr Pointer to MPI communication data.
92 * \param[out] xcoll Collective array of positions, idential on all nodes
93 * after this routine has been called.
94 * \param[in,out] shifts Collective array of shifts for xcoll, needed to make
95 * the group whole. This array remembers the shifts
96 * since the start of the simulation (where the group
97 * is whole) and must therefore not be changed outside
99 * \param[out] extra_shifts Extra shifts since last time step, only needed as
100 * buffer variable [0..nr].
101 * \param[in] bNS Neighborsearching/domain redecomposition has been
102 * performed at the begin of this time step such that
103 * the shifts have changed and need to be updated.
104 * \param[in] x_loc Pointer to the local atom positions this node has.
105 * \param[in] nr Total number of atoms in the group.
106 * \param[in] nr_loc Number of group atoms on the local node.
107 * \param[in] anrs_loc Array of the local atom indices.
108 * \param[in] coll_ind This array of size nr stores for each local atom where
109 * it belongs in the collective array so that the local
110 * contributions can be gmx_summed. It is provided by
111 * dd_make_local_group_indices.
112 * \param[in,out] xcoll_old Positions from the last time step, used to make the
114 * \param[in] box Simulation box matrix, needed to shift xcoll such that
115 * the group becomes whole.
117 extern void communicate_group_positions(t_commrec *cr, rvec *xcoll, ivec *shifts,
118 ivec *extra_shifts, const gmx_bool bNS,
119 rvec *x_loc, const int nr, const int nr_loc,
120 int *anrs_loc, int *coll_ind, rvec *xcoll_old,
124 /*! \brief Calculates the center of the positions x locally.
126 * Calculates the center of mass (if masses are given in the weight array) or
127 * the geometrical center (if NULL is passed as weight).
129 * \param[in] x Positions.
130 * \param[in] weight Can be NULL or an array of weights. If masses are
131 * given as weights, the COM is calculated.
132 * \param[in] nr Number of positions and weights if present.
133 * \param[out] center The (weighted) center of the positions.
136 extern void get_center(rvec x[], real weight[], const int nr, rvec center);
139 /*! \brief Calculates the sum of the positions x locally.
141 * Calculates the (weighted) sum of position vectors and returns the sum of
142 * weights, which is needed when local contributions shall be summed to a
143 * global weighted center.
145 * \param[in] x Array of positions.
146 * \param[in] weight Can be NULL or an array of weights.
147 * \param[in] nr Number of positions and weights if present.
148 * \param[out] dsumvec The (weighted) sum of the positions.
149 * \return Sum of weights.
152 extern double get_sum_of_positions(rvec x[], real weight[], const int nr, dvec dsumvec);
155 /*! \brief Calculates the global center of all local arrays x.
157 * Get the center from local positions [0..nr_loc], this involves communication.
158 * Not that the positions must already have the correct PBC representation. Use
159 * this routine if no collective coordinates are assembled from which the center
160 * could be calculated without communication.
162 * \param[in] cr Pointer to MPI communication data.
163 * \param[in] x_loc Array of local positions [0..nr_loc].
164 * \param[in] weight_loc Array of local weights, these are the masses if the
165 * center of mass is to be calculated.
166 * \param[in] nr_loc The number of positions on the local node.
167 * \param[in] nr_group The number of positions in the whole group. Since
168 * this is known anyway, we do not need to communicate
169 * and sum nr_loc if we pass it over.
170 * \param[out] center The (weighted) center of all x_loc from all the
173 extern void get_center_comm(t_commrec *cr, rvec x_loc[], real weight_loc[],
174 int nr_loc, int nr_group, rvec center);
177 /*! \brief Translate positions.
179 * Add a translation vector to the positions x.
181 * \param[in,out] x Array of positions.
182 * \param[in] nr Number of entries in the position array.
183 * \param[in] transvec Translation vector to be added to all positions.
186 extern void translate_x(rvec x[], const int nr, const rvec transvec);
189 /*! \brief Rotate positions.
191 * Rotate the positions with the rotation matrix.
193 * \param[in,out] x Array of positions.
194 * \param[in] nr Number of entries in the position array.
195 * \param[in] rmat Rotation matrix to operate on all positions.
198 extern void rotate_x(rvec x[], const int nr, matrix rmat);