Replace some [BR] in help texts with rst literals
[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, 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/legacyheaders/macros.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/selection/nbsearch.h"
50 #include "gromacs/selection/position.h"
51 #include "gromacs/utility/exceptions.h"
52
53 #include "selmethod.h"
54
55 /*! \internal
56  * \brief
57  * Data structure for distance-based selection method.
58  *
59  * The same data structure is used by all the distance-based methods.
60  *
61  * \ingroup module_selection
62  */
63 struct t_methoddata_distance
64 {
65     t_methoddata_distance() : cutoff(-1.0)
66     {
67     }
68
69     /** Cutoff distance. */
70     real                             cutoff;
71     /** Positions of the reference points. */
72     gmx_ana_pos_t                    p;
73     /** Neighborhood search data. */
74     gmx::AnalysisNeighborhood        nb;
75     /** Neighborhood search for an invididual frame. */
76     gmx::AnalysisNeighborhoodSearch  nbsearch;
77 };
78
79 /*! \brief
80  * Allocates data for distance-based selection methods.
81  *
82  * \param[in]     npar  Not used (should be 2).
83  * \param[in,out] param Method parameters (should point to one of the distance
84  *   parameter arrays).
85  * \returns       Pointer to the allocated data (\c t_methoddata_distance).
86  *
87  * Allocates memory for a \c t_methoddata_distance structure and
88  * initializes the parameter as follows:
89  *  - the first parameter defines the value for
90  *    \c t_methoddata_distance::cutoff.
91  *  - the second parameter defines the reference positions and the value is
92  *    stored in \c t_methoddata_distance::p.
93  */
94 static void *
95 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
111 init_common(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
112 /** Frees the data allocated for a distance-based selection method. */
113 static void
114 free_data_common(void *data);
115 /*! \brief
116  * Initializes the evaluation of a distance-based within selection method for a
117  * frame.
118  *
119  * \param[in]  top  Not used.
120  * \param[in]  fr   Current frame.
121  * \param[in]  pbc  PBC structure.
122  * \param      data Should point to a \c t_methoddata_distance.
123  * \returns    0 on success, a non-zero error code on error.
124  *
125  * Initializes the neighborhood search for the current frame.
126  */
127 static void
128 init_frame_common(t_topology *top, t_trxframe * fr, t_pbc *pbc, void *data);
129 /** Evaluates the \p distance selection method. */
130 static void
131 evaluate_distance(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
132                   gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
133 /** Evaluates the \p within selection method. */
134 static void
135 evaluate_within(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
136                 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
137
138 /** Parameters for the \p distance selection method. */
139 static gmx_ana_selparam_t smparams_distance[] = {
140     {"cutoff", {REAL_VALUE, 1, {NULL}}, NULL, SPAR_OPTIONAL},
141     {"from",   {POS_VALUE,  1, {NULL}}, NULL, SPAR_DYNAMIC},
142 };
143
144 /** Parameters for the \p mindistance selection method. */
145 static gmx_ana_selparam_t smparams_mindistance[] = {
146     {"cutoff", {REAL_VALUE, 1, {NULL}}, NULL, SPAR_OPTIONAL},
147     {"from",   {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
148 };
149
150 /** Parameters for the \p within selection method. */
151 static gmx_ana_selparam_t smparams_within[] = {
152     {NULL, {REAL_VALUE,  1, {NULL}}, NULL, 0},
153     {"of", {POS_VALUE,  -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
154 };
155
156 /** Help text for the distance selection methods. */
157 static const char *help_distance[] = {
158     "DISTANCE-BASED SELECTION KEYWORDS",
159     "",
160     "::",
161     "",
162     "  distance from POS [cutoff REAL]",
163     "  mindistance from POS_EXPR [cutoff REAL]",
164     "  within REAL of POS_EXPR",
165     "",
166     "[TT]distance[tt] and [TT]mindistance[tt] calculate the distance from the",
167     "given position(s), the only difference being in that [TT]distance[tt]",
168     "only accepts a single position, while any number of positions can be",
169     "given for [TT]mindistance[tt], which then calculates the distance to the",
170     "closest position.",
171     "[TT]within[tt] directly selects atoms that are within [TT]REAL[tt] of",
172     "[TT]POS_EXPR[tt].[PAR]",
173
174     "For the first two keywords, it is possible to specify a cutoff to speed",
175     "up the evaluation: all distances above the specified cutoff are",
176     "returned as equal to the cutoff.",
177 };
178
179 /** Selection method data for the \p distance method. */
180 gmx_ana_selmethod_t sm_distance = {
181     "distance", REAL_VALUE, SMETH_DYNAMIC,
182     asize(smparams_distance), smparams_distance,
183     &init_data_common,
184     NULL,
185     &init_common,
186     NULL,
187     &free_data_common,
188     &init_frame_common,
189     NULL,
190     &evaluate_distance,
191     {"distance from POS [cutoff REAL]", 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", REAL_VALUE, SMETH_DYNAMIC,
197     asize(smparams_mindistance), smparams_mindistance,
198     &init_data_common,
199     NULL,
200     &init_common,
201     NULL,
202     &free_data_common,
203     &init_frame_common,
204     NULL,
205     &evaluate_distance,
206     {"mindistance from POS_EXPR [cutoff REAL]", asize(help_distance), help_distance},
207 };
208
209 /** Selection method data for the \p within method. */
210 gmx_ana_selmethod_t sm_within = {
211     "within", GROUP_VALUE, SMETH_DYNAMIC,
212     asize(smparams_within), smparams_within,
213     &init_data_common,
214     NULL,
215     &init_common,
216     NULL,
217     &free_data_common,
218     &init_frame_common,
219     NULL,
220     &evaluate_within,
221     {"within REAL of POS_EXPR", asize(help_distance), help_distance},
222 };
223
224 static void *
225 init_data_common(int /* npar */, gmx_ana_selparam_t *param)
226 {
227     t_methoddata_distance *data = new t_methoddata_distance();
228     param[0].val.u.r = &data->cutoff;
229     param[1].val.u.p = &data->p;
230     return data;
231 }
232
233 static void
234 init_common(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
235 {
236     t_methoddata_distance *d = static_cast<t_methoddata_distance *>(data);
237
238     if ((param[0].flags & SPAR_SET) && d->cutoff <= 0)
239     {
240         GMX_THROW(gmx::InvalidInputError("Distance cutoff should be > 0"));
241     }
242     d->nb.setCutoff(d->cutoff);
243 }
244
245 /*!
246  * \param data Data to free (should point to a \c t_methoddata_distance).
247  *
248  * Frees the memory allocated for \c t_methoddata_distance::xref and
249  * \c t_methoddata_distance::nb.
250  */
251 static void
252 free_data_common(void *data)
253 {
254     delete static_cast<t_methoddata_distance *>(data);
255 }
256
257 static void
258 init_frame_common(t_topology * /* top */, t_trxframe * /* fr */, t_pbc *pbc, 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(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
275 evaluate_distance(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
276                   gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
277 {
278     t_methoddata_distance *d = static_cast<t_methoddata_distance *>(data);
279
280     out->nr = pos->count();
281     for (int i = 0; i < pos->count(); ++i)
282     {
283         out->u.r[i] = d->nbsearch.minimumDistance(pos->x[i]);
284     }
285 }
286
287 /*!
288  * See sel_updatefunc() for description of the parameters.
289  * \p data should point to a \c t_methoddata_distance.
290  *
291  * Finds the atoms that are closer than the defined cutoff to
292  * \c t_methoddata_distance::xref and puts them in \p out.g.
293  */
294 static void
295 evaluate_within(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
296                 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
297 {
298     t_methoddata_distance *d = static_cast<t_methoddata_distance *>(data);
299
300     out->u.g->isize = 0;
301     for (int b = 0; b < pos->count(); ++b)
302     {
303         if (d->nbsearch.isWithin(pos->x[b]))
304         {
305             gmx_ana_pos_add_to_group(out->u.g, pos, b);
306         }
307     }
308 }