2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013, by the GROMACS development team, led by
5 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
6 * others, as listed in the AUTHORS file in the top-level source
7 * 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 the \p same selection method.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
44 #include "gromacs/legacyheaders/macros.h"
45 #include "gromacs/legacyheaders/smalloc.h"
46 #include "gromacs/legacyheaders/string2.h"
48 #include "gromacs/selection/selmethod.h"
49 #include "gromacs/utility/exceptions.h"
52 #include "parsetree.h"
57 * Data structure for the \p same selection method.
59 * To avoid duplicate initialization code, the same data structure is used
60 * for matching both integer and string keywords; hence the unions.
62 * \ingroup module_selection
66 /** Value for each atom to match. */
74 * Number of values in the \p as array.
76 * For string values, this is actually the number of values in the
77 * \p as_s_sorted array.
80 /** Values to match against. */
88 * Separate array for sorted \p as.s array.
90 * The array of strings returned as the output value of a parameter should
91 * not be messed with to avoid memory corruption (the pointers in the array
92 * may be reused for several evaluations), so we keep our own copy for
96 /** Whether simple matching can be used. */
100 /** Allocates data for the \p same selection method. */
102 init_data_same(int npar, gmx_ana_selparam_t *param);
103 /** Initializes the \p same selection method. */
105 init_same(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
106 /** Frees the data allocated for the \p same selection method. */
108 free_data_same(void *data);
109 /** Initializes the evaluation of the \p same selection method for a frame. */
111 init_frame_same_int(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
112 /** Evaluates the \p same selection method. */
114 evaluate_same_int(t_topology *top, t_trxframe *fr, t_pbc *pbc,
115 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
116 /** Initializes the evaluation of the \p same selection method for a frame. */
118 init_frame_same_str(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
119 /** Evaluates the \p same selection method. */
121 evaluate_same_str(t_topology *top, t_trxframe *fr, t_pbc *pbc,
122 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
124 /** Parameters for the \p same selection method. */
125 static gmx_ana_selparam_t smparams_same_int[] = {
126 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_ATOMVAL},
127 {"as", {INT_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
130 /** Parameters for the \p same selection method. */
131 static gmx_ana_selparam_t smparams_same_str[] = {
132 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_ATOMVAL},
133 {"as", {STR_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
136 /** Help text for the \p same selection method. */
137 static const char *help_same[] = {
138 "EXTENDING SELECTIONS[PAR]",
140 "[TT]same KEYWORD as ATOM_EXPR[tt][PAR]",
142 "The keyword [TT]same[tt] can be used to select all atoms for which",
143 "the given [TT]KEYWORD[tt] matches any of the atoms in [TT]ATOM_EXPR[tt].",
144 "Keywords that evaluate to integer or string values are supported.",
147 /*! \internal \brief Selection method data for the \p same method. */
148 gmx_ana_selmethod_t sm_same = {
149 "same", GROUP_VALUE, 0,
150 asize(smparams_same_int), smparams_same_int,
156 &init_frame_same_int,
159 {"same KEYWORD as ATOM_EXPR", asize(help_same), help_same},
163 * Selection method data for the \p same method.
165 * This selection method is used for matching string keywords. The parser
166 * never sees this method; _gmx_selelem_custom_init_same() replaces sm_same
167 * with this method in cases where it is required.
169 static gmx_ana_selmethod_t sm_same_str = {
170 "same", GROUP_VALUE, SMETH_SINGLEVAL,
171 asize(smparams_same_str), smparams_same_str,
177 &init_frame_same_str,
180 {"same KEYWORD as ATOM_EXPR", asize(help_same), help_same},
184 * \param[in] npar Not used (should be 2).
185 * \param[in,out] param Method parameters (should point to a copy of
186 * ::smparams_same_int or ::smparams_same_str).
187 * \returns Pointer to the allocated data (\ref t_methoddata_same).
190 init_data_same(int npar, gmx_ana_selparam_t *param)
192 t_methoddata_same *data;
195 data->as_s_sorted = NULL;
196 param[1].nvalptr = &data->nas;
201 * \param[in,out] method The method to initialize.
202 * \param[in,out] params Pointer to the first parameter.
203 * \param[in] scanner Scanner data structure.
204 * \returns 0 on success, a non-zero error code on error.
206 * If \p *method is not a \c same method, this function returns zero
210 _gmx_selelem_custom_init_same(gmx_ana_selmethod_t **method,
211 const gmx::SelectionParserParameterListPointer ¶ms,
215 /* Do nothing if this is not a same method. */
216 if (!*method || (*method)->name != sm_same.name || params->empty())
221 const gmx::SelectionParserValueList &kwvalues = params->front().values();
222 if (kwvalues.size() != 1 || !kwvalues.front().hasExpressionValue()
223 || kwvalues.front().expr->type != SEL_EXPRESSION)
225 _gmx_selparser_error(scanner, "'same' should be followed by a single keyword");
228 gmx_ana_selmethod_t *kwmethod = kwvalues.front().expr->u.expr.method;
229 if (kwmethod->type == STR_VALUE)
231 *method = &sm_same_str;
234 /* We do custom processing for the "as" parameter. */
235 gmx::SelectionParserParameterList::iterator asparam = ++params->begin();
236 if (asparam != params->end() && asparam->name() == sm_same.param[1].name)
238 gmx::SelectionParserParameterList kwparams;
239 gmx::SelectionParserValueListPointer values(
240 new gmx::SelectionParserValueList(asparam->values()));
242 gmx::SelectionParserParameter::create(NULL, move(values)));
244 /* Create a second keyword evaluation element for the keyword given as
245 * the first parameter, evaluating the keyword in the group given by the
246 * second parameter. */
247 gmx::SelectionTreeElementPointer kwelem
248 = _gmx_sel_init_keyword_evaluator(kwmethod, kwparams, scanner);
249 // FIXME: Use exceptions.
254 /* Replace the second parameter with one with a value from \p kwelem. */
255 std::string pname = asparam->name();
256 *asparam = gmx::SelectionParserParameter::createFromExpression(pname, kwelem);
262 * \param top Not used.
263 * \param npar Not used (should be 2).
264 * \param param Initialized method parameters (should point to a copy of
265 * ::smparams_same_int or ::smparams_same_str).
266 * \param data Pointer to \ref t_methoddata_same to initialize.
267 * \returns 0 on success, -1 on failure.
270 init_same(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
272 t_methoddata_same *d = (t_methoddata_same *)data;
274 d->val.ptr = param[0].val.u.ptr;
275 d->as.ptr = param[1].val.u.ptr;
276 if (param[1].val.type == STR_VALUE)
278 snew(d->as_s_sorted, d->nas);
280 if (!(param[0].flags & SPAR_ATOMVAL))
282 GMX_THROW(gmx::InvalidInputError(
283 "The 'same' selection keyword combined with a "
284 "non-keyword does not make sense"));
289 * \param data Data to free (should point to a \ref t_methoddata_same).
292 free_data_same(void *data)
294 t_methoddata_same *d = (t_methoddata_same *)data;
296 sfree(d->as_s_sorted);
301 * Helper function for comparison of two integers.
304 cmp_int(const void *a, const void *b)
306 if (*(int *)a < *(int *)b)
310 if (*(int *)a > *(int *)b)
318 * \param[in] top Not used.
319 * \param[in] fr Current frame.
320 * \param[in] pbc PBC structure.
321 * \param data Should point to a \ref t_methoddata_same.
323 * Sorts the \c data->as.i array and removes identical values for faster and
327 init_frame_same_int(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
329 t_methoddata_same *d = (t_methoddata_same *)data;
332 /* Collapse adjacent values, and check whether the array is sorted. */
334 for (i = 1, j = 0; i < d->nas; ++i)
336 if (d->as.i[i] != d->as.i[j])
338 if (d->as.i[i] < d->as.i[j])
343 d->as.i[j] = d->as.i[i];
350 qsort(d->as.i, d->nas, sizeof(d->as.i[0]), &cmp_int);
351 /* More identical values may become adjacent after sorting. */
352 for (i = 1, j = 0; i < d->nas; ++i)
354 if (d->as.i[i] != d->as.i[j])
357 d->as.i[j] = d->as.i[i];
365 * See sel_updatefunc() for description of the parameters.
366 * \p data should point to a \c t_methoddata_same.
368 * Calculates which values in \c data->val.i can be found in \c data->as.i
369 * (assumed sorted), and writes the corresponding atoms to output.
370 * If \c data->val is sorted, uses a linear scan of both arrays, otherwise a
371 * binary search of \c data->as is performed for each block of values in
375 evaluate_same_int(t_topology *top, t_trxframe *fr, t_pbc *pbc,
376 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
378 t_methoddata_same *d = (t_methoddata_same *)data;
387 /* If we are sorted, we can do a simple linear scan. */
388 while (i < d->nas && d->as.i[i] < d->val.i[j])
395 /* If not, we must do a binary search of all the values. */
402 int itry = (i1 + i2) / 2;
403 if (d->as.i[itry] <= d->val.i[j])
412 i = (d->as.i[i1] == d->val.i[j] ? i1 : d->nas);
414 /* Check whether the value was found in the as list. */
415 if (i == d->nas || d->as.i[i] != d->val.i[j])
417 /* If not, skip all atoms with the same value. */
418 int tmpval = d->val.i[j];
420 while (j < g->isize && d->val.i[j] == tmpval)
427 /* Copy all the atoms with this value to the output. */
428 while (j < g->isize && d->val.i[j] == d->as.i[i])
430 out->u.g->index[out->u.g->isize++] = g->index[j];
434 if (j < g->isize && d->val.i[j] < d->val.i[j - 1])
442 * Helper function for comparison of two strings.
445 cmp_str(const void *a, const void *b)
447 return strcmp(*(char **)a, *(char **)b);
451 * \param[in] top Not used.
452 * \param[in] fr Current frame.
453 * \param[in] pbc PBC structure.
454 * \param data Should point to a \ref t_methoddata_same.
456 * Sorts the \c data->as.s array and removes identical values for faster and
460 init_frame_same_str(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
462 t_methoddata_same *d = (t_methoddata_same *)data;
465 /* Collapse adjacent values.
466 * For strings, it's unlikely that the values would be sorted originally,
467 * so set bSorted always to false. */
469 d->as_s_sorted[0] = d->as.s[0];
470 for (i = 1, j = 0; i < d->nas; ++i)
472 if (strcmp(d->as.s[i], d->as_s_sorted[j]) != 0)
475 d->as_s_sorted[j] = d->as.s[i];
480 qsort(d->as_s_sorted, d->nas, sizeof(d->as_s_sorted[0]), &cmp_str);
481 /* More identical values may become adjacent after sorting. */
482 for (i = 1, j = 0; i < d->nas; ++i)
484 if (strcmp(d->as_s_sorted[i], d->as_s_sorted[j]) != 0)
487 d->as_s_sorted[j] = d->as_s_sorted[i];
494 * See sel_updatefunc() for description of the parameters.
495 * \p data should point to a \c t_methoddata_same.
497 * Calculates which strings in \c data->val.s can be found in \c data->as.s
498 * (assumed sorted), and writes the corresponding atoms to output.
499 * A binary search of \c data->as is performed for each block of values in
503 evaluate_same_str(t_topology *top, t_trxframe *fr, t_pbc *pbc,
504 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
506 t_methoddata_same *d = (t_methoddata_same *)data;
513 /* Do a binary search of the strings. */
515 ptr = bsearch(&d->val.s[j], d->as_s_sorted, d->nas,
516 sizeof(d->as_s_sorted[0]), &cmp_str);
517 /* Check whether the value was found in the as list. */
520 /* If not, skip all atoms with the same value. */
521 const char *tmpval = d->val.s[j];
523 while (j < g->isize && strcmp(d->val.s[j], tmpval) == 0)
530 const char *tmpval = d->val.s[j];
531 /* Copy all the atoms with this value to the output. */
532 while (j < g->isize && strcmp(d->val.s[j], tmpval) == 0)
534 out->u.g->index[out->u.g->isize++] = g->index[j];