3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
9 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11 * Copyright (c) 2001-2009, The GROMACS development team,
12 * check out http://www.gromacs.org for more information.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
33 * Implements internal selection methods for numeric and string keyword
36 * \author Teemu Murtola <teemu.murtola@cbr.su.se>
37 * \ingroup module_selection
44 #ifdef HAVE_SYS_TYPES_H
45 #include <sys/types.h> /*old Mac needs types before regex.h*/
54 #include "gromacs/legacyheaders/macros.h"
55 #include "gromacs/legacyheaders/smalloc.h"
56 #include "gromacs/legacyheaders/string2.h"
58 #include "gromacs/selection/selmethod.h"
59 #include "gromacs/utility/exceptions.h"
60 #include "gromacs/utility/messagestringcollector.h"
61 #include "gromacs/utility/stringutil.h"
64 #include "parsetree.h"
68 /** Allocates data for integer keyword evaluation. */
70 init_data_kwint(int npar, gmx_ana_selparam_t *param);
71 /** Allocates data for real keyword evaluation. */
73 init_data_kwreal(int npar, gmx_ana_selparam_t *param);
74 /** Allocates data for string keyword evaluation. */
76 init_data_kwstr(int npar, gmx_ana_selparam_t *param);
77 /** Initializes data for integer keyword evaluation. */
79 init_kwint(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
80 /** Initializes data for real keyword evaluation. */
82 init_kwreal(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
83 /** Initializes data for string keyword evaluation. */
85 init_kwstr(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
86 /** Frees the memory allocated for string keyword evaluation. */
88 free_data_kwstr(void *data);
89 /** Evaluates integer selection keywords. */
91 evaluate_keyword_int(t_topology *top, t_trxframe *fr, t_pbc *pbc,
92 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
93 /** Evaluates real selection keywords. */
95 evaluate_keyword_real(t_topology *top, t_trxframe *fr, t_pbc *pbc,
96 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
97 /** Evaluates string selection keywords. */
99 evaluate_keyword_str(t_topology *top, t_trxframe *fr, t_pbc *pbc,
100 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
103 * Data structure for integer keyword expression evaluation.
105 typedef struct t_methoddata_kwint
107 /** Array of values for the keyword. */
109 /** Number of ranges in the \p r array. */
112 * Array of sorted integer ranges to match against.
114 * Each range is made of two integers, giving the endpoints (inclusive).
115 * This field stores the pointer to the ranges allocated by the
116 * parameter parser; see \ref SPAR_RANGES for more information.
119 } t_methoddata_kwint;
122 * Data structure for real keyword expression evaluation.
124 typedef struct t_methoddata_kwreal
126 /** Array of values for the keyword. */
128 /** Number of ranges in the \p r array. */
131 * Array of sorted ranges to match against.
133 * Each range is made of two values, giving the endpoints (inclusive).
134 * This field stores the pointer to the ranges allocated by the
135 * parameter parser; see \ref SPAR_RANGES for more information.
138 } t_methoddata_kwreal;
141 * Data structure for string keyword expression evaluation.
143 typedef struct t_methoddata_kwstr
145 /** Matching type for the strings. */
146 gmx::SelectionStringMatchType matchType;
147 /** Array of values for the keyword. */
149 /** Number of elements in the \p val array. */
152 * Array of strings/regular expressions to match against.
154 struct t_methoddata_kwstr_match {
155 /** true if the expression is a regular expression, false otherwise. */
157 /** The value to match against. */
160 /** Compiled regular expression if \p bRegExp is true. */
163 /** The string if \p bRegExp is false; */
167 /**< Array of strings/regular expressions to match against.*/
168 } t_methoddata_kwstr;
170 /** Parameters for integer keyword evaluation. */
171 static gmx_ana_selparam_t smparams_keyword_int[] = {
172 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
173 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
176 /** Parameters for real keyword evaluation. */
177 static gmx_ana_selparam_t smparams_keyword_real[] = {
178 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL | SPAR_DYNAMIC},
179 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
182 /** Parameters for string keyword evaluation. */
183 static gmx_ana_selparam_t smparams_keyword_str[] = {
184 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
185 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_VARNUM},
188 /** \internal Selection method data for integer keyword evaluation. */
189 gmx_ana_selmethod_t sm_keyword_int = {
190 "kw_int", GROUP_VALUE, SMETH_SINGLEVAL,
191 asize(smparams_keyword_int), smparams_keyword_int,
198 &evaluate_keyword_int,
203 /** \internal Selection method data for real keyword evaluation. */
204 gmx_ana_selmethod_t sm_keyword_real = {
205 "kw_real", GROUP_VALUE, SMETH_SINGLEVAL,
206 asize(smparams_keyword_real), smparams_keyword_real,
213 &evaluate_keyword_real,
218 /** \internal Selection method data for string keyword evaluation. */
219 gmx_ana_selmethod_t sm_keyword_str = {
220 "kw_str", GROUP_VALUE, SMETH_SINGLEVAL,
221 asize(smparams_keyword_str), smparams_keyword_str,
228 &evaluate_keyword_str,
233 /** Initializes keyword evaluation for an arbitrary group. */
235 init_kweval(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
236 /** Initializes output for keyword evaluation in an arbitrary group. */
238 init_output_kweval(t_topology *top, gmx_ana_selvalue_t *out, void *data);
239 /** Frees the data allocated for keyword evaluation in an arbitrary group. */
241 free_data_kweval(void *data);
242 /** Initializes frame evaluation for keyword evaluation in an arbitrary group. */
244 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
245 /** Evaluates keywords in an arbitrary group. */
247 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc,
248 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
251 * Data structure for keyword evaluation in arbitrary groups.
255 /** Wrapped keyword method for evaluating the values. */
256 gmx_ana_selmethod_t *kwmethod;
257 /** Method data for \p kwmethod. */
259 /** Group in which \p kwmethod should be evaluated. */
261 } t_methoddata_kweval;
263 /** Parameters for keyword evaluation in an arbitrary group. */
264 static gmx_ana_selparam_t smparams_kweval[] = {
265 {NULL, {GROUP_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
269 /********************************************************************
270 * INTEGER KEYWORD EVALUATION
271 ********************************************************************/
274 * \param[in] npar Not used.
275 * \param param Not used.
276 * \returns Pointer to the allocated data (\ref t_methoddata_kwint).
278 * Allocates memory for a \ref t_methoddata_kwint structure.
281 init_data_kwint(int npar, gmx_ana_selparam_t *param)
283 t_methoddata_kwint *data;
290 * \param[in] top Not used.
291 * \param[in] npar Not used (should be 2).
292 * \param[in] param Method parameters (should point to \ref smparams_keyword_int).
293 * \param[in] data Should point to \ref t_methoddata_kwint.
296 init_kwint(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
298 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
300 d->v = param[0].val.u.i;
301 d->n = param[1].val.nr;
302 d->r = param[1].val.u.i;
306 * See sel_updatefunc() for description of the parameters.
307 * \p data should point to a \c t_methoddata_kwint.
309 * Does a binary search to find which atoms match the ranges in the
310 * \c t_methoddata_kwint structure for this selection.
311 * Matching atoms are stored in \p out->u.g.
314 evaluate_keyword_int(t_topology *top, t_trxframe *fr, t_pbc *pbc,
315 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
317 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
318 int n, i, j, jmin, jmax;
323 for (i = 0; i < g->isize; ++i)
326 if (d->r[0] > val || d->r[2*n-1] < val)
332 while (jmax - jmin > 1)
334 j = jmin + (jmax - jmin) / 2;
342 if (val <= d->r[2*j+1])
349 if (val <= d->r[2*jmin+1])
351 out->u.g->index[out->u.g->isize++] = g->index[i];
357 /********************************************************************
358 * REAL KEYWORD EVALUATION
359 ********************************************************************/
362 * \param[in] npar Not used.
363 * \param param Not used.
364 * \returns Pointer to the allocated data (\ref t_methoddata_kwreal).
366 * Allocates memory for a \ref t_methoddata_kwreal structure.
369 init_data_kwreal(int npar, gmx_ana_selparam_t *param)
371 t_methoddata_kwreal *data;
378 * \param[in] top Not used.
379 * \param[in] npar Not used (should be 2).
380 * \param[in] param Method parameters (should point to \ref smparams_keyword_real).
381 * \param[in] data Should point to \ref t_methoddata_kwreal.
382 * \returns 0 (the initialization always succeeds).
385 init_kwreal(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
387 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
389 d->v = param[0].val.u.r;
390 d->n = param[1].val.nr;
391 d->r = param[1].val.u.r;
395 * See sel_updatefunc() for description of the parameters.
396 * \p data should point to a \c t_methoddata_kwreal.
398 * Does a binary search to find which atoms match the ranges in the
399 * \c t_methoddata_kwreal structure for this selection.
400 * Matching atoms are stored in \p out->u.g.
403 evaluate_keyword_real(t_topology *top, t_trxframe *fr, t_pbc *pbc,
404 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
406 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
407 int n, i, j, jmin, jmax;
412 for (i = 0; i < g->isize; ++i)
415 if (d->r[0] > val || d->r[2*n-1] < val)
421 while (jmax - jmin > 1)
423 j = jmin + (jmax - jmin) / 2;
431 if (val <= d->r[2*j+1])
438 if (val <= d->r[2*jmin+1])
440 out->u.g->index[out->u.g->isize++] = g->index[i];
446 /********************************************************************
447 * STRING KEYWORD EVALUATION
448 ********************************************************************/
451 * \param[in] npar Not used.
452 * \param param Not used.
453 * \returns Pointer to the allocated data (\ref t_methoddata_kwstr).
455 * Allocates memory for a \ref t_methoddata_kwstr structure.
458 init_data_kwstr(int npar, gmx_ana_selparam_t *param)
460 t_methoddata_kwstr *data;
463 data->matchType = gmx::eStringMatchType_Auto;
468 * \param[in,out] sel Selection element to initialize.
469 * \param[in] matchType Method to use to match string values.
471 * Sets the string matching method for string keyword matching.
474 _gmx_selelem_set_kwstr_match_type(const gmx::SelectionTreeElementPointer &sel,
475 gmx::SelectionStringMatchType matchType)
477 t_methoddata_kwstr *d = (t_methoddata_kwstr *)sel->u.expr.mdata;
479 if (sel->type != SEL_EXPRESSION || !sel->u.expr.method
480 || sel->u.expr.method->name != sm_keyword_str.name)
484 d->matchType = matchType;
488 * \param[in] top Not used.
489 * \param[in] npar Not used (should be 2).
490 * \param[in] param Method parameters (should point to \ref smparams_keyword_str).
491 * \param[in] data Should point to \ref t_methoddata_kwstr.
494 init_kwstr(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
496 t_methoddata_kwstr *d = (t_methoddata_kwstr *)data;
499 d->v = param[0].val.u.s;
500 d->n = param[1].val.nr;
501 /* Return if this is not the first time */
507 for (i = 0; i < d->n; ++i)
509 const char *s = param[1].val.u.s[i];
510 bool bRegExp = (d->matchType == gmx::eStringMatchType_RegularExpression);
511 if (d->matchType == gmx::eStringMatchType_Auto)
513 for (size_t j = 0; j < strlen(s); ++j)
515 if (ispunct(s[j]) && s[j] != '?' && s[j] != '*')
525 std::string buf(gmx::formatString("^%s$", s));
526 if (regcomp(&d->m[i].u.r, buf.c_str(), REG_EXTENDED | REG_NOSUB) != 0)
528 GMX_THROW(gmx::InvalidInputError(gmx::formatString(
529 "Error in regular expression \"%s\"", s)));
532 GMX_THROW(gmx::InvalidInputError(gmx::formatString(
533 "No regular expression support, cannot match \"%s\"", s)));
540 d->m[i].bRegExp = bRegExp;
545 * \param data Data to free (should point to a \ref t_methoddata_kwstr).
547 * Frees the memory allocated for t_methoddata_kwstr::val.
550 free_data_kwstr(void *data)
552 t_methoddata_kwstr *d = (t_methoddata_kwstr *)data;
555 for (i = 0; i < d->n; ++i)
560 /* This branch should only be taken if regular expressions
561 * are available, but the ifdef is still needed. */
562 regfree(&d->m[i].u.r);
571 * See sel_updatefunc() for description of the parameters.
572 * \p data should point to a \c t_methoddata_kwstr.
574 * Does a linear search to find which atoms match the strings in the
575 * \c t_methoddata_kwstr structure for this selection.
576 * Wildcards are allowed in the strings.
577 * Matching atoms are stored in \p out->u.g.
580 evaluate_keyword_str(t_topology *top, t_trxframe *fr, t_pbc *pbc,
581 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
583 t_methoddata_kwstr *d = (t_methoddata_kwstr *)data;
588 for (i = 0; i < g->isize; ++i)
591 for (j = 0; j < d->n && !bFound; ++j)
596 /* This branch should only be taken if regular expressions
597 * are available, but the ifdef is still needed. */
598 bFound = (regexec(&d->m[j].u.r, d->v[i], 0, NULL, 0) == 0);
601 else if (d->matchType == gmx::eStringMatchType_Exact)
603 bFound = (strcmp(d->m[j].u.s, d->v[i]) == 0);
607 bFound = (gmx_wcmatch(d->m[j].u.s, d->v[i]) == 0);
612 out->u.g->index[out->u.g->isize++] = g->index[i];
618 /********************************************************************
619 * KEYWORD EVALUATION FOR ARBITRARY GROUPS
620 ********************************************************************/
623 * \param[in] top Not used.
624 * \param[in] npar Not used.
625 * \param[in] param Not used.
626 * \param[in] data Should point to \ref t_methoddata_kweval.
627 * \returns 0 on success, a non-zero error code on return.
629 * Calls the initialization method of the wrapped keyword.
632 init_kweval(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
634 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
636 d->kwmethod->init(top, 0, NULL, d->kwmdata);
640 * \param[in] top Not used.
641 * \param[in,out] out Pointer to output data structure.
642 * \param[in,out] data Should point to \c t_methoddata_kweval.
643 * \returns 0 for success.
646 init_output_kweval(t_topology *top, gmx_ana_selvalue_t *out, void *data)
648 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
650 out->nr = d->g.isize;
651 _gmx_selvalue_reserve(out, out->nr);
655 * \param data Data to free (should point to a \c t_methoddata_kweval).
657 * Frees the memory allocated for all the members of \c t_methoddata_kweval.
660 free_data_kweval(void *data)
662 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
664 _gmx_selelem_free_method(d->kwmethod, d->kwmdata);
669 * \param[in] top Topology.
670 * \param[in] fr Current frame.
671 * \param[in] pbc PBC structure.
672 * \param data Should point to a \ref t_methoddata_kweval.
673 * \returns 0 on success, a non-zero error code on error.
675 * Creates a lookup structure that enables fast queries of whether a point
676 * is within the solid angle or not.
679 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
681 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
683 d->kwmethod->init_frame(top, fr, pbc, d->kwmdata);
687 * See sel_updatefunc() for description of the parameters.
688 * \p data should point to a \c t_methoddata_kweval.
690 * Calls the evaluation function of the wrapped keyword with the given
691 * parameters, with the exception of using \c t_methoddata_kweval::g for the
695 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc,
696 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
698 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
700 d->kwmethod->update(top, fr, pbc, &d->g, out, d->kwmdata);
704 * \param[in] method Keyword selection method to evaluate.
705 * \param[in] params Parameter that gives the group to evaluate \p method in.
706 * \param[in] scanner Scanner data structure.
707 * \returns Pointer to the created selection element (NULL on error).
709 * Creates a \ref SEL_EXPRESSION selection element that evaluates the keyword
710 * method given by \p method in the group given by \p param.
712 * The name of \p param should be empty.
714 gmx::SelectionTreeElementPointer
715 _gmx_sel_init_keyword_evaluator(gmx_ana_selmethod_t *method,
716 const gmx::SelectionParserParameterList ¶ms,
719 gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
721 sprintf(buf, "In evaluation of '%s'", method->name);
722 gmx::MessageStringContext context(errors, buf);
724 if ((method->flags & (SMETH_SINGLEVAL | SMETH_VARNUMVAL))
725 || method->outinit || method->pupdate)
727 GMX_THROW(gmx::InternalError(
728 "Unsupported keyword method for arbitrary group evaluation"));
731 gmx::SelectionTreeElementPointer sel(
732 new gmx::SelectionTreeElement(SEL_EXPRESSION));
733 _gmx_selelem_set_method(sel, method, scanner);
735 t_methoddata_kweval *data;
737 data->kwmethod = sel->u.expr.method;
738 data->kwmdata = sel->u.expr.mdata;
739 gmx_ana_index_clear(&data->g);
741 snew(sel->u.expr.method, 1);
742 memcpy(sel->u.expr.method, data->kwmethod, sizeof(gmx_ana_selmethod_t));
743 sel->u.expr.method->flags |= SMETH_VARNUMVAL;
744 sel->u.expr.method->init_data = NULL;
745 sel->u.expr.method->set_poscoll = NULL;
746 sel->u.expr.method->init = method->init ? &init_kweval : NULL;
747 sel->u.expr.method->outinit = &init_output_kweval;
748 sel->u.expr.method->free = &free_data_kweval;
749 sel->u.expr.method->init_frame = method->init_frame ? &init_frame_kweval : NULL;
750 sel->u.expr.method->update = &evaluate_kweval;
751 sel->u.expr.method->pupdate = NULL;
752 sel->u.expr.method->nparams = asize(smparams_kweval);
753 sel->u.expr.method->param = smparams_kweval;
754 _gmx_selelem_init_method_params(sel, scanner);
755 sel->u.expr.mdata = data;
757 sel->u.expr.method->param[0].val.u.g = &data->g;
759 if (!_gmx_sel_parse_params(params, sel->u.expr.method->nparams,
760 sel->u.expr.method->param, sel, scanner))
762 return gmx::SelectionTreeElementPointer();