Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / selection / compiler.cpp
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
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.
13
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.
18  *
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.
25  *
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.
28  *
29  * For more info, check our website at http://www.gromacs.org
30  */
31 /*! \internal \file
32  * \brief Selection compilation and optimization.
33  *
34  * \todo
35  * Better error handling and memory management in error situations.
36  * At least, the main compilation function leaves the selection collection in
37  * a bad state if an error occurs.
38  *
39  * \todo
40  * The memory usage could still be optimized.
41  * Use of memory pooling could still be extended, and a lot of redundant
42  * gmin/gmax data could be eliminated for complex arithmetic expressions.
43  *
44  * \author Teemu Murtola <teemu.murtola@cbr.su.se>
45  * \ingroup module_selection
46  */
47 /*! \internal
48  * \page page_module_selection_compiler Selection compilation
49  *
50  * The compiler takes the selection element tree from the selection parser
51  * (see \ref page_module_selection_parser) as input.
52  * The selection parser is quite independent of selection evaluation details,
53  * and the compiler processes the tree to conform to what the evaluation
54  * functions expect.
55  * For better control and optimization possibilities, the compilation is
56  * done on all selections simultaneously.
57  * Hence, all the selections should be parsed before the compiler can be
58  * called.
59  *
60  * The compiler initializes all fields in gmx::SelectionTreeElement not
61  * initialized by the parser: gmx::SelectionTreeElement::v (some fields have
62  * already been initialized by the parser),
63  * gmx::SelectionTreeElement::evaluate, and gmx::SelectionTreeElement::u
64  * (again, some elements have been initialized in the parser).
65  * The gmx::SelectionTreeElement::cdata field is used during the compilation to store
66  * internal data, but the data is freed when the compiler returns.
67  *
68  * In addition to initializing the elements, the compiler reorganizes the tree
69  * to simplify and optimize evaluation. The compiler also evaluates the static
70  * parts of the selection: in the end of the compilation, static parts have
71  * been replaced by the result of the evaluation.
72  *
73  * The compiler is invoked using gmx::SelectionCompiler.
74  * The gmx::SelectionCompiler::compile() method does the compilation in several
75  * passes over the gmx::SelectionTreeElement tree.
76  *  -# Defaults are set for the position type and flags of position calculation
77  *     methods that were not explicitly specified in the user input.
78  *  -# Subexpressions are extracted: a separate root is created for each
79  *     subexpression, and placed before the expression is first used.
80  *     Currently, only variables and expressions used to evaluate parameter
81  *     values are extracted, but common subexpression could also be detected
82  *     here.
83  *  -# A second pass (in fact, multiple passes because of interdependencies)
84  *     with simple reordering and initialization is done:
85  *    -# Boolean expressions are combined such that one element can evaluate,
86  *       e.g., "A and B and C". The subexpressions in boolean expression are
87  *       reordered such that static expressions come first without otherwise
88  *       altering the relative order of the expressions.
89  *    -# The compiler data structure is allocated for each element, and
90  *       the fields are initialized, with the exception of the contents of
91  *       \c gmax and \c gmin fields.  This is the part that needs multiple
92  *       passes, because some flags are set recursively based on which elements
93  *       refer to an element, and these flags need to be set to initialize
94  *       other fields.
95  *    -# The gmx::SelectionTreeElement::evaluate field is set to the correct
96  *       evaluation function from evaluate.h.
97  *    .
98  *  -# The evaluation function of all elements is replaced with the
99  *     analyze_static() function to be able to initialize the element before
100  *     the actual evaluation function is called.
101  *     The evaluation machinery is then called to initialize the whole tree,
102  *     while simultaneously evaluating the static expressions.
103  *     During the evaluation, track is kept of the smallest and largest
104  *     possible selections, and these are stored in the internal compiler
105  *     data structure for each element.
106  *     To be able to do this for all possible values of dynamical expressions,
107  *     special care needs to be taken with boolean expressions because they
108  *     are short-circuiting. This is done through the
109  *     \c SEL_CDATA_EVALMAX flag, which makes dynamic child expressions
110  *     of \c BOOL_OR expressions evaluate to empty groups, while subexpressions
111  *     of \c BOOL_AND are evaluated to largest possible groups.
112  *     Memory is also allocated to store the results of the evaluation.
113  *     For each element, analyze_static() calls the actual evaluation function
114  *     after the element has been properly initialized.
115  *  -# Another evaluation pass is done over subexpressions with more than
116  *     one reference to them. These cannot be completely processed during the
117  *     first pass, because it is not known whether later references require
118  *     additional evaluation of static expressions.
119  *  -# Unused subexpressions are removed. For efficiency reasons (and to avoid
120  *     some checks), this is actually done several times already earlier in
121  *     the compilation process.
122  *  -# Most of the processing is now done, and the next pass simply sets the
123  *     evaluation group of root elements to the largest selection as determined
124  *     in pass 4.  For root elements of subexpressions that should not be
125  *     evaluated before they are referred to, the evaluation group/function is
126  *     cleared.  At the same time, position calculation data is initialized for
127  *     for selection method elements that require it.  Compiler data is also
128  *     freed as it is no longer needed.
129  *  -# A final pass initializes the total masses and charges in the
130  *     \c gmx_ana_selection_t data structures.
131  *
132  * The actual evaluation of the selection is described in the documentation
133  * of the functions in evaluate.h.
134  *
135  * \todo
136  * Some combinations of method parameter flags are not yet properly treated by
137  * the compiler or the evaluation functions in evaluate.cpp. All the ones used by
138  * currently implemented methods should work, but new combinations might not.
139  *
140  *
141  * \section selcompiler_tree Element tree after compilation
142  *
143  * After the compilation, the selection element tree is suitable for
144  * gmx_ana_selcollection_evaluate().
145  * Enough memory has been allocated for gmx::SelectionTreeElement::v
146  * (and gmx::SelectionTreeElement::cgrp for \ref SEL_SUBEXPR elements) to allow
147  * the selection to be evaluated without allocating any memory.
148  *
149  *
150  * \subsection selcompiler_tree_root Root elements
151  *
152  * The top level of the tree consists of a chain of \ref SEL_ROOT elements.
153  * These are used for two purposes:
154  *  -# A selection that should be evaluated.
155  *     These elements appear in the same order as the selections in the input.
156  *     For these elements, gmx::SelectionTreeElement::v has been set to the
157  *     maximum possible group that the selection can evaluate to (only for
158  *     dynamic selections), and gmx::SelectionTreeElement::cgrp has been set to
159  *     use a NULL group for evaluation.
160  *  -# A subexpression that appears in one or more selections.
161  *     Each selection that gives a value for a method parameter is a
162  *     potential subexpression, as is any variable value.
163  *     Only subexpressions that require evaluation for each frame are left
164  *     after the selection is compiled.
165  *     Each subexpression appears in the chain before any references to it.
166  *     For these elements, gmx::SelectionTreeElement::cgrp has been set to the
167  *     group that should be used to evaluate the subexpression.
168  *     If gmx::SelectionTreeElement::cgrp is empty, the total evaluation group
169  *     is not known in advance or it is more efficient to evaluate the
170  *     subexpression only when it is referenced.  If this is the case,
171  *     gmx::SelectionTreeElement::evaluate is also NULL.
172  *
173  * The children of the \ref SEL_ROOT elements can be used to distinguish
174  * the two types of root elements from each other; the rules are the same
175  * as for the parsed tree (see \ref selparser_tree_root).
176  * Subexpressions are treated as if they had been provided through variables.
177  *
178  * Selection names are stored as after parsing (see \ref selparser_tree_root).
179  *
180  *
181  * \subsection selcompiler_tree_const Constant elements
182  *
183  * All (sub)selections that do not require particle positions have been
184  * replaced with \ref SEL_CONST elements.
185  * Constant elements from the parser are also retained if present in
186  * dynamic parts of the selections.
187  * Several constant elements with a NULL gmx::SelectionTreeElement::evaluate
188  * are left for debugging purposes; of these, only the ones for \ref BOOL_OR
189  * expressions are used during evaluation.
190  *
191  * The value is stored in gmx::SelectionTreeElement::v, and for group values
192  * with an evaluation function set, also in gmx::SelectionTreeElement::cgrp.
193  * For \ref GROUP_VALUE elements, unnecessary atoms (i.e., atoms that
194  * could never be selected) have been removed from the value.
195  *
196  * \ref SEL_CONST elements have no children.
197  *
198  *
199  * \subsection selcompiler_tree_method Method evaluation elements
200  *
201  * All selection methods that need to be evaluated dynamically are described
202  * by a \ref SEL_EXPRESSION element. The gmx::SelectionTreeElement::method and
203  * gmx::SelectionTreeElement::mdata fields have already been initialized by the parser,
204  * and the compiler only calls the initialization functions in the method
205  * data structure to do some additional initialization of these fields at
206  * appropriate points. If the gmx::SelectionTreeElement::pc data field has been
207  * created by the parser, the compiler initializes the data structure properly
208  * once the required positions are known. If the gmx::SelectionTreeElement::pc
209  * field is NULL after the parser, but the method provides only
210  * sel_updatefunc_pos(), an appropriate position calculation data structure is
211  * created.  If gmx::SelectionTreeElement::pc is not NULL,
212  * gmx::SelectionTreeElement::pos is also initialized to hold the positions
213  * calculated.
214  *
215  * Children of these elements are of type \ref SEL_SUBEXPRREF, and describe
216  * parameter values that need to be evaluated for each frame. See the next
217  * section for more details.
218  * \ref SEL_CONST children can also appear, and stand for parameters that get
219  * their value from a static expression. These elements are present only for
220  * debugging purposes: they always have a NULL evaluation function.
221  *
222  *
223  * \subsection selcompiler_tree_subexpr Subexpression elements
224  *
225  * As described in \ref selcompiler_tree_root, subexpressions are created
226  * for each variable and each expression that gives a value to a selection
227  * method parameter. As the only child of the \ref SEL_ROOT element,
228  * these elements have a \ref SEL_SUBEXPR element. The \ref SEL_SUBEXPR
229  * element has a single child, which evaluates the actual expression.
230  * After compilation, only subexpressions that require particle positions
231  * for evaluation are left.
232  * For non-variable subexpression, automatic names have been generated to
233  * help in debugging.
234  *
235  * For \ref SEL_SUBEXPR elements, memory has been allocated for
236  * gmx::SelectionTreeElement::cgrp to store the group for which the expression
237  * has been evaluated during the current frame.  This is only done if full
238  * subexpression evaluation by _gmx_sel_evaluate_subexpr() is needed; the other
239  * evaluation functions do not require this memory.
240  *
241  * \ref SEL_SUBEXPRREF elements are used to describe references to
242  * subexpressions. They have always a single child, which is the
243  * \ref SEL_SUBEXPR element being referenced.
244  *
245  * If a subexpression is used only once, the evaluation has been optimized by
246  * setting the child of the \ref SEL_SUBEXPR element to evaluate the value of
247  * \ref SEL_SUBEXPRREF directly (in the case of memory pooling, this is managed
248  * by the evaluation functions).  In such cases, the evaluation routines for the
249  * \ref SEL_SUBEXPRREF and \ref SEL_SUBEXPR elements only propagate some status
250  * information, but do not unnecessarily copy the values.
251  *
252  *
253  * \subsection selcompiler_tree_bool Boolean elements
254  *
255  * \ref SEL_BOOLEAN elements have been merged such that one element
256  * may carry out evaluation of more than one operation of the same type.
257  * The static parts of the expressions have been evaluated, and are placed
258  * in the first child. These are followed by the dynamic expressions, in the
259  * order provided by the user.
260  *
261  *
262  * \subsection selcompiler_tree_arith Arithmetic elements
263  *
264  * Constant and static expressions in \ref SEL_ARITHMETIC elements have been
265  * calculated.
266  * Currently, no other processing is done.
267  */
268 #include "compiler.h"
269
270 #ifdef HAVE_CONFIG_H
271 #include <config.h>
272 #endif
273
274 #include <algorithm>
275
276 #include <math.h>
277 #include <stdarg.h>
278
279 #include "smalloc.h"
280 #include "string2.h"
281 #include "vec.h"
282
283 #include "gromacs/selection/indexutil.h"
284 #include "gromacs/selection/poscalc.h"
285 #include "gromacs/selection/selection.h"
286 #include "gromacs/selection/selmethod.h"
287 #include "gromacs/utility/exceptions.h"
288 #include "gromacs/utility/stringutil.h"
289
290 #include "evaluate.h"
291 #include "keywords.h"
292 #include "mempool.h"
293 #include "selectioncollection-impl.h"
294 #include "selelem.h"
295
296 using std::min;
297 using gmx::SelectionTreeElement;
298 using gmx::SelectionTreeElementPointer;
299
300 /*! \internal \brief
301  * Compiler flags.
302  */
303 enum
304 {
305     /*! \brief
306      * Whether a subexpression needs to evaluated for all atoms.
307      *
308      * This flag is set for \ref SEL_SUBEXPR elements that are used to
309      * evaluate non-atom-valued selection method parameters, as well as
310      * those that are used directly as values of selections.
311      */
312     SEL_CDATA_FULLEVAL    =  1,
313     /*! \brief
314      * Whether the whole subexpression should be treated as static.
315      *
316      * This flag is always false if \ref SEL_DYNAMIC is set for the element,
317      * but it is also false for static elements within common subexpressions.
318      */
319     SEL_CDATA_STATIC      =  2,
320     /** Whether the subexpression will always be evaluated in the same group. */
321     SEL_CDATA_STATICEVAL  =  4,
322     /** Whether the compiler evaluation routine should return the maximal selection. */
323     SEL_CDATA_EVALMAX     =  8,
324     /** Whether memory has been allocated for \p gmin and \p gmax. */
325     SEL_CDATA_MINMAXALLOC = 16,
326     /** Whether to update \p gmin and \p gmax in static analysis. */
327     SEL_CDATA_DOMINMAX      = 128,
328     /** Whether the subexpression uses simple pass evaluation functions. */
329     SEL_CDATA_SIMPLESUBEXPR = 32,
330     /** Whether this expression is a part of a common subexpression. */
331     SEL_CDATA_COMMONSUBEXPR = 64
332 };
333
334 /*! \internal \brief
335  * Internal data structure used by the compiler.
336  */
337 typedef struct t_compiler_data
338 {
339     /** The real evaluation method. */
340     gmx::sel_evalfunc evaluate;
341     /** Number of references to a \ref SEL_SUBEXPR element. */
342     int              refcount;
343     /** Flags for specifying how to treat this element during compilation. */
344     int              flags;
345     /** Smallest selection that can be selected by the subexpression. */
346     gmx_ana_index_t *gmin;
347     /** Largest selection that can be selected by the subexpression. */
348     gmx_ana_index_t *gmax;
349 } t_compiler_data;
350
351
352 /********************************************************************
353  * COMPILER UTILITY FUNCTIONS
354  ********************************************************************/
355
356 /*! \brief
357  * Helper method for printing out debug information about a min/max group.
358  */
359 static void
360 print_group_info(FILE *fp, const char *name,
361                  const SelectionTreeElement &sel,
362                  gmx_ana_index_t *g)
363 {
364     fprintf(fp, " %s=", name);
365     if (!g)
366     {
367         fprintf(fp, "(null)");
368     }
369     else if (sel.cdata->flags & SEL_CDATA_MINMAXALLOC)
370     {
371         fprintf(fp, "(%d atoms, %p)", g->isize, (void*)g);
372     }
373     else if (sel.v.type == GROUP_VALUE && g == sel.v.u.g)
374     {
375         fprintf(fp, "(static, %p)", (void*)g);
376     }
377     else
378     {
379         fprintf(fp, "%p", (void*)g);
380     }
381 }
382
383 /*!
384  * \param[in] fp      File handle to receive the output.
385  * \param[in] sel     Selection element to print.
386  * \param[in] level   Indentation level, starting from zero.
387  */
388 void
389 _gmx_selelem_print_compiler_info(FILE *fp, const SelectionTreeElement &sel,
390                                  int level)
391 {
392     if (!sel.cdata)
393     {
394         return;
395     }
396     fprintf(fp, "%*c cdata: flg=", level*2+1, ' ');
397     if (sel.cdata->flags & SEL_CDATA_FULLEVAL)
398     {
399         fprintf(fp, "F");
400     }
401     if (!(sel.cdata->flags & SEL_CDATA_STATIC))
402     {
403         fprintf(fp, "D");
404     }
405     if (sel.cdata->flags & SEL_CDATA_STATICEVAL)
406     {
407         fprintf(fp, "S");
408     }
409     if (sel.cdata->flags & SEL_CDATA_EVALMAX)
410     {
411         fprintf(fp, "M");
412     }
413     if (sel.cdata->flags & SEL_CDATA_MINMAXALLOC)
414     {
415         fprintf(fp, "A");
416     }
417     if (sel.cdata->flags & SEL_CDATA_SIMPLESUBEXPR)
418     {
419         fprintf(fp, "Ss");
420     }
421     if (sel.cdata->flags & SEL_CDATA_COMMONSUBEXPR)
422     {
423         fprintf(fp, "Sc");
424     }
425     if (!sel.cdata->flags)
426     {
427         fprintf(fp, "0");
428     }
429     if (sel.cdata->refcount > 0)
430     {
431         fprintf(fp, " refc=%d", sel.cdata->refcount);
432     }
433     fprintf(fp, " eval=");
434     _gmx_sel_print_evalfunc_name(fp, sel.cdata->evaluate);
435     print_group_info(fp, "gmin", sel, sel.cdata->gmin);
436     print_group_info(fp, "gmax", sel, sel.cdata->gmax);
437     fprintf(fp, "\n");
438 }
439
440 namespace gmx
441 {
442
443 void SelectionTreeElement::freeCompilerData()
444 {
445     if (cdata)
446     {
447         evaluate = cdata->evaluate;
448         if (cdata->flags & SEL_CDATA_MINMAXALLOC)
449         {
450             cdata->gmin->name = NULL;
451             cdata->gmax->name = NULL;
452             gmx_ana_index_deinit(cdata->gmin);
453             gmx_ana_index_deinit(cdata->gmax);
454             sfree(cdata->gmin);
455             sfree(cdata->gmax);
456         }
457         sfree(cdata);
458     }
459     cdata = NULL;
460 }
461
462 } // namespace gmx
463
464 /*! \brief
465  * Allocates memory for storing the evaluated value of a selection element.
466  *
467  * \param     sel   Selection element to initialize
468  * \param[in] isize Maximum evaluation group size.
469  * \param[in] bChildEval true if children have already been processed.
470  *
471  * If called more than once, memory is (re)allocated to ensure that the
472  * maximum of the \p isize values can be stored.
473  */
474 static void
475 alloc_selection_data(const SelectionTreeElementPointer &sel,
476                      int isize, bool bChildEval)
477 {
478     int        nalloc;
479
480     if (sel->mempool)
481     {
482         return;
483     }
484     /* Find out the number of elements to allocate */
485     if (sel->flags & SEL_SINGLEVAL)
486     {
487         nalloc = 1;
488     }
489     else if (sel->flags & SEL_ATOMVAL)
490     {
491         nalloc = isize;
492     }
493     else /* sel->flags should contain SEL_VARNUMVAL */
494     {
495         if (!bChildEval)
496         {
497             return;
498         }
499         SelectionTreeElementPointer child = sel;
500         if (sel->type == SEL_SUBEXPRREF)
501         {
502             GMX_RELEASE_ASSERT(sel->child && sel->child->type == SEL_SUBEXPR,
503                 "Subexpression expected for subexpression reference");
504             child = sel->child->child;
505             GMX_RELEASE_ASSERT(child,
506                 "Subexpression elements should always have a child element");
507         }
508         nalloc = (sel->v.type == POS_VALUE) ? child->v.u.p->nr : child->v.nr;
509     }
510     /* For positions, we actually want to allocate just a single structure
511      * for nalloc positions. */
512     if (sel->v.type == POS_VALUE)
513     {
514         isize  = nalloc;
515         nalloc = 1;
516     }
517     /* Allocate memory for sel->v.u if needed */
518     if (sel->flags & SEL_ALLOCVAL)
519     {
520         _gmx_selvalue_reserve(&sel->v, nalloc);
521     }
522     /* Reserve memory inside group and position structures if
523      * SEL_ALLOCDATA is set. */
524     if (sel->flags & SEL_ALLOCDATA)
525     {
526         if (sel->v.type == GROUP_VALUE)
527         {
528             gmx_ana_index_reserve(sel->v.u.g, isize);
529         }
530         else if (sel->v.type == POS_VALUE)
531         {
532             gmx_ana_pos_reserve(sel->v.u.p, isize, 0);
533         }
534     }
535 }
536
537 /*! \brief
538  * Replace the evaluation function of each element in the subtree.
539  *
540  * \param     sel  Root of the selection subtree to process.
541  * \param[in] eval The new evaluation function.
542  */
543 static void
544 set_evaluation_function(const SelectionTreeElementPointer &sel,
545                         gmx::sel_evalfunc eval)
546 {
547     sel->evaluate = eval;
548     if (sel->type != SEL_SUBEXPRREF)
549     {
550         SelectionTreeElementPointer child = sel->child;
551         while (child)
552         {
553             set_evaluation_function(child, eval);
554             child = child->next;
555         }
556     }
557 }
558
559
560 /********************************************************************
561  * POSITION KEYWORD DEFAULT INITIALIZATION
562  ********************************************************************/
563
564 /*! \brief
565  * Initializes default values for position keyword evaluation.
566  *
567  * \param[in,out] root       Root of the element tree to initialize.
568  * \param[in]     spost      Default output position type.
569  * \param[in]     rpost      Default reference position type.
570  * \param[in]     sel        Selection that the element evaluates the positions
571  *      for, or NULL if the element is an internal element.
572  */
573 static void
574 init_pos_keyword_defaults(SelectionTreeElement *root,
575                           const char *spost, const char *rpost,
576                           const gmx::internal::SelectionData *sel)
577 {
578     /* Selections use largest static group by default, while
579      * reference positions use the whole residue/molecule. */
580     if (root->type == SEL_EXPRESSION)
581     {
582         bool bSelection = (sel != NULL);
583         int flags = bSelection ? POS_COMPLMAX : POS_COMPLWHOLE;
584         if (bSelection)
585         {
586             if (sel->hasFlag(gmx::efSelection_DynamicMask))
587             {
588                 flags |= POS_MASKONLY;
589             }
590             if (sel->hasFlag(gmx::efSelection_EvaluateVelocities))
591             {
592                 flags |= POS_VELOCITIES;
593             }
594             if (sel->hasFlag(gmx::efSelection_EvaluateForces))
595             {
596                 flags |= POS_FORCES;
597             }
598         }
599         _gmx_selelem_set_kwpos_type(root, bSelection ? spost : rpost);
600         _gmx_selelem_set_kwpos_flags(root, flags);
601     }
602     /* Change the defaults once we are no longer processing modifiers */
603     if (root->type != SEL_ROOT && root->type != SEL_MODIFIER
604         && root->type != SEL_SUBEXPRREF && root->type != SEL_SUBEXPR)
605     {
606         sel = NULL;
607     }
608     /* Recurse into children */
609     SelectionTreeElementPointer child = root->child;
610     while (child)
611     {
612         init_pos_keyword_defaults(child.get(), spost, rpost, sel);
613         child = child->next;
614     }
615 }
616
617
618 /********************************************************************
619  * SUBEXPRESSION PROCESSING
620  ********************************************************************/
621
622 /*! \brief
623  * Reverses the chain of selection elements starting at \p root.
624  *
625  * \param   root First selection in the whole selection chain.
626  * \returns The new first element for the chain.
627  */
628 static SelectionTreeElementPointer
629 reverse_selelem_chain(const SelectionTreeElementPointer &root)
630 {
631     SelectionTreeElementPointer prev;
632     SelectionTreeElementPointer item = root;
633     while (item)
634     {
635         SelectionTreeElementPointer next = item->next;
636         item->next = prev;
637         prev = item;
638         item = next;
639     }
640     return prev;
641 }
642
643 /*! \brief
644  * Removes subexpressions that don't have any references.
645  *
646  * \param     root First selection in the whole selection chain.
647  * \returns   The new first element for the chain.
648  *
649  * The elements are processed in reverse order to correctly detect
650  * subexpressions only referred to by other subexpressions.
651  */
652 static SelectionTreeElementPointer
653 remove_unused_subexpressions(SelectionTreeElementPointer root)
654 {
655     if (!root)
656     {
657         return SelectionTreeElementPointer();
658     }
659     root = reverse_selelem_chain(root);
660     while (root->child->type == SEL_SUBEXPR && root->child.unique())
661     {
662         // Frees the root element.
663         root = root->next;
664     }
665     SelectionTreeElementPointer prev = root;
666     SelectionTreeElementPointer item = root->next;
667     while (item)
668     {
669         SelectionTreeElementPointer next = item->next;
670         if (item->child->type == SEL_SUBEXPR && item->child.unique())
671         {
672             // Frees the current item when it goes out of scope.
673             prev->next = next;
674         }
675         else
676         {
677             prev = item;
678         }
679         item = next;
680     }
681     return reverse_selelem_chain(root);
682 }
683
684 /*! \brief
685  * Creates a name with a running number for a subexpression.
686  *
687  * \param[in,out] sel The subexpression to be named.
688  * \param[in]     i   Running number for the subexpression.
689  *
690  * The name of the selection becomes "SubExpr N", where N is \p i;
691  * Memory is allocated for the name and the name is stored both as the
692  * name of the subexpression element and as
693  * gmx::SelectionTreeElement::u::cgrp::name; the latter is freed by
694  * _gmx_selelem_free().
695  */
696 static void
697 create_subexpression_name(const SelectionTreeElementPointer &sel, int i)
698 {
699     std::string name(gmx::formatString("SubExpr %d", i));
700     sel->setName(name);
701     sel->u.cgrp.name = strdup(name.c_str());
702 }
703
704 /*! \brief
705  * Processes and extracts subexpressions from a given selection subtree.
706  *
707  * \param   sel      Root of the subtree to process.
708  * \param   subexprn Pointer to a subexpression counter.
709  * \returns Pointer to a chain of subselections, or NULL if none were found.
710  *
711  * This function finds recursively all \ref SEL_SUBEXPRREF elements below
712  * the given root element and ensures that their children are within
713  * \ref SEL_SUBEXPR elements. It also creates a chain of \ref SEL_ROOT elements
714  * that contain the subexpression as their children and returns the first
715  * of these root elements.
716  */
717 static SelectionTreeElementPointer
718 extract_item_subselections(const SelectionTreeElementPointer &sel,
719                            int *subexprn)
720 {
721     SelectionTreeElementPointer root;
722     SelectionTreeElementPointer subexpr;
723     SelectionTreeElementPointer child = sel->child;
724
725     while (child)
726     {
727         if (!root)
728         {
729             root = subexpr = extract_item_subselections(child, subexprn);
730         }
731         else
732         {
733             subexpr->next = extract_item_subselections(child, subexprn);
734         }
735         while (subexpr && subexpr->next)
736         {
737             subexpr = subexpr->next;
738         }
739         /* The latter check excludes variable references. */
740         if (child->type == SEL_SUBEXPRREF && child->child->type != SEL_SUBEXPR)
741         {
742             /* Create the root element for the subexpression */
743             if (!root)
744             {
745                 root.reset(new SelectionTreeElement(SEL_ROOT));
746                 subexpr = root;
747             }
748             else
749             {
750                 subexpr->next.reset(new SelectionTreeElement(SEL_ROOT));
751                 subexpr = subexpr->next;
752             }
753             /* Create the subexpression element and
754              * move the actual subexpression under the created element. */
755             subexpr->child.reset(new SelectionTreeElement(SEL_SUBEXPR));
756             _gmx_selelem_set_vtype(subexpr->child, child->v.type);
757             subexpr->child->child = child->child;
758             child->child          = subexpr->child;
759             create_subexpression_name(subexpr->child, ++*subexprn);
760             /* Set the flags for the created elements */
761             subexpr->flags          |= (child->flags & SEL_VALFLAGMASK);
762             subexpr->child->flags   |= (child->flags & SEL_VALFLAGMASK);
763         }
764         if (child->type == SEL_SUBEXPRREF)
765         {
766             child->setName(child->child->name());
767         }
768         child = child->next;
769     }
770
771     return root;
772 }
773
774 /*! \brief
775  * Extracts subexpressions of the selection chain.
776  * 
777  * \param   sel First selection in the whole selection chain.
778  * \returns The new first element for the chain.
779  *
780  * Finds all the subexpressions (and their subexpressions) in the
781  * selection chain starting from \p sel and creates \ref SEL_SUBEXPR
782  * elements for them.
783  * \ref SEL_ROOT elements are also created for each subexpression
784  * and inserted into the selection chain before the expressions that
785  * refer to them.
786  */
787 static SelectionTreeElementPointer
788 extract_subexpressions(SelectionTreeElementPointer sel)
789 {
790     SelectionTreeElementPointer root;
791     SelectionTreeElementPointer next = sel;
792     int subexprn = 0;
793     while (next)
794     {
795         SelectionTreeElementPointer item
796             = extract_item_subselections(next, &subexprn);
797         if (item)
798         {
799             if (!root)
800             {
801                 root = item;
802             }
803             else
804             {
805                 sel->next = item;
806             }
807             while (item->next)
808             {
809                 item = item->next;
810             }
811             item->next = next;
812         }
813         else if (!root)
814         {
815             root = next;
816         }
817         sel = next;
818         next = next->next;
819     }
820     return root;
821 }
822
823
824 /********************************************************************
825  * BOOLEAN OPERATION REORDERING
826  ********************************************************************/
827
828 /*! \brief
829  * Removes redundant boolean selection elements.
830  *
831  * \param  sel Root of the selection subtree to optimize.
832  *
833  * This function merges similar boolean operations (e.g., (A or B) or C becomes
834  * a single OR operation with three operands).
835  */
836 static void
837 optimize_boolean_expressions(const SelectionTreeElementPointer &sel)
838 {
839     /* Do recursively for children */
840     if (sel->type != SEL_SUBEXPRREF)
841     {
842         SelectionTreeElementPointer prev;
843         SelectionTreeElementPointer child = sel->child;
844         while (child)
845         {
846             optimize_boolean_expressions(child);
847             /* Remove double negations */
848             if (child->type == SEL_BOOLEAN && child->u.boolt == BOOL_NOT
849                 && child->child->type == SEL_BOOLEAN && child->child->u.boolt == BOOL_NOT)
850             {
851                 /* Move the doubly negated expression up two levels */
852                 if (!prev)
853                 {
854                     sel->child = child->child->child;
855                     prev       = sel->child;
856                 }
857                 else
858                 {
859                     prev->next = child->child->child;
860                     prev       = prev->next;
861                 }
862                 child->child->child->next = child->next;
863                 // Discards the two negations.
864                 child = prev;
865             }
866             prev  = child;
867             child = child->next;
868         }
869     }
870     if (sel->type != SEL_BOOLEAN || sel->u.boolt == BOOL_NOT)
871     {
872         return;
873     }
874     /* Merge subsequent binary operations */
875     SelectionTreeElementPointer prev;
876     SelectionTreeElementPointer child = sel->child;
877     while (child)
878     {
879         if (child->type == SEL_BOOLEAN && child->u.boolt == sel->u.boolt)
880         {
881             if (!prev)
882             {
883                 sel->child = child->child;
884                 prev       = sel->child;
885             }
886             else
887             {
888                 prev->next = child->child;
889             }
890             while (prev->next)
891             {
892                 prev = prev->next;
893             }
894             prev->next = child->next;
895             // Discards the old child.
896             child = prev->next;
897         }
898         else
899         {
900             prev = child;
901             child = child->next;
902         }
903     }
904 }
905
906 /*! \brief
907  * Reorders children of boolean expressions such that static selections
908  * come first.
909  *
910  * \param  sel Root of the selection subtree to reorder.
911  *
912  * The relative order of static expressions does not change.
913  * The same is true for the dynamic expressions.
914  */
915 static void
916 reorder_boolean_static_children(const SelectionTreeElementPointer &sel)
917 {
918     /* Do recursively for children */
919     if (sel->type != SEL_SUBEXPRREF)
920     {
921         SelectionTreeElementPointer child = sel->child;
922         while (child)
923         {
924             reorder_boolean_static_children(child);
925             child = child->next;
926         }
927     }
928
929     /* Reorder boolean expressions such that static selections come first */
930     if (sel->type == SEL_BOOLEAN && (sel->flags & SEL_DYNAMIC))
931     {
932         // Add a dummy head element that precedes the first child.
933         SelectionTreeElementPointer dummy(
934                 new SelectionTreeElement(SEL_BOOLEAN));
935         dummy->next = sel->child;
936         SelectionTreeElementPointer prev  = dummy;
937         SelectionTreeElementPointer child = dummy;
938         while (child->next)
939         {
940             /* child is the last handled static expression */
941             /* prev is the last handled non-static expression */
942             SelectionTreeElementPointer next = prev->next;
943             while (next && (next->flags & SEL_DYNAMIC))
944             {
945                 prev = next;
946                 next = next->next;
947             }
948             /* next is now the first static expression after child */
949             if (!next)
950             {
951                 break;
952             }
953             /* Reorder such that next comes after child */
954             if (prev != child)
955             {
956                 prev->next  = next->next;
957                 next->next  = child->next;
958                 child->next = next;
959             }
960             else
961             {
962                 prev = prev->next;
963             }
964             /* Advance child by one */
965             child = next;
966         }
967
968         sel->child = dummy->next;
969     }
970 }
971
972
973 /********************************************************************
974  * ARITHMETIC EXPRESSION PROCESSING
975  ********************************************************************/
976
977 /*! \brief
978  * Processes arithmetic expressions to simplify and speed up evaluation.
979  *
980  * \param  sel Root of the selection subtree to process.
981  *
982  * Currently, this function only converts integer constants to reals
983  * within arithmetic expressions.
984  */
985 static void
986 optimize_arithmetic_expressions(const SelectionTreeElementPointer &sel)
987 {
988     /* Do recursively for children. */
989     if (sel->type != SEL_SUBEXPRREF)
990     {
991         SelectionTreeElementPointer child = sel->child;
992         while (child)
993         {
994             optimize_arithmetic_expressions(child);
995             child = child->next;
996         }
997     }
998
999     if (sel->type != SEL_ARITHMETIC)
1000     {
1001         return;
1002     }
1003
1004     /* Convert integer constants to reals. */
1005     SelectionTreeElementPointer child = sel->child;
1006     while (child)
1007     {
1008         if (child->v.type == INT_VALUE)
1009         {
1010             real  *r;
1011
1012             if (child->type != SEL_CONST)
1013             {
1014                 GMX_THROW(gmx::InconsistentInputError("Non-constant integer expressions not implemented in arithmetic evaluation"));
1015             }
1016             snew(r, 1);
1017             r[0] = child->v.u.i[0];
1018             sfree(child->v.u.i);
1019             child->v.u.r = r;
1020             child->v.type = REAL_VALUE;
1021         }
1022         else if (child->v.type != REAL_VALUE)
1023         {
1024             GMX_THROW(gmx::InternalError("Non-numerical value in arithmetic expression"));
1025         }
1026         child = child->next;
1027     }
1028 }
1029
1030
1031 /********************************************************************
1032  * EVALUATION PREPARATION COMPILER
1033  ********************************************************************/
1034
1035 /*! \brief
1036  * Sets the evaluation functions for the selection (sub)tree.
1037  *
1038  * \param[in,out] sel Root of the selection subtree to process.
1039  *
1040  * This function sets the evaluation function
1041  * (gmx::SelectionTreeElement::evaluate) for the selection elements.
1042  */
1043 static void
1044 init_item_evalfunc(const SelectionTreeElementPointer &sel)
1045 {
1046     /* Process children. */
1047     if (sel->type != SEL_SUBEXPRREF)
1048     {
1049         SelectionTreeElementPointer child = sel->child;
1050         while (child)
1051         {
1052             init_item_evalfunc(child);
1053             child = child->next;
1054         }
1055     }
1056
1057     /* Set the evaluation function */
1058     switch (sel->type)
1059     {
1060         case SEL_CONST:
1061             if (sel->v.type == GROUP_VALUE)
1062             {
1063                 sel->evaluate = &_gmx_sel_evaluate_static;
1064             }
1065             break;
1066
1067         case SEL_EXPRESSION:
1068             if (!(sel->flags & SEL_DYNAMIC) && sel->u.expr.method
1069                 && sel->u.expr.method->init_frame)
1070             {
1071                 sel->flags |= SEL_INITFRAME;
1072             }
1073             sel->evaluate = &_gmx_sel_evaluate_method;
1074             break;
1075
1076         case SEL_ARITHMETIC:
1077             sel->evaluate = &_gmx_sel_evaluate_arithmetic;
1078             break;
1079
1080         case SEL_MODIFIER:
1081             if (sel->v.type != NO_VALUE)
1082             {
1083                 sel->evaluate = &_gmx_sel_evaluate_modifier;
1084             }
1085             break;
1086
1087         case SEL_BOOLEAN:
1088             switch (sel->u.boolt)
1089             {
1090                 case BOOL_NOT: sel->evaluate = &_gmx_sel_evaluate_not; break;
1091                 case BOOL_AND: sel->evaluate = &_gmx_sel_evaluate_and; break;
1092                 case BOOL_OR:  sel->evaluate = &_gmx_sel_evaluate_or;  break;
1093                 case BOOL_XOR:
1094                     GMX_THROW(gmx::NotImplementedError("xor expressions not implemented"));
1095             }
1096             break;
1097
1098         case SEL_ROOT:
1099             sel->evaluate = &_gmx_sel_evaluate_root;
1100             break;
1101
1102         case SEL_SUBEXPR:
1103             sel->evaluate = ((sel->cdata->flags & SEL_CDATA_SIMPLESUBEXPR)
1104                              ? &_gmx_sel_evaluate_subexpr_simple
1105                              : &_gmx_sel_evaluate_subexpr);
1106             break;
1107
1108         case SEL_SUBEXPRREF:
1109             sel->evaluate = ((sel->cdata->flags & SEL_CDATA_SIMPLESUBEXPR)
1110                              ? &_gmx_sel_evaluate_subexprref_simple
1111                              : &_gmx_sel_evaluate_subexprref);
1112             break;
1113
1114         case SEL_GROUPREF:
1115             GMX_THROW(gmx::APIError("Unresolved group reference in compilation"));
1116     }
1117     sel->cdata->evaluate = sel->evaluate;
1118 }
1119
1120 /*! \brief
1121  * Sets the memory pool for selection elements that can use it.
1122  *
1123  * \param     sel      Root of the selection subtree to process.
1124  * \param[in] mempool  Memory pool to use.
1125  */
1126 static void
1127 setup_memory_pooling(const SelectionTreeElementPointer &sel,
1128                      gmx_sel_mempool_t *mempool)
1129 {
1130     if (sel->type != SEL_SUBEXPRREF)
1131     {
1132         SelectionTreeElementPointer child = sel->child;
1133         while (child)
1134         {
1135             if ((sel->type == SEL_BOOLEAN && (child->flags & SEL_DYNAMIC))
1136                 || (sel->type == SEL_ARITHMETIC && child->type != SEL_CONST
1137                     && !(child->flags & SEL_SINGLEVAL))
1138                 || (sel->type == SEL_SUBEXPR
1139                     && !(sel->cdata->flags & SEL_CDATA_SIMPLESUBEXPR)))
1140             {
1141                 child->mempool = mempool;
1142                 if (child->type == SEL_SUBEXPRREF
1143                     && (child->child->cdata->flags & SEL_CDATA_SIMPLESUBEXPR))
1144                 {
1145                     child->child->child->mempool = mempool;
1146                 }
1147             }
1148             setup_memory_pooling(child, mempool);
1149             child = child->next;
1150         }
1151     }
1152 }
1153
1154 /*! \brief
1155  * Prepares the selection (sub)tree for evaluation.
1156  *
1157  * \param[in,out] sel Root of the selection subtree to prepare.
1158  *
1159  * It also allocates memory for the \p sel->v.u.g or \p sel->v.u.p
1160  * structure if required.
1161  */
1162 static void
1163 init_item_evaloutput(const SelectionTreeElementPointer &sel)
1164 {
1165     GMX_ASSERT(!(sel->child == NULL &&
1166                  (sel->type == SEL_SUBEXPRREF || sel->type == SEL_SUBEXPR)),
1167                "Subexpression elements should always have a child element");
1168
1169     /* Process children. */
1170     if (sel->type != SEL_SUBEXPRREF)
1171     {
1172         SelectionTreeElementPointer child = sel->child;
1173         while (child)
1174         {
1175             init_item_evaloutput(child);
1176             child = child->next;
1177         }
1178     }
1179
1180     if (sel->type == SEL_SUBEXPR
1181         && (sel->cdata->flags & SEL_CDATA_SIMPLESUBEXPR))
1182     {
1183         sel->flags &= ~(SEL_ALLOCVAL | SEL_ALLOCDATA);
1184         if (sel->v.type == GROUP_VALUE || sel->v.type == POS_VALUE)
1185         {
1186             _gmx_selvalue_setstore(&sel->v, sel->child->v.u.ptr);
1187         }
1188     }
1189     else if (sel->type == SEL_SUBEXPR
1190              && (sel->cdata->flags & SEL_CDATA_FULLEVAL))
1191     {
1192         sel->evaluate = &_gmx_sel_evaluate_subexpr_staticeval;
1193         sel->cdata->evaluate = sel->evaluate;
1194         sel->child->mempool = NULL;
1195         sel->flags &= ~(SEL_ALLOCVAL | SEL_ALLOCDATA);
1196         if (sel->v.type == GROUP_VALUE || sel->v.type == POS_VALUE)
1197         {
1198             _gmx_selvalue_setstore(&sel->v, sel->child->v.u.ptr);
1199         }
1200     }
1201     else if (sel->type == SEL_SUBEXPRREF
1202              && (sel->cdata->flags & SEL_CDATA_SIMPLESUBEXPR))
1203     {
1204         if (sel->v.u.ptr)
1205         {
1206             _gmx_selvalue_setstore(&sel->child->v, sel->v.u.ptr);
1207             sel->child->child->freeValues();
1208             sel->child->child->flags &= ~(SEL_ALLOCVAL | SEL_ALLOCDATA);
1209             sel->child->child->flags |= (sel->flags & SEL_ALLOCDATA);
1210             _gmx_selvalue_setstore(&sel->child->child->v, sel->v.u.ptr);
1211         }
1212         else if (sel->v.type == GROUP_VALUE || sel->v.type == POS_VALUE)
1213         {
1214             _gmx_selvalue_setstore(&sel->v, sel->child->child->v.u.ptr);
1215         }
1216         sel->flags &= ~(SEL_ALLOCVAL | SEL_ALLOCDATA);
1217     }
1218
1219     /* Make sure that the group/position structure is allocated. */
1220     if (!sel->v.u.ptr && (sel->flags & SEL_ALLOCVAL))
1221     {
1222         if (sel->v.type == GROUP_VALUE || sel->v.type == POS_VALUE)
1223         {
1224             _gmx_selvalue_reserve(&sel->v, 1);
1225             sel->v.nr = 1;
1226         }
1227     }
1228 }
1229
1230
1231 /********************************************************************
1232  * COMPILER DATA INITIALIZATION
1233  ********************************************************************/
1234
1235 /*! \brief
1236  * Allocates memory for the compiler data and initializes the structure.
1237  *
1238  * \param sel Root of the selection subtree to process.
1239  */
1240 static void
1241 init_item_compilerdata(const SelectionTreeElementPointer &sel)
1242 {
1243     /* Allocate the compiler data structure */
1244     snew(sel->cdata, 1);
1245
1246     /* Store the real evaluation method because the compiler will replace it */
1247     sel->cdata->evaluate = sel->evaluate;
1248
1249     /* This will be computed separately. */
1250     sel->cdata->refcount = 0;
1251
1252     /* Initialize the flags */
1253     sel->cdata->flags = SEL_CDATA_STATICEVAL;
1254     if (!(sel->flags & SEL_DYNAMIC))
1255     {
1256         sel->cdata->flags |= SEL_CDATA_STATIC;
1257     }
1258     if (sel->type == SEL_SUBEXPR)
1259     {
1260         sel->cdata->flags |= SEL_CDATA_EVALMAX;
1261     }
1262     /* Set the full evaluation flag for subexpressions that require it;
1263      * the subexpression has already been initialized, so we can simply
1264      * access its compilation flags.*/
1265     if (sel->type == SEL_EXPRESSION || sel->type == SEL_MODIFIER)
1266     {
1267         SelectionTreeElementPointer child = sel->child;
1268         while (child)
1269         {
1270             if (!(child->flags & SEL_ATOMVAL) && child->child)
1271             {
1272                 child->child->cdata->flags |= SEL_CDATA_FULLEVAL;
1273             }
1274             child = child->next;
1275         }
1276     }
1277     else if (sel->type == SEL_ROOT && sel->child->type == SEL_SUBEXPRREF)
1278     {
1279         sel->child->child->cdata->flags |= SEL_CDATA_FULLEVAL;
1280     }
1281
1282     /* Initialize children */
1283     if (sel->type != SEL_SUBEXPRREF)
1284     {
1285         SelectionTreeElementPointer child = sel->child;
1286         while (child)
1287         {
1288             init_item_compilerdata(child);
1289             child = child->next;
1290         }
1291     }
1292
1293     /* Determine whether we should evaluate the minimum or the maximum
1294      * for the children of this element. */
1295     if (sel->type == SEL_BOOLEAN)
1296     {
1297         bool  bEvalMax;
1298
1299         bEvalMax = (sel->u.boolt == BOOL_AND);
1300         SelectionTreeElementPointer child = sel->child;
1301         while (child)
1302         {
1303             if (bEvalMax)
1304             {
1305                 child->cdata->flags |= SEL_CDATA_EVALMAX;
1306             }
1307             else if (child->type == SEL_BOOLEAN && child->u.boolt == BOOL_NOT)
1308             {
1309                 child->child->cdata->flags |= SEL_CDATA_EVALMAX;
1310             }
1311             child = child->next;
1312         }
1313     }
1314     else if (sel->type == SEL_EXPRESSION || sel->type == SEL_MODIFIER
1315              || sel->type == SEL_SUBEXPR)
1316     {
1317         SelectionTreeElementPointer child = sel->child;
1318         while (child)
1319         {
1320             child->cdata->flags |= SEL_CDATA_EVALMAX;
1321             child = child->next;
1322         }
1323     }
1324 }
1325
1326 /*! \brief
1327  * Initializes the static evaluation flag for a selection subtree.
1328  *
1329  * \param[in,out] sel  Root of the selection subtree to process.
1330  *
1331  * Sets the \c bStaticEval in the compiler data structure:
1332  * for any element for which the evaluation group may depend on the trajectory
1333  * frame, the flag is cleared.
1334  *
1335  * reorder_boolean_static_children() should have been called.
1336  */
1337 static void
1338 init_item_staticeval(const SelectionTreeElementPointer &sel)
1339 {
1340     /* Subexpressions with full evaluation should always have bStaticEval,
1341      * so don't do anything if a reference to them is encountered. */
1342     if (sel->type == SEL_SUBEXPRREF
1343         && (sel->child->cdata->flags & SEL_CDATA_FULLEVAL))
1344     {
1345         return;
1346     }
1347
1348     /* Propagate the bStaticEval flag to children if it is not set */
1349     if (!(sel->cdata->flags & SEL_CDATA_STATICEVAL))
1350     {
1351         SelectionTreeElementPointer child = sel->child;
1352         while (child)
1353         {
1354             if ((sel->type != SEL_EXPRESSION && sel->type != SEL_MODIFIER)
1355                 || (child->flags & SEL_ATOMVAL))
1356             {
1357                 if (child->cdata->flags & SEL_CDATA_STATICEVAL)
1358                 {
1359                     child->cdata->flags &= ~SEL_CDATA_STATICEVAL;
1360                     init_item_staticeval(child);
1361                 }
1362             }
1363             child = child->next;
1364         }
1365     }
1366     else /* bStaticEval is set */
1367     {
1368         /* For boolean expressions, any expression after the first dynamic
1369          * expression should not have bStaticEval. */
1370         if (sel->type == SEL_BOOLEAN)
1371         {
1372             SelectionTreeElementPointer child = sel->child;
1373             while (child && !(child->flags & SEL_DYNAMIC))
1374             {
1375                 child = child->next;
1376             }
1377             if (child)
1378             {
1379                 child = child->next;
1380             }
1381             while (child)
1382             {
1383                 child->cdata->flags &= ~SEL_CDATA_STATICEVAL;
1384                 child = child->next;
1385             }
1386         }
1387
1388         /* Process the children */
1389         SelectionTreeElementPointer child = sel->child;
1390         while (child)
1391         {
1392             init_item_staticeval(child);
1393             child = child->next;
1394         }
1395     }
1396 }
1397
1398 /*! \brief
1399  * Compute reference counts for subexpressions.
1400  *
1401  * \param sel Root of the selection subtree to process.
1402  */
1403 static void
1404 init_item_subexpr_refcount(const SelectionTreeElementPointer &sel)
1405 {
1406     // Reset the counter when the subexpression is first encountered.
1407     if (sel->type == SEL_ROOT && sel->child->type == SEL_SUBEXPR
1408         && sel->child->cdata)
1409     {
1410         sel->child->cdata->refcount = 0;
1411     }
1412
1413     if (sel->type == SEL_SUBEXPRREF)
1414     {
1415         ++sel->child->cdata->refcount;
1416     }
1417     else
1418     {
1419         SelectionTreeElementPointer child = sel->child;
1420         while (child)
1421         {
1422             init_item_subexpr_refcount(child);
1423             child = child->next;
1424         }
1425     }
1426 }
1427
1428 /*! \brief
1429  * Initializes compiler flags for subexpressions.
1430  *
1431  * \param sel Root of the selection subtree to process.
1432  */
1433 static void
1434 init_item_subexpr_flags(const SelectionTreeElementPointer &sel)
1435 {
1436     if (sel->type == SEL_SUBEXPR)
1437     {
1438         if (sel->cdata->refcount == 1)
1439         {
1440             sel->cdata->flags |= SEL_CDATA_SIMPLESUBEXPR;
1441         }
1442         else if (!(sel->cdata->flags & SEL_CDATA_FULLEVAL))
1443         {
1444             sel->cdata->flags |= SEL_CDATA_COMMONSUBEXPR;
1445         }
1446     }
1447     else if (sel->type == SEL_SUBEXPRREF
1448              && (sel->child->cdata->flags & SEL_CDATA_SIMPLESUBEXPR))
1449     {
1450         sel->cdata->flags |= SEL_CDATA_SIMPLESUBEXPR;
1451     }
1452
1453     /* Process children, but only follow subexpression references if the
1454      * common subexpression flag needs to be propagated. */
1455     if (sel->type != SEL_SUBEXPRREF
1456         || ((sel->cdata->flags & SEL_CDATA_COMMONSUBEXPR)
1457             && !(sel->child->cdata->flags & SEL_CDATA_SIMPLESUBEXPR)))
1458     {
1459         SelectionTreeElementPointer child = sel->child;
1460         while (child)
1461         {
1462             if (!(child->cdata->flags & SEL_CDATA_COMMONSUBEXPR))
1463             {
1464                 if (sel->type != SEL_EXPRESSION || (child->flags & SEL_ATOMVAL))
1465                 {
1466                     child->cdata->flags |=
1467                         (sel->cdata->flags & SEL_CDATA_COMMONSUBEXPR);
1468                 }
1469                 init_item_subexpr_flags(child);
1470             }
1471             child = child->next;
1472         }
1473     }
1474 }
1475
1476 /*! \brief
1477  * Initializes the gmin and gmax fields of the compiler data structure.
1478  *
1479  * \param sel Root of the selection subtree to process.
1480  */
1481 static void
1482 init_item_minmax_groups(const SelectionTreeElementPointer &sel)
1483 {
1484     /* Process children. */
1485     if (sel->type != SEL_SUBEXPRREF)
1486     {
1487         SelectionTreeElementPointer child = sel->child;
1488         while (child)
1489         {
1490             init_item_minmax_groups(child);
1491             child = child->next;
1492         }
1493     }
1494
1495     /* Initialize the minimum and maximum evaluation groups. */
1496     if (sel->type != SEL_ROOT && sel->v.type != NO_VALUE)
1497     {
1498         if (sel->v.type == GROUP_VALUE
1499             && (sel->cdata->flags & SEL_CDATA_STATIC))
1500         {
1501             sel->cdata->gmin = sel->v.u.g;
1502             sel->cdata->gmax = sel->v.u.g;
1503         }
1504         else if (sel->type == SEL_SUBEXPR
1505                  && ((sel->cdata->flags & SEL_CDATA_SIMPLESUBEXPR)
1506                      || (sel->cdata->flags & SEL_CDATA_FULLEVAL)))
1507         {
1508             GMX_ASSERT(sel->child,
1509                        "Subexpression elements should always have a child element");
1510             sel->cdata->gmin = sel->child->cdata->gmin;
1511             sel->cdata->gmax = sel->child->cdata->gmax;
1512         }
1513         else
1514         {
1515             sel->cdata->flags |= SEL_CDATA_MINMAXALLOC | SEL_CDATA_DOMINMAX;
1516             snew(sel->cdata->gmin, 1);
1517             snew(sel->cdata->gmax, 1);
1518         }
1519     }
1520 }
1521
1522
1523 /********************************************************************
1524  * EVALUATION GROUP INITIALIZATION
1525  ********************************************************************/
1526
1527 /*! \brief
1528  * Initializes evaluation groups for root items.
1529  *
1530  * \param[in,out] sc   Selection collection data.
1531  *
1532  * The evaluation group of each \ref SEL_ROOT element corresponding to a
1533  * selection in \p sc is set to \p gall.  The same is done for \ref SEL_ROOT
1534  * elements corresponding to subexpressions that need full evaluation.
1535  */
1536 static void
1537 initialize_evalgrps(gmx_ana_selcollection_t *sc)
1538 {
1539     SelectionTreeElementPointer root = sc->root;
1540     while (root)
1541     {
1542         GMX_RELEASE_ASSERT(root->child,
1543                            "Root elements should always have a child");
1544         if (root->child->type != SEL_SUBEXPR
1545             || (root->child->cdata->flags & SEL_CDATA_FULLEVAL))
1546         {
1547             gmx_ana_index_set(&root->u.cgrp, sc->gall.isize, sc->gall.index,
1548                               root->u.cgrp.name, 0);
1549         }
1550         root = root->next;
1551     }
1552 }
1553
1554
1555 /********************************************************************
1556  * STATIC ANALYSIS
1557  ********************************************************************/
1558
1559 /*! \brief
1560  * Marks a subtree completely dynamic or undoes such a change.
1561  *
1562  * \param     sel      Selection subtree to mark.
1563  * \param[in] bDynamic If true, the \p bStatic flag of the whole
1564  *   selection subtree is cleared. If false, the flag is restored to
1565  *   using \ref SEL_DYNAMIC.
1566  *
1567  * Does not descend into parameters of methods unless the parameters
1568  * are evaluated for each atom.
1569  */
1570 static void
1571 mark_subexpr_dynamic(const SelectionTreeElementPointer &sel,
1572                      bool bDynamic)
1573 {
1574     if (!bDynamic && !(sel->flags & SEL_DYNAMIC))
1575     {
1576         sel->cdata->flags |= SEL_CDATA_STATIC;
1577     }
1578     else
1579     {
1580         sel->cdata->flags &= ~SEL_CDATA_STATIC;
1581     }
1582     SelectionTreeElementPointer child = sel->child;
1583     while (child)
1584     {
1585         if (sel->type != SEL_EXPRESSION || child->type != SEL_SUBEXPRREF
1586             || (child->u.param->flags & SPAR_ATOMVAL))
1587         {
1588             mark_subexpr_dynamic(child, bDynamic);
1589         }
1590         child = child->next;
1591     }
1592 }
1593
1594 /*! \brief
1595  * Frees memory for subexpressions that are no longer needed.
1596  *
1597  * \param     sel      Selection subtree to check.
1598  *
1599  * Checks whether the subtree rooted at \p sel refers to any \ref SEL_SUBEXPR
1600  * elements that are not referred to by anything else except their own root
1601  * element. If such elements are found, all memory allocated for them is freed
1602  * except the actual element. The element is left because otherwise a dangling
1603  * pointer would be left at the root element, which is not traversed by this
1604  * function. Later compilation passes remove the stub elements.
1605  */
1606 static void
1607 release_subexpr_memory(const SelectionTreeElementPointer &sel)
1608 {
1609     if (sel->type == SEL_SUBEXPRREF)
1610     {
1611         const SelectionTreeElementPointer &subexpr = sel->child;
1612         if (subexpr.use_count() == 2)
1613         {
1614             release_subexpr_memory(subexpr);
1615             // Free children.
1616             subexpr->child.reset();
1617             subexpr->freeValues();
1618             subexpr->freeExpressionData();
1619             subexpr->freeCompilerData();
1620         }
1621     }
1622     else
1623     {
1624         SelectionTreeElementPointer child = sel->child;
1625         while (child)
1626         {
1627             release_subexpr_memory(child);
1628             child = child->next;
1629         }
1630     }
1631 }
1632
1633 /*! \brief
1634  * Makes an evaluated selection element static.
1635  *
1636  * \param     sel   Selection element to make static.
1637  *
1638  * The evaluated value becomes the value of the static element.
1639  * The element type is changed to SEL_CONST and the children are
1640  * deleted.
1641  */
1642 static void
1643 make_static(const SelectionTreeElementPointer &sel)
1644 {
1645     /* If this is a subexpression reference and the data is stored in the
1646      * child, we transfer data ownership before doing anything else. */
1647     if (sel->type == SEL_SUBEXPRREF
1648         && (sel->cdata->flags & SEL_CDATA_SIMPLESUBEXPR))
1649     {
1650         if (sel->child->child->flags & SEL_ALLOCDATA)
1651         {
1652             sel->flags |= SEL_ALLOCDATA;
1653             sel->child->child->flags &= ~SEL_ALLOCDATA;
1654         }
1655         if (sel->child->child->flags & SEL_ALLOCVAL)
1656         {
1657             sel->flags |= SEL_ALLOCVAL;
1658             sel->v.nalloc = sel->child->child->v.nalloc;
1659             sel->child->child->flags &= ~SEL_ALLOCVAL;
1660             sel->child->child->v.nalloc = -1;
1661         }
1662     }
1663     /* Free the children. */
1664     release_subexpr_memory(sel);
1665     sel->child.reset();
1666     /* Free the expression data as it is no longer needed */
1667     sel->freeExpressionData();
1668     /* Make the item static */
1669     sel->type            = SEL_CONST;
1670     sel->evaluate        = NULL;
1671     sel->cdata->evaluate = NULL;
1672     /* Set the group value.
1673      * freeExpressionData() frees the cgrp group, so we can just override it.
1674      * */
1675     if (sel->v.type == GROUP_VALUE)
1676     {
1677         gmx_ana_index_set(&sel->u.cgrp, sel->v.u.g->isize, sel->v.u.g->index, NULL, 0);
1678     }
1679 }
1680
1681 /*! \brief
1682  * Evaluates a constant expression during analyze_static().
1683  *
1684  * \param[in]     data Evaluation data.
1685  * \param[in,out] sel Selection to process.
1686  * \param[in]     g   The evaluation group.
1687  * \returns       0 on success, a non-zero error code on error.
1688  */
1689 static void
1690 process_const(gmx_sel_evaluate_t *data,
1691               const SelectionTreeElementPointer &sel,
1692               gmx_ana_index_t *g)
1693 {
1694     if (sel->v.type == GROUP_VALUE)
1695     {
1696         if (sel->cdata->evaluate)
1697         {
1698             sel->cdata->evaluate(data, sel, g);
1699         }
1700     }
1701     /* Other constant expressions do not need evaluation */
1702 }
1703
1704 /*! \brief
1705  * Sets the parameter value pointer for \ref SEL_SUBEXPRREF params.
1706  *
1707  * \param[in,out] sel Selection to process.
1708  *
1709  * Copies the value pointer of \p sel to \c sel->u.param if one is present
1710  * and should receive the value from the compiler
1711  * (most parameter values are handled during parsing).
1712  * If \p sel is not of type \ref SEL_SUBEXPRREF, or if \c sel->u.param is NULL,
1713  * the function does nothing.
1714  * Also, if the \c sel->u.param does not have \ref SPAR_VARNUM or
1715  * \ref SPAR_ATOMVAL, the function returns immediately.
1716  */
1717 static void
1718 store_param_val(const SelectionTreeElementPointer &sel)
1719 {
1720     /* Return immediately if there is no parameter. */
1721     if (sel->type != SEL_SUBEXPRREF || !sel->u.param)
1722     {
1723         return;
1724     }
1725
1726     /* Or if the value does not need storing. */
1727     if (!(sel->u.param->flags & (SPAR_VARNUM | SPAR_ATOMVAL)))
1728     {
1729         return;
1730     }
1731
1732     if (sel->v.type == INT_VALUE || sel->v.type == REAL_VALUE
1733         || sel->v.type == STR_VALUE)
1734     {
1735         _gmx_selvalue_setstore(&sel->u.param->val, sel->v.u.ptr);
1736     }
1737 }
1738
1739 /*! \brief
1740  * Handles the initialization of a selection method during analyze_static() pass.
1741  *
1742  * \param[in,out] sel Selection element to process.
1743  * \param[in]     top Topology structure.
1744  * \param[in]     isize Size of the evaluation group for the element.
1745  * \returns       0 on success, a non-zero error code on return.
1746  *
1747  * Calls sel_initfunc() (and possibly sel_outinitfunc()) to initialize the
1748  * method.
1749  * If no \ref SPAR_ATOMVAL parameters are present, multiple initialization
1750  * is prevented by using \ref SEL_METHODINIT and \ref SEL_OUTINIT flags.
1751  */
1752 static void
1753 init_method(const SelectionTreeElementPointer &sel, t_topology *top, int isize)
1754 {
1755     /* Find out whether there are any atom-valued parameters */
1756     bool bAtomVal = false;
1757     SelectionTreeElementPointer child = sel->child;
1758     while (child)
1759     {
1760         if (child->flags & SEL_ATOMVAL)
1761         {
1762             bAtomVal = true;
1763         }
1764         child = child->next;
1765     }
1766
1767     /* Initialize the method */
1768     if (sel->u.expr.method->init
1769         && (bAtomVal || !(sel->flags & SEL_METHODINIT)))
1770     {
1771         sel->flags |= SEL_METHODINIT;
1772         sel->u.expr.method->init(top, sel->u.expr.method->nparams,
1773                 sel->u.expr.method->param, sel->u.expr.mdata);
1774     }
1775     if (bAtomVal || !(sel->flags & SEL_OUTINIT))
1776     {
1777         sel->flags |= SEL_OUTINIT;
1778         if (sel->u.expr.method->outinit)
1779         {
1780             sel->u.expr.method->outinit(top, &sel->v, sel->u.expr.mdata);
1781             if (sel->v.type != POS_VALUE && sel->v.type != GROUP_VALUE)
1782             {
1783                 alloc_selection_data(sel, isize, true);
1784             }
1785         }
1786         else
1787         {
1788             alloc_selection_data(sel, isize, true);
1789             if ((sel->flags & SEL_DYNAMIC)
1790                 && sel->v.type != GROUP_VALUE && sel->v.type != POS_VALUE)
1791             {
1792                 sel->v.nr = isize;
1793             }
1794             /* If the method is char-valued, pre-allocate the strings. */
1795             if (sel->u.expr.method->flags & SMETH_CHARVAL)
1796             {
1797                 int  i;
1798
1799                 /* A sanity check */
1800                 if (sel->v.type != STR_VALUE)
1801                 {
1802                     GMX_THROW(gmx::InternalError("Char-valued selection method in non-string element"));
1803                 }
1804                 sel->flags |= SEL_ALLOCDATA;
1805                 for (i = 0; i < isize; ++i)
1806                 {
1807                     if (sel->v.u.s[i] == NULL)
1808                     {
1809                         snew(sel->v.u.s[i], 2);
1810                     }
1811                 }
1812             }
1813         }
1814         /* Clear the values for dynamic output to avoid valgrind warnings. */
1815         if ((sel->flags & SEL_DYNAMIC) && sel->v.type == REAL_VALUE)
1816         {
1817             int i;
1818
1819             for (i = 0; i < sel->v.nr; ++i)
1820             {
1821                 sel->v.u.r[i] = 0.0;
1822             }
1823         }
1824     }
1825 }
1826
1827 /*! \brief
1828  * Evaluates the static part of a boolean expression.
1829  *
1830  * \param[in]     data Evaluation data.
1831  * \param[in,out] sel Boolean selection element whose children should be
1832  *   processed.
1833  * \param[in]     g   The evaluation group.
1834  * \returns       0 on success, a non-zero error code on error.
1835  *
1836  * reorder_item_static_children() should have been called.
1837  */
1838 static void
1839 evaluate_boolean_static_part(gmx_sel_evaluate_t *data,
1840                              const SelectionTreeElementPointer &sel,
1841                              gmx_ana_index_t *g)
1842 {
1843     /* Find the last static subexpression */
1844     SelectionTreeElementPointer child = sel->child;
1845     while (child->next && (child->next->cdata->flags & SEL_CDATA_STATIC))
1846     {
1847         child = child->next;
1848     }
1849     if (!(child->cdata->flags & SEL_CDATA_STATIC))
1850     {
1851         return;
1852     }
1853
1854     /* Evalute the static part if there is more than one expression */
1855     if (child != sel->child)
1856     {
1857         SelectionTreeElementPointer next = child->next;
1858         child->next.reset();
1859         sel->cdata->evaluate(data, sel, g);
1860         /* Replace the subexpressions with the result */
1861         child.reset(new SelectionTreeElement(SEL_CONST));
1862         child->flags      = SEL_FLAGSSET | SEL_SINGLEVAL | SEL_ALLOCVAL | SEL_ALLOCDATA;
1863         _gmx_selelem_set_vtype(child, GROUP_VALUE);
1864         child->evaluate   = NULL;
1865         _gmx_selvalue_reserve(&child->v, 1);
1866         gmx_ana_index_copy(child->v.u.g, sel->v.u.g, true);
1867         init_item_compilerdata(child);
1868         init_item_minmax_groups(child);
1869         child->cdata->flags &= ~SEL_CDATA_STATICEVAL;
1870         child->cdata->flags |= sel->cdata->flags & SEL_CDATA_STATICEVAL;
1871         child->next = next;
1872         // Frees the old static subexpressions.
1873         sel->child = child;
1874     }
1875     else if (child->evaluate)
1876     {
1877         child->evaluate(data, child, g);
1878     }
1879     /* Set the evaluation function for the constant element.
1880      * We never need to evaluate the element again during compilation,
1881      * but we may need to evaluate the static part again if the
1882      * expression is not an OR with a static evaluation group.
1883      * If we reach here with a NOT expression, the NOT expression
1884      * is also static, and will be made a constant later, so don't waste
1885      * time copying the group. */
1886     child->evaluate = NULL;
1887     if (sel->u.boolt == BOOL_NOT
1888         || ((sel->cdata->flags & SEL_CDATA_STATICEVAL)
1889             && sel->u.boolt == BOOL_OR))
1890     {
1891         child->cdata->evaluate = NULL;
1892     }
1893     else
1894     {
1895         child->cdata->evaluate = &_gmx_sel_evaluate_static;
1896         /* The cgrp has only been allocated if it originated from an
1897          * external index group. In that case, we need special handling
1898          * to preserve the name of the group and to not leak memory.
1899          * If cgrp has been set in make_static(), it is not allocated,
1900          * and hence we can overwrite it safely. */
1901         if (child->u.cgrp.nalloc_index > 0)
1902         {
1903             char *name = child->u.cgrp.name;
1904             gmx_ana_index_copy(&child->u.cgrp, child->v.u.g, false);
1905             gmx_ana_index_squeeze(&child->u.cgrp);
1906             child->u.cgrp.name = name;
1907         }
1908         else
1909         {
1910             gmx_ana_index_copy(&child->u.cgrp, child->v.u.g, true);
1911         }
1912     }
1913 }
1914
1915 /*! \brief
1916  * Evaluates the minimum and maximum groups for a boolean expression.
1917  *
1918  * \param[in]  sel  \ref SEL_BOOLEAN element currently being evaluated.
1919  * \param[in]  g    Group for which \p sel has been evaluated.
1920  * \param[out] gmin Largest subset of the possible values of \p sel.
1921  * \param[out] gmax Smallest superset of the possible values of \p sel.
1922  *
1923  * This is a helper function for analyze_static() that is called for
1924  * dynamic \ref SEL_BOOLEAN elements after they have been evaluated.
1925  * It uses the minimum and maximum groups of the children to calculate
1926  * the minimum and maximum groups for \p sel, and also updates the static
1927  * part of \p sel (which is in the first child) if the children give
1928  * cause for this.
1929  *
1930  * This function may allocate some extra memory for \p gmin and \p gmax,
1931  * but as these groups are freed at the end of analyze_static() (which is
1932  * reached shortly after this function returns), this should not be a major
1933  * problem.
1934  */
1935 static void
1936 evaluate_boolean_minmax_grps(const SelectionTreeElementPointer &sel,
1937                              gmx_ana_index_t *g,
1938                              gmx_ana_index_t *gmin, gmx_ana_index_t *gmax)
1939 {
1940     SelectionTreeElementPointer child;
1941
1942     switch (sel->u.boolt)
1943     {
1944         case BOOL_NOT:
1945             gmx_ana_index_reserve(gmin, g->isize);
1946             gmx_ana_index_reserve(gmax, g->isize);
1947             gmx_ana_index_difference(gmax, g, sel->child->cdata->gmin);
1948             gmx_ana_index_difference(gmin, g, sel->child->cdata->gmax);
1949             break;
1950
1951         case BOOL_AND:
1952             gmx_ana_index_copy(gmin, sel->child->cdata->gmin, true);
1953             gmx_ana_index_copy(gmax, sel->child->cdata->gmax, true);
1954             child = sel->child->next;
1955             while (child && gmax->isize > 0)
1956             {
1957                 gmx_ana_index_intersection(gmin, gmin, child->cdata->gmin);
1958                 gmx_ana_index_intersection(gmax, gmax, child->cdata->gmax);
1959                 child = child->next;
1960             }
1961             /* Update the static part if other expressions limit it */
1962             if ((sel->child->cdata->flags & SEL_CDATA_STATIC)
1963                 && sel->child->v.u.g->isize > gmax->isize)
1964             {
1965                 gmx_ana_index_copy(sel->child->v.u.g, gmax, false);
1966                 gmx_ana_index_squeeze(sel->child->v.u.g);
1967                 if (sel->child->u.cgrp.isize > 0)
1968                 {
1969                     gmx_ana_index_copy(&sel->child->u.cgrp, gmax, false);
1970                     gmx_ana_index_squeeze(&sel->child->u.cgrp);
1971                 }
1972             }
1973             break;
1974
1975         case BOOL_OR:
1976             /* We can assume here that the gmin of children do not overlap
1977              * because of the way _gmx_sel_evaluate_or() works. */
1978             gmx_ana_index_reserve(gmin, g->isize);
1979             gmx_ana_index_reserve(gmax, g->isize);
1980             gmx_ana_index_copy(gmin, sel->child->cdata->gmin, false);
1981             gmx_ana_index_copy(gmax, sel->child->cdata->gmax, false);
1982             child = sel->child->next;
1983             while (child && gmin->isize < g->isize)
1984             {
1985                 gmx_ana_index_merge(gmin, gmin, child->cdata->gmin);
1986                 gmx_ana_index_union(gmax, gmax, child->cdata->gmax);
1987                 child = child->next;
1988             }
1989             /* Update the static part if other expressions have static parts
1990              * that are not included. */
1991             if ((sel->child->cdata->flags & SEL_CDATA_STATIC)
1992                 && sel->child->v.u.g->isize < gmin->isize)
1993             {
1994                 GMX_RELEASE_ASSERT(sel->child->type == SEL_CONST,
1995                         "The first child should have already been evaluated "
1996                         "to a constant expression");
1997                 gmx_ana_index_reserve(sel->child->v.u.g, gmin->isize);
1998                 gmx_ana_index_copy(sel->child->v.u.g, gmin, false);
1999                 if (sel->child->u.cgrp.nalloc_index > 0)
2000                 {
2001                     /* Keep the name as in evaluate_boolean_static_part(). */
2002                     char *name = sel->child->u.cgrp.name;
2003                     gmx_ana_index_reserve(&sel->child->u.cgrp, gmin->isize);
2004                     gmx_ana_index_copy(&sel->child->u.cgrp, gmin, false);
2005                     sel->child->u.cgrp.name = name;
2006                 }
2007                 else
2008                 {
2009                     GMX_RELEASE_ASSERT(sel->child->u.cgrp.index == sel->child->v.u.g->index,
2010                             "If not allocated, the static group should equal the value");
2011                     sel->child->u.cgrp.isize = sel->child->v.u.g->isize;
2012                 }
2013             }
2014             break;
2015
2016         case BOOL_XOR: /* Should not be reached */
2017             GMX_THROW(gmx::NotImplementedError("xor expressions not implemented"));
2018             break;
2019     }
2020 }
2021
2022 /*! \brief
2023  * Evaluates the static parts of \p sel and analyzes the structure.
2024  * 
2025  * \param[in]     data Evaluation data.
2026  * \param[in,out] sel  Selection currently being evaluated.
2027  * \param[in]     g    Group for which \p sel should be evaluated.
2028  * \returns       0 on success, a non-zero error code on error.
2029  *
2030  * This function is used as the replacement for the
2031  * gmx::SelectionTreeElement::evaluate function pointer.
2032  * It does the single most complex task in the compiler: after all elements
2033  * have been processed, the \p gmin and \p gmax fields of \p t_compiler_data
2034  * have been properly initialized, enough memory has been allocated for
2035  * storing the value of each expression, and the static parts of the 
2036  * expressions have been evaluated.
2037  * The above is exactly true only for elements other than subexpressions:
2038  * another pass is required for subexpressions that are referred to more than
2039  * once and whose evaluation group is not known in advance.
2040  */
2041 static void
2042 analyze_static(gmx_sel_evaluate_t *data,
2043                const SelectionTreeElementPointer &sel,
2044                gmx_ana_index_t *g)
2045 {
2046     bool             bDoMinMax;
2047
2048     if (sel->type != SEL_ROOT && g)
2049     {
2050         alloc_selection_data(sel, g->isize, false);
2051     }
2052
2053     bDoMinMax = (sel->cdata->flags & SEL_CDATA_DOMINMAX);
2054     if (sel->type != SEL_SUBEXPR && bDoMinMax)
2055     {
2056         gmx_ana_index_deinit(sel->cdata->gmin);
2057         gmx_ana_index_deinit(sel->cdata->gmax);
2058     }
2059
2060     /* TODO: This switch is awfully long... */
2061     switch (sel->type)
2062     {
2063         case SEL_CONST:
2064             process_const(data, sel, g);
2065             break;
2066
2067         case SEL_EXPRESSION:
2068         case SEL_MODIFIER:
2069             GMX_ASSERT(g, "group cannot be null");
2070             _gmx_sel_evaluate_method_params(data, sel, g);
2071             init_method(sel, data->top, g->isize);
2072             if (!(sel->flags & SEL_DYNAMIC))
2073             {
2074                 sel->cdata->evaluate(data, sel, g);
2075                 if (sel->cdata->flags & SEL_CDATA_STATIC)
2076                 {
2077                     make_static(sel);
2078                 }
2079             }
2080             else
2081             {
2082                 /* Modifiers need to be evaluated even though they process
2083                  * positions to get the modified output groups from the
2084                  * maximum possible selections. */
2085                 if (sel->type == SEL_MODIFIER)
2086                 {
2087                     sel->cdata->evaluate(data, sel, g);
2088                 }
2089                 if (bDoMinMax)
2090                 {
2091                     gmx_ana_index_copy(sel->cdata->gmax, g, true);
2092                 }
2093             }
2094             break;
2095
2096         case SEL_BOOLEAN:
2097             if (!(sel->flags & SEL_DYNAMIC))
2098             {
2099                 sel->cdata->evaluate(data, sel, g);
2100                 if (sel->cdata->flags & SEL_CDATA_STATIC)
2101                 {
2102                     make_static(sel);
2103                 }
2104             }
2105             else
2106             {
2107                 /* Evalute the static part if there is more than one expression */
2108                 evaluate_boolean_static_part(data, sel, g);
2109
2110                 /* Evaluate the selection.
2111                  * If the type is boolean, we must explicitly handle the
2112                  * static part evaluated in evaluate_boolean_static_part()
2113                  * here because g may be larger. */
2114                 if (sel->u.boolt == BOOL_AND && sel->child->type == SEL_CONST)
2115                 {
2116                     sel->cdata->evaluate(data, sel, sel->child->v.u.g);
2117                 }
2118                 else
2119                 {
2120                     sel->cdata->evaluate(data, sel, g);
2121                 }
2122
2123                 /* Evaluate minimal and maximal selections */
2124                 evaluate_boolean_minmax_grps(sel, g, sel->cdata->gmin,
2125                                              sel->cdata->gmax);
2126             }
2127             break;
2128
2129         case SEL_ARITHMETIC:
2130             sel->cdata->evaluate(data, sel, g);
2131             if (!(sel->flags & SEL_DYNAMIC))
2132             {
2133                 if (sel->cdata->flags & SEL_CDATA_STATIC)
2134                 {
2135                     make_static(sel);
2136                 }
2137             }
2138             else if (bDoMinMax)
2139             {
2140                 gmx_ana_index_copy(sel->cdata->gmax, g, true);
2141             }
2142             break;
2143
2144         case SEL_ROOT:
2145             sel->cdata->evaluate(data, sel, g);
2146             break;
2147
2148         case SEL_SUBEXPR:
2149             if (sel->cdata->flags & (SEL_CDATA_SIMPLESUBEXPR | SEL_CDATA_FULLEVAL))
2150             {
2151                 sel->cdata->evaluate(data, sel, g);
2152                 _gmx_selvalue_setstore(&sel->v, sel->child->v.u.ptr);
2153             }
2154             else if (sel->u.cgrp.isize == 0)
2155             {
2156                 GMX_ASSERT(g, "group cannot be null");
2157                 gmx_ana_index_reserve(&sel->u.cgrp, g->isize);
2158                 sel->cdata->evaluate(data, sel, g);
2159                 if (bDoMinMax)
2160                 {
2161                     gmx_ana_index_copy(sel->cdata->gmin, sel->child->cdata->gmin, true);
2162                     gmx_ana_index_copy(sel->cdata->gmax, sel->child->cdata->gmax, true);
2163                 }
2164             }
2165             else
2166             {
2167                 int isize = gmx_ana_index_difference_size(g, &sel->u.cgrp);
2168                 if (isize > 0)
2169                 {
2170                     isize += sel->u.cgrp.isize;
2171                     gmx_ana_index_reserve(&sel->u.cgrp, isize);
2172                     alloc_selection_data(sel, isize, false);
2173                 }
2174                 sel->cdata->evaluate(data, sel, g);
2175                 if (isize > 0 && bDoMinMax)
2176                 {
2177                     gmx_ana_index_reserve(sel->cdata->gmin,
2178                                           sel->cdata->gmin->isize
2179                                           + sel->child->cdata->gmin->isize);
2180                     gmx_ana_index_reserve(sel->cdata->gmax,
2181                                           sel->cdata->gmax->isize
2182                                           + sel->child->cdata->gmax->isize);
2183                     gmx_ana_index_merge(sel->cdata->gmin, sel->cdata->gmin,
2184                                         sel->child->cdata->gmin);
2185                     gmx_ana_index_merge(sel->cdata->gmax, sel->cdata->gmax,
2186                                         sel->child->cdata->gmax);
2187                 }
2188             }
2189             break;
2190
2191         case SEL_SUBEXPRREF:
2192             if (!g && !(sel->cdata->flags & SEL_CDATA_SIMPLESUBEXPR))
2193             {
2194                 /* The subexpression should have been evaluated if g is NULL
2195                  * (i.e., this is a method parameter or a direct value of a
2196                  * selection). */
2197                 alloc_selection_data(sel, sel->child->cdata->gmax->isize, true);
2198             }
2199             sel->cdata->evaluate(data, sel, g);
2200             if ((sel->cdata->flags & SEL_CDATA_SIMPLESUBEXPR)
2201                 && (sel->child->child->flags & SEL_ALLOCVAL))
2202             {
2203                 _gmx_selvalue_setstore(&sel->v, sel->child->child->v.u.ptr);
2204             }
2205             /* Store the parameter value if required */
2206             store_param_val(sel);
2207             if (!(sel->flags & SEL_DYNAMIC))
2208             {
2209                 if (sel->cdata->flags & SEL_CDATA_STATIC)
2210                 {
2211                     make_static(sel);
2212                 }
2213             }
2214             else if (bDoMinMax)
2215             {
2216                 if ((sel->cdata->flags & SEL_CDATA_SIMPLESUBEXPR) || !g)
2217                 {
2218                     gmx_ana_index_copy(sel->cdata->gmin, sel->child->cdata->gmin, true);
2219                     gmx_ana_index_copy(sel->cdata->gmax, sel->child->cdata->gmax, true);
2220                 }
2221                 else
2222                 {
2223                     gmx_ana_index_reserve(sel->cdata->gmin,
2224                                           min(g->isize, sel->child->cdata->gmin->isize));
2225                     gmx_ana_index_reserve(sel->cdata->gmax,
2226                                           min(g->isize, sel->child->cdata->gmax->isize));
2227                     gmx_ana_index_intersection(sel->cdata->gmin,
2228                                                sel->child->cdata->gmin, g);
2229                     gmx_ana_index_intersection(sel->cdata->gmax,
2230                                                sel->child->cdata->gmax, g);
2231                 }
2232             }
2233             break;
2234
2235         case SEL_GROUPREF: /* Should not be reached */
2236             GMX_THROW(gmx::APIError("Unresolved group reference in compilation"));
2237     }
2238
2239     /* Update the minimal and maximal evaluation groups */
2240     if (bDoMinMax)
2241     {
2242         gmx_ana_index_squeeze(sel->cdata->gmin);
2243         gmx_ana_index_squeeze(sel->cdata->gmax);
2244         sfree(sel->cdata->gmin->name);
2245         sfree(sel->cdata->gmax->name);
2246         sel->cdata->gmin->name = NULL;
2247         sel->cdata->gmax->name = NULL;
2248     }
2249
2250     /* Replace the result of the evaluation */
2251     /* This is not necessary for subexpressions or for boolean negations
2252      * because the evaluation function already has done it properly. */
2253     if (sel->v.type == GROUP_VALUE && (sel->flags & SEL_DYNAMIC)
2254         && sel->type != SEL_SUBEXPR
2255         && !(sel->type == SEL_BOOLEAN && sel->u.boolt == BOOL_NOT))
2256     {
2257         if (sel->cdata->flags & SEL_CDATA_EVALMAX)
2258         {
2259             gmx_ana_index_copy(sel->v.u.g, sel->cdata->gmax, false);
2260         }
2261         else
2262         {
2263             gmx_ana_index_copy(sel->v.u.g, sel->cdata->gmin, false);
2264         }
2265     }
2266 }
2267
2268
2269 /********************************************************************
2270  * EVALUATION GROUP INITIALIZATION
2271  ********************************************************************/
2272
2273 /*! \brief
2274  * Initializes the evaluation group for a \ref SEL_ROOT element.
2275  *
2276  * \param     root Root element to initialize.
2277  * \param[in] gall Group of all atoms.
2278  *
2279  * Checks whether it is necessary to evaluate anything through the root
2280  * element, and either clears the evaluation function or initializes the
2281  * evaluation group.
2282  */
2283 static void
2284 init_root_item(const SelectionTreeElementPointer &root,
2285                gmx_ana_index_t *gall)
2286 {
2287     const SelectionTreeElementPointer &expr = root->child;
2288     /* Subexpressions with non-static evaluation group should not be
2289      * evaluated by the root, and neither should be single-reference
2290      * subexpressions that don't evaluate for all atoms. */
2291     if (expr->type == SEL_SUBEXPR
2292         && (!(root->child->cdata->flags & SEL_CDATA_STATICEVAL)
2293             || ((root->child->cdata->flags & SEL_CDATA_SIMPLESUBEXPR)
2294                 && !(root->child->cdata->flags & SEL_CDATA_FULLEVAL))))
2295     {
2296         root->evaluate = NULL;
2297         if (root->cdata)
2298         {
2299             root->cdata->evaluate = NULL;
2300         }
2301     }
2302
2303     /* Set the evaluation group */
2304     char *name = root->u.cgrp.name;
2305     if (root->evaluate)
2306     {
2307         /* Non-atom-valued non-group expressions don't care about the group, so
2308          * don't allocate any memory for it. */
2309         if ((expr->flags & SEL_VARNUMVAL)
2310             || ((expr->flags & SEL_SINGLEVAL) && expr->v.type != GROUP_VALUE))
2311         {
2312             gmx_ana_index_set(&root->u.cgrp, -1, NULL, NULL, 0);
2313         }
2314         else if (expr->cdata->gmax->isize == gall->isize)
2315         {
2316             /* Save some memory by only referring to the global group. */
2317             gmx_ana_index_set(&root->u.cgrp, gall->isize, gall->index, NULL, 0);
2318         }
2319         else
2320         {
2321             gmx_ana_index_copy(&root->u.cgrp, expr->cdata->gmax, true);
2322         }
2323         /* For selections, store the maximum group for
2324          * gmx_ana_selcollection_evaluate_fin() as the value of the root
2325          * element (unused otherwise). */
2326         if (expr->type != SEL_SUBEXPR && expr->v.u.p->g)
2327         {
2328             SelectionTreeElementPointer child = expr;
2329
2330             /* TODO: This code is copied from parsetree.c; it would be better
2331              * to have this hardcoded only in one place. */
2332             while (child->type == SEL_MODIFIER)
2333             {
2334                 child = child->child;
2335                 if (child->type == SEL_SUBEXPRREF)
2336                 {
2337                     child = child->child->child;
2338                 }
2339             }
2340             if (child->type == SEL_SUBEXPRREF)
2341             {
2342                 child = child->child->child;
2343             }
2344             if (child->child->flags & SEL_DYNAMIC)
2345             {
2346                 _gmx_selelem_set_vtype(root, GROUP_VALUE);
2347                 root->flags  |= (SEL_ALLOCVAL | SEL_ALLOCDATA);
2348                 _gmx_selvalue_reserve(&root->v, 1);
2349                 gmx_ana_index_copy(root->v.u.g, expr->v.u.p->g, true);
2350             }
2351         }
2352     }
2353     else
2354     {
2355         gmx_ana_index_clear(&root->u.cgrp);
2356     }
2357     root->u.cgrp.name = name;
2358 }
2359
2360
2361 /********************************************************************
2362  * FINAL SUBEXPRESSION OPTIMIZATION
2363  ********************************************************************/
2364
2365 /*! \brief
2366  * Optimizes subexpression evaluation.
2367  *
2368  * \param     sel Root of the selection subtree to process.
2369  *
2370  * Optimizes away some unnecessary evaluation of subexpressions that are only
2371  * referenced once.
2372  */
2373 static void
2374 postprocess_item_subexpressions(const SelectionTreeElementPointer &sel)
2375 {
2376     GMX_ASSERT(!(sel->child == NULL &&
2377                  (sel->type == SEL_SUBEXPRREF || sel->type == SEL_SUBEXPR)),
2378                "Subexpression elements should always have a child element");
2379
2380     /* Process children. */
2381     if (sel->type != SEL_SUBEXPRREF)
2382     {
2383         SelectionTreeElementPointer child = sel->child;
2384         while (child)
2385         {
2386             postprocess_item_subexpressions(child);
2387             child = child->next;
2388         }
2389     }
2390
2391     /* Replace the evaluation function of statically evaluated subexpressions
2392      * for which the static group was not known in advance. */
2393     if (sel->type == SEL_SUBEXPR && sel->cdata->refcount > 1
2394         && (sel->cdata->flags & SEL_CDATA_STATICEVAL)
2395         && !(sel->cdata->flags & SEL_CDATA_FULLEVAL))
2396     {
2397         char *name;
2398
2399         /* We need to free memory allocated for the group, because it is no
2400          * longer needed (and would be lost on next call to the evaluation
2401          * function). But we need to preserve the name. */
2402         name = sel->u.cgrp.name;
2403         gmx_ana_index_deinit(&sel->u.cgrp);
2404         sel->u.cgrp.name = name;
2405
2406         sel->evaluate = &_gmx_sel_evaluate_subexpr_staticeval;
2407         sel->cdata->evaluate = sel->evaluate;
2408
2409         sel->child->freeValues();
2410         sel->child->mempool = NULL;
2411         _gmx_selvalue_setstore(&sel->child->v, sel->v.u.ptr);
2412         sel->child->flags &= ~(SEL_ALLOCVAL | SEL_ALLOCDATA);
2413     }
2414
2415     /* Adjust memory allocation flags for subexpressions that are used only
2416      * once.  This is not strictly necessary, but we do it to have the memory
2417      * managed consistently for all types of subexpressions. */
2418     if (sel->type == SEL_SUBEXPRREF
2419         && (sel->cdata->flags & SEL_CDATA_SIMPLESUBEXPR))
2420     {
2421         if (sel->child->child->flags & SEL_ALLOCVAL)
2422         {
2423             sel->flags |= SEL_ALLOCVAL;
2424             sel->flags |= (sel->child->child->flags & SEL_ALLOCDATA);
2425             sel->v.nalloc = sel->child->child->v.nalloc;
2426             sel->child->child->flags &= ~(SEL_ALLOCVAL | SEL_ALLOCDATA);
2427             sel->child->child->v.nalloc = -1;
2428         }
2429     }
2430
2431     /* Do the same for subexpressions that are evaluated at once for all atoms. */
2432     if (sel->type == SEL_SUBEXPR
2433         && !(sel->cdata->flags & SEL_CDATA_SIMPLESUBEXPR)
2434         && (sel->cdata->flags & SEL_CDATA_FULLEVAL))
2435     {
2436         sel->flags |= SEL_ALLOCVAL;
2437         sel->flags |= (sel->child->flags & SEL_ALLOCDATA);
2438         sel->v.nalloc = sel->child->v.nalloc;
2439         sel->child->flags &= ~(SEL_ALLOCVAL | SEL_ALLOCDATA);
2440         sel->child->v.nalloc = -1;
2441     }
2442 }
2443
2444
2445 /********************************************************************
2446  * COM CALCULATION INITIALIZATION
2447  ********************************************************************/
2448
2449 /*! \brief
2450  * Initializes COM/COG calculation for method expressions that require it.
2451  *
2452  * \param     sel    Selection subtree to process.
2453  * \param[in,out] pcc   Position calculation collection to use.
2454  * \param[in] type   Default position calculation type.
2455  * \param[in] flags  Flags for default position calculation.
2456  *
2457  * Searches recursively through the selection tree for dynamic
2458  * \ref SEL_EXPRESSION elements that define the \c gmx_ana_selmethod_t::pupdate
2459  * function.
2460  * For each such element found, position calculation is initialized
2461  * for the maximal evaluation group.
2462  * The type of the calculation is determined by \p type and \p flags.
2463  * No calculation is initialized if \p type equals \ref POS_ATOM and
2464  * the method also defines the \c gmx_ana_selmethod_t::update method.
2465  */
2466 static void
2467 init_item_comg(const SelectionTreeElementPointer &sel,
2468                gmx::PositionCalculationCollection *pcc,
2469                e_poscalc_t type, int flags)
2470 {
2471     /* Initialize COM calculation for dynamic selections now that we know the maximal evaluation group */
2472     if (sel->type == SEL_EXPRESSION && sel->u.expr.method
2473         && sel->u.expr.method->pupdate)
2474     {
2475         if (!sel->u.expr.method->update || type != POS_ATOM)
2476         {
2477             /* Create a default calculation if one does not yet exist */
2478             int cflags = 0;
2479             if (!(sel->cdata->flags & SEL_CDATA_STATICEVAL))
2480             {
2481                 cflags |= POS_DYNAMIC;
2482             }
2483             if (!sel->u.expr.pc)
2484             {
2485                 cflags |= flags;
2486                 sel->u.expr.pc = pcc->createCalculation(type, cflags);
2487             }
2488             else
2489             {
2490                 gmx_ana_poscalc_set_flags(sel->u.expr.pc, cflags);
2491             }
2492             gmx_ana_poscalc_set_maxindex(sel->u.expr.pc, sel->cdata->gmax);
2493             snew(sel->u.expr.pos, 1);
2494             gmx_ana_poscalc_init_pos(sel->u.expr.pc, sel->u.expr.pos);
2495         }
2496     }
2497
2498     /* Call recursively for all children unless the children have already been processed */
2499     if (sel->type != SEL_SUBEXPRREF)
2500     {
2501         SelectionTreeElementPointer child = sel->child;
2502         while (child)
2503         {
2504             init_item_comg(child, pcc, type, flags);
2505             child = child->next;
2506         }
2507     }
2508 }
2509
2510
2511 /********************************************************************
2512  * COMPILER DATA FREEING
2513  ********************************************************************/
2514
2515 /*! \brief
2516  * Frees the allocated compiler data recursively.
2517  *
2518  * \param     sel Root of the selection subtree to process.
2519  *
2520  * Frees the data allocated for the compilation process.
2521  */
2522 static void
2523 free_item_compilerdata(const SelectionTreeElementPointer &sel)
2524 {
2525     /* Free compilation data */
2526     sel->freeCompilerData();
2527
2528     /* Call recursively for all children unless the children have already been processed */
2529     if (sel->type != SEL_SUBEXPRREF)
2530     {
2531         SelectionTreeElementPointer child = sel->child;
2532         while (child)
2533         {
2534             free_item_compilerdata(child);
2535             child = child->next;
2536         }
2537     }
2538 }
2539
2540
2541 /********************************************************************
2542  * MAIN COMPILATION FUNCTION
2543  ********************************************************************/
2544
2545 namespace gmx
2546 {
2547
2548 SelectionCompiler::SelectionCompiler()
2549 {
2550 }
2551
2552 /*!
2553  * \param[in,out] coll Selection collection to be compiled.
2554  * \returns       0 on successful compilation, a non-zero error code on error.
2555  *
2556  * Before compilation, the selection collection should have been initialized
2557  * with gmx_ana_selcollection_parse_*().
2558  * The compiled selection collection can be passed to
2559  * gmx_ana_selcollection_evaluate() to evaluate the selection for a frame.
2560  * If an error occurs, \p sc is cleared.
2561  *
2562  * The covered fraction information in \p sc is initialized to
2563  * \ref CFRAC_NONE.
2564  */
2565 void
2566 SelectionCompiler::compile(SelectionCollection *coll)
2567 {
2568     gmx_ana_selcollection_t *sc = &coll->impl_->sc_;
2569     gmx_sel_evaluate_t  evaldata;
2570     SelectionTreeElementPointer item;
2571     e_poscalc_t  post;
2572     size_t       i;
2573     int          flags;
2574     bool         bDebug = (coll->impl_->debugLevel_ >= 2
2575                            && coll->impl_->debugLevel_ != 3);
2576
2577     /* FIXME: Clean up the collection on exceptions */
2578
2579     sc->mempool = _gmx_sel_mempool_create();
2580     _gmx_sel_evaluate_init(&evaldata, sc->mempool, &sc->gall,
2581                            sc->top, NULL, NULL);
2582
2583     /* Clear the symbol table because it is not possible to parse anything
2584      * after compilation, and variable references in the symbol table can
2585      * also mess up the compilation and/or become invalid.
2586      */
2587     coll->impl_->clearSymbolTable();
2588
2589     /* Loop through selections and initialize position keyword defaults if no
2590      * other value has been provided.
2591      */
2592     for (i = 0; i < sc->sel.size(); ++i)
2593     {
2594         gmx::internal::SelectionData &sel = *sc->sel[i];
2595         init_pos_keyword_defaults(&sel.rootElement(),
2596                                   coll->impl_->spost_.c_str(),
2597                                   coll->impl_->rpost_.c_str(),
2598                                   &sel);
2599     }
2600
2601     /* Remove any unused variables. */
2602     sc->root = remove_unused_subexpressions(sc->root);
2603     /* Extract subexpressions into separate roots */
2604     sc->root = extract_subexpressions(sc->root);
2605
2606     /* Initialize the evaluation callbacks and process the tree structure
2607      * to conform to the expectations of the callback functions. */
2608     /* Also, initialize and allocate the compiler data structure */
2609     item = sc->root;
2610     while (item)
2611     {
2612         /* Process boolean and arithmetic expressions. */
2613         optimize_boolean_expressions(item);
2614         reorder_boolean_static_children(item);
2615         optimize_arithmetic_expressions(item);
2616         /* Initialize the compiler data */
2617         init_item_compilerdata(item);
2618         init_item_staticeval(item);
2619         init_item_subexpr_refcount(item);
2620         item = item->next;
2621     }
2622     /* Initialize subexpression flags.
2623      * Requires compiler flags for the full tree. */
2624     item = sc->root;
2625     while (item)
2626     {
2627         init_item_subexpr_flags(item);
2628         item = item->next;
2629     }
2630     /* Initialize evaluation.
2631      * Requires subexpression flags. */
2632     item = sc->root;
2633     while (item)
2634     {
2635         init_item_evalfunc(item);
2636         setup_memory_pooling(item, sc->mempool);
2637         init_item_evaloutput(item);
2638         item = item->next;
2639     }
2640     /* Initialize minimum/maximum index groups.
2641      * Requires evaluation output for the full tree. */
2642     item = sc->root;
2643     while (item)
2644     {
2645         init_item_minmax_groups(item);
2646         item = item->next;
2647     }
2648     /* Initialize the evaluation index groups */
2649     initialize_evalgrps(sc);
2650
2651     if (bDebug)
2652     {
2653         fprintf(stderr, "\nTree after initial compiler processing:\n");
2654         coll->printTree(stderr, false);
2655     }
2656
2657     /* Evaluate all static parts of the selection and analyze the tree
2658      * to allocate enough memory to store the value of each dynamic subtree. */
2659     item = sc->root;
2660     while (item)
2661     {
2662         if (item->child->cdata->flags & SEL_CDATA_COMMONSUBEXPR)
2663         {
2664             mark_subexpr_dynamic(item->child, true);
2665         }
2666         set_evaluation_function(item, &analyze_static);
2667         item->evaluate(&evaldata, item, NULL);
2668         item = item->next;
2669     }
2670
2671     /* At this point, static subexpressions no longer have references to them,
2672      * so they can be removed. */
2673     sc->root = remove_unused_subexpressions(sc->root);
2674     // Update the reference counts for consistency (only used for the
2675     // debugging output below).
2676     item = sc->root;
2677     while (item)
2678     {
2679         init_item_subexpr_refcount(item);
2680         item = item->next;
2681     }
2682
2683     if (bDebug)
2684     {
2685         fprintf(stderr, "\nTree after first analysis pass:\n");
2686         coll->printTree(stderr, false);
2687     }
2688
2689     /* Do a second pass to evaluate static parts of common subexpressions */
2690     item = sc->root;
2691     while (item)
2692     {
2693         if (item->child->cdata->flags & SEL_CDATA_COMMONSUBEXPR)
2694         {
2695             bool bMinMax = item->child->cdata->flags & SEL_CDATA_DOMINMAX;
2696
2697             mark_subexpr_dynamic(item->child, false);
2698             item->child->u.cgrp.isize = 0;
2699             /* We won't clear item->child->v.u.g here, because it may
2700              * be static, and hence actually point to item->child->cdata->gmax,
2701              * which is used below. We could also check whether this is the
2702              * case and only clear the group otherwise, but because the value
2703              * is actually overwritten immediately in the evaluate call, we
2704              * won't, because similar problems may arise if gmax handling ever
2705              * changes and the check were not updated.
2706              * For the same reason, we clear the min/max flag so that the
2707              * evaluation group doesn't get messed up. */
2708             set_evaluation_function(item, &analyze_static);
2709             item->child->cdata->flags &= ~SEL_CDATA_DOMINMAX;
2710             item->evaluate(&evaldata, item->child, item->child->cdata->gmax);
2711             if (bMinMax)
2712             {
2713                 item->child->cdata->flags |= SEL_CDATA_DOMINMAX;
2714             }
2715         }
2716         item = item->next;
2717     }
2718
2719     /* We need a yet another pass of subexpression removal to remove static
2720      * subexpressions referred to by common dynamic subexpressions. */
2721     sc->root = remove_unused_subexpressions(sc->root);
2722     // Update the reference counts, used by postprocess_item_subexpressions().
2723     item = sc->root;
2724     while (item)
2725     {
2726         init_item_subexpr_refcount(item);
2727         item = item->next;
2728     }
2729
2730     if (bDebug)
2731     {
2732         fprintf(stderr, "\nTree after second analysis pass:\n");
2733         coll->printTree(stderr, false);
2734     }
2735
2736     /* Initialize evaluation groups, position calculations for methods, perform
2737      * some final optimization, and free the memory allocated for the
2738      * compilation. */
2739     /* By default, use whole residues/molecules. */
2740     flags = POS_COMPLWHOLE;
2741     PositionCalculationCollection::typeFromEnum(coll->impl_->rpost_.c_str(),
2742                                                 &post, &flags);
2743     item = sc->root;
2744     while (item)
2745     {
2746         init_root_item(item, &sc->gall);
2747         postprocess_item_subexpressions(item);
2748         init_item_comg(item, &sc->pcc, post, flags);
2749         free_item_compilerdata(item);
2750         item = item->next;
2751     }
2752
2753     /* Allocate memory for the evaluation memory pool. */
2754     _gmx_sel_mempool_reserve(sc->mempool, 0);
2755
2756     /* Finish up by calculating total masses and charges. */
2757     for (i = 0; i < sc->sel.size(); ++i)
2758     {
2759         sc->sel[i]->initializeMassesAndCharges(sc->top);
2760     }
2761 }
2762
2763 } // namespace gmx