Title: | Occlusion Surface Using the Occluded Surface and Fibonacci Occluded Surface |
Version: | 2.0.1 |
Maintainer: | Herson Soares <hersonhebert@hotmail.com> |
Description: | The Occluded Surface (OS) algorithm is a widely used approach for analyzing atomic packing in biomolecules as described by Pattabiraman N, Ward KB, Fleming PJ (1995) <doi:10.1002/jmr.300080603>. Here, we introduce 'fibos', an 'R' and 'Python' package that extends the 'OS' methodology, as presented in Soares HHM, Romanelli JPR, Fleming PJ, da Silveira CH (2024) <doi:10.1101/2024.11.01.621530>. |
Depends: | R (≥ 4.1.0) |
License: | GPL-3 |
Encoding: | UTF-8 |
RoxygenNote: | 7.3.2 |
Imports: | fs, dplyr, readr, stringr, tidyr, reticulate, glue |
NeedsCompilation: | no |
Packaged: | 2025-06-03 23:47:24 UTC; hersonhebert |
Author: | Herson Soares [cre, aut], Carlos Silveira [aut], João Romanelli [aut], Patrick Fleming [aut], Posit Software, PBC [cph, fnd] |
Repository: | CRAN |
Date/Publication: | 2025-06-04 00:00:02 UTC |
Surface Calculation
Description
The function executes the implemented methods. Using this function, it is possible to calculate occluded areas through the traditional methodology, 'Occluded Surface', or by applying the 'Fibonacci OS' methodology. At the end of the method execution, the 'prot.srf' file is generated, and returned for the function. The data in this file refers to all contacts between atoms of molecules present in a protein's PDB.
Usage
execute(pdb, method, density_dots)
Arguments
pdb |
4-digit PDB id (will fetch it from the RCSB repository) or the path to a PDB local file. |
method |
Method to be used: 'OS' (classic) or 'FIBOS' (default). The classic 'OS' covers the surface radially with one of the axes as a reference when allocating the dots. In 'FIBOS', Fibonacci spirals were used to allocate the dots, which is known to produce lower axial anisotropy as well as more evenly spaced points on a sphere. |
density_dots |
Distribution density of atomic dots for surface occlusion calculation. |
Value
A data frame containing the protein surface data read from the generated SRF file.
Author(s)
Carlos Henrique da Silveira (carlos.silveira@unifei.edu.br)
Herson Hebert Mendes Soares (hersonhebert@hotmail.com)
Joao Paulo Roquim Romanelli (joaoromanelli@unifei.edu.br)
Patrick Fleming (Pat.Fleming@jhu.edu)
Install the 'Python' 'fibos' Module
Description
This function creates a 'Python' virtual environment and installs the 'Python' module 'fibos' required for the full functionality of this package. It handles different system configurations and ensures that the correct compiler paths are set.
Usage
fibos_config()
Note
This function will install external software (a 'Python' package) on your system. Administrator/sudo privileges might be required on some systems.
See Also
Examples
# Set up the 'Python' environment and install the required module
fibos_config()
Load Radii Values
Description
The 'get_radii' function is responsible for loading the atomic radii values used for surface-occlusion calculations. The values it returns are those currently employed in those calculations.
Usage
get_radii()
Value
A data frame containing the radii values.
Author(s)
Carlos Henrique da Silveira (carlos.silveira@unifei.edu.br)
Herson Hebert Mendes Soares (hersonhebert@hotmail.com)
Joao Paulo Roquim Romanelli (joaoromanelli@unifei.edu.br)
Patrick Fleming (Pat.Fleming@jhu.edu)
See Also
Examples
library(fibos)
fibos_config()
#Loads the radius values that have been configured for code execution.
radii = get_radii()
#Displays the first three lines.
radii |> utils::head(3) |> print()
Occluded Surface (OS)
Description
The 'Occluded Surface (OS)' algorithm is a widely used approach for analyzing atomic packing in biomolecules. Here, we introduce 'FIBOS', an 'R' and 'Python' package that extends the 'OS' methodology with enhancements. The homonymous function 'occluded_surface' calculates 'OS' per atom.
Usage
occluded_surface(pdb, method = "FIBOS", density_dots = 5)
Arguments
pdb |
4-digit PDB id (will fetch it from the RCSB repository) or the path to a PDB local file. |
method |
Method to be used: 'OS' (classic) or 'FIBOS'(default).The classic 'OS' covers the surface radially with one of the axes as a reference when allocating the dots. In 'FIBOS', Fibonacci spirals were used to allocate the dots, which is known to produce lower axial anisotropy as well as more evenly spaced points on a sphere. |
density_dots |
Distribution density of atomic dots for surface occlusion calculation. |
Details
'Occluded Surface (OS)' (Pattabiraman et al. 1995) method distributes dots (representing patches of area) across the atom surfaces. Each dot has a normal that extends until it reaches either a van der Waals surface of a neighboring atom (the dot is considered occluded) or covers a distance greater than the diameter of a water molecule (the dot is considered non-occluded and disregarded). Thus, with the summed areas of dots and the lengths of normals, it is possible to compose robust metrics capable of inferring the average packing density of atoms, residues, proteins, as well as any other group of biomolecules.
For more details, see (Fleming et al, 2000) and (Soares, et al, 2024)
Value
A table containing:
ATOM
the atomic contacts for each atom.
NUMBER OF POINTS
the number of dots (patches of area) on atomic surface.
AREA
the summed areas of dots.
RAYLENGTH
the average lengths of normals normalized by 2.8
\text{\AA}
(water diameter). So, raylen is a value between 0 and 1. A raylen close to 1 indicates worse packaging.DISTANCE
the average distances of contacts in (
\text{\AA}
).
Author(s)
Herson Soares, Joao Romanelli, Patrick Fleming, Carlos Silveira.
References
Fleming PJ, Richards FM (2000). "Protein packing: Dependence on protein size, secondary structure and amino acid composition." doi:10.1006/jmbi.2000.3750
Pattabiraman N, Ward KB, Fleming PJ (1995). "Occluded molecular surface: Analysis of protein packing." doi:10.1002/jmr.300080603
Soares HHM, Romanelli JPR, Fleming PJ, da Silveira CH (2024). "bioRxiv, 2024.11.01.621530." doi:10.1101/2024.11.01.621530
See Also
Examples
library(fibos)
#Configure the environment
fibos_config()
# Calculate FIBOS per atom and create .srf files in fibos_files folder
pdb_fibos <- occluded_surface("1ptx", method = "FIBOS", density_dots = 5.0)
# Calculate OSP metric per residue from .srf file in fibos_files folder
pdb_osp <- osp(fs::path("fibos_files","prot_1ptx.srf"))
Occluded Surface Packing (OSP)
Description
Implements the 'occluded surface' packing density metric (OSP) averaged by residue, as described in (Fleming and Richards 2000).
Usage
osp(file)
Arguments
file |
a SRF File (.srf) generated by 'occluded_surface' in fibos_files folder. |
Value
A table containing:
Resnum
residue id.
Resname
residue name.
OS
the summed areas of dots in residue.
`os*[1-raylen]`
'OS' areas weighted by (1-raylen). Raylen is the average lengths of normals normalized by 2.8
\text{\AA}
(water diameter). So, raylen is a value between 0 and 1. A raylen close to 1 indicates worse packaging, and the 'OS' will be reduced.OSP
average occluded surface packing value (OSP) by residue.
Author(s)
Herson Soares
Joao Romanelli
Patrick Fleming
Carlos Silveira.
References
Fleming PJ, Richards FM (2000). "Protein packing: Dependence on protein size, secondary structure and amino acid composition." doi:10.1006/jmbi.2000.3750
Pattabiraman N, Ward KB, Fleming PJ (1995). "Occluded molecular surface: Analysis of protein packing." doi:10.1002/jmr.300080603
Soares HHM, Romanelli JPR, Fleming PJ, da Silveira CH (2024). "bioRxiv, 2024.11.01.621530." doi:10.1101/2024.11.01.621530
See Also
Examples
library(fibos)
#Configure the Environment
fibos_config()
# Calculate FIBOS per atom and create .srf files in fibos_files folder
pdb_fibos <- occluded_surface("1ptx", method = "FIBOS", density_dots = 5.0)
# Calculate OSP metric per residue from .srf file in fibos_files folder
pdb_osp <- osp(fs::path("fibos_files","prot_1ptx.srf"))
Calculate OSP Metric
Description
Internal function to calculate the OSP metric from an SRF file. This function calls the 'Python' implementation of the OSP calculation.
Usage
osp_internal(file)
Arguments
file |
Path to the SRF file generated by the 'occluded_surface' function. |
Value
A data frame containing the OSP metric data read from the generated PAK file.
Author(s)
Carlos Henrique da Silveira (carlos.silveira@unifei.edu.br)
Herson Hebert Mendes Soares (hersonhebert@hotmail.com)
Joao Paulo Roquim Romanelli (joaoromanelli@unifei.edu.br)
Patrick Fleming (Pat.Fleming@jhu.edu)
Read OSP Values
Description
Internal function to read OSP values from a PAK file generated by the OSP calculation.
Usage
read_osp(prot_file)
Arguments
prot_file |
OSP File (.pak). |
Value
A data frame containing the OSP data.
Read Protein Surface File
Description
Reads and processes data from protein surface files generated by the 'occluded_surface' function. The results enable important analyses of protein behavior at the atomic level.
Arguments
file |
Path or name of protein surface file (.srf). |
Details
The function reads a surface file containing data calculated using the 'OS' or 'FIBOS' methodology. The file should be in the format generated by the 'occluded_surface' function.
Value
A data frame with the following columns:
- ATOM
atomic contact identifiers
- NUMBER_POINTS
number of points (integer) on atomic surface
- AREA
surface area (numeric)
- RAYLENGTH
average ray length (numeric)
- DISTANCE
average distance of contacts (numeric)
Author(s)
Carlos Henrique da Silveira (carlos.silveira@unifei.edu.br)
Herson Hebert Mendes Soares (hersonhebert@hotmail.com)
Joao Paulo Roquim Romanelli (joaoromanelli@unifei.edu.br)
Patrick Fleming (Pat.Fleming@jhu.edu)
Reset Radii Values
Description
This function reloads the 'OS' default atomic radii values.
Usage
reset_radii()
Author(s)
Carlos Henrique da Silveira (carlos.silveira@unifei.edu.br)
Herson Hebert Mendes Soares (hersonhebert@hotmail.com)
Joao Paulo Roquim Romanelli (joaoromanelli@unifei.edu.br)
Patrick Fleming (Pat.Fleming@jhu.edu)
See Also
Examples
library(fibos)
fibos_config()
#Loads the radius values that have been configured for code execution.
radii = get_radii()
#Displays the first three lines.
radii |> utils::head(3) |> print()
#Modifies the value of a specific radius.
radii$RAY[1] = 2.15
#Sets the radius value from a tibble.
set_radii(radii)
#Displays the first three lines.
get_radii() |> utils::head(3) |> print()
#Loads the default radius values.
reset_radii()
#Displays the first three lines.
get_radii() |> utils::head(3) |> print()
Change Radii Values
Description
This function enables modification of the radius values by passing a 'data.frame' as an argument.
Usage
set_radii(radii_values)
Arguments
radii_values |
A 'data.frame' containing atomic radii values. |
Author(s)
Carlos Henrique da Silveira (carlos.silveira@unifei.edu.br)
Herson Hebert Mendes Soares (hersonhebert@hotmail.com)
Joao Paulo Roquim Romanelli (joaoromanelli@unifei.edu.br)
Patrick Fleming (Pat.Fleming@jhu.edu)
See Also
Examples
library(fibos)
fibos_config()
#Loads the radius values that have been configured for code execution.
radii = get_radii()
#Displays the first three lines.
radii |> utils::head(3) |> print()
#Modifies the value of a specific radius.
radii$RAY[1] = 2.15
#Sets the radius value from a tibble.
set_radii(radii)
#Displays the first three lines.
get_radii() |> utils::head(3) |> print()