SolidAngleR: Computing normalized solid angles of polyhedral cones
Source:R/SolidAngleR-package.R
SolidAngleR-package.RdThe SolidAngleR package provides comprehensive methods for computing normalized solid angle measures of polyhedral cones in arbitrary dimensions. It implements state-of-the-art algorithms from computational geometry, offering multiple approaches tailored to different cone properties and dimensions.
Details
Overview
Solid angle measures quantify the "size" of a polyhedral cone as a fraction of the total space. This package provides closed formulas for 2D and 3D cones (Euler-Lagrange formula); hypergeometric series for positive definite cones (Ribando 2006); tridiagonal series for cones with special structure (Fitisone & Zhou 2023); decomposition methods for general cones (Fitisone & Zhou 2023); multivariate normal integration for ecological applications (Song et al. 2018); Monte Carlo estimation for complex geometries; geometric methods for circular cones and intersecting cones (Mazonka 2012); and spanning tree methods for extreme ray computation (Malekmohammadi & Mostafaee 2017).
Main functions
Unified interface
compute_solid_angleMain function with automatic method selection
diagnose_coneDiagnostic tool to suggest best method
Dimension-specific formulas
solid_angle_2d2D planar angle computation
solid_angle_3d3D Euler-Lagrange formula
spherical_triangle_areaSpherical triangle area
Series methods
hypergeometric_seriesRibando's hypergeometric series (Theorem 1.5)
tridiagonal_seriesOptimized series for tridiagonal cones (Theorem 4.1)
Geometric methods
solid_angle_coneRight circular cone
solid_angle_polyhedralVan Oosterom & Strackee formula
solid_angle_intersecting_conesIntersecting circular cones
plotIntersectingConesInteractive 3D visualization
Extreme rays
generate_spanning_treesSpanning tree generation
solid_angle_3d_from_raysSolid angle from extreme rays
Utilities
compute_associated_matrixAssociated matrix M_n(C)
is_positive_definitePositive definiteness check
create_tridiagonal_coneConstruct tridiagonal cone
angle_betweenAngle between vectors
lhuilier_angleL'Huilier's theorem for spherical triangles
Mathematical background
Normalized solid angle
The normalized solid angle measure \(\omega\) quantifies the fraction of \(n\)-dimensional space occupied by a polyhedral cone. In 2D, \(\omega = \theta/(2\pi)\) where \(\theta\) is the planar angle; in 3D, \(\omega = \Omega/(4\pi)\) where \(\Omega\) is the solid angle in steradians; and in n-D, \(\omega\) is the fraction of the unit sphere covered by the cone.
Key properties
The normalized solid angle \(\omega\) has range \(\omega \in [0, 1]\); for an orthant (axes-aligned cone), \(\omega = 1/2^n\); for full space, \(\omega = 1\); and for an empty cone, \(\omega = 0\).
Associated matrix
For a simplicial cone \(C\) generated by unit vectors \(v_1, \ldots, v_n\), the associated matrix \(M_n(C)\) is defined as:
$$M_n(C)_{ij} = \begin{cases} 1 & \text{if } i = j \\ -|v_i \cdot v_j| & \text{if } i \neq j \end{cases}$$
The hypergeometric series converges if and only if \(M_n(C)\) is positive definite (Theorem 1.5, Ribando 2006).
Method selection guide
| Cone property | Dimension | Method | Function |
| Any | 2 | Closed formula | solid_angle_2d() |
| Any | 3 | Euler-Lagrange | solid_angle_3d() |
| Positive definite | Any | Hypergeometric series | method = "series" |
| Tridiagonal V'V | Any | Tridiagonal series | method = "tridiagonal" |
| Any | Any | Decomposition | method = "decomposition" |
| Any | Any | Multivariate normal | Omega() |
Use method = "auto" in compute_solid_angle for automatic
selection, or diagnose_cone for detailed recommendations.
Typical workflows
Workflow 1: General cone computation
library(SolidAngleR)
# Define cone by generator matrix
V <- diag(4) # 4D orthant
# Automatic method selection
omega <- compute_solid_angle(V)
print(omega) # 0.0625 (1/16)Workflow 3: Geometric shapes
# Circular cone
theta <- pi/3 # 60-degree apex
omega <- solid_angle_cone(theta)
# Intersecting cones
omega_int <- solid_angle_intersecting_cones(
theta1 = pi/3,
theta2 = pi/4,
alpha = pi/6
)
# Visualize
fig <- plotIntersectingCones(
theta1 = pi/3,
theta2 = pi/4,
phi = pi/3
)
figWorkflow 4: Method comparison
V <- diag(3)
# Compare methods
omega_formula <- compute_solid_angle(V, method = "formula")
omega_series <- compute_solid_angle(V, method = "series")
omega_decomp <- compute_solid_angle(V, method = "decomposition")
omega_mvn <- Omega(V)$Omega
# All should give 0.125
c(omega_formula, omega_series, omega_decomp, omega_mvn)Applications
Computational geometry
Applications include polyhedral cone volume computation; spherical geometry and geodesic calculations; and solid angle subtended by complex 3D objects.
Note
This package is the original creation of the author in all conceptual, theoretical, and design aspects. Implementation was assisted by Anthropic's Claude Code v.2 (Opus 4.5) to streamline package development. All original ideas, creativity, and scientific contributions belong to the author, who maintains full responsibility for the package's correctness and reliability. While all code has been thoroughly reviewed and tested, users are encouraged to report any issues through the package's issue tracker.
References
Primary references
Ribando, J. M. (2006). Measuring solid angles beyond dimension three. Discrete & Computational Geometry, 36(3), 479-487. doi:10.1007/s00454-006-1253-4
Fitisone, A., & Zhou, Y. (2023). Solid angle measure of polyhedral cones. arXiv:2304.11102 (math.CO). https://arxiv.org/abs/2304.11102
Mazonka, O. (2012). Solid angle of conical surfaces, polyhedral cones, and intersecting spherical caps. arXiv:1205.1396v2 (math.MG). https://arxiv.org/abs/1205.1396
Malekmohammadi, N., & Mostafaee, A. (2017). Obtaining all extreme rays of a special cone using spanning trees in a complete digraph: application in DEA. Journal of the Operational Research Society, 69(3), 465-472. doi:10.1057/s41274-017-0265-9
Song, C., Rohr, R. P., & Saavedra, S. (2018). A guideline to study the feasibility domain of multi-trophic and changing ecological communities. Journal of Theoretical Biology, 450, 30-36. doi:10.1016/j.jtbi.2018.04.030
Additional references
Van Oosterom, A., & Strackee, J. (1983). The solid angle of a plane triangle. IEEE Transactions on Biomedical Engineering, BME-30(2), 125-126. doi:10.1109/TBME.1983.325207
Arun, I., & Venkatapathi, M. (2025). An O(n) algorithm for generating uniform random vectors in n-dimensional cones. Sankhya A: The Indian Journal of Statistics, 87(2), 327-348. doi:10.1007/s13171-025-00387-9
Aomoto, K. (1977). Analytic structure of Schlafli function. Nagoya Mathematical Journal, 68, 1-16. https://projecteuclid.org/euclid.nmj/1118786429
Girard, A. (1626). Tables des sinus, tangentes, & secantes. Leiden.
See also
Package website: https://robustecologies.github.io/SolidAngleR
GitHub repository: https://github.com/robustecologies/SolidAngleR
Bug reports: https://github.com/robustecologies/SolidAngleR/issues
Vignettes:
vignette("mathematical-theory", package = "SolidAngleR")vignette("mazonka-geometric-methods", package = "SolidAngleR")vignette("spanning-trees-extreme-rays", package = "SolidAngleR")
Author
Pablo Almaraz pablo.almaraz@csic.es (ORCID)
Examples
# Quick start examples
# ========================================================================== #
# Example 1: 2D cone ####
# ========================================================================== #
v1 <- c(1, 0)
v2 <- c(1, 1) / sqrt(2)
V_2d <- cbind(v1, v2)
omega_2d <- compute_solid_angle(V_2d)
print(omega_2d) # 0.125 (45° / 360°)
#> [1] 0.125
# ========================================================================== #
# Example 2: 3D orthogonal cone (octant) ####
# ========================================================================== #
V_3d <- diag(3)
omega_3d <- compute_solid_angle(V_3d)
print(omega_3d) # 0.125 (1/8 of space)
#> [1] 0.125
# ========================================================================== #
# Example 3: Diagnose cone properties ####
# ========================================================================== #
V <- diag(4)
diag_report <- diagnose_cone(V)
print(diag_report)
#> Cone diagnostic report
#> ======================
#>
#> Dimension: 4
#> Linearly independent: ✓ Yes
#> Associated matrix PD: ✓ Yes
#> Tridiagonal structure: ✓ Yes
#>
#> Eigenvalues: 1, 1, 1, 1
#> Minimum eigenvalue: 1
#>
#> Suggested method: tridiagonal series (most efficient) ✓
# ========================================================================== #
# Example 4: Circular cone ####
# ========================================================================== #
theta <- pi/3 # 60-degree apex
omega_cone <- solid_angle_cone(theta)
print(omega_cone)
#> [1] 0.25
# ========================================================================== #
# Example 5: Ecological feasibility domain ####
# ========================================================================== #
if (FALSE) { # \dontrun{
# 4-species community
S <- 4
alpha <- matrix(runif(S * S, -1, 0), S, S)
diag(alpha) <- -1
result <- Omega(alpha)
cat("Feasibility volume:", result$Omega, "\n")
cat("Per-species measure:", result$Omega_scaled, "\n")
} # }