molecular_simulations.analysis.sasa module¶
Solvent-Accessible Surface Area (SASA) analysis module.
This module computes SASA using the Shrake-Rupley algorithm, implemented as MDAnalysis AnalysisBase classes for ease of use and built-in parallelism. Adapted from BioPython and MDTraj implementations.
- class molecular_simulations.analysis.sasa.SASA(ag, probe_radius=1.4, n_points=256, **kwargs)[source]¶
Bases:
AnalysisBaseCompute solvent-accessible surface area using Shrake-Rupley algorithm.
Implements SASA calculation as an MDAnalysis AnalysisBase instance, supporting ease of deployment and built-in parallelism. Code is adapted from BioPython and MDTraj implementations as well as an unmerged PR from MDAnalysis.
- ag¶
The AtomGroup being analyzed.
- probe_radius¶
Probe radius for SASA calculation.
- n_points¶
Number of points on each atomic sphere.
- radii¶
Array of atomic radii (VDW + probe).
- max_radii¶
Maximum pairwise radius sum.
- results.sasa¶
Per-residue SASA values after run().
- Parameters:
ag (
AtomGroup) – AtomGroup for which to compute SASA.probe_radius (
float) – Probe radius in Angstroms. Defaults to 1.4 Å (standard water probe).n_points (
int) – Number of points on each sphere. Higher values are more accurate but slower. Defaults to 256.**kwargs – Additional arguments passed to AnalysisBase.
- Raises:
TypeError – If ag is an UpdatingAtomGroup.
ValueError – If the Universe has no ‘elements’ property.
Example
>>> u = mda.Universe("system.prmtop", "traj.dcd") >>> sasa = SASA(u.select_atoms("protein")) >>> sasa.run() >>> print(sasa.results.sasa)
Initialize the SASA analysis.
- Parameters:
- class molecular_simulations.analysis.sasa.RelativeSASA(ag, probe_radius=1.4, n_points=256, **kwargs)[source]¶
Bases:
SASACompute relative SASA for an AtomGroup.
Relative SASA is defined as measured SASA divided by maximum accessible surface area. For proteins, this is computed in a tripeptide context to avoid overestimating SASA of the amide linkage and its neighbors.
- results.sasa¶
Absolute SASA values.
- results.relative_area¶
Relative SASA values (0-1 scale).
- Parameters:
- Raises:
ValueError – If the Universe has no ‘bonds’ property.
Example
>>> u = mda.Universe("system.prmtop", "traj.dcd") >>> rsasa = RelativeSASA(u.select_atoms("protein")) >>> rsasa.run() >>> print(rsasa.results.relative_area)
Initialize the RelativeSASA analysis.
- Parameters: