Fixing copyright issues and code contributors
[alexxy/gromacs.git] / src / 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  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012,2013, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42
43 #include "groupcoord.h"
44 #include "mpelogging.h"
45 #include "network.h"
46 #include "pbc.h"
47 #include "vec.h"
48 #include "smalloc.h"
49 #include "gmx_ga2la.h"
50
51 #define MIN(a,b) (((a)<(b))?(a):(b))
52
53
54
55 /* Select the indices of the group's atoms which are local and store them in 
56  * anrs_loc[0..nr_loc]. The indices are saved in coll_ind[] for later reduction
57  * in communicate_group_positions()
58  */
59 extern void dd_make_local_group_indices(
60         gmx_ga2la_t   ga2la,
61         const int     nr,          /* IN:  Total number of atoms in the group */
62         int           anrs[],      /* IN:  Global atom numbers of the groups atoms */
63         int           *nr_loc,     /* OUT: Number of group atoms found locally */
64         int           *anrs_loc[], /* OUT: Local atom numbers of the group  */
65         int           *nalloc_loc, /* IN+OUT: Allocation size of anrs_loc */
66         int           coll_ind[])  /* OUT (opt): Where is this position found in the collective array? */
67 {
68     int  i,ii;
69     int  localnr;       
70
71     
72     /* Loop over all the atom indices of the group to check
73      * which ones are on the local node */
74     localnr = 0;
75     for(i=0; i<nr; i++)
76     {
77         if (ga2la_get_home(ga2la,anrs[i],&ii))
78         {
79             /* The atom with this index is a home atom */
80             if (localnr >= *nalloc_loc) /* Check whether memory suffices */
81             {
82                 *nalloc_loc = over_alloc_dd(localnr+1);
83                 /* We never need more memory than the number of atoms in the group */
84                 *nalloc_loc = MIN(*nalloc_loc, nr);
85                 srenew(*anrs_loc,*nalloc_loc);
86             }
87             /* Save the atoms index in the local atom numbers array */
88             (*anrs_loc)[localnr] = ii;
89
90             if (coll_ind != NULL)
91             {
92                 /* Keep track of where this local atom belongs in the collective index array.
93                  * This is needed when reducing the local arrays to a collective/global array
94                  * in communicate_group_positions */
95                 coll_ind[localnr] = i;
96             }
97
98             /* add one to the local atom count */
99             localnr++;
100         }
101     }
102  
103     /* Return the number of local atoms that were found */
104     *nr_loc = localnr;
105 }
106
107
108 static void get_shifts_group(
109         int    npbcdim, 
110         matrix box,
111         rvec   *xcoll,     /* IN:  Collective set of positions [0..nr] */
112         int    nr,         /* IN:  Total number of atoms in the group */
113         rvec   *xcoll_old, /* IN:  Positions from the last time step [0...nr] */
114         ivec   *shifts)    /* OUT: Shifts for xcoll */
115 {
116     int  i,m,d;
117     rvec dx;
118
119
120     /* Get the shifts such that each atom is within closest
121      * distance to its position at the last NS time step after shifting.
122      * If we start with a whole group, and always keep track of 
123      * shift changes, the group will stay whole this way */
124     for (i=0; i < nr; i++)
125         clear_ivec(shifts[i]);
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                     dx[d] += box[m][d];
139                 shifts[i][m]++;
140             }
141             while (dx[m] >= 0.5*box[m][m])
142             {
143                 for(d=0; d<DIM; d++)
144                     dx[d] -= box[m][d];
145                 shifts[i][m]--;
146             }
147         }
148     }
149 }
150
151
152 static void shift_positions_group(
153         matrix box, 
154         rvec   x[],      /* The positions [0..nr] */ 
155         ivec   *is,      /* The shifts [0..nr] */ 
156         int    nr)       /* The number of positions and shifts */
157 {
158     int      i,tx,ty,tz;
159
160
161     GMX_MPE_LOG(ev_shift_start);
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     } 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     GMX_MPE_LOG(ev_shift_finish);
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     GMX_MPE_LOG(ev_get_group_x_start);
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         copy_rvec(x_loc[anrs_loc[i]], xcoll[coll_ind[i]]);
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             copy_rvec(xcoll[i],xcoll_old[i]);
255     }
256     
257     GMX_MPE_LOG(ev_get_group_x_finish);
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         denom = weight_sum; /* Divide by the sum of weight */
308     else
309         denom = nr;        /* Divide by the number of atoms */
310         
311     dsvmul(1.0/denom, dcenter, dcenter);
312     
313     rcenter[XX] = dcenter[XX];
314     rcenter[YY] = dcenter[YY];
315     rcenter[ZZ] = dcenter[ZZ];
316 }
317
318
319 /* Get the center from local positions that already have the correct
320  * PBC representation */
321 extern void get_center_comm(
322         t_commrec *cr,
323         rvec x_loc[],       /* Local positions */
324         real weight_loc[],  /* Local masses or other weights */
325         int nr_loc,         /* Local number of atoms */
326         int nr_group,       /* Total number of atoms of the group */ 
327         rvec center)        /* Weighted center */
328 {
329     double weight_sum, denom;
330     dvec   dsumvec;
331     double buf[4];    
332     
333     
334     weight_sum = get_sum_of_positions(x_loc, weight_loc, nr_loc, dsumvec);
335     
336     /* Add the local contributions from all nodes. Put the sum vector and the 
337      * weight in a buffer array so that we get along with a single communication
338      * call. */
339     if (PAR(cr))
340     {
341         buf[0] = dsumvec[XX];
342         buf[1] = dsumvec[YY];
343         buf[2] = dsumvec[ZZ];
344         buf[3] = weight_sum;
345         
346         /* Communicate buffer */
347         gmx_sumd(4, buf, cr);
348         
349         dsumvec[XX] = buf[0];
350         dsumvec[YY] = buf[1];
351         dsumvec[ZZ] = buf[2];
352         weight_sum  = buf[3];
353     }
354     
355     if (weight_loc != NULL)
356         denom = 1.0/weight_sum; /* Divide by the sum of weight to get center of mass e.g. */
357     else
358         denom = 1.0/nr_group;   /* Divide by the number of atoms to get the geometrical center */
359         
360     center[XX] = dsumvec[XX]*denom;
361     center[YY] = dsumvec[YY]*denom;
362     center[ZZ] = dsumvec[ZZ]*denom;
363 }
364
365
366 /* Translate x with transvec */
367 extern void translate_x(rvec x[], const int nr, const rvec transvec)
368 {
369     int i;
370     
371     
372     for (i=0; i<nr; i++)
373         rvec_inc(x[i], transvec);
374 }
375
376
377 extern void rotate_x(rvec x[], const int nr, matrix rmat)
378 {
379     int i,j,k;
380     rvec x_old;
381
382     
383     /* Apply the rotation matrix */
384     for (i=0; i<nr; i++)
385     {
386         for (j=0; j<3; j++)
387             x_old[j] = x[i][j];
388         for (j=0; j<3; j++)
389         {
390             x[i][j] = 0;
391             for (k=0; k<3; k++)
392                 x[i][j] += rmat[k][j]*x_old[k];
393         }
394     }
395 }
396