A bit more modular selection position handling.
[gromacs/qmmm-gamess-us.git] / src / gmxlib / selection / parsetree.c
blob34f2182bd392947892add9071c0e29455ca22e2e
1 /*
3 * This source code is part of
5 * G R O M A C S
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
31 /*! \internal \file
32 * \brief Implementation of functions of in parsetree.h.
34 /*! \internal
35 * \page selparser Selection parsing
37 * The selection parser is implemented in the following files:
38 * - scanner.l:
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.
50 * - parser.y:
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
58 * in selection.h.
59 * - params.c:
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
75 * appropriately.
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:
80 * \ref selcompiler
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
90 * discussed below.
91 * Fields that are not initialized are set to zero, NULL, or other similar
92 * value.
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
105 * as in the input.
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
130 * groups,
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
166 * the variable.
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
181 * the expression.
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
195 * are possible.
197 #ifdef HAVE_CONFIG_H
198 #include <config.h>
199 #endif
201 #include <stdio.h>
203 #include <futil.h>
204 #include <smalloc.h>
205 #include <string2.h>
207 #include <poscalc.h>
208 #include <selection.h>
209 #include <selmethod.h>
211 #include "keywords.h"
212 #include "parsetree.h"
213 #include "selcollection.h"
214 #include "selelem.h"
215 #include "selhelp.h"
216 #include "symrec.h"
218 #include "scanner.h"
220 /* In parser.y */
221 /*! \brief
222 * Parser function generated by Bison.
225 _gmx_sel_yyparse(void *scanner);
228 * It is a simple wrapper for fprintf(stderr, ...).
230 void
231 _gmx_selparser_error(const char *fmt, ...)
233 va_list ap;
234 va_start(ap, fmt);
235 fprintf(stderr, "selection parser: ");
236 vfprintf(stderr, fmt, ap);
237 fprintf(stderr, "\n");
238 va_end(ap);
242 * \param[in] type Type for the new value.
243 * \returns Pointer to the newly allocated value.
245 t_selexpr_value *
246 _gmx_selexpr_create_value(e_selvalue_t type)
248 t_selexpr_value *value;
249 snew(value, 1);
250 value->type = type;
251 value->bExpr = FALSE;
252 value->next = NULL;
253 return value;
257 * \param[in] expr Expression for the value.
258 * \returns Pointer to the newly allocated value.
260 t_selexpr_value *
261 _gmx_selexpr_create_value_expr(t_selelem *expr)
263 t_selexpr_value *value;
264 snew(value, 1);
265 value->type = expr->v.type;
266 value->bExpr = TRUE;
267 value->u.expr = expr;
268 value->next = NULL;
269 return value;
273 * \param[in] name Name for the new parameter.
274 * \returns Pointer to the newly allocated parameter.
276 * No copy of \p name is made.
278 t_selexpr_param *
279 _gmx_selexpr_create_param(char *name)
281 t_selexpr_param *param;
282 snew(param, 1);
283 param->name = name;
284 param->next = NULL;
285 return param;
289 * \param value Pointer to the beginning of the value list to free.
291 * The expressions referenced by the values are also freed
292 * (to prevent this, set the expression to NULL before calling the function).
294 void
295 _gmx_selexpr_free_values(t_selexpr_value *value)
297 t_selexpr_value *old;
299 while (value)
301 if (value->bExpr)
303 if (value->u.expr)
305 _gmx_selelem_free(value->u.expr);
308 else if (value->type == STR_VALUE)
310 sfree(value->u.s);
312 old = value;
313 value = value->next;
314 sfree(old);
319 * \param param Pointer the the beginning of the parameter list to free.
321 * The values of the parameters are freed with free_selexpr_values().
323 void
324 _gmx_selexpr_free_params(t_selexpr_param *param)
326 t_selexpr_param *old;
328 while (param)
330 _gmx_selexpr_free_values(param->value);
331 old = param;
332 param = param->next;
333 sfree(old->name);
334 sfree(old);
339 * \param[in,out] sel Root of the selection element tree to initialize.
340 * \returns 0 on success, an error code on error.
342 * Propagates the \ref SEL_DYNAMIC flag from the children of \p sel to \p sel
343 * (if any child of \p sel is dynamic, \p sel is also marked as such).
344 * The \ref SEL_DYNAMIC flag is also set for \ref SEL_EXPRESSION elements with
345 * a dynamic method.
346 * Also, sets one of the \ref SEL_SINGLEVAL, \ref SEL_ATOMVAL, or
347 * \ref SEL_VARNUMVAL flags, either based on the children or on the type of
348 * the selection method.
349 * If the types of the children conflict, an error is returned.
351 * The flags of the children of \p sel are also updated if not done earlier.
352 * The flags are initialized only once for any element; if \ref SEL_FLAGSSET
353 * is set for an element, the function returns immediately, and the recursive
354 * operation does not descend beyond such elements.
357 _gmx_selelem_update_flags(t_selelem *sel)
359 t_selelem *child;
360 int rc;
361 bool bUseChildType;
363 /* Return if the flags have already been set */
364 if (sel->flags & SEL_FLAGSSET)
366 return 0;
368 /* Set the flags based on the current element type */
369 switch (sel->type)
371 case SEL_CONST:
372 sel->flags |= SEL_SINGLEVAL;
373 bUseChildType = FALSE;
374 break;
376 case SEL_EXPRESSION:
377 if (sel->u.expr.method->flags & SMETH_DYNAMIC)
379 sel->flags |= SEL_DYNAMIC;
381 if (sel->u.expr.method->flags & SMETH_SINGLEVAL)
383 sel->flags |= SEL_SINGLEVAL;
385 else if (sel->u.expr.method->flags & SMETH_VARNUMVAL)
387 sel->flags |= SEL_VARNUMVAL;
389 else
391 sel->flags |= SEL_ATOMVAL;
393 bUseChildType = FALSE;
394 break;
396 case SEL_MODIFIER:
397 if (sel->v.type != NO_VALUE)
399 sel->flags |= SEL_VARNUMVAL;
401 bUseChildType = FALSE;
402 break;
404 case SEL_ROOT:
405 bUseChildType = FALSE;
406 break;
408 default:
409 bUseChildType = TRUE;
410 break;
412 /* Loop through children to propagate their flags upwards */
413 child = sel->child;
414 while (child)
416 /* Update the child */
417 rc = _gmx_selelem_update_flags(child);
418 if (rc != 0)
420 return rc;
422 /* Propagate the dynamic flag */
423 sel->flags |= (child->flags & SEL_DYNAMIC);
424 /* Propagate the type flag if necessary and check for problems */
425 if (bUseChildType)
427 if ((sel->flags & SEL_VALTYPEMASK)
428 && !(sel->flags & child->flags & SEL_VALTYPEMASK))
430 _gmx_selparser_error("invalid combination of selection expressions");
431 return EINVAL;
433 sel->flags |= (child->flags & SEL_VALTYPEMASK);
436 child = child->next;
438 /* Mark that the flags are set */
439 sel->flags |= SEL_FLAGSSET;
440 /* For root elements, the type should be propagated here, after the
441 * children have been updated. */
442 if (sel->type == SEL_ROOT)
444 sel->flags |= (sel->child->flags & SEL_VALTYPEMASK);
446 return 0;
449 /*! \brief
450 * Initializes the method parameter data of \ref SEL_EXPRESSION and
451 * \ref SEL_MODIFIER elements.
453 * \param[in] sc Selection collection.
454 * \param[in,out] sel Selection element to initialize.
456 * A deep copy of the parameters is made to allow several
457 * expressions with the same method to coexist peacefully.
458 * Calls sel_datafunc() if one is specified for the method.
460 static void
461 init_method_params(gmx_ana_selcollection_t *sc, t_selelem *sel)
463 int nparams;
464 gmx_ana_selparam_t *orgparam;
465 gmx_ana_selparam_t *param;
466 int i;
467 void *mdata;
469 nparams = sel->u.expr.method->nparams;
470 orgparam = sel->u.expr.method->param;
471 snew(param, nparams);
472 memcpy(param, orgparam, nparams*sizeof(gmx_ana_selparam_t));
473 for (i = 0; i < nparams; ++i)
475 param[i].flags &= ~SPAR_SET;
476 _gmx_selvalue_setstore(&param[i].val, NULL);
477 if (param[i].flags & SPAR_VARNUM)
479 param[i].val.nr = -1;
481 /* Duplicate the enum value array if it is given statically */
482 if ((param[i].flags & SPAR_ENUMVAL) && orgparam[i].val.u.ptr != NULL)
484 int n;
486 /* Count the values */
487 n = 1;
488 while (orgparam[i].val.u.s[n] != NULL)
490 ++n;
492 _gmx_selvalue_reserve(&param[i].val, n+1);
493 memcpy(param[i].val.u.s, orgparam[i].val.u.s,
494 (n+1)*sizeof(param[i].val.u.s[0]));
497 mdata = NULL;
498 if (sel->u.expr.method->init_data)
500 mdata = sel->u.expr.method->init_data(nparams, param);
501 if (mdata == NULL)
503 gmx_fatal(FARGS, "Method data initialization failed");
506 if (sel->u.expr.method->set_poscoll)
508 sel->u.expr.method->set_poscoll(sc->pcc, mdata);
510 /* Store the values */
511 sel->u.expr.method->param = param;
512 sel->u.expr.mdata = mdata;
515 /*! \brief
516 * Initializes the method for a \ref SEL_EXPRESSION selection element.
518 * \param[in] sc Selection collection.
519 * \param[in,out] sel Selection element to initialize.
520 * \param[in] method Selection method to set.
522 * Makes a copy of \p method and stores it in \p sel->u.expr.method,
523 * and calls init_method_params();
525 static void
526 set_method(gmx_ana_selcollection_t *sc, t_selelem *sel,
527 gmx_ana_selmethod_t *method)
529 int i;
531 _gmx_selelem_set_vtype(sel, method->type);
532 sel->name = method->name;
533 snew(sel->u.expr.method, 1);
534 memcpy(sel->u.expr.method, method, sizeof(gmx_ana_selmethod_t));
535 init_method_params(sc, sel);
538 /*! \brief
539 * Initializes the reference position calculation for a \ref SEL_EXPRESSION
540 * element.
542 * \param[in,out] pcc Position calculation collection to use.
543 * \param[in,out] sel Selection element to initialize.
544 * \param[in] rpost Reference position type to use (NULL = default).
545 * \returns 0 on success, a non-zero error code on error.
547 static int
548 set_refpos_type(gmx_ana_poscalc_coll_t *pcc, t_selelem *sel, const char *rpost)
550 int rc;
552 if (!rpost)
554 return 0;
557 rc = 0;
558 if (sel->u.expr.method->pupdate)
560 /* By default, use whole residues/molecules. */
561 rc = gmx_ana_poscalc_create_enum(&sel->u.expr.pc, pcc, rpost,
562 POS_COMPLWHOLE);
564 else
566 _gmx_selparser_error("warning: '%d' modifier for '%s' ignored",
567 rpost, sel->u.expr.method->name);
569 return rc;
573 * \param[in] left Selection element for the left hand side.
574 * \param[in] right Selection element for the right hand side.
575 * \param[in] cmpop String representation of the comparison operator.
576 * \param[in] scanner Scanner data structure.
577 * \returns The created selection element.
579 * This function handles the creation of a \c t_selelem object for
580 * comparison expressions.
582 t_selelem *
583 _gmx_sel_init_comparison(t_selelem *left, t_selelem *right, char *cmpop,
584 yyscan_t scanner)
586 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
587 t_selelem *sel;
588 t_selexpr_param *params, *param;
589 char *name;
590 int rc;
592 sel = _gmx_selelem_create(SEL_EXPRESSION);
593 set_method(sc, sel, &sm_compare);
594 /* Create the parameter for the left expression */
595 name = left->v.type == INT_VALUE ? "int1" : "real1";
596 params = param = _gmx_selexpr_create_param(strdup(name));
597 param->nval = 1;
598 param->value = _gmx_selexpr_create_value_expr(left);
599 /* Create the parameter for the right expression */
600 name = right->v.type == INT_VALUE ? "int2" : "real2";
601 param = _gmx_selexpr_create_param(strdup(name));
602 param->nval = 1;
603 param->value = _gmx_selexpr_create_value_expr(right);
604 params->next = param;
605 /* Create the parameter for the operator */
606 param = _gmx_selexpr_create_param(strdup("op"));
607 param->nval = 1;
608 param->value = _gmx_selexpr_create_value(STR_VALUE);
609 param->value->u.s = cmpop;
610 params->next->next = param;
611 if (!_gmx_sel_parse_params(params, sel->u.expr.method->nparams,
612 sel->u.expr.method->param, sel, scanner))
614 _gmx_selparser_error("error in comparison initialization");
615 _gmx_selelem_free(sel);
616 return NULL;
619 return sel;
623 * \param[in] method Method to use.
624 * \param[in] args Pointer to the first argument.
625 * \param[in] rpost Reference position type to use (NULL = default).
626 * \param[in] scanner Scanner data structure.
627 * \returns The created selection element.
629 * This function handles the creation of a \c t_selelem object for
630 * selection methods that do not take parameters.
632 t_selelem *
633 _gmx_sel_init_keyword(gmx_ana_selmethod_t *method, t_selexpr_value *args,
634 const char *rpost, yyscan_t scanner)
636 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
637 t_selelem *root, *child;
638 t_selexpr_param *params, *param;
639 t_selexpr_value *arg;
640 int nargs;
641 int rc;
643 if (method->nparams > 0)
645 gmx_bug("internal error");
646 return NULL;
649 root = _gmx_selelem_create(SEL_EXPRESSION);
650 child = root;
651 set_method(sc, child, method);
653 /* Initialize the evaluation of keyword matching if values are provided */
654 if (args)
656 gmx_ana_selmethod_t *kwmethod;
657 switch (method->type)
659 case INT_VALUE: kwmethod = &sm_keyword_int; break;
660 case REAL_VALUE: kwmethod = &sm_keyword_real; break;
661 case STR_VALUE: kwmethod = &sm_keyword_str; break;
662 default:
663 _gmx_selparser_error("unknown type for keyword selection");
664 _gmx_selexpr_free_values(args);
665 goto on_error;
667 /* Count the arguments */
668 nargs = 0;
669 arg = args;
670 while (arg)
672 ++nargs;
673 arg = arg->next;
675 /* Initialize the selection element */
676 root = _gmx_selelem_create(SEL_EXPRESSION);
677 set_method(sc, root, kwmethod);
678 params = param = _gmx_selexpr_create_param(NULL);
679 param->nval = 1;
680 param->value = _gmx_selexpr_create_value_expr(child);
681 param = _gmx_selexpr_create_param(NULL);
682 param->nval = nargs;
683 param->value = args;
684 params->next = param;
685 if (!_gmx_sel_parse_params(params, root->u.expr.method->nparams,
686 root->u.expr.method->param, root, scanner))
688 _gmx_selparser_error("error in keyword selection initialization");
689 goto on_error;
692 rc = set_refpos_type(sc->pcc, child, rpost);
693 if (rc != 0)
695 goto on_error;
698 return root;
700 /* On error, free all memory and return NULL. */
701 on_error:
702 _gmx_selelem_free(root);
703 return NULL;
707 * \param[in] method Method to use for initialization.
708 * \param[in] params Pointer to the first parameter.
709 * \param[in] rpost Reference position type to use (NULL = default).
710 * \param[in] scanner Scanner data structure.
711 * \returns The created selection element.
713 * This function handles the creation of a \c t_selelem object for
714 * selection methods that take parameters.
716 t_selelem *
717 _gmx_sel_init_method(gmx_ana_selmethod_t *method, t_selexpr_param *params,
718 const char *rpost, yyscan_t scanner)
720 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
721 t_selelem *root;
722 int rc;
724 root = _gmx_selelem_create(SEL_EXPRESSION);
725 set_method(sc, root, method);
726 /* Process the parameters */
727 if (!_gmx_sel_parse_params(params, root->u.expr.method->nparams,
728 root->u.expr.method->param, root, scanner))
730 _gmx_selelem_free(root);
731 return NULL;
733 rc = set_refpos_type(sc->pcc, root, rpost);
734 if (rc != 0)
736 _gmx_selelem_free(root);
737 return NULL;
740 return root;
744 * \param[in] method Modifier to use for initialization.
745 * \param[in] params Pointer to the first parameter.
746 * \param[in] sel Selection element that the modifier should act on.
747 * \param[in] scanner Scanner data structure.
748 * \returns The created selection element.
750 * This function handles the creation of a \c t_selelem object for
751 * selection modifiers.
753 t_selelem *
754 _gmx_sel_init_modifier(gmx_ana_selmethod_t *method, t_selexpr_param *params,
755 t_selelem *sel, yyscan_t scanner)
757 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
758 t_selelem *root;
759 t_selelem *mod;
760 t_selexpr_param *vparam;
761 int i;
763 mod = _gmx_selelem_create(SEL_MODIFIER);
764 set_method(sc, mod, method);
765 if (method->type == NO_VALUE)
767 t_selelem *child;
769 child = sel;
770 while (child->next)
772 child = child->next;
774 child->next = mod;
775 root = sel;
777 else
779 vparam = _gmx_selexpr_create_param(NULL);
780 vparam->nval = 1;
781 vparam->value = _gmx_selexpr_create_value_expr(sel);
782 vparam->next = params;
783 params = vparam;
784 root = mod;
786 /* Process the parameters */
787 if (!_gmx_sel_parse_params(params, mod->u.expr.method->nparams,
788 mod->u.expr.method->param, mod, scanner))
790 if (mod->child != sel)
792 _gmx_selelem_free(sel);
794 _gmx_selelem_free(mod);
795 return NULL;
798 return root;
802 * \param[in] expr Input selection element for the position calculation.
803 * \param[in] type Reference position type or NULL for default.
804 * \param[in] scanner Scanner data structure.
805 * \returns The created selection element.
807 * This function handles the creation of a \c t_selelem object for
808 * evaluation of reference positions.
810 t_selelem *
811 _gmx_sel_init_position(t_selelem *expr, const char *type, yyscan_t scanner)
813 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
814 t_selelem *root;
815 t_selexpr_param *params;
817 root = _gmx_selelem_create(SEL_EXPRESSION);
818 set_method(sc, root, &sm_keyword_pos);
819 _gmx_selelem_set_kwpos_type(root, type);
820 /* Create the parameters for the parameter parser. */
821 params = _gmx_selexpr_create_param(NULL);
822 params->nval = 1;
823 params->value = _gmx_selexpr_create_value_expr(expr);
824 /* Parse the parameters. */
825 if (!_gmx_sel_parse_params(params, root->u.expr.method->nparams,
826 root->u.expr.method->param, root, scanner))
828 _gmx_selelem_free(root);
829 return NULL;
832 return root;
836 * \param[in] x,y,z Coordinates for the position.
837 * \returns The creates selection element.
839 t_selelem *
840 _gmx_sel_init_const_position(real x, real y, real z)
842 t_selelem *sel;
843 rvec pos;
845 sel = _gmx_selelem_create(SEL_CONST);
846 _gmx_selelem_set_vtype(sel, POS_VALUE);
847 _gmx_selvalue_reserve(&sel->v, 1);
848 pos[XX] = x;
849 pos[YY] = y;
850 pos[ZZ] = z;
851 gmx_ana_pos_init_const(sel->v.u.p, pos);
852 return sel;
856 * \param[in] name Name of an index group to search for.
857 * \param[in] scanner Scanner data structure.
858 * \returns The created constant selection element, or NULL if no matching
859 * index group found.
861 * See gmx_ana_indexgrps_find() for information on how \p name is matched
862 * against the index groups.
864 t_selelem *
865 _gmx_sel_init_group_by_name(const char *name, yyscan_t scanner)
867 gmx_ana_indexgrps_t *grps = _gmx_sel_lexer_indexgrps(scanner);
868 t_selelem *sel;
870 if (!grps)
872 return NULL;
874 sel = _gmx_selelem_create(SEL_CONST);
875 _gmx_selelem_set_vtype(sel, GROUP_VALUE);
876 /* FIXME: The constness should not be cast away */
877 if (!gmx_ana_indexgrps_find(&sel->u.cgrp, grps, (char *)name))
879 _gmx_selelem_free(sel);
880 return NULL;
882 sel->name = sel->u.cgrp.name;
883 return sel;
887 * \param[in] id Zero-based index number of the group to extract.
888 * \param[in] scanner Scanner data structure.
889 * \returns The created constant selection element, or NULL if no matching
890 * index group found.
892 t_selelem *
893 _gmx_sel_init_group_by_id(int id, yyscan_t scanner)
895 gmx_ana_indexgrps_t *grps = _gmx_sel_lexer_indexgrps(scanner);
896 t_selelem *sel;
898 if (!grps)
900 return NULL;
902 sel = _gmx_selelem_create(SEL_CONST);
903 _gmx_selelem_set_vtype(sel, GROUP_VALUE);
904 if (!gmx_ana_indexgrps_extract(&sel->u.cgrp, grps, id))
906 _gmx_selelem_free(sel);
907 return NULL;
909 sel->name = sel->u.cgrp.name;
910 return sel;
914 * \param[in,out] sel Value of the variable.
915 * \returns The created selection element that references \p sel.
917 * The reference count of \p sel is updated, but no other modifications are
918 * made.
920 t_selelem *
921 _gmx_sel_init_variable_ref(t_selelem *sel)
923 t_selelem *ref;
925 if (sel->v.type == POS_VALUE && sel->type == SEL_CONST)
927 ref = sel;
929 else
931 ref = _gmx_selelem_create(SEL_SUBEXPRREF);
932 _gmx_selelem_set_vtype(ref, sel->v.type);
933 ref->name = sel->name;
934 ref->child = sel;
936 sel->refcount++;
937 return ref;
940 /*! \brief
941 * Initializes default values for position keyword evaluation.
943 * \param[in,out] root Root of the element tree to initialize.
944 * \param[in] sc Selection collection to use defaults from.
945 * \param[in] bSelection Whether the element evaluates the positions for a
946 * selection.
948 static void
949 init_pos_keyword_defaults(t_selelem *root, gmx_ana_selcollection_t *sc, bool bSelection)
951 t_selelem *child;
952 int flags;
954 /* Selections use largest static group by default, while
955 * reference positions use the whole residue/molecule. */
956 if (root->type == SEL_EXPRESSION)
958 flags = bSelection ? POS_COMPLMAX : POS_COMPLWHOLE;
959 if (bSelection && sc->bMaskOnly)
961 flags |= POS_MASKONLY;
963 _gmx_selelem_set_kwpos_type(root, bSelection ? sc->spost : sc->rpost);
964 _gmx_selelem_set_kwpos_flags(root, flags);
966 /* Change the defaults once we are no longer processing modifiers */
967 if (root->type != SEL_ROOT && root->type != SEL_MODIFIER
968 && root->type != SEL_SUBEXPRREF && root->type != SEL_SUBEXPR)
970 bSelection = FALSE;
972 /* Recurse into children */
973 child = root->child;
974 while (child)
976 init_pos_keyword_defaults(child, sc, bSelection);
977 child = child->next;
982 * \param[in] name Name for the selection
983 * (if NULL, a default name is constructed).
984 * \param[in] sel The selection element that evaluates the selection.
985 * \param scanner Scanner data structure.
986 * \returns The created root selection element.
988 * This function handles the creation of root (\ref SEL_ROOT) \c t_selelem
989 * objects for selections.
991 t_selelem *
992 _gmx_sel_init_selection(char *name, t_selelem *sel, yyscan_t scanner)
994 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
995 t_selelem *root;
996 int rc;
998 if (sel->v.type != POS_VALUE)
1000 gmx_bug("each selection must evaluate to a position");
1001 /* FIXME: Better handling of this error */
1002 sfree(name);
1003 return NULL;
1006 root = _gmx_selelem_create(SEL_ROOT);
1007 root->child = sel;
1008 /* Assign the name (this is done here to free it automatically in the case
1009 * of an error below). */
1010 if (name)
1012 root->name = root->u.cgrp.name = name;
1014 /* Update the flags */
1015 rc = _gmx_selelem_update_flags(root);
1016 if (rc != 0)
1018 _gmx_selelem_free(root);
1019 return NULL;
1021 /* Initialize defaults for position keywords */
1022 init_pos_keyword_defaults(sel, sc, TRUE);
1024 /* If there is no name provided by the user, check whether the actual
1025 * selection given was from an external group, and if so, use the name
1026 * of the external group. */
1027 if (!root->name)
1029 t_selelem *child = root->child;
1030 while (child->type == SEL_MODIFIER)
1032 if (!child->child || child->child->type != SEL_SUBEXPRREF
1033 || !child->child->child)
1035 break;
1037 child = child->child->child;
1039 if (child->type == SEL_EXPRESSION
1040 && child->child && child->child->type == SEL_SUBEXPRREF
1041 && child->child->child
1042 && child->child->child->type == SEL_CONST
1043 && child->child->child->v.type == GROUP_VALUE)
1045 root->name = root->u.cgrp.name =
1046 strdup(child->child->child->u.cgrp.name);
1049 /* If there still is no name, use the selection string */
1050 if (!root->name)
1052 root->name = root->u.cgrp.name
1053 = strdup(_gmx_sel_lexer_pselstr(scanner));
1056 /* Print out some information if the parser is interactive */
1057 if (_gmx_sel_is_lexer_interactive(scanner))
1059 fprintf(stderr, "Selection '%s' parsed\n",
1060 _gmx_sel_lexer_pselstr(scanner));
1063 return root;
1068 * \param[in] name Name of the variable (should not be freed after this
1069 * function).
1070 * \param[in] expr The selection element that evaluates the variable.
1071 * \param scanner Scanner data structure.
1072 * \returns The created root selection element.
1074 * This function handles the creation of root \c t_selelem objects for
1075 * variable assignments. A \ref SEL_ROOT element and a \ref SEL_SUBEXPR
1076 * element are both created.
1078 t_selelem *
1079 _gmx_sel_assign_variable(char *name, t_selelem *expr, yyscan_t scanner)
1081 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
1082 const char *pselstr = _gmx_sel_lexer_pselstr(scanner);
1083 t_selelem *root = NULL;
1084 int rc;
1086 rc = _gmx_selelem_update_flags(expr);
1087 if (rc != 0)
1089 sfree(name);
1090 _gmx_selelem_free(expr);
1091 return NULL;
1093 /* Check if this is a constant non-group value */
1094 if (expr->type == SEL_CONST && expr->v.type != GROUP_VALUE)
1096 /* If so, just assign the constant value to the variable */
1097 if (!_gmx_sel_add_var_symbol(sc->symtab, name, expr))
1099 _gmx_selelem_free(expr);
1100 sfree(name);
1101 return NULL;
1103 _gmx_selelem_free(expr);
1104 sfree(name);
1105 goto finish;
1107 /* Check if we are assigning a variable to another variable */
1108 if (expr->type == SEL_SUBEXPRREF)
1110 /* If so, make a simple alias */
1111 if (!_gmx_sel_add_var_symbol(sc->symtab, name, expr->child))
1113 _gmx_selelem_free(expr);
1114 sfree(name);
1115 return NULL;
1117 _gmx_selelem_free(expr);
1118 sfree(name);
1119 goto finish;
1121 /* Create the root element */
1122 root = _gmx_selelem_create(SEL_ROOT);
1123 root->name = name;
1124 root->u.cgrp.name = name;
1125 /* Create the subexpression element */
1126 root->child = _gmx_selelem_create(SEL_SUBEXPR);
1127 _gmx_selelem_set_vtype(root->child, expr->v.type);
1128 root->child->name = name;
1129 root->child->child = expr;
1130 /* Update flags */
1131 rc = _gmx_selelem_update_flags(root);
1132 if (rc != 0)
1134 _gmx_selelem_free(root);
1135 return NULL;
1137 /* Add the variable to the symbol table */
1138 if (!_gmx_sel_add_var_symbol(sc->symtab, name, root->child))
1140 _gmx_selelem_free(root);
1141 return NULL;
1143 finish:
1144 srenew(sc->varstrs, sc->nvars + 1);
1145 sc->varstrs[sc->nvars] = strdup(pselstr);
1146 ++sc->nvars;
1147 if (_gmx_sel_is_lexer_interactive(scanner))
1149 fprintf(stderr, "Variable '%s' parsed\n", pselstr);
1151 return root;
1155 * \param sel Selection to append (can be NULL, in which
1156 * case nothing is done).
1157 * \param last Last selection, or NULL if not present or not known.
1158 * \param scanner Scanner data structure.
1159 * \returns The last selection after the append.
1161 * Appends \p sel after the last root element, and returns either \p sel
1162 * (if it was non-NULL) or the last element (if \p sel was NULL).
1164 t_selelem *
1165 _gmx_sel_append_selection(t_selelem *sel, t_selelem *last, yyscan_t scanner)
1167 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
1169 /* Append sel after last, or the last element of sc if last is NULL */
1170 if (last)
1172 last->next = sel;
1174 else
1176 if (sc->root)
1178 last = sc->root;
1179 while (last->next)
1181 last = last->next;
1183 last->next = sel;
1185 else
1187 sc->root = sel;
1190 /* Initialize a selection object if necessary */
1191 if (sel)
1193 last = sel;
1194 /* Add the new selection to the collection if it is not a variable. */
1195 if (sel->child->type != SEL_SUBEXPR)
1197 int i;
1199 sc->nr++;
1200 srenew(sc->sel, sc->nr);
1201 i = sc->nr - 1;
1202 snew(sc->sel[i], 1);
1203 sc->sel[i]->name = strdup(sel->name);
1204 sc->sel[i]->selstr = strdup(_gmx_sel_lexer_pselstr(scanner));
1206 if (sel->child->type == SEL_CONST)
1208 gmx_ana_pos_copy(&sc->sel[i]->p, sel->child->v.u.p, TRUE);
1209 sc->sel[i]->bDynamic = FALSE;
1211 else
1213 t_selelem *child;
1215 child = sel->child;
1216 child->flags &= ~SEL_ALLOCVAL;
1217 _gmx_selvalue_setstore(&child->v, &sc->sel[i]->p);
1218 /* We should also skip any modifiers to determine the dynamic
1219 * status. */
1220 while (child->type == SEL_MODIFIER)
1222 child = child->child;
1223 if (child->type == SEL_SUBEXPRREF)
1225 child = child->child->child;
1228 /* For variable references, we should skip the
1229 * SEL_SUBEXPRREF and SEL_SUBEXPR elements. */
1230 if (child->type == SEL_SUBEXPRREF)
1232 child = child->child->child;
1234 sc->sel[i]->bDynamic = (child->child->flags & SEL_DYNAMIC);
1236 /* The group will be set after compilation */
1237 sc->sel[i]->g = NULL;
1238 sc->sel[i]->selelem = sel;
1239 gmx_ana_selection_init_coverfrac(sc->sel[i], CFRAC_NONE);
1242 /* Clear the selection string now that we've saved it */
1243 _gmx_sel_lexer_clear_pselstr(scanner);
1244 return last;
1248 * \param[in] scanner Scanner data structure.
1249 * \returns TRUE if the parser should finish, FALSE if parsing should
1250 * continue.
1252 * This function is called always after _gmx_sel_append_selection() to
1253 * check whether a sufficient number of selections has already been provided.
1254 * This is used to terminate interactive parsers when the correct number of
1255 * selections has been provided.
1257 bool
1258 _gmx_sel_parser_should_finish(yyscan_t scanner)
1260 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
1261 return sc->nr == _gmx_sel_lexer_exp_selcount(scanner);
1265 * \param[in] scanner Scanner data structure.
1267 void
1268 _gmx_sel_handle_empty_cmd(yyscan_t scanner)
1270 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
1271 gmx_ana_indexgrps_t *grps = _gmx_sel_lexer_indexgrps(scanner);
1272 int i;
1274 if (!_gmx_sel_is_lexer_interactive(scanner))
1275 return;
1277 if (grps)
1279 fprintf(stderr, "Available index groups:\n");
1280 gmx_ana_indexgrps_print(_gmx_sel_lexer_indexgrps(scanner), 0);
1282 if (sc->nvars > 0 || sc->nr > 0)
1284 fprintf(stderr, "Currently provided selections:\n");
1285 for (i = 0; i < sc->nvars; ++i)
1287 fprintf(stderr, " %s\n", sc->varstrs[i]);
1289 for (i = 0; i < sc->nr; ++i)
1291 fprintf(stderr, " %2d. %s\n", i+1, sc->sel[i]->selstr);
1297 * \param[in] topic Topic for which help was requested, or NULL for general
1298 * help.
1299 * \param[in] scanner Scanner data structure.
1301 * \p topic is freed by this function.
1303 void
1304 _gmx_sel_handle_help_cmd(char *topic, yyscan_t scanner)
1306 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
1308 _gmx_sel_print_help(sc, topic);
1309 if (topic)
1311 sfree(topic);
1315 /*! \brief
1316 * Internal helper function used by gmx_ana_selcollection_parse_*() to do the actual work.
1318 * \param[in] maxnr Maximum number of selections to parse
1319 * (if -1, parse as many as provided by the user).
1320 * \param[in,out] scanner Scanner data structure.
1321 * \returns 0 on success, -1 on error.
1323 static int
1324 run_parser(int maxnr, yyscan_t scanner)
1326 gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
1327 bool bOk;
1328 int nr;
1330 nr = sc->nr;
1331 bOk = !_gmx_sel_yyparse(scanner);
1332 _gmx_sel_free_lexer(scanner);
1333 nr = sc->nr - nr;
1334 if (maxnr > 0 && nr != maxnr)
1336 return -1;
1338 return bOk ? 0 : -1;
1342 * \param[in,out] sc Selection collection to use for output.
1343 * \param[in] nr Number of selections to parse
1344 * (if -1, parse as many as provided by the user).
1345 * \param[in] grps External index groups (can be NULL).
1346 * \param[in] bInteractive Whether the parser should behave interactively.
1347 * \returns 0 on success, -1 on error.
1349 * The number of selections parsed can be accessed with
1350 * gmx_ana_selcollection_get_count() (note that if you call the parser
1351 * multiple times, this function returns the total count).
1354 gmx_ana_selcollection_parse_stdin(gmx_ana_selcollection_t *sc, int nr,
1355 gmx_ana_indexgrps_t *grps, bool bInteractive)
1357 yyscan_t scanner;
1358 int rc;
1360 rc = _gmx_sel_init_lexer(&scanner, sc, bInteractive, nr, grps);
1361 if (rc != 0)
1363 return rc;
1365 _gmx_sel_set_lex_input_file(scanner, stdin);
1366 return run_parser(nr, scanner);
1370 * \param[in,out] sc Selection collection to use for output.
1371 * \param[in] fnm Name of the file to parse selections from.
1372 * \param[in] grps External index groups (can be NULL).
1373 * \returns 0 on success, -1 on error.
1375 * The number of selections parsed can be accessed with
1376 * gmx_ana_selcollection_get_count() (note that if you call the parser
1377 * multiple times, this function returns the total count).
1380 gmx_ana_selcollection_parse_file(gmx_ana_selcollection_t *sc, const char *fnm,
1381 gmx_ana_indexgrps_t *grps)
1383 yyscan_t scanner;
1384 FILE *fp;
1385 int rc;
1387 rc = _gmx_sel_init_lexer(&scanner, sc, FALSE, -1, grps);
1388 if (rc != 0)
1390 return rc;
1392 fp = ffopen(fnm, "r");
1393 _gmx_sel_set_lex_input_file(scanner, fp);
1394 rc = run_parser(-1, scanner);
1395 fclose(fp);
1396 return rc;
1400 * \param[in,out] sc Selection collection to use for output.
1401 * \param[in] str String to parse selections from.
1402 * \param[in] grps External index groups (can be NULL).
1403 * \returns 0 on success, -1 on error.
1405 * The number of selections parsed can be accessed with
1406 * gmx_ana_selcollection_get_count() (note that if you call the parser
1407 * multiple times, this function returns the total count).
1410 gmx_ana_selcollection_parse_str(gmx_ana_selcollection_t *sc, const char *str,
1411 gmx_ana_indexgrps_t *grps)
1413 yyscan_t scanner;
1414 int rc;
1416 rc = _gmx_sel_init_lexer(&scanner, sc, FALSE, -1, grps);
1417 if (rc != 0)
1419 return rc;
1421 _gmx_sel_set_lex_input_str(scanner, str);
1422 return run_parser(-1, scanner);