Traj. analysis structure now passes oenv around.
[gromacs/adressmacs.git] / include / centerofmass.h
blobc64893d34168c3ded286d28b70a3cadb4f756fb7
1 /*
3 * This source code is part of
5 * G R O M A C S
7 * GROningen MAchine for Chemical Simulations
9 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11 * Copyright (c) 2001-2009, The GROMACS development team,
12 * check out http://www.gromacs.org for more information.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
31 /*! \file
32 * \brief API for calculation of centers of mass/geometry.
34 * This header defines a few functions that can be used to calculate
35 * centers of mass/geometry for a group of atoms.
36 * These routines can be used independently of the other parts of the
37 * library, but they are also used internally by the selection engine.
38 * In most cases, it should not be necessary to call these functions
39 * directly.
40 * Instead, one should write an analysis tool such that it gets all
41 * positions through selections.
43 * The functions in the header can be divided into a few groups based on the
44 * parameters they take. The simplest group of functions calculates the center
45 * of a single group of atoms:
46 * - gmx_calc_cog(): Calculates the center of geometry (COG) of a given
47 * group of atoms.
48 * - gmx_calc_com(): Calculates the center of mass (COM) of a given group
49 * of atoms.
50 * - gmx_calc_comg(): Calculates either the COM or COG, based on a
51 * boolean flag.
53 * A second set of routines is provided for calculating the centers for groups
54 * that wrap over periodic boundaries (gmx_calc_cog_pbc(), gmx_calc_com_pbc(),
55 * gmx_calc_comg_pbc()). These functions are slower, because they need to
56 * adjust the center iteratively.
58 * It is also possible to calculate centers for several groups of atoms in
59 * one call. The functions gmx_calc_cog_block(), gmx_calc_com_block() and
60 * gmx_calc_comg_block() take an index group and a partitioning of that index
61 * group (as a \c t_block structure), and calculate the centers for
62 * each group defined by the \c t_block structure separately.
64 * Finally, there is a function gmx_calc_comg_blocka() that takes both the
65 * index group and the partitioning as a single \c t_blocka structure.
67 #ifndef CENTEROFMASS_H
68 #define CENTEROFMASS_H
70 #include <typedefs.h>
72 #ifdef __cplusplus
73 extern "C"
75 #endif
77 //! Calculate a single center of geometry.
78 extern int
79 gmx_calc_cog(t_topology *top, rvec x[], int nrefat, atom_id index[], rvec xref);
80 //! Calculate a single center of mass.
81 extern int
82 gmx_calc_com(t_topology *top, rvec x[], int nrefat, atom_id index[], rvec xref);
83 //! Calculate a single center of mass/geometry.
84 extern int
85 gmx_calc_comg(t_topology *top, rvec x[], int nrefat, atom_id index[],
86 bool bMass, rvec xout);
88 //! Calculate a single center of geometry iteratively, taking PBC into account.
89 extern int
90 gmx_calc_cog_pbc(t_topology *top, rvec x[], t_pbc *pbc,
91 int nrefat, atom_id index[], rvec xout);
92 //! Calculate a single center of mass iteratively, taking PBC into account.
93 extern int
94 gmx_calc_com_pbc(t_topology *top, rvec x[], t_pbc *pbc,
95 int nrefat, atom_id index[], rvec xout);
96 //! Calculate a single center of mass/geometry iteratively with PBC.
97 extern int
98 gmx_calc_comg_pbc(t_topology *top, rvec x[], t_pbc *pbc,
99 int nrefat, atom_id index[], bool bMass, rvec xout);
101 //! Calculate centers of geometry for a blocked index.
102 extern int
103 gmx_calc_cog_block(t_topology *top, rvec x[], t_block *block,
104 atom_id index[], rvec xout[]);
105 //! Calculate centers of mass for a blocked index.
106 extern int
107 gmx_calc_com_block(t_topology *top, rvec x[], t_block *block,
108 atom_id index[], rvec xout[]);
109 //! Calculate centers of mass/geometry for a blocked index.
110 extern int
111 gmx_calc_comg_block(t_topology *top, rvec x[], t_block *block,
112 atom_id index[], bool bMass, rvec xout[]);
113 //! Calculate centers of mass/geometry for a set of blocks;
114 extern int
115 gmx_calc_comg_blocka(t_topology *top, rvec x[], t_blocka *block,
116 bool bMass, rvec xout[]);
118 #ifdef __cplusplus
120 #endif
122 #endif