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         clear_ivec(shifts[i]);
122
123     for (i=0; i<nr; i++)
124     {
125         /* The distance this atom moved since the last time step */
126         /* If this is more than just a bit, it has changed its home pbc box */
127         rvec_sub(xcoll[i],xcoll_old[i],dx);
128
129         for(m=npbcdim-1; m>=0; m--)
130         {
131             while (dx[m] < -0.5*box[m][m])
132             {
133                 for(d=0; d<DIM; d++)
134                     dx[d] += box[m][d];
135                 shifts[i][m]++;
136             }
137             while (dx[m] >= 0.5*box[m][m])
138             {
139                 for(d=0; d<DIM; d++)
140                     dx[d] -= box[m][d];
141                 shifts[i][m]--;
142             }
143         }
144     }
145 }
146
147
148 static void shift_positions_group(
149         matrix box, 
150         rvec   x[],      /* The positions [0..nr] */ 
151         ivec   *is,      /* The shifts [0..nr] */ 
152         int    nr)       /* The number of positions and shifts */
153 {
154     int      i,tx,ty,tz;
155
156
157     /* Loop over the group's atoms */
158     if(TRICLINIC(box)) 
159     {
160         for (i=0; i < nr; i++)
161         {
162             tx=is[i][XX];
163             ty=is[i][YY];
164             tz=is[i][ZZ];
165
166             x[i][XX]=x[i][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
167             x[i][YY]=x[i][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
168             x[i][ZZ]=x[i][ZZ]+tz*box[ZZ][ZZ];
169         }
170     } else
171     {
172         for (i=0; i < nr; i++)
173         {
174             tx=is[i][XX];
175             ty=is[i][YY];
176             tz=is[i][ZZ];
177
178             x[i][XX]=x[i][XX]+tx*box[XX][XX];
179             x[i][YY]=x[i][YY]+ty*box[YY][YY];
180             x[i][ZZ]=x[i][ZZ]+tz*box[ZZ][ZZ];
181         }
182     }    
183 }
184
185
186 /* Assemble the positions of the group such that every node has all of them. 
187  * The atom indices are retrieved from anrs_loc[0..nr_loc] 
188  * Note that coll_ind[i] = i is needed in the serial case */
189 extern void communicate_group_positions(
190         t_commrec  *cr, 
191         rvec       *xcoll,        /* OUT: Collective array of positions */
192         ivec       *shifts,       /* IN+OUT: Collective array of shifts for xcoll */
193         ivec       *extra_shifts, /* BUF: Extra shifts since last time step */
194         const gmx_bool bNS,           /* IN:  NS step, the shifts have changed */
195         rvec       *x_loc,        /* IN:  Local positions on this node */ 
196         const int  nr,            /* IN:  Total number of atoms in the group */
197         const int  nr_loc,        /* IN:  Local number of atoms in the group */
198         int        *anrs_loc,     /* IN:  Local atom numbers */
199         int        *coll_ind,     /* IN:  Collective index */
200         rvec       *xcoll_old,    /* IN+OUT: Positions from the last time step, used to make group whole */
201         matrix     box)
202 {
203     int i;
204
205
206     /* Zero out the groups' global position array */
207     clear_rvecs(nr, xcoll);
208
209     /* Put the local positions that this node has into the right place of 
210      * the collective array. Note that in the serial case, coll_ind[i] = i */
211     for (i=0; i<nr_loc; i++)
212         copy_rvec(x_loc[anrs_loc[i]], xcoll[coll_ind[i]]);
213
214     if (PAR(cr))
215     {
216         /* Add the arrays from all nodes together */
217         gmx_sum(nr*3, xcoll[0], cr);
218
219         /* To make the group whole, start with a whole group and each
220          * step move the assembled positions at closest distance to the positions 
221          * from the last step. First shift the positions with the saved shift 
222          * vectors (these are 0 when this routine is called for the first time!) */
223         shift_positions_group(box, xcoll, shifts, nr);
224         
225         /* Now check if some shifts changed since the last step.
226          * This only needs to be done when the shifts are expected to have changed,
227          * i.e. after neighboursearching */
228         if (bNS) 
229         {
230             get_shifts_group(3, box, xcoll, nr, xcoll_old, extra_shifts);
231             
232             /* Shift with the additional shifts such that we get a whole group now */
233             shift_positions_group(box, xcoll, extra_shifts, nr);
234             
235             /* Add the shift vectors together for the next time step */
236             for (i=0; i<nr; i++)
237             {
238                 shifts[i][XX] += extra_shifts[i][XX];
239                 shifts[i][YY] += extra_shifts[i][YY];
240                 shifts[i][ZZ] += extra_shifts[i][ZZ];
241             }
242             
243             /* Store current correctly-shifted positions for comparison in the next NS time step */
244             for (i=0; i<nr; i++)
245                 copy_rvec(xcoll[i],xcoll_old[i]);   
246         }
247     }
248 }
249
250
251 /* Determine the (weighted) sum vector from positions x */
252 extern double get_sum_of_positions(rvec x[], real weight[], const int nat, dvec dsumvec)
253 {
254     int i;
255     rvec x_weighted;
256     double weight_sum = 0.0;
257
258
259     /* Zero out the center */
260     clear_dvec(dsumvec);
261
262     /* Loop over all atoms and add their weighted position vectors */
263     if (weight != NULL)
264     {
265         for (i=0; i<nat; i++)
266         {
267             weight_sum += weight[i];
268             svmul(weight[i], x[i], x_weighted);
269             dsumvec[XX] += x_weighted[XX];
270             dsumvec[YY] += x_weighted[YY];
271             dsumvec[ZZ] += x_weighted[ZZ];
272         }
273     }
274     else
275     {
276         for (i=0; i<nat; i++)
277         {
278             dsumvec[XX] += x[i][XX];
279             dsumvec[YY] += x[i][YY];
280             dsumvec[ZZ] += x[i][ZZ];
281         }
282     }
283     return weight_sum;
284 }
285
286
287 /* Determine center of structure from collective positions x */
288 extern void get_center(rvec x[], real weight[], const int nr, rvec rcenter)
289 {
290     dvec   dcenter;
291     double weight_sum, denom;
292
293     
294     weight_sum = get_sum_of_positions(x, weight, nr, dcenter);
295     
296     if (weight != NULL)
297         denom = weight_sum; /* Divide by the sum of weight */
298     else
299         denom = nr;        /* Divide by the number of atoms */
300         
301     dsvmul(1.0/denom, dcenter, dcenter);
302     
303     rcenter[XX] = dcenter[XX];
304     rcenter[YY] = dcenter[YY];
305     rcenter[ZZ] = dcenter[ZZ];
306 }
307
308
309 /* Get the center from local positions that already have the correct
310  * PBC representation */
311 extern void get_center_comm(
312         t_commrec *cr,
313         rvec x_loc[],       /* Local positions */
314         real weight_loc[],  /* Local masses or other weights */
315         int nr_loc,         /* Local number of atoms */
316         int nr_group,       /* Total number of atoms of the group */ 
317         rvec center)        /* Weighted center */
318 {
319     double weight_sum, denom;
320     dvec   dsumvec;
321     double buf[4];    
322     
323     
324     weight_sum = get_sum_of_positions(x_loc, weight_loc, nr_loc, dsumvec);
325     
326     /* Add the local contributions from all nodes. Put the sum vector and the 
327      * weight in a buffer array so that we get along with a single communication
328      * call. */
329     if (PAR(cr))
330     {
331         buf[0] = dsumvec[XX];
332         buf[1] = dsumvec[YY];
333         buf[2] = dsumvec[ZZ];
334         buf[3] = weight_sum;
335         
336         /* Communicate buffer */
337         gmx_sumd(4, buf, cr);
338         
339         dsumvec[XX] = buf[0];
340         dsumvec[YY] = buf[1];
341         dsumvec[ZZ] = buf[2];
342         weight_sum  = buf[3];
343     }
344     
345     if (weight_loc != NULL)
346         denom = 1.0/weight_sum; /* Divide by the sum of weight to get center of mass e.g. */
347     else
348         denom = 1.0/nr_group;   /* Divide by the number of atoms to get the geometrical center */
349         
350     center[XX] = dsumvec[XX]*denom;
351     center[YY] = dsumvec[YY]*denom;
352     center[ZZ] = dsumvec[ZZ]*denom;
353 }
354
355
356 /* Translate x with transvec */
357 extern void translate_x(rvec x[], const int nr, const rvec transvec)
358 {
359     int i;
360     
361     
362     for (i=0; i<nr; i++)
363         rvec_inc(x[i], transvec);
364 }
365
366
367 extern void rotate_x(rvec x[], const int nr, matrix rmat)
368 {
369     int i,j,k;
370     rvec x_old;
371
372     
373     /* Apply the rotation matrix */
374     for (i=0; i<nr; i++)
375     {
376         for (j=0; j<3; j++)
377             x_old[j] = x[i][j];
378         for (j=0; j<3; j++)
379         {
380             x[i][j] = 0;
381             for (k=0; k<3; k++)
382                 x[i][j] += rmat[k][j]*x_old[k];
383         }
384     }
385 }
386