Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / selection / sm_distance.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2016,2017,2018,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
37  * Implements distance-based selection methods.
38  *
39  * This file implements the \p distance, \p mindistance and \p within
40  * selection methods.
41  *
42  * \author Teemu Murtola <teemu.murtola@gmail.com>
43  * \ingroup module_selection
44  */
45 #include "gmxpre.h"
46
47 #include "gromacs/math/vec.h"
48 #include "gromacs/selection/nbsearch.h"
49 #include "gromacs/utility/arraysize.h"
50 #include "gromacs/utility/exceptions.h"
51
52 #include "position.h"
53 #include "selmethod.h"
54 #include "selmethod_impl.h"
55
56 /*! \internal
57  * \brief
58  * Data structure for distance-based selection method.
59  *
60  * The same data structure is used by all the distance-based methods.
61  *
62  * \ingroup module_selection
63  */
64 struct t_methoddata_distance
65 {
66     t_methoddata_distance() : cutoff(-1.0) {}
67
68     /** Cutoff distance. */
69     real cutoff;
70     /** Positions of the reference points. */
71     gmx_ana_pos_t p;
72     /** Neighborhood search data. */
73     gmx::AnalysisNeighborhood nb;
74     /** Neighborhood search for an invididual frame. */
75     gmx::AnalysisNeighborhoodSearch nbsearch;
76 };
77
78 /*! \brief
79  * Allocates data for distance-based selection methods.
80  *
81  * \param[in]     npar  Not used (should be 2).
82  * \param[in,out] param Method parameters (should point to one of the distance
83  *   parameter arrays).
84  * \returns       Pointer to the allocated data (\c t_methoddata_distance).
85  *
86  * Allocates memory for a \c t_methoddata_distance structure and
87  * initializes the parameter as follows:
88  *  - the first parameter defines the value for
89  *    \c t_methoddata_distance::cutoff.
90  *  - the second parameter defines the reference positions and the value is
91  *    stored in \c t_methoddata_distance::p.
92  */
93 static void* init_data_common(int npar, gmx_ana_selparam_t* param);
94 /*! \brief
95  * Initializes a distance-based selection method.
96  *
97  * \param   top   Not used.
98  * \param   npar  Not used (should be 2).
99  * \param   param Method parameters (should point to one of the distance
100  *   parameter arrays).
101  * \param   data  Pointer to \c t_methoddata_distance to initialize.
102  * \returns 0 on success, a non-zero error code on failure.
103  *
104  * Initializes the neighborhood search data structure
105  * (\c t_methoddata_distance::nb).
106  * Also checks that the cutoff is valid.
107  */
108 static void init_common(const gmx_mtop_t* top, int npar, gmx_ana_selparam_t* param, void* data);
109 /** Frees the data allocated for a distance-based selection method. */
110 static void free_data_common(void* data);
111 /*! \brief
112  * Initializes the evaluation of a distance-based within selection method for a
113  * frame.
114  *
115  * \param[in]  context Evaluation context.
116  * \param      data    Should point to a \c t_methoddata_distance.
117  * \returns    0 on success, a non-zero error code on error.
118  *
119  * Initializes the neighborhood search for the current frame.
120  */
121 static void init_frame_common(const gmx::SelMethodEvalContext& context, void* data);
122 /** Evaluates the \p distance selection method. */
123 static void evaluate_distance(const gmx::SelMethodEvalContext& /*context*/,
124                               gmx_ana_pos_t*      pos,
125                               gmx_ana_selvalue_t* out,
126                               void*               data);
127 /** Evaluates the \p within selection method. */
128 static void evaluate_within(const gmx::SelMethodEvalContext& /*context*/,
129                             gmx_ana_pos_t*      pos,
130                             gmx_ana_selvalue_t* out,
131                             void*               data);
132
133 /** Parameters for the \p distance selection method. */
134 static gmx_ana_selparam_t smparams_distance[] = {
135     { "cutoff", { REAL_VALUE, 1, { nullptr } }, nullptr, SPAR_OPTIONAL },
136     { "from", { POS_VALUE, 1, { nullptr } }, nullptr, SPAR_DYNAMIC },
137 };
138
139 /** Parameters for the \p mindistance selection method. */
140 static gmx_ana_selparam_t smparams_mindistance[] = {
141     { "cutoff", { REAL_VALUE, 1, { nullptr } }, nullptr, SPAR_OPTIONAL },
142     { "from", { POS_VALUE, -1, { nullptr } }, nullptr, SPAR_DYNAMIC | SPAR_VARNUM },
143 };
144
145 /** Parameters for the \p within selection method. */
146 static gmx_ana_selparam_t smparams_within[] = {
147     { nullptr, { REAL_VALUE, 1, { nullptr } }, nullptr, 0 },
148     { "of", { POS_VALUE, -1, { nullptr } }, nullptr, SPAR_DYNAMIC | SPAR_VARNUM },
149 };
150
151 //! Help title for distance selection methods.
152 static const char helptitle_distance[] = "Selecting based on distance";
153 //! Help text for distance selection methods.
154 static const char* const help_distance[] = {
155     "::",
156     "",
157     "  distance from POS [cutoff REAL]",
158     "  mindistance from POS_EXPR [cutoff REAL]",
159     "  within REAL of POS_EXPR",
160     "",
161     "[TT]distance[tt] and [TT]mindistance[tt] calculate the distance from the",
162     "given position(s), the only difference being in that [TT]distance[tt]",
163     "only accepts a single position, while any number of positions can be",
164     "given for [TT]mindistance[tt], which then calculates the distance to the",
165     "closest position.",
166     "[TT]within[tt] directly selects atoms that are within [TT]REAL[tt] of",
167     "[TT]POS_EXPR[tt].[PAR]",
168
169     "For the first two keywords, it is possible to specify a cutoff to speed",
170     "up the evaluation: all distances above the specified cutoff are",
171     "returned as equal to the cutoff.",
172 };
173
174 /** Selection method data for the \p distance method. */
175 gmx_ana_selmethod_t sm_distance = {
176     "distance",
177     REAL_VALUE,
178     SMETH_DYNAMIC,
179     asize(smparams_distance),
180     smparams_distance,
181     &init_data_common,
182     nullptr,
183     &init_common,
184     nullptr,
185     &free_data_common,
186     &init_frame_common,
187     nullptr,
188     &evaluate_distance,
189     { "distance from POS [cutoff REAL]", helptitle_distance, asize(help_distance), help_distance },
190 };
191
192 /** Selection method data for the \p distance method. */
193 gmx_ana_selmethod_t sm_mindistance = {
194     "mindistance",
195     REAL_VALUE,
196     SMETH_DYNAMIC,
197     asize(smparams_mindistance),
198     smparams_mindistance,
199     &init_data_common,
200     nullptr,
201     &init_common,
202     nullptr,
203     &free_data_common,
204     &init_frame_common,
205     nullptr,
206     &evaluate_distance,
207     { "mindistance from POS_EXPR [cutoff REAL]", helptitle_distance, asize(help_distance), help_distance },
208 };
209
210 /** Selection method data for the \p within method. */
211 gmx_ana_selmethod_t sm_within = {
212     "within",
213     GROUP_VALUE,
214     SMETH_DYNAMIC,
215     asize(smparams_within),
216     smparams_within,
217     &init_data_common,
218     nullptr,
219     &init_common,
220     nullptr,
221     &free_data_common,
222     &init_frame_common,
223     nullptr,
224     &evaluate_within,
225     { "within REAL of POS_EXPR", helptitle_distance, asize(help_distance), help_distance },
226 };
227
228 static void* init_data_common(int /* npar */, gmx_ana_selparam_t* param)
229 {
230     t_methoddata_distance* data = new t_methoddata_distance();
231     param[0].val.u.r            = &data->cutoff;
232     param[1].val.u.p            = &data->p;
233     return data;
234 }
235
236 static void init_common(const gmx_mtop_t* /* top */, int /* npar */, gmx_ana_selparam_t* param, void* data)
237 {
238     t_methoddata_distance* d = static_cast<t_methoddata_distance*>(data);
239
240     if ((param[0].flags & SPAR_SET) && d->cutoff <= 0)
241     {
242         GMX_THROW(gmx::InvalidInputError("Distance cutoff should be > 0"));
243     }
244     d->nb.setCutoff(d->cutoff);
245 }
246
247 /*!
248  * \param data Data to free (should point to a \c t_methoddata_distance).
249  *
250  * Frees the memory allocated for \c t_methoddata_distance::xref and
251  * \c t_methoddata_distance::nb.
252  */
253 static void free_data_common(void* data)
254 {
255     delete static_cast<t_methoddata_distance*>(data);
256 }
257
258 static void init_frame_common(const gmx::SelMethodEvalContext& context, void* data)
259 {
260     t_methoddata_distance* d = static_cast<t_methoddata_distance*>(data);
261
262     d->nbsearch.reset();
263     gmx::AnalysisNeighborhoodPositions pos(d->p.x, d->p.count());
264     d->nbsearch = d->nb.initSearch(context.pbc, pos);
265 }
266
267 /*!
268  * See sel_updatefunc_pos() for description of the parameters.
269  * \p data should point to a \c t_methoddata_distance.
270  *
271  * Calculates the distance of each position from \c t_methoddata_distance::p
272  * and puts them in \p out->u.r.
273  */
274 static void evaluate_distance(const gmx::SelMethodEvalContext& /*context*/,
275                               gmx_ana_pos_t*      pos,
276                               gmx_ana_selvalue_t* out,
277                               void*               data)
278 {
279     t_methoddata_distance* d = static_cast<t_methoddata_distance*>(data);
280
281     out->nr = pos->count();
282     for (int i = 0; i < pos->count(); ++i)
283     {
284         out->u.r[i] = d->nbsearch.minimumDistance(pos->x[i]);
285     }
286 }
287
288 /*!
289  * See sel_updatefunc() for description of the parameters.
290  * \p data should point to a \c t_methoddata_distance.
291  *
292  * Finds the atoms that are closer than the defined cutoff to
293  * \c t_methoddata_distance::xref and puts them in \p out.g.
294  */
295 static void evaluate_within(const gmx::SelMethodEvalContext& /*context*/,
296                             gmx_ana_pos_t*      pos,
297                             gmx_ana_selvalue_t* out,
298                             void*               data)
299 {
300     t_methoddata_distance* d = static_cast<t_methoddata_distance*>(data);
301
302     out->u.g->isize = 0;
303     for (int b = 0; b < pos->count(); ++b)
304     {
305         if (d->nbsearch.isWithin(pos->x[b]))
306         {
307             gmx_ana_pos_add_to_group(out->u.g, pos, b);
308         }
309     }
310 }