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