Simplified uniform GPU selection in CMake
[gromacs.git] / src / gromacs / mdlib / groupcoord.cpp
blob34a7ca9156c69ee72bf1098a7b43d1a7f41d6bfa
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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, 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.
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.
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.
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.
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.
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.
38 #include "gmxpre.h"
40 #include "groupcoord.h"
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"
49 #define MIN(a, b) (((a) < (b)) ? (a) : (b))
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()
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? */
64 GMX_ASSERT(ga2la, "We need a valid ga2la object");
66 /* Loop over all the atom indices of the group to check
67 * which ones are on the local node */
68 int localnr = 0;
69 for (int i = 0; i < nr; i++)
71 if (const int* ii = ga2la->findHome(anrs[i]))
73 /* The atom with this index is a home atom */
74 if (localnr >= *nalloc_loc) /* Check whether memory suffices */
76 *nalloc_loc = over_alloc_dd(localnr + 1);
77 /* We never need more memory than the number of atoms in the group */
78 *nalloc_loc = MIN(*nalloc_loc, nr);
79 srenew(*anrs_loc, *nalloc_loc);
81 /* Save the atoms index in the local atom numbers array */
82 (*anrs_loc)[localnr] = *ii;
84 if (coll_ind != nullptr)
86 /* Keep track of where this local atom belongs in the collective index array.
87 * This is needed when reducing the local arrays to a collective/global array
88 * in communicate_group_positions */
89 coll_ind[localnr] = i;
92 /* add one to the local atom count */
93 localnr++;
97 /* Return the number of local atoms that were found */
98 *nr_loc = localnr;
102 static void get_shifts_group(int npbcdim,
103 const matrix box,
104 rvec* xcoll, /* IN: Collective set of positions [0..nr] */
105 int nr, /* IN: Total number of atoms in the group */
106 rvec* xcoll_old, /* IN: Positions from the last time step [0...nr] */
107 ivec* shifts) /* OUT: Shifts for xcoll */
109 int i, m, d;
110 rvec dx;
113 /* Get the shifts such that each atom is within closest
114 * distance to its position at the last NS time step after shifting.
115 * If we start with a whole group, and always keep track of
116 * shift changes, the group will stay whole this way */
117 for (i = 0; i < nr; i++)
119 clear_ivec(shifts[i]);
122 for (i = 0; i < nr; i++)
124 /* The distance this atom moved since the last time step */
125 /* If this is more than just a bit, it has changed its home pbc box */
126 rvec_sub(xcoll[i], xcoll_old[i], dx);
128 for (m = npbcdim - 1; m >= 0; m--)
130 while (dx[m] < -0.5 * box[m][m])
132 for (d = 0; d < DIM; d++)
134 dx[d] += box[m][d];
136 shifts[i][m]++;
138 while (dx[m] >= 0.5 * box[m][m])
140 for (d = 0; d < DIM; d++)
142 dx[d] -= box[m][d];
144 shifts[i][m]--;
151 static void shift_positions_group(const matrix box,
152 rvec x[], /* The positions [0..nr] */
153 ivec* is, /* The shifts [0..nr] */
154 int nr) /* The number of positions and shifts */
156 int i, tx, ty, tz;
159 /* Loop over the group's atoms */
160 if (TRICLINIC(box))
162 for (i = 0; i < nr; i++)
164 tx = is[i][XX];
165 ty = is[i][YY];
166 tz = is[i][ZZ];
168 x[i][XX] = x[i][XX] + tx * box[XX][XX] + ty * box[YY][XX] + tz * box[ZZ][XX];
169 x[i][YY] = x[i][YY] + ty * box[YY][YY] + tz * box[ZZ][YY];
170 x[i][ZZ] = x[i][ZZ] + tz * box[ZZ][ZZ];
173 else
175 for (i = 0; i < nr; i++)
177 tx = is[i][XX];
178 ty = is[i][YY];
179 tz = is[i][ZZ];
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];
189 /* Assemble the positions of the group such that every node has all of them.
190 * The atom indices are retrieved from anrs_loc[0..nr_loc]
191 * Note that coll_ind[i] = i is needed in the serial case */
192 extern void communicate_group_positions(const t_commrec* cr, /* Pointer to MPI communication data */
193 rvec* xcoll, /* Collective array of positions */
194 ivec* shifts, /* Collective array of shifts for xcoll (can be NULL) */
195 ivec* extra_shifts, /* (optional) Extra shifts since last time step */
196 const gmx_bool bNS, /* (optional) NS step, the shifts have changed */
197 const rvec* x_loc, /* Local positions on this node */
198 const int nr, /* Total number of atoms in the group */
199 const int nr_loc, /* Local number of atoms in the group */
200 const int* anrs_loc, /* Local atom numbers */
201 const int* coll_ind, /* Collective index */
202 rvec* xcoll_old, /* (optional) Positions from the last time
203 step, used to make group whole */
204 const matrix box) /* (optional) The box */
206 int i;
209 /* Zero out the groups' global position array */
210 clear_rvecs(nr, xcoll);
212 /* Put the local positions that this node has into the right place of
213 * the collective array. Note that in the serial case, coll_ind[i] = i */
214 for (i = 0; i < nr_loc; i++)
216 copy_rvec(x_loc[anrs_loc[i]], xcoll[coll_ind[i]]);
219 if (PAR(cr))
221 /* Add the arrays from all nodes together */
222 gmx_sum(nr * 3, xcoll[0], cr);
224 /* Now we have all the positions of the group in the xcoll array present on all
225 * nodes.
227 * The rest of the code is for making the group whole again in case atoms changed
228 * their PBC representation / crossed a box boundary. We only do that if the
229 * shifts array is allocated. */
230 if (nullptr != shifts)
232 /* To make the group whole, start with a whole group and each
233 * step move the assembled positions at closest distance to the positions
234 * from the last step. First shift the positions with the saved shift
235 * vectors (these are 0 when this routine is called for the first time!) */
236 shift_positions_group(box, xcoll, shifts, nr);
238 /* Now check if some shifts changed since the last step.
239 * This only needs to be done when the shifts are expected to have changed,
240 * i.e. after neighbor searching */
241 if (bNS)
243 get_shifts_group(3, box, xcoll, nr, xcoll_old, extra_shifts);
245 /* Shift with the additional shifts such that we get a whole group now */
246 shift_positions_group(box, xcoll, extra_shifts, nr);
248 /* Add the shift vectors together for the next time step */
249 for (i = 0; i < nr; i++)
251 shifts[i][XX] += extra_shifts[i][XX];
252 shifts[i][YY] += extra_shifts[i][YY];
253 shifts[i][ZZ] += extra_shifts[i][ZZ];
256 /* Store current correctly-shifted positions for comparison in the next NS time step */
257 for (i = 0; i < nr; i++)
259 copy_rvec(xcoll[i], xcoll_old[i]);
266 /* Determine the (weighted) sum vector from positions x */
267 extern double get_sum_of_positions(rvec x[], real weight[], const int nat, dvec dsumvec)
269 int i;
270 rvec x_weighted;
271 double weight_sum = 0.0;
274 /* Zero out the center */
275 clear_dvec(dsumvec);
277 /* Loop over all atoms and add their weighted position vectors */
278 if (weight != nullptr)
280 for (i = 0; i < nat; i++)
282 weight_sum += weight[i];
283 svmul(weight[i], x[i], x_weighted);
284 dsumvec[XX] += x_weighted[XX];
285 dsumvec[YY] += x_weighted[YY];
286 dsumvec[ZZ] += x_weighted[ZZ];
289 else
291 for (i = 0; i < nat; i++)
293 dsumvec[XX] += x[i][XX];
294 dsumvec[YY] += x[i][YY];
295 dsumvec[ZZ] += x[i][ZZ];
298 return weight_sum;
302 /* Determine center of structure from collective positions x */
303 extern void get_center(rvec x[], real weight[], const int nr, rvec rcenter)
305 dvec dcenter;
306 double weight_sum, denom;
309 weight_sum = get_sum_of_positions(x, weight, nr, dcenter);
311 if (weight != nullptr)
313 denom = weight_sum; /* Divide by the sum of weight */
315 else
317 denom = nr; /* Divide by the number of atoms */
319 dsvmul(1.0 / denom, dcenter, dcenter);
321 rcenter[XX] = dcenter[XX];
322 rcenter[YY] = dcenter[YY];
323 rcenter[ZZ] = dcenter[ZZ];
327 /* Get the center from local positions that already have the correct
328 * PBC representation */
329 extern void get_center_comm(const t_commrec* cr,
330 rvec x_loc[], /* Local positions */
331 real weight_loc[], /* Local masses or other weights */
332 int nr_loc, /* Local number of atoms */
333 int nr_group, /* Total number of atoms of the group */
334 rvec center) /* Weighted center */
336 double weight_sum, denom;
337 dvec dsumvec;
338 double buf[4];
341 weight_sum = get_sum_of_positions(x_loc, weight_loc, nr_loc, dsumvec);
343 /* Add the local contributions from all nodes. Put the sum vector and the
344 * weight in a buffer array so that we get along with a single communication
345 * call. */
346 if (PAR(cr))
348 buf[0] = dsumvec[XX];
349 buf[1] = dsumvec[YY];
350 buf[2] = dsumvec[ZZ];
351 buf[3] = weight_sum;
353 /* Communicate buffer */
354 gmx_sumd(4, buf, cr);
356 dsumvec[XX] = buf[0];
357 dsumvec[YY] = buf[1];
358 dsumvec[ZZ] = buf[2];
359 weight_sum = buf[3];
362 if (weight_loc != nullptr)
364 denom = 1.0 / weight_sum; /* Divide by the sum of weight to get center of mass e.g. */
366 else
368 denom = 1.0 / nr_group; /* Divide by the number of atoms to get the geometrical center */
370 center[XX] = dsumvec[XX] * denom;
371 center[YY] = dsumvec[YY] * denom;
372 center[ZZ] = dsumvec[ZZ] * denom;
376 /* Translate x with transvec */
377 extern void translate_x(rvec x[], const int nr, const rvec transvec)
379 int i;
382 for (i = 0; i < nr; i++)
384 rvec_inc(x[i], transvec);
389 extern void rotate_x(rvec x[], const int nr, matrix rmat)
391 int i, j, k;
392 rvec x_old;
395 /* Apply the rotation matrix */
396 for (i = 0; i < nr; i++)
398 for (j = 0; j < 3; j++)
400 x_old[j] = x[i][j];
402 for (j = 0; j < 3; j++)
404 x[i][j] = 0;
405 for (k = 0; k < 3; k++)
407 x[i][j] += rmat[k][j] * x_old[k];