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