pub enum RayCastAlgorithm<const N: usize> {
BracketAndBisect {
lambda: f32,
d_lambda: f32,
max_bisection_iterations: usize,
},
Analytical,
Newton {
max_iterations: usize,
nudge_distance: f32,
step_size: f32,
},
}Expand description
Ray-surface intersection algorithms for implicit surface ray casting.
This enum provides multiple algorithmic approaches for solving ray-surface intersection problems, each optimized for different scenarios in safety-critical applications where precision and reliability are paramount.
§Algorithm Selection Guidelines
-
Use
Analyticalfor:- Simple primitive shapes (spheres, planes, cylinders)
- Maximum accuracy and performance requirements
- Real-time applications with tight computational budgets
-
Use
Newtonfor:- Smooth, well-conditioned implicit surfaces
- When autodifferentiation is available
- Applications requiring fast convergence with few iterations
-
Use
BracketAndBisectfor:- Complex CSG constructions with sharp edges
- Robustness over speed requirements
- Surfaces with potential numerical instabilities
For theoretical background, see the Ray Casting chapter.
Variants§
BracketAndBisect
Robust numerical method using bracketing followed by bisection.
This algorithm first advances along each ray to bracket sign changes in the scalar field, then uses bisection to refine intersection points to high precision. It’s the most reliable method for complex CSG constructions.
§Parameters
lambda- Maximum ray extension distance for bracketing phased_lambda- Step size during bracketing searchmax_bisection_iterations- Maximum refinement iterations (typically 20-50)
§Convergence
Guaranteed to converge for continuous scalar fields with precision ε ≈ λ/2ⁿ where n is the number of bisection iterations.
Fields
Analytical
Direct analytical solutions for primitive geometric shapes.
Computes exact ray-surface intersections using closed-form mathematical formulas. This is the fastest and most accurate method when applicable, but is limited to simple primitives with known intersection equations.
§Supported Primitives
- Spheres: Quadratic equation solutions
- Planes: Linear ray-plane intersection
- Cylinders: Analytical cylinder intersection
- Cones: Quadratic cone intersection formulas
§Limitations
Cannot be used with:
- Arbitrary scalar field expressions
- CSG operations beyond simple unions
- Transformed or deformed primitives
Newton
Iterative Newton-Raphson root finding using automatic differentiation.
Uses gradient information to rapidly converge to surface intersections. Exhibits quadratic convergence near solutions but requires smooth, differentiable scalar fields for optimal performance.
§Parameters
max_iterations- Maximum Newton iterations (typically 5-20)nudge_distance- Perturbation for degenerate gradient casesstep_size- Relaxation factor for Newton steps (usually 0.5-1.0)
§Convergence Properties
- Quadratic: Error ∝ ε² per iteration near solutions
- Sensitive: May fail near sharp edges or gradient discontinuities
- Fast: Typically converges in 3-10 iterations for smooth surfaces
Implementations§
Source§impl<const N: usize> RayCastAlgorithm<N>
impl<const N: usize> RayCastAlgorithm<N>
Sourcepub fn solve<B: AutodiffBackend>(
&self,
rays: Rays<B, N>,
region: &Region<N, B>,
algebra: &impl CSGAlgebra<B>,
) -> RayCastResult<B, N>
pub fn solve<B: AutodiffBackend>( &self, rays: Rays<B, N>, region: &Region<N, B>, algebra: &impl CSGAlgebra<B>, ) -> RayCastResult<B, N>
Computes ray-surface intersections using the specified algorithm.
This method finds the points where rays intersect with the boundary of an implicit surface region. It supports multiple algorithmic approaches, each with different trade-offs between accuracy, performance, and applicability.
§Mathematical Problem
Given rays r(t) = o + td and a region R with characteristic function f(x), find parameter values t where f(r(t)) = 0, indicating surface intersection.
§Algorithm Variants
-
Analytical: Direct closed-form solutions for primitive shapes (spheres, planes, etc.)
- Fastest and most accurate when applicable
- Limited to simple geometric primitives with known ray intersection formulas
-
Newton: Iterative root-finding using gradient-based optimization
- Quadratic convergence near solutions
- Requires autodifferentiable scalar fields
- Best for smooth, well-conditioned surfaces
-
BracketAndBisect: Robust bracketing followed by bisection refinement
- Guaranteed convergence for continuous functions
- Works with any scalar field implementation
- Slower but most reliable for complex CSG constructions
§Parameters
rays- Collection of rays to cast against the surfaceregion- Implicit surface region to intersect withalgebra- CSG algebra for evaluating boolean operations in composite regions
§Returns
RayCastResult containing intersection points, distances, field values, and metadata
for each input ray. Rays that miss the surface are marked with sentinel values.
§Performance Notes
- Analytical: O(1) per ray for primitive surfaces
- Newton: O(k) where k is iteration count (typically 3-10)
- BracketAndBisect: O(log ε⁻¹) where ε is desired precision
For complex CSG regions, cost scales with the number of constituent half-spaces.
§Examples
use crater::analysis::prelude::*;
use crater::csg::prelude::*;
use burn::prelude::*;
use burn::backend::ndarray::{NdArrayDevice, NdArray};
use burn::backend::Autodiff;
// Create a sphere region
let sphere = Field3D::<Autodiff<NdArray>>::sphere(1.0, NdArrayDevice::Cpu).into_isosurface(0.0);
let region = Region::HalfSpace(sphere, Side::Negative);
let algebra = Algebra::<Autodiff<NdArray>>::default();
// Create rays pointing at the sphere
let origins = Tensor::from_data([[-2.0, 0.0, 0.0], [0.0, -2.0, 0.0]], &NdArrayDevice::Cpu);
let directions = Tensor::from_data([[1.0, 0.0, 0.0], [0.0, 1.0, 0.0]], &NdArrayDevice::Cpu);
let rays = Rays::<Autodiff<NdArray>, 3>::new(origins, directions);
// Use analytical algorithm for primitive shapes
let algorithm = RayCastAlgorithm::<3>::Analytical;
// let result = algorithm.solve(rays, ®ion, &algebra);
// Check intersection distances (should be ≈ 1.0 for unit sphere)
// let distances = result.distances();use crater::analysis::prelude::*;
use crater::csg::prelude::*;
use burn::prelude::*;
use burn::backend::wgpu::WgpuDevice;
// Complex CSG region requiring numerical methods
let sphere = Field3D::<burn::backend::wgpu::Wgpu>::sphere(1.0, WgpuDevice::default()).into_isosurface(0.0);
let cube = Field3D::<burn::backend::wgpu::Wgpu>::sphere(1.5, WgpuDevice::default()).into_isosurface(0.0);
let sphere_region = Region::HalfSpace(sphere, Side::Negative);
let cube_region = Region::HalfSpace(cube, Side::Negative);
let intersection = Region::Intersection(Box::new(sphere_region), Box::new(cube_region));
// Use bracket-and-bisect for robust intersection with CSG
let algorithm = RayCastAlgorithm::<3>::BracketAndBisect {
lambda: 10.0,
d_lambda: 0.01,
max_bisection_iterations: 50,
};
// let result = algorithm.solve(rays, &intersection, &algebra);Trait Implementations§
Source§impl<const N: usize> Clone for RayCastAlgorithm<N>
impl<const N: usize> Clone for RayCastAlgorithm<N>
Source§fn clone(&self) -> RayCastAlgorithm<N>
fn clone(&self) -> RayCastAlgorithm<N>
1.0.0 · Source§fn clone_from(&mut self, source: &Self)
fn clone_from(&mut self, source: &Self)
source. Read moreSource§impl<const N: usize> Debug for RayCastAlgorithm<N>
impl<const N: usize> Debug for RayCastAlgorithm<N>
impl<const N: usize> Copy for RayCastAlgorithm<N>
Auto Trait Implementations§
impl<const N: usize> Freeze for RayCastAlgorithm<N>
impl<const N: usize> RefUnwindSafe for RayCastAlgorithm<N>
impl<const N: usize> Send for RayCastAlgorithm<N>
impl<const N: usize> Sync for RayCastAlgorithm<N>
impl<const N: usize> Unpin for RayCastAlgorithm<N>
impl<const N: usize> UnwindSafe for RayCastAlgorithm<N>
Blanket Implementations§
Source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
Source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Source§impl<T> CloneToUninit for Twhere
T: Clone,
impl<T> CloneToUninit for Twhere
T: Clone,
Source§impl<T> IntoEither for T
impl<T> IntoEither for T
Source§fn into_either(self, into_left: bool) -> Either<Self, Self>
fn into_either(self, into_left: bool) -> Either<Self, Self>
self into a Left variant of Either<Self, Self>
if into_left is true.
Converts self into a Right variant of Either<Self, Self>
otherwise. Read moreSource§fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
self into a Left variant of Either<Self, Self>
if into_left(&self) returns true.
Converts self into a Right variant of Either<Self, Self>
otherwise. Read more