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