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"
54 #include "gromacs/selection/selmethod.h"
55 #include "gromacs/utility/cstringutil.h"
56 #include "gromacs/utility/exceptions.h"
57 #include "gromacs/utility/gmxregex.h"
58 #include "gromacs/utility/messagestringcollector.h"
59 #include "gromacs/utility/stringutil.h"
60 #include "gromacs/utility/smalloc.h"
63 #include "parsetree.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())
221 GMX_THROW(gmx::InvalidInputError(gmx::formatString(
222 "No regular expression support, "
223 "cannot match \"%s\"", str)));
225 regex_.reset(new gmx::Regex(str));
230 * Checks whether this item matches a string.
232 * \param[in] matchType String matching type.
233 * \param[in] value String to match.
234 * \returns true if this item matches \p value.
236 bool match(gmx::SelectionStringMatchType matchType,
237 const char *value) const
239 if (matchType == gmx::eStringMatchType_Exact)
241 return str_ == value;
245 return gmx::regexMatch(value, *regex_);
249 return gmx_wcmatch(str_.c_str(), value) == 0;
254 //! The raw string passed for the matcher.
256 //! Regular expression compiled from \p str_, if applicable.
257 boost::shared_ptr<gmx::Regex> regex_;
261 * Data structure for string keyword expression evaluation.
263 struct t_methoddata_kwstr
265 /** Matching type for the strings. */
266 gmx::SelectionStringMatchType matchType;
267 /** Array of values for the keyword. */
269 /** Array of strings/regular expressions to match against.*/
270 std::vector<StringKeywordMatchItem> matches;
275 /** Parameters for integer keyword evaluation. */
276 static gmx_ana_selparam_t smparams_keyword_int[] = {
277 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
278 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
281 /** Parameters for real keyword evaluation. */
282 static gmx_ana_selparam_t smparams_keyword_real[] = {
283 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL | SPAR_DYNAMIC},
284 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
287 /** Parameters for string keyword evaluation. */
288 static gmx_ana_selparam_t smparams_keyword_str[] = {
289 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
290 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_VARNUM},
293 /** Selection method data for integer keyword evaluation. */
294 gmx_ana_selmethod_t sm_keyword_int = {
295 "kw_int", GROUP_VALUE, SMETH_SINGLEVAL,
296 asize(smparams_keyword_int), smparams_keyword_int,
303 &evaluate_keyword_int,
308 /** Selection method data for real keyword evaluation. */
309 gmx_ana_selmethod_t sm_keyword_real = {
310 "kw_real", GROUP_VALUE, SMETH_SINGLEVAL,
311 asize(smparams_keyword_real), smparams_keyword_real,
318 &evaluate_keyword_real,
323 /** Selection method data for string keyword evaluation. */
324 gmx_ana_selmethod_t sm_keyword_str = {
325 "kw_str", GROUP_VALUE, SMETH_SINGLEVAL,
326 asize(smparams_keyword_str), smparams_keyword_str,
333 &evaluate_keyword_str,
339 * Initializes keyword evaluation for an arbitrary group.
341 * \param[in] top Not used.
342 * \param[in] npar Not used.
343 * \param[in] param Not used.
344 * \param[in] data Should point to \ref t_methoddata_kweval.
345 * \returns 0 on success, a non-zero error code on return.
347 * Calls the initialization method of the wrapped keyword.
350 init_kweval(t_topology *top, int npar, gmx_ana_selparam_t * param, void *data);
352 * Initializes output for keyword evaluation in an arbitrary group.
354 * \param[in] top Not used.
355 * \param[in,out] out Pointer to output data structure.
356 * \param[in,out] data Should point to \c t_methoddata_kweval.
357 * \returns 0 for success.
360 init_output_kweval(t_topology *top, gmx_ana_selvalue_t *out, void *data);
361 /** Frees the data allocated for keyword evaluation in an arbitrary group. */
363 free_data_kweval(void *data);
364 /** Initializes frame evaluation for keyword evaluation in an arbitrary group. */
366 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
367 /** Evaluates keywords in an arbitrary group. */
369 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, gmx_ana_index_t * /* g */, gmx_ana_selvalue_t *out, void *data);
372 * Data structure for keyword evaluation in arbitrary groups.
376 /** Wrapped keyword method for evaluating the values. */
377 gmx_ana_selmethod_t *kwmethod;
378 /** Method data for \p kwmethod. */
380 /** Group in which \p kwmethod should be evaluated. */
382 } t_methoddata_kweval;
384 /** Parameters for keyword evaluation in an arbitrary group. */
385 static gmx_ana_selparam_t smparams_kweval[] = {
386 {NULL, {GROUP_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
390 /********************************************************************
391 * INTEGER KEYWORD EVALUATION
392 ********************************************************************/
395 init_data_kwint(int /* npar */, gmx_ana_selparam_t * /* param */)
397 t_methoddata_kwint *data;
404 init_kwint(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
406 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
408 d->v = param[0].val.u.i;
409 d->n = param[1].val.nr;
410 d->r = param[1].val.u.i;
414 * See sel_updatefunc() for description of the parameters.
415 * \p data should point to a \c t_methoddata_kwint.
417 * Does a binary search to find which atoms match the ranges in the
418 * \c t_methoddata_kwint structure for this selection.
419 * Matching atoms are stored in \p out->u.g.
422 evaluate_keyword_int(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
423 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
425 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
426 int n, i, j, jmin, jmax;
431 for (i = 0; i < g->isize; ++i)
434 if (d->r[0] > val || d->r[2*n-1] < val)
440 while (jmax - jmin > 1)
442 j = jmin + (jmax - jmin) / 2;
450 if (val <= d->r[2*j+1])
457 if (val <= d->r[2*jmin+1])
459 out->u.g->index[out->u.g->isize++] = g->index[i];
465 /********************************************************************
466 * REAL KEYWORD EVALUATION
467 ********************************************************************/
470 init_data_kwreal(int /* npar */, gmx_ana_selparam_t * /* param */)
472 t_methoddata_kwreal *data;
479 init_kwreal(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
481 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
483 d->v = param[0].val.u.r;
484 d->n = param[1].val.nr;
485 d->r = param[1].val.u.r;
489 * See sel_updatefunc() for description of the parameters.
490 * \p data should point to a \c t_methoddata_kwreal.
492 * Does a binary search to find which atoms match the ranges in the
493 * \c t_methoddata_kwreal structure for this selection.
494 * Matching atoms are stored in \p out->u.g.
497 evaluate_keyword_real(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
498 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
500 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
501 int n, i, j, jmin, jmax;
506 for (i = 0; i < g->isize; ++i)
509 if (d->r[0] > val || d->r[2*n-1] < val)
515 while (jmax - jmin > 1)
517 j = jmin + (jmax - jmin) / 2;
525 if (val <= d->r[2*j+1])
532 if (val <= d->r[2*jmin+1])
534 out->u.g->index[out->u.g->isize++] = g->index[i];
540 /********************************************************************
541 * STRING KEYWORD EVALUATION
542 ********************************************************************/
545 init_data_kwstr(int /* npar */, gmx_ana_selparam_t * /* param */)
547 t_methoddata_kwstr *data = new t_methoddata_kwstr();
548 data->matchType = gmx::eStringMatchType_Auto;
553 * \param[in,out] sel Selection element to initialize.
554 * \param[in] matchType Method to use to match string values.
556 * Sets the string matching method for string keyword matching.
559 _gmx_selelem_set_kwstr_match_type(const gmx::SelectionTreeElementPointer &sel,
560 gmx::SelectionStringMatchType matchType)
562 t_methoddata_kwstr *d = static_cast<t_methoddata_kwstr *>(sel->u.expr.mdata);
564 if (sel->type != SEL_EXPRESSION || !sel->u.expr.method
565 || sel->u.expr.method->name != sm_keyword_str.name)
569 d->matchType = matchType;
573 init_kwstr(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
575 t_methoddata_kwstr *d = static_cast<t_methoddata_kwstr *>(data);
577 d->v = param[0].val.u.s;
578 /* Return if this is not the first time */
579 if (!d->matches.empty())
583 int n = param[1].val.nr;
584 d->matches.reserve(n);
585 for (int i = 0; i < n; ++i)
587 const char *s = param[1].val.u.s[i];
588 d->matches.push_back(StringKeywordMatchItem(d->matchType, s));
593 * \param data Data to free (should point to a t_methoddata_kwstr).
596 free_data_kwstr(void *data)
598 t_methoddata_kwstr *d = static_cast<t_methoddata_kwstr *>(data);
603 * See sel_updatefunc() for description of the parameters.
604 * \p data should point to a \c t_methoddata_kwstr.
606 * Does a linear search to find which atoms match the strings in the
607 * \c t_methoddata_kwstr structure for this selection.
608 * Wildcards are allowed in the strings.
609 * Matching atoms are stored in \p out->u.g.
612 evaluate_keyword_str(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
613 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
615 t_methoddata_kwstr *d = static_cast<t_methoddata_kwstr *>(data);
618 for (int i = 0; i < g->isize; ++i)
620 for (size_t j = 0; j < d->matches.size(); ++j)
622 if (d->matches[j].match(d->matchType, d->v[i]))
624 out->u.g->index[out->u.g->isize++] = g->index[i];
632 /********************************************************************
633 * KEYWORD EVALUATION FOR ARBITRARY GROUPS
634 ********************************************************************/
637 init_kweval(t_topology *top, int /* npar */, gmx_ana_selparam_t * /* param */, void *data)
639 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
641 d->kwmethod->init(top, 0, NULL, d->kwmdata);
645 init_output_kweval(t_topology * /* top */, gmx_ana_selvalue_t *out, void *data)
647 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
649 out->nr = d->g.isize;
650 _gmx_selvalue_reserve(out, out->nr);
654 * \param data Data to free (should point to a \c t_methoddata_kweval).
656 * Frees the memory allocated for all the members of \c t_methoddata_kweval.
659 free_data_kweval(void *data)
661 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
663 _gmx_selelem_free_method(d->kwmethod, d->kwmdata);
668 * \param[in] top Topology.
669 * \param[in] fr Current frame.
670 * \param[in] pbc PBC structure.
671 * \param data Should point to a \ref t_methoddata_kweval.
672 * \returns 0 on success, a non-zero error code on error.
674 * Creates a lookup structure that enables fast queries of whether a point
675 * is within the solid angle or not.
678 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
680 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
682 d->kwmethod->init_frame(top, fr, pbc, d->kwmdata);
686 * See sel_updatefunc() for description of the parameters.
687 * \p data should point to a \c t_methoddata_kweval.
689 * Calls the evaluation function of the wrapped keyword with the given
690 * parameters, with the exception of using \c t_methoddata_kweval::g for the
694 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc,
695 gmx_ana_index_t * /* g */, gmx_ana_selvalue_t *out, void *data)
697 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
699 d->kwmethod->update(top, fr, pbc, &d->g, out, d->kwmdata);
703 * \param[in] method Keyword selection method to evaluate.
704 * \param[in] params Parameter that gives the group to evaluate \p method in.
705 * \param[in] scanner Scanner data structure.
706 * \returns Pointer to the created selection element (NULL on error).
708 * Creates a \ref SEL_EXPRESSION selection element that evaluates the keyword
709 * method given by \p method in the group given by \p param.
711 * The name of \p param should be empty.
713 gmx::SelectionTreeElementPointer
714 _gmx_sel_init_keyword_evaluator(gmx_ana_selmethod_t *method,
715 const gmx::SelectionParserParameterList ¶ms,
718 gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
720 sprintf(buf, "In evaluation of '%s'", method->name);
721 gmx::MessageStringContext context(errors, buf);
723 if ((method->flags & (SMETH_SINGLEVAL | SMETH_VARNUMVAL))
724 || method->outinit || method->pupdate)
726 GMX_THROW(gmx::InternalError(
727 "Unsupported keyword method for arbitrary group evaluation"));
730 gmx::SelectionTreeElementPointer sel(
731 new gmx::SelectionTreeElement(SEL_EXPRESSION));
732 _gmx_selelem_set_method(sel, method, scanner);
734 t_methoddata_kweval *data;
736 data->kwmethod = sel->u.expr.method;
737 data->kwmdata = sel->u.expr.mdata;
738 gmx_ana_index_clear(&data->g);
740 snew(sel->u.expr.method, 1);
741 memcpy(sel->u.expr.method, data->kwmethod, sizeof(gmx_ana_selmethod_t));
742 sel->u.expr.method->flags |= SMETH_VARNUMVAL;
743 sel->u.expr.method->init_data = NULL;
744 sel->u.expr.method->set_poscoll = NULL;
745 sel->u.expr.method->init = method->init ? &init_kweval : NULL;
746 sel->u.expr.method->outinit = &init_output_kweval;
747 sel->u.expr.method->free = &free_data_kweval;
748 sel->u.expr.method->init_frame = method->init_frame ? &init_frame_kweval : NULL;
749 sel->u.expr.method->update = &evaluate_kweval;
750 sel->u.expr.method->pupdate = NULL;
751 sel->u.expr.method->nparams = asize(smparams_kweval);
752 sel->u.expr.method->param = smparams_kweval;
753 _gmx_selelem_init_method_params(sel, scanner);
754 sel->u.expr.mdata = data;
756 sel->u.expr.method->param[0].val.u.g = &data->g;
758 if (!_gmx_sel_parse_params(params, sel->u.expr.method->nparams,
759 sel->u.expr.method->param, sel, scanner))
761 return gmx::SelectionTreeElementPointer();