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 internal selection methods for numeric and string keyword
40 * \author Teemu Murtola <teemu.murtola@gmail.com>
41 * \ingroup module_selection
48 #include <boost/shared_ptr.hpp>
50 #include "gromacs/legacyheaders/macros.h"
51 #include "gromacs/legacyheaders/smalloc.h"
52 #include "gromacs/legacyheaders/string2.h"
54 #include "gromacs/selection/selmethod.h"
55 #include "gromacs/utility/exceptions.h"
56 #include "gromacs/utility/gmxregex.h"
57 #include "gromacs/utility/messagestringcollector.h"
58 #include "gromacs/utility/stringutil.h"
61 #include "parsetree.h"
66 * Allocates data for integer keyword evaluation.
68 * \param[in] npar Not used.
69 * \param param Not used.
70 * \returns Pointer to the allocated data (\ref t_methoddata_kwint).
72 * Allocates memory for a \ref t_methoddata_kwint structure.
75 init_data_kwint(int npar, gmx_ana_selparam_t * param);
77 * Allocates data for real keyword evaluation.
79 * \param[in] npar Not used.
80 * \param param Not used.
81 * \returns Pointer to the allocated data (\ref t_methoddata_kwreal).
83 * Allocates memory for a \ref t_methoddata_kwreal structure.
86 init_data_kwreal(int npar, gmx_ana_selparam_t * param);
88 * Allocates data for string keyword evaluation.
90 * \param[in] npar Not used.
91 * \param param Not used.
92 * \returns Pointer to the allocated data (t_methoddata_kwstr).
94 * Allocates memory for a t_methoddata_kwstr structure.
97 init_data_kwstr(int npar, gmx_ana_selparam_t * param);
98 /** /brief Initializes data for integer keyword evaluation.
100 * \param[in] top Not used.
101 * \param[in] npar Not used (should be 2).
102 * \param[in] param Method parameters (should point to \ref smparams_keyword_int).
103 * \param[in] data Should point to \ref t_methoddata_kwint.
106 init_kwint(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
108 * Initializes data for real keyword evaluation.
110 * \param[in] top Not used.
111 * \param[in] npar Not used (should be 2).
112 * \param[in] param Method parameters (should point to \ref smparams_keyword_real).
113 * \param[in] data Should point to \ref t_methoddata_kwreal.
114 * \returns 0 (the initialization always succeeds).
117 init_kwreal(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
119 * Initializes data for string keyword evaluation.
121 * \param[in] top Not used.
122 * \param[in] npar Not used (should be 2).
123 * \param[in] param Method parameters (should point to \ref smparams_keyword_str).
124 * \param[in] data Should point to t_methoddata_kwstr.
127 init_kwstr(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
128 /** Frees the memory allocated for string keyword evaluation. */
130 free_data_kwstr(void *data);
131 /** Evaluates integer selection keywords. */
133 evaluate_keyword_int(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
134 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
135 /** Evaluates real selection keywords. */
137 evaluate_keyword_real(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
138 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
139 /** Evaluates string selection keywords. */
141 evaluate_keyword_str(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
142 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
145 * Data structure for integer keyword expression evaluation.
147 typedef struct t_methoddata_kwint
149 /** Array of values for the keyword. */
151 /** Number of ranges in the \p r array. */
154 * Array of sorted integer ranges to match against.
156 * Each range is made of two integers, giving the endpoints (inclusive).
157 * This field stores the pointer to the ranges allocated by the
158 * parameter parser; see \ref SPAR_RANGES for more information.
161 } t_methoddata_kwint;
164 * Data structure for real keyword expression evaluation.
166 typedef struct t_methoddata_kwreal
168 /** Array of values for the keyword. */
170 /** Number of ranges in the \p r array. */
173 * Array of sorted ranges to match against.
175 * Each range is made of two values, giving the endpoints (inclusive).
176 * This field stores the pointer to the ranges allocated by the
177 * parameter parser; see \ref SPAR_RANGES for more information.
180 } t_methoddata_kwreal;
186 * Single item in the list of strings/regular expressions to match.
188 * \ingroup module_selection
190 class StringKeywordMatchItem
194 * Constructs a matcher from a string.
196 * \param[in] matchType String matching type.
197 * \param[in] str String to use for matching.
199 StringKeywordMatchItem(gmx::SelectionStringMatchType matchType,
203 bool bRegExp = (matchType == gmx::eStringMatchType_RegularExpression);
204 if (matchType == gmx::eStringMatchType_Auto)
206 for (size_t j = 0; j < std::strlen(str); ++j)
208 if (std::ispunct(str[j]) && str[j] != '?' && str[j] != '*')
217 if (!gmx::Regex::isSupported())
219 GMX_THROW(gmx::InvalidInputError(gmx::formatString(
220 "No regular expression support, "
221 "cannot match \"%s\"", str)));
223 regex_.reset(new gmx::Regex(str));
228 * Checks whether this item matches a string.
230 * \param[in] matchType String matching type.
231 * \param[in] value String to match.
232 * \returns true if this item matches \p value.
234 bool match(gmx::SelectionStringMatchType matchType,
235 const char *value) const
237 if (matchType == gmx::eStringMatchType_Exact)
239 return str_ == value;
243 return gmx::regexMatch(value, *regex_);
247 return gmx_wcmatch(str_.c_str(), value) == 0;
252 //! The raw string passed for the matcher.
254 //! Regular expression compiled from \p str_, if applicable.
255 boost::shared_ptr<gmx::Regex> regex_;
259 * Data structure for string keyword expression evaluation.
261 struct t_methoddata_kwstr
263 /** Matching type for the strings. */
264 gmx::SelectionStringMatchType matchType;
265 /** Array of values for the keyword. */
267 /** Array of strings/regular expressions to match against.*/
268 std::vector<StringKeywordMatchItem> matches;
273 /** Parameters for integer keyword evaluation. */
274 static gmx_ana_selparam_t smparams_keyword_int[] = {
275 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
276 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
279 /** Parameters for real keyword evaluation. */
280 static gmx_ana_selparam_t smparams_keyword_real[] = {
281 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL | SPAR_DYNAMIC},
282 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
285 /** Parameters for string keyword evaluation. */
286 static gmx_ana_selparam_t smparams_keyword_str[] = {
287 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
288 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_VARNUM},
291 /** Selection method data for integer keyword evaluation. */
292 gmx_ana_selmethod_t sm_keyword_int = {
293 "kw_int", GROUP_VALUE, SMETH_SINGLEVAL,
294 asize(smparams_keyword_int), smparams_keyword_int,
301 &evaluate_keyword_int,
306 /** Selection method data for real keyword evaluation. */
307 gmx_ana_selmethod_t sm_keyword_real = {
308 "kw_real", GROUP_VALUE, SMETH_SINGLEVAL,
309 asize(smparams_keyword_real), smparams_keyword_real,
316 &evaluate_keyword_real,
321 /** Selection method data for string keyword evaluation. */
322 gmx_ana_selmethod_t sm_keyword_str = {
323 "kw_str", GROUP_VALUE, SMETH_SINGLEVAL,
324 asize(smparams_keyword_str), smparams_keyword_str,
331 &evaluate_keyword_str,
337 * Initializes keyword evaluation for an arbitrary group.
339 * \param[in] top Not used.
340 * \param[in] npar Not used.
341 * \param[in] param Not used.
342 * \param[in] data Should point to \ref t_methoddata_kweval.
343 * \returns 0 on success, a non-zero error code on return.
345 * Calls the initialization method of the wrapped keyword.
348 init_kweval(t_topology *top, int npar, gmx_ana_selparam_t * param, void *data);
350 * Initializes output for keyword evaluation in an arbitrary group.
352 * \param[in] top Not used.
353 * \param[in,out] out Pointer to output data structure.
354 * \param[in,out] data Should point to \c t_methoddata_kweval.
355 * \returns 0 for success.
358 init_output_kweval(t_topology *top, gmx_ana_selvalue_t *out, void *data);
359 /** Frees the data allocated for keyword evaluation in an arbitrary group. */
361 free_data_kweval(void *data);
362 /** Initializes frame evaluation for keyword evaluation in an arbitrary group. */
364 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
365 /** Evaluates keywords in an arbitrary group. */
367 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, gmx_ana_index_t * /* g */, gmx_ana_selvalue_t *out, void *data);
370 * Data structure for keyword evaluation in arbitrary groups.
374 /** Wrapped keyword method for evaluating the values. */
375 gmx_ana_selmethod_t *kwmethod;
376 /** Method data for \p kwmethod. */
378 /** Group in which \p kwmethod should be evaluated. */
380 } t_methoddata_kweval;
382 /** Parameters for keyword evaluation in an arbitrary group. */
383 static gmx_ana_selparam_t smparams_kweval[] = {
384 {NULL, {GROUP_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
388 /********************************************************************
389 * INTEGER KEYWORD EVALUATION
390 ********************************************************************/
393 init_data_kwint(int /* npar */, gmx_ana_selparam_t * /* param */)
395 t_methoddata_kwint *data;
402 init_kwint(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
404 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
406 d->v = param[0].val.u.i;
407 d->n = param[1].val.nr;
408 d->r = param[1].val.u.i;
412 * See sel_updatefunc() for description of the parameters.
413 * \p data should point to a \c t_methoddata_kwint.
415 * Does a binary search to find which atoms match the ranges in the
416 * \c t_methoddata_kwint structure for this selection.
417 * Matching atoms are stored in \p out->u.g.
420 evaluate_keyword_int(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
421 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
423 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
424 int n, i, j, jmin, jmax;
429 for (i = 0; i < g->isize; ++i)
432 if (d->r[0] > val || d->r[2*n-1] < val)
438 while (jmax - jmin > 1)
440 j = jmin + (jmax - jmin) / 2;
448 if (val <= d->r[2*j+1])
455 if (val <= d->r[2*jmin+1])
457 out->u.g->index[out->u.g->isize++] = g->index[i];
463 /********************************************************************
464 * REAL KEYWORD EVALUATION
465 ********************************************************************/
468 init_data_kwreal(int /* npar */, gmx_ana_selparam_t * /* param */)
470 t_methoddata_kwreal *data;
477 init_kwreal(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
479 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
481 d->v = param[0].val.u.r;
482 d->n = param[1].val.nr;
483 d->r = param[1].val.u.r;
487 * See sel_updatefunc() for description of the parameters.
488 * \p data should point to a \c t_methoddata_kwreal.
490 * Does a binary search to find which atoms match the ranges in the
491 * \c t_methoddata_kwreal structure for this selection.
492 * Matching atoms are stored in \p out->u.g.
495 evaluate_keyword_real(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
496 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
498 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
499 int n, i, j, jmin, jmax;
504 for (i = 0; i < g->isize; ++i)
507 if (d->r[0] > val || d->r[2*n-1] < val)
513 while (jmax - jmin > 1)
515 j = jmin + (jmax - jmin) / 2;
523 if (val <= d->r[2*j+1])
530 if (val <= d->r[2*jmin+1])
532 out->u.g->index[out->u.g->isize++] = g->index[i];
538 /********************************************************************
539 * STRING KEYWORD EVALUATION
540 ********************************************************************/
543 init_data_kwstr(int /* npar */, gmx_ana_selparam_t * /* param */)
545 t_methoddata_kwstr *data = new t_methoddata_kwstr();
546 data->matchType = gmx::eStringMatchType_Auto;
551 * \param[in,out] sel Selection element to initialize.
552 * \param[in] matchType Method to use to match string values.
554 * Sets the string matching method for string keyword matching.
557 _gmx_selelem_set_kwstr_match_type(const gmx::SelectionTreeElementPointer &sel,
558 gmx::SelectionStringMatchType matchType)
560 t_methoddata_kwstr *d = static_cast<t_methoddata_kwstr *>(sel->u.expr.mdata);
562 if (sel->type != SEL_EXPRESSION || !sel->u.expr.method
563 || sel->u.expr.method->name != sm_keyword_str.name)
567 d->matchType = matchType;
571 init_kwstr(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
573 t_methoddata_kwstr *d = static_cast<t_methoddata_kwstr *>(data);
575 d->v = param[0].val.u.s;
576 /* Return if this is not the first time */
577 if (!d->matches.empty())
581 int n = param[1].val.nr;
582 d->matches.reserve(n);
583 for (int i = 0; i < n; ++i)
585 const char *s = param[1].val.u.s[i];
586 d->matches.push_back(StringKeywordMatchItem(d->matchType, s));
591 * \param data Data to free (should point to a t_methoddata_kwstr).
594 free_data_kwstr(void *data)
596 t_methoddata_kwstr *d = static_cast<t_methoddata_kwstr *>(data);
601 * See sel_updatefunc() for description of the parameters.
602 * \p data should point to a \c t_methoddata_kwstr.
604 * Does a linear search to find which atoms match the strings in the
605 * \c t_methoddata_kwstr structure for this selection.
606 * Wildcards are allowed in the strings.
607 * Matching atoms are stored in \p out->u.g.
610 evaluate_keyword_str(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
611 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
613 t_methoddata_kwstr *d = static_cast<t_methoddata_kwstr *>(data);
616 for (int i = 0; i < g->isize; ++i)
618 for (size_t j = 0; j < d->matches.size(); ++j)
620 if (d->matches[j].match(d->matchType, d->v[i]))
622 out->u.g->index[out->u.g->isize++] = g->index[i];
630 /********************************************************************
631 * KEYWORD EVALUATION FOR ARBITRARY GROUPS
632 ********************************************************************/
635 init_kweval(t_topology *top, int /* npar */, gmx_ana_selparam_t * /* param */, void *data)
637 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
639 d->kwmethod->init(top, 0, NULL, d->kwmdata);
643 init_output_kweval(t_topology * /* top */, gmx_ana_selvalue_t *out, void *data)
645 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
647 out->nr = d->g.isize;
648 _gmx_selvalue_reserve(out, out->nr);
652 * \param data Data to free (should point to a \c t_methoddata_kweval).
654 * Frees the memory allocated for all the members of \c t_methoddata_kweval.
657 free_data_kweval(void *data)
659 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
661 _gmx_selelem_free_method(d->kwmethod, d->kwmdata);
666 * \param[in] top Topology.
667 * \param[in] fr Current frame.
668 * \param[in] pbc PBC structure.
669 * \param data Should point to a \ref t_methoddata_kweval.
670 * \returns 0 on success, a non-zero error code on error.
672 * Creates a lookup structure that enables fast queries of whether a point
673 * is within the solid angle or not.
676 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
678 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
680 d->kwmethod->init_frame(top, fr, pbc, d->kwmdata);
684 * See sel_updatefunc() for description of the parameters.
685 * \p data should point to a \c t_methoddata_kweval.
687 * Calls the evaluation function of the wrapped keyword with the given
688 * parameters, with the exception of using \c t_methoddata_kweval::g for the
692 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc,
693 gmx_ana_index_t * /* g */, gmx_ana_selvalue_t *out, void *data)
695 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
697 d->kwmethod->update(top, fr, pbc, &d->g, out, d->kwmdata);
701 * \param[in] method Keyword selection method to evaluate.
702 * \param[in] params Parameter that gives the group to evaluate \p method in.
703 * \param[in] scanner Scanner data structure.
704 * \returns Pointer to the created selection element (NULL on error).
706 * Creates a \ref SEL_EXPRESSION selection element that evaluates the keyword
707 * method given by \p method in the group given by \p param.
709 * The name of \p param should be empty.
711 gmx::SelectionTreeElementPointer
712 _gmx_sel_init_keyword_evaluator(gmx_ana_selmethod_t *method,
713 const gmx::SelectionParserParameterList ¶ms,
716 gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
718 sprintf(buf, "In evaluation of '%s'", method->name);
719 gmx::MessageStringContext context(errors, buf);
721 if ((method->flags & (SMETH_SINGLEVAL | SMETH_VARNUMVAL))
722 || method->outinit || method->pupdate)
724 GMX_THROW(gmx::InternalError(
725 "Unsupported keyword method for arbitrary group evaluation"));
728 gmx::SelectionTreeElementPointer sel(
729 new gmx::SelectionTreeElement(SEL_EXPRESSION));
730 _gmx_selelem_set_method(sel, method, scanner);
732 t_methoddata_kweval *data;
734 data->kwmethod = sel->u.expr.method;
735 data->kwmdata = sel->u.expr.mdata;
736 gmx_ana_index_clear(&data->g);
738 snew(sel->u.expr.method, 1);
739 memcpy(sel->u.expr.method, data->kwmethod, sizeof(gmx_ana_selmethod_t));
740 sel->u.expr.method->flags |= SMETH_VARNUMVAL;
741 sel->u.expr.method->init_data = NULL;
742 sel->u.expr.method->set_poscoll = NULL;
743 sel->u.expr.method->init = method->init ? &init_kweval : NULL;
744 sel->u.expr.method->outinit = &init_output_kweval;
745 sel->u.expr.method->free = &free_data_kweval;
746 sel->u.expr.method->init_frame = method->init_frame ? &init_frame_kweval : NULL;
747 sel->u.expr.method->update = &evaluate_kweval;
748 sel->u.expr.method->pupdate = NULL;
749 sel->u.expr.method->nparams = asize(smparams_kweval);
750 sel->u.expr.method->param = smparams_kweval;
751 _gmx_selelem_init_method_params(sel, scanner);
752 sel->u.expr.mdata = data;
754 sel->u.expr.method->param[0].val.u.g = &data->g;
756 if (!_gmx_sel_parse_params(params, sel->u.expr.method->nparams,
757 sel->u.expr.method->param, sel, scanner))
759 return gmx::SelectionTreeElementPointer();