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