Fix UB when generating local indices for constraints
[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,2018 by the GROMACS development team.
7  * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source 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 #include "gmxpre.h"
39
40 #include "groupcoord.h"
41
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 /* 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 void dd_make_local_group_indices(const 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     GMX_ASSERT(ga2la, "We need a valid ga2la object");
65     GMX_RELEASE_ASSERT(anrs != *anrs_loc, "Can not update indices in-place");
66
67     /* Loop over all the atom indices of the group to check
68      * which ones are on the local node */
69     int localnr = 0;
70     for (int i = 0; i < nr; i++)
71     {
72         if (const int* ii = ga2la->findHome(anrs[i]))
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 != nullptr)
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(int          npbcdim,
104                              const matrix box,
105                              rvec*        xcoll, /* IN:  Collective set of positions [0..nr] */
106                              int          nr,    /* IN:  Total number of atoms in the group */
107                              rvec* xcoll_old, /* IN:  Positions from the last time step [0...nr] */
108                              ivec* shifts)    /* OUT: Shifts for xcoll */
109 {
110     int  i, m, d;
111     rvec dx;
112
113
114     /* Get the shifts such that each atom is within closest
115      * distance to its position at the last NS time step after shifting.
116      * If we start with a whole group, and always keep track of
117      * shift changes, the group will stay whole this way */
118     for (i = 0; i < nr; i++)
119     {
120         clear_ivec(shifts[i]);
121     }
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                 {
135                     dx[d] += box[m][d];
136                 }
137                 shifts[i][m]++;
138             }
139             while (dx[m] >= 0.5 * box[m][m])
140             {
141                 for (d = 0; d < DIM; d++)
142                 {
143                     dx[d] -= box[m][d];
144                 }
145                 shifts[i][m]--;
146             }
147         }
148     }
149 }
150
151
152 static void shift_positions_group(const matrix box,
153                                   rvec         x[], /* The positions [0..nr] */
154                                   ivec*        is,  /* The shifts [0..nr] */
155                                   int          nr)           /* The number of positions and shifts */
156 {
157     int i, tx, ty, tz;
158
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     }
174     else
175     {
176         for (i = 0; i < nr; i++)
177         {
178             tx = is[i][XX];
179             ty = is[i][YY];
180             tz = is[i][ZZ];
181
182             x[i][XX] = x[i][XX] + tx * box[XX][XX];
183             x[i][YY] = x[i][YY] + ty * box[YY][YY];
184             x[i][ZZ] = x[i][ZZ] + tz * box[ZZ][ZZ];
185         }
186     }
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(const t_commrec* cr, /* Pointer to MPI communication data */
194                                         rvec*            xcoll, /* Collective array of positions */
195                                         ivec* shifts, /* Collective array of shifts for xcoll (can be NULL) */
196                                         ivec* extra_shifts, /* (optional) Extra shifts since last time step */
197                                         const gmx_bool bNS, /* (optional) NS step, the shifts have changed */
198                                         const rvec* x_loc,  /* Local positions on this node */
199                                         const int   nr,     /* Total number of atoms in the group */
200                                         const int   nr_loc, /* Local number of atoms in the group */
201                                         const int*  anrs_loc, /* Local atom numbers */
202                                         const int*  coll_ind, /* Collective index */
203                                         rvec* xcoll_old,  /* (optional) Positions from the last time
204                                                              step,  used to make group whole */
205                                         const matrix box) /* (optional) The box */
206 {
207     int i;
208
209
210     /* Zero out the groups' global position array */
211     clear_rvecs(nr, xcoll);
212
213     /* Put the local positions that this node has into the right place of
214      * the collective array. Note that in the serial case, coll_ind[i] = i */
215     for (i = 0; i < nr_loc; i++)
216     {
217         copy_rvec(x_loc[anrs_loc[i]], xcoll[coll_ind[i]]);
218     }
219
220     if (PAR(cr))
221     {
222         /* Add the arrays from all nodes together */
223         gmx_sum(nr * 3, xcoll[0], cr);
224     }
225     /* Now we have all the positions of the group in the xcoll array present on all
226      * nodes.
227      *
228      * The rest of the code is for making the group whole again in case atoms changed
229      * their PBC representation / crossed a box boundary. We only do that if the
230      * shifts array is allocated. */
231     if (nullptr != shifts)
232     {
233         /* To make the group whole, start with a whole group and each
234          * step move the assembled positions at closest distance to the positions
235          * from the last step. First shift the positions with the saved shift
236          * vectors (these are 0 when this routine is called for the first time!) */
237         shift_positions_group(box, xcoll, shifts, nr);
238
239         /* Now check if some shifts changed since the last step.
240          * This only needs to be done when the shifts are expected to have changed,
241          * i.e. after neighbor searching */
242         if (bNS)
243         {
244             get_shifts_group(3, box, xcoll, nr, xcoll_old, extra_shifts);
245
246             /* Shift with the additional shifts such that we get a whole group now */
247             shift_positions_group(box, xcoll, extra_shifts, nr);
248
249             /* Add the shift vectors together for the next time step */
250             for (i = 0; i < nr; i++)
251             {
252                 shifts[i][XX] += extra_shifts[i][XX];
253                 shifts[i][YY] += extra_shifts[i][YY];
254                 shifts[i][ZZ] += extra_shifts[i][ZZ];
255             }
256
257             /* Store current correctly-shifted positions for comparison in the next NS time step */
258             for (i = 0; i < nr; i++)
259             {
260                 copy_rvec(xcoll[i], xcoll_old[i]);
261             }
262         }
263     }
264 }
265
266
267 /* Determine the (weighted) sum vector from positions x */
268 extern double get_sum_of_positions(rvec x[], real weight[], const int nat, dvec dsumvec)
269 {
270     int    i;
271     rvec   x_weighted;
272     double weight_sum = 0.0;
273
274
275     /* Zero out the center */
276     clear_dvec(dsumvec);
277
278     /* Loop over all atoms and add their weighted position vectors */
279     if (weight != nullptr)
280     {
281         for (i = 0; i < nat; i++)
282         {
283             weight_sum += weight[i];
284             svmul(weight[i], x[i], x_weighted);
285             dsumvec[XX] += x_weighted[XX];
286             dsumvec[YY] += x_weighted[YY];
287             dsumvec[ZZ] += x_weighted[ZZ];
288         }
289     }
290     else
291     {
292         for (i = 0; i < nat; i++)
293         {
294             dsumvec[XX] += x[i][XX];
295             dsumvec[YY] += x[i][YY];
296             dsumvec[ZZ] += x[i][ZZ];
297         }
298     }
299     return weight_sum;
300 }
301
302
303 /* Determine center of structure from collective positions x */
304 extern void get_center(rvec x[], real weight[], const int nr, rvec rcenter)
305 {
306     dvec   dcenter;
307     double weight_sum, denom;
308
309
310     weight_sum = get_sum_of_positions(x, weight, nr, dcenter);
311
312     if (weight != nullptr)
313     {
314         denom = weight_sum; /* Divide by the sum of weight */
315     }
316     else
317     {
318         denom = nr; /* Divide by the number of atoms */
319     }
320     dsvmul(1.0 / denom, dcenter, dcenter);
321
322     rcenter[XX] = dcenter[XX];
323     rcenter[YY] = dcenter[YY];
324     rcenter[ZZ] = dcenter[ZZ];
325 }
326
327
328 /* Get the center from local positions that already have the correct
329  * PBC representation */
330 extern void get_center_comm(const t_commrec* cr,
331                             rvec             x_loc[],      /* Local positions */
332                             real             weight_loc[], /* Local masses or other weights */
333                             int              nr_loc,       /* Local number of atoms */
334                             int              nr_group,     /* Total number of atoms of the group */
335                             rvec             center)                   /* Weighted center */
336 {
337     double weight_sum, denom;
338     dvec   dsumvec;
339     double buf[4];
340
341
342     weight_sum = get_sum_of_positions(x_loc, weight_loc, nr_loc, dsumvec);
343
344     /* Add the local contributions from all nodes. Put the sum vector and the
345      * weight in a buffer array so that we get along with a single communication
346      * call. */
347     if (PAR(cr))
348     {
349         buf[0] = dsumvec[XX];
350         buf[1] = dsumvec[YY];
351         buf[2] = dsumvec[ZZ];
352         buf[3] = weight_sum;
353
354         /* Communicate buffer */
355         gmx_sumd(4, buf, cr);
356
357         dsumvec[XX] = buf[0];
358         dsumvec[YY] = buf[1];
359         dsumvec[ZZ] = buf[2];
360         weight_sum  = buf[3];
361     }
362
363     if (weight_loc != nullptr)
364     {
365         denom = 1.0 / weight_sum; /* Divide by the sum of weight to get center of mass e.g. */
366     }
367     else
368     {
369         denom = 1.0 / nr_group; /* Divide by the number of atoms to get the geometrical center */
370     }
371     center[XX] = dsumvec[XX] * denom;
372     center[YY] = dsumvec[YY] * denom;
373     center[ZZ] = dsumvec[ZZ] * denom;
374 }
375
376
377 /* Translate x with transvec */
378 extern void translate_x(rvec x[], const int nr, const rvec transvec)
379 {
380     int i;
381
382
383     for (i = 0; i < nr; i++)
384     {
385         rvec_inc(x[i], transvec);
386     }
387 }
388
389
390 extern void rotate_x(rvec x[], const int nr, matrix rmat)
391 {
392     int  i, j, k;
393     rvec x_old;
394
395
396     /* Apply the rotation matrix */
397     for (i = 0; i < nr; i++)
398     {
399         for (j = 0; j < 3; j++)
400         {
401             x_old[j] = x[i][j];
402         }
403         for (j = 0; j < 3; j++)
404         {
405             x[i][j] = 0;
406             for (k = 0; k < 3; k++)
407             {
408                 x[i][j] += rmat[k][j] * x_old[k];
409             }
410         }
411     }
412 }