Gromacs  2020.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Typedefs | Enumerations | Functions
#include "gmxpre.h"
#include "compiler.h"
#include <cmath>
#include <cstdarg>
#include <algorithm>
#include "gromacs/math/vec.h"
#include "gromacs/selection/indexutil.h"
#include "gromacs/selection/selection.h"
#include "gromacs/utility/exceptions.h"
#include "gromacs/utility/gmxassert.h"
#include "gromacs/utility/smalloc.h"
#include "gromacs/utility/stringutil.h"
#include "evaluate.h"
#include "keywords.h"
#include "mempool.h"
#include "poscalc.h"
#include "selectioncollection_impl.h"
#include "selelem.h"
#include "selmethod.h"
+ Include dependency graph for compiler.cpp:

Description

Selection compilation and optimization.

Todo:
Better error handling and memory management in error situations. At least, the main compilation function leaves the selection collection in a bad state if an error occurs.
Todo:
The memory usage could still be optimized. Use of memory pooling could still be extended, and a lot of redundant gmin/gmax data could be eliminated for complex arithmetic expressions.
Author
Teemu Murtola teemu.nosp@m..mur.nosp@m.tola@.nosp@m.gmai.nosp@m.l.com

Classes

struct  t_compiler_data
 Internal data structure used by the compiler. More...
 

Typedefs

typedef struct t_compiler_data t_compiler_data
 Internal data structure used by the compiler. More...
 

Enumerations

enum  {
  SEL_CDATA_FULLEVAL = 1, SEL_CDATA_STATIC = 2, SEL_CDATA_STATICEVAL = 4, SEL_CDATA_EVALMAX = 8,
  SEL_CDATA_MINMAXALLOC = 16, SEL_CDATA_DOMINMAX = 256, SEL_CDATA_SIMPLESUBEXPR = 32, SEL_CDATA_STATICMULTIEVALSUBEXPR = 64,
  SEL_CDATA_COMMONSUBEXPR = 128
}
 Compiler flags. More...
 

Functions

static void print_group_info (FILE *fp, const char *name, const SelectionTreeElement &sel, gmx_ana_index_t *g)
 Helper method for printing out debug information about a min/max group.
 
void _gmx_selelem_print_compiler_info (FILE *fp, const SelectionTreeElement &sel, int level)
 Prints a human-readable version of the internal compiler data structure. More...
 
static void alloc_selection_data (const SelectionTreeElementPointer &sel, int isize, bool bChildEval)
 Allocates memory for storing the evaluated value of a selection element. More...
 
static void alloc_selection_pos_data (const SelectionTreeElementPointer &sel)
 Allocates memory for storing the evaluated value of a selection element. More...
 
static void set_evaluation_function (const SelectionTreeElementPointer &sel, gmx::sel_evalfunc eval)
 Replace the evaluation function of each element in the subtree. More...
 
static void init_pos_keyword_defaults (SelectionTreeElement *root, const char *spost, const char *rpost, const gmx::internal::SelectionData *sel)
 Initializes default values for position keyword evaluation. More...
 
static SelectionTreeElementPointer reverse_selelem_chain (const SelectionTreeElementPointer &root)
 Reverses the chain of selection elements starting at root. More...
 
static SelectionTreeElementPointer remove_unused_subexpressions (SelectionTreeElementPointer root)
 Removes subexpressions that don't have any references. More...
 
static void create_subexpression_name (const SelectionTreeElementPointer &sel, int i)
 Creates a name with a running number for a subexpression. More...
 
static SelectionTreeElementPointer extract_item_subselections (const SelectionTreeElementPointer &sel, int *subexprn)
 Processes and extracts subexpressions from a given selection subtree. More...
 
static SelectionTreeElementPointer extract_subexpressions (SelectionTreeElementPointer sel)
 Extracts subexpressions of the selection chain. More...
 
static void optimize_boolean_expressions (const SelectionTreeElementPointer &sel)
 Removes redundant boolean selection elements. More...
 
static void reorder_boolean_static_children (const SelectionTreeElementPointer &sel)
 Reorders children of boolean expressions such that static selections come first. More...
 
static void optimize_arithmetic_expressions (const SelectionTreeElementPointer &sel)
 Processes arithmetic expressions to simplify and speed up evaluation. More...
 
static void init_item_evalfunc (const SelectionTreeElementPointer &sel)
 Sets the evaluation functions for the selection (sub)tree. More...
 
static void setup_memory_pooling (const SelectionTreeElementPointer &sel, gmx_sel_mempool_t *mempool)
 Sets the memory pool for selection elements that can use it. More...
 
static void init_item_evaloutput (const SelectionTreeElementPointer &sel)
 Prepares the selection (sub)tree for evaluation. More...
 
static void init_item_compilerdata (const SelectionTreeElementPointer &sel)
 Allocates memory for the compiler data and initializes the structure. More...
 
static void init_item_staticeval (const SelectionTreeElementPointer &sel)
 Initializes the static evaluation flag for a selection subtree. More...
 
static void init_item_subexpr_refcount (const SelectionTreeElementPointer &sel)
 Compute reference counts for subexpressions. More...
 
static void init_item_subexpr_flags (const SelectionTreeElementPointer &sel)
 Initializes compiler flags for subexpressions. More...
 
static void init_item_minmax_groups (const SelectionTreeElementPointer &sel)
 Initializes the gmin and gmax fields of the compiler data structure. More...
 
static void initialize_evalgrps (gmx_ana_selcollection_t *sc)
 Initializes evaluation groups for root items. More...
 
static void mark_subexpr_dynamic (const SelectionTreeElementPointer &sel, bool bDynamic)
 Marks a subtree completely dynamic or undoes such a change. More...
 
static void release_subexpr_memory (const SelectionTreeElementPointer &sel)
 Frees memory for subexpressions that are no longer needed. More...
 
static void make_static (const SelectionTreeElementPointer &sel)
 Makes an evaluated selection element static. More...
 
static void process_const (gmx_sel_evaluate_t *data, const SelectionTreeElementPointer &sel, gmx_ana_index_t *g)
 Evaluates a constant expression during analyze_static(). More...
 
static void store_param_val (const SelectionTreeElementPointer &sel)
 Sets the parameter value pointer for SEL_SUBEXPRREF params. More...
 
static void init_method (const SelectionTreeElementPointer &sel, const gmx_mtop_t *top, int isize)
 Handles the initialization of a selection method during analyze_static() pass. More...
 
static void evaluate_boolean_static_part (gmx_sel_evaluate_t *data, const SelectionTreeElementPointer &sel, gmx_ana_index_t *g)
 Evaluates the static part of a boolean expression. More...
 
static void evaluate_boolean_minmax_grps (const SelectionTreeElementPointer &sel, gmx_ana_index_t *g, gmx_ana_index_t *gmin, gmx_ana_index_t *gmax)
 Evaluates the minimum and maximum groups for a boolean expression. More...
 
static void analyze_static (gmx_sel_evaluate_t *data, const SelectionTreeElementPointer &sel, gmx_ana_index_t *g)
 Evaluates the static parts of sel and analyzes the structure. More...
 
static void init_root_item (const SelectionTreeElementPointer &root, gmx_ana_index_t *gall)
 Initializes the evaluation group for a SEL_ROOT element. More...
 
static void init_required_atoms (const SelectionTreeElementPointer &sel, gmx_ana_index_t *requiredAtoms)
 Finds the highest atom index required to evaluate a selection subtree. More...
 
static void postprocess_item_subexpressions (const SelectionTreeElementPointer &sel)
 Optimizes subexpression evaluation. More...
 
static void init_item_comg (const SelectionTreeElementPointer &sel, gmx::PositionCalculationCollection *pcc, e_poscalc_t type, int flags)
 Initializes COM/COG calculation for method expressions that require it. More...
 
static void free_item_compilerdata (const SelectionTreeElementPointer &sel)
 Frees the allocated compiler data recursively. More...
 

Typedef Documentation

Internal data structure used by the compiler.

Enumeration Type Documentation

anonymous enum

Compiler flags.

Enumerator
SEL_CDATA_FULLEVAL 

Whether a subexpression needs to evaluated for all atoms.

This flag is set for SEL_SUBEXPR elements that are used to evaluate non-atom-valued selection method parameters, as well as those that are used directly as values of selections.

SEL_CDATA_STATIC 

Whether the whole subexpression should be treated as static.

This flag is always false if SEL_DYNAMIC is set for the element, but it is also false for static elements within common subexpressions.

SEL_CDATA_STATICEVAL 

Whether the subexpression will always be evaluated in the same group.

SEL_CDATA_EVALMAX 

Whether the compiler evaluation routine should return the maximal selection.

SEL_CDATA_MINMAXALLOC 

Whether memory has been allocated for gmin and gmax.

SEL_CDATA_DOMINMAX 

Whether to update gmin and gmax in static analysis.

SEL_CDATA_SIMPLESUBEXPR 

Whether the subexpression uses simple pass evaluation functions.

SEL_CDATA_STATICMULTIEVALSUBEXPR 

Whether a static subexpression needs to support multiple evaluations.

This flag may only be set on SEL_SUBEXPR elements that also have SEL_CDATA_SIMPLESUBEXPR.

SEL_CDATA_COMMONSUBEXPR 

Whether this expression is a part of a common subexpression.

Function Documentation

void _gmx_selelem_print_compiler_info ( FILE *  fp,
const SelectionTreeElement sel,
int  level 
)

Prints a human-readable version of the internal compiler data structure.

Parameters
[in]fpFile handle to receive the output.
[in]selSelection element to print.
[in]levelIndentation level, starting from zero.
static void alloc_selection_data ( const SelectionTreeElementPointer &  sel,
int  isize,
bool  bChildEval 
)
static

Allocates memory for storing the evaluated value of a selection element.

Parameters
selSelection element to initialize
[in]isizeMaximum evaluation group size.
[in]bChildEvaltrue if children have already been processed.

If called more than once, memory is (re)allocated to ensure that the maximum of the isize values can be stored.

Allocation of POS_VALUE selection elements is a special case, and is handled by alloc_selection_pos_data().

static void alloc_selection_pos_data ( const SelectionTreeElementPointer &  sel)
static

Allocates memory for storing the evaluated value of a selection element.

Parameters
selSelection element to initialize.

Allocation of POS_VALUE selection elements is a special case, and is handled by this function instead of by alloc_selection_data().

static void analyze_static ( gmx_sel_evaluate_t data,
const SelectionTreeElementPointer &  sel,
gmx_ana_index_t g 
)
static

Evaluates the static parts of sel and analyzes the structure.

Parameters
[in]dataEvaluation data.
[in,out]selSelection currently being evaluated.
[in]gGroup for which sel should be evaluated.
Returns
0 on success, a non-zero error code on error.

This function is used as the replacement for the gmx::SelectionTreeElement::evaluate function pointer. It does the single most complex task in the compiler: after all elements have been processed, the gmin and gmax fields of t_compiler_data have been properly initialized, enough memory has been allocated for storing the value of each expression, and the static parts of the expressions have been evaluated. The above is exactly true only for elements other than subexpressions: another pass is required for subexpressions that are referred to more than once and whose evaluation group is not known in advance.

static void create_subexpression_name ( const SelectionTreeElementPointer &  sel,
int  i 
)
static

Creates a name with a running number for a subexpression.

Parameters
[in,out]selThe subexpression to be named.
[in]iRunning number for the subexpression.

The name of the selection becomes "SubExpr N", where N is i; Memory is allocated for the name and the name is stored both as the name of the subexpression element and as gmx::SelectionTreeElement::u::cgrp::name; the latter is freed by _gmx_selelem_free().

static void evaluate_boolean_minmax_grps ( const SelectionTreeElementPointer &  sel,
gmx_ana_index_t g,
gmx_ana_index_t gmin,
gmx_ana_index_t gmax 
)
static

Evaluates the minimum and maximum groups for a boolean expression.

Parameters
[in]selSEL_BOOLEAN element currently being evaluated.
[in]gGroup for which sel has been evaluated.
[out]gminLargest subset of the possible values of sel.
[out]gmaxSmallest superset of the possible values of sel.

This is a helper function for analyze_static() that is called for dynamic SEL_BOOLEAN elements after they have been evaluated. It uses the minimum and maximum groups of the children to calculate the minimum and maximum groups for sel, and also updates the static part of sel (which is in the first child) if the children give cause for this.

This function may allocate some extra memory for gmin and gmax, but as these groups are freed at the end of analyze_static() (which is reached shortly after this function returns), this should not be a major problem.

static void evaluate_boolean_static_part ( gmx_sel_evaluate_t data,
const SelectionTreeElementPointer &  sel,
gmx_ana_index_t g 
)
static

Evaluates the static part of a boolean expression.

Parameters
[in]dataEvaluation data.
[in,out]selBoolean selection element whose children should be processed.
[in]gThe evaluation group.
Returns
0 on success, a non-zero error code on error.

reorder_item_static_children() should have been called.

static SelectionTreeElementPointer extract_item_subselections ( const SelectionTreeElementPointer &  sel,
int *  subexprn 
)
static

Processes and extracts subexpressions from a given selection subtree.

Parameters
selRoot of the subtree to process.
subexprnPointer to a subexpression counter.
Returns
Pointer to a chain of subselections, or NULL if none were found.

This function finds recursively all SEL_SUBEXPRREF elements below the given root element and ensures that their children are within SEL_SUBEXPR elements. It also creates a chain of SEL_ROOT elements that contain the subexpression as their children and returns the first of these root elements.

static SelectionTreeElementPointer extract_subexpressions ( SelectionTreeElementPointer  sel)
static

Extracts subexpressions of the selection chain.

Parameters
selFirst selection in the whole selection chain.
Returns
The new first element for the chain.

Finds all the subexpressions (and their subexpressions) in the selection chain starting from sel and creates SEL_SUBEXPR elements for them. SEL_ROOT elements are also created for each subexpression and inserted into the selection chain before the expressions that refer to them.

static void free_item_compilerdata ( const SelectionTreeElementPointer &  sel)
static

Frees the allocated compiler data recursively.

Parameters
selRoot of the selection subtree to process.

Frees the data allocated for the compilation process.

static void init_item_comg ( const SelectionTreeElementPointer &  sel,
gmx::PositionCalculationCollection pcc,
e_poscalc_t  type,
int  flags 
)
static

Initializes COM/COG calculation for method expressions that require it.

Parameters
selSelection subtree to process.
[in,out]pccPosition calculation collection to use.
[in]typeDefault position calculation type.
[in]flagsFlags for default position calculation.

Searches recursively through the selection tree for dynamic SEL_EXPRESSION elements that define the gmx_ana_selmethod_t::pupdate function. For each such element found, position calculation is initialized for the maximal evaluation group. The type of the calculation is determined by type and flags. No calculation is initialized if type equals POS_ATOM and the method also defines the gmx_ana_selmethod_t::update method.

static void init_item_compilerdata ( const SelectionTreeElementPointer &  sel)
static

Allocates memory for the compiler data and initializes the structure.

Parameters
selRoot of the selection subtree to process.
static void init_item_evalfunc ( const SelectionTreeElementPointer &  sel)
static

Sets the evaluation functions for the selection (sub)tree.

Parameters
[in,out]selRoot of the selection subtree to process.

This function sets the evaluation function (gmx::SelectionTreeElement::evaluate) for the selection elements.

static void init_item_evaloutput ( const SelectionTreeElementPointer &  sel)
static

Prepares the selection (sub)tree for evaluation.

Parameters
[in,out]selRoot of the selection subtree to prepare.

It also allocates memory for the sel->v.u.g or sel->v.u.p structure if required.

static void init_item_minmax_groups ( const SelectionTreeElementPointer &  sel)
static

Initializes the gmin and gmax fields of the compiler data structure.

Parameters
selRoot of the selection subtree to process.
static void init_item_staticeval ( const SelectionTreeElementPointer &  sel)
static

Initializes the static evaluation flag for a selection subtree.

Parameters
[in,out]selRoot of the selection subtree to process.

Sets the bStaticEval in the compiler data structure: for any element for which the evaluation group may depend on the trajectory frame, the flag is cleared.

reorder_boolean_static_children() should have been called.

static void init_item_subexpr_flags ( const SelectionTreeElementPointer &  sel)
static

Initializes compiler flags for subexpressions.

Parameters
selRoot of the selection subtree to process.
static void init_item_subexpr_refcount ( const SelectionTreeElementPointer &  sel)
static

Compute reference counts for subexpressions.

Parameters
selRoot of the selection subtree to process.
static void init_method ( const SelectionTreeElementPointer &  sel,
const gmx_mtop_t *  top,
int  isize 
)
static

Handles the initialization of a selection method during analyze_static() pass.

Parameters
[in,out]selSelection element to process.
[in]topTopology structure.
[in]isizeSize of the evaluation group for the element.
Returns
0 on success, a non-zero error code on return.

Calls sel_initfunc() (and possibly sel_outinitfunc()) to initialize the method. If no SPAR_ATOMVAL parameters are present, multiple initialization is prevented by using SEL_METHODINIT and SEL_OUTINIT flags.

static void init_pos_keyword_defaults ( SelectionTreeElement root,
const char *  spost,
const char *  rpost,
const gmx::internal::SelectionData sel 
)
static

Initializes default values for position keyword evaluation.

Parameters
[in,out]rootRoot of the element tree to initialize.
[in]spostDefault output position type.
[in]rpostDefault reference position type.
[in]selSelection that the element evaluates the positions for, or NULL if the element is an internal element.
static void init_required_atoms ( const SelectionTreeElementPointer &  sel,
gmx_ana_index_t requiredAtoms 
)
static

Finds the highest atom index required to evaluate a selection subtree.

Parameters
[in]selRoot of the selection subtree to process.
[in,out]requiredAtomsThe atoms required to evaluate the subtree. The existing group is only added to, so multiple calls with the same parameter will compute the union over several subtrees.

For evaluation that starts from a SEL_ROOT element with a fixed group, children will never extend the evaluation group except for method parameter evaluation (which have their own root element), so it is sufficient to check the root. However, children of SEL_EXPRESSION elements (i.e., the method parameters) may have been independently evaluated to a static group that no longer has a separate root, so those need to be checked as well.

Position calculations are not considered here, but are analyzed through the position calculation collection in the main compilation method.

static void init_root_item ( const SelectionTreeElementPointer &  root,
gmx_ana_index_t gall 
)
static

Initializes the evaluation group for a SEL_ROOT element.

Parameters
rootRoot element to initialize.
[in]gallGroup of all atoms.

Checks whether it is necessary to evaluate anything through the root element, and either clears the evaluation function or initializes the evaluation group.

static void initialize_evalgrps ( gmx_ana_selcollection_t sc)
static

Initializes evaluation groups for root items.

Parameters
[in,out]scSelection collection data.

The evaluation group of each SEL_ROOT element corresponding to a selection in sc is set to NULL. The evaluation group for SEL_ROOT elements corresponding to subexpressions that need full evaluation is set to sc->gall.

static void make_static ( const SelectionTreeElementPointer &  sel)
static

Makes an evaluated selection element static.

Parameters
selSelection element to make static.

The evaluated value becomes the value of the static element. The element type is changed to SEL_CONST and the children are deleted.

static void mark_subexpr_dynamic ( const SelectionTreeElementPointer &  sel,
bool  bDynamic 
)
static

Marks a subtree completely dynamic or undoes such a change.

Parameters
selSelection subtree to mark.
[in]bDynamicIf true, the bStatic flag of the whole selection subtree is cleared. If false, the flag is restored to using SEL_DYNAMIC.

Does not descend into parameters of methods unless the parameters are evaluated for each atom.

static void optimize_arithmetic_expressions ( const SelectionTreeElementPointer &  sel)
static

Processes arithmetic expressions to simplify and speed up evaluation.

Parameters
selRoot of the selection subtree to process.

Currently, this function only converts integer constants to reals within arithmetic expressions.

static void optimize_boolean_expressions ( const SelectionTreeElementPointer &  sel)
static

Removes redundant boolean selection elements.

Parameters
selRoot of the selection subtree to optimize.

This function merges similar boolean operations (e.g., (A or B) or C becomes a single OR operation with three operands).

static void postprocess_item_subexpressions ( const SelectionTreeElementPointer &  sel)
static

Optimizes subexpression evaluation.

Parameters
selRoot of the selection subtree to process.

Optimizes away some unnecessary evaluation of subexpressions that are only referenced once.

static void process_const ( gmx_sel_evaluate_t data,
const SelectionTreeElementPointer &  sel,
gmx_ana_index_t g 
)
static

Evaluates a constant expression during analyze_static().

Parameters
[in]dataEvaluation data.
[in,out]selSelection to process.
[in]gThe evaluation group.
Returns
0 on success, a non-zero error code on error.
static void release_subexpr_memory ( const SelectionTreeElementPointer &  sel)
static

Frees memory for subexpressions that are no longer needed.

Parameters
selSelection subtree to check.

Checks whether the subtree rooted at sel refers to any SEL_SUBEXPR elements that are not referred to by anything else except their own root element. If such elements are found, all memory allocated for them is freed except the actual element. The element is left because otherwise a dangling pointer would be left at the root element, which is not traversed by this function. Later compilation passes remove the stub elements.

static SelectionTreeElementPointer remove_unused_subexpressions ( SelectionTreeElementPointer  root)
static

Removes subexpressions that don't have any references.

Parameters
rootFirst selection in the whole selection chain.
Returns
The new first element for the chain.

The elements are processed in reverse order to correctly detect subexpressions only referred to by other subexpressions.

static void reorder_boolean_static_children ( const SelectionTreeElementPointer &  sel)
static

Reorders children of boolean expressions such that static selections come first.

Parameters
selRoot of the selection subtree to reorder.

The relative order of static expressions does not change. The same is true for the dynamic expressions.

static SelectionTreeElementPointer reverse_selelem_chain ( const SelectionTreeElementPointer &  root)
static

Reverses the chain of selection elements starting at root.

Parameters
rootFirst selection in the whole selection chain.
Returns
The new first element for the chain.
static void set_evaluation_function ( const SelectionTreeElementPointer &  sel,
gmx::sel_evalfunc  eval 
)
static

Replace the evaluation function of each element in the subtree.

Parameters
selRoot of the selection subtree to process.
[in]evalThe new evaluation function.
static void setup_memory_pooling ( const SelectionTreeElementPointer &  sel,
gmx_sel_mempool_t mempool 
)
static

Sets the memory pool for selection elements that can use it.

Parameters
selRoot of the selection subtree to process.
[in]mempoolMemory pool to use.
static void store_param_val ( const SelectionTreeElementPointer &  sel)
static

Sets the parameter value pointer for SEL_SUBEXPRREF params.

Parameters
[in,out]selSelection to process.

Copies the value pointer of sel to sel->u.param if one is present and should receive the value from the compiler (most parameter values are handled during parsing). If sel is not of type SEL_SUBEXPRREF, or if sel->u.param is NULL, the function does nothing. Also, if the sel->u.param does not have SPAR_VARNUM or SPAR_ATOMVAL, the function returns immediately.