2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2009, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief Implementations of internal selection methods for numeric and
40 * string keyword evaluation.
47 #ifdef HAVE_SYS_TYPES_H
48 #include <sys/types.h> /*old Mac needs types before regex.h*/
55 #include <gmx_fatal.h>
60 #include <selmethod.h>
63 #include "parsetree.h"
66 /** Allocates data for integer keyword evaluation. */
68 init_data_kwint(int npar, gmx_ana_selparam_t *param);
69 /** Allocates data for real keyword evaluation. */
71 init_data_kwreal(int npar, gmx_ana_selparam_t *param);
72 /** Allocates data for string keyword evaluation. */
74 init_data_kwstr(int npar, gmx_ana_selparam_t *param);
75 /** Initializes data for integer keyword evaluation. */
77 init_kwint(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
78 /** Initializes data for real keyword evaluation. */
80 init_kwreal(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
81 /** Initializes data for string keyword evaluation. */
83 init_kwstr(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
84 /** Frees the memory allocated for string keyword evaluation. */
86 free_data_kwstr(void *data);
87 /** Evaluates integer selection keywords. */
89 evaluate_keyword_int(t_topology *top, t_trxframe *fr, t_pbc *pbc,
90 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
91 /** Evaluates real selection keywords. */
93 evaluate_keyword_real(t_topology *top, t_trxframe *fr, t_pbc *pbc,
94 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
95 /** Evaluates string selection keywords. */
97 evaluate_keyword_str(t_topology *top, t_trxframe *fr, t_pbc *pbc,
98 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
101 * Data structure for integer keyword expression evaluation.
103 typedef struct t_methoddata_kwint
105 /** Array of values for the keyword. */
107 /** Number of ranges in the \p r array. */
110 * Array of sorted integer ranges to match against.
112 * Each range is made of two integers, giving the endpoints (inclusive).
113 * This field stores the pointer to the ranges allocated by the
114 * parameter parser; see \ref SPAR_RANGES for more information.
117 } t_methoddata_kwint;
120 * Data structure for real keyword expression evaluation.
122 typedef struct t_methoddata_kwreal
124 /** Array of values for the keyword. */
126 /** Number of ranges in the \p r array. */
129 * Array of sorted ranges to match against.
131 * Each range is made of two values, giving the endpoints (inclusive).
132 * This field stores the pointer to the ranges allocated by the
133 * parameter parser; see \ref SPAR_RANGES for more information.
136 } t_methoddata_kwreal;
139 * Data structure for string keyword expression evaluation.
141 typedef struct t_methoddata_kwstr
143 /** Array of values for the keyword. */
145 /** Number of elements in the \p val array. */
148 * Array of strings/regular expressions to match against.
150 struct t_methoddata_kwstr_match {
151 /** TRUE if the expression is a regular expression, FALSE otherwise. */
153 /** The value to match against. */
156 /** Compiled regular expression if \p bRegExp is TRUE. */
159 /** The string if \p bRegExp is FALSE; */
163 } t_methoddata_kwstr;
165 /** Parameters for integer keyword evaluation. */
166 static gmx_ana_selparam_t smparams_keyword_int[] = {
167 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
168 {NULL, {INT_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
171 /** Parameters for real keyword evaluation. */
172 static gmx_ana_selparam_t smparams_keyword_real[] = {
173 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL | SPAR_DYNAMIC},
174 {NULL, {REAL_VALUE, -1, {NULL}}, NULL, SPAR_RANGES | SPAR_VARNUM},
177 /** Parameters for string keyword evaluation. */
178 static gmx_ana_selparam_t smparams_keyword_str[] = {
179 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_ATOMVAL},
180 {NULL, {STR_VALUE, -1, {NULL}}, NULL, SPAR_VARNUM},
183 /** \internal Selection method data for integer keyword evaluation. */
184 gmx_ana_selmethod_t sm_keyword_int = {
185 "kw_int", GROUP_VALUE, SMETH_SINGLEVAL,
186 asize(smparams_keyword_int), smparams_keyword_int,
193 &evaluate_keyword_int,
198 /** \internal Selection method data for real keyword evaluation. */
199 gmx_ana_selmethod_t sm_keyword_real = {
200 "kw_real", GROUP_VALUE, SMETH_SINGLEVAL,
201 asize(smparams_keyword_real), smparams_keyword_real,
208 &evaluate_keyword_real,
213 /** \internal Selection method data for string keyword evaluation. */
214 gmx_ana_selmethod_t sm_keyword_str = {
215 "kw_str", GROUP_VALUE, SMETH_SINGLEVAL,
216 asize(smparams_keyword_str), smparams_keyword_str,
223 &evaluate_keyword_str,
228 /** Initializes keyword evaluation for an arbitrary group. */
230 init_kweval(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
231 /** Initializes output for keyword evaluation in an arbitrary group. */
233 init_output_kweval(t_topology *top, gmx_ana_selvalue_t *out, void *data);
234 /** Frees the data allocated for keyword evaluation in an arbitrary group. */
236 free_data_kweval(void *data);
237 /** Initializes frame evaluation for keyword evaluation in an arbitrary group. */
239 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
240 /** Evaluates keywords in an arbitrary group. */
242 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc,
243 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data);
246 * Data structure for keyword evaluation in arbitrary groups.
250 /** Wrapped keyword method for evaluating the values. */
251 gmx_ana_selmethod_t *kwmethod;
252 /** Method data for \p kwmethod. */
254 /** Group in which \p kwmethod should be evaluated. */
256 } t_methoddata_kweval;
258 /** Parameters for keyword evaluation in an arbitrary group. */
259 static gmx_ana_selparam_t smparams_kweval[] = {
260 {NULL, {GROUP_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
264 /********************************************************************
265 * INTEGER KEYWORD EVALUATION
266 ********************************************************************/
269 * \param[in] npar Not used.
270 * \param param Not used.
271 * \returns Pointer to the allocated data (\ref t_methoddata_kwint).
273 * Allocates memory for a \ref t_methoddata_kwint structure.
276 init_data_kwint(int npar, gmx_ana_selparam_t *param)
278 t_methoddata_kwint *data;
285 * \param[in] top Not used.
286 * \param[in] npar Not used (should be 2).
287 * \param[in] param Method parameters (should point to \ref smparams_keyword_int).
288 * \param[in] data Should point to \ref t_methoddata_kwint.
289 * \returns 0 (the initialization always succeeds).
292 init_kwint(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
294 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
296 d->v = param[0].val.u.i;
297 d->n = param[1].val.nr;
298 d->r = param[1].val.u.i;
303 * See sel_updatefunc() for description of the parameters.
304 * \p data should point to a \c t_methoddata_kwint.
306 * Does a binary search to find which atoms match the ranges in the
307 * \c t_methoddata_kwint structure for this selection.
308 * Matching atoms are stored in \p out->u.g.
311 evaluate_keyword_int(t_topology *top, t_trxframe *fr, t_pbc *pbc,
312 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
314 t_methoddata_kwint *d = (t_methoddata_kwint *)data;
315 int n, i, j, jmin, jmax;
320 for (i = 0; i < g->isize; ++i)
323 if (d->r[0] > val || d->r[2*n-1] < val)
329 while (jmax - jmin > 1)
331 j = jmin + (jmax - jmin) / 2;
339 if (val <= d->r[2*j+1])
346 if (val <= d->r[2*jmin+1])
348 out->u.g->index[out->u.g->isize++] = g->index[i];
355 /********************************************************************
356 * REAL KEYWORD EVALUATION
357 ********************************************************************/
360 * \param[in] npar Not used.
361 * \param param Not used.
362 * \returns Pointer to the allocated data (\ref t_methoddata_kwreal).
364 * Allocates memory for a \ref t_methoddata_kwreal structure.
367 init_data_kwreal(int npar, gmx_ana_selparam_t *param)
369 t_methoddata_kwreal *data;
376 * \param[in] top Not used.
377 * \param[in] npar Not used (should be 2).
378 * \param[in] param Method parameters (should point to \ref smparams_keyword_real).
379 * \param[in] data Should point to \ref t_methoddata_kwreal.
380 * \returns 0 (the initialization always succeeds).
383 init_kwreal(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
385 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
387 d->v = param[0].val.u.r;
388 d->n = param[1].val.nr;
389 d->r = param[1].val.u.r;
394 * See sel_updatefunc() for description of the parameters.
395 * \p data should point to a \c t_methoddata_kwreal.
397 * Does a binary search to find which atoms match the ranges in the
398 * \c t_methoddata_kwreal structure for this selection.
399 * Matching atoms are stored in \p out->u.g.
402 evaluate_keyword_real(t_topology *top, t_trxframe *fr, t_pbc *pbc,
403 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
405 t_methoddata_kwreal *d = (t_methoddata_kwreal *)data;
406 int n, i, j, jmin, jmax;
411 for (i = 0; i < g->isize; ++i)
414 if (d->r[0] > val || d->r[2*n-1] < val)
420 while (jmax - jmin > 1)
422 j = jmin + (jmax - jmin) / 2;
430 if (val <= d->r[2*j+1])
437 if (val <= d->r[2*jmin+1])
439 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;
467 * \param[in] top Not used.
468 * \param[in] npar Not used (should be 2).
469 * \param[in] param Method parameters (should point to \ref smparams_keyword_str).
470 * \param[in] data Should point to \ref t_methoddata_kwstr.
471 * \returns 0 (the initialization always succeeds).
474 init_kwstr(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
476 t_methoddata_kwstr *d = (t_methoddata_kwstr *)data;
483 d->v = param[0].val.u.s;
484 d->n = param[1].val.nr;
485 /* Return if this is not the first time */
491 for (i = 0; i < d->n; ++i)
493 s = param[1].val.u.s[i];
495 for (j = 0; j < strlen(s); ++j)
497 if (ispunct(s[j]) && s[j] != '?' && s[j] != '*')
506 snew(buf, strlen(s) + 3);
507 sprintf(buf, "^%s$", s);
508 if (regcomp(&d->m[i].u.r, buf, REG_EXTENDED | REG_NOSUB))
511 fprintf(stderr, "WARNING: error in regular expression,\n"
512 " will match '%s' as a simple string\n", s);
517 fprintf(stderr, "WARNING: no regular expressions support,\n"
518 " will match '%s' as a simple string\n", s);
525 d->m[i].bRegExp = bRegExp;
531 * \param data Data to free (should point to a \ref t_methoddata_kwstr).
533 * Frees the memory allocated for t_methoddata_kwstr::val.
536 free_data_kwstr(void *data)
538 t_methoddata_kwstr *d = (t_methoddata_kwstr *)data;
541 for (i = 0; i < d->n; ++i)
546 /* This branch should only be taken if regular expressions
547 * are available, but the ifdef is still needed. */
548 regfree(&d->m[i].u.r);
556 * See sel_updatefunc() for description of the parameters.
557 * \p data should point to a \c t_methoddata_kwstr.
559 * Does a linear search to find which atoms match the strings in the
560 * \c t_methoddata_kwstr structure for this selection.
561 * Wildcards are allowed in the strings.
562 * Matching atoms are stored in \p out->u.g.
565 evaluate_keyword_str(t_topology *top, t_trxframe *fr, t_pbc *pbc,
566 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
568 t_methoddata_kwstr *d = (t_methoddata_kwstr *)data;
573 for (i = 0; i < g->isize; ++i)
576 for (j = 0; j < d->n && !bFound; ++j)
581 /* This branch should only be taken if regular expressions
582 * are available, but the ifdef is still needed. */
583 if (!regexec(&d->m[j].u.r, d->v[i], 0, NULL, 0))
591 if (gmx_wcmatch(d->m[j].u.s, d->v[i]) == 0)
599 out->u.g->index[out->u.g->isize++] = g->index[i];
606 /********************************************************************
607 * KEYWORD EVALUATION FOR ARBITRARY GROUPS
608 ********************************************************************/
611 * \param[in] top Not used.
612 * \param[in] npar Not used.
613 * \param[in] param Not used.
614 * \param[in] data Should point to \ref t_methoddata_kweval.
615 * \returns 0 on success, a non-zero error code on return.
617 * Calls the initialization method of the wrapped keyword.
620 init_kweval(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
622 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
624 return d->kwmethod->init(top, 0, NULL, d->kwmdata);
628 * \param[in] top Not used.
629 * \param[in,out] out Pointer to output data structure.
630 * \param[in,out] data Should point to \c t_methoddata_kweval.
631 * \returns 0 for success.
634 init_output_kweval(t_topology *top, gmx_ana_selvalue_t *out, void *data)
636 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
638 out->nr = d->g.isize;
639 _gmx_selvalue_reserve(out, out->nr);
644 * \param data Data to free (should point to a \c t_methoddata_kweval).
646 * Frees the memory allocated for all the members of \c t_methoddata_kweval.
649 free_data_kweval(void *data)
651 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
653 _gmx_selelem_free_method(d->kwmethod, d->kwmdata);
657 * \param[in] top Topology.
658 * \param[in] fr Current frame.
659 * \param[in] pbc PBC structure.
660 * \param data Should point to a \ref t_methoddata_kweval.
661 * \returns 0 on success, a non-zero error code on error.
663 * Creates a lookup structure that enables fast queries of whether a point
664 * is within the solid angle or not.
667 init_frame_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
669 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
671 return d->kwmethod->init_frame(top, fr, pbc, d->kwmdata);
675 * See sel_updatefunc() for description of the parameters.
676 * \p data should point to a \c t_methoddata_kweval.
678 * Calls the evaluation function of the wrapped keyword with the given
679 * parameters, with the exception of using \c t_methoddata_kweval::g for the
683 evaluate_kweval(t_topology *top, t_trxframe *fr, t_pbc *pbc,
684 gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
686 t_methoddata_kweval *d = (t_methoddata_kweval *)data;
688 return d->kwmethod->update(top, fr, pbc, &d->g, out, d->kwmdata);
692 * \param[out] selp Pointer to receive a pointer to the created selection
693 * element (set to NULL on error).
694 * \param[in] method Keyword selection method to evaluate.
695 * \param[in] param Parameter that gives the group to evaluate \p method in.
696 * \param[in] scanner Scanner data structure.
697 * \returns 0 on success, non-zero error code on error.
699 * Creates a \ref SEL_EXPRESSION selection element (pointer put in \c *selp)
700 * that evaluates the keyword method given by \p method in the group given by
704 _gmx_sel_init_keyword_evaluator(t_selelem **selp, gmx_ana_selmethod_t *method,
705 t_selexpr_param *param, void *scanner)
708 t_methoddata_kweval *data;
710 if ((method->flags & (SMETH_SINGLEVAL | SMETH_VARNUMVAL))
711 || method->outinit || method->pupdate)
713 _gmx_selexpr_free_params(param);
714 gmx_incons("unsupported keyword method for arbitrary group evaluation");
719 sel = _gmx_selelem_create(SEL_EXPRESSION);
720 _gmx_selelem_set_method(sel, method, scanner);
723 data->kwmethod = sel->u.expr.method;
724 data->kwmdata = sel->u.expr.mdata;
725 gmx_ana_index_clear(&data->g);
727 snew(sel->u.expr.method, 1);
728 memcpy(sel->u.expr.method, data->kwmethod, sizeof(gmx_ana_selmethod_t));
729 sel->u.expr.method->flags |= SMETH_VARNUMVAL;
730 sel->u.expr.method->init_data = NULL;
731 sel->u.expr.method->set_poscoll = NULL;
732 sel->u.expr.method->init = method->init ? &init_kweval : NULL;
733 sel->u.expr.method->outinit = &init_output_kweval;
734 sel->u.expr.method->free = &free_data_kweval;
735 sel->u.expr.method->init_frame = method->init_frame ? &init_frame_kweval : NULL;
736 sel->u.expr.method->update = &evaluate_kweval;
737 sel->u.expr.method->pupdate = NULL;
738 sel->u.expr.method->nparams = asize(smparams_kweval);
739 sel->u.expr.method->param = smparams_kweval;
740 _gmx_selelem_init_method_params(sel, scanner);
741 sel->u.expr.mdata = data;
743 sel->u.expr.method->param[0].val.u.g = &data->g;
747 if (!_gmx_sel_parse_params(param, sel->u.expr.method->nparams,
748 sel->u.expr.method->param, sel, scanner))
750 _gmx_selelem_free(sel);