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
32 * \brief Implementation of distance-based selection methods.
34 * This file implements the \p distance, \p mindistance and \p within
48 #include <selmethod.h>
51 * Data structure for distance-based selection method.
53 * The same data structure is used by all the distance-based methods.
57 /** Cutoff distance. */
59 /** Positions of the reference points. */
61 /** Neighborhood search data. */
62 gmx_ana_nbsearch_t
*nb
;
63 } t_methoddata_distance
;
65 /** Allocates data for distance-based selection methods. */
67 init_data_common(int npar
, gmx_ana_selparam_t
*param
);
68 /** Initializes a distance-based selection method. */
70 init_common(t_topology
*top
, int npar
, gmx_ana_selparam_t
*param
, void *data
);
71 /** Frees the data allocated for a distance-based selection method. */
73 free_data_common(void *data
);
74 /** Initializes the evaluation of a distance-based within selection method for a frame. */
76 init_frame_common(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
, void *data
);
77 /** Evaluates the \p distance selection method. */
79 evaluate_distance(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
,
80 gmx_ana_pos_t
*pos
, gmx_ana_selvalue_t
*out
, void *data
);
81 /** Evaluates the \p within selection method. */
83 evaluate_within(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
,
84 gmx_ana_pos_t
*pos
, gmx_ana_selvalue_t
*out
, void *data
);
86 /** Parameters for the \p distance selection method. */
87 static gmx_ana_selparam_t smparams_distance
[] = {
88 {"cutoff", {REAL_VALUE
, 1, {NULL
}}, NULL
, SPAR_OPTIONAL
},
89 {"from", {POS_VALUE
, 1, {NULL
}}, NULL
, SPAR_DYNAMIC
},
92 /** Parameters for the \p mindistance selection method. */
93 static gmx_ana_selparam_t smparams_mindistance
[] = {
94 {"cutoff", {REAL_VALUE
, 1, {NULL
}}, NULL
, SPAR_OPTIONAL
},
95 {"from", {POS_VALUE
, -1, {NULL
}}, NULL
, SPAR_DYNAMIC
| SPAR_VARNUM
},
98 /** Parameters for the \p within selection method. */
99 static gmx_ana_selparam_t smparams_within
[] = {
100 {NULL
, {REAL_VALUE
, 1, {NULL
}}, NULL
, 0},
101 {"of", {POS_VALUE
, -1, {NULL
}}, NULL
, SPAR_DYNAMIC
| SPAR_VARNUM
},
104 /** Help text for the distance selection methods. */
105 static const char *help_distance
[] = {
106 "DISTANCE-BASED SELECTION KEYWORDS[PAR]",
108 "[TT]distance from POS [cutoff REAL][tt][BR]",
109 "[TT]mindistance from POS_EXPR [cutoff REAL][tt][BR]",
110 "[TT]within REAL of POS_EXPR[tt][PAR]",
112 "[TT]distance[tt] and [TT]mindistance[tt] calculate the distance from the",
113 "given position(s), the only difference being in that [TT]distance[tt]",
114 "only accepts a single position, while any number of positions can be",
115 "given for [TT]mindistance[tt], which then calculates the distance to the",
117 "[TT]within[tt] directly selects atoms that are within [TT]REAL[tt] of",
118 "[TT]POS_EXPR[tt].[PAR]",
120 "For the first two keywords, it is possible to specify a cutoff to speed",
121 "up the evaluation: all distances above the specified cutoff are",
122 "returned as equal to the cutoff.",
123 "Currently, this does nothing, but in the future, it allows the use of",
124 "grid-based neighborhood search techniques.",
127 /** \internal Selection method data for the \p distance method. */
128 gmx_ana_selmethod_t sm_distance
= {
129 "distance", REAL_VALUE
, SMETH_DYNAMIC
,
130 asize(smparams_distance
), smparams_distance
,
139 {"distance from POS [cutoff REAL]", asize(help_distance
), help_distance
},
142 /** \internal Selection method data for the \p distance method. */
143 gmx_ana_selmethod_t sm_mindistance
= {
144 "mindistance", REAL_VALUE
, SMETH_DYNAMIC
,
145 asize(smparams_mindistance
), smparams_mindistance
,
154 {"mindistance from POS_EXPR [cutoff REAL]", asize(help_distance
), help_distance
},
157 /** \internal Selection method data for the \p within method. */
158 gmx_ana_selmethod_t sm_within
= {
159 "within", GROUP_VALUE
, SMETH_DYNAMIC
,
160 asize(smparams_within
), smparams_within
,
169 {"within REAL of POS_EXPR", asize(help_distance
), help_distance
},
173 * \param[in] npar Not used (should be 2).
174 * \param[in,out] param Method parameters (should point to one of the distance
176 * \returns Pointer to the allocated data (\c t_methoddata_distance).
178 * Allocates memory for a \c t_methoddata_distance structure and
179 * initializes the parameter as follows:
180 * - the first parameter defines the value for
181 * \c t_methoddata_distance::cutoff.
182 * - the second parameter defines the reference positions and the value is
183 * stored in \c t_methoddata_distance::p.
186 init_data_common(int npar
, gmx_ana_selparam_t
*param
)
188 t_methoddata_distance
*data
;
192 param
[0].val
.u
.r
= &data
->cutoff
;
193 param
[1].val
.u
.p
= &data
->p
;
198 * \param top Not used.
199 * \param npar Not used (should be 2).
200 * \param param Method parameters (should point to one of the distance
202 * \param data Pointer to \c t_methoddata_distance to initialize.
203 * \returns 0 on success, a non-zero error code on failure.
205 * Initializes the neighborhood search data structure
206 * (\c t_methoddata_distance::nb).
207 * Also checks that the cutoff is valid.
210 init_common(t_topology
*top
, int npar
, gmx_ana_selparam_t
*param
, void *data
)
212 t_methoddata_distance
*d
= (t_methoddata_distance
*)data
;
214 if ((param
[0].flags
& SPAR_SET
) && d
->cutoff
<= 0)
216 fprintf(stderr
, "error: distance cutoff should be > 0");
219 return gmx_ana_nbsearch_create(&d
->nb
, d
->cutoff
, d
->p
.nr
);
223 * \param data Data to free (should point to a \c t_methoddata_distance).
225 * Frees the memory allocated for \c t_methoddata_distance::xref and
226 * \c t_methoddata_distance::nb.
229 free_data_common(void *data
)
231 t_methoddata_distance
*d
= (t_methoddata_distance
*)data
;
233 gmx_ana_nbsearch_free(d
->nb
);
237 * \param[in] top Not used.
238 * \param[in] fr Current frame.
239 * \param[in] pbc PBC structure.
240 * \param data Should point to a \c t_methoddata_distance.
241 * \returns 0 on success, a non-zero error code on error.
243 * Initializes the neighborhood search for the current frame.
246 init_frame_common(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
, void *data
)
248 t_methoddata_distance
*d
= (t_methoddata_distance
*)data
;
250 return gmx_ana_nbsearch_pos_init(d
->nb
, pbc
, &d
->p
);
254 * See sel_updatefunc_pos() for description of the parameters.
255 * \p data should point to a \c t_methoddata_distance.
257 * Calculates the distance of each position from \c t_methoddata_distance::p
258 * and puts them in \p out->u.r.
261 evaluate_distance(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
,
262 gmx_ana_pos_t
*pos
, gmx_ana_selvalue_t
*out
, void *data
)
264 t_methoddata_distance
*d
= (t_methoddata_distance
*)data
;
268 out
->nr
= pos
->g
->isize
;
269 for (b
= 0; b
< pos
->nr
; ++b
)
271 n
= gmx_ana_nbsearch_pos_mindist(d
->nb
, pos
, b
);
272 for (i
= pos
->m
.mapb
.index
[b
]; i
< pos
->m
.mapb
.index
[b
+1]; ++i
)
281 * See sel_updatefunc() for description of the parameters.
282 * \p data should point to a \c t_methoddata_distance.
284 * Finds the atoms that are closer than the defined cutoff to
285 * \c t_methoddata_distance::xref and puts them in \p out.g.
288 evaluate_within(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
,
289 gmx_ana_pos_t
*pos
, gmx_ana_selvalue_t
*out
, void *data
)
291 t_methoddata_distance
*d
= (t_methoddata_distance
*)data
;
295 for (b
= 0; b
< pos
->nr
; ++b
)
297 if (gmx_ana_nbsearch_pos_is_within(d
->nb
, pos
, b
))
299 gmx_ana_pos_append(NULL
, out
->u
.g
, pos
, b
, 0);