Separate management of GPU contexts from modules
[gromacs.git] / src / gromacs / selection / centerofmass.cpp
blob0639890a0db148b09f437eda7c01997ed4340c81
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2016, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
35 /*! \internal \file
36 * \brief
37 * Implements functions in centerofmass.h.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
42 #include "gmxpre.h"
44 #include "centerofmass.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/pbcutil/pbc.h"
48 #include "gromacs/topology/block.h"
49 #include "gromacs/topology/mtop_lookup.h"
50 #include "gromacs/topology/topology.h"
51 #include "gromacs/utility/gmxassert.h"
53 void
54 gmx_calc_cog(const gmx_mtop_t * /* top */, rvec x[], int nrefat, const int index[], rvec xout)
56 int m, ai;
58 clear_rvec(xout);
59 for (m = 0; m < nrefat; ++m)
61 ai = index[m];
62 rvec_inc(xout, x[ai]);
64 svmul(1.0/nrefat, xout, xout);
67 /*!
68 * \param[in] top Topology structure with masses.
69 * \param[in] x Position vectors of all atoms.
70 * \param[in] nrefat Number of atoms in the index.
71 * \param[in] index Indices of atoms.
72 * \param[out] xout COM position for the indexed atoms.
74 * Works exactly as gmx_calc_cog() with the exception that a center of
75 * mass are calculated, and hence a topology with masses is required.
77 void
78 gmx_calc_com(const gmx_mtop_t *top, rvec x[], int nrefat, const int index[], rvec xout)
80 GMX_RELEASE_ASSERT(gmx_mtop_has_masses(top),
81 "No masses available while mass weighting was requested");
82 clear_rvec(xout);
83 real mtot = 0;
84 int molb = 0;
85 for (int m = 0; m < nrefat; ++m)
87 const int ai = index[m];
88 const real mass = mtopGetAtomMass(top, ai, &molb);
89 for (int j = 0; j < DIM; ++j)
91 xout[j] += mass * x[ai][j];
93 mtot += mass;
95 svmul(1.0/mtot, xout, xout);
98 /*!
99 * \param[in] top Topology structure with masses.
100 * \param[in] f Forces on all atoms.
101 * \param[in] nrefat Number of atoms in the index.
102 * \param[in] index Indices of atoms.
103 * \param[out] fout Force on the COG position for the indexed atoms.
105 void
106 gmx_calc_cog_f(const gmx_mtop_t *top, rvec f[], int nrefat, const int index[], rvec fout)
108 GMX_RELEASE_ASSERT(gmx_mtop_has_masses(top),
109 "No masses available while mass weighting was requested");
110 clear_rvec(fout);
111 real mtot = 0;
112 int molb = 0;
113 for (int m = 0; m < nrefat; ++m)
115 const int ai = index[m];
116 const real mass = mtopGetAtomMass(top, ai, &molb);
117 for (int j = 0; j < DIM; ++j)
119 fout[j] += f[ai][j] / mass;
121 mtot += mass;
123 svmul(mtot / nrefat, fout, fout);
126 void
127 gmx_calc_com_f(const gmx_mtop_t * /* top */, rvec f[], int nrefat, const int index[], rvec fout)
129 clear_rvec(fout);
130 for (int m = 0; m < nrefat; ++m)
132 const int ai = index[m];
133 rvec_inc(fout, f[ai]);
138 * \param[in] top Topology structure with masses
139 * (can be NULL if \p bMASS==false).
140 * \param[in] x Position vectors of all atoms.
141 * \param[in] nrefat Number of atoms in the index.
142 * \param[in] index Indices of atoms.
143 * \param[in] bMass If true, mass weighting is used.
144 * \param[out] xout COM/COG position for the indexed atoms.
146 * Calls either gmx_calc_com() or gmx_calc_cog() depending on the value of
147 * \p bMass.
148 * Other parameters are passed unmodified to these functions.
150 void
151 gmx_calc_comg(const gmx_mtop_t *top, rvec x[], int nrefat, const int index[],
152 bool bMass, rvec xout)
154 if (bMass)
156 gmx_calc_com(top, x, nrefat, index, xout);
158 else
160 gmx_calc_cog(top, x, nrefat, index, xout);
165 * \param[in] top Topology structure with masses
166 * (can be NULL if \p bMASS==true).
167 * \param[in] f Forces on all atoms.
168 * \param[in] nrefat Number of atoms in the index.
169 * \param[in] index Indices of atoms.
170 * \param[in] bMass If true, force on COM is calculated.
171 * \param[out] fout Force on the COM/COG position for the indexed atoms.
173 * Calls either gmx_calc_cog_f() or gmx_calc_com_f() depending on the value of
174 * \p bMass.
175 * Other parameters are passed unmodified to these functions.
177 void
178 gmx_calc_comg_f(const gmx_mtop_t *top, rvec f[], int nrefat, const int index[],
179 bool bMass, rvec fout)
181 if (bMass)
183 gmx_calc_com_f(top, f, nrefat, index, fout);
185 else
187 gmx_calc_cog_f(top, f, nrefat, index, fout);
193 * \param[in] top Topology structure (unused, can be NULL).
194 * \param[in] x Position vectors of all atoms.
195 * \param[in] pbc Periodic boundary conditions structure.
196 * \param[in] nrefat Number of atoms in the index.
197 * \param[in] index Indices of atoms.
198 * \param[out] xout COG position for the indexed atoms.
200 * Works exactly as gmx_calc_com_pbc(), but calculates the center of geometry.
202 void
203 gmx_calc_cog_pbc(const gmx_mtop_t *top, rvec x[], const t_pbc *pbc,
204 int nrefat, const int index[], rvec xout)
206 const real tol = 1e-4;
207 bool bChanged;
208 int m, j, ai, iter;
209 rvec dx, xtest;
211 /* First simple calculation */
212 gmx_calc_cog(top, x, nrefat, index, xout);
213 /* Now check if any atom is more than half the box from the COM */
214 if (pbc)
216 iter = 0;
219 bChanged = false;
220 for (m = 0; m < nrefat; ++m)
222 ai = index[m];
223 pbc_dx(pbc, x[ai], xout, dx);
224 rvec_add(xout, dx, xtest);
225 for (j = 0; j < DIM; ++j)
227 if (fabs(xtest[j] - x[ai][j]) > tol)
229 /* Here we have used the wrong image for contributing to the COM */
230 xout[j] += (xtest[j] - x[ai][j]) / nrefat;
231 x[ai][j] = xtest[j];
232 bChanged = true;
236 iter++;
238 while (bChanged);
243 * \param[in] top Topology structure with masses.
244 * \param[in] x Position vectors of all atoms.
245 * \param[in] pbc Periodic boundary conditions structure.
246 * \param[in] nrefat Number of atoms in the index.
247 * \param[in] index Indices of atoms.
248 * \param[out] xout COM position for the indexed atoms.
250 * Works as gmx_calc_com(), but takes into account periodic boundary
251 * conditions: If any atom is more than half the box from the COM,
252 * it is wrapped around and a new COM is calculated. This is repeated
253 * until no atoms violate the condition.
255 * Modified from src/tools/gmx_sorient.c in Gromacs distribution.
257 void
258 gmx_calc_com_pbc(const gmx_mtop_t *top, rvec x[], const t_pbc *pbc,
259 int nrefat, const int index[], rvec xout)
261 GMX_RELEASE_ASSERT(gmx_mtop_has_masses(top),
262 "No masses available while mass weighting was requested");
263 /* First simple calculation */
264 clear_rvec(xout);
265 real mtot = 0;
266 int molb = 0;
267 for (int m = 0; m < nrefat; ++m)
269 const int ai = index[m];
270 const real mass = mtopGetAtomMass(top, ai, &molb);
271 for (int j = 0; j < DIM; ++j)
273 xout[j] += mass * x[ai][j];
275 mtot += mass;
277 svmul(1.0/mtot, xout, xout);
278 /* Now check if any atom is more than half the box from the COM */
279 if (pbc)
281 const real tol = 1e-4;
282 bool bChanged;
285 bChanged = false;
286 molb = 0;
287 for (int m = 0; m < nrefat; ++m)
289 rvec dx, xtest;
290 const int ai = index[m];
291 const real mass = mtopGetAtomMass(top, ai, &molb) / mtot;
292 pbc_dx(pbc, x[ai], xout, dx);
293 rvec_add(xout, dx, xtest);
294 for (int j = 0; j < DIM; ++j)
296 if (fabs(xtest[j] - x[ai][j]) > tol)
298 /* Here we have used the wrong image for contributing to the COM */
299 xout[j] += mass * (xtest[j] - x[ai][j]);
300 x[ai][j] = xtest[j];
301 bChanged = true;
306 while (bChanged);
311 * \param[in] top Topology structure with masses
312 * (can be NULL if \p bMASS==false).
313 * \param[in] x Position vectors of all atoms.
314 * \param[in] pbc Periodic boundary conditions structure.
315 * \param[in] nrefat Number of atoms in the index.
316 * \param[in] index Indices of atoms.
317 * \param[in] bMass If true, mass weighting is used.
318 * \param[out] xout COM/COG position for the indexed atoms.
320 * Calls either gmx_calc_com() or gmx_calc_cog() depending on the value of
321 * \p bMass.
322 * Other parameters are passed unmodified to these functions.
324 void
325 gmx_calc_comg_pbc(const gmx_mtop_t *top, rvec x[], const t_pbc *pbc,
326 int nrefat, const int index[], bool bMass, rvec xout)
328 if (bMass)
330 gmx_calc_com_pbc(top, x, pbc, nrefat, index, xout);
332 else
334 gmx_calc_cog_pbc(top, x, pbc, nrefat, index, xout);
339 void
340 gmx_calc_cog_block(const gmx_mtop_t * /* top */, rvec x[], const t_block *block, const int index[],
341 rvec xout[])
343 int b, i, ai;
344 rvec xb;
346 for (b = 0; b < block->nr; ++b)
348 clear_rvec(xb);
349 for (i = block->index[b]; i < block->index[b+1]; ++i)
351 ai = index[i];
352 rvec_inc(xb, x[ai]);
354 svmul(1.0/(block->index[b+1] - block->index[b]), xb, xout[b]);
359 * \param[in] top Topology structure with masses.
360 * \param[in] x Position vectors of all atoms.
361 * \param[in] block t_block structure that divides \p index into blocks.
362 * \param[in] index Indices of atoms.
363 * \param[out] xout \p block->nr COM positions.
365 * Works exactly as gmx_calc_cog_block() with the exception that centers of
366 * mass are calculated, and hence a topology with masses is required.
368 void
369 gmx_calc_com_block(const gmx_mtop_t *top, rvec x[], const t_block *block, const int index[],
370 rvec xout[])
372 GMX_RELEASE_ASSERT(gmx_mtop_has_masses(top),
373 "No masses available while mass weighting was requested");
374 int molb = 0;
375 for (int b = 0; b < block->nr; ++b)
377 rvec xb;
378 clear_rvec(xb);
379 real mtot = 0;
380 for (int i = block->index[b]; i < block->index[b+1]; ++i)
382 const int ai = index[i];
383 const real mass = mtopGetAtomMass(top, ai, &molb);
384 for (int d = 0; d < DIM; ++d)
386 xb[d] += mass * x[ai][d];
388 mtot += mass;
390 svmul(1.0/mtot, xb, xout[b]);
395 * \param[in] top Topology structure with masses.
396 * \param[in] f Forces on all atoms.
397 * \param[in] block t_block structure that divides \p index into blocks.
398 * \param[in] index Indices of atoms.
399 * \param[out] fout \p block->nr Forces on COG positions.
401 void
402 gmx_calc_cog_f_block(const gmx_mtop_t *top, rvec f[], const t_block *block, const int index[],
403 rvec fout[])
405 GMX_RELEASE_ASSERT(gmx_mtop_has_masses(top),
406 "No masses available while mass weighting was requested");
407 int molb = 0;
408 for (int b = 0; b < block->nr; ++b)
410 rvec fb;
411 clear_rvec(fb);
412 real mtot = 0;
413 for (int i = block->index[b]; i < block->index[b+1]; ++i)
415 const int ai = index[i];
416 const real mass = mtopGetAtomMass(top, ai, &molb);
417 for (int d = 0; d < DIM; ++d)
419 fb[d] += f[ai][d] / mass;
421 mtot += mass;
423 svmul(mtot / (block->index[b+1] - block->index[b]), fb, fout[b]);
427 void
428 gmx_calc_com_f_block(const gmx_mtop_t * /* top */, rvec f[], const t_block *block, const int index[],
429 rvec fout[])
431 for (int b = 0; b < block->nr; ++b)
433 rvec fb;
434 clear_rvec(fb);
435 for (int i = block->index[b]; i < block->index[b+1]; ++i)
437 const int ai = index[i];
438 rvec_inc(fb, f[ai]);
440 copy_rvec(fb, fout[b]);
445 * \param[in] top Topology structure with masses
446 * (can be NULL if \p bMASS==false).
447 * \param[in] x Position vectors of all atoms.
448 * \param[in] block t_block structure that divides \p index into blocks.
449 * \param[in] index Indices of atoms.
450 * \param[in] bMass If true, mass weighting is used.
451 * \param[out] xout \p block->nr COM/COG positions.
453 * Calls either gmx_calc_com_block() or gmx_calc_cog_block() depending on the
454 * value of \p bMass.
455 * Other parameters are passed unmodified to these functions.
457 void
458 gmx_calc_comg_block(const gmx_mtop_t *top, rvec x[], const t_block *block, const int index[],
459 bool bMass, rvec xout[])
461 if (bMass)
463 gmx_calc_com_block(top, x, block, index, xout);
465 else
467 gmx_calc_cog_block(top, x, block, index, xout);
472 * \param[in] top Topology structure with masses
473 * (can be NULL if \p bMASS==false).
474 * \param[in] f Forces on all atoms.
475 * \param[in] block t_block structure that divides \p index into blocks.
476 * \param[in] index Indices of atoms.
477 * \param[in] bMass If true, force on COM is calculated.
478 * \param[out] fout \p block->nr forces on the COM/COG positions.
480 * Calls either gmx_calc_com_f_block() or gmx_calc_cog_f_block() depending on
481 * the value of \p bMass.
482 * Other parameters are passed unmodified to these functions.
484 void
485 gmx_calc_comg_f_block(const gmx_mtop_t *top, rvec f[], const t_block *block, const int index[],
486 bool bMass, rvec fout[])
488 if (bMass)
490 gmx_calc_com_f_block(top, f, block, index, fout);
492 else
494 gmx_calc_cog_f_block(top, f, block, index, fout);
499 * \param[in] top Topology structure with masses
500 * (can be NULL if \p bMASS==false).
501 * \param[in] x Position vectors of all atoms.
502 * \param[in] block Blocks for calculation.
503 * \param[in] bMass If true, mass weighting is used.
504 * \param[out] xout \p block->nr COM/COG positions.
506 * Calls gmx_calc_comg_block(), converting the \p t_blocka structure into
507 * a \p t_block and an index. Other parameters are passed unmodified.
509 * \attention
510 * This function assumes that a pointer to \c t_blocka can be safely typecast
511 * into \c t_block such that the index fields can still be referenced.
512 * With the present Gromacs defitions of these types, this is the case,
513 * but if the layout of these structures is changed, this may lead to strange
514 * crashes.
516 void
517 gmx_calc_comg_blocka(const gmx_mtop_t *top, rvec x[], const t_blocka *block,
518 bool bMass, rvec xout[])
520 /* TODO: It would probably be better to do this without the type cast */
521 gmx_calc_comg_block(top, x, reinterpret_cast<const t_block *>(block), block->a, bMass, xout);
525 * \param[in] top Topology structure with masses
526 * (can be NULL if \p bMASS==true).
527 * \param[in] f Forces on all atoms.
528 * \param[in] block Blocks for calculation.
529 * \param[in] bMass If true, force on COM is calculated.
530 * \param[out] fout \p block->nr forces on the COM/COG positions.
532 * Calls gmx_calc_comg_f_block(), converting the \p t_blocka structure into
533 * a \p t_block and an index. Other parameters are passed unmodified.
535 * \attention
536 * This function assumes that a pointer to \c t_blocka can be safely typecast
537 * into \c t_block such that the index fields can still be referenced.
538 * With the present Gromacs defitions of these types, this is the case,
539 * but if the layout of these structures is changed, this may lead to strange
540 * crashes.
542 void
543 gmx_calc_comg_f_blocka(const gmx_mtop_t *top, rvec f[], const t_blocka *block,
544 bool bMass, rvec fout[])
546 /* TODO: It would probably be better to do this without the type cast */
547 gmx_calc_comg_f_block(top, f, reinterpret_cast<const t_block *>(block), block->a, bMass, fout);