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
33 * Implements the \p same selection method.
35 * \author Teemu Murtola <teemu.murtola@cbr.su.se>
36 * \ingroup module_selection
48 #include "gromacs/fatalerror/exceptions.h"
49 #include "gromacs/selection/selmethod.h"
52 #include "parsetree.h"
56 * Data structure for the \p same selection method.
58 * To avoid duplicate initialization code, the same data structure is used
59 * for matching both integer and string keywords; hence the unions.
63 /** Value for each atom to match. */
71 * Number of values in the \p as array.
73 * For string values, this is actually the number of values in the
74 * \p as_s_sorted array.
77 /** Values to match against. */
85 * Separate array for sorted \p as.s array.
87 * The array of strings returned as the output value of a parameter should
88 * not be messed with to avoid memory corruption (the pointers in the array
89 * may be reused for several evaluations), so we keep our own copy for
93 /** Whether simple matching can be used. */
97 /** Allocates data for the \p same selection method. */
99 init_data_same(int npar, gmx_ana_selparam_t *param);
100 /** Initializes the \p same selection method. */
102 init_same(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
103 /** Frees the data allocated for the \p same selection method. */
105 free_data_same(void *data);
106 /** Initializes the evaluation of the \p same selection method for a frame. */
108 init_frame_same_int(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
109 /** Evaluates the \p same selection method. */
111 evaluate_same_int(t_topology *top, t_trxframe *fr, t_pbc *pbc,
112 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
113 /** Initializes the evaluation of the \p same selection method for a frame. */
115 init_frame_same_str(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
116 /** Evaluates the \p same selection method. */
118 evaluate_same_str(t_topology *top, t_trxframe *fr, t_pbc *pbc,
119 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
121 /** Parameters for the \p same selection method. */
122 static gmx_ana_selparam_t smparams_same_int[] = {
123 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_ATOMVAL},
124 {"as", {INT_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
127 /** Parameters for the \p same selection method. */
128 static gmx_ana_selparam_t smparams_same_str[] = {
129 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_ATOMVAL},
130 {"as", {STR_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
133 /** Help text for the \p same selection method. */
134 static const char *help_same[] = {
135 "EXTENDING SELECTIONS[PAR]",
137 "[TT]same KEYWORD as ATOM_EXPR[tt][PAR]",
139 "The keyword [TT]same[tt] can be used to select all atoms for which",
140 "the given [TT]KEYWORD[tt] matches any of the atoms in [TT]ATOM_EXPR[tt].",
141 "Keywords that evaluate to integer or string values are supported.",
144 /*! \internal \brief Selection method data for the \p same method. */
145 gmx_ana_selmethod_t sm_same = {
146 "same", GROUP_VALUE, 0,
147 asize(smparams_same_int), smparams_same_int,
153 &init_frame_same_int,
156 {"same KEYWORD as ATOM_EXPR", asize(help_same), help_same},
160 * Selection method data for the \p same method.
162 * This selection method is used for matching string keywords. The parser
163 * never sees this method; _gmx_selelem_custom_init_same() replaces sm_same
164 * with this method in cases where it is required.
166 static gmx_ana_selmethod_t sm_same_str = {
167 "same", GROUP_VALUE, SMETH_SINGLEVAL,
168 asize(smparams_same_str), smparams_same_str,
174 &init_frame_same_str,
177 {"same KEYWORD as ATOM_EXPR", asize(help_same), help_same},
181 * \param[in] npar Not used (should be 2).
182 * \param[in,out] param Method parameters (should point to
183 * \ref smparams_same).
184 * \returns Pointer to the allocated data (\ref t_methoddata_same).
187 init_data_same(int npar, gmx_ana_selparam_t *param)
189 t_methoddata_same *data;
192 data->as_s_sorted = NULL;
193 param[1].nvalptr = &data->nas;
198 * \param[in,out] method The method to initialize.
199 * \param[in,out] params Pointer to the first parameter.
200 * \param[in] scanner Scanner data structure.
201 * \returns 0 on success, a non-zero error code on error.
203 * If \p *method is not a \c same method, this function returns zero
207 _gmx_selelem_custom_init_same(gmx_ana_selmethod_t **method,
208 t_selexpr_param *params,
211 gmx_ana_selmethod_t *kwmethod;
213 t_selexpr_param *param;
217 /* Do nothing if this is not a same method. */
218 if (!*method || (*method)->name != sm_same.name)
223 if (params->nval != 1 || !params->value->bExpr
224 || params->value->u.expr->type != SEL_EXPRESSION)
226 _gmx_selparser_error(scanner, "'same' should be followed by a single keyword");
229 kwmethod = params->value->u.expr->u.expr.method;
231 if (kwmethod->type == STR_VALUE)
233 *method = &sm_same_str;
236 /* We do custom processing with the second parameter, so remove it from
237 * the params list, but save the name for later. */
238 param = params->next;
242 /* Create a second keyword evaluation element for the keyword given as
243 * the first parameter, evaluating the keyword in the group given by the
244 * second parameter. */
245 rc = _gmx_sel_init_keyword_evaluator(&kwelem, kwmethod, param, scanner);
251 /* Replace the second parameter with one with a value from \p kwelem. */
252 param = _gmx_selexpr_create_param(pname);
254 param->value = _gmx_selexpr_create_value_expr(kwelem);
255 params->next = param;
260 * \param top Not used.
261 * \param npar Not used (should be 2).
262 * \param param Initialized method parameters (should point to a copy of
263 * \ref smparams_same).
264 * \param data Pointer to \ref t_methoddata_same to initialize.
265 * \returns 0 on success, -1 on failure.
268 init_same(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
270 t_methoddata_same *d = (t_methoddata_same *)data;
272 d->val.ptr = param[0].val.u.ptr;
273 d->as.ptr = param[1].val.u.ptr;
274 if (param[1].val.type == STR_VALUE)
276 snew(d->as_s_sorted, d->nas);
278 if (!(param[0].flags & SPAR_ATOMVAL))
280 GMX_THROW(gmx::InvalidInputError(
281 "The 'same' selection keyword combined with a "
282 "non-keyword does not make sense"));
287 * \param data Data to free (should point to a \ref t_methoddata_same).
290 free_data_same(void *data)
292 t_methoddata_same *d = (t_methoddata_same *)data;
294 sfree(d->as_s_sorted);
298 * Helper function for comparison of two integers.
301 cmp_int(const void *a, const void *b)
303 if (*(int *)a < *(int *)b)
307 if (*(int *)a > *(int *)b)
315 * \param[in] top Not used.
316 * \param[in] fr Current frame.
317 * \param[in] pbc PBC structure.
318 * \param data Should point to a \ref t_methoddata_same.
320 * Sorts the \c data->as.i array and removes identical values for faster and
324 init_frame_same_int(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
326 t_methoddata_same *d = (t_methoddata_same *)data;
329 /* Collapse adjacent values, and check whether the array is sorted. */
331 for (i = 1, j = 0; i < d->nas; ++i)
333 if (d->as.i[i] != d->as.i[j])
335 if (d->as.i[i] < d->as.i[j])
340 d->as.i[j] = d->as.i[i];
347 qsort(d->as.i, d->nas, sizeof(d->as.i[0]), &cmp_int);
348 /* More identical values may become adjacent after sorting. */
349 for (i = 1, j = 0; i < d->nas; ++i)
351 if (d->as.i[i] != d->as.i[j])
354 d->as.i[j] = d->as.i[i];
362 * See sel_updatefunc() for description of the parameters.
363 * \p data should point to a \c t_methoddata_same.
365 * Calculates which values in \c data->val.i can be found in \c data->as.i
366 * (assumed sorted), and writes the corresponding atoms to output.
367 * If \c data->val is sorted, uses a linear scan of both arrays, otherwise a
368 * binary search of \c data->as is performed for each block of values in
372 evaluate_same_int(t_topology *top, t_trxframe *fr, t_pbc *pbc,
373 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
375 t_methoddata_same *d = (t_methoddata_same *)data;
384 /* If we are sorted, we can do a simple linear scan. */
385 while (i < d->nas && d->as.i[i] < d->val.i[j]) ++i;
389 /* If not, we must do a binary search of all the values. */
396 int itry = (i1 + i2) / 2;
397 if (d->as.i[itry] <= d->val.i[j])
406 i = (d->as.i[i1] == d->val.i[j] ? i1 : d->nas);
408 /* Check whether the value was found in the as list. */
409 if (i == d->nas || d->as.i[i] != d->val.i[j])
411 /* If not, skip all atoms with the same value. */
412 int tmpval = d->val.i[j];
414 while (j < g->isize && d->val.i[j] == tmpval) ++j;
418 /* Copy all the atoms with this value to the output. */
419 while (j < g->isize && d->val.i[j] == d->as.i[i])
421 out->u.g->index[out->u.g->isize++] = g->index[j];
425 if (j < g->isize && d->val.i[j] < d->val.i[j - 1])
433 * Helper function for comparison of two strings.
436 cmp_str(const void *a, const void *b)
438 return strcmp(*(char **)a, *(char **)b);
442 * \param[in] top Not used.
443 * \param[in] fr Current frame.
444 * \param[in] pbc PBC structure.
445 * \param data Should point to a \ref t_methoddata_same.
447 * Sorts the \c data->as.s array and removes identical values for faster and
451 init_frame_same_str(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
453 t_methoddata_same *d = (t_methoddata_same *)data;
456 /* Collapse adjacent values.
457 * For strings, it's unlikely that the values would be sorted originally,
458 * so set bSorted always to false. */
460 d->as_s_sorted[0] = d->as.s[0];
461 for (i = 1, j = 0; i < d->nas; ++i)
463 if (strcmp(d->as.s[i], d->as_s_sorted[j]) != 0)
466 d->as_s_sorted[j] = d->as.s[i];
471 qsort(d->as_s_sorted, d->nas, sizeof(d->as_s_sorted[0]), &cmp_str);
472 /* More identical values may become adjacent after sorting. */
473 for (i = 1, j = 0; i < d->nas; ++i)
475 if (strcmp(d->as_s_sorted[i], d->as_s_sorted[j]) != 0)
478 d->as_s_sorted[j] = d->as_s_sorted[i];
485 * See sel_updatefunc() for description of the parameters.
486 * \p data should point to a \c t_methoddata_same.
488 * Calculates which strings in \c data->val.s can be found in \c data->as.s
489 * (assumed sorted), and writes the corresponding atoms to output.
490 * A binary search of \c data->as is performed for each block of values in
494 evaluate_same_str(t_topology *top, t_trxframe *fr, t_pbc *pbc,
495 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
497 t_methoddata_same *d = (t_methoddata_same *)data;
504 /* Do a binary search of the strings. */
506 ptr = bsearch(&d->val.s[j], d->as_s_sorted, d->nas,
507 sizeof(d->as_s_sorted[0]), &cmp_str);
508 /* Check whether the value was found in the as list. */
511 /* If not, skip all atoms with the same value. */
512 const char *tmpval = d->val.s[j];
514 while (j < g->isize && strcmp(d->val.s[j], tmpval) == 0) ++j;
518 const char *tmpval = d->val.s[j];
519 /* Copy all the atoms with this value to the output. */
520 while (j < g->isize && strcmp(d->val.s[j], tmpval) == 0)
522 out->u.g->index[out->u.g->isize++] = g->index[j];