2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015, 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/legacyheaders/macros.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/selection/nbsearch.h"
50 #include "gromacs/selection/position.h"
51 #include "gromacs/utility/exceptions.h"
53 #include "selmethod.h"
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. */
71 /** Positions of the reference points. */
73 /** Neighborhood search data. */
74 gmx::AnalysisNeighborhood nb
;
75 /** Neighborhood search for an invididual frame. */
76 gmx::AnalysisNeighborhoodSearch nbsearch
;
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
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.
95 init_data_common(int npar
, gmx_ana_selparam_t
*param
);
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
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.
111 init_common(t_topology
*top
, int npar
, gmx_ana_selparam_t
*param
, void *data
);
112 /** Frees the data allocated for a distance-based selection method. */
114 free_data_common(void *data
);
116 * Initializes the evaluation of a distance-based within selection method for a
119 * \param[in] top Not used.
120 * \param[in] fr Current frame.
121 * \param[in] pbc PBC structure.
122 * \param data Should point to a \c t_methoddata_distance.
123 * \returns 0 on success, a non-zero error code on error.
125 * Initializes the neighborhood search for the current frame.
128 init_frame_common(t_topology
*top
, t_trxframe
* fr
, t_pbc
*pbc
, void *data
);
129 /** Evaluates the \p distance selection method. */
131 evaluate_distance(t_topology
* /* top */, t_trxframe
* /* fr */, t_pbc
* /* pbc */,
132 gmx_ana_pos_t
*pos
, gmx_ana_selvalue_t
*out
, void *data
);
133 /** Evaluates the \p within selection method. */
135 evaluate_within(t_topology
* /* top */, t_trxframe
* /* fr */, t_pbc
* /* pbc */,
136 gmx_ana_pos_t
*pos
, gmx_ana_selvalue_t
*out
, void *data
);
138 /** Parameters for the \p distance selection method. */
139 static gmx_ana_selparam_t smparams_distance
[] = {
140 {"cutoff", {REAL_VALUE
, 1, {NULL
}}, NULL
, SPAR_OPTIONAL
},
141 {"from", {POS_VALUE
, 1, {NULL
}}, NULL
, SPAR_DYNAMIC
},
144 /** Parameters for the \p mindistance selection method. */
145 static gmx_ana_selparam_t smparams_mindistance
[] = {
146 {"cutoff", {REAL_VALUE
, 1, {NULL
}}, NULL
, SPAR_OPTIONAL
},
147 {"from", {POS_VALUE
, -1, {NULL
}}, NULL
, SPAR_DYNAMIC
| SPAR_VARNUM
},
150 /** Parameters for the \p within selection method. */
151 static gmx_ana_selparam_t smparams_within
[] = {
152 {NULL
, {REAL_VALUE
, 1, {NULL
}}, NULL
, 0},
153 {"of", {POS_VALUE
, -1, {NULL
}}, NULL
, SPAR_DYNAMIC
| SPAR_VARNUM
},
156 /** Help text for the distance selection methods. */
157 static const char *help_distance
[] = {
158 "DISTANCE-BASED SELECTION KEYWORDS",
162 " distance from POS [cutoff REAL]",
163 " mindistance from POS_EXPR [cutoff REAL]",
164 " within REAL of POS_EXPR",
166 "[TT]distance[tt] and [TT]mindistance[tt] calculate the distance from the",
167 "given position(s), the only difference being in that [TT]distance[tt]",
168 "only accepts a single position, while any number of positions can be",
169 "given for [TT]mindistance[tt], which then calculates the distance to the",
171 "[TT]within[tt] directly selects atoms that are within [TT]REAL[tt] of",
172 "[TT]POS_EXPR[tt].[PAR]",
174 "For the first two keywords, it is possible to specify a cutoff to speed",
175 "up the evaluation: all distances above the specified cutoff are",
176 "returned as equal to the cutoff.",
179 /** Selection method data for the \p distance method. */
180 gmx_ana_selmethod_t sm_distance
= {
181 "distance", REAL_VALUE
, SMETH_DYNAMIC
,
182 asize(smparams_distance
), smparams_distance
,
191 {"distance from POS [cutoff REAL]", 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]", 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
,
221 {"within REAL of POS_EXPR", asize(help_distance
), help_distance
},
225 init_data_common(int /* npar */, gmx_ana_selparam_t
*param
)
227 t_methoddata_distance
*data
= new t_methoddata_distance();
228 param
[0].val
.u
.r
= &data
->cutoff
;
229 param
[1].val
.u
.p
= &data
->p
;
234 init_common(t_topology
* /* top */, int /* npar */, gmx_ana_selparam_t
*param
, void *data
)
236 t_methoddata_distance
*d
= static_cast<t_methoddata_distance
*>(data
);
238 if ((param
[0].flags
& SPAR_SET
) && d
->cutoff
<= 0)
240 GMX_THROW(gmx::InvalidInputError("Distance cutoff should be > 0"));
242 d
->nb
.setCutoff(d
->cutoff
);
246 * \param data Data to free (should point to a \c t_methoddata_distance).
248 * Frees the memory allocated for \c t_methoddata_distance::xref and
249 * \c t_methoddata_distance::nb.
252 free_data_common(void *data
)
254 delete static_cast<t_methoddata_distance
*>(data
);
258 init_frame_common(t_topology
* /* top */, t_trxframe
* /* fr */, t_pbc
*pbc
, void *data
)
260 t_methoddata_distance
*d
= static_cast<t_methoddata_distance
*>(data
);
263 gmx::AnalysisNeighborhoodPositions
pos(d
->p
.x
, d
->p
.count());
264 d
->nbsearch
= d
->nb
.initSearch(pbc
, pos
);
268 * See sel_updatefunc_pos() for description of the parameters.
269 * \p data should point to a \c t_methoddata_distance.
271 * Calculates the distance of each position from \c t_methoddata_distance::p
272 * and puts them in \p out->u.r.
275 evaluate_distance(t_topology
* /* top */, t_trxframe
* /* fr */, t_pbc
* /* pbc */,
276 gmx_ana_pos_t
*pos
, gmx_ana_selvalue_t
*out
, void *data
)
278 t_methoddata_distance
*d
= static_cast<t_methoddata_distance
*>(data
);
280 out
->nr
= pos
->count();
281 for (int i
= 0; i
< pos
->count(); ++i
)
283 out
->u
.r
[i
] = d
->nbsearch
.minimumDistance(pos
->x
[i
]);
288 * See sel_updatefunc() for description of the parameters.
289 * \p data should point to a \c t_methoddata_distance.
291 * Finds the atoms that are closer than the defined cutoff to
292 * \c t_methoddata_distance::xref and puts them in \p out.g.
295 evaluate_within(t_topology
* /* top */, t_trxframe
* /* fr */, t_pbc
* /* pbc */,
296 gmx_ana_pos_t
*pos
, gmx_ana_selvalue_t
*out
, void *data
)
298 t_methoddata_distance
*d
= static_cast<t_methoddata_distance
*>(data
);
301 for (int b
= 0; b
< pos
->count(); ++b
)
303 if (d
->nbsearch
.isWithin(pos
->x
[b
]))
305 gmx_ana_pos_add_to_group(out
->u
.g
, pos
, b
);