Gromacs
2018.8
|
#include <gromacs/trajectoryanalysis/modules/surfacearea.h>
Computes surface areas for a group of atoms/spheres.
This class provides a surface area/volume calculator.
The algorithm is based on representing each atom/sphere surface as a set of dots, and determining which dots are on the surface (not covered by any other atom/sphere). The dots are distributed evenly using an icosahedron- or a dodecahedron-based method (see the original reference cited in the code). The area is then estimated from the area represented by each dot. The volume is calculated by selecting a fixed point and integrating over the surface dots, summing up the cones whose tip is at the fixed point and base at the surface points.
The default dot density per sphere is 32, which gives quite inaccurate areas and volumes, but a reasonable number of surface points. According to original documentation of the method, a density of 600-700 dots gives an accuracy of 1.5 A^2 per atom.
Public Member Functions | |
SurfaceAreaCalculator () | |
Initializes a surface area calculator. More... | |
void | setDotCount (int dotCount) |
Sets the number of surface dots per sphere to use. More... | |
void | setRadii (const ArrayRef< const real > &radius) |
Sets the radii of spheres to use in the calculation. More... | |
void | setCalculateVolume (bool bVolume) |
Requests calculation of volume. More... | |
void | setCalculateAtomArea (bool bAtomArea) |
Requests output of per-atom areas. More... | |
void | setCalculateSurfaceDots (bool bDots) |
Requests output of all surface dots. More... | |
void | calculate (const rvec *x, const t_pbc *pbc, int nat, int index[], int flags, real *area, real *volume, real **at_area, real **lidots, int *n_dots) const |
Calculates the surface area for a set of positions. More... | |
gmx::SurfaceAreaCalculator::SurfaceAreaCalculator | ( | ) |
Initializes a surface area calculator.
std::bad_alloc | if out of memory. |
void gmx::SurfaceAreaCalculator::calculate | ( | const rvec * | x, |
const t_pbc * | pbc, | ||
int | nat, | ||
int | index[], | ||
int | flags, | ||
real * | area, | ||
real * | volume, | ||
real ** | at_area, | ||
real ** | lidots, | ||
int * | n_dots | ||
) | const |
Calculates the surface area for a set of positions.
[in] | x | Atom positions (sphere centers). |
[in] | pbc | PBC information (if NULL , calculation is done without PBC). |
[in] | nat | Number of atoms to calculate. |
[in] | index | Atom indices to include in the calculation. |
[in] | flags | Additional flags for the calculation. |
[out] | area | Total surface area (must be non-NULL ). |
[out] | volume | Total volume (can be NULL ). |
[out] | at_area | Surface area for each atom in index (nat values) (can be NULL ). |
[out] | lidots | Surface dots as x,y,z triplets (3*lidots values) (can be NULL ). |
[out] | n_dots | Number of surface dots in lidots (can be NULL ). |
Calculates the surface area of spheres centered at x[index[0]]
, ..., x[index[nat-1]]
, with radii radii[index[0]]
, ..., where radii
is the array passed to setRadii().
If flags
is 0, the calculation is done for the items specified with setCalculateVolume(), setCalculateAtomArea(), and setCalculateSurfaceDots(). Flags can specify FLAG_VOLUME, FLAG_ATOM_AREA, and/or FLAG_DOTS to request additional output for this particular calculation. If any output is NULL
, that output is not calculated, irrespective of the calculation mode set.
void gmx::SurfaceAreaCalculator::setCalculateAtomArea | ( | bool | bAtomArea | ) |
Requests output of per-atom areas.
If not called, and FLAG_ATOM_AREA is not passed to calculate(), the atom area output is not produced.
Does not throw.
void gmx::SurfaceAreaCalculator::setCalculateSurfaceDots | ( | bool | bDots | ) |
Requests output of all surface dots.
If not called, and FLAG_DOTS is not passed to calculate(), the surface dot output is not produced.
Does not throw.
void gmx::SurfaceAreaCalculator::setCalculateVolume | ( | bool | bVolume | ) |
Requests calculation of volume.
If not called, and FLAG_VOLUME is not passed to calculate(), the volume output is not produced.
Does not throw.
void gmx::SurfaceAreaCalculator::setDotCount | ( | int | dotCount | ) |
Sets the number of surface dots per sphere to use.
This function must be called before calculate() to set the desired accuracy/computational cost.
Sets the radii of spheres to use in the calculation.
[in] | radius | Radius for each atom/sphere. |
This function must be called before calculate() to set the radii for the spheres. All calculations must use the same set of radii to share the same grid search. These radii are used as-is, without adding any probe radius. The passed array must remain valid for the lifetime of this object.
Does not throw.