Fixed include statements such that double precision version of genborn.c
[gromacs/rigid-bodies.git] / src / mdlib / groupcoord.c
blob1449e01c8fded3ffae067985fa159dde7a9dc9e0
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 #include "groupcoord.h"
37 #include "mpelogging.h"
38 #include "network.h"
39 #include "pbc.h"
40 #include "vec.h"
41 #include "smalloc.h"
42 #include "gmx_ga2la.h"
44 #define MIN(a,b) (((a)<(b))?(a):(b))
48 /* Select the indices of the group's atoms which are local and store them in
49 * anrs_loc[0..nr_loc]. The indices are saved in coll_ind[] for later reduction
50 * in communicate_group_positions()
52 extern void dd_make_local_group_indices(
53 gmx_ga2la_t ga2la,
54 const int nr, /* IN: Total number of atoms in the group */
55 int anrs[], /* IN: Global atom numbers of the groups atoms */
56 int *nr_loc, /* OUT: Number of group atoms found locally */
57 int *anrs_loc[], /* OUT: Local atom numbers of the group */
58 int *nalloc_loc, /* IN+OUT: Allocation size of anrs_loc */
59 int coll_ind[]) /* OUT (opt): Where is this position found in the collective array? */
61 int i,ii;
62 int localnr;
65 /* Loop over all the atom indices of the group to check
66 * which ones are on the local node */
67 localnr = 0;
68 for(i=0; i<nr; i++)
70 if (ga2la_get_home(ga2la,anrs[i],&ii))
72 /* The atom with this index is a home atom */
73 if (localnr >= *nalloc_loc) /* Check whether memory suffices */
75 *nalloc_loc = over_alloc_dd(localnr+1);
76 /* We never need more memory than the number of atoms in the group */
77 *nalloc_loc = MIN(*nalloc_loc, nr);
78 srenew(*anrs_loc,*nalloc_loc);
80 /* Save the atoms index in the local atom numbers array */
81 (*anrs_loc)[localnr] = ii;
83 if (coll_ind != NULL)
85 /* Keep track of where this local atom belongs in the collective index array.
86 * This is needed when reducing the local arrays to a collective/global array
87 * in communicate_group_positions */
88 coll_ind[localnr] = i;
91 /* add one to the local atom count */
92 localnr++;
96 /* Return the number of local atoms that were found */
97 *nr_loc = localnr;
101 static void get_shifts_group(
102 int npbcdim,
103 matrix box,
104 rvec *xcoll, /* IN: Collective set of positions [0..nr] */
105 int nr, /* IN: Total number of atoms in the group */
106 rvec *xcoll_old, /* IN: Positions from the last time step [0...nr] */
107 ivec *shifts) /* OUT: Shifts for xcoll */
109 int i,m,d;
110 rvec dx;
113 /* Get the shifts such that each atom is within closest
114 * distance to its position at the last NS time step after shifting.
115 * If we start with a whole group, and always keep track of
116 * shift changes, the group will stay whole this way */
117 for (i=0; i < nr; i++)
118 clear_ivec(shifts[i]);
120 for (i=0; i<nr; i++)
122 /* The distance this atom moved since the last time step */
123 /* If this is more than just a bit, it has changed its home pbc box */
124 rvec_sub(xcoll[i],xcoll_old[i],dx);
126 for(m=npbcdim-1; m>=0; m--)
128 while (dx[m] < -0.5*box[m][m])
130 for(d=0; d<DIM; d++)
131 dx[d] += box[m][d];
132 shifts[i][m]++;
134 while (dx[m] >= 0.5*box[m][m])
136 for(d=0; d<DIM; d++)
137 dx[d] -= box[m][d];
138 shifts[i][m]--;
145 static void shift_positions_group(
146 matrix box,
147 rvec x[], /* The positions [0..nr] */
148 ivec *is, /* The shifts [0..nr] */
149 int nr) /* The number of positions and shifts */
151 int i,tx,ty,tz;
154 GMX_MPE_LOG(ev_shift_start);
156 /* Loop over the group's atoms */
157 if(TRICLINIC(box))
159 for (i=0; i < nr; i++)
161 tx=is[i][XX];
162 ty=is[i][YY];
163 tz=is[i][ZZ];
165 x[i][XX]=x[i][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
166 x[i][YY]=x[i][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
167 x[i][ZZ]=x[i][ZZ]+tz*box[ZZ][ZZ];
169 } else
171 for (i=0; i < nr; i++)
173 tx=is[i][XX];
174 ty=is[i][YY];
175 tz=is[i][ZZ];
177 x[i][XX]=x[i][XX]+tx*box[XX][XX];
178 x[i][YY]=x[i][YY]+ty*box[YY][YY];
179 x[i][ZZ]=x[i][ZZ]+tz*box[ZZ][ZZ];
182 GMX_MPE_LOG(ev_shift_finish);
186 /* Assemble the positions of the group such that every node has all of them.
187 * The atom indices are retrieved from anrs_loc[0..nr_loc]
188 * Note that coll_ind[i] = i is needed in the serial case */
189 extern void communicate_group_positions(
190 t_commrec *cr,
191 rvec *xcoll, /* OUT: Collective array of positions */
192 ivec *shifts, /* IN+OUT: Collective array of shifts for xcoll */
193 ivec *extra_shifts, /* BUF: Extra shifts since last time step */
194 const bool bNS, /* IN: NS step, the shifts have changed */
195 rvec *x_loc, /* IN: Local positions on this node */
196 const int nr, /* IN: Total number of atoms in the group */
197 const int nr_loc, /* IN: Local number of atoms in the group */
198 int *anrs_loc, /* IN: Local atom numbers */
199 int *coll_ind, /* IN: Collective index */
200 rvec *xcoll_old, /* IN+OUT: Positions from the last time step, used to make group whole */
201 matrix box)
203 int i;
206 GMX_MPE_LOG(ev_get_group_x_start);
208 /* Zero out the groups' global position array */
209 clear_rvecs(nr, xcoll);
211 /* Put the local positions that this node has into the right place of
212 * the collective array. Note that in the serial case, coll_ind[i] = i */
213 for (i=0; i<nr_loc; i++)
214 copy_rvec(x_loc[anrs_loc[i]], xcoll[coll_ind[i]]);
216 if (PAR(cr))
218 /* Add the arrays from all nodes together */
219 gmx_sum(nr*3, xcoll[0], cr);
221 /* To make the group whole, start with a whole group and each
222 * step move the assembled positions at closest distance to the positions
223 * from the last step. First shift the positions with the saved shift
224 * vectors (these are 0 when this routine is called for the first time!) */
225 shift_positions_group(box, xcoll, shifts, nr);
227 /* Now check if some shifts changed since the last step.
228 * This only needs to be done when the shifts are expected to have changed,
229 * i.e. after neighboursearching */
230 if (bNS)
232 get_shifts_group(3, box, xcoll, nr, xcoll_old, extra_shifts);
234 /* Shift with the additional shifts such that we get a whole group now */
235 shift_positions_group(box, xcoll, extra_shifts, nr);
237 /* Add the shift vectors together for the next time step */
238 for (i=0; i<nr; i++)
240 shifts[i][XX] += extra_shifts[i][XX];
241 shifts[i][YY] += extra_shifts[i][YY];
242 shifts[i][ZZ] += extra_shifts[i][ZZ];
245 /* Store current correctly-shifted positions for comparison in the next NS time step */
246 for (i=0; i<nr; i++)
247 copy_rvec(xcoll[i],xcoll_old[i]);
251 GMX_MPE_LOG(ev_get_group_x_finish);
255 /* Determine the (weighted) sum vector from positions x */
256 extern double get_sum_of_positions(rvec x[], real weight[], const int nat, dvec dsumvec)
258 int i;
259 rvec x_weighted;
260 double weight_sum = 0.0;
263 /* Zero out the center */
264 clear_dvec(dsumvec);
266 /* Loop over all atoms and add their weighted position vectors */
267 if (weight != NULL)
269 for (i=0; i<nat; i++)
271 weight_sum += weight[i];
272 svmul(weight[i], x[i], x_weighted);
273 dsumvec[XX] += x_weighted[XX];
274 dsumvec[YY] += x_weighted[YY];
275 dsumvec[ZZ] += x_weighted[ZZ];
278 else
280 for (i=0; i<nat; i++)
282 dsumvec[XX] += x[i][XX];
283 dsumvec[YY] += x[i][YY];
284 dsumvec[ZZ] += x[i][ZZ];
287 return weight_sum;
291 /* Determine center of structure from collective positions x */
292 extern void get_center(rvec x[], real weight[], const int nr, rvec rcenter)
294 dvec dcenter;
295 double weight_sum, denom;
298 weight_sum = get_sum_of_positions(x, weight, nr, dcenter);
300 if (weight != NULL)
301 denom = weight_sum; /* Divide by the sum of weight */
302 else
303 denom = nr; /* Divide by the number of atoms */
305 dsvmul(1.0/denom, dcenter, dcenter);
307 rcenter[XX] = dcenter[XX];
308 rcenter[YY] = dcenter[YY];
309 rcenter[ZZ] = dcenter[ZZ];
313 /* Get the center from local positions that already have the correct
314 * PBC representation */
315 extern void get_center_comm(
316 t_commrec *cr,
317 rvec x_loc[], /* Local positions */
318 real weight_loc[], /* Local masses or other weights */
319 int nr_loc, /* Local number of atoms */
320 int nr_group, /* Total number of atoms of the group */
321 rvec center) /* Weighted center */
323 double weight_sum, denom;
324 dvec dsumvec;
325 double buf[4];
328 weight_sum = get_sum_of_positions(x_loc, weight_loc, nr_loc, dsumvec);
330 /* Add the local contributions from all nodes. Put the sum vector and the
331 * weight in a buffer array so that we get along with a single communication
332 * call. */
333 if (PAR(cr))
335 buf[0] = dsumvec[XX];
336 buf[1] = dsumvec[YY];
337 buf[2] = dsumvec[ZZ];
338 buf[3] = weight_sum;
340 /* Communicate buffer */
341 gmx_sumd(4, buf, cr);
343 dsumvec[XX] = buf[0];
344 dsumvec[YY] = buf[1];
345 dsumvec[ZZ] = buf[2];
346 weight_sum = buf[3];
349 if (weight_loc != NULL)
350 denom = 1.0/weight_sum; /* Divide by the sum of weight to get center of mass e.g. */
351 else
352 denom = 1.0/nr_group; /* Divide by the number of atoms to get the geometrical center */
354 center[XX] = dsumvec[XX]*denom;
355 center[YY] = dsumvec[YY]*denom;
356 center[ZZ] = dsumvec[ZZ]*denom;
360 /* Translate x with transvec */
361 extern void translate_x(rvec x[], const int nr, const rvec transvec)
363 int i;
366 for (i=0; i<nr; i++)
367 rvec_inc(x[i], transvec);
371 extern void rotate_x(rvec x[], const int nr, matrix rmat)
373 int i,j,k;
374 rvec x_old;
377 /* Apply the rotation matrix */
378 for (i=0; i<nr; i++)
380 for (j=0; j<3; j++)
381 x_old[j] = x[i][j];
382 for (j=0; j<3; j++)
384 x[i][j] = 0;
385 for (k=0; k<3; k++)
386 x[i][j] += rmat[k][j]*x_old[k];