2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 * Implements distance-based selection methods.
39 * This file implements the \p distance, \p mindistance and \p within
42 * \author Teemu Murtola <teemu.murtola@gmail.com>
43 * \ingroup module_selection
45 #include "gromacs/legacyheaders/macros.h"
46 #include "gromacs/legacyheaders/pbc.h"
47 #include "gromacs/legacyheaders/vec.h"
49 #include "gromacs/selection/nbsearch.h"
50 #include "gromacs/selection/position.h"
51 #include "gromacs/selection/selmethod.h"
52 #include "gromacs/utility/exceptions.h"
55 * Data structure for distance-based selection method.
57 * The same data structure is used by all the distance-based methods.
59 struct t_methoddata_distance
61 t_methoddata_distance() : cutoff(-1.0)
63 gmx_ana_pos_clear(&p);
66 /** Cutoff distance. */
68 /** Positions of the reference points. */
70 /** Neighborhood search data. */
71 gmx::AnalysisNeighborhood nb;
72 /** Neighborhood search for an invididual frame. */
73 gmx::AnalysisNeighborhoodSearch nbsearch;
76 /** Allocates data for distance-based selection methods. */
78 init_data_common(int npar, gmx_ana_selparam_t *param);
79 /** Initializes a distance-based selection method. */
81 init_common(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
82 /** Frees the data allocated for a distance-based selection method. */
84 free_data_common(void *data);
85 /** Initializes the evaluation of a distance-based within selection method for a frame. */
87 init_frame_common(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
88 /** Evaluates the \p distance selection method. */
90 evaluate_distance(t_topology *top, t_trxframe *fr, t_pbc *pbc,
91 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
92 /** Evaluates the \p within selection method. */
94 evaluate_within(t_topology *top, t_trxframe *fr, t_pbc *pbc,
95 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
97 /** Parameters for the \p distance selection method. */
98 static gmx_ana_selparam_t smparams_distance[] = {
99 {"cutoff", {REAL_VALUE, 1, {NULL}}, NULL, SPAR_OPTIONAL},
100 {"from", {POS_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
103 /** Parameters for the \p mindistance selection method. */
104 static gmx_ana_selparam_t smparams_mindistance[] = {
105 {"cutoff", {REAL_VALUE, 1, {NULL}}, NULL, SPAR_OPTIONAL},
106 {"from", {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
109 /** Parameters for the \p within selection method. */
110 static gmx_ana_selparam_t smparams_within[] = {
111 {NULL, {REAL_VALUE, 1, {NULL}}, NULL, 0},
112 {"of", {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
115 /** Help text for the distance selection methods. */
116 static const char *help_distance[] = {
117 "DISTANCE-BASED SELECTION KEYWORDS[PAR]",
119 "[TT]distance from POS [cutoff REAL][tt][BR]",
120 "[TT]mindistance from POS_EXPR [cutoff REAL][tt][BR]",
121 "[TT]within REAL of POS_EXPR[tt][PAR]",
123 "[TT]distance[tt] and [TT]mindistance[tt] calculate the distance from the",
124 "given position(s), the only difference being in that [TT]distance[tt]",
125 "only accepts a single position, while any number of positions can be",
126 "given for [TT]mindistance[tt], which then calculates the distance to the",
128 "[TT]within[tt] directly selects atoms that are within [TT]REAL[tt] of",
129 "[TT]POS_EXPR[tt].[PAR]",
131 "For the first two keywords, it is possible to specify a cutoff to speed",
132 "up the evaluation: all distances above the specified cutoff are",
133 "returned as equal to the cutoff.",
136 /** \internal Selection method data for the \p distance method. */
137 gmx_ana_selmethod_t sm_distance = {
138 "distance", REAL_VALUE, SMETH_DYNAMIC,
139 asize(smparams_distance), smparams_distance,
148 {"distance from POS [cutoff REAL]", asize(help_distance), help_distance},
151 /** \internal Selection method data for the \p distance method. */
152 gmx_ana_selmethod_t sm_mindistance = {
153 "mindistance", REAL_VALUE, SMETH_DYNAMIC,
154 asize(smparams_mindistance), smparams_mindistance,
163 {"mindistance from POS_EXPR [cutoff REAL]", asize(help_distance), help_distance},
166 /** \internal Selection method data for the \p within method. */
167 gmx_ana_selmethod_t sm_within = {
168 "within", GROUP_VALUE, SMETH_DYNAMIC,
169 asize(smparams_within), smparams_within,
178 {"within REAL of POS_EXPR", asize(help_distance), help_distance},
182 * \param[in] npar Not used (should be 2).
183 * \param[in,out] param Method parameters (should point to one of the distance
185 * \returns Pointer to the allocated data (\c t_methoddata_distance).
187 * Allocates memory for a \c t_methoddata_distance structure and
188 * initializes the parameter as follows:
189 * - the first parameter defines the value for
190 * \c t_methoddata_distance::cutoff.
191 * - the second parameter defines the reference positions and the value is
192 * stored in \c t_methoddata_distance::p.
195 init_data_common(int npar, gmx_ana_selparam_t *param)
197 t_methoddata_distance *data = new t_methoddata_distance();
198 param[0].val.u.r = &data->cutoff;
199 param[1].val.u.p = &data->p;
204 * \param top Not used.
205 * \param npar Not used (should be 2).
206 * \param param Method parameters (should point to one of the distance
208 * \param data Pointer to \c t_methoddata_distance to initialize.
209 * \returns 0 on success, a non-zero error code on failure.
211 * Initializes the neighborhood search data structure
212 * (\c t_methoddata_distance::nb).
213 * Also checks that the cutoff is valid.
216 init_common(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
218 t_methoddata_distance *d = (t_methoddata_distance *)data;
220 if ((param[0].flags & SPAR_SET) && d->cutoff <= 0)
222 GMX_THROW(gmx::InvalidInputError("Distance cutoff should be > 0"));
224 d->nb.setCutoff(d->cutoff);
228 * \param data Data to free (should point to a \c t_methoddata_distance).
230 * Frees the memory allocated for \c t_methoddata_distance::xref and
231 * \c t_methoddata_distance::nb.
234 free_data_common(void *data)
236 delete static_cast<t_methoddata_distance *>(data);
240 * \param[in] top Not used.
241 * \param[in] fr Current frame.
242 * \param[in] pbc PBC structure.
243 * \param data Should point to a \c t_methoddata_distance.
244 * \returns 0 on success, a non-zero error code on error.
246 * Initializes the neighborhood search for the current frame.
249 init_frame_common(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
251 t_methoddata_distance *d = (t_methoddata_distance *)data;
254 gmx::AnalysisNeighborhoodPositions pos(d->p.x, d->p.count());
255 d->nbsearch = d->nb.initSearch(pbc, pos);
259 * See sel_updatefunc_pos() for description of the parameters.
260 * \p data should point to a \c t_methoddata_distance.
262 * Calculates the distance of each position from \c t_methoddata_distance::p
263 * and puts them in \p out->u.r.
266 evaluate_distance(t_topology *top, t_trxframe *fr, t_pbc *pbc,
267 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
269 t_methoddata_distance *d = (t_methoddata_distance *)data;
271 out->nr = pos->m.mapb.nra;
272 for (int b = 0; b < pos->count(); ++b)
274 real dist = d->nbsearch.minimumDistance(pos->x[b]);
275 for (int i = pos->m.mapb.index[b]; i < pos->m.mapb.index[b+1]; ++i)
283 * See sel_updatefunc() for description of the parameters.
284 * \p data should point to a \c t_methoddata_distance.
286 * Finds the atoms that are closer than the defined cutoff to
287 * \c t_methoddata_distance::xref and puts them in \p out.g.
290 evaluate_within(t_topology *top, t_trxframe *fr, t_pbc *pbc,
291 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
293 t_methoddata_distance *d = (t_methoddata_distance *)data;
296 for (int b = 0; b < pos->count(); ++b)
298 if (d->nbsearch.isWithin(pos->x[b]))
300 gmx_ana_pos_add_to_group(out->u.g, pos, b);