Skip to main content

Core API

Basic SASA calculation functions.

calculate_sasa

def calculate_sasa(
coords: NDArray[np.float64],
radii: NDArray[np.float64],
*,
algorithm: Literal["sr", "lr"] = "sr",
n_points: int = 100,
n_slices: int = 20,
probe_radius: float = 1.4,
n_threads: int = 0,
use_bitmask: bool = False,
) -> SasaResult

Calculate Solvent Accessible Surface Area.

Parameters:

ParameterTypeDefaultDescription
coordsNDArray[float64]requiredAtom coordinates as (N, 3) array
radiiNDArray[float64]requiredAtom radii as (N,) array
algorithm"sr" or "lr""sr"Shrake-Rupley or Lee-Richards
n_pointsint100Test points per atom (SR only)
n_slicesint20Slices per atom (LR only)
probe_radiusfloat1.4Water probe radius in Å
n_threadsint0Number of threads (0 = auto)
use_bitmaskboolFalseUse bitmask LUT optimization (SR only, n_points must be 1..1024)

Returns: SasaResult

Raises: ValueError for invalid input

SasaResult

@dataclass
class SasaResult:
total_area: float # Total SASA in Ų
atom_areas: NDArray[float64] # Per-atom SASA values

Batch API

For trajectory analysis with multiple frames.

calculate_sasa_batch

def calculate_sasa_batch(
coordinates: NDArray[np.floating],
radii: NDArray[np.floating],
*,
algorithm: Literal["sr", "lr"] = "sr",
n_points: int = 100,
n_slices: int = 20,
probe_radius: float = 1.4,
n_threads: int = 0,
precision: Literal["f64", "f32"] = "f64",
use_bitmask: bool = False,
) -> BatchSasaResult

Calculate SASA for multiple frames in batch.

Parameters:

ParameterTypeDefaultDescription
coordinatesNDArray[floating]requiredCoordinates as (n_frames, n_atoms, 3) array
radiiNDArray[floating]requiredAtom radii as (n_atoms,) array
algorithm"sr" or "lr""sr"Shrake-Rupley or Lee-Richards
n_pointsint100Test points per atom (SR only)
n_slicesint20Slices per atom (LR only)
probe_radiusfloat1.4Water probe radius in Å
n_threadsint0Number of threads (0 = auto)
precision"f64" or "f32""f64"Floating-point precision
use_bitmaskboolFalseUse bitmask LUT optimization (SR only, n_points must be 1..1024)

Returns: BatchSasaResult

BatchSasaResult

@dataclass
class BatchSasaResult:
atom_areas: NDArray[float32] # Per-atom SASA, shape (n_frames, n_atoms)

# Properties
n_frames: int # Number of frames
n_atoms: int # Number of atoms
total_areas: NDArray[float32] # Total SASA per frame, shape (n_frames,)

Precision Parameter

The precision parameter controls the floating-point precision used internally:

PrecisionDescriptionUse Case
"f64" (default)64-bit double precisionHigher accuracy, recommended
"f32"32-bit single precisionComparison with RustSASA/mdsasa-bolt

Note: The difference between f64 and f32 is typically < 0.01% for SASA calculations. f64 is recommended unless you need exact comparison with f32-based implementations.

Example:

import numpy as np
from zsasa import calculate_sasa_batch

# Multiple frames
n_frames = 100
n_atoms = 1000
coords = np.random.rand(n_frames, n_atoms, 3).astype(np.float32) * 50
radii = np.full(n_atoms, 1.5, dtype=np.float32)

# f64 precision (default)
result = calculate_sasa_batch(coords, radii, n_threads=4)
print(f"Shape: {result.atom_areas.shape}") # (100, 1000)

# f32 precision (for comparison with RustSASA)
result_f32 = calculate_sasa_batch(coords, radii, precision="f32")

Directory Batch Processing

Process all supported structure files in a directory at once.

process_directory

def process_directory(
input_dir: str | Path,
*,
output_dir: str | Path | None = None,
algorithm: Literal["sr", "lr"] = "sr",
n_points: int = 100,
n_slices: int = 20,
probe_radius: float = 1.4,
n_threads: int = 0,
classifier: ClassifierType | None = ClassifierType.PROTOR,
include_hydrogens: bool = False,
include_hetatm: bool = False,
) -> BatchDirResult

Process all supported structure files (.pdb, .cif, .mmcif, .ent, .json, and .gz variants) in a directory.

Parameters:

ParameterTypeDefaultDescription
input_dirstr | PathrequiredDirectory containing structure files
output_dirstr | Path | NoneNoneOutput directory for per-file results
algorithm"sr" or "lr""sr"Shrake-Rupley or Lee-Richards
n_pointsint100Test points per atom (SR only)
n_slicesint20Slices per atom (LR only)
probe_radiusfloat1.4Water probe radius in Å
n_threadsint0Number of threads (0 = auto)
classifierClassifierType | NonePROTORClassifier for radii. None = use input radii
include_hydrogensboolFalseInclude hydrogen atoms
include_hetatmboolFalseInclude HETATM records

Returns: BatchDirResult

Raises:

ExceptionCause
ValueErrorInvalid parameters (algorithm, probe_radius, n_points, etc.)
FileNotFoundErrorDirectory does not exist or is not readable
MemoryErrorOut of memory
RuntimeErrorCalculation or other processing error

BatchDirResult

@dataclass
class BatchDirResult:
total_files: int # Number of supported files found
successful: int # Number successfully processed
failed: int # Number that failed
filenames: list[str] # Per-file basenames
n_atoms: list[int] # Per-file atom counts
total_sasa: list[float] # Per-file SASA in Ų (NaN if failed)
status: list[int] # Per-file status (1=ok, 0=failed)

Example:

from zsasa import process_directory

# Process all structure files in a directory
result = process_directory("path/to/structures/")
print(f"Processed {result.successful}/{result.total_files} files")

# Per-file results
for i in range(result.total_files):
if result.status[i] == 1:
print(f" {result.filenames[i]}: {result.n_atoms[i]} atoms, "
f"SASA={result.total_sasa[i]:.1f} Ų")

# With output files and LR algorithm
result = process_directory(
"input/",
output_dir="output/",
algorithm="lr",
n_threads=4,
)