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