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
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/messagestringcollector.h"
59 #include "gromacs/utility/smalloc.h"
60 #include "gromacs/utility/stringutil.h"
63 #include "parsetree.h"
66 #include "selmethod.h"
69 * Allocates data for integer keyword evaluation.
71 * \param[in] npar Not used.
72 * \param param Not used.
73 * \returns Pointer to the allocated data (\ref t_methoddata_kwint).
75 * Allocates memory for a \ref t_methoddata_kwint structure.
78 init_data_kwint(int npar, gmx_ana_selparam_t * param);
80 * Allocates data for real keyword evaluation.
82 * \param[in] npar Not used.
83 * \param param Not used.
84 * \returns Pointer to the allocated data (\ref t_methoddata_kwreal).
86 * Allocates memory for a \ref t_methoddata_kwreal structure.
89 init_data_kwreal(int npar, gmx_ana_selparam_t * param);
91 * Allocates data for string keyword evaluation.
93 * \param[in] npar Not used.
94 * \param param Not used.
95 * \returns Pointer to the allocated data (t_methoddata_kwstr).
97 * Allocates memory for a t_methoddata_kwstr structure.
100 init_data_kwstr(int npar, gmx_ana_selparam_t * param);
101 /** /brief Initializes data for integer keyword evaluation.
103 * \param[in] top Not used.
104 * \param[in] npar Not used (should be 2).
105 * \param[in] param Method parameters (should point to \ref smparams_keyword_int).
106 * \param[in] data Should point to \ref t_methoddata_kwint.
109 init_kwint(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
111 * Initializes data for real keyword evaluation.
113 * \param[in] top Not used.
114 * \param[in] npar Not used (should be 2).
115 * \param[in] param Method parameters (should point to \ref smparams_keyword_real).
116 * \param[in] data Should point to \ref t_methoddata_kwreal.
117 * \returns 0 (the initialization always succeeds).
120 init_kwreal(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
122 * Initializes data for string keyword evaluation.
124 * \param[in] top Not used.
125 * \param[in] npar Not used (should be 2).
126 * \param[in] param Method parameters (should point to \ref smparams_keyword_str).
127 * \param[in] data Should point to t_methoddata_kwstr.
130 init_kwstr(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
131 /** Frees the memory allocated for string keyword evaluation. */
133 free_data_kwstr(void *data);
134 /** Evaluates integer selection keywords. */
136 evaluate_keyword_int(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
137 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
138 /** Evaluates real selection keywords. */
140 evaluate_keyword_real(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
141 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
142 /** Evaluates string selection keywords. */
144 evaluate_keyword_str(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
145 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
148 * Data structure for integer keyword expression evaluation.
150 typedef struct t_methoddata_kwint
152 /** Array of values for the keyword. */
154 /** Number of ranges in the \p r array. */
157 * Array of sorted integer ranges to match against.
159 * Each range is made of two integers, giving the endpoints (inclusive).
160 * This field stores the pointer to the ranges allocated by the
161 * parameter parser; see \ref SPAR_RANGES for more information.
164 } t_methoddata_kwint;
167 * Data structure for real keyword expression evaluation.
169 typedef struct t_methoddata_kwreal
171 /** Array of values for the keyword. */
173 /** Number of ranges in the \p r array. */
176 * Array of sorted ranges to match against.
178 * Each range is made of two values, giving the endpoints (inclusive).
179 * This field stores the pointer to the ranges allocated by the
180 * parameter parser; see \ref SPAR_RANGES for more information.
183 } t_methoddata_kwreal;
189 * Single item in the list of strings/regular expressions to match.
191 * \ingroup module_selection
193 class StringKeywordMatchItem
197 * Constructs a matcher from a string.
199 * \param[in] matchType String matching type.
200 * \param[in] str String to use for matching.
202 StringKeywordMatchItem(gmx::SelectionStringMatchType matchType,
206 bool bRegExp = (matchType == gmx::eStringMatchType_RegularExpression);
207 if (matchType == gmx::eStringMatchType_Auto)
209 for (size_t j = 0; j < std::strlen(str); ++j)
211 if (std::ispunct(str[j]) && str[j] != '?' && str[j] != '*')
220 if (!gmx::Regex::isSupported())
223 = gmx::formatString("No regular expression support, "
224 "cannot match \"%s\"", str);
225 GMX_THROW(gmx::InvalidInputError(message));
227 regex_.reset(new gmx::Regex(str));
232 * Checks whether this item matches a string.
234 * \param[in] matchType String matching type.
235 * \param[in] value String to match.
236 * \returns true if this item matches \p value.
238 bool match(gmx::SelectionStringMatchType matchType,
239 const char *value) const
241 if (matchType == gmx::eStringMatchType_Exact)
243 return str_ == value;
247 return gmx::regexMatch(value, *regex_);
251 return gmx_wcmatch(str_.c_str(), value) == 0;
256 //! The raw string passed for the matcher.
258 //! Regular expression compiled from \p str_, if applicable.
259 boost::shared_ptr<gmx::Regex> regex_;
263 * Data structure for string keyword expression evaluation.
265 struct t_methoddata_kwstr
267 /** Matching type for the strings. */
268 gmx::SelectionStringMatchType matchType;
269 /** Array of values for the keyword. */
271 /** Array of strings/regular expressions to match against.*/
272 std::vector<StringKeywordMatchItem> matches;
277 /** Parameters for integer keyword evaluation. */
278 static gmx_ana_selparam_t smparams_keyword_int[] = {
279 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
280 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
283 /** Parameters for real keyword evaluation. */
284 static gmx_ana_selparam_t smparams_keyword_real[] = {
285 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL | SPAR_DYNAMIC},
286 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
289 /** Parameters for string keyword evaluation. */
290 static gmx_ana_selparam_t smparams_keyword_str[] = {
291 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
292 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_VARNUM},
295 /** Selection method data for integer keyword evaluation. */
296 gmx_ana_selmethod_t sm_keyword_int = {
297 "kw_int", GROUP_VALUE, SMETH_SINGLEVAL,
298 asize(smparams_keyword_int), smparams_keyword_int,
305 &evaluate_keyword_int,
310 /** Selection method data for real keyword evaluation. */
311 gmx_ana_selmethod_t sm_keyword_real = {
312 "kw_real", GROUP_VALUE, SMETH_SINGLEVAL,
313 asize(smparams_keyword_real), smparams_keyword_real,
320 &evaluate_keyword_real,
325 /** Selection method data for string keyword evaluation. */
326 gmx_ana_selmethod_t sm_keyword_str = {
327 "kw_str", GROUP_VALUE, SMETH_SINGLEVAL,
328 asize(smparams_keyword_str), smparams_keyword_str,
335 &evaluate_keyword_str,
341 * Initializes keyword evaluation for an arbitrary group.
343 * \param[in] top Not used.
344 * \param[in] npar Not used.
345 * \param[in] param Not used.
346 * \param[in] data Should point to \ref t_methoddata_kweval.
347 * \returns 0 on success, a non-zero error code on return.
349 * Calls the initialization method of the wrapped keyword.
352 init_kweval(t_topology *top, int npar, gmx_ana_selparam_t * param, void *data);
354 * Initializes output for keyword evaluation in an arbitrary group.
356 * \param[in] top Not used.
357 * \param[in,out] out Pointer to output data structure.
358 * \param[in,out] data Should point to \c t_methoddata_kweval.
359 * \returns 0 for success.
362 init_output_kweval(t_topology *top, gmx_ana_selvalue_t *out, void *data);
363 /** Frees the data allocated for keyword evaluation in an arbitrary group. */
365 free_data_kweval(void *data);
366 /** Initializes frame evaluation for keyword evaluation in an arbitrary group. */
368 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
370 * Evaluates keywords in an arbitrary group.
372 * See sel_updatefunc() for description of the parameters.
373 * \p data should point to a \c t_methoddata_kweval.
375 * Calls the evaluation function of the wrapped keyword with the given
376 * parameters, with the exception of using \c t_methoddata_kweval::g for the
380 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, gmx_ana_index_t *g,
381 gmx_ana_selvalue_t *out, void *data);
383 * Evaluates keywords in an arbitrary set of positions.
385 * See sel_updatefunc() for description of the parameters.
386 * \p data should point to a \c t_methoddata_kweval.
388 * Calls the evaluation function of the wrapped keyword with the given
389 * parameters, with the exception of using \c t_methoddata_kweval::p for the
390 * evaluation positions.
393 evaluate_kweval_pos(t_topology *top, t_trxframe *fr, t_pbc *pbc, gmx_ana_index_t *g,
394 gmx_ana_selvalue_t *out, void *data);
397 * Data structure for keyword evaluation in arbitrary groups.
399 struct t_methoddata_kweval
401 //! Initialize new keyword evaluator for the given keyword.
402 t_methoddata_kweval(gmx_ana_selmethod_t *method, void *data)
403 : kwmethod(method), kwmdata(data)
405 gmx_ana_index_clear(&g);
407 ~t_methoddata_kweval()
409 _gmx_selelem_free_method(kwmethod, kwmdata);
412 //! Wrapped keyword method for evaluating the values.
413 gmx_ana_selmethod_t *kwmethod;
414 //! Method data for \p kwmethod.
416 //! Group in which \p kwmethod should be evaluated.
418 //! Positions for which \p kwmethod should be evaluated.
422 /** Parameters for keyword evaluation in an arbitrary group. */
423 static gmx_ana_selparam_t smparams_kweval_group[] = {
424 {NULL, {GROUP_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
426 /** Parameters for keyword evaluation from positions. */
427 static gmx_ana_selparam_t smparams_kweval_pos[] = {
428 {NULL, {POS_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
432 /********************************************************************
433 * INTEGER KEYWORD EVALUATION
434 ********************************************************************/
437 init_data_kwint(int /* npar */, gmx_ana_selparam_t * /* param */)
439 t_methoddata_kwint *data;
446 init_kwint(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
448 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
450 d->v = param[0].val.u.i;
451 d->n = param[1].val.nr;
452 d->r = param[1].val.u.i;
456 * See sel_updatefunc() for description of the parameters.
457 * \p data should point to a \c t_methoddata_kwint.
459 * Does a binary search to find which atoms match the ranges in the
460 * \c t_methoddata_kwint structure for this selection.
461 * Matching atoms are stored in \p out->u.g.
464 evaluate_keyword_int(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
465 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
467 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
468 int n, i, j, jmin, jmax;
473 for (i = 0; i < g->isize; ++i)
476 if (d->r[0] > val || d->r[2*n-1] < val)
482 while (jmax - jmin > 1)
484 j = jmin + (jmax - jmin) / 2;
492 if (val <= d->r[2*j+1])
499 if (val <= d->r[2*jmin+1])
501 out->u.g->index[out->u.g->isize++] = g->index[i];
507 /********************************************************************
508 * REAL KEYWORD EVALUATION
509 ********************************************************************/
512 init_data_kwreal(int /* npar */, gmx_ana_selparam_t * /* param */)
514 t_methoddata_kwreal *data;
521 init_kwreal(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
523 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
525 d->v = param[0].val.u.r;
526 d->n = param[1].val.nr;
527 d->r = param[1].val.u.r;
531 * See sel_updatefunc() for description of the parameters.
532 * \p data should point to a \c t_methoddata_kwreal.
534 * Does a binary search to find which atoms match the ranges in the
535 * \c t_methoddata_kwreal structure for this selection.
536 * Matching atoms are stored in \p out->u.g.
539 evaluate_keyword_real(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
540 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
542 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
543 int n, i, j, jmin, jmax;
548 for (i = 0; i < g->isize; ++i)
551 if (d->r[0] > val || d->r[2*n-1] < val)
557 while (jmax - jmin > 1)
559 j = jmin + (jmax - jmin) / 2;
567 if (val <= d->r[2*j+1])
574 if (val <= d->r[2*jmin+1])
576 out->u.g->index[out->u.g->isize++] = g->index[i];
582 /********************************************************************
583 * STRING KEYWORD EVALUATION
584 ********************************************************************/
587 init_data_kwstr(int /* npar */, gmx_ana_selparam_t * /* param */)
589 t_methoddata_kwstr *data = new t_methoddata_kwstr();
590 data->matchType = gmx::eStringMatchType_Auto;
595 * \param[in,out] sel Selection element to initialize.
596 * \param[in] matchType Method to use to match string values.
598 * Sets the string matching method for string keyword matching.
601 _gmx_selelem_set_kwstr_match_type(const gmx::SelectionTreeElementPointer &sel,
602 gmx::SelectionStringMatchType matchType)
604 t_methoddata_kwstr *d = static_cast<t_methoddata_kwstr *>(sel->u.expr.mdata);
606 if (sel->type != SEL_EXPRESSION || !sel->u.expr.method
607 || sel->u.expr.method->name != sm_keyword_str.name)
611 d->matchType = matchType;
615 init_kwstr(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
617 t_methoddata_kwstr *d = static_cast<t_methoddata_kwstr *>(data);
619 d->v = param[0].val.u.s;
620 /* Return if this is not the first time */
621 if (!d->matches.empty())
625 int n = param[1].val.nr;
626 d->matches.reserve(n);
627 for (int i = 0; i < n; ++i)
629 const char *s = param[1].val.u.s[i];
630 d->matches.push_back(StringKeywordMatchItem(d->matchType, s));
635 * \param data Data to free (should point to a t_methoddata_kwstr).
638 free_data_kwstr(void *data)
640 t_methoddata_kwstr *d = static_cast<t_methoddata_kwstr *>(data);
645 * See sel_updatefunc() for description of the parameters.
646 * \p data should point to a \c t_methoddata_kwstr.
648 * Does a linear search to find which atoms match the strings in the
649 * \c t_methoddata_kwstr structure for this selection.
650 * Wildcards are allowed in the strings.
651 * Matching atoms are stored in \p out->u.g.
654 evaluate_keyword_str(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
655 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
657 t_methoddata_kwstr *d = static_cast<t_methoddata_kwstr *>(data);
660 for (int i = 0; i < g->isize; ++i)
662 for (size_t j = 0; j < d->matches.size(); ++j)
664 if (d->matches[j].match(d->matchType, d->v[i]))
666 out->u.g->index[out->u.g->isize++] = g->index[i];
674 /********************************************************************
675 * KEYWORD EVALUATION FOR ARBITRARY GROUPS
676 ********************************************************************/
679 init_kweval(t_topology *top, int /* npar */, gmx_ana_selparam_t * /* param */, void *data)
681 t_methoddata_kweval *d = static_cast<t_methoddata_kweval *>(data);
683 d->kwmethod->init(top, 0, NULL, d->kwmdata);
687 init_output_kweval(t_topology * /* top */, gmx_ana_selvalue_t *out, void *data)
689 t_methoddata_kweval *d = static_cast<t_methoddata_kweval *>(data);
691 out->nr = d->g.isize;
692 _gmx_selvalue_reserve(out, out->nr);
696 * \param data Data to free (should point to a \c t_methoddata_kweval).
698 * Frees the memory allocated for all the members of \c t_methoddata_kweval.
701 free_data_kweval(void *data)
703 t_methoddata_kweval *d = static_cast<t_methoddata_kweval *>(data);
709 * \param[in] top Topology.
710 * \param[in] fr Current frame.
711 * \param[in] pbc PBC structure.
712 * \param data Should point to a \ref t_methoddata_kweval.
713 * \returns 0 on success, a non-zero error code on error.
715 * Creates a lookup structure that enables fast queries of whether a point
716 * is within the solid angle or not.
719 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
721 t_methoddata_kweval *d = static_cast<t_methoddata_kweval *>(data);
723 d->kwmethod->init_frame(top, fr, pbc, d->kwmdata);
727 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc,
728 gmx_ana_index_t * /* g */, gmx_ana_selvalue_t *out, void *data)
730 t_methoddata_kweval *d = static_cast<t_methoddata_kweval *>(data);
732 d->kwmethod->update(top, fr, pbc, &d->g, out, d->kwmdata);
736 evaluate_kweval_pos(t_topology *top, t_trxframe *fr, t_pbc *pbc,
737 gmx_ana_index_t * /* g */, gmx_ana_selvalue_t *out, void *data)
739 t_methoddata_kweval *d = static_cast<t_methoddata_kweval *>(data);
741 d->kwmethod->pupdate(top, fr, pbc, &d->p, out, d->kwmdata);
745 * Initializes a selection element for evaluating a keyword in a given group.
747 * \param[in] method Keyword selection method to evaluate.
748 * \param[in] params Parameters to pass to initialization (the child group).
749 * \param[in] scanner Scanner data structure.
750 * \returns Pointer to the created selection element.
752 * Implements _gmx_sel_init_keyword_evaluator() for \ref GROUP_VALUE input
755 static gmx::SelectionTreeElementPointer
756 init_evaluator_group(gmx_ana_selmethod_t *method,
757 const gmx::SelectionParserParameterList ¶ms,
760 if ((method->flags & (SMETH_SINGLEVAL | SMETH_VARNUMVAL))
761 || method->outinit || method->pupdate)
764 = gmx::formatString("Keyword '%s' cannot be evaluated in this context",
766 GMX_THROW(gmx::InvalidInputError(message));
769 // TODO: For same ... as ..., some other location could be more intuitive.
770 gmx::SelectionTreeElementPointer sel(
771 new gmx::SelectionTreeElement(
772 SEL_EXPRESSION, _gmx_sel_lexer_get_current_location(scanner)));
773 _gmx_selelem_set_method(sel, method, scanner);
775 t_methoddata_kweval *data
776 = new t_methoddata_kweval(sel->u.expr.method, sel->u.expr.mdata);
778 snew(sel->u.expr.method, 1);
779 sel->u.expr.method->name = data->kwmethod->name;
780 sel->u.expr.method->type = data->kwmethod->type;
781 sel->u.expr.method->flags = data->kwmethod->flags | SMETH_VARNUMVAL;
782 sel->u.expr.method->init_data = NULL;
783 sel->u.expr.method->set_poscoll = NULL;
784 sel->u.expr.method->init = method->init ? &init_kweval : NULL;
785 sel->u.expr.method->outinit = &init_output_kweval;
786 sel->u.expr.method->free = &free_data_kweval;
787 sel->u.expr.method->init_frame = method->init_frame ? &init_frame_kweval : NULL;
788 sel->u.expr.method->update = &evaluate_kweval;
789 sel->u.expr.method->pupdate = NULL;
790 sel->u.expr.method->nparams = asize(smparams_kweval_group);
791 sel->u.expr.method->param = smparams_kweval_group;
792 _gmx_selelem_init_method_params(sel, scanner);
793 sel->u.expr.mdata = data;
795 sel->u.expr.method->param[0].val.u.g = &data->g;
797 _gmx_sel_parse_params(params, sel->u.expr.method->nparams,
798 sel->u.expr.method->param, sel, scanner);
803 * Initializes a selection element for evaluating a keyword in given positions.
805 * \param[in] method Keyword selection method to evaluate.
806 * \param[in] params Parameters to pass to initialization (the child positions).
807 * \param[in] scanner Scanner data structure.
808 * \returns Pointer to the created selection element.
810 * Implements _gmx_sel_init_keyword_evaluator() for \ref POS_VALUE input
813 static gmx::SelectionTreeElementPointer
814 init_evaluator_pos(gmx_ana_selmethod_t *method,
815 const gmx::SelectionParserParameterList ¶ms,
818 if ((method->flags & (SMETH_SINGLEVAL | SMETH_VARNUMVAL))
819 || method->outinit || method->pupdate == NULL)
822 = gmx::formatString("Keyword '%s' cannot be evaluated in this context",
824 GMX_THROW(gmx::InvalidInputError(message));
827 gmx::SelectionTreeElementPointer sel(
828 new gmx::SelectionTreeElement(
829 SEL_EXPRESSION, _gmx_sel_lexer_get_current_location(scanner)));
830 _gmx_selelem_set_method(sel, method, scanner);
832 t_methoddata_kweval *data
833 = new t_methoddata_kweval(sel->u.expr.method, sel->u.expr.mdata);
835 snew(sel->u.expr.method, 1);
836 sel->u.expr.method->name = data->kwmethod->name;
837 sel->u.expr.method->type = data->kwmethod->type;
838 sel->u.expr.method->flags = data->kwmethod->flags | SMETH_SINGLEVAL;
839 sel->u.expr.method->init_data = NULL;
840 sel->u.expr.method->set_poscoll = NULL;
841 sel->u.expr.method->init = method->init ? &init_kweval : NULL;
842 sel->u.expr.method->outinit = NULL;
843 sel->u.expr.method->free = &free_data_kweval;
844 sel->u.expr.method->init_frame = method->init_frame ? &init_frame_kweval : NULL;
845 sel->u.expr.method->update = &evaluate_kweval_pos;
846 sel->u.expr.method->pupdate = NULL;
847 sel->u.expr.method->nparams = asize(smparams_kweval_pos);
848 sel->u.expr.method->param = smparams_kweval_pos;
849 _gmx_selelem_init_method_params(sel, scanner);
850 sel->u.expr.mdata = data;
852 sel->u.expr.method->param[0].val.u.p = &data->p;
854 _gmx_sel_parse_params(params, sel->u.expr.method->nparams,
855 sel->u.expr.method->param, sel, scanner);
859 gmx::SelectionTreeElementPointer
860 _gmx_sel_init_keyword_evaluator(gmx_ana_selmethod_t *method,
861 const gmx::SelectionTreeElementPointer &child,
864 gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
866 sprintf(buf, "In evaluation of '%s'", method->name);
867 gmx::MessageStringContext context(errors, buf);
869 gmx::SelectionParserParameterList params;
871 gmx::SelectionParserParameter::createFromExpression(NULL, child));
872 if (child->v.type == GROUP_VALUE)
874 return init_evaluator_group(method, params, scanner);
876 else if (child->v.type == POS_VALUE)
878 return init_evaluator_pos(method, params, scanner);
882 std::string text(_gmx_sel_lexer_get_text(scanner, child->location()));
884 = gmx::formatString("Expression '%s' cannot be used to evaluate keywords",
886 GMX_THROW(gmx::InvalidInputError(message));