A freeze group can now be allowed to move rigidly in some dimension by using "freezed...
[gromacs/rigid-bodies.git] / src / mdlib / groupcoord.h
blob79d107dcc358ebf52aaa96d55b201b1d50044d29
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2008, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
36 /*! \file groupcoord.h
38 * @brief Assemble atom positions for comparison with a reference set.
40 * This file contains functions to assemble the positions of a subset of the
41 * atoms and to do operations on it like determining the center of mass, or
42 * doing translations and rotations. These functions are useful when
43 * a subset of the positions needs to be compared to some set of reference
44 * positions, as e.g. done for essential dynamics.
48 #ifdef HAVE_CONFIG_H
49 #include <config.h>
50 #endif
52 #include <stdio.h>
53 #include "typedefs.h"
56 /*! \brief Select local atoms of a group.
58 * Selects the indices of local atoms of a group and stores them in anrs_loc[0..nr_loc].
59 * If you need the positions of the group's atoms on all nodes, provide a coll_ind[0..nr]
60 * array and pass it on to communicate_group_positions. Thus the collective array
61 * will always have the same atom order (ascending indices).
63 * \param ga2la[in] Global to local atom index conversion data.
64 * \param nr[in] The total number of atoms that the group contains.
65 * \param anrs[in] The global atom number of the group's atoms.
66 * \param nr_loc[out] The number of group atoms present on the local node.
67 * \param anrs_loc[out] The local atom numbers of the group.
68 * \param nalloc_loc[inout] Local allocation size of anrs_loc array.
69 * \param coll_ind[opt] If not NULL this array must be of size nr. It stores
70 * for each local atom where it belongs in the global
71 * (collective) array such that it can be gmx_summed
72 * in the communicate_group_positions routine.
74 extern void dd_make_local_group_indices(gmx_ga2la_t ga2la,
75 const int nr, int anrs[], int *nr_loc,
76 int *anrs_loc[], int *nalloc_loc,
77 int coll_ind[]);
80 /*! \brief Assemble local positions into a collective array present on all nodes.
82 * Communicate the positions of the group's atoms such that every node has all of
83 * them. Unless running on huge number of cores, this is not a big performance impact
84 * as long as the collective subset [0..nr] is kept small. The atom indices are
85 * retrieved from anrs_loc[0..nr_loc]. If you call the routine for the serial case,
86 * provide an array coll_ind[i] = i for i in 1..nr.
88 * \param cr[in] Pointer to MPI communication data.
89 * \param xcoll[out] Collective array of positions, idential on all nodes
90 * after this routine has been called.
91 * \param shifts[inout] Collective array of shifts for xcoll, needed to make
92 * the group whole. This array remembers the shifts
93 * since the start of the simulation (where the group
94 * is whole) and must therefore not be changed outside
95 * of this routine!
96 * \param extra_shifts[buf] Extra shifts since last time step, only needed as
97 * buffer variable [0..nr].
98 * \param bNS[in] Neighborsearching/domain redecomposition has been
99 * performed at the begin of this time step such that
100 * the shifts have changed and need to be updated.
101 * \param x_loc[in] Pointer to the local atom positions this node has.
102 * \param nr[in] Total number of atoms in the group.
103 * \param nr_loc[in] Number of group atoms on the local node.
104 * \param anrs_loc[in] Array of the local atom indices.
105 * \param coll_ind[in] This array of size nr stores for each local atom where
106 * it belongs in the collective array so that the local
107 * contributions can be gmx_summed. It is provided by
108 * dd_make_local_group_indices.
109 * \param xcoll_old[inout] Positions from the last time step, used to make the
110 * group whole.
111 * \param box[in] Simulation box matrix, needed to shift xcoll such that
112 * the group becomes whole.
114 extern void communicate_group_positions(t_commrec *cr, rvec *xcoll, ivec *shifts,
115 ivec *extra_shifts, const gmx_bool bNS,
116 rvec *x_loc, const int nr, const int nr_loc,
117 int *anrs_loc, int *coll_ind, rvec *xcoll_old,
118 matrix box);
121 /*! \brief Calculates the center of the positions x locally.
123 * Calculates the center of mass (if masses are given in the weight array) or
124 * the geometrical center (if NULL is passed as weight).
126 * \param x[in] Positions.
127 * \param weight[in] Can be NULL or an array of weights. If masses are
128 * given as weights, the COM is calculated.
129 * \param nr[in] Number of positions and weights if present.
130 * \param center[out] The (weighted) center of the positions.
133 extern void get_center(rvec x[], real weight[], const int nr, rvec center);
136 /*! \brief Calculates the sum of the positions x locally.
138 * Calculates the (weighted) sum of position vectors and returns the sum of
139 * weights, which is needed when local contributions shall be summed to a
140 * global weighted center.
142 * \param x[in] Array of positions.
143 * \param weight[in] Can be NULL or an array of weights.
144 * \param nr[in] Number of positions and weights if present.
145 * \param dsumvec[out] The (weighted) sum of the positions.
146 * \return Sum of weights.
149 extern double get_sum_of_positions(rvec x[], real weight[], const int nr, dvec dsumvec);
152 /*! \brief Calculates the global center of all local arrays x.
154 * Get the center from local positions [0..nr_loc], this involves communication.
155 * Not that the positions must already have the correct PBC representation. Use
156 * this routine if no collective coordinates are assembled from which the center
157 * could be calculated without communication.
159 * \param cr[in] Pointer to MPI communication data.
160 * \param x_loc[in] Array of local positions [0..nr_loc].
161 * \param weight_loc[in] Array of local weights, these are the masses if the
162 * center of mass is to be calculated.
163 * \param nr_loc[in] The number of positions on the local node.
164 * \param nr_group[in] The number of positions in the whole group. Since
165 * this is known anyway, we do not need to communicate
166 * and sum nr_loc if we pass it over.
167 * \param center[out] The (weighted) center of all x_loc from all the
168 * nodes.
170 extern void get_center_comm(t_commrec *cr, rvec x_loc[], real weight_loc[],
171 int nr_loc, int nr_group, rvec center);
174 /*! \brief Translate positions.
176 * Add a translation vector to the positions x.
178 * \param x[inout] Array of positions.
179 * \param nr[in] Number of entries in the position array.
180 * \param transvec[in] Translation vector to be added to all positions.
183 extern void translate_x(rvec x[], const int nr, const rvec transvec);
186 /*! \brief Rotate positions.
188 * Rotate the positions with the rotation matrix.
190 * \param x[inout] Array of positions.
191 * \param nr[in] Number of entries in the position array.
192 * \param rmat[in] Rotation matrix to operate on all positions.
195 extern void rotate_x(rvec x[], const int nr, matrix rmat);