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