Use gmx_mtop_t in selections, part 3
[gromacs.git] / src / gromacs / selection / sm_distance.cpp
blob738657c92f051fd9cf11fa8da9407d0d0bae0c12
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2016, 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.
35 /*! \internal \file
36 * \brief
37 * Implements distance-based selection methods.
39 * This file implements the \p distance, \p mindistance and \p within
40 * selection methods.
42 * \author Teemu Murtola <teemu.murtola@gmail.com>
43 * \ingroup module_selection
45 #include "gmxpre.h"
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"
55 /*! \internal
56 * \brief
57 * Data structure for distance-based selection method.
59 * The same data structure is used by all the distance-based methods.
61 * \ingroup module_selection
63 struct t_methoddata_distance
65 t_methoddata_distance() : cutoff(-1.0)
69 /** Cutoff distance. */
70 real cutoff;
71 /** Positions of the reference points. */
72 gmx_ana_pos_t p;
73 /** Neighborhood search data. */
74 gmx::AnalysisNeighborhood nb;
75 /** Neighborhood search for an invididual frame. */
76 gmx::AnalysisNeighborhoodSearch nbsearch;
79 /*! \brief
80 * Allocates data for distance-based selection methods.
82 * \param[in] npar Not used (should be 2).
83 * \param[in,out] param Method parameters (should point to one of the distance
84 * parameter arrays).
85 * \returns Pointer to the allocated data (\c t_methoddata_distance).
87 * Allocates memory for a \c t_methoddata_distance structure and
88 * initializes the parameter as follows:
89 * - the first parameter defines the value for
90 * \c t_methoddata_distance::cutoff.
91 * - the second parameter defines the reference positions and the value is
92 * stored in \c t_methoddata_distance::p.
94 static void *
95 init_data_common(int npar, gmx_ana_selparam_t *param);
96 /*! \brief
97 * Initializes a distance-based selection method.
99 * \param top Not used.
100 * \param npar Not used (should be 2).
101 * \param param Method parameters (should point to one of the distance
102 * parameter arrays).
103 * \param data Pointer to \c t_methoddata_distance to initialize.
104 * \returns 0 on success, a non-zero error code on failure.
106 * Initializes the neighborhood search data structure
107 * (\c t_methoddata_distance::nb).
108 * Also checks that the cutoff is valid.
110 static void
111 init_common(const gmx_mtop_t *top, int npar, gmx_ana_selparam_t *param, void *data);
112 /** Frees the data allocated for a distance-based selection method. */
113 static void
114 free_data_common(void *data);
115 /*! \brief
116 * Initializes the evaluation of a distance-based within selection method for a
117 * frame.
119 * \param[in] context Evaluation context.
120 * \param data Should point to a \c t_methoddata_distance.
121 * \returns 0 on success, a non-zero error code on error.
123 * Initializes the neighborhood search for the current frame.
125 static void
126 init_frame_common(const gmx::SelMethodEvalContext &context, void *data);
127 /** Evaluates the \p distance selection method. */
128 static void
129 evaluate_distance(const gmx::SelMethodEvalContext & /*context*/,
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(const gmx::SelMethodEvalContext & /*context*/,
134 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
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},
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},
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},
154 //! Help title for distance selection methods.
155 static const char helptitle_distance[] = "Selecting based on distance";
156 //! Help text for distance selection methods.
157 static const char *const help_distance[] = {
158 "::",
160 " distance from POS [cutoff REAL]",
161 " mindistance from POS_EXPR [cutoff REAL]",
162 " within REAL of POS_EXPR",
164 "[TT]distance[tt] and [TT]mindistance[tt] calculate the distance from the",
165 "given position(s), the only difference being in that [TT]distance[tt]",
166 "only accepts a single position, while any number of positions can be",
167 "given for [TT]mindistance[tt], which then calculates the distance to the",
168 "closest position.",
169 "[TT]within[tt] directly selects atoms that are within [TT]REAL[tt] of",
170 "[TT]POS_EXPR[tt].[PAR]",
172 "For the first two keywords, it is possible to specify a cutoff to speed",
173 "up the evaluation: all distances above the specified cutoff are",
174 "returned as equal to the cutoff.",
177 /** Selection method data for the \p distance method. */
178 gmx_ana_selmethod_t sm_distance = {
179 "distance", REAL_VALUE, SMETH_DYNAMIC,
180 asize(smparams_distance), smparams_distance,
181 &init_data_common,
182 NULL,
183 &init_common,
184 NULL,
185 &free_data_common,
186 &init_frame_common,
187 NULL,
188 &evaluate_distance,
189 {"distance from POS [cutoff REAL]",
190 helptitle_distance, asize(help_distance), help_distance},
193 /** Selection method data for the \p distance method. */
194 gmx_ana_selmethod_t sm_mindistance = {
195 "mindistance", REAL_VALUE, SMETH_DYNAMIC,
196 asize(smparams_mindistance), smparams_mindistance,
197 &init_data_common,
198 NULL,
199 &init_common,
200 NULL,
201 &free_data_common,
202 &init_frame_common,
203 NULL,
204 &evaluate_distance,
205 {"mindistance from POS_EXPR [cutoff REAL]",
206 helptitle_distance, asize(help_distance), help_distance},
209 /** Selection method data for the \p within method. */
210 gmx_ana_selmethod_t sm_within = {
211 "within", GROUP_VALUE, SMETH_DYNAMIC,
212 asize(smparams_within), smparams_within,
213 &init_data_common,
214 NULL,
215 &init_common,
216 NULL,
217 &free_data_common,
218 &init_frame_common,
219 NULL,
220 &evaluate_within,
221 {"within REAL of POS_EXPR",
222 helptitle_distance, asize(help_distance), help_distance},
225 static void *
226 init_data_common(int /* npar */, gmx_ana_selparam_t *param)
228 t_methoddata_distance *data = new t_methoddata_distance();
229 param[0].val.u.r = &data->cutoff;
230 param[1].val.u.p = &data->p;
231 return data;
234 static void
235 init_common(const gmx_mtop_t * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
237 t_methoddata_distance *d = static_cast<t_methoddata_distance *>(data);
239 if ((param[0].flags & SPAR_SET) && d->cutoff <= 0)
241 GMX_THROW(gmx::InvalidInputError("Distance cutoff should be > 0"));
243 d->nb.setCutoff(d->cutoff);
247 * \param data Data to free (should point to a \c t_methoddata_distance).
249 * Frees the memory allocated for \c t_methoddata_distance::xref and
250 * \c t_methoddata_distance::nb.
252 static void
253 free_data_common(void *data)
255 delete static_cast<t_methoddata_distance *>(data);
258 static void
259 init_frame_common(const gmx::SelMethodEvalContext &context, void *data)
261 t_methoddata_distance *d = static_cast<t_methoddata_distance *>(data);
263 d->nbsearch.reset();
264 gmx::AnalysisNeighborhoodPositions pos(d->p.x, d->p.count());
265 d->nbsearch = d->nb.initSearch(context.pbc, pos);
269 * See sel_updatefunc_pos() for description of the parameters.
270 * \p data should point to a \c t_methoddata_distance.
272 * Calculates the distance of each position from \c t_methoddata_distance::p
273 * and puts them in \p out->u.r.
275 static void
276 evaluate_distance(const gmx::SelMethodEvalContext & /*context*/,
277 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
279 t_methoddata_distance *d = static_cast<t_methoddata_distance *>(data);
281 out->nr = pos->count();
282 for (int i = 0; i < pos->count(); ++i)
284 out->u.r[i] = d->nbsearch.minimumDistance(pos->x[i]);
289 * See sel_updatefunc() for description of the parameters.
290 * \p data should point to a \c t_methoddata_distance.
292 * Finds the atoms that are closer than the defined cutoff to
293 * \c t_methoddata_distance::xref and puts them in \p out.g.
295 static void
296 evaluate_within(const gmx::SelMethodEvalContext & /*context*/,
297 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
299 t_methoddata_distance *d = static_cast<t_methoddata_distance *>(data);
301 out->u.g->isize = 0;
302 for (int b = 0; b < pos->count(); ++b)
304 if (d->nbsearch.isWithin(pos->x[b]))
306 gmx_ana_pos_add_to_group(out->u.g, pos, b);