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
32 * \brief Implementations of internal selection methods for numeric and
33 * string keyword evaluation.
40 #ifdef HAVE_SYS_TYPES_H
41 #include <sys/types.h> /*old Mac needs types before regex.h*/
48 #include <gmx_fatal.h>
53 #include <selmethod.h>
56 #include "parsetree.h"
59 /** Allocates data for integer keyword evaluation. */
61 init_data_kwint(int npar, gmx_ana_selparam_t *param);
62 /** Allocates data for real keyword evaluation. */
64 init_data_kwreal(int npar, gmx_ana_selparam_t *param);
65 /** Allocates data for string keyword evaluation. */
67 init_data_kwstr(int npar, gmx_ana_selparam_t *param);
68 /** Initializes data for integer keyword evaluation. */
70 init_kwint(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
71 /** Initializes data for real keyword evaluation. */
73 init_kwreal(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
74 /** Initializes data for string keyword evaluation. */
76 init_kwstr(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
77 /** Frees the memory allocated for string keyword evaluation. */
79 free_data_kwstr(void *data);
80 /** Evaluates integer selection keywords. */
82 evaluate_keyword_int(t_topology *top, t_trxframe *fr, t_pbc *pbc,
83 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
84 /** Evaluates real selection keywords. */
86 evaluate_keyword_real(t_topology *top, t_trxframe *fr, t_pbc *pbc,
87 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
88 /** Evaluates string selection keywords. */
90 evaluate_keyword_str(t_topology *top, t_trxframe *fr, t_pbc *pbc,
91 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
94 * Data structure for integer keyword expression evaluation.
96 typedef struct t_methoddata_kwint
98 /** Array of values for the keyword. */
100 /** Number of ranges in the \p r array. */
103 * Array of sorted integer ranges to match against.
105 * Each range is made of two integers, giving the endpoints (inclusive).
106 * This field stores the pointer to the ranges allocated by the
107 * parameter parser; see \ref SPAR_RANGES for more information.
110 } t_methoddata_kwint;
113 * Data structure for real keyword expression evaluation.
115 typedef struct t_methoddata_kwreal
117 /** Array of values for the keyword. */
119 /** Number of ranges in the \p r array. */
122 * Array of sorted ranges to match against.
124 * Each range is made of two values, giving the endpoints (inclusive).
125 * This field stores the pointer to the ranges allocated by the
126 * parameter parser; see \ref SPAR_RANGES for more information.
129 } t_methoddata_kwreal;
132 * Data structure for string keyword expression evaluation.
134 typedef struct t_methoddata_kwstr
136 /** Array of values for the keyword. */
138 /** Number of elements in the \p val array. */
141 * Array of strings/regular expressions to match against.
143 struct t_methoddata_kwstr_match {
144 /** TRUE if the expression is a regular expression, FALSE otherwise. */
146 /** The value to match against. */
149 /** Compiled regular expression if \p bRegExp is TRUE. */
152 /** The string if \p bRegExp is FALSE; */
156 } t_methoddata_kwstr;
158 /** Parameters for integer keyword evaluation. */
159 static gmx_ana_selparam_t smparams_keyword_int[] = {
160 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
161 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
164 /** Parameters for real keyword evaluation. */
165 static gmx_ana_selparam_t smparams_keyword_real[] = {
166 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL | SPAR_DYNAMIC},
167 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
170 /** Parameters for string keyword evaluation. */
171 static gmx_ana_selparam_t smparams_keyword_str[] = {
172 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
173 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_VARNUM},
176 /** \internal Selection method data for integer keyword evaluation. */
177 gmx_ana_selmethod_t sm_keyword_int = {
178 "kw_int", GROUP_VALUE, SMETH_SINGLEVAL,
179 asize(smparams_keyword_int), smparams_keyword_int,
186 &evaluate_keyword_int,
191 /** \internal Selection method data for real keyword evaluation. */
192 gmx_ana_selmethod_t sm_keyword_real = {
193 "kw_real", GROUP_VALUE, SMETH_SINGLEVAL,
194 asize(smparams_keyword_real), smparams_keyword_real,
201 &evaluate_keyword_real,
206 /** \internal Selection method data for string keyword evaluation. */
207 gmx_ana_selmethod_t sm_keyword_str = {
208 "kw_str", GROUP_VALUE, SMETH_SINGLEVAL,
209 asize(smparams_keyword_str), smparams_keyword_str,
216 &evaluate_keyword_str,
221 /** Initializes keyword evaluation for an arbitrary group. */
223 init_kweval(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
224 /** Initializes output for keyword evaluation in an arbitrary group. */
226 init_output_kweval(t_topology *top, gmx_ana_selvalue_t *out, void *data);
227 /** Frees the data allocated for keyword evaluation in an arbitrary group. */
229 free_data_kweval(void *data);
230 /** Initializes frame evaluation for keyword evaluation in an arbitrary group. */
232 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
233 /** Evaluates keywords in an arbitrary group. */
235 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc,
236 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
239 * Data structure for keyword evaluation in arbitrary groups.
243 /** Wrapped keyword method for evaluating the values. */
244 gmx_ana_selmethod_t *kwmethod;
245 /** Method data for \p kwmethod. */
247 /** Group in which \p kwmethod should be evaluated. */
249 } t_methoddata_kweval;
251 /** Parameters for keyword evaluation in an arbitrary group. */
252 static gmx_ana_selparam_t smparams_kweval[] = {
253 {NULL, {GROUP_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
257 /********************************************************************
258 * INTEGER KEYWORD EVALUATION
259 ********************************************************************/
262 * \param[in] npar Not used.
263 * \param param Not used.
264 * \returns Pointer to the allocated data (\ref t_methoddata_kwint).
266 * Allocates memory for a \ref t_methoddata_kwint structure.
269 init_data_kwint(int npar, gmx_ana_selparam_t *param)
271 t_methoddata_kwint *data;
278 * \param[in] top Not used.
279 * \param[in] npar Not used (should be 2).
280 * \param[in] param Method parameters (should point to \ref smparams_keyword_int).
281 * \param[in] data Should point to \ref t_methoddata_kwint.
282 * \returns 0 (the initialization always succeeds).
285 init_kwint(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
287 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
289 d->v = param[0].val.u.i;
290 d->n = param[1].val.nr;
291 d->r = param[1].val.u.i;
296 * See sel_updatefunc() for description of the parameters.
297 * \p data should point to a \c t_methoddata_kwint.
299 * Does a binary search to find which atoms match the ranges in the
300 * \c t_methoddata_kwint structure for this selection.
301 * Matching atoms are stored in \p out->u.g.
304 evaluate_keyword_int(t_topology *top, t_trxframe *fr, t_pbc *pbc,
305 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
307 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
308 int n, i, j, jmin, jmax;
313 for (i = 0; i < g->isize; ++i)
316 if (d->r[0] > val || d->r[2*n-1] < val)
322 while (jmax - jmin > 1)
324 j = jmin + (jmax - jmin) / 2;
332 if (val <= d->r[2*j+1])
339 if (val <= d->r[2*jmin+1])
341 out->u.g->index[out->u.g->isize++] = g->index[i];
348 /********************************************************************
349 * REAL KEYWORD EVALUATION
350 ********************************************************************/
353 * \param[in] npar Not used.
354 * \param param Not used.
355 * \returns Pointer to the allocated data (\ref t_methoddata_kwreal).
357 * Allocates memory for a \ref t_methoddata_kwreal structure.
360 init_data_kwreal(int npar, gmx_ana_selparam_t *param)
362 t_methoddata_kwreal *data;
369 * \param[in] top Not used.
370 * \param[in] npar Not used (should be 2).
371 * \param[in] param Method parameters (should point to \ref smparams_keyword_real).
372 * \param[in] data Should point to \ref t_methoddata_kwreal.
373 * \returns 0 (the initialization always succeeds).
376 init_kwreal(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
378 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
380 d->v = param[0].val.u.r;
381 d->n = param[1].val.nr;
382 d->r = param[1].val.u.r;
387 * See sel_updatefunc() for description of the parameters.
388 * \p data should point to a \c t_methoddata_kwreal.
390 * Does a binary search to find which atoms match the ranges in the
391 * \c t_methoddata_kwreal structure for this selection.
392 * Matching atoms are stored in \p out->u.g.
395 evaluate_keyword_real(t_topology *top, t_trxframe *fr, t_pbc *pbc,
396 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
398 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
399 int n, i, j, jmin, jmax;
404 for (i = 0; i < g->isize; ++i)
407 if (d->r[0] > val || d->r[2*n-1] < val)
413 while (jmax - jmin > 1)
415 j = jmin + (jmax - jmin) / 2;
423 if (val <= d->r[2*j+1])
430 if (val <= d->r[2*jmin+1])
432 out->u.g->index[out->u.g->isize++] = g->index[i];
439 /********************************************************************
440 * STRING KEYWORD EVALUATION
441 ********************************************************************/
444 * \param[in] npar Not used.
445 * \param param Not used.
446 * \returns Pointer to the allocated data (\ref t_methoddata_kwstr).
448 * Allocates memory for a \ref t_methoddata_kwstr structure.
451 init_data_kwstr(int npar, gmx_ana_selparam_t *param)
453 t_methoddata_kwstr *data;
460 * \param[in] top Not used.
461 * \param[in] npar Not used (should be 2).
462 * \param[in] param Method parameters (should point to \ref smparams_keyword_str).
463 * \param[in] data Should point to \ref t_methoddata_kwstr.
464 * \returns 0 (the initialization always succeeds).
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] != '*')
499 snew(buf, strlen(s) + 3);
500 sprintf(buf, "^%s$", s);
501 if (regcomp(&d->m[i].u.r, buf, REG_EXTENDED | REG_NOSUB))
504 fprintf(stderr, "WARNING: error in regular expression,\n"
505 " will match '%s' as a simple string\n", s);
510 fprintf(stderr, "WARNING: no regular expressions support,\n"
511 " will match '%s' as a simple string\n", s);
518 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];
599 /********************************************************************
600 * KEYWORD EVALUATION FOR ARBITRARY GROUPS
601 ********************************************************************/
604 * \param[in] top Not used.
605 * \param[in] npar Not used.
606 * \param[in] param Not used.
607 * \param[in] data Should point to \ref t_methoddata_kweval.
608 * \returns 0 on success, a non-zero error code on return.
610 * Calls the initialization method of the wrapped keyword.
613 init_kweval(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
615 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
617 return d->kwmethod->init(top, 0, NULL, d->kwmdata);
621 * \param[in] top Not used.
622 * \param[in,out] out Pointer to output data structure.
623 * \param[in,out] data Should point to \c t_methoddata_kweval.
624 * \returns 0 for success.
627 init_output_kweval(t_topology *top, gmx_ana_selvalue_t *out, void *data)
629 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
631 out->nr = d->g.isize;
636 * \param data Data to free (should point to a \c t_methoddata_kweval).
638 * Frees the memory allocated for all the members of \c t_methoddata_kweval.
641 free_data_kweval(void *data)
643 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
645 _gmx_selelem_free_method(d->kwmethod, d->kwmdata);
649 * \param[in] top Topology.
650 * \param[in] fr Current frame.
651 * \param[in] pbc PBC structure.
652 * \param data Should point to a \ref t_methoddata_kweval.
653 * \returns 0 on success, a non-zero error code on error.
655 * Creates a lookup structure that enables fast queries of whether a point
656 * is within the solid angle or not.
659 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
661 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
663 return d->kwmethod->init_frame(top, fr, pbc, d->kwmdata);
667 * See sel_updatefunc() for description of the parameters.
668 * \p data should point to a \c t_methoddata_kweval.
670 * Calls the evaluation function of the wrapped keyword with the given
671 * parameters, with the exception of using \c t_methoddata_kweval::g for the
675 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc,
676 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
678 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
680 return d->kwmethod->update(top, fr, pbc, &d->g, out, d->kwmdata);
684 * \param[out] selp Pointer to receive a pointer to the created selection
685 * element (set to NULL on error).
686 * \param[in] method Keyword selection method to evaluate.
687 * \param[in] param Parameter that gives the group to evaluate \p method in.
688 * \param[in] scanner Scanner data structure.
689 * \returns 0 on success, non-zero error code on error.
691 * Creates a \ref SEL_EXPRESSION selection element (pointer put in \c *selp)
692 * that evaluates the keyword method given by \p method in the group given by
696 _gmx_sel_init_keyword_evaluator(t_selelem **selp, gmx_ana_selmethod_t *method,
697 t_selexpr_param *param, void *scanner)
700 t_methoddata_kweval *data;
702 if ((method->flags & (SMETH_SINGLEVAL | SMETH_VARNUMVAL))
703 || method->outinit || method->pupdate)
705 _gmx_selexpr_free_params(param);
706 gmx_incons("unsupported keyword method for arbitrary group evaluation");
711 sel = _gmx_selelem_create(SEL_EXPRESSION);
712 _gmx_selelem_set_method(sel, method, scanner);
715 data->kwmethod = sel->u.expr.method;
716 data->kwmdata = sel->u.expr.mdata;
717 gmx_ana_index_clear(&data->g);
719 snew(sel->u.expr.method, 1);
720 memcpy(sel->u.expr.method, data->kwmethod, sizeof(gmx_ana_selmethod_t));
721 sel->u.expr.method->flags |= SMETH_VARNUMVAL;
722 sel->u.expr.method->init_data = NULL;
723 sel->u.expr.method->set_poscoll = NULL;
724 sel->u.expr.method->init = method->init ? &init_kweval : NULL;
725 sel->u.expr.method->outinit = &init_output_kweval;
726 sel->u.expr.method->free = &free_data_kweval;
727 sel->u.expr.method->init_frame = method->init_frame ? &init_frame_kweval : NULL;
728 sel->u.expr.method->update = &evaluate_kweval;
729 sel->u.expr.method->pupdate = NULL;
730 sel->u.expr.method->nparams = asize(smparams_kweval);
731 sel->u.expr.method->param = smparams_kweval;
732 _gmx_selelem_init_method_params(sel, scanner);
733 sel->u.expr.mdata = data;
735 sel->u.expr.method->param[0].val.u.g = &data->g;
739 if (!_gmx_sel_parse_params(param, sel->u.expr.method->nparams,
740 sel->u.expr.method->param, sel, scanner))
742 _gmx_selelem_free(sel);