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