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