Split lines with many copyright years
[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 by the GROMACS development team.
5  * Copyright (c) 2015,2016,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \internal \file
37  * \brief API for calculation of centers of mass/geometry.
38  *
39  * This header defines a few functions that can be used to calculate
40  * centers of mass/geometry for a group of atoms.
41  * These routines can be used independently of the other parts of the
42  * library, but they are also used internally by the selection engine.
43  * In most cases, it should not be necessary to call these functions
44  * directly.
45  * Instead, one should write an analysis tool such that it gets all
46  * positions through selections.
47  *
48  * The functions in the header can be divided into a few groups based on the
49  * parameters they take. The simplest group of functions calculates the center
50  * of a single group of atoms:
51  *  - gmx_calc_cog(): Calculates the center of geometry (COG) of a given
52  *    group of atoms.
53  *  - gmx_calc_com(): Calculates the center of mass (COM) of a given group
54  *    of atoms.
55  *  - gmx_calc_comg(): Calculates either the COM or COG, based on a
56  *    boolean flag.
57  *
58  * A second set of routines is provided for calculating the centers for groups
59  * that wrap over periodic boundaries (gmx_calc_cog_pbc(), gmx_calc_com_pbc(),
60  * gmx_calc_comg_pbc()). These functions are slower, because they need to
61  * adjust the center iteratively.
62  *
63  * It is also possible to calculate centers for several groups of atoms in
64  * one call. The functions gmx_calc_cog_block(), gmx_calc_com_block() and
65  * gmx_calc_comg_block() take an index group and a partitioning of that index
66  * group (as a \c t_block structure), and calculate the centers for
67  * each group defined by the \c t_block structure separately.
68  *
69  * Finally, there is a function gmx_calc_comg_blocka() that takes both the
70  * index group and the partitioning as a single \c t_blocka structure.
71  *
72  * \author Teemu Murtola <teemu.murtola@gmail.com>
73  * \ingroup module_selection
74  */
75 #ifndef GMX_SELECTION_CENTEROFMASS_H
76 #define GMX_SELECTION_CENTEROFMASS_H
77
78 #include "gromacs/math/vectypes.h"
79
80 struct gmx_mtop_t;
81 struct t_block;
82 struct t_blocka;
83 struct t_pbc;
84
85 /*! \brief
86  * Calculate a single center of geometry.
87  *
88  * \param[in]  top    Topology structure (unused, can be NULL).
89  * \param[in]  x      Position vectors of all atoms.
90  * \param[in]  nrefat Number of atoms in the index.
91  * \param[in]  index  Indices of atoms.
92  * \param[out] xout   COG position for the indexed atoms.
93  */
94 void gmx_calc_cog(const gmx_mtop_t* top, rvec x[], int nrefat, const int index[], rvec xout);
95 /** Calculate a single center of mass. */
96 void gmx_calc_com(const gmx_mtop_t* top, rvec x[], int nrefat, const int index[], rvec xout);
97 /** Calculate force on a single center of geometry. */
98 void gmx_calc_cog_f(const gmx_mtop_t* top, rvec f[], int nrefat, const int index[], rvec fout);
99 /*! \brief
100  * Calculate force on a single center of mass.
101  *
102  * \param[in]  top    Topology structure (unused, can be NULL).
103  * \param[in]  f      Forces on all atoms.
104  * \param[in]  nrefat Number of atoms in the index.
105  * \param[in]  index  Indices of atoms.
106  * \param[out] fout   Force on the COM position for the indexed atoms.
107  */
108 void gmx_calc_com_f(const gmx_mtop_t* top, rvec f[], int nrefat, const int index[], rvec fout);
109 /** Calculate a single center of mass/geometry. */
110 void gmx_calc_comg(const gmx_mtop_t* top, rvec x[], int nrefat, const int index[], bool bMass, rvec xout);
111 /** Calculate force on a single center of mass/geometry. */
112 void gmx_calc_comg_f(const gmx_mtop_t* top, rvec f[], int nrefat, const int index[], bool bMass, rvec fout);
113
114 /** Calculate a single center of geometry iteratively, taking PBC into account. */
115 void gmx_calc_cog_pbc(const gmx_mtop_t* top, rvec x[], const t_pbc* pbc, int nrefat, const int index[], rvec xout);
116 /** Calculate a single center of mass iteratively, taking PBC into account. */
117 void gmx_calc_com_pbc(const gmx_mtop_t* top, rvec x[], const t_pbc* pbc, int nrefat, const int index[], rvec xout);
118 /** Calculate a single center of mass/geometry iteratively with PBC. */
119 void gmx_calc_comg_pbc(const gmx_mtop_t* top,
120                        rvec              x[],
121                        const t_pbc*      pbc,
122                        int               nrefat,
123                        const int         index[],
124                        bool              bMass,
125                        rvec              xout);
126
127 /*! \brief
128  * Calculate centers of geometry for a blocked index.
129  *
130  * \param[in]  top   Topology structure (unused, can be NULL).
131  * \param[in]  x     Position vectors of all atoms.
132  * \param[in]  block t_block structure that divides \p index into blocks.
133  * \param[in]  index Indices of atoms.
134  * \param[out] xout  \p block->nr COG positions.
135  */
136 void gmx_calc_cog_block(const gmx_mtop_t* top, rvec x[], const t_block* block, const int index[], rvec xout[]);
137 /** Calculate centers of mass for a blocked index. */
138 void gmx_calc_com_block(const gmx_mtop_t* top, rvec x[], const t_block* block, const int index[], rvec xout[]);
139 /** Calculate forces on centers of geometry for a blocked index. */
140 void gmx_calc_cog_f_block(const gmx_mtop_t* top, rvec f[], const t_block* block, const int index[], rvec fout[]);
141 /*! \brief
142  * Calculate forces on centers of mass for a blocked index.
143  *
144  * \param[in]  top   Topology structure (unused, can be NULL).
145  * \param[in]  f     Forces on all atoms.
146  * \param[in]  block t_block structure that divides \p index into blocks.
147  * \param[in]  index Indices of atoms.
148  * \param[out] fout  \p block->nr Forces on COM positions.
149  */
150 void gmx_calc_com_f_block(const gmx_mtop_t* top, rvec f[], const t_block* block, const int index[], rvec fout[]);
151 /** Calculate centers of mass/geometry for a blocked index. */
152 void gmx_calc_comg_block(const gmx_mtop_t* top,
153                          rvec              x[],
154                          const t_block*    block,
155                          const int         index[],
156                          bool              bMass,
157                          rvec              xout[]);
158 /** Calculate forces on centers of mass/geometry for a blocked index. */
159 void gmx_calc_comg_f_block(const gmx_mtop_t* top,
160                            rvec              f[],
161                            const t_block*    block,
162                            const int         index[],
163                            bool              bMass,
164                            rvec              fout[]);
165 /** Calculate centers of mass/geometry for a set of blocks; */
166 void gmx_calc_comg_blocka(const gmx_mtop_t* top, rvec x[], const t_blocka* block, bool bMass, rvec xout[]);
167 /** Calculate forces on centers of mass/geometry for a set of blocks; */
168 void gmx_calc_comg_f_blocka(const gmx_mtop_t* top, rvec x[], const t_blocka* block, bool bMass, rvec xout[]);
169
170 #endif