Merge branch 'release-4-6'
[alexxy/gromacs.git] / src / gromacs / mdlib / groupcoord.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 4.5
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.
14
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.
19  *
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.
26  *
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.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Groningen Machine for Chemical Simulation
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39
40 #include "groupcoord.h"
41 #include "network.h"
42 #include "pbc.h"
43 #include "vec.h"
44 #include "smalloc.h"
45 #include "gmx_ga2la.h"
46
47 #define MIN(a, b) (((a) < (b)) ? (a) : (b))
48
49
50
51 /* Select the indices of the group's atoms which are local and store them in
52  * anrs_loc[0..nr_loc]. The indices are saved in coll_ind[] for later reduction
53  * in communicate_group_positions()
54  */
55 extern void dd_make_local_group_indices(
56         gmx_ga2la_t    ga2la,
57         const int      nr,         /* IN:  Total number of atoms in the group */
58         int            anrs[],     /* IN:  Global atom numbers of the groups atoms */
59         int           *nr_loc,     /* OUT: Number of group atoms found locally */
60         int           *anrs_loc[], /* OUT: Local atom numbers of the group  */
61         int           *nalloc_loc, /* IN+OUT: Allocation size of anrs_loc */
62         int            coll_ind[]) /* OUT (opt): Where is this position found in the collective array? */
63 {
64     int  i, ii;
65     int  localnr;
66
67
68     /* Loop over all the atom indices of the group to check
69      * which ones are on the local node */
70     localnr = 0;
71     for (i = 0; i < nr; i++)
72     {
73         if (ga2la_get_home(ga2la, anrs[i], &ii))
74         {
75             /* The atom with this index is a home atom */
76             if (localnr >= *nalloc_loc) /* Check whether memory suffices */
77             {
78                 *nalloc_loc = over_alloc_dd(localnr+1);
79                 /* We never need more memory than the number of atoms in the group */
80                 *nalloc_loc = MIN(*nalloc_loc, nr);
81                 srenew(*anrs_loc, *nalloc_loc);
82             }
83             /* Save the atoms index in the local atom numbers array */
84             (*anrs_loc)[localnr] = ii;
85
86             if (coll_ind != NULL)
87             {
88                 /* Keep track of where this local atom belongs in the collective index array.
89                  * This is needed when reducing the local arrays to a collective/global array
90                  * in communicate_group_positions */
91                 coll_ind[localnr] = i;
92             }
93
94             /* add one to the local atom count */
95             localnr++;
96         }
97     }
98
99     /* Return the number of local atoms that were found */
100     *nr_loc = localnr;
101 }
102
103
104 static void get_shifts_group(
105         int     npbcdim,
106         matrix  box,
107         rvec   *xcoll,     /* IN:  Collective set of positions [0..nr] */
108         int     nr,        /* IN:  Total number of atoms in the group */
109         rvec   *xcoll_old, /* IN:  Positions from the last time step [0...nr] */
110         ivec   *shifts)    /* OUT: Shifts for xcoll */
111 {
112     int  i, m, d;
113     rvec dx;
114
115
116     /* Get the shifts such that each atom is within closest
117      * distance to its position at the last NS time step after shifting.
118      * If we start with a whole group, and always keep track of
119      * shift changes, the group will stay whole this way */
120     for (i = 0; i < nr; i++)
121     {
122         clear_ivec(shifts[i]);
123     }
124
125     for (i = 0; i < nr; i++)
126     {
127         /* The distance this atom moved since the last time step */
128         /* If this is more than just a bit, it has changed its home pbc box */
129         rvec_sub(xcoll[i], xcoll_old[i], dx);
130
131         for (m = npbcdim-1; m >= 0; m--)
132         {
133             while (dx[m] < -0.5*box[m][m])
134             {
135                 for (d = 0; d < DIM; d++)
136                 {
137                     dx[d] += box[m][d];
138                 }
139                 shifts[i][m]++;
140             }
141             while (dx[m] >= 0.5*box[m][m])
142             {
143                 for (d = 0; d < DIM; d++)
144                 {
145                     dx[d] -= box[m][d];
146                 }
147                 shifts[i][m]--;
148             }
149         }
150     }
151 }
152
153
154 static void shift_positions_group(
155         matrix  box,
156         rvec    x[],     /* The positions [0..nr] */
157         ivec   *is,      /* The shifts [0..nr] */
158         int     nr)      /* The number of positions and shifts */
159 {
160     int      i, tx, ty, tz;
161
162
163     /* Loop over the group's atoms */
164     if (TRICLINIC(box))
165     {
166         for (i = 0; i < nr; i++)
167         {
168             tx = is[i][XX];
169             ty = is[i][YY];
170             tz = is[i][ZZ];
171
172             x[i][XX] = x[i][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
173             x[i][YY] = x[i][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
174             x[i][ZZ] = x[i][ZZ]+tz*box[ZZ][ZZ];
175         }
176     }
177     else
178     {
179         for (i = 0; i < nr; i++)
180         {
181             tx = is[i][XX];
182             ty = is[i][YY];
183             tz = is[i][ZZ];
184
185             x[i][XX] = x[i][XX]+tx*box[XX][XX];
186             x[i][YY] = x[i][YY]+ty*box[YY][YY];
187             x[i][ZZ] = x[i][ZZ]+tz*box[ZZ][ZZ];
188         }
189     }
190 }
191
192
193 /* Assemble the positions of the group such that every node has all of them.
194  * The atom indices are retrieved from anrs_loc[0..nr_loc]
195  * Note that coll_ind[i] = i is needed in the serial case */
196 extern void communicate_group_positions(
197         t_commrec     *cr,
198         rvec          *xcoll,        /* OUT: Collective array of positions */
199         ivec          *shifts,       /* IN+OUT: Collective array of shifts for xcoll */
200         ivec          *extra_shifts, /* BUF: Extra shifts since last time step */
201         const gmx_bool bNS,          /* IN:  NS step, the shifts have changed */
202         rvec          *x_loc,        /* IN:  Local positions on this node */
203         const int      nr,           /* IN:  Total number of atoms in the group */
204         const int      nr_loc,       /* IN:  Local number of atoms in the group */
205         int           *anrs_loc,     /* IN:  Local atom numbers */
206         int           *coll_ind,     /* IN:  Collective index */
207         rvec          *xcoll_old,    /* IN+OUT: Positions from the last time step, used to make group whole */
208         matrix         box)
209 {
210     int i;
211
212
213     /* Zero out the groups' global position array */
214     clear_rvecs(nr, xcoll);
215
216     /* Put the local positions that this node has into the right place of
217      * the collective array. Note that in the serial case, coll_ind[i] = i */
218     for (i = 0; i < nr_loc; i++)
219     {
220         copy_rvec(x_loc[anrs_loc[i]], xcoll[coll_ind[i]]);
221     }
222
223     if (PAR(cr))
224     {
225         /* Add the arrays from all nodes together */
226         gmx_sum(nr*3, xcoll[0], cr);
227     }
228     /* To make the group whole, start with a whole group and each
229      * step move the assembled positions at closest distance to the positions
230      * from the last step. First shift the positions with the saved shift
231      * vectors (these are 0 when this routine is called for the first time!) */
232     shift_positions_group(box, xcoll, shifts, nr);
233
234     /* Now check if some shifts changed since the last step.
235      * This only needs to be done when the shifts are expected to have changed,
236      * i.e. after neighboursearching */
237     if (bNS)
238     {
239         get_shifts_group(3, box, xcoll, nr, xcoll_old, extra_shifts);
240
241         /* Shift with the additional shifts such that we get a whole group now */
242         shift_positions_group(box, xcoll, extra_shifts, nr);
243
244         /* Add the shift vectors together for the next time step */
245         for (i = 0; i < nr; i++)
246         {
247             shifts[i][XX] += extra_shifts[i][XX];
248             shifts[i][YY] += extra_shifts[i][YY];
249             shifts[i][ZZ] += extra_shifts[i][ZZ];
250         }
251
252         /* Store current correctly-shifted positions for comparison in the next NS time step */
253         for (i = 0; i < nr; i++)
254         {
255             copy_rvec(xcoll[i], xcoll_old[i]);
256         }
257     }
258 }
259
260
261 /* Determine the (weighted) sum vector from positions x */
262 extern double get_sum_of_positions(rvec x[], real weight[], const int nat, dvec dsumvec)
263 {
264     int    i;
265     rvec   x_weighted;
266     double weight_sum = 0.0;
267
268
269     /* Zero out the center */
270     clear_dvec(dsumvec);
271
272     /* Loop over all atoms and add their weighted position vectors */
273     if (weight != NULL)
274     {
275         for (i = 0; i < nat; i++)
276         {
277             weight_sum += weight[i];
278             svmul(weight[i], x[i], x_weighted);
279             dsumvec[XX] += x_weighted[XX];
280             dsumvec[YY] += x_weighted[YY];
281             dsumvec[ZZ] += x_weighted[ZZ];
282         }
283     }
284     else
285     {
286         for (i = 0; i < nat; i++)
287         {
288             dsumvec[XX] += x[i][XX];
289             dsumvec[YY] += x[i][YY];
290             dsumvec[ZZ] += x[i][ZZ];
291         }
292     }
293     return weight_sum;
294 }
295
296
297 /* Determine center of structure from collective positions x */
298 extern void get_center(rvec x[], real weight[], const int nr, rvec rcenter)
299 {
300     dvec   dcenter;
301     double weight_sum, denom;
302
303
304     weight_sum = get_sum_of_positions(x, weight, nr, dcenter);
305
306     if (weight != NULL)
307     {
308         denom = weight_sum; /* Divide by the sum of weight */
309     }
310     else
311     {
312         denom = nr;        /* Divide by the number of atoms */
313
314     }
315     dsvmul(1.0/denom, dcenter, dcenter);
316
317     rcenter[XX] = dcenter[XX];
318     rcenter[YY] = dcenter[YY];
319     rcenter[ZZ] = dcenter[ZZ];
320 }
321
322
323 /* Get the center from local positions that already have the correct
324  * PBC representation */
325 extern void get_center_comm(
326         t_commrec *cr,
327         rvec       x_loc[],      /* Local positions */
328         real       weight_loc[], /* Local masses or other weights */
329         int        nr_loc,       /* Local number of atoms */
330         int        nr_group,     /* Total number of atoms of the group */
331         rvec       center)       /* Weighted center */
332 {
333     double weight_sum, denom;
334     dvec   dsumvec;
335     double buf[4];
336
337
338     weight_sum = get_sum_of_positions(x_loc, weight_loc, nr_loc, dsumvec);
339
340     /* Add the local contributions from all nodes. Put the sum vector and the
341      * weight in a buffer array so that we get along with a single communication
342      * call. */
343     if (PAR(cr))
344     {
345         buf[0] = dsumvec[XX];
346         buf[1] = dsumvec[YY];
347         buf[2] = dsumvec[ZZ];
348         buf[3] = weight_sum;
349
350         /* Communicate buffer */
351         gmx_sumd(4, buf, cr);
352
353         dsumvec[XX] = buf[0];
354         dsumvec[YY] = buf[1];
355         dsumvec[ZZ] = buf[2];
356         weight_sum  = buf[3];
357     }
358
359     if (weight_loc != NULL)
360     {
361         denom = 1.0/weight_sum; /* Divide by the sum of weight to get center of mass e.g. */
362     }
363     else
364     {
365         denom = 1.0/nr_group;   /* Divide by the number of atoms to get the geometrical center */
366
367     }
368     center[XX] = dsumvec[XX]*denom;
369     center[YY] = dsumvec[YY]*denom;
370     center[ZZ] = dsumvec[ZZ]*denom;
371 }
372
373
374 /* Translate x with transvec */
375 extern void translate_x(rvec x[], const int nr, const rvec transvec)
376 {
377     int i;
378
379
380     for (i = 0; i < nr; i++)
381     {
382         rvec_inc(x[i], transvec);
383     }
384 }
385
386
387 extern void rotate_x(rvec x[], const int nr, matrix rmat)
388 {
389     int  i, j, k;
390     rvec x_old;
391
392
393     /* Apply the rotation matrix */
394     for (i = 0; i < nr; i++)
395     {
396         for (j = 0; j < 3; j++)
397         {
398             x_old[j] = x[i][j];
399         }
400         for (j = 0; j < 3; j++)
401         {
402             x[i][j] = 0;
403             for (k = 0; k < 3; k++)
404             {
405                 x[i][j] += rmat[k][j]*x_old[k];
406             }
407         }
408     }
409 }