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 Implementation of functions of in parsetree.h.
35 * \page selparser Selection parsing
37 * The selection parser is implemented in the following files:
39 * Tokenizer implemented using Flex, splits the input into tokens
40 * (scanner.c and scanner_flex.h are generated from this file).
41 * - scanner.h, scanner_internal.h, scanner_internal.c:
42 * Helper functions for scanner.l and for interfacing between
43 * scanner.l and parser.y. Functions in scanner_internal.h are only
44 * used from scanner.l, while scanner.h is used from the parser.
45 * - symrec.h, symrec.c:
46 * Functions used by the tokenizer to handle the symbol table, i.e.,
47 * the recognized keywords. Some basic keywords are hardcoded into
48 * scanner.l, but all method and variable references go through the
49 * symbol table, as do position evaluation keywords.
51 * Semantic rules for parsing the grammar
52 * (parser.c and parser.h are generated from this file by Bison).
53 * - parsetree.h, parsetree.c:
54 * Functions called from actions in parser.y to construct the
55 * evaluation elements corresponding to different grammar elements.
56 * parsetree.c also defines the external interface of the parser,
57 * i.e., the \c gmx_ana_selcollection_parse_*() functions declared
60 * Defines a function that processes the parameters of selection
61 * methods and initializes the children of the method element.
63 * The basic control flow in the parser is as follows: when a
64 * fmx_ana_selcollection_parse_*() gets called, it performs some
65 * initialization, and then calls the _gmx_sel_yyparse() function generated
66 * by Bison. This function then calls _gmx_sel_yylex() to repeatedly read
67 * tokens from the input (more complex tasks related to token recognition
68 * and bookkeeping are done by functions in scanner_internal.c) and uses the
69 * grammar rules to decide what to do with them. Whenever a grammar rule
70 * matches, a corresponding function in parsetree.c is called to construct
71 * either a temporary representation for the object or a \c t_selelem object
72 * (some simple rules are handled internally in parser.y).
73 * When a complete selection has been parsed, the functions in parsetree.c
74 * also take care of updating the \c gmx_ana_selcollection_t structure
77 * The rest of this page describes the resulting \c t_selelem object tree.
78 * Before the selections can be evaluated, this tree needs to be passed to
79 * the selection compiler, which is described on a separate page:
83 * \section selparser_tree Element tree constructed by the parser
85 * The parser initializes the following fields in all selection elements:
86 * \c t_selelem::name, \c t_selelem::type, \c t_selelem::v\c .type,
87 * \c t_selelem::flags, \c t_selelem::child, \c t_selelem::next, and
88 * \c t_selelem::refcount.
89 * Some other fields are also initialized for particular element types as
91 * Fields that are not initialized are set to zero, NULL, or other similar
95 * \subsection selparser_tree_root Root elements
97 * The parser creates a \ref SEL_ROOT selection element for each variable
98 * assignment and each selection. However, there are two exceptions that do
99 * not result in a \ref SEL_ROOT element (in these cases, only the symbol
100 * table is modified):
101 * - Variable assignments that assign a variable to another variable.
102 * - Variable assignments that assign a non-group constant.
104 * The \ref SEL_ROOT elements are linked together in a chain in the same order
107 * The children of the \ref SEL_ROOT elements can be used to distinguish
108 * the two types of root elements from each other:
109 * - For variable assignments, the first and only child is always
110 * a \ref SEL_SUBEXPR element.
111 * - For selections, the first child is a \ref SEL_EXPRESSION or a
112 * \ref SEL_MODIFIER element that evaluates the final positions (if the
113 * selection defines a constant position, the child is a \ref SEL_CONST).
114 * The rest of the children are \ref SEL_MODIFIER elements with
115 * \ref NO_VALUE, in the order given by the user.
117 * The name of the selection/variable is stored in \c t_selelem::cgrp\c .name.
118 * It is set to either the name provided by the user or the selection string
119 * for selections not explicitly named by the user.
120 * \ref SEL_ROOT or \ref SEL_SUBEXPR elements do not appear anywhere else.
123 * \subsection selparser_tree_const Constant elements
125 * \ref SEL_CONST elements are created for every constant that is required
126 * for later evaluation.
127 * Currently, \ref SEL_CONST elements can be present for
128 * - selections that consist of a constant position,
129 * - \ref GROUP_VALUE method parameters if provided using external index
132 * For group-valued elements, the value is stored in \c t_selelem::cgrp;
133 * other types of values are stored in \c t_selelem::v.
134 * Constants that appear as parameters for selection methods are not present
135 * in the selection tree unless they have \ref GROUP_VALUE.
136 * \ref SEL_CONST elements have no children.
139 * \subsection selparser_tree_method Method evaluation elements
141 * \ref SEL_EXPRESSION and \ref SEL_MODIFIER elements are treated very
142 * similarly. The \c gmx_ana_selmethod_t structure corresponding to the
143 * evaluation method is in \c t_selelem::method, and the method data in
144 * \c t_selelem::mdata has been allocated using sel_datafunc().
145 * If a non-standard reference position type was set, \c t_selelem::pc has
146 * also been created, but only the type has been set.
147 * All children of these elements are of the type \ref SEL_SUBEXPRREF, and
148 * each describes a selection that needs to be evaluated to obtain a value
149 * for one parameter of the method.
150 * No children are present for parameters that were given a constant
151 * non-\ref GROUP_VALUE value.
152 * The children are sorted in the order in which the parameters appear in the
153 * \ref gmx_ana_selmethod_t structure.
155 * In addition to actual selection keywords, \ref SEL_EXPRESSION elements
156 * are used internally to implement numerical comparisons (e.g., "x < 5")
157 * and keyword matching (e.g., "resnr 1 to 3" or "name CA").
160 * \subsection selparser_tree_subexpr Subexpression elements
162 * \ref SEL_SUBEXPR elements only appear for variables, as described above.
163 * \c t_selelem::name points to the name of the variable (from the
164 * \ref SEL_ROOT element).
165 * The element always has exactly one child, which represents the value of
167 * \ref SEL_SUBEXPR element is the only element type that can have
168 * \c t_selelem::refcount different from 1.
170 * \ref SEL_SUBEXPRREF elements are used for two purposes:
171 * - Variable references that need to be evaluated (i.e., there is a
172 * \ref SEL_SUBEXPR element for the variable) are represented using
173 * \ref SEL_SUBEXPRREF elements.
174 * In this case, \c t_selelem::param is NULL, and the first and only
175 * child of the element is the \ref SEL_SUBEXPR element of the variable.
176 * Such references can appear anywhere where the variable value
177 * (the child of the \ref SEL_SUBEXPR element) would be valid.
178 * - Children of \ref SEL_EXPRESSION and \ref SEL_MODIFIER elements are
179 * always of this type. For these elements, \c t_selelem::param is
180 * initialized to point to the parameter that receives the value from
182 * Each such element has exactly one child, which can be of any type;
183 * the \ref SEL_SUBEXPR element of a variable is used if the value comes
184 * from a variable, otherwise the child type is not \ref SEL_SUBEXPR.
187 * \subsection selparser_tree_bool Boolean elements
189 * One \ref SEL_BOOLEAN element is created for each boolean keyword in the
190 * input, and the tree structure represents the evaluation order.
191 * The \c t_selelem::boolt type gives the type of the operation.
192 * Each element has exactly two children (one for \ref BOOL_NOT elements),
193 * which are in the order given in the input.
194 * The children always have \ref GROUP_VALUE, but different element types
198 * \subsection selparser_tree_arith Arithmetic elements
200 * One \ref SEL_ARITHMETIC element is created for each arithmetic operation in
201 * the input, and the tree structure represents the evaluation order.
202 * The \c t_selelem::optype type gives the name of the operation.
203 * Each element has exactly two children (one for unary negation elements),
204 * which are in the order given in the input.
216 #include <gmx_fatal.h>
219 #include <selection.h>
220 #include <selmethod.h>
222 #include "keywords.h"
223 #include "parsetree.h"
224 #include "selcollection.h"
233 * Parser function generated by Bison.
236 _gmx_sel_yybparse(void *scanner
);
239 * It is a simple wrapper for fprintf(stderr, ...).
242 _gmx_selparser_error(const char *fmt
, ...)
246 fprintf(stderr
, "selection parser: ");
247 vfprintf(stderr
, fmt
, ap
);
248 fprintf(stderr
, "\n");
253 * \param[in] type Type for the new value.
254 * \returns Pointer to the newly allocated value.
257 _gmx_selexpr_create_value(e_selvalue_t type
)
259 t_selexpr_value
*value
;
262 value
->bExpr
= FALSE
;
268 * \param[in] expr Expression for the value.
269 * \returns Pointer to the newly allocated value.
272 _gmx_selexpr_create_value_expr(t_selelem
*expr
)
274 t_selexpr_value
*value
;
276 value
->type
= expr
->v
.type
;
278 value
->u
.expr
= expr
;
284 * \param[in] name Name for the new parameter.
285 * \returns Pointer to the newly allocated parameter.
287 * No copy of \p name is made.
290 _gmx_selexpr_create_param(char *name
)
292 t_selexpr_param
*param
;
300 * \param value Pointer to the beginning of the value list to free.
302 * The expressions referenced by the values are also freed
303 * (to prevent this, set the expression to NULL before calling the function).
306 _gmx_selexpr_free_values(t_selexpr_value
*value
)
308 t_selexpr_value
*old
;
316 _gmx_selelem_free(value
->u
.expr
);
319 else if (value
->type
== STR_VALUE
)
330 * \param param Pointer the the beginning of the parameter list to free.
332 * The values of the parameters are freed with free_selexpr_values().
335 _gmx_selexpr_free_params(t_selexpr_param
*param
)
337 t_selexpr_param
*old
;
341 _gmx_selexpr_free_values(param
->value
);
350 * \param[in,out] sel Root of the selection element tree to initialize.
351 * \returns 0 on success, an error code on error.
353 * Propagates the \ref SEL_DYNAMIC flag from the children of \p sel to \p sel
354 * (if any child of \p sel is dynamic, \p sel is also marked as such).
355 * The \ref SEL_DYNAMIC flag is also set for \ref SEL_EXPRESSION elements with
357 * Also, sets one of the \ref SEL_SINGLEVAL, \ref SEL_ATOMVAL, or
358 * \ref SEL_VARNUMVAL flags, either based on the children or on the type of
359 * the selection method.
360 * If the types of the children conflict, an error is returned.
362 * The flags of the children of \p sel are also updated if not done earlier.
363 * The flags are initialized only once for any element; if \ref SEL_FLAGSSET
364 * is set for an element, the function returns immediately, and the recursive
365 * operation does not descend beyond such elements.
368 _gmx_selelem_update_flags(t_selelem
*sel
)
372 bool bUseChildType
=FALSE
;
373 bool bOnlySingleChildren
;
375 /* Return if the flags have already been set */
376 if (sel
->flags
& SEL_FLAGSSET
)
380 /* Set the flags based on the current element type */
384 sel
->flags
|= SEL_SINGLEVAL
;
385 bUseChildType
= FALSE
;
389 if (sel
->u
.expr
.method
->flags
& SMETH_DYNAMIC
)
391 sel
->flags
|= SEL_DYNAMIC
;
393 if (sel
->u
.expr
.method
->flags
& SMETH_SINGLEVAL
)
395 sel
->flags
|= SEL_SINGLEVAL
;
397 else if (sel
->u
.expr
.method
->flags
& SMETH_VARNUMVAL
)
399 sel
->flags
|= SEL_VARNUMVAL
;
403 sel
->flags
|= SEL_ATOMVAL
;
405 bUseChildType
= FALSE
;
409 sel
->flags
|= SEL_ATOMVAL
;
410 bUseChildType
= FALSE
;
414 if (sel
->v
.type
!= NO_VALUE
)
416 sel
->flags
|= SEL_VARNUMVAL
;
418 bUseChildType
= FALSE
;
422 bUseChildType
= FALSE
;
428 bUseChildType
= TRUE
;
431 /* Loop through children to propagate their flags upwards */
432 bOnlySingleChildren
= TRUE
;
436 /* Update the child */
437 rc
= _gmx_selelem_update_flags(child
);
442 /* Propagate the dynamic flag */
443 sel
->flags
|= (child
->flags
& SEL_DYNAMIC
);
444 /* Propagate the type flag if necessary and check for problems */
447 if ((sel
->flags
& SEL_VALTYPEMASK
)
448 && !(sel
->flags
& child
->flags
& SEL_VALTYPEMASK
))
450 _gmx_selparser_error("invalid combination of selection expressions");
453 sel
->flags
|= (child
->flags
& SEL_VALTYPEMASK
);
455 if (!(child
->flags
& SEL_SINGLEVAL
))
457 bOnlySingleChildren
= FALSE
;
462 /* For arithmetic expressions consisting only of single values,
463 * the result is also a single value. */
464 if (sel
->type
== SEL_ARITHMETIC
&& bOnlySingleChildren
)
466 sel
->flags
= (sel
->flags
& ~SEL_VALTYPEMASK
) | SEL_SINGLEVAL
;
468 /* For root elements, the type should be propagated here, after the
469 * children have been updated. */
470 if (sel
->type
== SEL_ROOT
)
472 sel
->flags
|= (sel
->child
->flags
& SEL_VALTYPEMASK
);
474 /* Mark that the flags are set */
475 sel
->flags
|= SEL_FLAGSSET
;
480 * \param[in,out] sel Selection element to initialize.
481 * \param[in] scanner Scanner data structure.
483 * A deep copy of the parameters is made to allow several
484 * expressions with the same method to coexist peacefully.
485 * Calls sel_datafunc() if one is specified for the method.
488 _gmx_selelem_init_method_params(t_selelem
*sel
, yyscan_t scanner
)
491 gmx_ana_selparam_t
*orgparam
;
492 gmx_ana_selparam_t
*param
;
496 nparams
= sel
->u
.expr
.method
->nparams
;
497 orgparam
= sel
->u
.expr
.method
->param
;
498 snew(param
, nparams
);
499 memcpy(param
, orgparam
, nparams
*sizeof(gmx_ana_selparam_t
));
500 for (i
= 0; i
< nparams
; ++i
)
502 param
[i
].flags
&= ~SPAR_SET
;
503 _gmx_selvalue_setstore(¶m
[i
].val
, NULL
);
504 if (param
[i
].flags
& SPAR_VARNUM
)
506 param
[i
].val
.nr
= -1;
508 /* Duplicate the enum value array if it is given statically */
509 if ((param
[i
].flags
& SPAR_ENUMVAL
) && orgparam
[i
].val
.u
.ptr
!= NULL
)
513 /* Count the values */
515 while (orgparam
[i
].val
.u
.s
[n
] != NULL
)
519 _gmx_selvalue_reserve(¶m
[i
].val
, n
+1);
520 memcpy(param
[i
].val
.u
.s
, orgparam
[i
].val
.u
.s
,
521 (n
+1)*sizeof(param
[i
].val
.u
.s
[0]));
525 if (sel
->u
.expr
.method
->init_data
)
527 mdata
= sel
->u
.expr
.method
->init_data(nparams
, param
);
530 gmx_fatal(FARGS
, "Method data initialization failed");
533 if (sel
->u
.expr
.method
->set_poscoll
)
535 gmx_ana_selcollection_t
*sc
= _gmx_sel_lexer_selcollection(scanner
);
537 sel
->u
.expr
.method
->set_poscoll(sc
->pcc
, mdata
);
539 /* Store the values */
540 sel
->u
.expr
.method
->param
= param
;
541 sel
->u
.expr
.mdata
= mdata
;
545 * \param[in,out] sel Selection element to initialize.
546 * \param[in] method Selection method to set.
547 * \param[in] scanner Scanner data structure.
549 * Makes a copy of \p method and stores it in \p sel->u.expr.method,
550 * and calls _gmx_selelem_init_method_params();
553 _gmx_selelem_set_method(t_selelem
*sel
, gmx_ana_selmethod_t
*method
,
558 _gmx_selelem_set_vtype(sel
, method
->type
);
559 sel
->name
= method
->name
;
560 snew(sel
->u
.expr
.method
, 1);
561 memcpy(sel
->u
.expr
.method
, method
, sizeof(gmx_ana_selmethod_t
));
562 _gmx_selelem_init_method_params(sel
, scanner
);
566 * Initializes the reference position calculation for a \ref SEL_EXPRESSION
569 * \param[in,out] pcc Position calculation collection to use.
570 * \param[in,out] sel Selection element to initialize.
571 * \param[in] rpost Reference position type to use (NULL = default).
572 * \returns 0 on success, a non-zero error code on error.
575 set_refpos_type(gmx_ana_poscalc_coll_t
*pcc
, t_selelem
*sel
, const char *rpost
)
585 if (sel
->u
.expr
.method
->pupdate
)
587 /* By default, use whole residues/molecules. */
588 rc
= gmx_ana_poscalc_create_enum(&sel
->u
.expr
.pc
, pcc
, rpost
,
593 _gmx_selparser_error("warning: '%s' modifier for '%s' ignored",
594 rpost
, sel
->u
.expr
.method
->name
);
600 * \param[in] left Selection element for the left hand side.
601 * \param[in] right Selection element for the right hand side.
602 * \param[in] op String representation of the operator.
603 * \param[in] scanner Scanner data structure.
604 * \returns The created selection element.
606 * This function handles the creation of a \c t_selelem object for
607 * arithmetic expressions.
610 _gmx_sel_init_arithmetic(t_selelem
*left
, t_selelem
*right
, char op
,
618 sel
= _gmx_selelem_create(SEL_ARITHMETIC
);
619 sel
->v
.type
= REAL_VALUE
;
622 case '+': sel
->u
.arith
.type
= ARITH_PLUS
; break;
623 case '-': sel
->u
.arith
.type
= (right
? ARITH_MINUS
: ARITH_NEG
); break;
624 case '*': sel
->u
.arith
.type
= ARITH_MULT
; break;
625 case '/': sel
->u
.arith
.type
= ARITH_DIV
; break;
626 case '^': sel
->u
.arith
.type
= ARITH_EXP
; break;
628 sel
->u
.arith
.opstr
= strdup(buf
);
629 sel
->name
= sel
->u
.arith
.opstr
;
631 sel
->child
->next
= right
;
636 * \param[in] left Selection element for the left hand side.
637 * \param[in] right Selection element for the right hand side.
638 * \param[in] cmpop String representation of the comparison operator.
639 * \param[in] scanner Scanner data structure.
640 * \returns The created selection element.
642 * This function handles the creation of a \c t_selelem object for
643 * comparison expressions.
646 _gmx_sel_init_comparison(t_selelem
*left
, t_selelem
*right
, char *cmpop
,
650 t_selexpr_param
*params
, *param
;
654 sel
= _gmx_selelem_create(SEL_EXPRESSION
);
655 _gmx_selelem_set_method(sel
, &sm_compare
, scanner
);
656 /* Create the parameter for the left expression */
657 name
= left
->v
.type
== INT_VALUE
? "int1" : "real1";
658 params
= param
= _gmx_selexpr_create_param(strdup(name
));
660 param
->value
= _gmx_selexpr_create_value_expr(left
);
661 /* Create the parameter for the right expression */
662 name
= right
->v
.type
== INT_VALUE
? "int2" : "real2";
663 param
= _gmx_selexpr_create_param(strdup(name
));
665 param
->value
= _gmx_selexpr_create_value_expr(right
);
666 params
->next
= param
;
667 /* Create the parameter for the operator */
668 param
= _gmx_selexpr_create_param(strdup("op"));
670 param
->value
= _gmx_selexpr_create_value(STR_VALUE
);
671 param
->value
->u
.s
= cmpop
;
672 params
->next
->next
= param
;
673 if (!_gmx_sel_parse_params(params
, sel
->u
.expr
.method
->nparams
,
674 sel
->u
.expr
.method
->param
, sel
, scanner
))
676 _gmx_selparser_error("error in comparison initialization");
677 _gmx_selelem_free(sel
);
685 * \param[in] method Method to use.
686 * \param[in] args Pointer to the first argument.
687 * \param[in] rpost Reference position type to use (NULL = default).
688 * \param[in] scanner Scanner data structure.
689 * \returns The created selection element.
691 * This function handles the creation of a \c t_selelem object for
692 * selection methods that do not take parameters.
695 _gmx_sel_init_keyword(gmx_ana_selmethod_t
*method
, t_selexpr_value
*args
,
696 const char *rpost
, yyscan_t scanner
)
698 gmx_ana_selcollection_t
*sc
= _gmx_sel_lexer_selcollection(scanner
);
699 t_selelem
*root
, *child
;
700 t_selexpr_param
*params
, *param
;
701 t_selexpr_value
*arg
;
705 if (method
->nparams
> 0)
707 gmx_bug("internal error");
711 root
= _gmx_selelem_create(SEL_EXPRESSION
);
713 _gmx_selelem_set_method(child
, method
, scanner
);
715 /* Initialize the evaluation of keyword matching if values are provided */
718 gmx_ana_selmethod_t
*kwmethod
;
719 switch (method
->type
)
721 case INT_VALUE
: kwmethod
= &sm_keyword_int
; break;
722 case REAL_VALUE
: kwmethod
= &sm_keyword_real
; break;
723 case STR_VALUE
: kwmethod
= &sm_keyword_str
; break;
725 _gmx_selparser_error("unknown type for keyword selection");
726 _gmx_selexpr_free_values(args
);
729 /* Count the arguments */
737 /* Initialize the selection element */
738 root
= _gmx_selelem_create(SEL_EXPRESSION
);
739 _gmx_selelem_set_method(root
, kwmethod
, scanner
);
740 params
= param
= _gmx_selexpr_create_param(NULL
);
742 param
->value
= _gmx_selexpr_create_value_expr(child
);
743 param
= _gmx_selexpr_create_param(NULL
);
746 params
->next
= param
;
747 if (!_gmx_sel_parse_params(params
, root
->u
.expr
.method
->nparams
,
748 root
->u
.expr
.method
->param
, root
, scanner
))
750 _gmx_selparser_error("error in keyword selection initialization");
754 rc
= set_refpos_type(sc
->pcc
, child
, rpost
);
762 /* On error, free all memory and return NULL. */
764 _gmx_selelem_free(root
);
769 * \param[in] method Method to use for initialization.
770 * \param[in] params Pointer to the first parameter.
771 * \param[in] rpost Reference position type to use (NULL = default).
772 * \param[in] scanner Scanner data structure.
773 * \returns The created selection element.
775 * This function handles the creation of a \c t_selelem object for
776 * selection methods that take parameters.
778 * Part of the behavior of the \c same selection keyword is hardcoded into
779 * this function (or rather, into _gmx_selelem_custom_init_same()) to allow the
780 * use of any keyword in \c "same KEYWORD as" without requiring special
781 * handling somewhere else (or sacrificing the simple syntax).
784 _gmx_sel_init_method(gmx_ana_selmethod_t
*method
, t_selexpr_param
*params
,
785 const char *rpost
, yyscan_t scanner
)
787 gmx_ana_selcollection_t
*sc
= _gmx_sel_lexer_selcollection(scanner
);
791 _gmx_sel_finish_method(scanner
);
792 /* The "same" keyword needs some custom massaging of the parameters. */
793 rc
= _gmx_selelem_custom_init_same(&method
, params
, scanner
);
796 _gmx_selexpr_free_params(params
);
799 root
= _gmx_selelem_create(SEL_EXPRESSION
);
800 _gmx_selelem_set_method(root
, method
, scanner
);
801 /* Process the parameters */
802 if (!_gmx_sel_parse_params(params
, root
->u
.expr
.method
->nparams
,
803 root
->u
.expr
.method
->param
, root
, scanner
))
805 _gmx_selelem_free(root
);
808 rc
= set_refpos_type(sc
->pcc
, root
, rpost
);
811 _gmx_selelem_free(root
);
819 * \param[in] method Modifier to use for initialization.
820 * \param[in] params Pointer to the first parameter.
821 * \param[in] sel Selection element that the modifier should act on.
822 * \param[in] scanner Scanner data structure.
823 * \returns The created selection element.
825 * This function handles the creation of a \c t_selelem object for
826 * selection modifiers.
829 _gmx_sel_init_modifier(gmx_ana_selmethod_t
*method
, t_selexpr_param
*params
,
830 t_selelem
*sel
, yyscan_t scanner
)
834 t_selexpr_param
*vparam
;
837 _gmx_sel_finish_method(scanner
);
838 mod
= _gmx_selelem_create(SEL_MODIFIER
);
839 _gmx_selelem_set_method(mod
, method
, scanner
);
840 if (method
->type
== NO_VALUE
)
854 vparam
= _gmx_selexpr_create_param(NULL
);
856 vparam
->value
= _gmx_selexpr_create_value_expr(sel
);
857 vparam
->next
= params
;
861 /* Process the parameters */
862 if (!_gmx_sel_parse_params(params
, mod
->u
.expr
.method
->nparams
,
863 mod
->u
.expr
.method
->param
, mod
, scanner
))
865 if (mod
->child
!= sel
)
867 _gmx_selelem_free(sel
);
869 _gmx_selelem_free(mod
);
877 * \param[in] expr Input selection element for the position calculation.
878 * \param[in] type Reference position type or NULL for default.
879 * \param[in] scanner Scanner data structure.
880 * \returns The created selection element.
882 * This function handles the creation of a \c t_selelem object for
883 * evaluation of reference positions.
886 _gmx_sel_init_position(t_selelem
*expr
, const char *type
, yyscan_t scanner
)
889 t_selexpr_param
*params
;
891 root
= _gmx_selelem_create(SEL_EXPRESSION
);
892 _gmx_selelem_set_method(root
, &sm_keyword_pos
, scanner
);
893 _gmx_selelem_set_kwpos_type(root
, type
);
894 /* Create the parameters for the parameter parser. */
895 params
= _gmx_selexpr_create_param(NULL
);
897 params
->value
= _gmx_selexpr_create_value_expr(expr
);
898 /* Parse the parameters. */
899 if (!_gmx_sel_parse_params(params
, root
->u
.expr
.method
->nparams
,
900 root
->u
.expr
.method
->param
, root
, scanner
))
902 _gmx_selelem_free(root
);
910 * \param[in] x,y,z Coordinates for the position.
911 * \returns The creates selection element.
914 _gmx_sel_init_const_position(real x
, real y
, real z
)
919 sel
= _gmx_selelem_create(SEL_CONST
);
920 _gmx_selelem_set_vtype(sel
, POS_VALUE
);
921 _gmx_selvalue_reserve(&sel
->v
, 1);
925 gmx_ana_pos_init_const(sel
->v
.u
.p
, pos
);
930 * \param[in] name Name of an index group to search for.
931 * \param[in] scanner Scanner data structure.
932 * \returns The created constant selection element, or NULL if no matching
935 * See gmx_ana_indexgrps_find() for information on how \p name is matched
936 * against the index groups.
939 _gmx_sel_init_group_by_name(const char *name
, yyscan_t scanner
)
941 gmx_ana_indexgrps_t
*grps
= _gmx_sel_lexer_indexgrps(scanner
);
948 sel
= _gmx_selelem_create(SEL_CONST
);
949 _gmx_selelem_set_vtype(sel
, GROUP_VALUE
);
950 /* FIXME: The constness should not be cast away */
951 if (!gmx_ana_indexgrps_find(&sel
->u
.cgrp
, grps
, (char *)name
))
953 _gmx_selelem_free(sel
);
956 sel
->name
= sel
->u
.cgrp
.name
;
961 * \param[in] id Zero-based index number of the group to extract.
962 * \param[in] scanner Scanner data structure.
963 * \returns The created constant selection element, or NULL if no matching
967 _gmx_sel_init_group_by_id(int id
, yyscan_t scanner
)
969 gmx_ana_indexgrps_t
*grps
= _gmx_sel_lexer_indexgrps(scanner
);
976 sel
= _gmx_selelem_create(SEL_CONST
);
977 _gmx_selelem_set_vtype(sel
, GROUP_VALUE
);
978 if (!gmx_ana_indexgrps_extract(&sel
->u
.cgrp
, grps
, id
))
980 _gmx_selelem_free(sel
);
983 sel
->name
= sel
->u
.cgrp
.name
;
988 * \param[in,out] sel Value of the variable.
989 * \returns The created selection element that references \p sel.
991 * The reference count of \p sel is updated, but no other modifications are
995 _gmx_sel_init_variable_ref(t_selelem
*sel
)
999 if (sel
->v
.type
== POS_VALUE
&& sel
->type
== SEL_CONST
)
1005 ref
= _gmx_selelem_create(SEL_SUBEXPRREF
);
1006 _gmx_selelem_set_vtype(ref
, sel
->v
.type
);
1007 ref
->name
= sel
->name
;
1015 * Initializes default values for position keyword evaluation.
1017 * \param[in,out] root Root of the element tree to initialize.
1018 * \param[in] sc Selection collection to use defaults from.
1019 * \param[in] bSelection Whether the element evaluates the positions for a
1023 init_pos_keyword_defaults(t_selelem
*root
, gmx_ana_selcollection_t
*sc
, bool bSelection
)
1028 /* Selections use largest static group by default, while
1029 * reference positions use the whole residue/molecule. */
1030 if (root
->type
== SEL_EXPRESSION
)
1032 flags
= bSelection
? POS_COMPLMAX
: POS_COMPLWHOLE
;
1033 if (bSelection
&& sc
->bMaskOnly
)
1035 flags
|= POS_MASKONLY
;
1037 if (bSelection
&& sc
->bVelocities
)
1039 flags
|= POS_VELOCITIES
;
1041 if (bSelection
&& sc
->bForces
)
1043 flags
|= POS_FORCES
;
1045 _gmx_selelem_set_kwpos_type(root
, bSelection
? sc
->spost
: sc
->rpost
);
1046 _gmx_selelem_set_kwpos_flags(root
, flags
);
1048 /* Change the defaults once we are no longer processing modifiers */
1049 if (root
->type
!= SEL_ROOT
&& root
->type
!= SEL_MODIFIER
1050 && root
->type
!= SEL_SUBEXPRREF
&& root
->type
!= SEL_SUBEXPR
)
1054 /* Recurse into children */
1055 child
= root
->child
;
1058 init_pos_keyword_defaults(child
, sc
, bSelection
);
1059 child
= child
->next
;
1064 * \param[in] name Name for the selection
1065 * (if NULL, a default name is constructed).
1066 * \param[in] sel The selection element that evaluates the selection.
1067 * \param scanner Scanner data structure.
1068 * \returns The created root selection element.
1070 * This function handles the creation of root (\ref SEL_ROOT) \c t_selelem
1071 * objects for selections.
1074 _gmx_sel_init_selection(char *name
, t_selelem
*sel
, yyscan_t scanner
)
1076 gmx_ana_selcollection_t
*sc
= _gmx_sel_lexer_selcollection(scanner
);
1080 if (sel
->v
.type
!= POS_VALUE
)
1082 gmx_bug("each selection must evaluate to a position");
1083 /* FIXME: Better handling of this error */
1088 root
= _gmx_selelem_create(SEL_ROOT
);
1090 /* Assign the name (this is done here to free it automatically in the case
1091 * of an error below). */
1094 root
->name
= root
->u
.cgrp
.name
= name
;
1096 /* Update the flags */
1097 rc
= _gmx_selelem_update_flags(root
);
1100 _gmx_selelem_free(root
);
1103 /* Initialize defaults for position keywords */
1104 init_pos_keyword_defaults(sel
, sc
, TRUE
);
1106 /* If there is no name provided by the user, check whether the actual
1107 * selection given was from an external group, and if so, use the name
1108 * of the external group. */
1111 t_selelem
*child
= root
->child
;
1112 while (child
->type
== SEL_MODIFIER
)
1114 if (!child
->child
|| child
->child
->type
!= SEL_SUBEXPRREF
1115 || !child
->child
->child
)
1119 child
= child
->child
->child
;
1121 if (child
->type
== SEL_EXPRESSION
1122 && child
->child
&& child
->child
->type
== SEL_SUBEXPRREF
1123 && child
->child
->child
1124 && child
->child
->child
->type
== SEL_CONST
1125 && child
->child
->child
->v
.type
== GROUP_VALUE
)
1127 root
->name
= root
->u
.cgrp
.name
=
1128 strdup(child
->child
->child
->u
.cgrp
.name
);
1131 /* If there still is no name, use the selection string */
1134 root
->name
= root
->u
.cgrp
.name
1135 = strdup(_gmx_sel_lexer_pselstr(scanner
));
1138 /* Print out some information if the parser is interactive */
1139 if (_gmx_sel_is_lexer_interactive(scanner
))
1141 fprintf(stderr
, "Selection '%s' parsed\n",
1142 _gmx_sel_lexer_pselstr(scanner
));
1150 * \param[in] name Name of the variable (should not be freed after this
1152 * \param[in] expr The selection element that evaluates the variable.
1153 * \param scanner Scanner data structure.
1154 * \returns The created root selection element.
1156 * This function handles the creation of root \c t_selelem objects for
1157 * variable assignments. A \ref SEL_ROOT element and a \ref SEL_SUBEXPR
1158 * element are both created.
1161 _gmx_sel_assign_variable(char *name
, t_selelem
*expr
, yyscan_t scanner
)
1163 gmx_ana_selcollection_t
*sc
= _gmx_sel_lexer_selcollection(scanner
);
1164 const char *pselstr
= _gmx_sel_lexer_pselstr(scanner
);
1165 t_selelem
*root
= NULL
;
1168 rc
= _gmx_selelem_update_flags(expr
);
1172 _gmx_selelem_free(expr
);
1175 /* Check if this is a constant non-group value */
1176 if (expr
->type
== SEL_CONST
&& expr
->v
.type
!= GROUP_VALUE
)
1178 /* If so, just assign the constant value to the variable */
1179 if (!_gmx_sel_add_var_symbol(sc
->symtab
, name
, expr
))
1181 _gmx_selelem_free(expr
);
1185 _gmx_selelem_free(expr
);
1189 /* Check if we are assigning a variable to another variable */
1190 if (expr
->type
== SEL_SUBEXPRREF
)
1192 /* If so, make a simple alias */
1193 if (!_gmx_sel_add_var_symbol(sc
->symtab
, name
, expr
->child
))
1195 _gmx_selelem_free(expr
);
1199 _gmx_selelem_free(expr
);
1203 /* Create the root element */
1204 root
= _gmx_selelem_create(SEL_ROOT
);
1206 root
->u
.cgrp
.name
= name
;
1207 /* Create the subexpression element */
1208 root
->child
= _gmx_selelem_create(SEL_SUBEXPR
);
1209 _gmx_selelem_set_vtype(root
->child
, expr
->v
.type
);
1210 root
->child
->name
= name
;
1211 root
->child
->child
= expr
;
1213 rc
= _gmx_selelem_update_flags(root
);
1216 _gmx_selelem_free(root
);
1219 /* Add the variable to the symbol table */
1220 if (!_gmx_sel_add_var_symbol(sc
->symtab
, name
, root
->child
))
1222 _gmx_selelem_free(root
);
1226 srenew(sc
->varstrs
, sc
->nvars
+ 1);
1227 sc
->varstrs
[sc
->nvars
] = strdup(pselstr
);
1229 if (_gmx_sel_is_lexer_interactive(scanner
))
1231 fprintf(stderr
, "Variable '%s' parsed\n", pselstr
);
1237 * \param sel Selection to append (can be NULL, in which
1238 * case nothing is done).
1239 * \param last Last selection, or NULL if not present or not known.
1240 * \param scanner Scanner data structure.
1241 * \returns The last selection after the append.
1243 * Appends \p sel after the last root element, and returns either \p sel
1244 * (if it was non-NULL) or the last element (if \p sel was NULL).
1247 _gmx_sel_append_selection(t_selelem
*sel
, t_selelem
*last
, yyscan_t scanner
)
1249 gmx_ana_selcollection_t
*sc
= _gmx_sel_lexer_selcollection(scanner
);
1251 /* Append sel after last, or the last element of sc if last is NULL */
1272 /* Initialize a selection object if necessary */
1276 /* Add the new selection to the collection if it is not a variable. */
1277 if (sel
->child
->type
!= SEL_SUBEXPR
)
1282 srenew(sc
->sel
, sc
->nr
);
1284 snew(sc
->sel
[i
], 1);
1285 sc
->sel
[i
]->name
= strdup(sel
->name
);
1286 sc
->sel
[i
]->selstr
= strdup(_gmx_sel_lexer_pselstr(scanner
));
1288 if (sel
->child
->type
== SEL_CONST
)
1290 gmx_ana_pos_copy(&sc
->sel
[i
]->p
, sel
->child
->v
.u
.p
, TRUE
);
1291 sc
->sel
[i
]->bDynamic
= FALSE
;
1298 child
->flags
&= ~SEL_ALLOCVAL
;
1299 _gmx_selvalue_setstore(&child
->v
, &sc
->sel
[i
]->p
);
1300 /* We should also skip any modifiers to determine the dynamic
1302 while (child
->type
== SEL_MODIFIER
)
1304 child
= child
->child
;
1305 if (child
->type
== SEL_SUBEXPRREF
)
1307 child
= child
->child
->child
;
1310 /* For variable references, we should skip the
1311 * SEL_SUBEXPRREF and SEL_SUBEXPR elements. */
1312 if (child
->type
== SEL_SUBEXPRREF
)
1314 child
= child
->child
->child
;
1316 sc
->sel
[i
]->bDynamic
= (child
->child
->flags
& SEL_DYNAMIC
);
1318 /* The group will be set after compilation */
1319 sc
->sel
[i
]->g
= NULL
;
1320 sc
->sel
[i
]->selelem
= sel
;
1321 gmx_ana_selection_init_coverfrac(sc
->sel
[i
], CFRAC_NONE
);
1324 /* Clear the selection string now that we've saved it */
1325 _gmx_sel_lexer_clear_pselstr(scanner
);
1330 * \param[in] scanner Scanner data structure.
1331 * \returns TRUE if the parser should finish, FALSE if parsing should
1334 * This function is called always after _gmx_sel_append_selection() to
1335 * check whether a sufficient number of selections has already been provided.
1336 * This is used to terminate interactive parsers when the correct number of
1337 * selections has been provided.
1340 _gmx_sel_parser_should_finish(yyscan_t scanner
)
1342 gmx_ana_selcollection_t
*sc
= _gmx_sel_lexer_selcollection(scanner
);
1343 return sc
->nr
== _gmx_sel_lexer_exp_selcount(scanner
);
1347 * \param[in] scanner Scanner data structure.
1350 _gmx_sel_handle_empty_cmd(yyscan_t scanner
)
1352 gmx_ana_selcollection_t
*sc
= _gmx_sel_lexer_selcollection(scanner
);
1353 gmx_ana_indexgrps_t
*grps
= _gmx_sel_lexer_indexgrps(scanner
);
1356 if (!_gmx_sel_is_lexer_interactive(scanner
))
1361 fprintf(stderr
, "Available index groups:\n");
1362 gmx_ana_indexgrps_print(_gmx_sel_lexer_indexgrps(scanner
), 0);
1364 if (sc
->nvars
> 0 || sc
->nr
> 0)
1366 fprintf(stderr
, "Currently provided selections:\n");
1367 for (i
= 0; i
< sc
->nvars
; ++i
)
1369 fprintf(stderr
, " %s\n", sc
->varstrs
[i
]);
1371 for (i
= 0; i
< sc
->nr
; ++i
)
1373 fprintf(stderr
, " %2d. %s\n", i
+1, sc
->sel
[i
]->selstr
);
1379 * \param[in] topic Topic for which help was requested, or NULL for general
1381 * \param[in] scanner Scanner data structure.
1383 * \p topic is freed by this function.
1386 _gmx_sel_handle_help_cmd(char *topic
, yyscan_t scanner
)
1388 gmx_ana_selcollection_t
*sc
= _gmx_sel_lexer_selcollection(scanner
);
1390 _gmx_sel_print_help(sc
, topic
);
1398 * Internal helper function used by gmx_ana_selcollection_parse_*() to do the actual work.
1400 * \param[in] maxnr Maximum number of selections to parse
1401 * (if -1, parse as many as provided by the user).
1402 * \param[in,out] scanner Scanner data structure.
1403 * \returns 0 on success, -1 on error.
1406 run_parser(int maxnr
, yyscan_t scanner
)
1408 gmx_ana_selcollection_t
*sc
= _gmx_sel_lexer_selcollection(scanner
);
1413 bOk
= !_gmx_sel_yybparse(scanner
);
1414 _gmx_sel_free_lexer(scanner
);
1416 if (maxnr
> 0 && nr
!= maxnr
)
1420 return bOk
? 0 : -1;
1424 * \param[in,out] sc Selection collection to use for output.
1425 * \param[in] nr Number of selections to parse
1426 * (if -1, parse as many as provided by the user).
1427 * \param[in] grps External index groups (can be NULL).
1428 * \param[in] bInteractive Whether the parser should behave interactively.
1429 * \returns 0 on success, -1 on error.
1431 * The number of selections parsed can be accessed with
1432 * gmx_ana_selcollection_get_count() (note that if you call the parser
1433 * multiple times, this function returns the total count).
1436 gmx_ana_selcollection_parse_stdin(gmx_ana_selcollection_t
*sc
, int nr
,
1437 gmx_ana_indexgrps_t
*grps
, bool bInteractive
)
1442 rc
= _gmx_sel_init_lexer(&scanner
, sc
, bInteractive
, nr
, grps
);
1447 /* We don't set the lexer input here, which causes it to use a special
1448 * internal implementation for reading from stdin. */
1449 return run_parser(nr
, scanner
);
1453 * \param[in,out] sc Selection collection to use for output.
1454 * \param[in] fnm Name of the file to parse selections from.
1455 * \param[in] grps External index groups (can be NULL).
1456 * \returns 0 on success, -1 on error.
1458 * The number of selections parsed can be accessed with
1459 * gmx_ana_selcollection_get_count() (note that if you call the parser
1460 * multiple times, this function returns the total count).
1463 gmx_ana_selcollection_parse_file(gmx_ana_selcollection_t
*sc
, const char *fnm
,
1464 gmx_ana_indexgrps_t
*grps
)
1470 rc
= _gmx_sel_init_lexer(&scanner
, sc
, FALSE
, -1, grps
);
1475 fp
= ffopen(fnm
, "r");
1476 _gmx_sel_set_lex_input_file(scanner
, fp
);
1477 rc
= run_parser(-1, scanner
);
1483 * \param[in,out] sc Selection collection to use for output.
1484 * \param[in] str String to parse selections from.
1485 * \param[in] grps External index groups (can be NULL).
1486 * \returns 0 on success, -1 on error.
1488 * The number of selections parsed can be accessed with
1489 * gmx_ana_selcollection_get_count() (note that if you call the parser
1490 * multiple times, this function returns the total count).
1493 gmx_ana_selcollection_parse_str(gmx_ana_selcollection_t
*sc
, const char *str
,
1494 gmx_ana_indexgrps_t
*grps
)
1499 rc
= _gmx_sel_init_lexer(&scanner
, sc
, FALSE
, -1, grps
);
1504 _gmx_sel_set_lex_input_str(scanner
, str
);
1505 return run_parser(-1, scanner
);