3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
9 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11 * Copyright (c) 2001-2009, The GROMACS development team,
12 * check out http://www.gromacs.org for more information.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
33 * Implements distance-based selection methods.
35 * This file implements the \p distance, \p mindistance and \p within
38 * \author Teemu Murtola <teemu.murtola@cbr.su.se>
39 * \ingroup module_selection
41 #include "gromacs/legacyheaders/macros.h"
42 #include "gromacs/legacyheaders/pbc.h"
43 #include "gromacs/legacyheaders/smalloc.h"
44 #include "gromacs/legacyheaders/vec.h"
46 #include "gromacs/selection/nbsearch.h"
47 #include "gromacs/selection/position.h"
48 #include "gromacs/selection/selmethod.h"
49 #include "gromacs/utility/exceptions.h"
52 * Data structure for distance-based selection method.
54 * The same data structure is used by all the distance-based methods.
58 /** Cutoff distance. */
60 /** Positions of the reference points. */
62 /** Neighborhood search data. */
63 gmx_ana_nbsearch_t *nb;
64 } t_methoddata_distance;
66 /** Allocates data for distance-based selection methods. */
68 init_data_common(int npar, gmx_ana_selparam_t *param);
69 /** Initializes a distance-based selection method. */
71 init_common(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
72 /** Frees the data allocated for a distance-based selection method. */
74 free_data_common(void *data);
75 /** Initializes the evaluation of a distance-based within selection method for a frame. */
77 init_frame_common(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
78 /** Evaluates the \p distance selection method. */
80 evaluate_distance(t_topology *top, t_trxframe *fr, t_pbc *pbc,
81 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
82 /** Evaluates the \p within selection method. */
84 evaluate_within(t_topology *top, t_trxframe *fr, t_pbc *pbc,
85 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
87 /** Parameters for the \p distance selection method. */
88 static gmx_ana_selparam_t smparams_distance[] = {
89 {"cutoff", {REAL_VALUE, 1, {NULL}}, NULL, SPAR_OPTIONAL},
90 {"from", {POS_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
93 /** Parameters for the \p mindistance selection method. */
94 static gmx_ana_selparam_t smparams_mindistance[] = {
95 {"cutoff", {REAL_VALUE, 1, {NULL}}, NULL, SPAR_OPTIONAL},
96 {"from", {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
99 /** Parameters for the \p within selection method. */
100 static gmx_ana_selparam_t smparams_within[] = {
101 {NULL, {REAL_VALUE, 1, {NULL}}, NULL, 0},
102 {"of", {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
105 /** Help text for the distance selection methods. */
106 static const char *help_distance[] = {
107 "DISTANCE-BASED SELECTION KEYWORDS[PAR]",
109 "[TT]distance from POS [cutoff REAL][tt][BR]",
110 "[TT]mindistance from POS_EXPR [cutoff REAL][tt][BR]",
111 "[TT]within REAL of POS_EXPR[tt][PAR]",
113 "[TT]distance[tt] and [TT]mindistance[tt] calculate the distance from the",
114 "given position(s), the only difference being in that [TT]distance[tt]",
115 "only accepts a single position, while any number of positions can be",
116 "given for [TT]mindistance[tt], which then calculates the distance to the",
118 "[TT]within[tt] directly selects atoms that are within [TT]REAL[tt] of",
119 "[TT]POS_EXPR[tt].[PAR]",
121 "For the first two keywords, it is possible to specify a cutoff to speed",
122 "up the evaluation: all distances above the specified cutoff are",
123 "returned as equal to the cutoff.",
126 /** \internal Selection method data for the \p distance method. */
127 gmx_ana_selmethod_t sm_distance = {
128 "distance", REAL_VALUE, SMETH_DYNAMIC,
129 asize(smparams_distance), smparams_distance,
138 {"distance from POS [cutoff REAL]", asize(help_distance), help_distance},
141 /** \internal Selection method data for the \p distance method. */
142 gmx_ana_selmethod_t sm_mindistance = {
143 "mindistance", REAL_VALUE, SMETH_DYNAMIC,
144 asize(smparams_mindistance), smparams_mindistance,
153 {"mindistance from POS_EXPR [cutoff REAL]", asize(help_distance), help_distance},
156 /** \internal Selection method data for the \p within method. */
157 gmx_ana_selmethod_t sm_within = {
158 "within", GROUP_VALUE, SMETH_DYNAMIC,
159 asize(smparams_within), smparams_within,
168 {"within REAL of POS_EXPR", asize(help_distance), help_distance},
172 * \param[in] npar Not used (should be 2).
173 * \param[in,out] param Method parameters (should point to one of the distance
175 * \returns Pointer to the allocated data (\c t_methoddata_distance).
177 * Allocates memory for a \c t_methoddata_distance structure and
178 * initializes the parameter as follows:
179 * - the first parameter defines the value for
180 * \c t_methoddata_distance::cutoff.
181 * - the second parameter defines the reference positions and the value is
182 * stored in \c t_methoddata_distance::p.
185 init_data_common(int npar, gmx_ana_selparam_t *param)
187 t_methoddata_distance *data;
191 param[0].val.u.r = &data->cutoff;
192 param[1].val.u.p = &data->p;
197 * \param top Not used.
198 * \param npar Not used (should be 2).
199 * \param param Method parameters (should point to one of the distance
201 * \param data Pointer to \c t_methoddata_distance to initialize.
202 * \returns 0 on success, a non-zero error code on failure.
204 * Initializes the neighborhood search data structure
205 * (\c t_methoddata_distance::nb).
206 * Also checks that the cutoff is valid.
209 init_common(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
211 t_methoddata_distance *d = (t_methoddata_distance *)data;
213 if ((param[0].flags & SPAR_SET) && d->cutoff <= 0)
215 GMX_THROW(gmx::InvalidInputError("Distance cutoff should be > 0"));
217 d->nb = gmx_ana_nbsearch_create(d->cutoff, d->p.nr);
221 * \param data Data to free (should point to a \c t_methoddata_distance).
223 * Frees the memory allocated for \c t_methoddata_distance::xref and
224 * \c t_methoddata_distance::nb.
227 free_data_common(void *data)
229 t_methoddata_distance *d = (t_methoddata_distance *)data;
233 gmx_ana_nbsearch_free(d->nb);
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.
245 * Initializes the neighborhood search for the current frame.
248 init_frame_common(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
250 t_methoddata_distance *d = (t_methoddata_distance *)data;
252 gmx_ana_nbsearch_pos_init(d->nb, pbc, &d->p);
256 * See sel_updatefunc_pos() for description of the parameters.
257 * \p data should point to a \c t_methoddata_distance.
259 * Calculates the distance of each position from \c t_methoddata_distance::p
260 * and puts them in \p out->u.r.
263 evaluate_distance(t_topology *top, t_trxframe *fr, t_pbc *pbc,
264 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
266 t_methoddata_distance *d = (t_methoddata_distance *)data;
270 out->nr = pos->g->isize;
271 for (b = 0; b < pos->nr; ++b)
273 n = gmx_ana_nbsearch_pos_mindist(d->nb, pos, b);
274 for (i = pos->m.mapb.index[b]; i < pos->m.mapb.index[b+1]; ++i)
282 * See sel_updatefunc() for description of the parameters.
283 * \p data should point to a \c t_methoddata_distance.
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.
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)
292 t_methoddata_distance *d = (t_methoddata_distance *)data;
296 for (b = 0; b < pos->nr; ++b)
298 if (gmx_ana_nbsearch_pos_is_within(d->nb, pos, b))
300 gmx_ana_pos_append(NULL, out->u.g, pos, b, 0);