2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2016,2017, 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
47 #include "gromacs/utility/arraysize.h"
48 #include "gromacs/utility/exceptions.h"
49 #include "gromacs/utility/smalloc.h"
52 #include "parsetree.h"
54 #include "selmethod.h"
58 * Data structure for the \p same selection method.
60 * To avoid duplicate initialization code, the same data structure is used
61 * for matching both integer and string keywords; hence the unions.
63 * \ingroup module_selection
67 /** Value for each atom to match. */
75 * Number of values in the \p as array.
77 * For string values, this is actually the number of values in the
78 * \p as_s_sorted array.
81 /** Values to match against. */
89 * Separate array for sorted \p as.s array.
91 * The array of strings returned as the output value of a parameter should
92 * not be messed with to avoid memory corruption (the pointers in the array
93 * may be reused for several evaluations), so we keep our own copy for
97 /** Whether simple matching can be used. */
102 * Allocates data for the \p same selection method.
104 * \param[in] npar Not used (should be 2).
105 * \param[in,out] param Method parameters (should point to a copy of
106 * ::smparams_same_int or ::smparams_same_str).
107 * \returns Pointer to the allocated data (\ref t_methoddata_same).
110 init_data_same(int npar
, gmx_ana_selparam_t
*param
);
112 * Initializes the \p same selection method.
114 * \param top Not used.
115 * \param npar Not used (should be 2).
116 * \param param Initialized method parameters (should point to a copy of
117 * ::smparams_same_int or ::smparams_same_str).
118 * \param data Pointer to \ref t_methoddata_same to initialize.
119 * \returns 0 on success, -1 on failure.
122 init_same(const gmx_mtop_t
*top
, int npar
, gmx_ana_selparam_t
*param
, void *data
);
123 /** Frees the data allocated for the \p same selection method. */
125 free_data_same(void *data
);
127 * Initializes the evaluation of the \p same selection method for a frame.
129 * \param[in] context Not used.
130 * \param data Should point to a \ref t_methoddata_same.
132 * Sorts the \c data->as.i array and removes identical values for faster and
136 init_frame_same_int(const gmx::SelMethodEvalContext
&context
, void *data
);
137 /** Evaluates the \p same selection method. */
139 evaluate_same_int(const gmx::SelMethodEvalContext
&context
,
140 gmx_ana_index_t
*g
, gmx_ana_selvalue_t
*out
, void *data
);
142 * Initializes the evaluation of the \p same selection method for a frame.
144 * \param[in] context Not used.
145 * \param data Should point to a \ref t_methoddata_same.
147 * Sorts the \c data->as.s array and removes identical values for faster and
151 init_frame_same_str(const gmx::SelMethodEvalContext
&context
, void *data
);
152 /** Evaluates the \p same selection method. */
154 evaluate_same_str(const gmx::SelMethodEvalContext
&context
,
155 gmx_ana_index_t
*g
, gmx_ana_selvalue_t
*out
, void *data
);
157 /** Parameters for the \p same selection method. */
158 static gmx_ana_selparam_t smparams_same_int
[] = {
159 {nullptr, {INT_VALUE
, -1, {nullptr}}, nullptr, SPAR_DYNAMIC
| SPAR_ATOMVAL
},
160 {"as", {INT_VALUE
, -1, {nullptr}}, nullptr, SPAR_DYNAMIC
| SPAR_VARNUM
},
163 /** Parameters for the \p same selection method. */
164 static gmx_ana_selparam_t smparams_same_str
[] = {
165 {nullptr, {STR_VALUE
, -1, {nullptr}}, nullptr, SPAR_DYNAMIC
| SPAR_ATOMVAL
},
166 {"as", {STR_VALUE
, -1, {nullptr}}, nullptr, SPAR_DYNAMIC
| SPAR_VARNUM
},
169 /** Help text for the \p same selection method. */
170 static const char *const help_same
[] = {
173 " same KEYWORD as ATOM_EXPR",
176 "The keyword [TT]same[tt] can be used to select all atoms for which",
177 "the given [TT]KEYWORD[tt] matches any of the atoms in [TT]ATOM_EXPR[tt].",
178 "Keywords that evaluate to integer or string values are supported.",
181 /*! \internal \brief Selection method data for the \p same method. */
182 gmx_ana_selmethod_t sm_same
= {
183 "same", GROUP_VALUE
, 0,
184 asize(smparams_same_int
), smparams_same_int
,
190 &init_frame_same_int
,
193 {"same KEYWORD as ATOM_EXPR",
194 "Extending selections", asize(help_same
), help_same
},
198 * Selection method data for the \p same method.
200 * This selection method is used for matching string keywords. The parser
201 * never sees this method; _gmx_selelem_custom_init_same() replaces sm_same
202 * with this method in cases where it is required.
204 static gmx_ana_selmethod_t sm_same_str
= {
205 "same", GROUP_VALUE
, SMETH_SINGLEVAL
,
206 asize(smparams_same_str
), smparams_same_str
,
212 &init_frame_same_str
,
215 {nullptr, nullptr, 0, nullptr},
219 init_data_same(int /* npar */, gmx_ana_selparam_t
*param
)
221 t_methoddata_same
*data
;
224 data
->as_s_sorted
= nullptr;
225 param
[1].nvalptr
= &data
->nas
;
230 * \param[in,out] method The method to initialize.
231 * \param[in,out] params Pointer to the first parameter.
232 * \param[in] scanner Scanner data structure.
234 * If \p *method is not a \c same method, this function returns
238 _gmx_selelem_custom_init_same(gmx_ana_selmethod_t
**method
,
239 const gmx::SelectionParserParameterListPointer
¶ms
,
243 /* Do nothing if this is not a same method. */
244 if (!*method
|| (*method
)->name
!= sm_same
.name
|| params
->empty())
249 const gmx::SelectionParserValueList
&kwvalues
= params
->front().values();
250 if (kwvalues
.size() != 1 || !kwvalues
.front().hasExpressionValue()
251 || kwvalues
.front().expr
->type
!= SEL_EXPRESSION
)
253 GMX_THROW(gmx::InvalidInputError("'same' should be followed by a single keyword"));
255 gmx_ana_selmethod_t
*kwmethod
= kwvalues
.front().expr
->u
.expr
.method
;
256 if (kwmethod
->type
== STR_VALUE
)
258 *method
= &sm_same_str
;
261 /* We do custom processing for the "as" parameter. */
262 gmx::SelectionParserParameterList::iterator asparam
= ++params
->begin();
263 if (asparam
!= params
->end() && asparam
->name() == sm_same
.param
[1].name
)
265 const gmx::SelectionParserValueList
&asvalues
= asparam
->values();
266 if (asvalues
.size() != 1 || !asvalues
.front().hasExpressionValue())
268 // TODO: Think about providing more informative context.
269 GMX_THROW(gmx::InvalidInputError("'same ... as' should be followed by a single expression"));
271 const gmx::SelectionTreeElementPointer
&child
= asvalues
.front().expr
;
272 /* Create a second keyword evaluation element for the keyword given as
273 * the first parameter, evaluating the keyword in the group given by the
274 * second parameter. */
275 gmx::SelectionTreeElementPointer kwelem
276 = _gmx_sel_init_keyword_evaluator(kwmethod
, child
, scanner
);
277 /* Replace the second parameter with one with a value from \p kwelem. */
278 std::string pname
= asparam
->name();
279 *asparam
= gmx::SelectionParserParameter::createFromExpression(pname
, kwelem
);
284 init_same(const gmx_mtop_t
* /* top */, int /* npar */, gmx_ana_selparam_t
*param
, void *data
)
286 t_methoddata_same
*d
= (t_methoddata_same
*)data
;
288 d
->val
.ptr
= param
[0].val
.u
.ptr
;
289 d
->as
.ptr
= param
[1].val
.u
.ptr
;
290 if (param
[1].val
.type
== STR_VALUE
)
292 snew(d
->as_s_sorted
, d
->nas
);
294 if (!(param
[0].flags
& SPAR_ATOMVAL
))
296 GMX_THROW(gmx::InvalidInputError(
297 "The 'same' selection keyword combined with a "
298 "non-keyword does not make sense"));
303 * \param data Data to free (should point to a \ref t_methoddata_same).
306 free_data_same(void *data
)
308 t_methoddata_same
*d
= (t_methoddata_same
*)data
;
310 sfree(d
->as_s_sorted
);
315 * Helper function for comparison of two integers.
318 cmp_int(const void *a
, const void *b
)
320 if (*(int *)a
< *(int *)b
)
324 if (*(int *)a
> *(int *)b
)
332 init_frame_same_int(const gmx::SelMethodEvalContext
& /*context*/, void *data
)
334 t_methoddata_same
*d
= (t_methoddata_same
*)data
;
337 /* Collapse adjacent values, and check whether the array is sorted. */
343 for (i
= 1, j
= 0; i
< d
->nas
; ++i
)
345 if (d
->as
.i
[i
] != d
->as
.i
[j
])
347 if (d
->as
.i
[i
] < d
->as
.i
[j
])
352 d
->as
.i
[j
] = d
->as
.i
[i
];
359 qsort(d
->as
.i
, d
->nas
, sizeof(d
->as
.i
[0]), &cmp_int
);
360 /* More identical values may become adjacent after sorting. */
361 for (i
= 1, j
= 0; i
< d
->nas
; ++i
)
363 if (d
->as
.i
[i
] != d
->as
.i
[j
])
366 d
->as
.i
[j
] = d
->as
.i
[i
];
374 * See sel_updatefunc() for description of the parameters.
375 * \p data should point to a \c t_methoddata_same.
377 * Calculates which values in \c data->val.i can be found in \c data->as.i
378 * (assumed sorted), and writes the corresponding atoms to output.
379 * If \c data->val is sorted, uses a linear scan of both arrays, otherwise a
380 * binary search of \c data->as is performed for each block of values in
384 evaluate_same_int(const gmx::SelMethodEvalContext
& /*context*/,
385 gmx_ana_index_t
*g
, gmx_ana_selvalue_t
*out
, void *data
)
387 t_methoddata_same
*d
= (t_methoddata_same
*)data
;
396 /* If we are sorted, we can do a simple linear scan. */
397 while (i
< d
->nas
&& d
->as
.i
[i
] < d
->val
.i
[j
])
404 /* If not, we must do a binary search of all the values. */
411 int itry
= (i1
+ i2
) / 2;
412 if (d
->as
.i
[itry
] <= d
->val
.i
[j
])
421 i
= (d
->as
.i
[i1
] == d
->val
.i
[j
] ? i1
: d
->nas
);
423 /* Check whether the value was found in the as list. */
424 if (i
== d
->nas
|| d
->as
.i
[i
] != d
->val
.i
[j
])
426 /* If not, skip all atoms with the same value. */
427 int tmpval
= d
->val
.i
[j
];
429 while (j
< g
->isize
&& d
->val
.i
[j
] == tmpval
)
436 /* Copy all the atoms with this value to the output. */
437 while (j
< g
->isize
&& d
->val
.i
[j
] == d
->as
.i
[i
])
439 out
->u
.g
->index
[out
->u
.g
->isize
++] = g
->index
[j
];
443 if (j
< g
->isize
&& d
->val
.i
[j
] < d
->val
.i
[j
- 1])
451 * Helper function for comparison of two strings.
454 cmp_str(const void *a
, const void *b
)
456 return strcmp(*(char **)a
, *(char **)b
);
460 init_frame_same_str(const gmx::SelMethodEvalContext
& /*context*/, 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. */
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(const gmx::SelMethodEvalContext
& /*context*/,
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
];