2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2016,2017,2018,2019, 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.
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
47 #include "gromacs/math/vec.h"
48 #include "gromacs/selection/nbsearch.h"
49 #include "gromacs/selection/position.h"
50 #include "gromacs/utility/arraysize.h"
51 #include "gromacs/utility/exceptions.h"
53 #include "selmethod.h"
54 #include "selmethod_impl.h"
58 * Data structure for distance-based selection method.
60 * The same data structure is used by all the distance-based methods.
62 * \ingroup module_selection
64 struct t_methoddata_distance
66 t_methoddata_distance() : cutoff(-1.0)
70 /** Cutoff distance. */
72 /** Positions of the reference points. */
74 /** Neighborhood search data. */
75 gmx::AnalysisNeighborhood nb;
76 /** Neighborhood search for an invididual frame. */
77 gmx::AnalysisNeighborhoodSearch nbsearch;
81 * Allocates data for distance-based selection methods.
83 * \param[in] npar Not used (should be 2).
84 * \param[in,out] param Method parameters (should point to one of the distance
86 * \returns Pointer to the allocated data (\c t_methoddata_distance).
88 * Allocates memory for a \c t_methoddata_distance structure and
89 * initializes the parameter as follows:
90 * - the first parameter defines the value for
91 * \c t_methoddata_distance::cutoff.
92 * - the second parameter defines the reference positions and the value is
93 * stored in \c t_methoddata_distance::p.
96 init_data_common(int npar, gmx_ana_selparam_t *param);
98 * Initializes a distance-based selection method.
100 * \param top Not used.
101 * \param npar Not used (should be 2).
102 * \param param Method parameters (should point to one of the distance
104 * \param data Pointer to \c t_methoddata_distance to initialize.
105 * \returns 0 on success, a non-zero error code on failure.
107 * Initializes the neighborhood search data structure
108 * (\c t_methoddata_distance::nb).
109 * Also checks that the cutoff is valid.
112 init_common(const gmx_mtop_t *top, int npar, gmx_ana_selparam_t *param, void *data);
113 /** Frees the data allocated for a distance-based selection method. */
115 free_data_common(void *data);
117 * Initializes the evaluation of a distance-based within selection method for a
120 * \param[in] context Evaluation context.
121 * \param data Should point to a \c t_methoddata_distance.
122 * \returns 0 on success, a non-zero error code on error.
124 * Initializes the neighborhood search for the current frame.
127 init_frame_common(const gmx::SelMethodEvalContext &context, void *data);
128 /** Evaluates the \p distance selection method. */
130 evaluate_distance(const gmx::SelMethodEvalContext & /*context*/,
131 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
132 /** Evaluates the \p within selection method. */
134 evaluate_within(const gmx::SelMethodEvalContext & /*context*/,
135 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
137 /** Parameters for the \p distance selection method. */
138 static gmx_ana_selparam_t smparams_distance[] = {
139 {"cutoff", {REAL_VALUE, 1, {nullptr}}, nullptr, SPAR_OPTIONAL},
140 {"from", {POS_VALUE, 1, {nullptr}}, nullptr, SPAR_DYNAMIC},
143 /** Parameters for the \p mindistance selection method. */
144 static gmx_ana_selparam_t smparams_mindistance[] = {
145 {"cutoff", {REAL_VALUE, 1, {nullptr}}, nullptr, SPAR_OPTIONAL},
146 {"from", {POS_VALUE, -1, {nullptr}}, nullptr, SPAR_DYNAMIC | SPAR_VARNUM},
149 /** Parameters for the \p within selection method. */
150 static gmx_ana_selparam_t smparams_within[] = {
151 {nullptr, {REAL_VALUE, 1, {nullptr}}, nullptr, 0},
152 {"of", {POS_VALUE, -1, {nullptr}}, nullptr, SPAR_DYNAMIC | SPAR_VARNUM},
155 //! Help title for distance selection methods.
156 static const char helptitle_distance[] = "Selecting based on distance";
157 //! Help text for distance selection methods.
158 static const char *const help_distance[] = {
161 " distance from POS [cutoff REAL]",
162 " mindistance from POS_EXPR [cutoff REAL]",
163 " within REAL of POS_EXPR",
165 "[TT]distance[tt] and [TT]mindistance[tt] calculate the distance from the",
166 "given position(s), the only difference being in that [TT]distance[tt]",
167 "only accepts a single position, while any number of positions can be",
168 "given for [TT]mindistance[tt], which then calculates the distance to the",
170 "[TT]within[tt] directly selects atoms that are within [TT]REAL[tt] of",
171 "[TT]POS_EXPR[tt].[PAR]",
173 "For the first two keywords, it is possible to specify a cutoff to speed",
174 "up the evaluation: all distances above the specified cutoff are",
175 "returned as equal to the cutoff.",
178 /** Selection method data for the \p distance method. */
179 gmx_ana_selmethod_t sm_distance = {
180 "distance", REAL_VALUE, SMETH_DYNAMIC,
181 asize(smparams_distance), smparams_distance,
190 {"distance from POS [cutoff REAL]",
191 helptitle_distance, asize(help_distance), help_distance},
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,
206 {"mindistance from POS_EXPR [cutoff REAL]",
207 helptitle_distance, asize(help_distance), help_distance},
210 /** Selection method data for the \p within method. */
211 gmx_ana_selmethod_t sm_within = {
212 "within", GROUP_VALUE, SMETH_DYNAMIC,
213 asize(smparams_within), smparams_within,
222 {"within REAL of POS_EXPR",
223 helptitle_distance, asize(help_distance), help_distance},
227 init_data_common(int /* npar */, gmx_ana_selparam_t *param)
229 t_methoddata_distance *data = new t_methoddata_distance();
230 param[0].val.u.r = &data->cutoff;
231 param[1].val.u.p = &data->p;
236 init_common(const gmx_mtop_t * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
238 t_methoddata_distance *d = static_cast<t_methoddata_distance *>(data);
240 if ((param[0].flags & SPAR_SET) && d->cutoff <= 0)
242 GMX_THROW(gmx::InvalidInputError("Distance cutoff should be > 0"));
244 d->nb.setCutoff(d->cutoff);
248 * \param data Data to free (should point to a \c t_methoddata_distance).
250 * Frees the memory allocated for \c t_methoddata_distance::xref and
251 * \c t_methoddata_distance::nb.
254 free_data_common(void *data)
256 delete static_cast<t_methoddata_distance *>(data);
260 init_frame_common(const gmx::SelMethodEvalContext &context, void *data)
262 t_methoddata_distance *d = static_cast<t_methoddata_distance *>(data);
265 gmx::AnalysisNeighborhoodPositions pos(d->p.x, d->p.count());
266 d->nbsearch = d->nb.initSearch(context.pbc, pos);
270 * See sel_updatefunc_pos() for description of the parameters.
271 * \p data should point to a \c t_methoddata_distance.
273 * Calculates the distance of each position from \c t_methoddata_distance::p
274 * and puts them in \p out->u.r.
277 evaluate_distance(const gmx::SelMethodEvalContext & /*context*/,
278 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
280 t_methoddata_distance *d = static_cast<t_methoddata_distance *>(data);
282 out->nr = pos->count();
283 for (int i = 0; i < pos->count(); ++i)
285 out->u.r[i] = d->nbsearch.minimumDistance(pos->x[i]);
290 * See sel_updatefunc() for description of the parameters.
291 * \p data should point to a \c t_methoddata_distance.
293 * Finds the atoms that are closer than the defined cutoff to
294 * \c t_methoddata_distance::xref and puts them in \p out.g.
297 evaluate_within(const gmx::SelMethodEvalContext & /*context*/,
298 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
300 t_methoddata_distance *d = static_cast<t_methoddata_distance *>(data);
303 for (int b = 0; b < pos->count(); ++b)
305 if (d->nbsearch.isWithin(pos->x[b]))
307 gmx_ana_pos_add_to_group(out->u.g, pos, b);