2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014, 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.
37 * Implements functions in centerofmass.h.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
42 #include "gromacs/selection/centerofmass.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/pbcutil/pbc.h"
48 #include "gromacs/topology/block.h"
49 #include "gromacs/topology/topology.h"
52 gmx_calc_cog(t_topology * /* top */, rvec x[], int nrefat, atom_id index[], rvec xout)
57 for (m = 0; m < nrefat; ++m)
60 rvec_inc(xout, x[ai]);
62 svmul(1.0/nrefat, xout, xout);
67 * \param[in] top Topology structure with masses.
68 * \param[in] x Position vectors of all atoms.
69 * \param[in] nrefat Number of atoms in the index.
70 * \param[in] index Indices of atoms.
71 * \param[out] xout COM position for the indexed atoms.
72 * \returns 0 on success, EINVAL if \p top is NULL.
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.
78 gmx_calc_com(t_topology *top, rvec x[], int nrefat, atom_id index[], rvec xout)
85 gmx_incons("no masses available while mass weighting was requested");
90 for (m = 0; m < nrefat; ++m)
93 mass = top->atoms.atom[ai].m;
94 for (j = 0; j < DIM; ++j)
96 xout[j] += mass * x[ai][j];
100 svmul(1.0/mtot, xout, xout);
105 * \param[in] top Topology structure with masses.
106 * \param[in] f Forces on all atoms.
107 * \param[in] nrefat Number of atoms in the index.
108 * \param[in] index Indices of atoms.
109 * \param[out] fout Force on the COG position for the indexed atoms.
110 * \returns 0 on success, EINVAL if \p top is NULL.
113 gmx_calc_cog_f(t_topology *top, rvec f[], int nrefat, atom_id index[], rvec fout)
120 gmx_incons("no masses available while mass weighting was needed");
125 for (m = 0; m < nrefat; ++m)
128 mass = top->atoms.atom[ai].m;
129 for (j = 0; j < DIM; ++j)
131 fout[j] += f[ai][j] / mass;
135 svmul(mtot / nrefat, fout, fout);
140 gmx_calc_com_f(t_topology * /* top */, rvec f[], int nrefat, atom_id index[], rvec fout)
143 for (int m = 0; m < nrefat; ++m)
145 const int ai = index[m];
146 rvec_inc(fout, f[ai]);
152 * \param[in] top Topology structure with masses
153 * (can be NULL if \p bMASS==false).
154 * \param[in] x Position vectors of all atoms.
155 * \param[in] nrefat Number of atoms in the index.
156 * \param[in] index Indices of atoms.
157 * \param[in] bMass If true, mass weighting is used.
158 * \param[out] xout COM/COG position for the indexed atoms.
159 * \returns 0 on success, EINVAL if \p top is NULL and \p bMass is true.
161 * Calls either gmx_calc_com() or gmx_calc_cog() depending on the value of
163 * Other parameters are passed unmodified to these functions.
166 gmx_calc_comg(t_topology *top, rvec x[], int nrefat, atom_id index[],
167 bool bMass, rvec xout)
171 return gmx_calc_com(top, x, nrefat, index, xout);
175 return gmx_calc_cog(top, x, nrefat, index, xout);
180 * \param[in] top Topology structure with masses
181 * (can be NULL if \p bMASS==true).
182 * \param[in] f Forces on all atoms.
183 * \param[in] nrefat Number of atoms in the index.
184 * \param[in] index Indices of atoms.
185 * \param[in] bMass If true, force on COM is calculated.
186 * \param[out] fout Force on the COM/COG position for the indexed atoms.
187 * \returns 0 on success, EINVAL if \p top is NULL and \p bMass is false.
189 * Calls either gmx_calc_cog_f() or gmx_calc_com_f() depending on the value of
191 * Other parameters are passed unmodified to these functions.
194 gmx_calc_comg_f(t_topology *top, rvec f[], int nrefat, atom_id index[],
195 bool bMass, rvec fout)
199 return gmx_calc_com_f(top, f, nrefat, index, fout);
203 return gmx_calc_cog_f(top, f, nrefat, index, fout);
209 * \param[in] top Topology structure (unused, can be NULL).
210 * \param[in] x Position vectors of all atoms.
211 * \param[in] pbc Periodic boundary conditions structure.
212 * \param[in] nrefat Number of atoms in the index.
213 * \param[in] index Indices of atoms.
214 * \param[out] xout COG position for the indexed atoms.
215 * \returns 0 on success.
217 * Works exactly as gmx_calc_com_pbc(), but calculates the center of geometry.
220 gmx_calc_cog_pbc(t_topology *top, rvec x[], t_pbc *pbc,
221 int nrefat, atom_id index[], rvec xout)
223 const real tol = 1e-4;
228 /* First simple calculation */
229 gmx_calc_cog(top, x, nrefat, index, xout);
230 /* Now check if any atom is more than half the box from the COM */
237 for (m = 0; m < nrefat; ++m)
240 pbc_dx(pbc, x[ai], xout, dx);
241 rvec_add(xout, dx, xtest);
242 for (j = 0; j < DIM; ++j)
244 if (fabs(xtest[j] - x[ai][j]) > tol)
246 /* Here we have used the wrong image for contributing to the COM */
247 xout[j] += (xtest[j] - x[ai][j]) / nrefat;
261 * \param[in] top Topology structure with masses.
262 * \param[in] x Position vectors of all atoms.
263 * \param[in] pbc Periodic boundary conditions structure.
264 * \param[in] nrefat Number of atoms in the index.
265 * \param[in] index Indices of atoms.
266 * \param[out] xout COM position for the indexed atoms.
267 * \returns 0 on success, EINVAL if \p top is NULL.
269 * Works as gmx_calc_com(), but takes into account periodic boundary
270 * conditions: If any atom is more than half the box from the COM,
271 * it is wrapped around and a new COM is calculated. This is repeated
272 * until no atoms violate the condition.
274 * Modified from src/tools/gmx_sorient.c in Gromacs distribution.
277 gmx_calc_com_pbc(t_topology *top, rvec x[], t_pbc *pbc,
278 int nrefat, atom_id index[], rvec xout)
280 const real tol = 1e-4;
288 gmx_incons("no masses available while mass weighting was requested");
291 /* First simple calculation */
294 for (m = 0; m < nrefat; ++m)
297 mass = top->atoms.atom[ai].m;
298 for (j = 0; j < DIM; ++j)
300 xout[j] += mass * x[ai][j];
304 svmul(1.0/mtot, xout, xout);
305 /* Now check if any atom is more than half the box from the COM */
312 for (m = 0; m < nrefat; ++m)
315 mass = top->atoms.atom[ai].m / mtot;
316 pbc_dx(pbc, x[ai], xout, dx);
317 rvec_add(xout, dx, xtest);
318 for (j = 0; j < DIM; ++j)
320 if (fabs(xtest[j] - x[ai][j]) > tol)
322 /* Here we have used the wrong image for contributing to the COM */
323 xout[j] += mass * (xtest[j] - x[ai][j]);
337 * \param[in] top Topology structure with masses
338 * (can be NULL if \p bMASS==false).
339 * \param[in] x Position vectors of all atoms.
340 * \param[in] pbc Periodic boundary conditions structure.
341 * \param[in] nrefat Number of atoms in the index.
342 * \param[in] index Indices of atoms.
343 * \param[in] bMass If true, mass weighting is used.
344 * \param[out] xout COM/COG position for the indexed atoms.
345 * \returns 0 on success, EINVAL if \p top is NULL and \p bMass is true.
347 * Calls either gmx_calc_com() or gmx_calc_cog() depending on the value of
349 * Other parameters are passed unmodified to these functions.
352 gmx_calc_comg_pbc(t_topology *top, rvec x[], t_pbc *pbc,
353 int nrefat, atom_id index[], bool bMass, rvec xout)
357 return gmx_calc_com_pbc(top, x, pbc, nrefat, index, xout);
361 return gmx_calc_cog_pbc(top, x, pbc, nrefat, index, xout);
367 gmx_calc_cog_block(t_topology * /* top */, rvec x[], t_block *block, atom_id index[],
373 for (b = 0; b < block->nr; ++b)
376 for (i = block->index[b]; i < block->index[b+1]; ++i)
381 svmul(1.0/(block->index[b+1] - block->index[b]), xb, xout[b]);
387 * \param[in] top Topology structure with masses.
388 * \param[in] x Position vectors of all atoms.
389 * \param[in] block t_block structure that divides \p index into blocks.
390 * \param[in] index Indices of atoms.
391 * \param[out] xout \p block->nr COM positions.
392 * \returns 0 on success, EINVAL if \p top is NULL.
394 * Works exactly as gmx_calc_cog_block() with the exception that centers of
395 * mass are calculated, and hence a topology with masses is required.
398 gmx_calc_com_block(t_topology *top, rvec x[], t_block *block, atom_id index[],
407 gmx_incons("no masses available while mass weighting was requested");
410 for (b = 0; b < block->nr; ++b)
414 for (i = block->index[b]; i < block->index[b+1]; ++i)
417 mass = top->atoms.atom[ai].m;
418 for (d = 0; d < DIM; ++d)
420 xb[d] += mass * x[ai][d];
424 svmul(1.0/mtot, xb, xout[b]);
430 * \param[in] top Topology structure with masses.
431 * \param[in] f Forces on all atoms.
432 * \param[in] block t_block structure that divides \p index into blocks.
433 * \param[in] index Indices of atoms.
434 * \param[out] fout \p block->nr Forces on COG positions.
435 * \returns 0 on success, EINVAL if \p top is NULL.
438 gmx_calc_cog_f_block(t_topology *top, rvec f[], t_block *block, atom_id index[],
447 gmx_incons("no masses available while mass weighting was needed");
450 for (b = 0; b < block->nr; ++b)
454 for (i = block->index[b]; i < block->index[b+1]; ++i)
457 mass = top->atoms.atom[ai].m;
458 for (d = 0; d < DIM; ++d)
460 fb[d] += f[ai][d] / mass;
464 svmul(mtot / (block->index[b+1] - block->index[b]), fb, fout[b]);
470 gmx_calc_com_f_block(t_topology * /* top */, rvec f[], t_block *block, atom_id index[],
473 for (int b = 0; b < block->nr; ++b)
477 for (int i = block->index[b]; i < block->index[b+1]; ++i)
479 const int ai = index[i];
482 copy_rvec(fb, fout[b]);
488 * \param[in] top Topology structure with masses
489 * (can be NULL if \p bMASS==false).
490 * \param[in] x Position vectors of all atoms.
491 * \param[in] block t_block structure that divides \p index into blocks.
492 * \param[in] index Indices of atoms.
493 * \param[in] bMass If true, mass weighting is used.
494 * \param[out] xout \p block->nr COM/COG positions.
495 * \returns 0 on success, EINVAL if \p top is NULL and \p bMass is true.
497 * Calls either gmx_calc_com_block() or gmx_calc_cog_block() depending on the
499 * Other parameters are passed unmodified to these functions.
502 gmx_calc_comg_block(t_topology *top, rvec x[], t_block *block, atom_id index[],
503 bool bMass, rvec xout[])
507 return gmx_calc_com_block(top, x, block, index, xout);
511 return gmx_calc_cog_block(top, x, block, index, xout);
516 * \param[in] top Topology structure with masses
517 * (can be NULL if \p bMASS==false).
518 * \param[in] f Forces on all atoms.
519 * \param[in] block t_block structure that divides \p index into blocks.
520 * \param[in] index Indices of atoms.
521 * \param[in] bMass If true, force on COM is calculated.
522 * \param[out] fout \p block->nr forces on the COM/COG positions.
523 * \returns 0 on success, EINVAL if \p top is NULL and \p bMass is true.
525 * Calls either gmx_calc_com_f_block() or gmx_calc_cog_f_block() depending on
526 * the value of \p bMass.
527 * Other parameters are passed unmodified to these functions.
530 gmx_calc_comg_f_block(t_topology *top, rvec f[], t_block *block, atom_id index[],
531 bool bMass, rvec fout[])
535 return gmx_calc_com_f_block(top, f, block, index, fout);
539 return gmx_calc_cog_f_block(top, f, block, index, fout);
544 * \param[in] top Topology structure with masses
545 * (can be NULL if \p bMASS==false).
546 * \param[in] x Position vectors of all atoms.
547 * \param[in] block Blocks for calculation.
548 * \param[in] bMass If true, mass weighting is used.
549 * \param[out] xout \p block->nr COM/COG positions.
550 * \returns 0 on success, EINVAL if \p top is NULL and \p bMass is true.
552 * Calls gmx_calc_comg_block(), converting the \p t_blocka structure into
553 * a \p t_block and an index. Other parameters are passed unmodified.
556 * This function assumes that a pointer to \c t_blocka can be safely typecast
557 * into \c t_block such that the index fields can still be referenced.
558 * With the present Gromacs defitions of these types, this is the case,
559 * but if the layout of these structures is changed, this may lead to strange
563 gmx_calc_comg_blocka(t_topology *top, rvec x[], t_blocka *block,
564 bool bMass, rvec xout[])
566 /* TODO: It would probably be better to do this without the type cast */
567 return gmx_calc_comg_block(top, x, (t_block *)block, block->a, bMass, xout);
571 * \param[in] top Topology structure with masses
572 * (can be NULL if \p bMASS==true).
573 * \param[in] f Forces on all atoms.
574 * \param[in] block Blocks for calculation.
575 * \param[in] bMass If true, force on COM is calculated.
576 * \param[out] fout \p block->nr forces on the COM/COG positions.
577 * \returns 0 on success, EINVAL if \p top is NULL and \p bMass is false.
579 * Calls gmx_calc_comg_f_block(), converting the \p t_blocka structure into
580 * a \p t_block and an index. Other parameters are passed unmodified.
583 * This function assumes that a pointer to \c t_blocka can be safely typecast
584 * into \c t_block such that the index fields can still be referenced.
585 * With the present Gromacs defitions of these types, this is the case,
586 * but if the layout of these structures is changed, this may lead to strange
590 gmx_calc_comg_f_blocka(t_topology *top, rvec f[], t_blocka *block,
591 bool bMass, rvec fout[])
593 /* TODO: It would probably be better to do this without the type cast */
594 return gmx_calc_comg_f_block(top, f, (t_block *)block, block->a, bMass, fout);