Gromacs  2022.2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Typedefs | Functions
#include "selparam.h"
#include "selvalue.h"
+ Include dependency graph for selmethod.h:
+ This graph shows which files directly or indirectly include this file:

Description

API for handling selection methods.

There should be no need to use the data structures or call the functions in this file directly unless implementing a custom selection method.

Instructions for implementing custom selection methods can be found on a separate page: Custom selection methods

Author
Teemu Murtola teemu.nosp@m..mur.nosp@m.tola@.nosp@m.gmai.nosp@m.l.com

Classes

struct  gmx::SelMethodEvalContext
 Evaluation context for selection methods. More...
 
struct  gmx_ana_selmethod_help_t
 Help information for a selection method. More...
 
struct  gmx_ana_selmethod_t
 Describes a selection method. More...
 

Macros

Selection method flags

#define SMETH_REQTOP   1
 If set, the method requires topology information.
 
#define SMETH_REQMASS   2
 If set, the method requires atom masses.
 
#define SMETH_DYNAMIC   4
 If set, the method can only be evaluated dynamically.
 
#define SMETH_SINGLEVAL   8
 If set, the method evaluates to a single value. More...
 
#define SMETH_VARNUMVAL   16
 If set, the method evaluates to an arbitrary number of values. More...
 
#define SMETH_CHARVAL   64
 If set, the method evaluates to single-character strings. More...
 
#define SMETH_ALLOW_UNSORTED   128
 If set, the method accepts unsorted atoms in its input parameters. More...
 
#define SMETH_MODIFIER   256
 If set, the method is a selection modifier. More...
 

Typedefs

typedef void *(* sel_datafunc )(int npar, gmx_ana_selparam_t *param)
 Allocates and initializes internal data and parameter values. More...
 
typedef void(* sel_posfunc )(gmx::PositionCalculationCollection *pcc, void *data)
 Sets the position calculation collection for the method. More...
 
typedef void(* sel_initfunc )(const gmx_mtop_t *top, int npar, gmx_ana_selparam_t *param, void *data)
 Does initialization based on topology and/or parameter values. More...
 
typedef void(* sel_outinitfunc )(const gmx_mtop_t *top, gmx_ana_selvalue_t *out, void *data)
 Initializes output data structure. More...
 
typedef void(* sel_freefunc )(void *data)
 Frees the internal data. More...
 
typedef void(* sel_framefunc )(const gmx::SelMethodEvalContext &context, void *data)
 Initializes the evaluation for a new frame. More...
 
typedef void(* sel_updatefunc )(const gmx::SelMethodEvalContext &context, gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)
 Evaluates a selection method. More...
 
typedef void(* sel_updatefunc_pos )(const gmx::SelMethodEvalContext &context, gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
 Evaluates a selection method using positions. More...
 

Functions

int gmx_ana_selmethod_register (gmx::SelectionParserSymbolTable *symtab, const char *name, gmx_ana_selmethod_t *method)
 Registers a selection method. More...
 
int gmx_ana_selmethod_register_defaults (gmx::SelectionParserSymbolTable *symtab)
 Registers all selection methods in the library. More...
 
gmx_ana_selparam_tgmx_ana_selmethod_find_param (const char *name, gmx_ana_selmethod_t *method)
 Finds a parameter from a selection method by name. More...
 

Macro Definition Documentation

#define SMETH_ALLOW_UNSORTED   128

If set, the method accepts unsorted atoms in its input parameters.

Currently, the support for this functionality is fairly limited, and only static index group references can actually contain unsorted atoms. But to make this single case work, the position evaluation must support unsorted atoms as well.

#define SMETH_CHARVAL   64

If set, the method evaluates to single-character strings.

This flag can only be set for STR_VALUE methods. If it is set, the selection engine automatically allocates and frees the required strings. The evaluation function should store the character values as the first character in the strings in the output data structure and should not change the string pointers.

#define SMETH_MODIFIER   256

If set, the method is a selection modifier.

The method type should be GROUP_VALUE or NO_VALUE . Cannot be combined with SMETH_SINGLEVAL or SMETH_VARNUMVAL .

#define SMETH_SINGLEVAL   8

If set, the method evaluates to a single value.

The default is that the method evaluates to a value for each input atom. Cannot be combined with SMETH_VARNUMVAL.

#define SMETH_VARNUMVAL   16

If set, the method evaluates to an arbitrary number of values.

The default is that the method evaluates to a value for each input atom. Cannot be combined with SMETH_SINGLEVAL or with GROUP_VALUE.

Typedef Documentation

typedef void*(* sel_datafunc)(int npar, gmx_ana_selparam_t *param)

Allocates and initializes internal data and parameter values.

Parameters
[in]nparNumber of parameters in param.
[in,out]paramPointer to (a copy of) the method's gmx_ana_selmethod_t::param.
Returns
Pointer to method-specific data structure. This pointer will be passed as the last parameter of all other function calls.
Exceptions
unspecifiedAny errors should be indicated by throwing an exception.

Should allocate and initialize any internal data required by the method. Should also initialize the value pointers (gmx_ana_selparam_t::val) in param to point to variables within the internal data structure, with the exception of parameters that specify the SPAR_VARNUM or the SPAR_ATOMVAL flag (these should be handled in sel_initfunc()). However, parameters with a position value should be initialized. It is also possible to initialize SPAR_ENUMVAL statically outside this function (see SPAR_ENUMVAL). The gmx_ana_selparam_t::nvalptr should also be initialized for non-position-valued parameters that have both SPAR_VARNUM and SPAR_DYNAMIC set (it can also be initialized for other parameters if desired, but the same information will be available through other means). For optional parameters, the default values can (and should) be initialized here, as the parameter values are not changed if the parameter is not provided.

For boolean parameters (type equals NO_VALUE), the default value should be set here. The user can override the value by giving the parameter either as 'NAME'/'noNAME', or as 'NAME on/off/yes/no'.

If the method takes any parameters, this function must be provided.

typedef void(* sel_framefunc)(const gmx::SelMethodEvalContext &context, void *data)

Initializes the evaluation for a new frame.

Parameters
[in]contextEvaluation context.
dataInternal data structure from sel_datafunc().
Returns
0 on success, a non-zero error code on failure.

This function should be implemented if the selection method needs to do some preprocessing for each frame, and the preprocessing does not depend on the evaluation group. Because sel_updatefunc_* can be called more than once for a frame, it is inefficient do the preprocessing there. It is ensured that this function will be called before sel_updatefunc_* for each frame, and that it will be called at most once for each frame. For static methods, it is called once, with fr and pbc set to NULL.

typedef void(* sel_freefunc)(void *data)

Frees the internal data.

Parameters
[in]dataInternal data structure from sel_datafunc().

This function should be provided if the internal data structure contains dynamically allocated data, and should free any such data. The data structure itself should also be freed. For convenience, if there is no dynamically allocated data within the structure and the structure is allocated using malloc()/snew(), this function is not needed: the selection engine automatically frees the structure using sfree(). Any memory pointers received as values of parameters are managed externally, and should not be freed. Pointers set as the value pointer of SPAR_ENUMVAL parameters should not be freed.

typedef void(* sel_initfunc)(const gmx_mtop_t *top, int npar, gmx_ana_selparam_t *param, void *data)

Does initialization based on topology and/or parameter values.

Parameters
[in]topTopology structure (can be NULL if SMETH_REQTOP or SMETH_REQMASS is not set).
[in]nparNumber of parameters in param.
[in]paramPointer to (an initialized copy of) the method's gmx_ana_selmethod_t::param.
dataInternal data structure from sel_datafunc().
Returns
0 on success, a non-zero error code on failure.

This function is called after the parameters have been processed: the values of the parameters are stored at the locations set in sel_datafunc(). The flags SPAR_DYNAMIC and SPAR_ATOMVAL are cleared before calling the function if the value is static or single-valued, respectively. If a parameter had the SPAR_VARNUM or SPAR_ATOMVAL flag (and is not POS_VALUE), a pointer to the memory allocated for the values is found in gmx_ana_selparam_t::val. The pointer should be stored by this function, otherwise the values cannot be accessed. For SPAR_VARNUM parameters, the number of values can be accessed through gmx_ana_selparam_t::val. For parameters with SPAR_DYNAMIC, the number is the maximum number of values (the actual number can be accessed in sel_framefunc() and in the update callback through the value pointed by gmx_ana_selparam_t::nvalptr). For SPAR_ATOMVAL parameters, gmx_ana_selparam_t::val::nr is set to 1 if a single value was provided, otherwise it is set to the maximum number of values possibly passed to the method. The value pointed by gmx_ana_selparam_t::nvalptr always contains the same value as gmx_ana_selparam_t::val::nr.

For dynamic GROUP_VALUE parameters (SPAR_DYNAMIC set), the value will be the largest possible selection that may occur during the evaluation. For other types of dynamic parameters, the values are undefined.

If the method takes any parameters with the SPAR_VARNUM or SPAR_ATOMVAL flags, this function must be provided, except if these parameters all have POS_VALUE.

This function may be called multiple times for the same method if the method takes parameters with SPAR_ATOMVAL set.

typedef void(* sel_outinitfunc)(const gmx_mtop_t *top, gmx_ana_selvalue_t *out, void *data)

Initializes output data structure.

Parameters
[in]topTopology structure (can be NULL if SMETH_REQTOP or SMETH_REQMASS is not set).
[in,out]outOutput data structure.
[in]dataInternal data structure from sel_datafunc().
Returns
0 on success, an error code on error.

This function is called immediately after sel_initfunc().

If the method evaluates to a position (POS_VALUE), this function should be provided, and it should initialize the gmx_ana_pos_t data structure pointed by out.p (the pointer is guaranteed to be non-NULL). The out.p->g pointer should be initialized to the group that is used to evaluate positions in sel_updatefunc() or sel_updatefunc_pos().

The function should also be provided for non-position-valued SMETH_VARNUMVAL methods. For these methods, it suffices to set the out->nr field to reflect the maximum number of values returned by the method.

Currently, this function is not needed for other types of methods.

This function may be called multiple times for the same method if the method takes parameters with SPAR_ATOMVAL set.

typedef void(* sel_posfunc)(gmx::PositionCalculationCollection *pcc, void *data)

Sets the position calculation collection for the method.

Parameters
[in]pccPosition calculation collection that the method should use for position calculations.
dataInternal data structure from sel_datafunc().

This function should be provided if the method uses the routines from poscalc.h for calculating positions. The pointer pcc should then be stored and used for initialization for any position calculation structures.

typedef void(* sel_updatefunc)(const gmx::SelMethodEvalContext &context, gmx_ana_index_t *g, gmx_ana_selvalue_t *out, void *data)

Evaluates a selection method.

Parameters
[in]contextEvaluation context.
[in]gIndex group for which the method should be evaluated.
[out]outOutput data structure.
dataInternal data structure from sel_datafunc().
Returns
0 on success, a non-zero error code on error.

This function should evaluate the method for each atom included in g, and write the output to out. The pointer in the union out->u that corresponds to the type of the method should be used. Enough memory has been allocated to store the output values. The number of values in out should also be updated if necessary. However, POS_VALUE or GROUP_VALUE methods should not touch out->nr (it should be 1 anyways).

For STR_VALUE methods, the pointers stored in out->s are discarded without freeing; it is the responsibility of this function to provide pointers that can be discarded without memory leaks.

If the method accesses fr outside the index group specified in g or what it receives from its parameters, it must check that fr actually contains such an atom in case the fr has been loaded from a trajectory that only contains a subset of the system.

typedef void(* sel_updatefunc_pos)(const gmx::SelMethodEvalContext &context, gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)

Evaluates a selection method using positions.

Parameters
[in]contextEvaluation context.
[in]posPositions for which the method should be evaluated.
[out]outOutput data structure.
dataInternal data structure from sel_datafunc().
Returns
0 on success, a non-zero error code on error.

This function should evaluate the method for each position in pos, and write the output values to out. The pointer in the union out->u that corresponds to the type of the method should be used. Enough memory has been allocated to store the output values. The number of values in out should also be updated if necessary. POS_VALUE or GROUP_VALUE methods should not touch out->nr (it should be 1 anyways). For other types of methods, the number of output values should equal the number of positions in pos.

For STR_VALUE methods, the pointers stored in out->s are discarded without freeing; it is the responsibility of this function to provide pointers that can be discarded without memory leaks.

If the method accesses fr outside the atoms referenced in pos or what it receives from its parameters, it must check that fr actually contains such an atom in case the fr has been loaded from a trajectory that only contains a subset of the system.

Function Documentation

gmx_ana_selparam_t* gmx_ana_selmethod_find_param ( const char *  name,
gmx_ana_selmethod_t method 
)

Finds a parameter from a selection method by name.

Parameters
[in]nameName of the parameter to search.
[in]methodMethod to search for the parameter.
Returns
Pointer to the parameter in the method->param array, or NULL if no parameter with name name was found.

This is a simple wrapper for gmx_ana_selparam_find().

int gmx_ana_selmethod_register ( gmx::SelectionParserSymbolTable symtab,
const char *  name,
gmx_ana_selmethod_t method 
)

Registers a selection method.

Parameters
[in,out]symtabSymbol table to register the method to.
[in]nameName under which the method should be registered.
[in]methodMethod to register.
Returns
0 on success, -1 if there was something wrong with the method.

name does not need to match the name of the method, and the same method can be registered multiple times under different names. If name equals some previously registered name, an error message is printed and the method is not registered.

The function also performs some sanity checking on the input method, and refuses to register it if there are problems. Some problems only generate warnings. All problems are described to stderr.

int gmx_ana_selmethod_register_defaults ( gmx::SelectionParserSymbolTable symtab)

Registers all selection methods in the library.

Parameters
[in,out]symtabSymbol table to register the methods to.
Returns
0 on success, -1 if any of the default methods could not be registered.