Gromacs  2024.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Typedefs | Functions | Variables
#include "gmxpre.h"
#include <cmath>
#include <algorithm>
#include "gromacs/math/functions.h"
#include "gromacs/math/units.h"
#include "gromacs/math/utilities.h"
#include "gromacs/math/vec.h"
#include "gromacs/pbcutil/pbc.h"
#include "gromacs/selection/indexutil.h"
#include "gromacs/selection/position.h"
#include "gromacs/selection/selection.h"
#include "gromacs/utility/arraysize.h"
#include "gromacs/utility/exceptions.h"
#include "gromacs/utility/smalloc.h"
#include "selelem.h"
#include "selmethod.h"
#include "selmethod_impl.h"
+ Include dependency graph for sm_insolidangle.cpp:

Description

Implements the insolidangle selection method.

Todo:
The implementation could be optimized quite a bit.
Todo:
Move the covered fraction stuff somewhere else and make it more generic (along the lines it is handled in selection.h and trajana.h in the old C API).
Author
Teemu Murtola teemu.nosp@m..mur.nosp@m.tola@.nosp@m.gmai.nosp@m.l.com

Classes

struct  t_partition_item
 Internal data structure for the insolidangle selection method. More...
 
struct  partition
 Internal data structure for the insolidangle selection method. More...
 
struct  spheresurfacebin
 Internal data structure for the insolidangle selection method. More...
 
struct  methoddata_insolidangle
 Data structure for the insolidangle selection method. More...
 

Typedefs

typedef struct partition t_partition
 Internal data structure for the insolidangle selection method. More...
 
typedef struct spheresurfacebin t_spheresurfacebin
 Internal data structure for the insolidangle selection method. More...
 
typedef struct
methoddata_insolidangle 
t_methoddata_insolidangle
 Data structure for the insolidangle selection method. More...
 

Functions

static void * init_data_insolidangle (int npar, gmx_ana_selparam_t *param)
 Allocates data for the insolidangle selection method. More...
 
static void init_insolidangle (const gmx_mtop_t *top, int npar, gmx_ana_selparam_t *param, void *data)
 Initializes the insolidangle selection method. More...
 
static void free_data_insolidangle (void *data)
 Frees the data allocated for the insolidangle selection method. More...
 
static void init_frame_insolidangle (const gmx::SelMethodEvalContext &context, void *data)
 Initializes the evaluation of the insolidangle selection method for a frame. More...
 
static bool accept_insolidangle (rvec x, const t_pbc *pbc, void *data)
 Internal helper function for evaluate_insolidangle(). More...
 
static void evaluate_insolidangle (const gmx::SelMethodEvalContext &context, gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
 Evaluates the insolidangle selection method. More...
 
static real sph_distc (rvec x1, rvec x2)
 Calculates the distance between unit vectors. More...
 
static int find_partition_bin (t_partition *p, real value)
 Does a binary search on a t_partition to find a bin for a value. More...
 
static int find_surface_bin (t_methoddata_insolidangle *surf, rvec x)
 Finds a bin that corresponds to a location on the unit sphere surface. More...
 
static void clear_surface_points (t_methoddata_insolidangle *surf)
 Clears/initializes the bins on the unit sphere surface. More...
 
static void free_surface_points (t_methoddata_insolidangle *surf)
 Frees memory allocated for storing the reference points in the surface bins. More...
 
static void add_surface_point (t_methoddata_insolidangle *surf, int tbin, int pbin, rvec x)
 Adds a reference point to a given bin. More...
 
static void mark_surface_covered (t_methoddata_insolidangle *surf, int tbin, int pbin)
 Marks a bin as completely covered. More...
 
static void update_surface_bin (t_methoddata_insolidangle *surf, int tbin, real phi, real pdelta1, real pdelta2, real pdeltamax, rvec x)
 Helper function for store_surface_point() to update a single zenith angle bin. More...
 
static void store_surface_point (t_methoddata_insolidangle *surf, rvec x)
 Adds a single reference point and updates the surface bins. More...
 
static void optimize_surface_points (t_methoddata_insolidangle *surf)
 Optimizes the surface bins for faster searching. More...
 
static real estimate_covered_fraction (t_methoddata_insolidangle *surf)
 Estimates the area covered by the reference cones. More...
 
static bool is_surface_covered (t_methoddata_insolidangle *surf, rvec x)
 Checks whether a point lies within a solid angle. More...
 
bool _gmx_selelem_can_estimate_cover (const gmx::SelectionTreeElement &sel)
 Returns true if the covered fraction of the selection can be calculated. More...
 
real _gmx_selelem_estimate_coverfrac (const gmx::SelectionTreeElement &sel)
 Returns the covered fraction of the selection for the current frame. More...
 

Variables

static gmx_ana_selparam_t smparams_insolidangle []
 Parameters for the insolidangle selection method. More...
 
static const char *const help_insolidangle []
 Help text for the insolidangle selection method. More...
 
gmx_ana_selmethod_t sm_insolidangle
 Selection method data for the insolidangle method. More...
 

Function Documentation

bool _gmx_selelem_can_estimate_cover ( const gmx::SelectionTreeElement sel)

Returns true if the covered fraction of the selection can be calculated.

Parameters
[in]selSelection element to query.
Returns
true if the covered fraction can be estimated for sel with _gmx_selelem_estimate_coverfrac(), false otherwise.
real _gmx_selelem_estimate_coverfrac ( const gmx::SelectionTreeElement sel)

Returns the covered fraction of the selection for the current frame.

Parameters
[in]selSelection for which the fraction should be calculated.
Returns
Fraction of angles covered by the selection (between zero and one).

The return value is undefined if _gmx_selelem_can_estimate_cover() returns false. Should be called after gmx_ana_evaluate_selections() has been called for the frame.

static bool accept_insolidangle ( rvec  x,
const t_pbc pbc,
void *  data 
)
static

Internal helper function for evaluate_insolidangle().

Parameters
[in]xTest point.
[in]pbcPBC data (if NULL, no PBC are used).
[in]dataPointer to a t_methoddata_insolidangle data structure.
Returns
true if x is within the solid angle, false otherwise.
static void add_surface_point ( t_methoddata_insolidangle surf,
int  tbin,
int  pbin,
rvec  x 
)
static

Adds a reference point to a given bin.

Parameters
[in,out]surfSurface data structure.
[in]tbinBin number in the zenith angle direction.
[in]pbinBin number in the azimuthal angle direction.
[in]xPoint to store.
static void clear_surface_points ( t_methoddata_insolidangle surf)
static

Clears/initializes the bins on the unit sphere surface.

Parameters
[in,out]surfSurface data structure.

Clears the reference points from the bins and (re)initializes the edges of the azimuthal bins.

static real estimate_covered_fraction ( t_methoddata_insolidangle surf)
static

Estimates the area covered by the reference cones.

Parameters
[in]surfSurface data structure.
Returns
An estimate for the area covered by the reference points.
static void evaluate_insolidangle ( const gmx::SelMethodEvalContext context,
gmx_ana_pos_t pos,
gmx_ana_selvalue_t out,
void *  data 
)
static

Evaluates the insolidangle selection method.

See sel_updatefunc() for description of the parameters. data should point to a t_methoddata_insolidangle.

Calculates which atoms in g are within the solid angle spanned by t_methoddata_insolidangle::span and centered at t_methoddata_insolidangle::center, and stores the result in out->u.g.

static int find_partition_bin ( t_partition p,
real  value 
)
static

Does a binary search on a t_partition to find a bin for a value.

Parameters
[in]pPartition to search.
[in]valueValue to search for.
Returns
The partition index in p that contains value.

If value is outside the range of p, the first/last index is returned. Otherwise, the return value i satisfies p->p[i].left<=value and p->p[i+1].left>value

static int find_surface_bin ( t_methoddata_insolidangle surf,
rvec  x 
)
static

Finds a bin that corresponds to a location on the unit sphere surface.

Parameters
[in]surfSurface data structure to search.
[in]xUnit vector to find.
Returns
The bin index that contains x.

The return value is an index to the surf->bin array.

static void free_data_insolidangle ( void *  data)
static

Frees the data allocated for the insolidangle selection method.

Parameters
dataData to free (should point to a t_methoddata_insolidangle).

Frees the memory allocated for t_methoddata_insolidangle::center and t_methoddata_insolidangle::span, as well as the memory for the internal bin structure.

static void free_surface_points ( t_methoddata_insolidangle surf)
static

Frees memory allocated for storing the reference points in the surface bins.

Parameters
[in,out]surfSurface data structure.
static void * init_data_insolidangle ( int  npar,
gmx_ana_selparam_t param 
)
static

Allocates data for the insolidangle selection method.

Parameters
[in]nparNot used (should be 3).
[in,out]paramMethod parameters (should point to smparams_insolidangle).
Returns
Pointer to the allocated data (t_methoddata_insolidangle).

Allocates memory for a t_methoddata_insolidangle structure and initializes the parameter as follows:

static void init_frame_insolidangle ( const gmx::SelMethodEvalContext context,
void *  data 
)
static

Initializes the evaluation of the insolidangle selection method for a frame.

Parameters
[in]contextEvaluation context.
dataShould point to a t_methoddata_insolidangle.

Creates a lookup structure that enables fast queries of whether a point is within the solid angle or not.

static void init_insolidangle ( const gmx_mtop_t *  top,
int  npar,
gmx_ana_selparam_t param,
void *  data 
)
static

Initializes the insolidangle selection method.

Parameters
topNot used.
nparNot used.
paramNot used.
dataPointer to t_methoddata_insolidangle to initialize.
Returns
0 on success, -1 on failure.

Converts t_methoddata_insolidangle::angcut to radians and allocates and allocates memory for the bins used during the evaluation.

static bool is_surface_covered ( t_methoddata_insolidangle surf,
rvec  x 
)
static

Checks whether a point lies within a solid angle.

Parameters
[in]surfSurface data structure to search.
[in]xUnit vector to check.
Returns
true if x is within the solid angle, false otherwise.
static void mark_surface_covered ( t_methoddata_insolidangle surf,
int  tbin,
int  pbin 
)
static

Marks a bin as completely covered.

Parameters
[in,out]surfSurface data structure.
[in]tbinBin number in the zenith angle direction.
[in]pbinBin number in the azimuthal angle direction.
static void optimize_surface_points ( t_methoddata_insolidangle surf)
static

Optimizes the surface bins for faster searching.

Parameters
[in,out]surfSurface data structure.

Currently, this function does nothing.

static real sph_distc ( rvec  x1,
rvec  x2 
)
static

Calculates the distance between unit vectors.

Parameters
[in]x1Unit vector 1.
[in]x2Unit vector 2.
Returns
Minus the dot product of x1 and x2.

This function is used internally to calculate the distance between the unit vectors x1 and x2 to find out whether x2 is within the cone centered at x1. Currently, the cosine of the angle is used for efficiency, and the minus is there to make it behave like a normal distance (larger values mean longer distances).

static void store_surface_point ( t_methoddata_insolidangle surf,
rvec  x 
)
static

Adds a single reference point and updates the surface bins.

Parameters
[in,out]surfSurface data structure.
[in]xPoint to store (should have unit length).

Finds all the bins covered by the cone centered at x and calls update_surface_bin() to update them.

static void update_surface_bin ( t_methoddata_insolidangle surf,
int  tbin,
real  phi,
real  pdelta1,
real  pdelta2,
real  pdeltamax,
rvec  x 
)
static

Helper function for store_surface_point() to update a single zenith angle bin.

Parameters
[in,out]surfSurface data structure.
[in]tbinBin number in the zenith angle direction.
[in]phiAzimuthal angle of x.
[in]pdelta1Width of the cone at the lower edge of tbin.
[in]pdelta2Width of the cone at the uppper edge of tbin.
[in]pdeltamaxMax. width of the cone inside tbin.
[in]xPoint to store (should have unit length).

Variable Documentation

const char* const help_insolidangle[]
static
Initial value:
= {
"::",
"",
" insolidangle center POS span POS_EXPR [cutoff REAL]",
"",
"This keyword selects atoms that are within [TT]REAL[tt] degrees",
"(default=5) of any position in [TT]POS_EXPR[tt] as seen from [TT]POS[tt]",
"a position expression that evaluates to a single position), i.e., atoms",
"in the solid angle spanned by the positions in [TT]POS_EXPR[tt] and",
"centered at [TT]POS[tt].[PAR]",
"Technically, the solid angle is constructed as a union of small cones",
"whose tip is at [TT]POS[tt] and the axis goes through a point in",
"[TT]POS_EXPR[tt]. There is such a cone for each position in",
"[TT]POS_EXPR[tt], and point is in the solid angle if it lies within any",
"of these cones. The cutoff determines the width of the cones.",
}

Help text for the insolidangle selection method.

gmx_ana_selmethod_t sm_insolidangle
Initial value:
= {
"insolidangle",
4 ,
nullptr,
nullptr,
nullptr,
{ "insolidangle center POS span POS_EXPR [cutoff REAL]",
"Selecting atoms in a solid angle",
}
One group of atoms.
Definition: selvalue.h:58
static void * init_data_insolidangle(int npar, gmx_ana_selparam_t *param)
Allocates data for the insolidangle selection method.
Definition: sm_insolidangle.cpp:362
static void init_insolidangle(const gmx_mtop_t *top, int npar, gmx_ana_selparam_t *param, void *data)
Initializes the insolidangle selection method.
Definition: sm_insolidangle.cpp:384
static gmx_ana_selparam_t smparams_insolidangle[]
Parameters for the insolidangle selection method.
Definition: sm_insolidangle.cpp:316
static void evaluate_insolidangle(const gmx::SelMethodEvalContext &context, gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
Evaluates the insolidangle selection method.
Definition: sm_insolidangle.cpp:499
constexpr int asize(T(&)[N])
Calculates the number of elements in a static array at compile time.
Definition: arraysize.h:50
static void free_data_insolidangle(void *data)
Frees the data allocated for the insolidangle selection method.
Definition: sm_insolidangle.cpp:425
static const char *const help_insolidangle[]
Help text for the insolidangle selection method.
Definition: sm_insolidangle.cpp:323
static void init_frame_insolidangle(const gmx::SelMethodEvalContext &context, void *data)
Initializes the evaluation of the insolidangle selection method for a frame.
Definition: sm_insolidangle.cpp:443

Selection method data for the insolidangle method.

gmx_ana_selparam_t smparams_insolidangle[]
static
Initial value:
= {
{ "center", { POS_VALUE, 1, { nullptr } }, nullptr, 4 },
{ "span", { POS_VALUE, -1, { nullptr } }, nullptr, 4 | 16 },
{ "cutoff", { REAL_VALUE, 1, { nullptr } }, nullptr, 2 },
}
One or more real values.
Definition: selvalue.h:55
One or more position values.
Definition: selvalue.h:57

Parameters for the insolidangle selection method.