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