2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015, 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
50 #include <boost/shared_ptr.hpp>
52 #include "gromacs/legacyheaders/macros.h"
53 #include "gromacs/selection/position.h"
54 #include "gromacs/utility/cstringutil.h"
55 #include "gromacs/utility/exceptions.h"
56 #include "gromacs/utility/gmxassert.h"
57 #include "gromacs/utility/gmxregex.h"
58 #include "gromacs/utility/smalloc.h"
59 #include "gromacs/utility/stringutil.h"
62 #include "parsetree.h"
65 #include "selmethod.h"
68 * Allocates data for integer keyword evaluation.
70 * \param[in] npar Not used.
71 * \param param Not used.
72 * \returns Pointer to the allocated data (\ref t_methoddata_kwint).
74 * Allocates memory for a \ref t_methoddata_kwint structure.
77 init_data_kwint(int npar, gmx_ana_selparam_t * param);
79 * Allocates data for real keyword evaluation.
81 * \param[in] npar Not used.
82 * \param param Not used.
83 * \returns Pointer to the allocated data (\ref t_methoddata_kwreal).
85 * Allocates memory for a \ref t_methoddata_kwreal structure.
88 init_data_kwreal(int npar, gmx_ana_selparam_t * param);
90 * Allocates data for string keyword evaluation.
92 * \param[in] npar Not used.
93 * \param param Not used.
94 * \returns Pointer to the allocated data (t_methoddata_kwstr).
96 * Allocates memory for a t_methoddata_kwstr structure.
99 init_data_kwstr(int npar, gmx_ana_selparam_t * param);
100 /** /brief Initializes data for integer keyword evaluation.
102 * \param[in] top Not used.
103 * \param[in] npar Not used (should be 2).
104 * \param[in] param Method parameters (should point to \ref smparams_keyword_int).
105 * \param[in] data Should point to \ref t_methoddata_kwint.
108 init_kwint(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
110 * Initializes data for real keyword evaluation.
112 * \param[in] top Not used.
113 * \param[in] npar Not used (should be 2).
114 * \param[in] param Method parameters (should point to \ref smparams_keyword_real).
115 * \param[in] data Should point to \ref t_methoddata_kwreal.
116 * \returns 0 (the initialization always succeeds).
119 init_kwreal(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
121 * Initializes data for string keyword evaluation.
123 * \param[in] top Not used.
124 * \param[in] npar Not used (should be 2).
125 * \param[in] param Method parameters (should point to \ref smparams_keyword_str).
126 * \param[in] data Should point to t_methoddata_kwstr.
129 init_kwstr(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
130 /** Frees the memory allocated for string keyword evaluation. */
132 free_data_kwstr(void *data);
133 /** Evaluates integer selection keywords. */
135 evaluate_keyword_int(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
136 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
137 /** Evaluates real selection keywords. */
139 evaluate_keyword_real(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
140 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
141 /** Evaluates string selection keywords. */
143 evaluate_keyword_str(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
144 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
147 * Data structure for integer keyword expression evaluation.
149 typedef struct t_methoddata_kwint
151 /** Array of values for the keyword. */
153 /** Number of ranges in the \p r array. */
156 * Array of sorted integer ranges to match against.
158 * Each range is made of two integers, giving the endpoints (inclusive).
159 * This field stores the pointer to the ranges allocated by the
160 * parameter parser; see \ref SPAR_RANGES for more information.
163 } t_methoddata_kwint;
166 * Data structure for real keyword expression evaluation.
168 typedef struct t_methoddata_kwreal
170 /** Array of values for the keyword. */
172 /** Number of ranges in the \p r array. */
175 * Array of sorted ranges to match against.
177 * Each range is made of two values, giving the endpoints (inclusive).
178 * This field stores the pointer to the ranges allocated by the
179 * parameter parser; see \ref SPAR_RANGES for more information.
182 } t_methoddata_kwreal;
188 * Single item in the list of strings/regular expressions to match.
190 * \ingroup module_selection
192 class StringKeywordMatchItem
196 * Constructs a matcher from a string.
198 * \param[in] matchType String matching type.
199 * \param[in] str String to use for matching.
201 StringKeywordMatchItem(gmx::SelectionStringMatchType matchType,
205 bool bRegExp = (matchType == gmx::eStringMatchType_RegularExpression);
206 if (matchType == gmx::eStringMatchType_Auto)
208 for (size_t j = 0; j < std::strlen(str); ++j)
210 if (std::ispunct(str[j]) && str[j] != '?' && str[j] != '*')
219 if (!gmx::Regex::isSupported())
222 = gmx::formatString("No regular expression support, "
223 "cannot match \"%s\"", str);
224 GMX_THROW(gmx::InvalidInputError(message));
226 regex_.reset(new gmx::Regex(str));
231 * Checks whether this item matches a string.
233 * \param[in] matchType String matching type.
234 * \param[in] value String to match.
235 * \returns true if this item matches \p value.
237 bool match(gmx::SelectionStringMatchType matchType,
238 const char *value) const
240 if (matchType == gmx::eStringMatchType_Exact)
242 return str_ == value;
246 return gmx::regexMatch(value, *regex_);
250 return gmx_wcmatch(str_.c_str(), value) == 0;
255 //! The raw string passed for the matcher.
257 //! Regular expression compiled from \p str_, if applicable.
258 boost::shared_ptr<gmx::Regex> regex_;
262 * Data structure for string keyword expression evaluation.
264 struct t_methoddata_kwstr
266 /** Matching type for the strings. */
267 gmx::SelectionStringMatchType matchType;
268 /** Array of values for the keyword. */
270 /** Array of strings/regular expressions to match against.*/
271 std::vector<StringKeywordMatchItem> matches;
276 /** Parameters for integer keyword evaluation. */
277 static gmx_ana_selparam_t smparams_keyword_int[] = {
278 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
279 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
282 /** Parameters for real keyword evaluation. */
283 static gmx_ana_selparam_t smparams_keyword_real[] = {
284 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL | SPAR_DYNAMIC},
285 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
288 /** Parameters for string keyword evaluation. */
289 static gmx_ana_selparam_t smparams_keyword_str[] = {
290 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
291 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_VARNUM},
294 /** Selection method data for integer keyword evaluation. */
295 gmx_ana_selmethod_t sm_keyword_int = {
296 "kw_int", GROUP_VALUE, SMETH_SINGLEVAL,
297 asize(smparams_keyword_int), smparams_keyword_int,
304 &evaluate_keyword_int,
309 /** Selection method data for real keyword evaluation. */
310 gmx_ana_selmethod_t sm_keyword_real = {
311 "kw_real", GROUP_VALUE, SMETH_SINGLEVAL,
312 asize(smparams_keyword_real), smparams_keyword_real,
319 &evaluate_keyword_real,
324 /** Selection method data for string keyword evaluation. */
325 gmx_ana_selmethod_t sm_keyword_str = {
326 "kw_str", GROUP_VALUE, SMETH_SINGLEVAL,
327 asize(smparams_keyword_str), smparams_keyword_str,
334 &evaluate_keyword_str,
340 * Initializes keyword evaluation for an arbitrary group.
342 * \param[in] top Not used.
343 * \param[in] npar Not used.
344 * \param[in] param Not used.
345 * \param[in] data Should point to \ref t_methoddata_kweval.
346 * \returns 0 on success, a non-zero error code on return.
348 * Calls the initialization method of the wrapped keyword.
351 init_kweval(t_topology *top, int npar, gmx_ana_selparam_t * param, void *data);
353 * Initializes output for keyword evaluation in an arbitrary group.
355 * \param[in] top Not used.
356 * \param[in,out] out Pointer to output data structure.
357 * \param[in,out] data Should point to \c t_methoddata_kweval.
358 * \returns 0 for success.
361 init_output_kweval(t_topology *top, gmx_ana_selvalue_t *out, void *data);
362 /** Frees the data allocated for keyword evaluation in an arbitrary group. */
364 free_data_kweval(void *data);
365 /** Initializes frame evaluation for keyword evaluation in an arbitrary group. */
367 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
369 * Evaluates keywords in an arbitrary group.
371 * See sel_updatefunc() for description of the parameters.
372 * \p data should point to a \c t_methoddata_kweval.
374 * Calls the evaluation function of the wrapped keyword with the given
375 * parameters, with the exception of using \c t_methoddata_kweval::g for the
379 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, gmx_ana_index_t *g,
380 gmx_ana_selvalue_t *out, void *data);
382 * Evaluates keywords in an arbitrary set of positions.
384 * See sel_updatefunc() for description of the parameters.
385 * \p data should point to a \c t_methoddata_kweval.
387 * Calls the evaluation function of the wrapped keyword with the given
388 * parameters, with the exception of using \c t_methoddata_kweval::p for the
389 * evaluation positions.
392 evaluate_kweval_pos(t_topology *top, t_trxframe *fr, t_pbc *pbc, gmx_ana_index_t *g,
393 gmx_ana_selvalue_t *out, void *data);
396 * Data structure for keyword evaluation in arbitrary groups.
398 struct t_methoddata_kweval
400 //! Initialize new keyword evaluator for the given keyword.
401 t_methoddata_kweval(gmx_ana_selmethod_t *method, void *data)
402 : kwmethod(method), kwmdata(data)
404 gmx_ana_index_clear(&g);
406 ~t_methoddata_kweval()
408 _gmx_selelem_free_method(kwmethod, kwmdata);
411 //! Wrapped keyword method for evaluating the values.
412 gmx_ana_selmethod_t *kwmethod;
413 //! Method data for \p kwmethod.
415 //! Group in which \p kwmethod should be evaluated.
417 //! Positions for which \p kwmethod should be evaluated.
421 /** Parameters for keyword evaluation in an arbitrary group. */
422 static gmx_ana_selparam_t smparams_kweval_group[] = {
423 {NULL, {GROUP_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
425 /** Parameters for keyword evaluation from positions. */
426 static gmx_ana_selparam_t smparams_kweval_pos[] = {
427 {NULL, {POS_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
431 /********************************************************************
432 * INTEGER KEYWORD EVALUATION
433 ********************************************************************/
436 init_data_kwint(int /* npar */, gmx_ana_selparam_t * /* param */)
438 t_methoddata_kwint *data;
445 init_kwint(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
447 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
449 d->v = param[0].val.u.i;
450 d->n = param[1].val.nr;
451 d->r = param[1].val.u.i;
455 * See sel_updatefunc() for description of the parameters.
456 * \p data should point to a \c t_methoddata_kwint.
458 * Does a binary search to find which atoms match the ranges in the
459 * \c t_methoddata_kwint structure for this selection.
460 * Matching atoms are stored in \p out->u.g.
463 evaluate_keyword_int(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
464 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
466 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
467 int n, i, j, jmin, jmax;
472 for (i = 0; i < g->isize; ++i)
475 if (d->r[0] > val || d->r[2*n-1] < val)
481 while (jmax - jmin > 1)
483 j = jmin + (jmax - jmin) / 2;
491 if (val <= d->r[2*j+1])
498 if (val <= d->r[2*jmin+1])
500 out->u.g->index[out->u.g->isize++] = g->index[i];
506 /********************************************************************
507 * REAL KEYWORD EVALUATION
508 ********************************************************************/
511 init_data_kwreal(int /* npar */, gmx_ana_selparam_t * /* param */)
513 t_methoddata_kwreal *data;
520 init_kwreal(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
522 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
524 d->v = param[0].val.u.r;
525 d->n = param[1].val.nr;
526 d->r = param[1].val.u.r;
530 * See sel_updatefunc() for description of the parameters.
531 * \p data should point to a \c t_methoddata_kwreal.
533 * Does a binary search to find which atoms match the ranges in the
534 * \c t_methoddata_kwreal structure for this selection.
535 * Matching atoms are stored in \p out->u.g.
538 evaluate_keyword_real(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
539 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
541 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
542 int n, i, j, jmin, jmax;
547 for (i = 0; i < g->isize; ++i)
550 if (d->r[0] > val || d->r[2*n-1] < val)
556 while (jmax - jmin > 1)
558 j = jmin + (jmax - jmin) / 2;
566 if (val <= d->r[2*j+1])
573 if (val <= d->r[2*jmin+1])
575 out->u.g->index[out->u.g->isize++] = g->index[i];
581 /********************************************************************
582 * STRING KEYWORD EVALUATION
583 ********************************************************************/
586 init_data_kwstr(int /* npar */, gmx_ana_selparam_t * /* param */)
588 t_methoddata_kwstr *data = new t_methoddata_kwstr();
589 data->matchType = gmx::eStringMatchType_Auto;
594 * \param[in,out] sel Selection element to initialize.
595 * \param[in] matchType Method to use to match string values.
597 * Sets the string matching method for string keyword matching.
600 _gmx_selelem_set_kwstr_match_type(const gmx::SelectionTreeElementPointer &sel,
601 gmx::SelectionStringMatchType matchType)
603 t_methoddata_kwstr *d = static_cast<t_methoddata_kwstr *>(sel->u.expr.mdata);
605 if (sel->type != SEL_EXPRESSION || !sel->u.expr.method
606 || sel->u.expr.method->name != sm_keyword_str.name)
610 d->matchType = matchType;
614 init_kwstr(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
616 t_methoddata_kwstr *d = static_cast<t_methoddata_kwstr *>(data);
618 d->v = param[0].val.u.s;
619 /* Return if this is not the first time */
620 if (!d->matches.empty())
624 int n = param[1].val.nr;
625 d->matches.reserve(n);
626 for (int i = 0; i < n; ++i)
628 const char *s = param[1].val.u.s[i];
629 d->matches.push_back(StringKeywordMatchItem(d->matchType, s));
634 * \param data Data to free (should point to a t_methoddata_kwstr).
637 free_data_kwstr(void *data)
639 t_methoddata_kwstr *d = static_cast<t_methoddata_kwstr *>(data);
644 * See sel_updatefunc() for description of the parameters.
645 * \p data should point to a \c t_methoddata_kwstr.
647 * Does a linear search to find which atoms match the strings in the
648 * \c t_methoddata_kwstr structure for this selection.
649 * Wildcards are allowed in the strings.
650 * Matching atoms are stored in \p out->u.g.
653 evaluate_keyword_str(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
654 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
656 t_methoddata_kwstr *d = static_cast<t_methoddata_kwstr *>(data);
659 for (int i = 0; i < g->isize; ++i)
661 for (size_t j = 0; j < d->matches.size(); ++j)
663 if (d->matches[j].match(d->matchType, d->v[i]))
665 out->u.g->index[out->u.g->isize++] = g->index[i];
673 /********************************************************************
674 * KEYWORD EVALUATION FOR ARBITRARY GROUPS
675 ********************************************************************/
678 init_kweval(t_topology *top, int /* npar */, gmx_ana_selparam_t * /* param */, void *data)
680 t_methoddata_kweval *d = static_cast<t_methoddata_kweval *>(data);
682 d->kwmethod->init(top, 0, NULL, d->kwmdata);
686 init_output_kweval(t_topology * /* top */, gmx_ana_selvalue_t *out, void *data)
688 t_methoddata_kweval *d = static_cast<t_methoddata_kweval *>(data);
690 out->nr = d->g.isize;
691 _gmx_selvalue_reserve(out, out->nr);
695 * \param data Data to free (should point to a \c t_methoddata_kweval).
697 * Frees the memory allocated for all the members of \c t_methoddata_kweval.
700 free_data_kweval(void *data)
702 t_methoddata_kweval *d = static_cast<t_methoddata_kweval *>(data);
708 * \param[in] top Topology.
709 * \param[in] fr Current frame.
710 * \param[in] pbc PBC structure.
711 * \param data Should point to a \ref t_methoddata_kweval.
712 * \returns 0 on success, a non-zero error code on error.
714 * Creates a lookup structure that enables fast queries of whether a point
715 * is within the solid angle or not.
718 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
720 t_methoddata_kweval *d = static_cast<t_methoddata_kweval *>(data);
722 d->kwmethod->init_frame(top, fr, pbc, d->kwmdata);
726 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc,
727 gmx_ana_index_t * /* g */, gmx_ana_selvalue_t *out, void *data)
729 t_methoddata_kweval *d = static_cast<t_methoddata_kweval *>(data);
731 d->kwmethod->update(top, fr, pbc, &d->g, out, d->kwmdata);
735 evaluate_kweval_pos(t_topology *top, t_trxframe *fr, t_pbc *pbc,
736 gmx_ana_index_t * /* g */, gmx_ana_selvalue_t *out, void *data)
738 t_methoddata_kweval *d = static_cast<t_methoddata_kweval *>(data);
740 d->kwmethod->pupdate(top, fr, pbc, &d->p, out, d->kwmdata);
744 * Initializes a selection element for evaluating a keyword in a given group.
746 * \param[in] method Keyword selection method to evaluate.
747 * \param[in] params Parameters to pass to initialization (the child group).
748 * \param[in] scanner Scanner data structure.
749 * \returns Pointer to the created selection element.
751 * Implements _gmx_sel_init_keyword_evaluator() for \ref GROUP_VALUE input
754 static gmx::SelectionTreeElementPointer
755 init_evaluator_group(gmx_ana_selmethod_t *method,
756 const gmx::SelectionParserParameterList ¶ms,
759 if ((method->flags & (SMETH_SINGLEVAL | SMETH_VARNUMVAL))
760 || method->outinit || method->pupdate)
763 = gmx::formatString("Keyword '%s' cannot be evaluated in this context",
765 GMX_THROW(gmx::InvalidInputError(message));
768 // TODO: For same ... as ..., some other location could be more intuitive.
769 gmx::SelectionTreeElementPointer sel(
770 new gmx::SelectionTreeElement(
771 SEL_EXPRESSION, _gmx_sel_lexer_get_current_location(scanner)));
772 _gmx_selelem_set_method(sel, method, scanner);
774 t_methoddata_kweval *data
775 = new t_methoddata_kweval(sel->u.expr.method, sel->u.expr.mdata);
777 snew(sel->u.expr.method, 1);
778 sel->u.expr.method->name = data->kwmethod->name;
779 sel->u.expr.method->type = data->kwmethod->type;
780 sel->u.expr.method->flags = data->kwmethod->flags | SMETH_VARNUMVAL;
781 sel->u.expr.method->init_data = NULL;
782 sel->u.expr.method->set_poscoll = NULL;
783 sel->u.expr.method->init = method->init ? &init_kweval : NULL;
784 sel->u.expr.method->outinit = &init_output_kweval;
785 sel->u.expr.method->free = &free_data_kweval;
786 sel->u.expr.method->init_frame = method->init_frame ? &init_frame_kweval : NULL;
787 sel->u.expr.method->update = &evaluate_kweval;
788 sel->u.expr.method->pupdate = NULL;
789 sel->u.expr.method->nparams = asize(smparams_kweval_group);
790 sel->u.expr.method->param = smparams_kweval_group;
791 _gmx_selelem_init_method_params(sel, scanner);
792 sel->u.expr.mdata = data;
794 sel->u.expr.method->param[0].val.u.g = &data->g;
796 _gmx_sel_parse_params(params, sel->u.expr.method->nparams,
797 sel->u.expr.method->param, sel, scanner);
802 * Initializes a selection element for evaluating a keyword in given positions.
804 * \param[in] method Keyword selection method to evaluate.
805 * \param[in] params Parameters to pass to initialization (the child positions).
806 * \param[in] scanner Scanner data structure.
807 * \returns Pointer to the created selection element.
809 * Implements _gmx_sel_init_keyword_evaluator() for \ref POS_VALUE input
812 static gmx::SelectionTreeElementPointer
813 init_evaluator_pos(gmx_ana_selmethod_t *method,
814 const gmx::SelectionParserParameterList ¶ms,
817 if ((method->flags & (SMETH_SINGLEVAL | SMETH_VARNUMVAL))
818 || method->outinit || method->pupdate == NULL)
821 = gmx::formatString("Keyword '%s' cannot be evaluated in this context",
823 GMX_THROW(gmx::InvalidInputError(message));
826 gmx::SelectionTreeElementPointer sel(
827 new gmx::SelectionTreeElement(
828 SEL_EXPRESSION, _gmx_sel_lexer_get_current_location(scanner)));
829 _gmx_selelem_set_method(sel, method, scanner);
831 t_methoddata_kweval *data
832 = new t_methoddata_kweval(sel->u.expr.method, sel->u.expr.mdata);
834 snew(sel->u.expr.method, 1);
835 sel->u.expr.method->name = data->kwmethod->name;
836 sel->u.expr.method->type = data->kwmethod->type;
837 sel->u.expr.method->flags = data->kwmethod->flags | SMETH_SINGLEVAL;
838 sel->u.expr.method->init_data = NULL;
839 sel->u.expr.method->set_poscoll = NULL;
840 sel->u.expr.method->init = method->init ? &init_kweval : NULL;
841 sel->u.expr.method->outinit = NULL;
842 sel->u.expr.method->free = &free_data_kweval;
843 sel->u.expr.method->init_frame = method->init_frame ? &init_frame_kweval : NULL;
844 sel->u.expr.method->update = &evaluate_kweval_pos;
845 sel->u.expr.method->pupdate = NULL;
846 sel->u.expr.method->nparams = asize(smparams_kweval_pos);
847 sel->u.expr.method->param = smparams_kweval_pos;
848 _gmx_selelem_init_method_params(sel, scanner);
849 sel->u.expr.mdata = data;
851 sel->u.expr.method->param[0].val.u.p = &data->p;
853 _gmx_sel_parse_params(params, sel->u.expr.method->nparams,
854 sel->u.expr.method->param, sel, scanner);
858 gmx::SelectionTreeElementPointer
859 _gmx_sel_init_keyword_evaluator(gmx_ana_selmethod_t *method,
860 const gmx::SelectionTreeElementPointer &child,
863 gmx::SelectionParserParameterList params;
865 gmx::SelectionParserParameter::createFromExpression(NULL, child));
866 if (child->v.type == GROUP_VALUE)
868 return init_evaluator_group(method, params, scanner);
870 else if (child->v.type == POS_VALUE)
872 return init_evaluator_pos(method, params, scanner);
876 std::string text(_gmx_sel_lexer_get_text(scanner, child->location()));
878 = gmx::formatString("Expression '%s' cannot be used to evaluate keywords",
880 GMX_THROW(gmx::InvalidInputError(message));