Single-Hit Poincare Sections

The singlehit module provides single-hit Poincare section detection for individual trajectories.

Backend Classes

class hiten.algorithms.poincare.singlehit.backend._SingleHitBackend[source]

Bases: _ReturnMapBackend

Concrete backend for single-hit Poincare section crossing search.

This class implements the generic surface-of-section crossing search for single-hit Poincare sections. It extends the abstract base class to provide a complete implementation using numerical integration and root finding.

The backend uses a two-stage approach: 1. Coarse integration to get near the section 2. Fine root finding to locate the exact crossing point

Notes

This backend is optimized for single-hit computations where only the first intersection with the section is needed. It uses efficient root finding to locate the exact crossing point after coarse integration.

The backend is stateless - all dynamic system data must be passed to the step_to_section method as arguments.

All time units are nondimensional unless otherwise specified.

run(seeds, *, dynsys, surface, dt=0.01, t_guess=None, forward=1, method='adaptive', order=8)[source]

Find the next crossing for every seed.

This method implements the core functionality of the single-hit backend, finding the first intersection of each seed trajectory with the Poincare section.

Parameters:
  • seeds (ndarray, shape (m, 6)) – Array of initial state vectors [x, y, z, vx, vy, vz] in nondimensional units.

  • dt (float, default=1e-2) – Integration time step (nondimensional units). Used for Runge-Kutta methods, ignored for adaptive methods.

  • t_guess (float, optional) – Initial guess for the crossing time. If None, uses a default value based on the orbital period.

  • dynsys (_DynamicalSystemProtocol)

  • surface (_SurfaceEvent)

  • forward (int)

  • method (Literal['fixed', 'symplectic', 'adaptive'])

  • order (int)

Returns:

  • points (ndarray, shape (k, 2)) – Array of 2D intersection points in the section plane. Only includes trajectories that successfully cross the section.

  • states (ndarray, shape (k, 6)) – Array of full state vectors at the intersection points. Shape matches points array.

Return type:

tuple[ndarray, ndarray]

Notes

This method processes each seed individually, finding the first intersection with the section. Trajectories that don’t cross the section are excluded from the results.

The method uses a two-stage approach: 1. Coarse integration to get near the section 2. Fine root finding to locate the exact crossing

The 2D projection uses the first two coordinates as a fallback projection method.

Functions

hiten.algorithms.poincare.singlehit.backend.find_crossing()[source]

Find a single crossing for a given state and surface.

Parameters:
  • dynsys (_DynamicalSystemProtocol) – The dynamical system providing the equations of motion.

  • state0 (array_like, shape (6,)) – Initial state vector [x, y, z, vx, vy, vz] in nondimensional units.

  • surface (_SurfaceEvent) – The Poincare section surface definition.

  • **kwargs – Additional keyword arguments passed to the backend constructor.

  • forward (int)

Returns:

Tuple of (points, states) arrays from the backend’s step_to_section method.

Return type:

tuple[numpy.ndarray, numpy.ndarray]

Notes

This is a convenience function that creates a single-hit backend and finds the crossing for a single state vector. It’s useful for simple crossing computations without needing to create a backend instance explicitly.

hiten.algorithms.poincare.singlehit.backend._plane_crossing_factory()[source]

Factory function for creating plane crossing functions.

Parameters:
  • coord (str) – Coordinate identifier for the plane (e.g., ‘x’, ‘y’, ‘z’).

  • value (float, default=0.0) – Plane offset value (nondimensional units).

  • direction ({1, -1, None}, optional) – Crossing direction filter.

Returns:

A function that finds crossings for the specified plane.

Return type:

callable

Notes

This factory function creates specialized crossing functions for specific coordinate planes. The returned function takes a dynamical system and initial state and returns the crossing time and state.

The returned function signature is:

_section_crossing(*, dynsys, x0, forward=1, **kwargs) -> (time, state)