RayCastAlgorithm

Enum RayCastAlgorithm 

Source
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 Analytical for:

    • Simple primitive shapes (spheres, planes, cylinders)
    • Maximum accuracy and performance requirements
    • Real-time applications with tight computational budgets
  • Use Newton for:

    • Smooth, well-conditioned implicit surfaces
    • When autodifferentiation is available
    • Applications requiring fast convergence with few iterations
  • Use BracketAndBisect for:

    • 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 phase
  • d_lambda - Step size during bracketing search
  • max_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

§lambda: f32

Maximum distance to extend rays during bracketing

§d_lambda: f32

Step size for bracketing phase

§max_bisection_iterations: usize

Maximum iterations for bisection refinement

§

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 cases
  • step_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

Fields

§max_iterations: usize

Maximum number of Newton iterations

§nudge_distance: f32

Distance to perturb points with zero gradients

§step_size: f32

Step size relaxation factor (0.0 to 1.0)

Implementations§

Source§

impl<const N: usize> RayCastAlgorithm<N>

Source

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 surface
  • region - Implicit surface region to intersect with
  • algebra - 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, &region, &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>

Source§

fn clone(&self) -> RayCastAlgorithm<N>

Returns a duplicate of the value. Read more
1.0.0 · Source§

fn clone_from(&mut self, source: &Self)

Performs copy-assignment from source. Read more
Source§

impl<const N: usize> Debug for RayCastAlgorithm<N>

Source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result

Formats the value using the given formatter. Read more
Source§

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> Any for T
where T: 'static + ?Sized,

Source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
Source§

impl<T> Borrow<T> for T
where T: ?Sized,

Source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
Source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

Source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
Source§

impl<T> CloneToUninit for T
where T: Clone,

Source§

unsafe fn clone_to_uninit(&self, dest: *mut u8)

🔬This is a nightly-only experimental API. (clone_to_uninit)
Performs copy-assignment from self to dest. Read more
Source§

impl<T> Downcast<T> for T

Source§

fn downcast(&self) -> &T

Source§

impl<T> DynClone for T
where T: Clone,

Source§

fn __clone_box(&self, _: Private) -> *mut ()

Source§

impl<T> From<T> for T

Source§

fn from(t: T) -> T

Returns the argument unchanged.

Source§

impl<T, U> Into<U> for T
where U: From<T>,

Source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

Source§

impl<T> IntoEither for T

Source§

fn into_either(self, into_left: bool) -> Either<Self, Self>

Converts 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 more
Source§

fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
where F: FnOnce(&Self) -> bool,

Converts 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
Source§

impl<T> Pointable for T

Source§

const ALIGN: usize

The alignment of pointer.
Source§

type Init = T

The type for initializers.
Source§

unsafe fn init(init: <T as Pointable>::Init) -> usize

Initializes a with the given initializer. Read more
Source§

unsafe fn deref<'a>(ptr: usize) -> &'a T

Dereferences the given pointer. Read more
Source§

unsafe fn deref_mut<'a>(ptr: usize) -> &'a mut T

Mutably dereferences the given pointer. Read more
Source§

unsafe fn drop(ptr: usize)

Drops the object pointed to by the given pointer. Read more
Source§

impl<T> ToOwned for T
where T: Clone,

Source§

type Owned = T

The resulting type after obtaining ownership.
Source§

fn to_owned(&self) -> T

Creates owned data from borrowed data, usually by cloning. Read more
Source§

fn clone_into(&self, target: &mut T)

Uses borrowed data to replace owned data, usually by cloning. Read more
Source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

Source§

type Error = Infallible

The type returned in the event of a conversion error.
Source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
Source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

Source§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
Source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.
Source§

impl<T> Upcast<T> for T

Source§

fn upcast(&self) -> Option<&T>

Source§

impl<V, T> VZip<V> for T
where V: MultiLane<T>,

Source§

fn vzip(self) -> V

Source§

impl<T> WasmNotSend for T
where T: Send,

Source§

impl<T> WasmNotSendSync for T

Source§

impl<T> WasmNotSync for T
where T: Sync,