Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / selection / centerofmass.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009,2010,2011,2013,2014,2015,2016,2019, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief API for calculation of centers of mass/geometry.
37  *
38  * This header defines a few functions that can be used to calculate
39  * centers of mass/geometry for a group of atoms.
40  * These routines can be used independently of the other parts of the
41  * library, but they are also used internally by the selection engine.
42  * In most cases, it should not be necessary to call these functions
43  * directly.
44  * Instead, one should write an analysis tool such that it gets all
45  * positions through selections.
46  *
47  * The functions in the header can be divided into a few groups based on the
48  * parameters they take. The simplest group of functions calculates the center
49  * of a single group of atoms:
50  *  - gmx_calc_cog(): Calculates the center of geometry (COG) of a given
51  *    group of atoms.
52  *  - gmx_calc_com(): Calculates the center of mass (COM) of a given group
53  *    of atoms.
54  *  - gmx_calc_comg(): Calculates either the COM or COG, based on a
55  *    boolean flag.
56  *
57  * A second set of routines is provided for calculating the centers for groups
58  * that wrap over periodic boundaries (gmx_calc_cog_pbc(), gmx_calc_com_pbc(),
59  * gmx_calc_comg_pbc()). These functions are slower, because they need to
60  * adjust the center iteratively.
61  *
62  * It is also possible to calculate centers for several groups of atoms in
63  * one call. The functions gmx_calc_cog_block(), gmx_calc_com_block() and
64  * gmx_calc_comg_block() take an index group and a partitioning of that index
65  * group (as a \c t_block structure), and calculate the centers for
66  * each group defined by the \c t_block structure separately.
67  *
68  * Finally, there is a function gmx_calc_comg_blocka() that takes both the
69  * index group and the partitioning as a single \c t_blocka structure.
70  *
71  * \author Teemu Murtola <teemu.murtola@gmail.com>
72  * \ingroup module_selection
73  */
74 #ifndef GMX_SELECTION_CENTEROFMASS_H
75 #define GMX_SELECTION_CENTEROFMASS_H
76
77 #include "gromacs/math/vectypes.h"
78
79 struct gmx_mtop_t;
80 struct t_block;
81 struct t_blocka;
82 struct t_pbc;
83
84 /*! \brief
85  * Calculate a single center of geometry.
86  *
87  * \param[in]  top    Topology structure (unused, can be NULL).
88  * \param[in]  x      Position vectors of all atoms.
89  * \param[in]  nrefat Number of atoms in the index.
90  * \param[in]  index  Indices of atoms.
91  * \param[out] xout   COG position for the indexed atoms.
92  */
93 void gmx_calc_cog(const gmx_mtop_t* top, rvec x[], int nrefat, const int index[], rvec xout);
94 /** Calculate a single center of mass. */
95 void gmx_calc_com(const gmx_mtop_t* top, rvec x[], int nrefat, const int index[], rvec xout);
96 /** Calculate force on a single center of geometry. */
97 void gmx_calc_cog_f(const gmx_mtop_t* top, rvec f[], int nrefat, const int index[], rvec fout);
98 /*! \brief
99  * Calculate force on a single center of mass.
100  *
101  * \param[in]  top    Topology structure (unused, can be NULL).
102  * \param[in]  f      Forces on all atoms.
103  * \param[in]  nrefat Number of atoms in the index.
104  * \param[in]  index  Indices of atoms.
105  * \param[out] fout   Force on the COM position for the indexed atoms.
106  */
107 void gmx_calc_com_f(const gmx_mtop_t* top, rvec f[], int nrefat, const int index[], rvec fout);
108 /** Calculate a single center of mass/geometry. */
109 void gmx_calc_comg(const gmx_mtop_t* top, rvec x[], int nrefat, const int index[], bool bMass, rvec xout);
110 /** Calculate force on a single center of mass/geometry. */
111 void gmx_calc_comg_f(const gmx_mtop_t* top, rvec f[], int nrefat, const int index[], bool bMass, rvec fout);
112
113 /** Calculate a single center of geometry iteratively, taking PBC into account. */
114 void gmx_calc_cog_pbc(const gmx_mtop_t* top, rvec x[], const t_pbc* pbc, int nrefat, const int index[], rvec xout);
115 /** Calculate a single center of mass iteratively, taking PBC into account. */
116 void gmx_calc_com_pbc(const gmx_mtop_t* top, rvec x[], const t_pbc* pbc, int nrefat, const int index[], rvec xout);
117 /** Calculate a single center of mass/geometry iteratively with PBC. */
118 void gmx_calc_comg_pbc(const gmx_mtop_t* top,
119                        rvec              x[],
120                        const t_pbc*      pbc,
121                        int               nrefat,
122                        const int         index[],
123                        bool              bMass,
124                        rvec              xout);
125
126 /*! \brief
127  * Calculate centers of geometry for a blocked index.
128  *
129  * \param[in]  top   Topology structure (unused, can be NULL).
130  * \param[in]  x     Position vectors of all atoms.
131  * \param[in]  block t_block structure that divides \p index into blocks.
132  * \param[in]  index Indices of atoms.
133  * \param[out] xout  \p block->nr COG positions.
134  */
135 void gmx_calc_cog_block(const gmx_mtop_t* top, rvec x[], const t_block* block, const int index[], rvec xout[]);
136 /** Calculate centers of mass for a blocked index. */
137 void gmx_calc_com_block(const gmx_mtop_t* top, rvec x[], const t_block* block, const int index[], rvec xout[]);
138 /** Calculate forces on centers of geometry for a blocked index. */
139 void gmx_calc_cog_f_block(const gmx_mtop_t* top, rvec f[], const t_block* block, const int index[], rvec fout[]);
140 /*! \brief
141  * Calculate forces on centers of mass for a blocked index.
142  *
143  * \param[in]  top   Topology structure (unused, can be NULL).
144  * \param[in]  f     Forces on all atoms.
145  * \param[in]  block t_block structure that divides \p index into blocks.
146  * \param[in]  index Indices of atoms.
147  * \param[out] fout  \p block->nr Forces on COM positions.
148  */
149 void gmx_calc_com_f_block(const gmx_mtop_t* top, rvec f[], const t_block* block, const int index[], rvec fout[]);
150 /** Calculate centers of mass/geometry for a blocked index. */
151 void gmx_calc_comg_block(const gmx_mtop_t* top,
152                          rvec              x[],
153                          const t_block*    block,
154                          const int         index[],
155                          bool              bMass,
156                          rvec              xout[]);
157 /** Calculate forces on centers of mass/geometry for a blocked index. */
158 void gmx_calc_comg_f_block(const gmx_mtop_t* top,
159                            rvec              f[],
160                            const t_block*    block,
161                            const int         index[],
162                            bool              bMass,
163                            rvec              fout[]);
164 /** Calculate centers of mass/geometry for a set of blocks; */
165 void gmx_calc_comg_blocka(const gmx_mtop_t* top, rvec x[], const t_blocka* block, bool bMass, rvec xout[]);
166 /** Calculate forces on centers of mass/geometry for a set of blocks; */
167 void gmx_calc_comg_f_blocka(const gmx_mtop_t* top, rvec x[], const t_blocka* block, bool bMass, rvec xout[]);
168
169 #endif