2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014, 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 the \p same selection method.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
46 #include "gromacs/legacyheaders/macros.h"
47 #include "gromacs/utility/exceptions.h"
48 #include "gromacs/utility/smalloc.h"
51 #include "parsetree.h"
53 #include "selmethod.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. */
101 * Allocates data for the \p same selection method.
103 * \param[in] npar Not used (should be 2).
104 * \param[in,out] param Method parameters (should point to a copy of
105 * ::smparams_same_int or ::smparams_same_str).
106 * \returns Pointer to the allocated data (\ref t_methoddata_same).
109 init_data_same(int npar, gmx_ana_selparam_t *param);
111 * Initializes the \p same selection method.
113 * \param top Not used.
114 * \param npar Not used (should be 2).
115 * \param param Initialized method parameters (should point to a copy of
116 * ::smparams_same_int or ::smparams_same_str).
117 * \param data Pointer to \ref t_methoddata_same to initialize.
118 * \returns 0 on success, -1 on failure.
121 init_same(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
122 /** Frees the data allocated for the \p same selection method. */
124 free_data_same(void *data);
126 * Initializes the evaluation of the \p same selection method for a frame.
128 * \param[in] top Not used.
129 * \param[in] fr Current frame.
130 * \param[in] pbc PBC structure.
131 * \param data Should point to a \ref t_methoddata_same.
133 * Sorts the \c data->as.i array and removes identical values for faster and
137 init_frame_same_int(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
138 /** Evaluates the \p same selection method. */
140 evaluate_same_int(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
141 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
143 * Initializes the evaluation of the \p same selection method for a frame.
145 * \param[in] top Not used.
146 * \param[in] fr Current frame.
147 * \param[in] pbc PBC structure.
148 * \param data Should point to a \ref t_methoddata_same.
150 * Sorts the \c data->as.s array and removes identical values for faster and
154 init_frame_same_str(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
155 /** Evaluates the \p same selection method. */
157 evaluate_same_str(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
158 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
160 /** Parameters for the \p same selection method. */
161 static gmx_ana_selparam_t smparams_same_int[] = {
162 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_ATOMVAL},
163 {"as", {INT_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
166 /** Parameters for the \p same selection method. */
167 static gmx_ana_selparam_t smparams_same_str[] = {
168 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_ATOMVAL},
169 {"as", {STR_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
172 /** Help text for the \p same selection method. */
173 static const char *help_same[] = {
174 "EXTENDING SELECTIONS[PAR]",
176 "[TT]same KEYWORD as ATOM_EXPR[tt][PAR]",
178 "The keyword [TT]same[tt] can be used to select all atoms for which",
179 "the given [TT]KEYWORD[tt] matches any of the atoms in [TT]ATOM_EXPR[tt].",
180 "Keywords that evaluate to integer or string values are supported.",
183 /*! \internal \brief Selection method data for the \p same method. */
184 gmx_ana_selmethod_t sm_same = {
185 "same", GROUP_VALUE, 0,
186 asize(smparams_same_int), smparams_same_int,
192 &init_frame_same_int,
195 {"same KEYWORD as ATOM_EXPR", asize(help_same), help_same},
199 * Selection method data for the \p same method.
201 * This selection method is used for matching string keywords. The parser
202 * never sees this method; _gmx_selelem_custom_init_same() replaces sm_same
203 * with this method in cases where it is required.
205 static gmx_ana_selmethod_t sm_same_str = {
206 "same", GROUP_VALUE, SMETH_SINGLEVAL,
207 asize(smparams_same_str), smparams_same_str,
213 &init_frame_same_str,
216 {"same KEYWORD as ATOM_EXPR", asize(help_same), help_same},
220 init_data_same(int /* npar */, gmx_ana_selparam_t *param)
222 t_methoddata_same *data;
225 data->as_s_sorted = NULL;
226 param[1].nvalptr = &data->nas;
231 * \param[in,out] method The method to initialize.
232 * \param[in,out] params Pointer to the first parameter.
233 * \param[in] scanner Scanner data structure.
234 * \returns 0 on success, a non-zero error code on error.
236 * If \p *method is not a \c same method, this function returns zero
240 _gmx_selelem_custom_init_same(gmx_ana_selmethod_t **method,
241 const gmx::SelectionParserParameterListPointer ¶ms,
245 /* Do nothing if this is not a same method. */
246 if (!*method || (*method)->name != sm_same.name || params->empty())
251 const gmx::SelectionParserValueList &kwvalues = params->front().values();
252 if (kwvalues.size() != 1 || !kwvalues.front().hasExpressionValue()
253 || kwvalues.front().expr->type != SEL_EXPRESSION)
255 _gmx_selparser_error(scanner, "'same' should be followed by a single keyword");
258 gmx_ana_selmethod_t *kwmethod = kwvalues.front().expr->u.expr.method;
259 if (kwmethod->type == STR_VALUE)
261 *method = &sm_same_str;
264 /* We do custom processing for the "as" parameter. */
265 gmx::SelectionParserParameterList::iterator asparam = ++params->begin();
266 if (asparam != params->end() && asparam->name() == sm_same.param[1].name)
268 gmx::SelectionParserParameterList kwparams;
269 gmx::SelectionParserValueListPointer values(
270 new gmx::SelectionParserValueList(asparam->values()));
272 gmx::SelectionParserParameter::create(NULL, move(values)));
274 /* Create a second keyword evaluation element for the keyword given as
275 * the first parameter, evaluating the keyword in the group given by the
276 * second parameter. */
277 gmx::SelectionTreeElementPointer kwelem
278 = _gmx_sel_init_keyword_evaluator(kwmethod, kwparams, scanner);
279 // FIXME: Use exceptions.
284 /* Replace the second parameter with one with a value from \p kwelem. */
285 std::string pname = asparam->name();
286 *asparam = gmx::SelectionParserParameter::createFromExpression(pname, kwelem);
292 init_same(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
294 t_methoddata_same *d = (t_methoddata_same *)data;
296 d->val.ptr = param[0].val.u.ptr;
297 d->as.ptr = param[1].val.u.ptr;
298 if (param[1].val.type == STR_VALUE)
300 snew(d->as_s_sorted, d->nas);
302 if (!(param[0].flags & SPAR_ATOMVAL))
304 GMX_THROW(gmx::InvalidInputError(
305 "The 'same' selection keyword combined with a "
306 "non-keyword does not make sense"));
311 * \param data Data to free (should point to a \ref t_methoddata_same).
314 free_data_same(void *data)
316 t_methoddata_same *d = (t_methoddata_same *)data;
318 sfree(d->as_s_sorted);
323 * Helper function for comparison of two integers.
326 cmp_int(const void *a, const void *b)
328 if (*(int *)a < *(int *)b)
332 if (*(int *)a > *(int *)b)
340 init_frame_same_int(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */, void *data)
342 t_methoddata_same *d = (t_methoddata_same *)data;
345 /* Collapse adjacent values, and check whether the array is sorted. */
347 for (i = 1, j = 0; i < d->nas; ++i)
349 if (d->as.i[i] != d->as.i[j])
351 if (d->as.i[i] < d->as.i[j])
356 d->as.i[j] = d->as.i[i];
363 qsort(d->as.i, d->nas, sizeof(d->as.i[0]), &cmp_int);
364 /* More identical values may become adjacent after sorting. */
365 for (i = 1, j = 0; i < d->nas; ++i)
367 if (d->as.i[i] != d->as.i[j])
370 d->as.i[j] = d->as.i[i];
378 * See sel_updatefunc() for description of the parameters.
379 * \p data should point to a \c t_methoddata_same.
381 * Calculates which values in \c data->val.i can be found in \c data->as.i
382 * (assumed sorted), and writes the corresponding atoms to output.
383 * If \c data->val is sorted, uses a linear scan of both arrays, otherwise a
384 * binary search of \c data->as is performed for each block of values in
388 evaluate_same_int(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
389 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
391 t_methoddata_same *d = (t_methoddata_same *)data;
400 /* If we are sorted, we can do a simple linear scan. */
401 while (i < d->nas && d->as.i[i] < d->val.i[j])
408 /* If not, we must do a binary search of all the values. */
415 int itry = (i1 + i2) / 2;
416 if (d->as.i[itry] <= d->val.i[j])
425 i = (d->as.i[i1] == d->val.i[j] ? i1 : d->nas);
427 /* Check whether the value was found in the as list. */
428 if (i == d->nas || d->as.i[i] != d->val.i[j])
430 /* If not, skip all atoms with the same value. */
431 int tmpval = d->val.i[j];
433 while (j < g->isize && d->val.i[j] == tmpval)
440 /* Copy all the atoms with this value to the output. */
441 while (j < g->isize && d->val.i[j] == d->as.i[i])
443 out->u.g->index[out->u.g->isize++] = g->index[j];
447 if (j < g->isize && d->val.i[j] < d->val.i[j - 1])
455 * Helper function for comparison of two strings.
458 cmp_str(const void *a, const void *b)
460 return strcmp(*(char **)a, *(char **)b);
464 init_frame_same_str(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */, void *data)
466 t_methoddata_same *d = (t_methoddata_same *)data;
469 /* Collapse adjacent values.
470 * For strings, it's unlikely that the values would be sorted originally,
471 * so set bSorted always to false. */
473 d->as_s_sorted[0] = d->as.s[0];
474 for (i = 1, j = 0; i < d->nas; ++i)
476 if (strcmp(d->as.s[i], d->as_s_sorted[j]) != 0)
479 d->as_s_sorted[j] = d->as.s[i];
484 qsort(d->as_s_sorted, d->nas, sizeof(d->as_s_sorted[0]), &cmp_str);
485 /* More identical values may become adjacent after sorting. */
486 for (i = 1, j = 0; i < d->nas; ++i)
488 if (strcmp(d->as_s_sorted[i], d->as_s_sorted[j]) != 0)
491 d->as_s_sorted[j] = d->as_s_sorted[i];
498 * See sel_updatefunc() for description of the parameters.
499 * \p data should point to a \c t_methoddata_same.
501 * Calculates which strings in \c data->val.s can be found in \c data->as.s
502 * (assumed sorted), and writes the corresponding atoms to output.
503 * A binary search of \c data->as is performed for each block of values in
507 evaluate_same_str(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
508 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
510 t_methoddata_same *d = (t_methoddata_same *)data;
517 /* Do a binary search of the strings. */
519 ptr = bsearch(&d->val.s[j], d->as_s_sorted, d->nas,
520 sizeof(d->as_s_sorted[0]), &cmp_str);
521 /* Check whether the value was found in the as list. */
524 /* If not, skip all atoms with the same value. */
525 const char *tmpval = d->val.s[j];
527 while (j < g->isize && strcmp(d->val.s[j], tmpval) == 0)
534 const char *tmpval = d->val.s[j];
535 /* Copy all the atoms with this value to the output. */
536 while (j < g->isize && strcmp(d->val.s[j], tmpval) == 0)
538 out->u.g->index[out->u.g->isize++] = g->index[j];