Rosetta
Classes | Typedefs | Enumerations | Functions | Variables
core::scoring::sasa Namespace Reference

Classes

class  LeGrandSasa
 LeGrand SASA approximation method Used by SasaCalc but can be used by itself. Virt atms are skipped as radii=0. More...
 
class  SasaCalc
 Main interface to sasa calculations outside of pose metrics. Virt atms are skipped as radii=0. More...
 
class  SasaMethod
 Abstract base class for SasaMethods. Feel free to edit as needed. Virt atms are skipped as radii=0. More...
 

Typedefs

typedef utility::pointer::shared_ptr< LeGrandSasaLeGrandSasaOP
 
typedef utility::pointer::shared_ptr< LeGrandSasa const > LeGrandSasaCOP
 
typedef utility::pointer::shared_ptr< SasaCalcSasaCalcOP
 
typedef utility::pointer::shared_ptr< SasaCalc const > SasaCalcCOP
 
typedef utility::pointer::shared_ptr< SasaMethodSasaMethodOP
 
typedef utility::pointer::shared_ptr< SasaMethod const > SasaMethodCOP
 

Enumerations

enum class  SasaMethodHPMode {
  ALL_SASA = 1 , POLAR_SASA , HYDROPHOBIC_SASA , INVALID_MODE ,
  END_OF_LIST = INVALID_MODE
}
 The selection mode. Allows selecting polar SASA only, hydrophobic SASA only, polar SASA only, etc. More...
 
enum  SasaRadii {
  LJ = 1 , legacy , naccess , reduce ,
  chothia =naccess , SasaRadii_total = reduce
}
 Type of Radii to use. More...
 
enum  SasaMethodEnum { LeGrand = 1 , SasaMethodType_total = LeGrand }
 

Functions

SasaMethodOP create_sasa_method (SasaMethodEnum method, core::Real probe_radius, SasaRadii radii_set)
 Very (very) basic implementation until I understand the regular implementation used by constraints/features/etc. Also used for me to debug everything else before creating the real factory. More...
 
utility::vector1< Realper_res_sc_sasa (const pose::Pose &pose)
 Absolute per residue sidechain SASA. More...
 
utility::vector1< Realrel_per_res_sc_sasa (const pose::Pose &pose)
 Relative per residue sidechain SASA. More...
 
bool is_res_exposed (const pose::Pose &pose, core::Size resnum)
 Is residue exposed? More...
 
std::pair< Real, Realget_sc_bb_sasa (const pose::Pose &pose, const id::AtomID_Map< Real > &atom_sasa)
 Calculate the sidechain and backbone sasa from atom sasa. More...
 
std::pair< utility::vector1< Real >, utility::vector1< Real > > get_sc_bb_sasa_per_res (const pose::Pose &pose, const id::AtomID_Map< Real > &atom_sasa)
 
SasaMethodEnum get_sasa_method_from_string (std::string method)
 Gets sasa enum from string passed by options system. More...
 
std::string get_sasa_radii_parameter_name (SasaRadii radii_set)
 Get string name of SASA radii used to obtain extra parameter index from atom_type_set. More...
 
SasaRadii get_sasa_radii_set_from_string (std::string radii_set)
 Gets sasa radii enum from string passed by options system. More...
 
ObjexxFCL::FArray2D_int angles (num_phi, num_theta)
 
ObjexxFCL::FArray2D_ubyte masks (num_bytes, num_overlaps *num_orientations)
 
ObjexxFCL::FArray2D_ubyte const & get_legrand_sasa_masks ()
 Returns const access to the masks FArray, which contains the information in the SASA database file sampling/SASA-masks.dat. Adding this in so that the values in the SASA database files can be used in SASA-based scores. (ronj) More...
 
ObjexxFCL::FArray2D_int const & get_legrand_sasa_angles ()
 Returns const access to the angles FArray, which contains the information in the SASA database file sampling/SASA-angles.dat. Adding this in so that the values in the SASA database files can be used in SASA-based scores. (ronj) More...
 
void input_legrand_sasa_dats ()
 Reads in the SASA database files sampling/SASA-angles.dat and sampling/SASA-masks.dat into FArrays above. More...
 
void get_legrand_atomic_overlap (Real const radius_a, Real const radius_b, Real const distance_ijxyz, int &degree_of_overlap)
 
void get_legrand_orientation (Vector const &a_xyz, Vector const &b_xyz, int &phi_index, int &theta_index, Real distance_ijxyz)
 Gets the orientation of a to b (i to j, see below). Does this by calculating two angles, aphi and theta. (j) More...
 
void get_legrand_2way_orientation (Vector const &a_xyz, Vector const &b_xyz, int &phi_a2b_index, int &theta_a2b_index, int &phi_b2a_index, int &theta_b2a_index, Real distance_ijxyz)
 Gets the orientation of a to b (i to j, see below). Does this by calculating two angles, aphi and theta. (j) More...
 

Variables

int const num_bytes = 21
 
int const num_phi = 64
 
int const num_theta = 64
 
int const num_overlaps = 100
 
int const num_orientations = 162
 

Typedef Documentation

◆ LeGrandSasaCOP

typedef utility::pointer::shared_ptr< LeGrandSasa const> core::scoring::sasa::LeGrandSasaCOP

◆ LeGrandSasaOP

typedef utility::pointer::shared_ptr< LeGrandSasa> core::scoring::sasa::LeGrandSasaOP

◆ SasaCalcCOP

typedef utility::pointer::shared_ptr< SasaCalc const> core::scoring::sasa::SasaCalcCOP

◆ SasaCalcOP

typedef utility::pointer::shared_ptr< SasaCalc> core::scoring::sasa::SasaCalcOP

◆ SasaMethodCOP

typedef utility::pointer::shared_ptr< SasaMethod const> core::scoring::sasa::SasaMethodCOP

◆ SasaMethodOP

typedef utility::pointer::shared_ptr< SasaMethod> core::scoring::sasa::SasaMethodOP

Enumeration Type Documentation

◆ SasaMethodEnum

Enumerator
LeGrand 
SasaMethodType_total 

◆ SasaMethodHPMode

The selection mode. Allows selecting polar SASA only, hydrophobic SASA only, polar SASA only, etc.

If you add to this list, update the SasaMethod::sasa_metric_name_from_mode() function!

Author
Vikram K. Mulligan (vmull.nosp@m.igan.nosp@m.@flat.nosp@m.iron.nosp@m.insti.nosp@m.tute.nosp@m..org).
Enumerator
ALL_SASA 
POLAR_SASA 
HYDROPHOBIC_SASA 
INVALID_MODE 
END_OF_LIST 

◆ SasaRadii

Type of Radii to use.

LJ: Refers to Leonard Jones radii - Rosetta uses radii at the minimum of the potential (sigma2). Legacy: Refers to radii optimized for a no longer in use term, but some protocols have been optimized to use it. naccess: Refers to radii used in the program naccess. Originally derived from Chothia. Do not use for all-atom SASA as hydrogens are implicitly included. 'The Nature of the Accessible and Buried Surfaces in Proteins' J. Mol. Biol. (1976) 105, 1-14 reduce: Radii used by the program reduce. Hydrogens are explicitly included in the radii.

Enumerator
LJ 
legacy 
naccess 
reduce 
chothia 
SasaRadii_total 

Function Documentation

◆ angles()

ObjexxFCL::FArray2D_int core::scoring::sasa::angles ( num_phi  ,
num_theta   
)

◆ create_sasa_method()

SasaMethodOP core::scoring::sasa::create_sasa_method ( SasaMethodEnum  ,
core::Real  probe_radius,
SasaRadii  radii_set 
)

Very (very) basic implementation until I understand the regular implementation used by constraints/features/etc. Also used for me to debug everything else before creating the real factory.

Referenced by core::scoring::sasa::SasaCalc::setup_sasa_method().

◆ get_legrand_2way_orientation()

void core::scoring::sasa::get_legrand_2way_orientation ( Vector const &  a_xyz,
Vector const &  b_xyz,
int &  phi_a2b_index,
int &  theta_a2b_index,
int &  phi_b2a_index,
int &  theta_b2a_index,
Real  distance_ijxyz 
)

Gets the orientation of a to b (i to j, see below). Does this by calculating two angles, aphi and theta. (j)

ronj This function is the same as the function above but get the orientation of a to b simultaneously with the ronj orientation of b to a. The same result could be achieved by making two separate get_legrand_2way_orientation() calls ronj but this method does it more efficiently by avoiding an atan2 and acos call. Instead, once you compute the ronj phi and theta for a on b, you can add/subtrate pi factors to get the phi and theta for b on a. ronj Still not sure how this method returns the correct values, though.

This function is the same as the function above but get the orientation of a to b simultaneously with the orientation of b to a. The same result could be achieved by making two separate get_2way_orientation() calls but this method does it more efficiently by avoiding an atan2 and acos call. Instead, once you compute the phi and theta for a on b, you can add/subtrate pi factors to get the phi and theta for b on a. Still not sure how this method returns the correct values, though. (ronj)

References num_phi, and num_theta.

Referenced by core::pack::interaction_graph::RotamerDots::get_atom_atom_coverage(), protocols::vardist_solaccess::VarSolDRotamerDots::get_atom_overlap_masks(), core::pack::interaction_graph::InvRotamerDots::overlap_exposed(), and core::pack::interaction_graph::InvRotamerDots::write_circle_intersection_mask_to_kinemage().

◆ get_legrand_atomic_overlap()

void core::scoring::sasa::get_legrand_atomic_overlap ( Real const  radius_a,
Real const  radius_b,
Real const  distance_ijxyz,
int &  degree_of_overlap 
)

getting overlap from a to b (or i to j, as the atoms are referred to in calc_per_atom_sasa below). this returns the degree of overlap between two atoms adapted from erics code in area.c GetD2 and returns value from 1 to 100. This calculation is based on the law of cosines. See LeGrand and Merz, Journal of Computational Chemistry 14(3):349-52 (1993). Note that equation (4) is wrong, the denominator should be 2*ri*riq instead of 2*ri*rq (j)

The function gets passed in the sasa radius of atom i (plus the probe radius), the sasa radius of atom j (plus the probe radius), the distance between the atom centers, and a reference to the degree of overlap (represented as an int). The degree of overlap that's returned can be thought of as how much of atom a is covered by atom b. A value of 100 means that atom a is completely covered up by atom b. A value of 1 means that not much of the surface of 'a' is covered up by 'b'. The law of cosines relates the cosine of one angle of a triangle to the lengths of its sides. More specifically, c^2 = a^2 + b^2 - 2*a*b*cos theta, where theta is the angle between sides a and b. For the function we want to compute the angle of the cone of intersection between spheres 'a' and 'b'. Let the radius of atom a be ri, and the radius of atom b be rq, and the distance between atom centers be riq. Let the angle between ri and riq be theta_iq. The cosine of theta_iq will be equivalent to ( ri^2 + riq^2 - rq^2 ) / 2 * ri * riq

Referenced by core::pack::interaction_graph::RotamerDots::get_atom_atom_coverage(), core::pack::interaction_graph::HPatchNode< V, E, G >::initialize_atom_atom_overlap_cache(), core::pack::interaction_graph::HPatchBackgroundNode< V, E, G >::initialize_atom_atom_overlaps(), core::pack::interaction_graph::HPatchInteractionGraph< V, E, G >::initialize_bg_bg_atom_atom_overlaps(), protocols::vardist_solaccess::VarSolDRotamerDots::overlap_atoms(), core::pack::interaction_graph::InvRotamerDots::overlap_exposed(), and core::pack::interaction_graph::InvRotamerDots::write_circle_intersection_mask_to_kinemage().

◆ get_legrand_orientation()

void core::scoring::sasa::get_legrand_orientation ( Vector const &  a_xyz,
Vector const &  b_xyz,
int &  phi_index,
int &  theta_index,
Real  distance_ijxyz 
)

Gets the orientation of a to b (i to j, see below). Does this by calculating two angles, aphi and theta. (j)

ronj This function is used to get two indexes (phi and theta) which are used to get the index of a dot on the ronj surface of the 'a' sphere. When calculating how much surface area sphere b covers on a, we can get the degree ronj of overlap from the function above, but it's not necessarily the case that the vector that connects the center ronj of atom 'a' and atom 'b' goes through one of the predetermined dot locations on the surface of 'a'. In fact, ronj it's very unlikely that the vector goes through a predetermined dot. Instead, what is done is the actual point ronj of intersection (the outermost point of a on the line from the center of 'a' to center of 'b') is converted ronj to spherical polar coordinates. Then, the values are used to find the location of the closest predetermined ronj point on the surface of 'a' using a lookup table. So what this function needs to do is convert the ronj cartesian coordinate of the actual point of intersection into polar coordinates. ronj ronj To get the spherical, polar coordinates of a cartesian point x,y,z, use these equations: ronj r = sqrt( x^2 + y^2 + z^2 ) ronj theta = arccos( z / r ) ronj phi = arctan( y / x )

ronj Then, once we have the true phi and theta, we need to translate this into an index (or offset) for the correct ronj value in the database file. There are 64 phi angle bin and 64 theta bins in the database file sampling/SASA-angles.dat. ronj We need to convert the phi and theta into indexes for this file by multiplying them by num_phi / 2*pi. ronj Note: I think phi and theta have been reversed in the function below. The code below uses the following: ronj phi = arccos( z ) ronj theta = arctan( y / x )

ronj After a couple of weeks trying to write tests for this function, I have been unsuccessful in figuring out why ronj it's doing what it does. Despite using the wrong equations, it seems to work. Comparing the total residue ronj SASA values calculated by calc_per_atom_sasa() below results in a correlation of 0.98 against what the program ronj NACCESS finds for the same residues. This test was done on a small 110aa protein. I also looked at the per-atom ronj total SASA and the correlation for all atoms (mini v. NACCESS) was approximately 0.94. I'm using exactly the same ronj van der Waals radii for both programs so I feel like the correlations should be 1.0. Explanations for the ronj differences can be 1) this method is doing something wrong in calculating the closest surface point, 2) this ronj method is correct but the masks that are in the database are not aligned to the surface points correctly, 3) the ronj differences are solely due to the different way that the two program calculate surface area.

This function is used to get two indexes (phi and theta) which are used to get the index of a dot on the surface of the 'a' sphere. When calculating how much surface area sphere b covers on a, we can get the degree of overlap from the function above, but it's not necessarily the case that the vector that connects the center of atom 'a' and atom 'b' goes through one of the predetermined dot locations on the surface of 'a'. In fact, it's very unlikely that the vector goes through a predetermined dot. Instead, what is done is the actual point of intersection (the outermost point of a on the line from the center of 'a' to center of 'b') is converted to spherical polar coordinates. Then, the values are used to find the location of the closest predetermined point on the surface of 'a' using a lookup table. So what this function needs to do is convert the cartesian coordinate of the actual point of intersection into polar coordinates.

To get the spherical, polar coordinates of a cartesian point x,y,z, use these equations: r = sqrt( x^2 + y^2 + z^2 ) theta = arccos( z / r ) phi = arctan( y / x )

Then, once we have the true phi and theta, we need to translate this into an index (or offset) for the correct value in the database file. There are 64 phi angle bin and 64 theta bins in the database file sampling/SASA-angles.dat. We need to convert the phi and theta into indexes for this file by multiplying them by num_phi / 2*pi. Note: I think phi and theta have been reversed in the function below. The code below uses the following: phi = arccos( z ) theta = arctan( y / x )

After a couple of weeks trying to write tests for this function, I have been unsuccessful in figuring out why it's doing what it does. Despite using the wrong equations, it seems to work. Comparing the total residue SASA values calculated by calc_per_atom_sasa() below results in a correlation of 0.98 against what the program NACCESS finds for the same residues. This test was done on a small 110aa protein. I also looked at the per-atom total SASA and the correlation for all atoms (mini v. NACCESS) was approximately 0.94. I'm using exactly the same van der Waals radii for both programs so I feel like the correlations should be 1.0. Explanations for the differences can be 1) this method is doing something wrong in calculating the closest surface point, 2) this method is correct but the masks that are in the database are not aligned to the surface points correctly, 3) the differences are solely due to the different way that the two program calculate surface area. (ronj)

References num_phi, and num_theta.

◆ get_legrand_sasa_angles()

ObjexxFCL::FArray2D_int const & core::scoring::sasa::get_legrand_sasa_angles ( )

Returns const access to the angles FArray, which contains the information in the SASA database file sampling/SASA-angles.dat. Adding this in so that the values in the SASA database files can be used in SASA-based scores. (ronj)

Return total SASA for side chains. @Note: Uses Naccess radii by default for all atom. Not for general use. Values returned used to generate SVM for ProQ/ProQ2.

Returns const access to the angles FArray, which contains the information in the SASA database file sampling/SASA-angles.dat. Adding this in so that the values in the SASA database files can be used in SASA-based scores. (ronj)

References angles(), and input_legrand_sasa_dats().

Referenced by core::pack::interaction_graph::RotamerDots::initialize_sasa_arrays(), and protocols::vardist_solaccess::VarSolDistSasaCalculator::initialize_sasa_arrays().

◆ get_legrand_sasa_masks()

ObjexxFCL::FArray2D_ubyte const & core::scoring::sasa::get_legrand_sasa_masks ( )

Returns const access to the masks FArray, which contains the information in the SASA database file sampling/SASA-masks.dat. Adding this in so that the values in the SASA database files can be used in SASA-based scores. (ronj)

References input_legrand_sasa_dats(), and masks().

Referenced by core::pack::interaction_graph::RotamerDots::initialize_sasa_arrays(), and protocols::vardist_solaccess::VarSolDistSasaCalculator::initialize_sasa_arrays().

◆ get_sasa_method_from_string()

SasaMethodEnum core::scoring::sasa::get_sasa_method_from_string ( std::string  method)

Gets sasa enum from string passed by options system.

Enum Management. Can go in separate file if it gets large.

References LeGrand.

Referenced by core::scoring::sasa::SasaCalc::set_defaults().

◆ get_sasa_radii_parameter_name()

std::string core::scoring::sasa::get_sasa_radii_parameter_name ( SasaRadii  radii_set)

Get string name of SASA radii used to obtain extra parameter index from atom_type_set.

References chothia, legacy, LJ, naccess, reduce, and SasaRadii_total.

Referenced by core::scoring::sasa::LeGrandSasa::calculate().

◆ get_sasa_radii_set_from_string()

SasaRadii core::scoring::sasa::get_sasa_radii_set_from_string ( std::string  radii_set)

Gets sasa radii enum from string passed by options system.

References chothia, legacy, LJ, naccess, and reduce.

Referenced by core::scoring::sasa::SasaCalc::set_defaults().

◆ get_sc_bb_sasa()

std::pair< Real, Real > core::scoring::sasa::get_sc_bb_sasa ( const pose::Pose pose,
const id::AtomID_Map< Real > &  atom_sasa 
)

◆ get_sc_bb_sasa_per_res()

std::pair< utility::vector1< Real >, utility::vector1< Real > > core::scoring::sasa::get_sc_bb_sasa_per_res ( const pose::Pose pose,
const id::AtomID_Map< Real > &  atom_sasa 
)

◆ input_legrand_sasa_dats()

void core::scoring::sasa::input_legrand_sasa_dats ( )

Reads in the SASA database files sampling/SASA-angles.dat and sampling/SASA-masks.dat into FArrays above.

References angles(), core::init::init(), masks(), num_bytes, num_orientations, num_overlaps, num_phi, and num_theta.

Referenced by get_legrand_sasa_angles(), and get_legrand_sasa_masks().

◆ is_res_exposed()

bool core::scoring::sasa::is_res_exposed ( const pose::Pose pose,
core::Size  resnum 
)

Is residue exposed?

Added by JKLeman (julia.nosp@m..koe.nosp@m.hler1.nosp@m.982@.nosp@m.gmail.nosp@m..com) Uses the function rel_per_res_sc_sasa above THIS IS EXPENSIVE, BE AWARE!!! IF YOU NEED TO RUN IT OVER THE ENTIRE PROTEIN, USE THE rel_per_res_sc_sasa FUNCTION INSTEAD!

References rel_per_res_sc_sasa().

◆ masks()

ObjexxFCL::FArray2D_ubyte core::scoring::sasa::masks ( num_bytes  ,
num_overlaps num_orientations 
)

◆ per_res_sc_sasa()

utility::vector1< Real > core::scoring::sasa::per_res_sc_sasa ( const pose::Pose pose)

Absolute per residue sidechain SASA.

Convenience Functions

Added by JKLeman (julia.nosp@m..koe.nosp@m.hler1.nosp@m.982@.nosp@m.gmail.nosp@m..com) GXG tripeptide values for sidechain SASA are taken from http://www.proteinsandproteomics.org/content/free/tables_1/table08.pdf

References core::scoring::sasa::SasaCalc::calculate(), and core::scoring::sasa::SasaCalc::get_residue_sasa_sc().

Referenced by rel_per_res_sc_sasa().

◆ rel_per_res_sc_sasa()

utility::vector1< Real > core::scoring::sasa::rel_per_res_sc_sasa ( const pose::Pose pose)

Variable Documentation

◆ num_bytes

int const core::scoring::sasa::num_bytes = 21

Referenced by input_legrand_sasa_dats().

◆ num_orientations

int const core::scoring::sasa::num_orientations = 162

Referenced by input_legrand_sasa_dats().

◆ num_overlaps

int const core::scoring::sasa::num_overlaps = 100

Referenced by input_legrand_sasa_dats().

◆ num_phi

int const core::scoring::sasa::num_phi = 64

◆ num_theta

int const core::scoring::sasa::num_theta = 64