sasa

Bindings for core::scoring::sasa namespace

class pyrosetta.rosetta.core.scoring.sasa.LeGrandSasa

Bases: SasaMethod

LeGrand SASA approximation method

Used by SasaCalc but can be used by itself. Virt atms are skipped as radii=0

LeGrand S, Merz KM. Rapid approximation to molecular surface area via the use of Boolean logic and look-up tables.

J Comput Chem 1993;14:349-352.

Fortran Implementation: Jerry Tsai C++ Translation: Jeff Gray Cleanup/Bugfixes/OOP: Jared Adolf-Bryfogle

assign(self: pyrosetta.rosetta.core.scoring.sasa.LeGrandSasa, : pyrosetta.rosetta.core.scoring.sasa.LeGrandSasa) pyrosetta.rosetta.core.scoring.sasa.LeGrandSasa

C++: core::scoring::sasa::LeGrandSasa::operator=(const class core::scoring::sasa::LeGrandSasa &) –> class core::scoring::sasa::LeGrandSasa &

calculate(self: pyrosetta.rosetta.core.scoring.sasa.LeGrandSasa, pose: pyrosetta.rosetta.core.pose.Pose, atom_subset: pyrosetta.rosetta.core.id.AtomID_Map_bool_t, atom_sasa: pyrosetta.rosetta.core.id.AtomID_Map_double_t, rsd_sasa: pyrosetta.rosetta.utility.vector1_double) float

Calculate Sasa. Atoms not calculated have -1 sasa. This is carried over for compatability purposes.

C++: core::scoring::sasa::LeGrandSasa::calculate(const class core::pose::Pose &, const class core::id::AtomID_Map<bool> &, class core::id::AtomID_Map<double> &, class utility::vector1<double, class std::allocator<double> > &) –> double

get_2way_orientation(self: pyrosetta.rosetta.core.scoring.sasa.LeGrandSasa, a_xyz: pyrosetta.rosetta.numeric.xyzVector_double_t, b_xyz: pyrosetta.rosetta.numeric.xyzVector_double_t, phi_a2b_index: int, theta_a2b_index: int, phi_b2a_index: int, theta_b2a_index: int, distance_ijxyz: float) None

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

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)

C++: core::scoring::sasa::LeGrandSasa::get_2way_orientation(const class numeric::xyzVector<double> &, const class numeric::xyzVector<double> &, int &, int &, int &, int &, double) const –> void

get_angles(self: pyrosetta.rosetta.core.scoring.sasa.LeGrandSasa) pyrosetta.rosetta.ObjexxFCL.FArray2D_int_t
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)

C++: core::scoring::sasa::LeGrandSasa::get_angles() const –> const class ObjexxFCL::FArray2D<int> &

get_masks(self: pyrosetta.rosetta.core.scoring.sasa.LeGrandSasa) pyrosetta.rosetta.ObjexxFCL.FArray2D_ObjexxFCL_ubyte_t
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)

C++: core::scoring::sasa::LeGrandSasa::get_masks() const –> const class ObjexxFCL::FArray2D<class ObjexxFCL::ubyte> &

get_name(self: pyrosetta.rosetta.core.scoring.sasa.LeGrandSasa) str

C++: core::scoring::sasa::LeGrandSasa::get_name() const –> std::string

get_orientation(self: pyrosetta.rosetta.core.scoring.sasa.LeGrandSasa, a_xyz: pyrosetta.rosetta.numeric.xyzVector_double_t, b_xyz: pyrosetta.rosetta.numeric.xyzVector_double_t, phi_index: int, theta_index: int, distance_ijxyz: float) None
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)

C++: core::scoring::sasa::LeGrandSasa::get_orientation(const class numeric::xyzVector<double> &, const class numeric::xyzVector<double> &, int &, int &, double) const –> void

get_overlap(self: pyrosetta.rosetta.core.scoring.sasa.LeGrandSasa, radius_a: float, radius_b: float, distance_ijxyz: float, degree_of_overlap: int) None
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

C++: core::scoring::sasa::LeGrandSasa::get_overlap(const double, const double, const double, int &) const –> void

static is_polar_atom(rsd: pyrosetta.rosetta.core.conformation.Residue, atom_index: int) bool

Classify an atom on a residue as “polar” for the purposes of SASA

C++: core::scoring::sasa::SasaMethod::is_polar_atom(const class core::conformation::Residue &, const unsigned long) –> bool

static list_sasa_method_hp_modes() str

Construct a comma-separeted string listing all of the sasa metric modes.

Vikram K. Mulligan (vmulligan.org).

C++: core::scoring::sasa::SasaMethod::list_sasa_method_hp_modes() –> std::string

sasa_method_hp_mode(self: pyrosetta.rosetta.core.scoring.sasa.SasaMethod) pyrosetta.rosetta.core.scoring.sasa.SasaMethodHPMode

Get whether we’re counting all SASA (default), polar SASA, or hydrophobic SASA.

Vikram K. Mulligan (vmulligan.org).

C++: core::scoring::sasa::SasaMethod::sasa_method_hp_mode() const –> enum core::scoring::sasa::SasaMethodHPMode

static sasa_metric_mode_from_name(mode_name: str) pyrosetta.rosetta.core.scoring.sasa.SasaMethodHPMode

Given the name of the SasaMethodHPMode, get the mode.

SasaMethodHPMode::INVALID_MODE if the string can’t be parsed; the correct SasaMethodHPMode if it can.

Vikram K. Mulligan (vmulligan.org).

C++: core::scoring::sasa::SasaMethod::sasa_metric_mode_from_name(const std::string &) –> enum core::scoring::sasa::SasaMethodHPMode

static sasa_metric_name_from_mode(mode: pyrosetta.rosetta.core.scoring.sasa.SasaMethodHPMode) str

Given the SasaMethodHPMode, get the name.

Vikram K. Mulligan (vmulligan.org).

C++: core::scoring::sasa::SasaMethod::sasa_metric_name_from_mode(const enum core::scoring::sasa::SasaMethodHPMode) –> std::string

set_include_probe_radius_in_calc(self: pyrosetta.rosetta.core.scoring.sasa.SasaMethod, include_probe_radius: bool) None

Include the probe radius in calc. Typical for SASA.

C++: core::scoring::sasa::SasaMethod::set_include_probe_radius_in_calc(bool) –> void

set_probe_radius(self: pyrosetta.rosetta.core.scoring.sasa.SasaMethod, probe_radius: float) None

Set the probe radius. Typical value is that of water at 1.4 A

C++: core::scoring::sasa::SasaMethod::set_probe_radius(double) –> void

set_radii_set(self: pyrosetta.rosetta.core.scoring.sasa.SasaMethod, radii_set: pyrosetta.rosetta.core.scoring.sasa.SasaRadii) None

Set the radii type.

C++: core::scoring::sasa::SasaMethod::set_radii_set(enum core::scoring::sasa::SasaRadii) –> void

set_sasa_method_hp_mode(self: pyrosetta.rosetta.core.scoring.sasa.SasaMethod, mode_in: pyrosetta.rosetta.core.scoring.sasa.SasaMethodHPMode) None

Set whether we’re counting all SASA (default), polar SASA, or hydrophobic SASA.

Vikram K. Mulligan (vmulligan.org).

C++: core::scoring::sasa::SasaMethod::set_sasa_method_hp_mode(const enum core::scoring::sasa::SasaMethodHPMode) –> void

set_use_big_polar_hydrogen(self: pyrosetta.rosetta.core.scoring.sasa.SasaMethod, big_polar_h: bool) None

Legacy option to increase polar hydrogen radii to 1.08A. Supported for now.

C++: core::scoring::sasa::SasaMethod::set_use_big_polar_hydrogen(bool) –> void

static skip_atom(rsd: pyrosetta.rosetta.core.conformation.Residue, atom_index: int, hp_mode: pyrosetta.rosetta.core.scoring.sasa.SasaMethodHPMode) bool
Given a residue, an atom index, and a SasaMethodHPMode, determine whether the atom is one to skip (returns true)

or count (returns false).

Vikram K. Mulligan (vmulligan.org).

C++: core::scoring::sasa::SasaMethod::skip_atom(const class core::conformation::Residue &, const unsigned long, const enum core::scoring::sasa::SasaMethodHPMode) –> bool

class pyrosetta.rosetta.core.scoring.sasa.SasaCalc

Bases: pybind11_object

Main interface to sasa calculations outside of pose metrics. Virt atms are skipped as radii=0

assign(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, : pyrosetta.rosetta.core.scoring.sasa.SasaCalc) pyrosetta.rosetta.core.scoring.sasa.SasaCalc

C++: core::scoring::sasa::SasaCalc::operator=(const class core::scoring::sasa::SasaCalc &) –> class core::scoring::sasa::SasaCalc &

calculate(*args, **kwargs)

Overloaded function.

  1. calculate(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, pose: pyrosetta.rosetta.core.pose.Pose) -> float

Calculate Sasa. Atoms not calculated have -1 sasa in AtomID_Map. This is carried over for compatability purposes.

C++: core::scoring::sasa::SasaCalc::calculate(const class core::pose::Pose &) –> double

  1. calculate(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, pose: pyrosetta.rosetta.core.pose.Pose, atom_sasa: pyrosetta.rosetta.core.id.AtomID_Map_double_t) -> float

// Legacy-style interfaces //////////

C++: core::scoring::sasa::SasaCalc::calculate(const class core::pose::Pose &, class core::id::AtomID_Map<double> &) –> double

  1. calculate(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, pose: pyrosetta.rosetta.core.pose.Pose, rsd_sasa: pyrosetta.rosetta.utility.vector1_double, rsd_hsasa: pyrosetta.rosetta.utility.vector1_double) -> float

C++: core::scoring::sasa::SasaCalc::calculate(const class core::pose::Pose &, class utility::vector1<double, class std::allocator<double> > &, class utility::vector1<double, class std::allocator<double> > &) –> double

  1. calculate(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, pose: pyrosetta.rosetta.core.pose.Pose, atom_sasa: pyrosetta.rosetta.core.id.AtomID_Map_double_t, rsd_sasa: pyrosetta.rosetta.utility.vector1_double, rsd_hsasa: pyrosetta.rosetta.utility.vector1_double) -> float

C++: core::scoring::sasa::SasaCalc::calculate(const class core::pose::Pose &, class core::id::AtomID_Map<double> &, class utility::vector1<double, class std::allocator<double> > &, class utility::vector1<double, class std::allocator<double> > &) –> double

  1. calculate(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, pose: pyrosetta.rosetta.core.pose.Pose, atom_sasa: pyrosetta.rosetta.core.id.AtomID_Map_double_t, rsd_sasa: pyrosetta.rosetta.utility.vector1_double, rsd_hsasa: pyrosetta.rosetta.utility.vector1_double, rsd_rel_hsasa: pyrosetta.rosetta.utility.vector1_double) -> float

C++: core::scoring::sasa::SasaCalc::calculate(const class core::pose::Pose &, class core::id::AtomID_Map<double> &, class utility::vector1<double, class std::allocator<double> > &, class utility::vector1<double, class std::allocator<double> > &, class utility::vector1<double, class std::allocator<double> > &) –> double

fill_data(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, total_hsasa: float, total_rel_hsasa: float, atom_sasa: pyrosetta.rosetta.core.id.AtomID_Map_double_t, rsd_sasa: pyrosetta.rosetta.utility.vector1_double, rsd_hsasa: pyrosetta.rosetta.utility.vector1_double, rel_hsasa: pyrosetta.rosetta.utility.vector1_double) None

Convenience function to fill most commonly used data.

C++: core::scoring::sasa::SasaCalc::fill_data(double &, double &, class core::id::AtomID_Map<double> &, class utility::vector1<double, class std::allocator<double> > &, class utility::vector1<double, class std::allocator<double> > &, class utility::vector1<double, class std::allocator<double> > &) –> void

get_atom_sasa(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) pyrosetta.rosetta.core.id.AtomID_Map_double_t

////Per Atom

C++: core::scoring::sasa::SasaCalc::get_atom_sasa() const –> class core::id::AtomID_Map<double>

get_rel_hphobic_sasa_by_charge(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) pyrosetta.rosetta.utility.vector1_double

C++: core::scoring::sasa::SasaCalc::get_rel_hphobic_sasa_by_charge() const –> class utility::vector1<double, class std::allocator<double> >

get_residue_hsasa(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) pyrosetta.rosetta.utility.vector1_double

C++: core::scoring::sasa::SasaCalc::get_residue_hsasa() const –> class utility::vector1<double, class std::allocator<double> >

get_residue_hsasa_bb(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) pyrosetta.rosetta.utility.vector1_double

C++: core::scoring::sasa::SasaCalc::get_residue_hsasa_bb() const –> class utility::vector1<double, class std::allocator<double> >

get_residue_hsasa_sc(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) pyrosetta.rosetta.utility.vector1_double

C++: core::scoring::sasa::SasaCalc::get_residue_hsasa_sc() const –> class utility::vector1<double, class std::allocator<double> >

get_residue_sasa(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) pyrosetta.rosetta.utility.vector1_double

/////Per Residue

C++: core::scoring::sasa::SasaCalc::get_residue_sasa() const –> class utility::vector1<double, class std::allocator<double> >

get_residue_sasa_bb(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) pyrosetta.rosetta.utility.vector1_double

C++: core::scoring::sasa::SasaCalc::get_residue_sasa_bb() const –> class utility::vector1<double, class std::allocator<double> >

get_residue_sasa_sc(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) pyrosetta.rosetta.utility.vector1_double

C++: core::scoring::sasa::SasaCalc::get_residue_sasa_sc() const –> class utility::vector1<double, class std::allocator<double> >

get_total_hsasa(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) float

C++: core::scoring::sasa::SasaCalc::get_total_hsasa() const –> double

get_total_hsasa_bb(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) float

C++: core::scoring::sasa::SasaCalc::get_total_hsasa_bb() const –> double

get_total_hsasa_sc(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) float

C++: core::scoring::sasa::SasaCalc::get_total_hsasa_sc() const –> double

get_total_rel_hsasa(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) float

C++: core::scoring::sasa::SasaCalc::get_total_rel_hsasa() const –> double

get_total_sasa(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) float

//////Totals

C++: core::scoring::sasa::SasaCalc::get_total_sasa() const –> double

get_total_sasa_bb(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) float

C++: core::scoring::sasa::SasaCalc::get_total_sasa_bb() const –> double

get_total_sasa_sc(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) float

C++: core::scoring::sasa::SasaCalc::get_total_sasa_sc() const –> double

sasa_method_hp_mode(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) pyrosetta.rosetta.core.scoring.sasa.SasaMethodHPMode

Get whether we’re counting all SASA (default), polar SASA, or hydrophobic SASA.

Vikram K. Mulligan (vmulligan.org).

C++: core::scoring::sasa::SasaCalc::sasa_method_hp_mode() const –> enum core::scoring::sasa::SasaMethodHPMode

set_calculation_method(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, method: pyrosetta.rosetta.core.scoring.sasa.SasaMethodEnum) None

C++: core::scoring::sasa::SasaCalc::set_calculation_method(enum core::scoring::sasa::SasaMethodEnum) –> void

set_defaults(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc) None

//////// Common Options ///////////

C++: core::scoring::sasa::SasaCalc::set_defaults() –> void

set_exclude_polar_atoms_by_charge(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, exclude_polar_all: bool) None
Polar carbons and other atoms should not be included in hydrophobic hSASA - though historically they were.

.4 is a relative number. This makes sure that carbonyl and carboxyl carbons are marked as polar, while others (protein-based) are non-polar

C++: core::scoring::sasa::SasaCalc::set_exclude_polar_atoms_by_charge(bool) –> void

set_explicit_hydrogen_included_radii_set(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, radii_set: pyrosetta.rosetta.core.scoring.sasa.SasaRadii) None

Radii set to use when including hydrogens (LJ/reduce)

C++: core::scoring::sasa::SasaCalc::set_explicit_hydrogen_included_radii_set(enum core::scoring::sasa::SasaRadii) –> void

set_implicit_hydrogen_included_radii_set(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, radii_set: pyrosetta.rosetta.core.scoring.sasa.SasaRadii) None
Radii set to use when not including hydrogens (naccess/chothia, legacy)

Do not use legacy unless you know what you are doing and why.

C++: core::scoring::sasa::SasaCalc::set_implicit_hydrogen_included_radii_set(enum core::scoring::sasa::SasaRadii) –> void

set_include_carbon_sulfer_only_in_hydrophobic_calc(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, include_c_s_only: bool) None

Typically, only carbon or sulfers are included in the calculation. If you are using ligands - this may not be good enough.

C++: core::scoring::sasa::SasaCalc::set_include_carbon_sulfer_only_in_hydrophobic_calc(bool) –> void

set_include_hydrogens_explicitly(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, include_hydrogens: bool) None

Include hydrogens explicitly

C++: core::scoring::sasa::SasaCalc::set_include_hydrogens_explicitly(bool) –> void

set_include_probe_radius_in_atom_radii(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, include_probe_radius: bool) None

This is typically done. Disabling it is more akin to obtaining the Surface Area than the SASA

C++: core::scoring::sasa::SasaCalc::set_include_probe_radius_in_atom_radii(bool) –> void

set_polar_charge_cutoff(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, cutoff: float) None

C++: core::scoring::sasa::SasaCalc::set_polar_charge_cutoff(double) –> void

set_probe_radius(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, probe_radius: float) None

Probe radius of 1.4 (water) is typically used

C++: core::scoring::sasa::SasaCalc::set_probe_radius(double) –> void

set_sasa_method_hp_mode(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, mode_in: pyrosetta.rosetta.core.scoring.sasa.SasaMethodHPMode) None

Set whether we’re counting all SASA (default), polar SASA, or hydrophobic SASA.

Vikram K. Mulligan (vmulligan.org).

C++: core::scoring::sasa::SasaCalc::set_sasa_method_hp_mode(const enum core::scoring::sasa::SasaMethodHPMode) –> void

set_use_big_polar_hydrogen(self: pyrosetta.rosetta.core.scoring.sasa.SasaCalc, big_polar_h: bool) None

Not for general use. Used to calculate unsatisfied buried polars with legacy radii (which implicitly had included hydrogens.)

C++: core::scoring::sasa::SasaCalc::set_use_big_polar_hydrogen(bool) –> void

class pyrosetta.rosetta.core.scoring.sasa.SasaMethod

Bases: pybind11_object

Abstract base class for SasaMethods. Feel free to edit as needed. Virt atms are skipped as radii=0

assign(self: pyrosetta.rosetta.core.scoring.sasa.SasaMethod, : pyrosetta.rosetta.core.scoring.sasa.SasaMethod) pyrosetta.rosetta.core.scoring.sasa.SasaMethod

C++: core::scoring::sasa::SasaMethod::operator=(const class core::scoring::sasa::SasaMethod &) –> class core::scoring::sasa::SasaMethod &

calculate(self: pyrosetta.rosetta.core.scoring.sasa.SasaMethod, pose: pyrosetta.rosetta.core.pose.Pose, atom_subset: pyrosetta.rosetta.core.id.AtomID_Map_bool_t, atom_sasa: pyrosetta.rosetta.core.id.AtomID_Map_double_t, rsd_sasa: pyrosetta.rosetta.utility.vector1_double) float

Calculate Sasa. Atoms not calculated have -1 sasa in AtomID_Map. This is carried over for compatability purposes.

C++: core::scoring::sasa::SasaMethod::calculate(const class core::pose::Pose &, const class core::id::AtomID_Map<bool> &, class core::id::AtomID_Map<double> &, class utility::vector1<double, class std::allocator<double> > &) –> double

get_name(self: pyrosetta.rosetta.core.scoring.sasa.SasaMethod) str

C++: core::scoring::sasa::SasaMethod::get_name() const –> std::string

static is_polar_atom(rsd: pyrosetta.rosetta.core.conformation.Residue, atom_index: int) bool

Classify an atom on a residue as “polar” for the purposes of SASA

C++: core::scoring::sasa::SasaMethod::is_polar_atom(const class core::conformation::Residue &, const unsigned long) –> bool

static list_sasa_method_hp_modes() str

Construct a comma-separeted string listing all of the sasa metric modes.

Vikram K. Mulligan (vmulligan.org).

C++: core::scoring::sasa::SasaMethod::list_sasa_method_hp_modes() –> std::string

sasa_method_hp_mode(self: pyrosetta.rosetta.core.scoring.sasa.SasaMethod) pyrosetta.rosetta.core.scoring.sasa.SasaMethodHPMode

Get whether we’re counting all SASA (default), polar SASA, or hydrophobic SASA.

Vikram K. Mulligan (vmulligan.org).

C++: core::scoring::sasa::SasaMethod::sasa_method_hp_mode() const –> enum core::scoring::sasa::SasaMethodHPMode

static sasa_metric_mode_from_name(mode_name: str) pyrosetta.rosetta.core.scoring.sasa.SasaMethodHPMode

Given the name of the SasaMethodHPMode, get the mode.

SasaMethodHPMode::INVALID_MODE if the string can’t be parsed; the correct SasaMethodHPMode if it can.

Vikram K. Mulligan (vmulligan.org).

C++: core::scoring::sasa::SasaMethod::sasa_metric_mode_from_name(const std::string &) –> enum core::scoring::sasa::SasaMethodHPMode

static sasa_metric_name_from_mode(mode: pyrosetta.rosetta.core.scoring.sasa.SasaMethodHPMode) str

Given the SasaMethodHPMode, get the name.

Vikram K. Mulligan (vmulligan.org).

C++: core::scoring::sasa::SasaMethod::sasa_metric_name_from_mode(const enum core::scoring::sasa::SasaMethodHPMode) –> std::string

set_include_probe_radius_in_calc(self: pyrosetta.rosetta.core.scoring.sasa.SasaMethod, include_probe_radius: bool) None

Include the probe radius in calc. Typical for SASA.

C++: core::scoring::sasa::SasaMethod::set_include_probe_radius_in_calc(bool) –> void

set_probe_radius(self: pyrosetta.rosetta.core.scoring.sasa.SasaMethod, probe_radius: float) None

Set the probe radius. Typical value is that of water at 1.4 A

C++: core::scoring::sasa::SasaMethod::set_probe_radius(double) –> void

set_radii_set(self: pyrosetta.rosetta.core.scoring.sasa.SasaMethod, radii_set: pyrosetta.rosetta.core.scoring.sasa.SasaRadii) None

Set the radii type.

C++: core::scoring::sasa::SasaMethod::set_radii_set(enum core::scoring::sasa::SasaRadii) –> void

set_sasa_method_hp_mode(self: pyrosetta.rosetta.core.scoring.sasa.SasaMethod, mode_in: pyrosetta.rosetta.core.scoring.sasa.SasaMethodHPMode) None

Set whether we’re counting all SASA (default), polar SASA, or hydrophobic SASA.

Vikram K. Mulligan (vmulligan.org).

C++: core::scoring::sasa::SasaMethod::set_sasa_method_hp_mode(const enum core::scoring::sasa::SasaMethodHPMode) –> void

set_use_big_polar_hydrogen(self: pyrosetta.rosetta.core.scoring.sasa.SasaMethod, big_polar_h: bool) None

Legacy option to increase polar hydrogen radii to 1.08A. Supported for now.

C++: core::scoring::sasa::SasaMethod::set_use_big_polar_hydrogen(bool) –> void

static skip_atom(rsd: pyrosetta.rosetta.core.conformation.Residue, atom_index: int, hp_mode: pyrosetta.rosetta.core.scoring.sasa.SasaMethodHPMode) bool
Given a residue, an atom index, and a SasaMethodHPMode, determine whether the atom is one to skip (returns true)

or count (returns false).

Vikram K. Mulligan (vmulligan.org).

C++: core::scoring::sasa::SasaMethod::skip_atom(const class core::conformation::Residue &, const unsigned long, const enum core::scoring::sasa::SasaMethodHPMode) –> bool

class pyrosetta.rosetta.core.scoring.sasa.SasaMethodEnum

Bases: pybind11_object

Members:

LeGrand

SasaMethodType_total

LeGrand = <SasaMethodEnum.LeGrand: 1>
SasaMethodType_total = <SasaMethodEnum.LeGrand: 1>
property name
property value
class pyrosetta.rosetta.core.scoring.sasa.SasaMethodHPMode

Bases: pybind11_object

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!

Vikram K. Mulligan (vmulligan.org).

Members:

ALL_SASA

POLAR_SASA

HYDROPHOBIC_SASA

INVALID_MODE

END_OF_LIST

ALL_SASA = <SasaMethodHPMode.ALL_SASA: 1>
END_OF_LIST = <SasaMethodHPMode.INVALID_MODE: 4>
HYDROPHOBIC_SASA = <SasaMethodHPMode.HYDROPHOBIC_SASA: 3>
INVALID_MODE = <SasaMethodHPMode.INVALID_MODE: 4>
POLAR_SASA = <SasaMethodHPMode.POLAR_SASA: 2>
property name
property value
class pyrosetta.rosetta.core.scoring.sasa.SasaRadii

Bases: pybind11_object

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.

Members:

LJ

legacy

naccess

reduce

chothia

SasaRadii_total

LJ = <SasaRadii.LJ: 1>
SasaRadii_total = <SasaRadii.reduce: 4>
chothia = <SasaRadii.naccess: 3>
legacy = <SasaRadii.legacy: 2>
naccess = <SasaRadii.naccess: 3>
property name
reduce = <SasaRadii.reduce: 4>
property value
pyrosetta.rosetta.core.scoring.sasa.create_sasa_method(method: pyrosetta.rosetta.core.scoring.sasa.SasaMethodEnum, probe_radius: float, radii_set: pyrosetta.rosetta.core.scoring.sasa.SasaRadii) pyrosetta.rosetta.core.scoring.sasa.SasaMethod
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.

C++: core::scoring::sasa::create_sasa_method(enum core::scoring::sasa::SasaMethodEnum, double, enum core::scoring::sasa::SasaRadii) –> class std::shared_ptr<class core::scoring::sasa::SasaMethod>

pyrosetta.rosetta.core.scoring.sasa.get_legrand_2way_orientation(a_xyz: pyrosetta.rosetta.numeric.xyzVector_double_t, b_xyz: pyrosetta.rosetta.numeric.xyzVector_double_t, phi_a2b_index: int, theta_a2b_index: int, phi_b2a_index: int, theta_b2a_index: int, distance_ijxyz: float) None

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

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)

C++: core::scoring::sasa::get_legrand_2way_orientation(const class numeric::xyzVector<double> &, const class numeric::xyzVector<double> &, int &, int &, int &, int &, double) –> void

pyrosetta.rosetta.core.scoring.sasa.get_legrand_atomic_overlap(radius_a: float, radius_b: float, distance_ijxyz: float, degree_of_overlap: int) None
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

C++: core::scoring::sasa::get_legrand_atomic_overlap(const double, const double, const double, int &) –> void

pyrosetta.rosetta.core.scoring.sasa.get_legrand_orientation(a_xyz: pyrosetta.rosetta.numeric.xyzVector_double_t, b_xyz: pyrosetta.rosetta.numeric.xyzVector_double_t, phi_index: int, theta_index: int, distance_ijxyz: float) None
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)

C++: core::scoring::sasa::get_legrand_orientation(const class numeric::xyzVector<double> &, const class numeric::xyzVector<double> &, int &, int &, double) –> void

pyrosetta.rosetta.core.scoring.sasa.get_legrand_sasa_angles() pyrosetta.rosetta.ObjexxFCL.FArray2D_int_t
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)

C++: core::scoring::sasa::get_legrand_sasa_angles() –> const class ObjexxFCL::FArray2D<int> &

pyrosetta.rosetta.core.scoring.sasa.get_legrand_sasa_masks() pyrosetta.rosetta.ObjexxFCL.FArray2D_ObjexxFCL_ubyte_t
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)

C++: core::scoring::sasa::get_legrand_sasa_masks() –> const class ObjexxFCL::FArray2D<class ObjexxFCL::ubyte> &

pyrosetta.rosetta.core.scoring.sasa.get_sasa_method_from_string(method: str) pyrosetta.rosetta.core.scoring.sasa.SasaMethodEnum

Gets sasa enum from string passed by options system.

C++: core::scoring::sasa::get_sasa_method_from_string(std::string) –> enum core::scoring::sasa::SasaMethodEnum

pyrosetta.rosetta.core.scoring.sasa.get_sasa_radii_parameter_name(radii_set: pyrosetta.rosetta.core.scoring.sasa.SasaRadii) str

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

C++: core::scoring::sasa::get_sasa_radii_parameter_name(enum core::scoring::sasa::SasaRadii) –> std::string

pyrosetta.rosetta.core.scoring.sasa.get_sasa_radii_set_from_string(radii_set: str) pyrosetta.rosetta.core.scoring.sasa.SasaRadii

Gets sasa radii enum from string passed by options system.

C++: core::scoring::sasa::get_sasa_radii_set_from_string(std::string) –> enum core::scoring::sasa::SasaRadii

pyrosetta.rosetta.core.scoring.sasa.get_sc_bb_sasa(pose: pyrosetta.rosetta.core.pose.Pose, atom_sasa: pyrosetta.rosetta.core.id.AtomID_Map_double_t) Tuple[float, float]

Calculate the sidechain and backbone sasa from atom sasa

C++: core::scoring::sasa::get_sc_bb_sasa(const class core::pose::Pose &, const class core::id::AtomID_Map<double> &) –> struct std::pair<double, double>

pyrosetta.rosetta.core.scoring.sasa.get_sc_bb_sasa_per_res(pose: pyrosetta.rosetta.core.pose.Pose, atom_sasa: pyrosetta.rosetta.core.id.AtomID_Map_double_t) Tuple[pyrosetta.rosetta.utility.vector1_double, pyrosetta.rosetta.utility.vector1_double]

C++: core::scoring::sasa::get_sc_bb_sasa_per_res(const class core::pose::Pose &, const class core::id::AtomID_Map<double> &) –> struct std::pair<class utility::vector1<double, class std::allocator<double> >, class utility::vector1<double, class std::allocator<double> > >

pyrosetta.rosetta.core.scoring.sasa.input_legrand_sasa_dats() None

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

C++: core::scoring::sasa::input_legrand_sasa_dats() –> void

pyrosetta.rosetta.core.scoring.sasa.is_res_exposed(pose: pyrosetta.rosetta.core.pose.Pose, resnum: int) bool

Is residue exposed?

Added by JKLeman (julia.koehler1982.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!

C++: core::scoring::sasa::is_res_exposed(const class core::pose::Pose &, unsigned long) –> bool

pyrosetta.rosetta.core.scoring.sasa.per_res_sc_sasa(pose: pyrosetta.rosetta.core.pose.Pose) pyrosetta.rosetta.utility.vector1_double

Absolute per residue sidechain SASA

Added by JKLeman (julia.koehler1982.com)

GXG tripeptide values for sidechain SASA are taken from http://www.proteinsandproteomics.org/content/free/tables_1/table08.pdf

C++: core::scoring::sasa::per_res_sc_sasa(const class core::pose::Pose &) –> class utility::vector1<double, class std::allocator<double> >

pyrosetta.rosetta.core.scoring.sasa.rel_per_res_sc_sasa(pose: pyrosetta.rosetta.core.pose.Pose) pyrosetta.rosetta.utility.vector1_double

Relative per residue sidechain SASA

Added by JKLeman (julia.koehler1982.com)

GXG tripeptide values for sidechain SASA are taken from http://www.proteinsandproteomics.org/content/free/tables_1/table08.pdf

C++: core::scoring::sasa::rel_per_res_sc_sasa(const class core::pose::Pose &) –> class utility::vector1<double, class std::allocator<double> >