Remove PrivateImplPointer in favour of std::unique_ptr
[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,2021, 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 <memory>
41
42 #include "gromacs/math/vectypes.h"
43 #include "gromacs/utility/arrayref.h"
44 #include "gromacs/utility/real.h"
45
46 struct t_pbc;
47
48 #define FLAG_DOTS 01
49 #define FLAG_VOLUME 02
50 #define FLAG_ATOM_AREA 04
51
52 namespace gmx
53 {
54
55 /*! \internal
56  * \brief
57  * Computes surface areas for a group of atoms/spheres.
58  *
59  * This class provides a surface area/volume calculator.
60  *
61  * The algorithm is based on representing each atom/sphere surface as a set of
62  * dots, and determining which dots are on the surface (not covered by any
63  * other atom/sphere).  The dots are distributed evenly using an icosahedron- or
64  * a dodecahedron-based method (see the original reference cited in the code).
65  * The area is then estimated from the area represented by each dot.
66  * The volume is calculated by selecting a fixed point and integrating over the
67  * surface dots, summing up the cones whose tip is at the fixed point and base
68  * at the surface points.
69  *
70  * The default dot density per sphere is 32, which gives quite inaccurate
71  * areas and volumes, but a reasonable number of surface points.  According to
72  * original documentation of the method, a density of 600-700 dots gives an
73  * accuracy of 1.5 A^2 per atom.
74  *
75  * \ingroup module_trajectoryanalysis
76  */
77 class SurfaceAreaCalculator
78 {
79 public:
80     /*! \brief
81      * Initializes a surface area calculator.
82      *
83      * \throws std::bad_alloc if out of memory.
84      */
85     SurfaceAreaCalculator();
86     ~SurfaceAreaCalculator();
87
88     /*! \brief
89      * Sets the number of surface dots per sphere to use.
90      *
91      * This function must be called before calculate() to set the desired
92      * accuracy/computational cost.
93      */
94     void setDotCount(int dotCount);
95     /*! \brief
96      * Sets the radii of spheres to use in the calculation.
97      *
98      * \param[in]  radius  Radius for each atom/sphere.
99      *
100      * This function must be called before calculate() to set the radii for
101      * the spheres.  All calculations must use the same set of radii to
102      * share the same grid search.
103      * These radii are used as-is, without adding any probe radius.
104      * The passed array must remain valid for the lifetime of this object.
105      *
106      * Does not throw.
107      */
108     void setRadii(const ArrayRef<const real>& radius);
109
110     /*! \brief
111      * Requests calculation of volume.
112      *
113      * If not called, and FLAG_VOLUME is not passed to calculate(), the
114      * volume output is not produced.
115      *
116      * Does not throw.
117      */
118     void setCalculateVolume(bool bVolume);
119     /*! \brief
120      * Requests output of per-atom areas.
121      *
122      * If not called, and FLAG_ATOM_AREA is not passed to calculate(), the
123      * atom area output is not produced.
124      *
125      * Does not throw.
126      */
127     void setCalculateAtomArea(bool bAtomArea);
128     /*! \brief
129      * Requests output of all surface dots.
130      *
131      * If not called, and FLAG_DOTS is not passed to calculate(), the
132      * surface dot output is not produced.
133      *
134      * Does not throw.
135      */
136     void setCalculateSurfaceDots(bool bDots);
137
138     /*! \brief
139      * Calculates the surface area for a set of positions.
140      *
141      * \param[in]  x       Atom positions (sphere centers).
142      * \param[in]  pbc     PBC information (if `NULL`, calculation is done
143      *     without PBC).
144      * \param[in]  nat     Number of atoms to calculate.
145      * \param[in]  index   Atom indices to include in the calculation.
146      * \param[in]  flags   Additional flags for the calculation.
147      * \param[out] area    Total surface area (must be non-`NULL`).
148      * \param[out] volume  Total volume (can be `NULL`).
149      * \param[out] at_area Surface area for each atom in \p index
150      *     (\p nat values) (can be `NULL`).
151      * \param[out] lidots  Surface dots as x,y,z triplets (`3*lidots` values)
152      *     (can be `NULL`).
153      * \param[out] n_dots Number of surface dots in \p lidots
154      *     (can be `NULL`).
155      *
156      * Calculates the surface area of spheres centered at `x[index[0]]`,
157      * ..., `x[index[nat-1]]`, with radii `radii[index[0]]`, ..., where
158      * `radii` is the array passed to setRadii().
159      *
160      * If \p flags is 0, the calculation is done for the items specified
161      * with setCalculateVolume(), setCalculateAtomArea(), and
162      * setCalculateSurfaceDots().  Flags can specify FLAG_VOLUME,
163      * FLAG_ATOM_AREA, and/or FLAG_DOTS to request additional output for
164      * this particular calculation.  If any output is `NULL`, that output
165      * is not calculated, irrespective of the calculation mode set.
166      *
167      * \todo
168      * Make the output options more C++-like, in particular for the array
169      * outputs.
170      */
171     void calculate(const rvec*  x,
172                    const t_pbc* pbc,
173                    int          nat,
174                    int          index[],
175                    int          flags,
176                    real*        area,
177                    real*        volume,
178                    real**       at_area,
179                    real**       lidots,
180                    int*         n_dots) const;
181
182 private:
183     class Impl;
184
185     std::unique_ptr<Impl> impl_;
186 };
187
188 } // namespace gmx
189
190 #endif