Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / trajectoryanalysis / modules / surfacearea.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2017,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #ifndef GMX_TRAJECTORYANALYSIS_SURFACEAREA_H
38 #define GMX_TRAJECTORYANALYSIS_SURFACEAREA_H
39
40 #include "gromacs/math/vectypes.h"
41 #include "gromacs/utility/arrayref.h"
42 #include "gromacs/utility/classhelpers.h"
43 #include "gromacs/utility/real.h"
44
45 struct t_pbc;
46
47 #define FLAG_DOTS 01
48 #define FLAG_VOLUME 02
49 #define FLAG_ATOM_AREA 04
50
51 namespace gmx
52 {
53
54 /*! \internal
55  * \brief
56  * Computes surface areas for a group of atoms/spheres.
57  *
58  * This class provides a surface area/volume calculator.
59  *
60  * The algorithm is based on representing each atom/sphere surface as a set of
61  * dots, and determining which dots are on the surface (not covered by any
62  * other atom/sphere).  The dots are distributed evenly using an icosahedron- or
63  * a dodecahedron-based method (see the original reference cited in the code).
64  * The area is then estimated from the area represented by each dot.
65  * The volume is calculated by selecting a fixed point and integrating over the
66  * surface dots, summing up the cones whose tip is at the fixed point and base
67  * at the surface points.
68  *
69  * The default dot density per sphere is 32, which gives quite inaccurate
70  * areas and volumes, but a reasonable number of surface points.  According to
71  * original documentation of the method, a density of 600-700 dots gives an
72  * accuracy of 1.5 A^2 per atom.
73  *
74  * \ingroup module_trajectoryanalysis
75  */
76 class SurfaceAreaCalculator
77 {
78 public:
79     /*! \brief
80      * Initializes a surface area calculator.
81      *
82      * \throws std::bad_alloc if out of memory.
83      */
84     SurfaceAreaCalculator();
85     ~SurfaceAreaCalculator();
86
87     /*! \brief
88      * Sets the number of surface dots per sphere to use.
89      *
90      * This function must be called before calculate() to set the desired
91      * accuracy/computational cost.
92      */
93     void setDotCount(int dotCount);
94     /*! \brief
95      * Sets the radii of spheres to use in the calculation.
96      *
97      * \param[in]  radius  Radius for each atom/sphere.
98      *
99      * This function must be called before calculate() to set the radii for
100      * the spheres.  All calculations must use the same set of radii to
101      * share the same grid search.
102      * These radii are used as-is, without adding any probe radius.
103      * The passed array must remain valid for the lifetime of this object.
104      *
105      * Does not throw.
106      */
107     void setRadii(const ArrayRef<const real>& radius);
108
109     /*! \brief
110      * Requests calculation of volume.
111      *
112      * If not called, and FLAG_VOLUME is not passed to calculate(), the
113      * volume output is not produced.
114      *
115      * Does not throw.
116      */
117     void setCalculateVolume(bool bVolume);
118     /*! \brief
119      * Requests output of per-atom areas.
120      *
121      * If not called, and FLAG_ATOM_AREA is not passed to calculate(), the
122      * atom area output is not produced.
123      *
124      * Does not throw.
125      */
126     void setCalculateAtomArea(bool bAtomArea);
127     /*! \brief
128      * Requests output of all surface dots.
129      *
130      * If not called, and FLAG_DOTS is not passed to calculate(), the
131      * surface dot output is not produced.
132      *
133      * Does not throw.
134      */
135     void setCalculateSurfaceDots(bool bDots);
136
137     /*! \brief
138      * Calculates the surface area for a set of positions.
139      *
140      * \param[in]  x       Atom positions (sphere centers).
141      * \param[in]  pbc     PBC information (if `NULL`, calculation is done
142      *     without PBC).
143      * \param[in]  nat     Number of atoms to calculate.
144      * \param[in]  index   Atom indices to include in the calculation.
145      * \param[in]  flags   Additional flags for the calculation.
146      * \param[out] area    Total surface area (must be non-`NULL`).
147      * \param[out] volume  Total volume (can be `NULL`).
148      * \param[out] at_area Surface area for each atom in \p index
149      *     (\p nat values) (can be `NULL`).
150      * \param[out] lidots  Surface dots as x,y,z triplets (`3*lidots` values)
151      *     (can be `NULL`).
152      * \param[out] n_dots Number of surface dots in \p lidots
153      *     (can be `NULL`).
154      *
155      * Calculates the surface area of spheres centered at `x[index[0]]`,
156      * ..., `x[index[nat-1]]`, with radii `radii[index[0]]`, ..., where
157      * `radii` is the array passed to setRadii().
158      *
159      * If \p flags is 0, the calculation is done for the items specified
160      * with setCalculateVolume(), setCalculateAtomArea(), and
161      * setCalculateSurfaceDots().  Flags can specify FLAG_VOLUME,
162      * FLAG_ATOM_AREA, and/or FLAG_DOTS to request additional output for
163      * this particular calculation.  If any output is `NULL`, that output
164      * is not calculated, irrespective of the calculation mode set.
165      *
166      * \todo
167      * Make the output options more C++-like, in particular for the array
168      * outputs.
169      */
170     void calculate(const rvec*  x,
171                    const t_pbc* pbc,
172                    int          nat,
173                    int          index[],
174                    int          flags,
175                    real*        area,
176                    real*        volume,
177                    real**       at_area,
178                    real**       lidots,
179                    int*         n_dots) const;
180
181 private:
182     class Impl;
183
184     PrivateImplPointer<Impl> impl_;
185 };
186
187 } // namespace gmx
188
189 #endif