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*/
56 #include "gromacs/fatalerror/errorcodes.h"
57 #include "gromacs/fatalerror/messagestringcollector.h"
58 #include "gromacs/selection/selmethod.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;
138 * Data structure for string keyword expression evaluation.
140 typedef struct t_methoddata_kwstr
142 /** Array of values for the keyword. */
144 /** Number of elements in the \p val array. */
147 * Array of strings/regular expressions to match against.
149 struct t_methoddata_kwstr_match {
150 /** true if the expression is a regular expression, false otherwise. */
152 /** The value to match against. */
155 /** Compiled regular expression if \p bRegExp is true. */
158 /** The string if \p bRegExp is false; */
162 } t_methoddata_kwstr;
164 /** Parameters for integer keyword evaluation. */
165 static gmx_ana_selparam_t smparams_keyword_int[] = {
166 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
167 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
170 /** Parameters for real keyword evaluation. */
171 static gmx_ana_selparam_t smparams_keyword_real[] = {
172 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL | SPAR_DYNAMIC},
173 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
176 /** Parameters for string keyword evaluation. */
177 static gmx_ana_selparam_t smparams_keyword_str[] = {
178 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
179 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_VARNUM},
182 /** \internal Selection method data for integer keyword evaluation. */
183 gmx_ana_selmethod_t sm_keyword_int = {
184 "kw_int", GROUP_VALUE, SMETH_SINGLEVAL,
185 asize(smparams_keyword_int), smparams_keyword_int,
192 &evaluate_keyword_int,
197 /** \internal Selection method data for real keyword evaluation. */
198 gmx_ana_selmethod_t sm_keyword_real = {
199 "kw_real", GROUP_VALUE, SMETH_SINGLEVAL,
200 asize(smparams_keyword_real), smparams_keyword_real,
207 &evaluate_keyword_real,
212 /** \internal Selection method data for string keyword evaluation. */
213 gmx_ana_selmethod_t sm_keyword_str = {
214 "kw_str", GROUP_VALUE, SMETH_SINGLEVAL,
215 asize(smparams_keyword_str), smparams_keyword_str,
222 &evaluate_keyword_str,
227 /** Initializes keyword evaluation for an arbitrary group. */
229 init_kweval(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
230 /** Initializes output for keyword evaluation in an arbitrary group. */
232 init_output_kweval(t_topology *top, gmx_ana_selvalue_t *out, void *data);
233 /** Frees the data allocated for keyword evaluation in an arbitrary group. */
235 free_data_kweval(void *data);
236 /** Initializes frame evaluation for keyword evaluation in an arbitrary group. */
238 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
239 /** Evaluates keywords in an arbitrary group. */
241 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc,
242 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
245 * Data structure for keyword evaluation in arbitrary groups.
249 /** Wrapped keyword method for evaluating the values. */
250 gmx_ana_selmethod_t *kwmethod;
251 /** Method data for \p kwmethod. */
253 /** Group in which \p kwmethod should be evaluated. */
255 } t_methoddata_kweval;
257 /** Parameters for keyword evaluation in an arbitrary group. */
258 static gmx_ana_selparam_t smparams_kweval[] = {
259 {NULL, {GROUP_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
263 /********************************************************************
264 * INTEGER KEYWORD EVALUATION
265 ********************************************************************/
268 * \param[in] npar Not used.
269 * \param param Not used.
270 * \returns Pointer to the allocated data (\ref t_methoddata_kwint).
272 * Allocates memory for a \ref t_methoddata_kwint structure.
275 init_data_kwint(int npar, gmx_ana_selparam_t *param)
277 t_methoddata_kwint *data;
284 * \param[in] top Not used.
285 * \param[in] npar Not used (should be 2).
286 * \param[in] param Method parameters (should point to \ref smparams_keyword_int).
287 * \param[in] data Should point to \ref t_methoddata_kwint.
290 init_kwint(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
292 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
294 d->v = param[0].val.u.i;
295 d->n = param[1].val.nr;
296 d->r = param[1].val.u.i;
300 * See sel_updatefunc() for description of the parameters.
301 * \p data should point to a \c t_methoddata_kwint.
303 * Does a binary search to find which atoms match the ranges in the
304 * \c t_methoddata_kwint structure for this selection.
305 * Matching atoms are stored in \p out->u.g.
308 evaluate_keyword_int(t_topology *top, t_trxframe *fr, t_pbc *pbc,
309 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
311 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
312 int n, i, j, jmin, jmax;
317 for (i = 0; i < g->isize; ++i)
320 if (d->r[0] > val || d->r[2*n-1] < val)
326 while (jmax - jmin > 1)
328 j = jmin + (jmax - jmin) / 2;
336 if (val <= d->r[2*j+1])
343 if (val <= d->r[2*jmin+1])
345 out->u.g->index[out->u.g->isize++] = g->index[i];
351 /********************************************************************
352 * REAL KEYWORD EVALUATION
353 ********************************************************************/
356 * \param[in] npar Not used.
357 * \param param Not used.
358 * \returns Pointer to the allocated data (\ref t_methoddata_kwreal).
360 * Allocates memory for a \ref t_methoddata_kwreal structure.
363 init_data_kwreal(int npar, gmx_ana_selparam_t *param)
365 t_methoddata_kwreal *data;
372 * \param[in] top Not used.
373 * \param[in] npar Not used (should be 2).
374 * \param[in] param Method parameters (should point to \ref smparams_keyword_real).
375 * \param[in] data Should point to \ref t_methoddata_kwreal.
376 * \returns 0 (the initialization always succeeds).
379 init_kwreal(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
381 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
383 d->v = param[0].val.u.r;
384 d->n = param[1].val.nr;
385 d->r = param[1].val.u.r;
389 * See sel_updatefunc() for description of the parameters.
390 * \p data should point to a \c t_methoddata_kwreal.
392 * Does a binary search to find which atoms match the ranges in the
393 * \c t_methoddata_kwreal structure for this selection.
394 * Matching atoms are stored in \p out->u.g.
397 evaluate_keyword_real(t_topology *top, t_trxframe *fr, t_pbc *pbc,
398 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
400 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
401 int n, i, j, jmin, jmax;
406 for (i = 0; i < g->isize; ++i)
409 if (d->r[0] > val || d->r[2*n-1] < val)
415 while (jmax - jmin > 1)
417 j = jmin + (jmax - jmin) / 2;
425 if (val <= d->r[2*j+1])
432 if (val <= d->r[2*jmin+1])
434 out->u.g->index[out->u.g->isize++] = g->index[i];
440 /********************************************************************
441 * STRING KEYWORD EVALUATION
442 ********************************************************************/
445 * \param[in] npar Not used.
446 * \param param Not used.
447 * \returns Pointer to the allocated data (\ref t_methoddata_kwstr).
449 * Allocates memory for a \ref t_methoddata_kwstr structure.
452 init_data_kwstr(int npar, gmx_ana_selparam_t *param)
454 t_methoddata_kwstr *data;
461 * \param[in] top Not used.
462 * \param[in] npar Not used (should be 2).
463 * \param[in] param Method parameters (should point to \ref smparams_keyword_str).
464 * \param[in] data Should point to \ref t_methoddata_kwstr.
467 init_kwstr(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
469 t_methoddata_kwstr *d = (t_methoddata_kwstr *)data;
476 d->v = param[0].val.u.s;
477 d->n = param[1].val.nr;
478 /* Return if this is not the first time */
484 for (i = 0; i < d->n; ++i)
486 s = param[1].val.u.s[i];
488 for (j = 0; j < strlen(s); ++j)
490 if (ispunct(s[j]) && s[j] != '?' && s[j] != '*')
498 // TODO: Get rid of these prints to stderr
500 snew(buf, strlen(s) + 3);
501 sprintf(buf, "^%s$", s);
502 if (regcomp(&d->m[i].u.r, buf, REG_EXTENDED | REG_NOSUB))
505 fprintf(stderr, "WARNING: error in regular expression,\n"
506 " will match '%s' as a simple string\n", s);
511 fprintf(stderr, "WARNING: no regular expressions support,\n"
512 " will match '%s' as a simple string\n", s);
519 d->m[i].bRegExp = bRegExp;
524 * \param data Data to free (should point to a \ref t_methoddata_kwstr).
526 * Frees the memory allocated for t_methoddata_kwstr::val.
529 free_data_kwstr(void *data)
531 t_methoddata_kwstr *d = (t_methoddata_kwstr *)data;
534 for (i = 0; i < d->n; ++i)
539 /* This branch should only be taken if regular expressions
540 * are available, but the ifdef is still needed. */
541 regfree(&d->m[i].u.r);
549 * See sel_updatefunc() for description of the parameters.
550 * \p data should point to a \c t_methoddata_kwstr.
552 * Does a linear search to find which atoms match the strings in the
553 * \c t_methoddata_kwstr structure for this selection.
554 * Wildcards are allowed in the strings.
555 * Matching atoms are stored in \p out->u.g.
558 evaluate_keyword_str(t_topology *top, t_trxframe *fr, t_pbc *pbc,
559 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
561 t_methoddata_kwstr *d = (t_methoddata_kwstr *)data;
566 for (i = 0; i < g->isize; ++i)
569 for (j = 0; j < d->n && !bFound; ++j)
574 /* This branch should only be taken if regular expressions
575 * are available, but the ifdef is still needed. */
576 if (!regexec(&d->m[j].u.r, d->v[i], 0, NULL, 0))
584 if (gmx_wcmatch(d->m[j].u.s, d->v[i]) == 0)
592 out->u.g->index[out->u.g->isize++] = g->index[i];
598 /********************************************************************
599 * KEYWORD EVALUATION FOR ARBITRARY GROUPS
600 ********************************************************************/
603 * \param[in] top Not used.
604 * \param[in] npar Not used.
605 * \param[in] param Not used.
606 * \param[in] data Should point to \ref t_methoddata_kweval.
607 * \returns 0 on success, a non-zero error code on return.
609 * Calls the initialization method of the wrapped keyword.
612 init_kweval(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
614 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
616 d->kwmethod->init(top, 0, NULL, d->kwmdata);
620 * \param[in] top Not used.
621 * \param[in,out] out Pointer to output data structure.
622 * \param[in,out] data Should point to \c t_methoddata_kweval.
623 * \returns 0 for success.
626 init_output_kweval(t_topology *top, gmx_ana_selvalue_t *out, void *data)
628 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
630 out->nr = d->g.isize;
634 * \param data Data to free (should point to a \c t_methoddata_kweval).
636 * Frees the memory allocated for all the members of \c t_methoddata_kweval.
639 free_data_kweval(void *data)
641 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
643 _gmx_selelem_free_method(d->kwmethod, d->kwmdata);
647 * \param[in] top Topology.
648 * \param[in] fr Current frame.
649 * \param[in] pbc PBC structure.
650 * \param data Should point to a \ref t_methoddata_kweval.
651 * \returns 0 on success, a non-zero error code on error.
653 * Creates a lookup structure that enables fast queries of whether a point
654 * is within the solid angle or not.
657 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
659 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
661 d->kwmethod->init_frame(top, fr, pbc, d->kwmdata);
665 * See sel_updatefunc() for description of the parameters.
666 * \p data should point to a \c t_methoddata_kweval.
668 * Calls the evaluation function of the wrapped keyword with the given
669 * parameters, with the exception of using \c t_methoddata_kweval::g for the
673 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc,
674 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
676 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
678 d->kwmethod->update(top, fr, pbc, &d->g, out, d->kwmdata);
682 * \param[out] selp Pointer to receive a pointer to the created selection
683 * element (set to NULL on error).
684 * \param[in] method Keyword selection method to evaluate.
685 * \param[in] param Parameter that gives the group to evaluate \p method in.
686 * \param[in] scanner Scanner data structure.
687 * \returns 0 on success, non-zero error code on error.
689 * Creates a \ref SEL_EXPRESSION selection element (pointer put in \c *selp)
690 * that evaluates the keyword method given by \p method in the group given by
694 _gmx_sel_init_keyword_evaluator(t_selelem **selp, gmx_ana_selmethod_t *method,
695 t_selexpr_param *param, void *scanner)
698 t_methoddata_kweval *data;
700 gmx::MessageStringCollector *errors = _gmx_sel_lexer_error_reporter(scanner);
702 sprintf(buf, "In evaluation of '%s'", method->name);
703 gmx::MessageStringContext context(errors, buf);
705 if ((method->flags & (SMETH_SINGLEVAL | SMETH_VARNUMVAL))
706 || method->outinit || method->pupdate)
708 _gmx_selexpr_free_params(param);
709 GMX_ERROR(gmx::eeInternalError,
710 "Unsupported keyword method for arbitrary group evaluation");
714 sel = _gmx_selelem_create(SEL_EXPRESSION);
715 _gmx_selelem_set_method(sel, method, scanner);
718 data->kwmethod = sel->u.expr.method;
719 data->kwmdata = sel->u.expr.mdata;
720 gmx_ana_index_clear(&data->g);
722 snew(sel->u.expr.method, 1);
723 memcpy(sel->u.expr.method, data->kwmethod, sizeof(gmx_ana_selmethod_t));
724 sel->u.expr.method->flags |= SMETH_VARNUMVAL;
725 sel->u.expr.method->init_data = NULL;
726 sel->u.expr.method->set_poscoll = NULL;
727 sel->u.expr.method->init = method->init ? &init_kweval : NULL;
728 sel->u.expr.method->outinit = &init_output_kweval;
729 sel->u.expr.method->free = &free_data_kweval;
730 sel->u.expr.method->init_frame = method->init_frame ? &init_frame_kweval : NULL;
731 sel->u.expr.method->update = &evaluate_kweval;
732 sel->u.expr.method->pupdate = NULL;
733 sel->u.expr.method->nparams = asize(smparams_kweval);
734 sel->u.expr.method->param = smparams_kweval;
735 _gmx_selelem_init_method_params(sel, scanner);
736 sel->u.expr.mdata = data;
738 sel->u.expr.method->param[0].val.u.g = &data->g;
742 if (!_gmx_sel_parse_params(param, sel->u.expr.method->nparams,
743 sel->u.expr.method->param, sel, scanner))
745 _gmx_selelem_free(sel);