Bplane#

class Bplane(epoch: EpochArray | float, bplane_name: str, bplane_spice_id: int | str, sc_name: str, sc_spice_id: int | str, target_name: str, target_spice_id: int | str, fk_file: str | Path, v_hat: bool = True, new_bplane: bool = True, scenario_kernel_dir: str | Path | None = None)#

Bases: object

B-plane and B-vector definition and targeting.

The B-plane is a plane passing through the centre of a target body, perpendicular to the incoming hyperbolic asymptote of the approach trajectory. The B-vector is the vector from the target-body centre to the point where the asymptote pierces that plane, and its components B·T and B·R are the primary targeting parameters used in planetary-flyby and encounter design.

B-plane geometry

The B-plane is defined at a reference epoch by three right-handed orthonormal vectors:

\[\begin{split}\hat{S} &= \frac{\mathbf{v}_\infty}{\|\mathbf{v}_\infty\|} \quad\text{(incoming asymptote direction)}\\ \hat{T} &= \frac{\hat{k} \times \hat{S}}{\|\hat{k} \times \hat{S}\|} \quad\text{(in B-plane, perpendicular to S)}\\ \hat{R} &= \hat{S} \times \hat{T}\end{split}\]

where \(\hat{k}\) is the ecliptic north pole and \(\mathbf{v}_\infty\) is the hyperbolic excess velocity vector.

The B-vector is the perpendicular from the target-body centre to the incoming asymptote:

\[\mathbf{B} = B_T\,\hat{T} + B_R\,\hat{R}, \qquad |\mathbf{B}| = \sqrt{B_T^2 + B_R^2}\]

For a hyperbolic flyby with gravitational parameter \(\mu\), the periapsis distance is

\[r_p = \frac{\mu}{v_\infty^2} \left(\sqrt{1 + \frac{v_\infty^2\,|\mathbf{B}|^2}{\mu^2}} - 1\right)\]

Targeting

target() minimizes \(\|\mathbf{b}_{target} - \mathbf{b}_{ref}\|\) by solving for the \(\Delta v\) correction at a specified maneuver epoch, where \(\mathbf{b}_{target} = [B\cdot T_{target},\;B\cdot R_{target},\;\mathrm{LTOF}_{target}]\).

Parameters:
  • epoch (EpochArray | datetime.datetime | float) – Reference epoch for the B-plane frame. Floats are interpreted as J2000 ephemeris seconds (ET); datetime objects are converted automatically.

  • bplane_name (str) – Name of the B-plane SPICE frame (e.g. 'BPLANE_EUROPA').

  • bplane_spice_id (int or str) – SPICE numeric frame ID for the B-plane frame.

  • sc_name (str) – Human-readable spacecraft name (used for labelling only).

  • sc_spice_id (int or str) – SPICE body ID (name string or integer) for the spacecraft.

  • target_name (str) – Human-readable target-body name; must match a key in scarabaeus.Constants that exposes a mu entry.

  • target_spice_id (int or str) – SPICE body ID (name string or integer) for the target body.

  • fk_file (str or pathlib.Path) – Path to the frames-kernel (.tf) file. If new_bplane is True the file is created here; otherwise it is loaded from this path.

  • v_hat (bool, optional) – When True (default) the secondary axis of the B-plane frame is aligned with the +Y direction; False flips it to −Y.

  • new_bplane (bool, optional) – True (default) creates and furnishes a new FK file. False loads an existing kernel from fk_file.

  • scenario_kernel_dir (str or pathlib.Path, optional) – Directory into which the FK file is written. Defaults to <cwd>/data/Kernels.

bt_ref, br_ref

Reference B-plane targeting parameters B·T and B·R [km].

Type:

ArrayWUnits

ltof_ref#

Linearised time-of-flight from the reference epoch to close approach [s].

Type:

ArrayWUnits

time_close_approach_ref#

Estimated close-approach epoch expressed as ET seconds.

Type:

ArrayWUnits

bpar_ref#

Concatenated reference vector [B·T, B·R, ltof].

Type:

ArrayWUnits

See also

scarabaeus.Trajectory

Trajectory container used as input to compute the B-plane.

scarabaeus.EpochArray

Epoch representation for the B-plane reference epoch.

Attributes:
bplane_name

SPICE frame name of the B-plane.

bplane_spice_id

SPICE numeric frame ID of the B-plane.

br_ref

Reference B·R component of the B-vector [km].

bt_ref

Reference B·T component of the B-vector [km].

epoch

Reference epoch as an EpochArray.

fk_file

Absolute path to the B-plane frames kernel file.

ltof_ref

Reference linearised time-of-flight [s].

new_bplane

Whether a new FK file was created by this instance.

rot_to_bplane

Rotation matrix (DCM) from J2000 to the B-plane frame.

rot_to_j2000

Rotation matrix (DCM) from the B-plane frame to J2000.

rotmat_to_bplane

Alias for rot_to_bplane.

rotmat_to_j2000

Alias for rot_to_j2000.

sc_name

Human-readable spacecraft name.

sc_spice_id

Resolved SPICE body identifier for the spacecraft.

target_name

Human-readable target-body name.

target_spice_id

Resolved SPICE body identifier for the target.

time_close_approach_ref

Estimated close-approach epoch [ET seconds].

v_hat

Secondary-axis orientation flag (True → +Y, False → −Y).

Methods

compute_hyperbolic_parameters(sc_rel_state_j2000)

Compute classical hyperbolic B-plane parameters from a Cartesian state.

compute_jacobian(sc_rel_state_j2000)

Finite-difference Jacobian of the B-plane parameters w.r.t.

compute_parameters(sc_rel_state_j2000)

Compute the linear B-plane parameters from a spacecraft state.

linear_targeting(sc_rel_state_j2000[, err, ...])

Compute a trajectory-correction manoeuvre (TCM) via linear targeting.

nonlinear_targeting(sc_rel_state_j2000[, err])

Compute a TCM via nonlinear optimisation (Nelder-Mead).

piecewise_targeting(trajectory, tcm_epochs)

Compute TCMs at multiple epochs along a trajectory and (optionally) generate a corrected SPK ephemeris file.

propagate_covariance_to_close_approach(...)

Propagate a state covariance to the close-approach epoch.

target_covariance(sc_cov, sc_rel_state_j2000)

Map a spacecraft state covariance to B-plane parameter covariance.

target_spacecraft_covariance(sc_targ_cov, ...)

Map a joint spacecraft-and-target covariance to B-plane parameter covariance.

compute_hyperbolic_parameters(sc_rel_state_j2000: ArrayWUnits) ArrayWUnits#

Compute classical hyperbolic B-plane parameters from a Cartesian state.

This routine evaluates the full conic geometry of the approach trajectory and derives the B-vector components B·T and B·R from the asymptote direction and the orbit’s semi-latus rectum. For elliptic trajectories (eccentricity < 1) the B-plane is undefined and the corresponding output entries are set to NaN.

Parameters:

sc_rel_state_j2000 (ArrayWUnits, shape (6,) or (N, 6)) – Relative state(s) in J2000 with units convertible to km / km·s⁻¹.

Returns:

[B·T [km], B·R [km], ltof [s]].

Return type:

ArrayWUnits, shape (3,) or (N, 3)

Notes

The B-plane convention follows Kizner (1961):

  • S – unit vector along the incoming asymptote.

  • T – perpendicular to S lying in the target equatorial plane (constructed as T = [-Sy, Sx, 0] / |[-Sy, Sx, 0]|).

  • R = S × T.

  • B = semi-minor axis × (S × ĥ).

Warning

The T-vector construction becomes degenerate when S is nearly parallel to the z-axis (i.e. the approach asymptote nearly coincides with the north pole direction). - The computation of beta = arccos(1/ecc) is only valid for hyperbolic or parabolic cases. - The construction of Tvec as [-Sy, Sx, 0] becomes degenerate when the projection of Svec onto the xy-plane is near zero. - If ecc == 1 exactly, the semimajor/seminor-axis expressions may become numerically delicate. - gm is retrieved from constants using

self.target_name.

compute_jacobian(sc_rel_state_j2000: ArrayWUnits) ArrayWUnits#

Finite-difference Jacobian of the B-plane parameters w.r.t. velocity.

Computes the 3×3 matrix d[B·T, B·R, ltof] / d[vx, vy, vz] via forward finite differences with perturbation size chosen as sqrt(eps_machine) * max(1, |v_i|).

Parameters:

sc_rel_state_j2000 (ArrayWUnits, shape (6,)) – Current relative state in J2000.

Returns:

jacobian – Jacobian matrix with units [km/(km/s), km/(km/s), s/(km/s)] i.e. [s, s, s²/km] per column.

Return type:

ArrayWUnits, shape (3, 3)

compute_parameters(sc_rel_state_j2000: ArrayWUnits) ArrayWUnits#

Compute the linear B-plane parameters from a spacecraft state.

Projects the relative state onto the B-plane and returns the intersection coordinates together with the linearised time-of-flight.

Parameters:

sc_rel_state_j2000 (ArrayWUnits, shape (6,)) – Relative Cartesian state [x, y, z, vx, vy, vz] of the spacecraft with respect to the target body, expressed in J2000.

Returns:

bpar[B·T, B·R, ltof] where B·T and B·R are in km and ltof is in seconds.

Return type:

ArrayWUnits, shape (3,)

linear_targeting(sc_rel_state_j2000: ArrayWUnits, err: float = 1e-08, max_iterations: int = 100) ArrayWUnits#

Compute a trajectory-correction manoeuvre (TCM) via linear targeting.

Iterates a Newton-style correction loop using the local Jacobian until the residual in B-plane parameter space drops below err.

Parameters:
  • sc_rel_state_j2000 (ArrayWUnits, shape (6,)) – Perturbed relative state in J2000 at the manoeuvre epoch.

  • err (float, optional) – Convergence threshold on the L2 norm of the B-plane residual.

  • max_iterations (int, optional) – Maximum number of Newton iterations before raising an error.

Returns:

DV_tot – Velocity increment vector [ΔVx, ΔVy, ΔVz] in km/s.

Return type:

ArrayWUnits, shape (3,)

Raises:

RuntimeError – If the iteration does not converge within max_iterations steps.

nonlinear_targeting(sc_rel_state_j2000: ArrayWUnits, err: float = 1e-10) ArrayWUnits#

Compute a TCM via nonlinear optimisation (Nelder-Mead).

Minimises the squared L2 norm of the B-plane residual ‖[B·T, B·R, ltof]_ref [B·T, B·R, ltof]‖² over the three velocity components at the manoeuvre epoch.

Parameters:
  • sc_rel_state_j2000 (ArrayWUnits, shape (6,)) – Perturbed relative state in J2000 at the manoeuvre epoch.

  • err (float, optional) – Tolerance passed to scipy.optimize.minimize.

Returns:

DV_tot – Velocity increment vector [ΔVx, ΔVy, ΔVz] in km/s.

Return type:

ArrayWUnits, shape (3,)

piecewise_targeting(trajectory: Trajectory, tcm_epochs: list, method: str = 'linear', err: float = 1e-08, max_iterations: int = 100, spk_output_path: str | Path | None = None) list[ArrayWUnits]#

Compute TCMs at multiple epochs along a trajectory and (optionally) generate a corrected SPK ephemeris file.

The algorithm proceeds segment-by-segment:

  1. Evaluate the current relative state at tcm_epoch[i].

  2. Solve for the ΔV that corrects the B-plane parameters to the reference values using the chosen method.

  3. Apply the ΔV to the state, propagate to the next TCM epoch (or to close approach for the last segment), and write the resulting arc to the SPK file.

Parameters:
  • trajectory (Trajectory) –

    A scarabaeus.Trajectory object that exposes at minimum:

    • get_state(epoch)ArrayWUnits relative state in J2000.

    • propagate(state, epoch_start, epoch_end) → updated Trajectory segment.

    • write_spk(filepath) → generate a SPICE SPK kernel.

  • tcm_epochs (list of EpochArray | datetime | float) – Ordered sequence of epochs at which TCMs are to be computed. Must be in chronological order.

  • method ({'linear', 'nonlinear'}, optional) – Targeting algorithm. 'linear' uses the Newton iteration (linear_targeting()); 'nonlinear' uses Nelder-Mead (nonlinear_targeting()).

  • err (float, optional) – Convergence tolerance forwarded to the chosen targeting method.

  • max_iterations (int, optional) – Maximum Newton iterations (linear method only).

  • spk_output_path (str or pathlib.Path, optional) – If provided, a SPICE SPK file covering the full corrected trajectory is written to this path after all TCMs have been applied.

Returns:

dv_list – ΔV vectors [ΔVx, ΔVy, ΔVz] [km/s] at each TCM epoch.

Return type:

list of ArrayWUnits

Raises:

ValueError – If method is not 'linear' or 'nonlinear'.

propagate_covariance_to_close_approach(sc_cov: ndarray, sc_rel_state_j2000: ArrayWUnits, process_noise: TypeAliasForwardRef('numpy.ndarray') | None = None) ArrayWUnits#

Propagate a state covariance to the close-approach epoch.

Note

Not yet implemented. A correct implementation requires a state transition matrix (STM) Φ from the manoeuvre epoch to close approach, so that:

P_ca = Φ @ P @ Φ.T + Q

where Q is the process-noise covariance. The mapped B-plane covariance is then obtained via _map_covariance() evaluated at the close-approach state. Implement once the Trajectory class exposes get_stm(epoch_start, epoch_end) -> ndarray.

Raises:

NotImplementedError – Always.

target_covariance(sc_cov: ndarray, sc_rel_state_j2000: ArrayWUnits) ArrayWUnits#

Map a spacecraft state covariance to B-plane parameter covariance.

Uses the analytical linearisation of the B-plane projection to propagate a 6×6 Cartesian covariance (in J2000) into the 3×3 space [B·T, B·R, ltof].

Assumes no cross-correlation between the B-plane position components and the linearised time-of-flight.

Parameters:
  • sc_cov (array-like, shape (6, 6)) – Full covariance matrix of the spacecraft state in J2000 (position–velocity).

  • sc_rel_state_j2000 (ArrayWUnits, shape (6,)) – Relative Cartesian state in J2000 at the close-approach epoch.

Returns:

P_bplane – B-plane parameter covariance with units [[km², km², km·s], [km², km², km·s], [km·s, km·s, s²]].

Return type:

ArrayWUnits, shape (3, 3)

target_spacecraft_covariance(sc_targ_cov: ndarray, sc_rel_state_j2000: ArrayWUnits) ArrayWUnits#

Map a joint spacecraft-and-target covariance to B-plane parameter covariance.

Accounts for uncertainty in both the spacecraft and the target body by first forming the relative covariance via the 6×12 subtraction matrix, then propagating as in target_covariance().

Parameters:
  • sc_targ_cov (array-like, shape (12, 12)) – Joint covariance matrix of the spacecraft and target states in J2000, ordered as [sc_pos, sc_vel, tgt_pos, tgt_vel].

  • sc_rel_state_j2000 (ArrayWUnits, shape (6,)) – Relative Cartesian state sc target in J2000 at close approach.

Returns:

P_bplane – B-plane parameter covariance (same units as target_covariance()).

Return type:

ArrayWUnits, shape (3, 3)

property bplane_name: str#

SPICE frame name of the B-plane.

property bplane_spice_id: int#

SPICE numeric frame ID of the B-plane.

property br_ref: ArrayWUnits#

Reference B·R component of the B-vector [km].

property bt_ref: ArrayWUnits#

Reference B·T component of the B-vector [km].

property epoch: EpochArray#

Reference epoch as an EpochArray.

property fk_file: Path#

Absolute path to the B-plane frames kernel file.

property new_bplane: bool#

Whether a new FK file was created by this instance.

property rot_to_bplane: ArrayWUnits#

Rotation matrix (DCM) from J2000 to the B-plane frame.

property rot_to_j2000: ArrayWUnits#

Rotation matrix (DCM) from the B-plane frame to J2000.

property rotmat_to_bplane: ArrayWUnits#

Alias for rot_to_bplane.

property rotmat_to_j2000: ArrayWUnits#

Alias for rot_to_j2000.

property sc_name: str#

Human-readable spacecraft name.

property sc_spice_id: int | str#

Resolved SPICE body identifier for the spacecraft.

property target_name: str#

Human-readable target-body name.

property target_spice_id: int | str#

Resolved SPICE body identifier for the target.

property v_hat: bool#

Secondary-axis orientation flag (True → +Y, False → −Y).