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