7. SRIF#
- class SRIF(trajectory: Trajectory, initial_covariance: CovarianceMatrix, measurement_data: list, dataset_names: list)#
Bases:
FilterOD
Implements the Square Root Information Filter (SRIF) for spacecraft orbit determination.
The SRIF is a variant of the Kalman filter that operates on the square root of the information matrix rather than the covariance matrix. This formulation improves numerical stability and allows more accurate handling of poorly conditioned systems.
The filter follows a predict-update cycle, where:
The prediction step propagates the state and information matrix forward in time.
The update step refines the state using measurement updates via an orthogonal transformation (Cholesky decomposition).
The advantage of the SRIF over the standard Kalman filter is its numerical stability, especially when dealing with high-dimensional state vectors and long-duration orbit determination problems.
- 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 SRIF is particularly useful for handling ill-conditioned sequential estimation problems.
The filter maintains positive definiteness of covariance matrices, avoiding numerical issues that arise in the standard Kalman filter.
References
Tapley, B. D., Schutz, B. E., & Born, G. H. (2004). Statistical orbit determination. Elsevier.
Methods
compute
([meas_corr, print_output, ...])Performs the state estimation using the Square Root Information Filter.
compute_batch
([meas_corr, ...])Computes the batch least-squares solution using the Square Root Information Filter (SRIF).
initialize_process_noise
(...[, state_definition])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.
smoother
()Computes the smoothed state deviation and covariance using backward recursion.
unwrap_residuals
(postfit_residuals, ...)Organizes postfit residuals into a dictionary categorized by dataset.
Validates consistency in state vector definitions across sequential estimation legs.
- compute(meas_corr: list = None, print_output: bool = True, underweighting_factor: float = 1.0) SolutionOD #
Performs the state estimation using the Square Root Information Filter.
Iterates through the measurement epochs, performing prediction and update steps. The predicted state is corrected based on new measurements using an orthogonal transformation (QR decomposition).
- Parameters:
meas_corr (
list
, optional) – List of correction factors for measurement covariance. Defaults toNone
.print_output (
bool
, optional) – IfTrue
, prints the progress of the filter. Defaults toTrue
.underweighting_factor (
float
, optional) – A scaling factor to modify the measurement covariance (used for robustness). Defaults to1.0
.
- Returns:
sol – An object containing the estimated state deviation history and covariance evolution.
- Return type:
- compute_batch(meas_corr: list = None, map_back_sequence: bool = True, print_output: bool = True) SolutionOD #
Computes the batch least-squares solution using the Square Root Information Filter (SRIF).
Performs batch estimation of the spacecraft state by processing all available measurements in a least-squares sense. The formulation relies on the square root of the information matrix, leveraging Cholesky decomposition to ensure numerical stability.
In sequence mode, maps the final leg deviation and covariance back through previous trajectory legs, ensuring continuity across estimation segments.
- Parameters:
meas_corr (
list
, optional) – List of correction factors applied to the measurement covariance. Defaults toNone
.map_back_sequence (
bool
, optional) – If T`rue, maps the final leg deviation and covariance back to previous legs. Defaults toTrue
.print_output (
bool
, optional) – IfTrue
, prints progress information during execution. Defaults toTrue
.
- Returns:
sol – An object containing the computed state deviation history, covariance evolution, and postfit residuals.
- Return type:
- initialize_process_noise(continuous_time_process_covariance: CovarianceMatrix, state_definition=None) None #
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:
- 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) – IfTrue
, expands the covariance matrix; otherwise, reduces it. Defaults toTrue
.
- Returns:
updated_covar – Updated covariance matrix reflecting the new state definition.
- Return type:
- 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) – IfTrue
, expands the deviation vector; otherwise, reduces it. Defaults toTrue
.
- Returns:
updated_dev – Updated deviation vector reflecting the new state definition.
- Return type:
- 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) – IfTrue
, expands the STM; otherwise, reduces it. Defaults toTrue
.
- Returns:
updated_stm – Updated State Transition Matrix reflecting the new state definition.
- Return type:
- 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:
- Returns:
info_matrix, obs_matrix – A tuple with the following values corresponding to their respective indices:
[0]
= info_matrixnumpy.ndarrayThe computed information matrix.
[1]
= obs_matrixnumpy.ndarrayThe observability matrix constructed from measurement partials.
- Return type:
tuple[numpy.ndarray
,numpy.ndarray]
- smoother()#
Computes the smoothed state deviation and covariance using backward recursion.
The smoother refines the estimated state by incorporating future measurements through a backward pass over the measurement sequence. It provides the best accurate trajectory reconstruction.
The smoothing process applies orthogonal transformations to extract the corrected state deviation and covariance at each epoch.
- Returns:
A tuple with the following values corresponding to their respective indices:
[0]
= smoothed_state_devslist[numpy.ndarray]The refined state deviations at each measurement epoch.
[1]
= smoothed_covars: list[numpy.ndarray]The corresponding smoothed covariance matrices.
- Return type:
smoothed_state_devs
,smooth_covars = tuple[list[numpy.ndarray]
,list[numpy.ndarray]]