2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012, by the GROMACS development team, led by
5 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
6 * others, as listed in the AUTHORS file in the top-level source
7 * 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"
65 /** Allocates data for integer keyword evaluation. */
67 init_data_kwint(int npar, gmx_ana_selparam_t *param);
68 /** Allocates data for real keyword evaluation. */
70 init_data_kwreal(int npar, gmx_ana_selparam_t *param);
71 /** Allocates data for string keyword evaluation. */
73 init_data_kwstr(int npar, gmx_ana_selparam_t *param);
74 /** Initializes data for integer keyword evaluation. */
76 init_kwint(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
77 /** Initializes data for real keyword evaluation. */
79 init_kwreal(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
80 /** Initializes data for string keyword evaluation. */
82 init_kwstr(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
83 /** Frees the memory allocated for string keyword evaluation. */
85 free_data_kwstr(void *data);
86 /** Evaluates integer selection keywords. */
88 evaluate_keyword_int(t_topology *top, t_trxframe *fr, t_pbc *pbc,
89 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
90 /** Evaluates real selection keywords. */
92 evaluate_keyword_real(t_topology *top, t_trxframe *fr, t_pbc *pbc,
93 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
94 /** Evaluates string selection keywords. */
96 evaluate_keyword_str(t_topology *top, t_trxframe *fr, t_pbc *pbc,
97 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
100 * Data structure for integer keyword expression evaluation.
102 typedef struct t_methoddata_kwint
104 /** Array of values for the keyword. */
106 /** Number of ranges in the \p r array. */
109 * Array of sorted integer ranges to match against.
111 * Each range is made of two integers, giving the endpoints (inclusive).
112 * This field stores the pointer to the ranges allocated by the
113 * parameter parser; see \ref SPAR_RANGES for more information.
116 } t_methoddata_kwint;
119 * Data structure for real keyword expression evaluation.
121 typedef struct t_methoddata_kwreal
123 /** Array of values for the keyword. */
125 /** Number of ranges in the \p r array. */
128 * Array of sorted ranges to match against.
130 * Each range is made of two values, giving the endpoints (inclusive).
131 * This field stores the pointer to the ranges allocated by the
132 * parameter parser; see \ref SPAR_RANGES for more information.
135 } t_methoddata_kwreal;
141 * Single item in the list of strings/regular expressions to match.
143 * \ingroup module_selection
145 class StringKeywordMatchItem
149 * Constructs a matcher from a string.
151 * \param[in] matchType String matching type.
152 * \param[in] str String to use for matching.
154 StringKeywordMatchItem(gmx::SelectionStringMatchType matchType,
158 bool bRegExp = (matchType == gmx::eStringMatchType_RegularExpression);
159 if (matchType == gmx::eStringMatchType_Auto)
161 for (size_t j = 0; j < std::strlen(str); ++j)
163 if (std::ispunct(str[j]) && str[j] != '?' && str[j] != '*')
172 if (!gmx::Regex::isSupported())
174 GMX_THROW(gmx::InvalidInputError(gmx::formatString(
175 "No regular expression support, "
176 "cannot match \"%s\"", str)));
178 regex_.reset(new gmx::Regex(str));
183 * Checks whether this item matches a string.
185 * \param[in] matchType String matching type.
186 * \param[in] value String to match.
187 * \returns true if this item matches \p value.
189 bool match(gmx::SelectionStringMatchType matchType,
190 const char *value) const
192 if (matchType == gmx::eStringMatchType_Exact)
194 return str_ == value;
198 return gmx::regexMatch(value, *regex_);
202 return gmx_wcmatch(str_.c_str(), value) == 0;
207 //! The raw string passed for the matcher.
209 //! Regular expression compiled from \p str_, if applicable.
210 boost::shared_ptr<gmx::Regex> regex_;
214 * Data structure for string keyword expression evaluation.
216 struct t_methoddata_kwstr
218 /** Matching type for the strings. */
219 gmx::SelectionStringMatchType matchType;
220 /** Array of values for the keyword. */
222 /** Array of strings/regular expressions to match against.*/
223 std::vector<StringKeywordMatchItem> matches;
228 /** Parameters for integer keyword evaluation. */
229 static gmx_ana_selparam_t smparams_keyword_int[] = {
230 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
231 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
234 /** Parameters for real keyword evaluation. */
235 static gmx_ana_selparam_t smparams_keyword_real[] = {
236 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL | SPAR_DYNAMIC},
237 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
240 /** Parameters for string keyword evaluation. */
241 static gmx_ana_selparam_t smparams_keyword_str[] = {
242 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
243 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_VARNUM},
246 /** \internal Selection method data for integer keyword evaluation. */
247 gmx_ana_selmethod_t sm_keyword_int = {
248 "kw_int", GROUP_VALUE, SMETH_SINGLEVAL,
249 asize(smparams_keyword_int), smparams_keyword_int,
256 &evaluate_keyword_int,
261 /** \internal Selection method data for real keyword evaluation. */
262 gmx_ana_selmethod_t sm_keyword_real = {
263 "kw_real", GROUP_VALUE, SMETH_SINGLEVAL,
264 asize(smparams_keyword_real), smparams_keyword_real,
271 &evaluate_keyword_real,
276 /** \internal Selection method data for string keyword evaluation. */
277 gmx_ana_selmethod_t sm_keyword_str = {
278 "kw_str", GROUP_VALUE, SMETH_SINGLEVAL,
279 asize(smparams_keyword_str), smparams_keyword_str,
286 &evaluate_keyword_str,
291 /** Initializes keyword evaluation for an arbitrary group. */
293 init_kweval(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
294 /** Initializes output for keyword evaluation in an arbitrary group. */
296 init_output_kweval(t_topology *top, gmx_ana_selvalue_t *out, void *data);
297 /** Frees the data allocated for keyword evaluation in an arbitrary group. */
299 free_data_kweval(void *data);
300 /** Initializes frame evaluation for keyword evaluation in an arbitrary group. */
302 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
303 /** Evaluates keywords in an arbitrary group. */
305 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc,
306 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
309 * Data structure for keyword evaluation in arbitrary groups.
313 /** Wrapped keyword method for evaluating the values. */
314 gmx_ana_selmethod_t *kwmethod;
315 /** Method data for \p kwmethod. */
317 /** Group in which \p kwmethod should be evaluated. */
319 } t_methoddata_kweval;
321 /** Parameters for keyword evaluation in an arbitrary group. */
322 static gmx_ana_selparam_t smparams_kweval[] = {
323 {NULL, {GROUP_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
327 /********************************************************************
328 * INTEGER KEYWORD EVALUATION
329 ********************************************************************/
332 * \param[in] npar Not used.
333 * \param param Not used.
334 * \returns Pointer to the allocated data (\ref t_methoddata_kwint).
336 * Allocates memory for a \ref t_methoddata_kwint structure.
339 init_data_kwint(int npar, gmx_ana_selparam_t *param)
341 t_methoddata_kwint *data;
348 * \param[in] top Not used.
349 * \param[in] npar Not used (should be 2).
350 * \param[in] param Method parameters (should point to \ref smparams_keyword_int).
351 * \param[in] data Should point to \ref t_methoddata_kwint.
354 init_kwint(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
356 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
358 d->v = param[0].val.u.i;
359 d->n = param[1].val.nr;
360 d->r = param[1].val.u.i;
364 * See sel_updatefunc() for description of the parameters.
365 * \p data should point to a \c t_methoddata_kwint.
367 * Does a binary search to find which atoms match the ranges in the
368 * \c t_methoddata_kwint structure for this selection.
369 * Matching atoms are stored in \p out->u.g.
372 evaluate_keyword_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_kwint *d = (t_methoddata_kwint *)data;
376 int n, i, j, jmin, jmax;
381 for (i = 0; i < g->isize; ++i)
384 if (d->r[0] > val || d->r[2*n-1] < val)
390 while (jmax - jmin > 1)
392 j = jmin + (jmax - jmin) / 2;
400 if (val <= d->r[2*j+1])
407 if (val <= d->r[2*jmin+1])
409 out->u.g->index[out->u.g->isize++] = g->index[i];
415 /********************************************************************
416 * REAL KEYWORD EVALUATION
417 ********************************************************************/
420 * \param[in] npar Not used.
421 * \param param Not used.
422 * \returns Pointer to the allocated data (\ref t_methoddata_kwreal).
424 * Allocates memory for a \ref t_methoddata_kwreal structure.
427 init_data_kwreal(int npar, gmx_ana_selparam_t *param)
429 t_methoddata_kwreal *data;
436 * \param[in] top Not used.
437 * \param[in] npar Not used (should be 2).
438 * \param[in] param Method parameters (should point to \ref smparams_keyword_real).
439 * \param[in] data Should point to \ref t_methoddata_kwreal.
440 * \returns 0 (the initialization always succeeds).
443 init_kwreal(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
445 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
447 d->v = param[0].val.u.r;
448 d->n = param[1].val.nr;
449 d->r = param[1].val.u.r;
453 * See sel_updatefunc() for description of the parameters.
454 * \p data should point to a \c t_methoddata_kwreal.
456 * Does a binary search to find which atoms match the ranges in the
457 * \c t_methoddata_kwreal structure for this selection.
458 * Matching atoms are stored in \p out->u.g.
461 evaluate_keyword_real(t_topology *top, t_trxframe *fr, t_pbc *pbc,
462 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
464 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
465 int n, i, j, jmin, jmax;
470 for (i = 0; i < g->isize; ++i)
473 if (d->r[0] > val || d->r[2*n-1] < val)
479 while (jmax - jmin > 1)
481 j = jmin + (jmax - jmin) / 2;
489 if (val <= d->r[2*j+1])
496 if (val <= d->r[2*jmin+1])
498 out->u.g->index[out->u.g->isize++] = g->index[i];
504 /********************************************************************
505 * STRING KEYWORD EVALUATION
506 ********************************************************************/
509 * \param[in] npar Not used.
510 * \param param Not used.
511 * \returns Pointer to the allocated data (t_methoddata_kwstr).
513 * Allocates memory for a t_methoddata_kwstr structure.
516 init_data_kwstr(int npar, gmx_ana_selparam_t *param)
518 t_methoddata_kwstr *data = new t_methoddata_kwstr();
519 data->matchType = gmx::eStringMatchType_Auto;
524 * \param[in,out] sel Selection element to initialize.
525 * \param[in] matchType Method to use to match string values.
527 * Sets the string matching method for string keyword matching.
530 _gmx_selelem_set_kwstr_match_type(const gmx::SelectionTreeElementPointer &sel,
531 gmx::SelectionStringMatchType matchType)
533 t_methoddata_kwstr *d = (t_methoddata_kwstr *)sel->u.expr.mdata;
535 if (sel->type != SEL_EXPRESSION || !sel->u.expr.method
536 || sel->u.expr.method->name != sm_keyword_str.name)
540 d->matchType = matchType;
544 * \param[in] top Not used.
545 * \param[in] npar Not used (should be 2).
546 * \param[in] param Method parameters (should point to \ref smparams_keyword_str).
547 * \param[in] data Should point to t_methoddata_kwstr.
550 init_kwstr(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
552 t_methoddata_kwstr *d = (t_methoddata_kwstr *)data;
554 d->v = param[0].val.u.s;
555 /* Return if this is not the first time */
556 if (!d->matches.empty())
560 int n = param[1].val.nr;
561 d->matches.reserve(n);
562 for (int i = 0; i < n; ++i)
564 const char *s = param[1].val.u.s[i];
565 d->matches.push_back(StringKeywordMatchItem(d->matchType, s));
570 * \param data Data to free (should point to a t_methoddata_kwstr).
573 free_data_kwstr(void *data)
575 t_methoddata_kwstr *d = (t_methoddata_kwstr *)data;
580 * See sel_updatefunc() for description of the parameters.
581 * \p data should point to a \c t_methoddata_kwstr.
583 * Does a linear search to find which atoms match the strings in the
584 * \c t_methoddata_kwstr structure for this selection.
585 * Wildcards are allowed in the strings.
586 * Matching atoms are stored in \p out->u.g.
589 evaluate_keyword_str(t_topology *top, t_trxframe *fr, t_pbc *pbc,
590 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
592 t_methoddata_kwstr *d = (t_methoddata_kwstr *)data;
595 for (int i = 0; i < g->isize; ++i)
597 for (size_t j = 0; j < d->matches.size(); ++j)
599 if (d->matches[j].match(d->matchType, d->v[i]))
601 out->u.g->index[out->u.g->isize++] = g->index[i];
609 /********************************************************************
610 * KEYWORD EVALUATION FOR ARBITRARY GROUPS
611 ********************************************************************/
614 * \param[in] top Not used.
615 * \param[in] npar Not used.
616 * \param[in] param Not used.
617 * \param[in] data Should point to \ref t_methoddata_kweval.
618 * \returns 0 on success, a non-zero error code on return.
620 * Calls the initialization method of the wrapped keyword.
623 init_kweval(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
625 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
627 d->kwmethod->init(top, 0, NULL, d->kwmdata);
631 * \param[in] top Not used.
632 * \param[in,out] out Pointer to output data structure.
633 * \param[in,out] data Should point to \c t_methoddata_kweval.
634 * \returns 0 for success.
637 init_output_kweval(t_topology *top, gmx_ana_selvalue_t *out, void *data)
639 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
641 out->nr = d->g.isize;
642 _gmx_selvalue_reserve(out, out->nr);
646 * \param data Data to free (should point to a \c t_methoddata_kweval).
648 * Frees the memory allocated for all the members of \c t_methoddata_kweval.
651 free_data_kweval(void *data)
653 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
655 _gmx_selelem_free_method(d->kwmethod, d->kwmdata);
660 * \param[in] top Topology.
661 * \param[in] fr Current frame.
662 * \param[in] pbc PBC structure.
663 * \param data Should point to a \ref t_methoddata_kweval.
664 * \returns 0 on success, a non-zero error code on error.
666 * Creates a lookup structure that enables fast queries of whether a point
667 * is within the solid angle or not.
670 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
672 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
674 d->kwmethod->init_frame(top, fr, pbc, d->kwmdata);
678 * See sel_updatefunc() for description of the parameters.
679 * \p data should point to a \c t_methoddata_kweval.
681 * Calls the evaluation function of the wrapped keyword with the given
682 * parameters, with the exception of using \c t_methoddata_kweval::g for the
686 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc,
687 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
689 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
691 d->kwmethod->update(top, fr, pbc, &d->g, out, d->kwmdata);
695 * \param[in] method Keyword selection method to evaluate.
696 * \param[in] params Parameter that gives the group to evaluate \p method in.
697 * \param[in] scanner Scanner data structure.
698 * \returns Pointer to the created selection element (NULL on error).
700 * Creates a \ref SEL_EXPRESSION selection element that evaluates the keyword
701 * method given by \p method in the group given by \p param.
703 * The name of \p param should be empty.
705 gmx::SelectionTreeElementPointer
706 _gmx_sel_init_keyword_evaluator(gmx_ana_selmethod_t *method,
707 const gmx::SelectionParserParameterList ¶ms,
710 gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
712 sprintf(buf, "In evaluation of '%s'", method->name);
713 gmx::MessageStringContext context(errors, buf);
715 if ((method->flags & (SMETH_SINGLEVAL | SMETH_VARNUMVAL))
716 || method->outinit || method->pupdate)
718 GMX_THROW(gmx::InternalError(
719 "Unsupported keyword method for arbitrary group evaluation"));
722 gmx::SelectionTreeElementPointer sel(
723 new gmx::SelectionTreeElement(SEL_EXPRESSION));
724 _gmx_selelem_set_method(sel, method, scanner);
726 t_methoddata_kweval *data;
728 data->kwmethod = sel->u.expr.method;
729 data->kwmdata = sel->u.expr.mdata;
730 gmx_ana_index_clear(&data->g);
732 snew(sel->u.expr.method, 1);
733 memcpy(sel->u.expr.method, data->kwmethod, sizeof(gmx_ana_selmethod_t));
734 sel->u.expr.method->flags |= SMETH_VARNUMVAL;
735 sel->u.expr.method->init_data = NULL;
736 sel->u.expr.method->set_poscoll = NULL;
737 sel->u.expr.method->init = method->init ? &init_kweval : NULL;
738 sel->u.expr.method->outinit = &init_output_kweval;
739 sel->u.expr.method->free = &free_data_kweval;
740 sel->u.expr.method->init_frame = method->init_frame ? &init_frame_kweval : NULL;
741 sel->u.expr.method->update = &evaluate_kweval;
742 sel->u.expr.method->pupdate = NULL;
743 sel->u.expr.method->nparams = asize(smparams_kweval);
744 sel->u.expr.method->param = smparams_kweval;
745 _gmx_selelem_init_method_params(sel, scanner);
746 sel->u.expr.mdata = data;
748 sel->u.expr.method->param[0].val.u.g = &data->g;
750 if (!_gmx_sel_parse_params(params, sel->u.expr.method->nparams,
751 sel->u.expr.method->param, sel, scanner))
753 return gmx::SelectionTreeElementPointer();