3. LKF#

class LKF(trajectory: Trajectory, initial_covariance: CovarianceMatrix, measurement_data: list, dataset_names: list)#

Bases: FilterOD

Implements a Linearized Kalman Filter (LKF) for spacecraft orbit determination.

Applies the Linearized Kalman Filter (LKF) method to estimate a spacecraft’s state using sequential measurements while propagating uncertainty through a dynamic model. The LKF is particularly useful when the state evolution and measurement models are nonlinear, but can be linearized using a State Transition Matrix (STM).

The filter follows a predict-update cycle, where:

  • The prediction step propagates the state and covariance forward in time.

  • The update step refines the state using measurement updates via the Kalman gain.

The LKF is recursive, meaning it updates the estimate as new measurements become available, making it computationally efficient compared to batch methods.

Parameters:
  • trajectory (Trajectory) –

    The spacecraft trajectory object containing:

    • The nominal trajectory state history.

    • State Transition Matrices (STMs) for covariance propagation.

  • initial_covariance (CovarianceMatrix) – The initial covariance matrix representing uncertainty in the estimated state.

  • measurement_data (list) –

    A list of measurement observations, where each entry includes:

    • Measurement time.

    • Residuals between observed and predicted values.

    • Measurement covariance matrix.

    • Measurement Jacobian (sensitivity matrix).

  • dataset_names (list) – Names of the datasets used in the estimation process.

See also

scarabaeus.FilterOD

Base class for orbit determination filters.

scarabaeus.SolutionOD

Stores estimated state history and covariance evolution.

scarabaeus.CovarianceMatrix

Defines measurement noise covariance.

scarabaeus.ProcessNoiseCovariance

Represents process noise in the system dynamics.

Notes

The LKF assumes Gaussian noise models for both system dynamics and measurement errors.

Process noise can be included to account for unmodeled dynamics in state propagation.

The filter can operate across multiple trajectory legs, adjusting the state size at event epochs.

Unlike batch least-squares, the Kalman filter updates estimates sequentially, allowing real-time implementation.

References

Tapley, B. D., Schutz, B. E., & Born, G. H. (2004). Statistical orbit determination. Elsevier.

Methods

compute([meas_corr, printOutput, ...])

Performs the state estimation using the Linearized Kalman Filter.

initialize_process_noise(...)

Initializes the process noise for the orbit determination.

map_covariance_definition_node(...[, ...])

Adjusts the covariance matrix to match a new state vector definition.

map_deviation_definition_node(...[, ...])

Adjusts the deviation vector to align with a new state vector definition.

map_stm_definition_node(counter_events, old_STM)

Adjusts the State Transition Matrix (STM) to match a new state definition.

observability_test([print_obs, threshold, ...])

Performs an observability analysis on the system by computing the information matrix.

unwrap_residuals(postfit_residuals, ...)

Organizes postfit residuals into a dictionary categorized by dataset.

validate_sequence_batch()

Validates consistency in state vector definitions across sequential estimation legs.

compute(meas_corr: list = None, printOutput: bool = True, underweighting_factor: float = 1.0) SolutionOD#

Performs the state estimation using the Linearized Kalman Filter.

Iterates through the measurement epochs, performing prediction and update steps. The predicted state is corrected based on new measurements using the Kalman gain.

Parameters:
  • meas_corr (list, optional) – List of correction factors for measurement covariance. Defaults to None.

  • printOutput (bool, optional) – If True, prints the progress of the filter. Defaults to True.

  • underweighting_factor (float, optional) – A scaling factor to modify the measurement covariance (used for robustness). Defaults to 1.0.

  • Returns

    solSolutionOD

    An object containing the estimated state deviation history and covariance evolution.

initialize_process_noise(continuous_time_process_covariance: CovarianceMatrix)#

Initializes the process noise for the orbit determination.

Sets up the continuous-time process noise covariance and computes its discrete-time equivalent for use in the filter updates.

Parameters:

continuous_time_process_covariance (CovarianceMatrix) – The process noise covariance matrix defined in continuous time.

Return type:

None

map_covariance_definition_node(counter_events: int, old_covariance: ndarray, flag_augment: bool = True) ndarray#

Adjusts the covariance matrix to match a new state vector definition.

Parameters:
  • counter_events (int) – Index of the current sequence event for state transition.

  • old_covariance (numpy.ndarray) – Previous covariance matrix to be updated.

  • flag_augment (bool, optional) – If True, expands the covariance matrix; otherwise, reduces it. Defaults to True.

Returns:

updated_covar – Updated covariance matrix reflecting the new state definition.

Return type:

numpy.ndarray

map_deviation_definition_node(counter_events: int, old_deviation: ndarray, flag_augment: bool = True) ndarray#

Adjusts the deviation vector to align with a new state vector definition.

Parameters:
  • counter_events (int) – Index of the current sequence event for state transition.

  • old_deviation (numpy.ndarray) – Previous state deviation vector.

  • flag_augment (bool, optional) – If True, expands the deviation vector; otherwise, reduces it. Defaults to True.

Returns:

updated_dev – Updated deviation vector reflecting the new state definition.

Return type:

numpy.ndarray

map_stm_definition_node(counter_events: int, old_STM: ndarray, flag_augment: bool = True) ndarray#

Adjusts the State Transition Matrix (STM) to match a new state definition.

Parameters:
  • counter_events (int) – Index of the current sequence event for state transition.

  • old_STM (numpy.ndarray) – Previous STM to be updated.

  • flag_augment (bool, optional) – If True, expands the STM; otherwise, reduces it. Defaults to True.

Returns:

updated_stm – Updated State Transition Matrix reflecting the new state definition.

Return type:

numpy.ndarray

observability_test(print_obs=False, threshold=0.1, meas_corr=None)#

Performs an observability analysis on the system by computing the information matrix.

Checks whether the system’s observability matrix has full rank and determines the positive definiteness of the information matrix.

Parameters:
  • print_obs (bool, optional) – If True, prints details of the observability analysis.

  • threshold (float, optional) – Minimum eigenvalue threshold for checking positive definiteness.

  • meas_corr (list, optional) – Measurement covariance correction factors.

Returns:

info_matrix, obs_matrix – A tuple with the following values corresponding to their respective indices:

  • [0] = info_matrixnumpy.ndarray

    The computed information matrix.

  • [1] = obs_matrixnumpy.ndarray

    The observability matrix constructed from measurement partials.

Return type:

tuple[numpy.ndarray, numpy.ndarray]

unwrap_residuals(postfit_residuals: list, dataset_names: list[str]) dict#

Organizes postfit residuals into a dictionary categorized by dataset.

Parameters:
  • postfit_residuals (list) – List of residuals and associated metadata.

  • dataset_names (list) – List of dataset names corresponding to different measurements.

Returns:

unwrapped_residuals – Dictionary mapping dataset names to corresponding residuals and sigmas.

Return type:

dict

validate_sequence_batch()#

Validates consistency in state vector definitions across sequential estimation legs.

Ensures that each state vector in a sequence contains all necessary elements from previous legs.

Return type:

None