estimation#

This module contains the functionality for managing the inputs and outputs of an estimation.

Functions#

simulate_observations(simulation_settings, ...)

Function to simulate observations.

compute_target_angles_and_range(bodies, ...)

Function to compute the azimuth angle, elevation angle and range at a ground station.

propagate_covariance(initial_covariance, ...)

Function to propagate system covariance through time.

propagate_formal_errors(initial_covariance, ...)

Function to propagate system formal errors through time.

estimation_convergence_checker([...])

Function for creating an EstimationConvergenceChecker object.

simulate_observations(simulation_settings: List[tudat::simulation_setup::ObservationSimulationSettings<double>], observation_simulators: List[tudatpy.kernel.numerical_simulation.estimation.ObservationSimulator], bodies: tudatpy.kernel.numerical_simulation.environment.SystemOfBodies) tudat::observation_models::ObservationCollection<double, double, 0>#

Function to simulate observations.

Function to simulate observations from set observation simulators and observation simulator settings. Automatically iterates over all provided observation simulators, generating the full set of simulated observations.

Parameters:
  • observation_to_simulate (List[ ObservationSimulationSettings ]) – List of settings objects, each object providing the observation time settings for simulating one type of observable and link end set.

  • observation_simulators (List[ ObservationSimulator ]) – List of ObservationSimulator objects, each object hosting the functionality for simulating one type of observable and link end set.

  • bodies (SystemOfBodies) – Object consolidating all bodies and environment models, including ground station models, that constitute the physical environment.

Returns:

Object collecting all products of the observation simulation.

Return type:

ObservationCollection

compute_target_angles_and_range(bodies: tudatpy.kernel.numerical_simulation.environment.SystemOfBodies, station_id: Tuple[str, str], target_body: str, observation_times: List[float], is_station_transmitting: bool) Dict[float, numpy.ndarray[numpy.float64[m, 1]]]#

Function to compute the azimuth angle, elevation angle and range at a ground station.

Function to compute the azimuth angle, elevation angle and range at a ground station. This functions is provided as a function of convenience, to prevent users having to manually define the relevant settings for this often-needed functionality. This function takes an observing station and a target body as input, and provides the observed angles and current range (without correction for aberrations, with correction for light time) as observed at that station

Parameters:
  • bodies (SystemOfBodies) – System of bodies that defines the full physical environment

  • station_id (tuple[ str, str]) – Identifier for the observing station, as a pair of strings: the body name and the station name.

  • target_body (str) – Name of body which is observed by ground station

  • observation_times (list[float]) – List of times at which the ground station observations are to be analyzed

  • is_station_transmitting (Bool) – Boolean defining whether the observation times define times at which the station is transmitting to, or receiving from, the ground station. This has an impact on the whether the light-time is computed forward or backward in time from the ground station to the target

Returns:

Dictionary with the required output. Key defines the observation time, the value is an array of size three containing entry 0 - elevation angle, entry 1 - azimuth angle, entry 2 - range

Return type:

dict[float,numpy.ndarray[numpy.float64[3, 1]]]

propagate_covariance(initial_covariance: numpy.ndarray[numpy.float64[m, n]], state_transition_interface: tudatpy.kernel.numerical_simulation.estimation.CombinedStateTransitionAndSensitivityMatrixInterface, output_times: List[float]) Dict[float, numpy.ndarray[numpy.float64[m, n]]]#

Function to propagate system covariance through time.

Function to propagate the covariance of a given system through time. The system dynamics and numerical settings of the propagation are prescribed by the state_transition_interface parameter.

Parameters:
  • initial_covariance (numpy.ndarray[numpy.float64[m, n]]) – System covariance matrix (symmetric and positive semi-definite) at initial time. Dimensions have to be consistent with estimatable parameters in the system (specified by state_transition_interface)

  • state_transition_interface (CombinedStateTransitionAndSensitivityMatrixInterface) – Interface to the variational equations of the system dynamics, handling the propagation of the covariance matrix through time.

  • output_times (List[ float ]) – Times at which the propagated covariance matrix shall be reported. Note that this argument has no impact on the integration time-steps of the covariance propagation, which always adheres to the integrator settings that the state_transition_interface links to. Output times which do not coincide with integration time steps are calculated via interpolation.

Returns:

Dictionary reporting the propagated covariances at each output time.

Return type:

Dict[ float, numpy.ndarray[numpy.float64[m, n]] ]

propagate_formal_errors(initial_covariance: numpy.ndarray[numpy.float64[m, n]], state_transition_interface: tudatpy.kernel.numerical_simulation.estimation.CombinedStateTransitionAndSensitivityMatrixInterface, output_times: List[float]) Dict[float, numpy.ndarray[numpy.float64[m, 1]]]#

Function to propagate system formal errors through time.

Function to propagate the formal errors of a given system through time. Note that in practice the entire covariance matrix is propagated, but only the formal errors (variances) are reported at the output times. The system dynamics and numerical settings of the propagation are prescribed by the state_transition_interface parameter.

Parameters:
  • initial_covariance (numpy.ndarray[numpy.float64[m, n]]) – System covariance matrix (symmetric and positive semi-definite) at initial time. Dimensions have to be consistent with estimatable parameters in the system (specified by state_transition_interface)

  • state_transition_interface (CombinedStateTransitionAndSensitivityMatrixInterface) – Interface to the variational equations of the system dynamics, handling the propagation of the covariance matrix through time.

  • output_times (List[ float ]) – Times at which the propagated covariance matrix shall be reported. Note that this argument has no impact on the integration time-steps of the covariance propagation, which always adheres to the integrator settings that the state_transition_interface links to. Output times which do not coincide with integration time steps are calculated via interpolation.

Returns:

Dictionary reporting the propagated formal errors at each output time.

Return type:

Dict[ float, numpy.ndarray[numpy.float64[m, 1]] ]

estimation_convergence_checker(maximum_iterations: int = 5, minimum_residual_change: float = 0.0, minimum_residual: float = 0.0, number_of_iterations_without_improvement: int = 2) tudatpy.kernel.numerical_simulation.estimation.EstimationConvergenceChecker#

Function for creating an EstimationConvergenceChecker object.

Function for creating an EstimationConvergenceChecker object, which is required for defining the convergence criteria of an estimation.

Parameters:
  • maximum_iterations (int, default = 5) – Maximum number of allowed iterations for estimation.

  • minimum_residual_change (float, default = 0.0) – Minimum required change in residual between two iterations.

  • minimum_residual (float, default = 0.0) – Minimum value of observation residual below which estimation is converged.

  • number_of_iterations_without_improvement (int, default = 2) – Number of iterations without reduction of residual.

Returns:

Instance of the EstimationConvergenceChecker class, defining the convergence criteria for an estimation.

Return type:

EstimationConvergenceChecker

Classes#

EstimatableParameterSet

Class containing a consolidated set of estimatable parameters.

ObservationViabilityCalculator

Template class for observation viability calculators.

ObservationSimulator

Class hosting the functionality for simulating observations.

ObservationCollection

Class collecting all observations and associated data for use in an estimation.

SingleObservationSet

Class collecting a single set of observations and associated data, of a given observable type, link ends, and ancilliary data.

CombinedStateTransitionAndSensitivityMatrixInterface

Class establishing an interface with the simulation's State Transition and Sensitivity Matrices.

EstimationConvergenceChecker

Class defining the convergence criteria for an estimation.

CovarianceAnalysisInput

Class for defining all specific inputs to a covariance analysis.

EstimationInput

Class for defining all inputs to the estimation.

CovarianceAnalysisOutput

Class collecting all outputs from the covariance analysis process.

EstimationOutput

Class collecting all outputs from the iterative estimation process.

class EstimatableParameterSet#

Class containing a consolidated set of estimatable parameters.

Class containing a consolidated set of estimatable parameters, linked to the environment and acceleration settings of the simulation. The user typically creates instances of this class via the create_parameters_to_estimate() factory function.

indices_for_parameter_type(self: tudatpy.kernel.numerical_simulation.estimation.EstimatableParameterSet, parameter_type: Tuple[tudat::estimatable_parameters::EstimatebleParametersEnum, Tuple[str, str]]) List[Tuple[int, int]]#

Function to retrieve the indices of a given type of parameter.

Function to retrieve the index of all parameters of a given type from the parameter set. This function can be very useful, since the order of parameters within the parameter set does not necessarily correspond to the order in which the elements were added to the set.

Parameters:

parameter_type (Tuple[ EstimatableParameterTypes, Tuple[str, str] ]) – help

Returns:

help

Return type:

List[ Tuple[int, int] ]

property constraints_size#

read-only

Total size of linear constraint that is to be applied during estimation.

Type:

int

property initial_multi_arc_states_size#

read-only

Amount of initial state parameters in the set, which are treated in a multi-arc fashion.

Type:

int

property initial_single_arc_states_size#

read-only

Amount of initial state parameters in the set, which are treated in a single-arc fashion.

Type:

int

property initial_states_size#

read-only

Amount of initial state parameters contained in the set.

Type:

int

property parameter_set_size#

read-only

Size of the parameter set, i.e. amount of estimatable parameters contained in the set.

Type:

int

property parameter_vector#

Vector containing the parameter values of all parameters in the set.

Type:

numpy.ndarray[numpy.float64[m, 1]]

class ObservationViabilityCalculator#

Template class for observation viability calculators.

Template class for classes which conducts viability calculations on simulated observations. Instances of the applicable ObservationViabilityCalculators are automatically created from the given ObservationSimulationSettings objects during the simulation of observations (simulate_observations()). The user typically does not interact directly with this class.

is_observation_viable(self: tudatpy.kernel.numerical_simulation.estimation.ObservationViabilityCalculator, link_end_states: List[numpy.ndarray[numpy.float64[6, 1]]], link_end_times: List[float]) bool#

Function to check whether an observation is viable.

Function to check whether an observation is viable. The calculation is performed based on the given times and link end states. Note, that this function is called automatically during the simulation of observations. Direct calls to this function are generally not required.

Parameters:
  • link_end_states (List[ numpy.ndarray[numpy.float64[6, 1]] ]) – Vector of states of the link ends involved in the observation.

  • link_end_times (List[float]) – Vector of times at the link ends involved in the observation.

Returns:

True if observation is viable, false if not.

Return type:

bool

class ObservationSimulator#

Class hosting the functionality for simulating observations.

Class hosting the functionality for simulating a given observable over a defined link geometry. Instances of this class are automatically created from the given ObservationSettings objects upon instantiation of the Estimator class.

class ObservationCollection#

Class collecting all observations and associated data for use in an estimation.

Class containing the full set of observations and associated data, typically for input into the estimation. When using simulated data, this class is instantiated via a call to the simulate_observations() function. More information is provided on the user guide

No documentation found.

Function to get all observation sets for a given observable type and link definition.

Parameters:
  • observable_type (ObservableType) – Observable type of which observations are to be simulated.

  • link_ends (LinkDefinition) – Link ends for which observations are to be simulated.

Returns:

List of observation sets for given observable type and link definition.

Return type:

list[ SingleObservationSet ]

read-only

Vector containing concatenated indices identifying the link ends. Each set of link ends is assigned a unique integer identifier (for a given instance of this class). The definition of a given integer identifier with the link ends is given by this class’ link_definition_ids() function. See user guide for details on storage order of the present vector.

Type:

numpy.ndarray[ int ]

property concatenated_observations#

read-only

Vector containing concatenated observable values. See user guide for details on storage order

Type:

numpy.ndarray[numpy.float64[m, 1]]

property concatenated_times#

read-only

Vector containing concatenated observation times. See user guide for details on storage order

Type:

numpy.ndarray[numpy.float64[m, 1]]

read-only

Dictionaty mapping a link end integer identifier to the specific link ends

Type:

dict[ int, dict[ LinkEndType, LinkEndId ] ]

No documentation found.

No documentation found.

property observable_type_start_index_and_size#

read-only

Dictionary defining per obervable type (dict key), the index in the full observation vector (concatenated_observations()) where the given observable type starts, and the number of subsequent entries in this vector containing a value of an observable of this type

Type:

dict[ ObservableType, [ int, int ] ]

property observation_set_start_index_and_size#

read-only

The nested dictionary/list returned by this property mirrors the structure of the sorted_observation_sets() property of this class. The present function provides the start index and size of the observables in the full observation vector that come from the correspoding SingleObservationSet in the sorted_observation_sets() Consequently, the present property returns a nested dictionary defining per obervable type, link end identifier, and SingleObservationSet index (for the given observable type and link end identifier), where the observables in the given SingleObservationSet starts, and the number of subsequent entries in this vector containing data from it.

Type:

dict[ ObservableType, dict[ int, list[ int, int ] ] ]

property observation_vector_size#

read-only

Length of the total vector of observations

Type:

int

property sorted_observation_sets#

read-only

The nested dictionary/list contains the list of SingleObservationSet objects, in the same method as they are stored internally in the present class. Specifics on the storage order are given in the user guide

Type:

dict[ ObservableType, dict[ int, list[ SingleObservationSet ] ] ]

class SingleObservationSet#

Class collecting a single set of observations and associated data, of a given observable type, link ends, and ancilliary data.

property ancilliary_settings#

read-only

Ancilliary settings all stored observations

Type:

ObservationAncilliarySimulationSettings

property concatenated_observations#

read-only

Concatenated vector of all stored observations

Type:

numpy.ndarray[numpy.float64[m, 1]]

read-only

Definition of the link ends for which the object stores observations

Type:

LinkDefinition

property list_of_observations#

read-only

List of separate stored observations. Each entry of this list is a vector containing a single observation. In cases where the observation is single-valued (range, Doppler), the vector is size 1, but for multi-valued observations such as angular position, each vector in the list will have size >1

Type:

list[ numpy.ndarray[numpy.float64[m, 1]] ]

property observable_type#

read-only

Type of observable for which the object stores observations

Type:

ObservableType

property observation_times#

read-only

Reference time for each of the observations in list_of_observations

Type:

list[ float]

property observations_history#

read-only

Dictionary of observations sorted by time. Created by making a dictionaty with observation_times as keys and list_of_observations as values

Type:

dict[ float, numpy.ndarray[numpy.float64[m, 1]] ]

read-only

Reference link end for all stored observations

Type:

LinkEndType

class CombinedStateTransitionAndSensitivityMatrixInterface#

Class establishing an interface with the simulation’s State Transition and Sensitivity Matrices.

Class establishing an interface to the State Transition and Sensitivity Matrices. Instances of this class are instantiated automatically upon creation of Estimator objects, using the simulation information in the observation, propagation and integration settings that the Estimator instance is linked to.

full_state_transition_sensitivity_at_epoch(self: tudatpy.kernel.numerical_simulation.estimation.CombinedStateTransitionAndSensitivityMatrixInterface, time: float, add_central_body_dependency: bool = True, arc_defining_bodies: List[str] = []) numpy.ndarray[numpy.float64[m, n]]#
Parameters:

time (float) – Time at which full concatenated state transition and sensitivity matrix are to be retrieved.

Returns:

Full concatenated state transition and sensitivity matrix at a given time.

Return type:

numpy.ndarray[numpy.float64[m, n]]

state_transition_sensitivity_at_epoch(self: tudatpy.kernel.numerical_simulation.estimation.CombinedStateTransitionAndSensitivityMatrixInterface, time: float, add_central_body_dependency: bool = True, arc_defining_bodies: List[str] = []) numpy.ndarray[numpy.float64[m, n]]#

Function to get the concatenated state transition and sensitivity matrix at a given time.

Function to get the concatenated state transition and sensitivity matrix at a given time. Entries corresponding to parameters which are not active at the current arc are omitted.

Parameters:

time (float) – Time at which concatenated state transition and sensitivity matrix are to be retrieved.

Returns:

Concatenated state transition and sensitivity matrix at a given time.

Return type:

numpy.ndarray[numpy.float64[m, n]]

property full_parameter_size#

read-only

Full amount of parameters w.r.t. which partials have been set up via State Transition and Sensitivity Matrices.

Type:

int

property sensitivity_size#

read-only

Number of columns in the sensitivity matrix.

Type:

int

property state_transition_size#

read-only

Size of the (square) state transition matrix.

Type:

int

class EstimationConvergenceChecker#

Class defining the convergence criteria for an estimation.

Class defining the convergence criteria for an estimation. The user typically creates instances of this class via the estimation_convergence_checker() factory function.

class CovarianceAnalysisInput#

Class for defining all specific inputs to a covariance analysis.

__init__(self: tudatpy.kernel.numerical_simulation.estimation.CovarianceAnalysisInput, observations_and_times: tudatpy.kernel.numerical_simulation.estimation.ObservationCollection, inverse_apriori_covariance: numpy.ndarray[numpy.float64[m, n]] = array([], shape=(0, 0), dtype=float64), consider_covariance: numpy.ndarray[numpy.float64[m, n]] = array([], shape=(0, 0), dtype=float64)) None#

Class constructor.

Constructor through which the user can create instances of this class. Note that the weight are all initiated as 1.0, and the default settings of define_covariance_settings are used.

Parameters:
  • observations_and_times (ObservationCollection) – Total data structure of observations and associated times/link ends/type/etc.

  • inverse_apriori_covariance (numpy.ndarray[numpy.float64[m, n]], default = [ ]) – A priori covariance matrix (unnormalized) of estimated parameters. This should be either a size 0x0 matrix (no a priori information), or a square matrix with the same size as the number of parameters that are considered

Returns:

Instance of the CovarianceAnalysisInput class, defining the data and other settings to be used for the covariance analysis.

Return type:

CovarianceAnalysisInput

define_covariance_settings(self: tudatpy.kernel.numerical_simulation.estimation.CovarianceAnalysisInput, reintegrate_equations_on_first_iteration: bool = True, reintegrate_variational_equations: bool = True, save_design_matrix: bool = True, print_output_to_terminal: bool = True, limit_condition_number_for_warning: float = 100000000.0) None#

Function to define specific settings for covariance analysis process

Function to define specific settings for covariance analysis process

Parameters:
  • reintegrate_equations (bool, default = True) – Boolean denoting whether the dynamics and variational equations are to be reintegrated or if existing values are to be used to perform first iteration.

  • reintegrate_variational_equations (bool, default = True) – Boolean denoting whether the variational equations are to be reintegrated during estimation (if this is set to False, and reintegrate_equations to true, only the dynamics are re-integrated)

  • save_design_matrix (bool, default = True) – Boolean denoting whether to save the partials matrix (also called design matrix) \(\mathbf{H}\) in the output. Setting this to false makes the \(\mathbf{H}\) matrix unavailable to the user, with the advantage of lower RAM usage.

  • print_output_to_terminal (bool, default = True) – Boolean denoting whether to print covariance-analysis-specific output to the terminal when running the estimation.

Returns:

Function modifies the object in-place.

Return type:

None

No documentation found.

No documentation found.

set_constant_single_observable_vector_weight(self: tudatpy.kernel.numerical_simulation.estimation.CovarianceAnalysisInput, observable_type: tudat::observation_models::ObservableType, weight: numpy.ndarray[numpy.float64[m, 1]]) None#

No documentation found.

set_constant_single_observable_weight(self: tudatpy.kernel.numerical_simulation.estimation.CovarianceAnalysisInput, observable_type: tudat::observation_models::ObservableType, weight: float) None#

No documentation found.

set_constant_vector_weight_per_observable(self: tudatpy.kernel.numerical_simulation.estimation.CovarianceAnalysisInput, weight_per_observable: Dict[tudat::observation_models::ObservableType, numpy.ndarray[numpy.float64[m, 1]]]) None#

No documentation found.

set_constant_weight(self: tudatpy.kernel.numerical_simulation.estimation.CovarianceAnalysisInput, weight: float) None#

Function to set a constant weight matrix for all observables.

Function to set a constant weight matrix for all observables. The weights are applied to all observations managed by the given PodInput object.

Parameters:

constant_weight (float) – Constant weight factor that is to be applied to all observations.

Returns:

Function modifies the object in-place.

Return type:

None

set_constant_weight_per_observable(self: tudatpy.kernel.numerical_simulation.estimation.CovarianceAnalysisInput, weight_per_observable: Dict[tudat::observation_models::ObservableType, float]) None#

Function to set a constant weight matrix for a given type of observable.

Function to set a constant weight matrix for a given type of observable. The weights are applied to all observations of the observable type specified by the weight_per_observable parameter.

Parameters:

constant_weight (Dict[ ObservableType, float ]) – Constant weight factor that is to be applied to all observations.

Returns:

Function modifies the object in-place.

Return type:

None

No documentation found.

property weight_matrix_diagonal#

read-only

Complete diagonal of the weights matrix that is to be used

Type:

numpy.ndarray[numpy.float64[n, 1]]

class EstimationInput#

Class for defining all inputs to the estimation.

__init__(self: tudatpy.kernel.numerical_simulation.estimation.EstimationInput, observations_and_times: tudatpy.kernel.numerical_simulation.estimation.ObservationCollection, inverse_apriori_covariance: numpy.ndarray[numpy.float64[m, n]] = array([], shape=(0, 0), dtype=float64), convergence_checker: tudatpy.kernel.numerical_simulation.estimation.EstimationConvergenceChecker = <tudatpy.kernel.numerical_simulation.estimation.EstimationConvergenceChecker object at 0x7f10b921eab0>, consider_covariance: numpy.ndarray[numpy.float64[m, n]] = array([], shape=(0, 0), dtype=float64), consider_parameters_deviations: numpy.ndarray[numpy.float64[m, 1]] = array([], dtype=float64), apply_final_parameter_correction: bool = True) None#

Class constructor.

Constructor through which the user can create instances of this class.

Parameters:
  • observations_and_times (ObservationCollection) – Total data structure of observations and associated times/link ends/type/etc.

  • inverse_apriori_covariance (numpy.ndarray[numpy.float64[m, n]], default = [ ]) – A priori covariance matrix (unnormalized) of estimated parameters. This should be either a size 0x0 matrix (no a priori information), or a square matrix with the same size as the number of parameters that are considered

  • convergence_checker (EstimationConvergenceChecker, default = estimation_convergence_checker()) – Object defining when the estimation is converged.

Returns:

Instance of the EstimationInput class, defining the data and other settings to be used for the estimation.

Return type:

EstimationInput

define_estimation_settings(self: tudatpy.kernel.numerical_simulation.estimation.EstimationInput, reintegrate_equations_on_first_iteration: bool = True, reintegrate_variational_equations: bool = True, save_design_matrix: bool = True, print_output_to_terminal: bool = True, save_residuals_and_parameters_per_iteration: bool = True, save_state_history_per_iteration: bool = False, limit_condition_number_for_warning: float = 100000000.0, condition_number_warning_each_iteration: bool = True) None#

Function to define specific settings for the estimation process

Function to define specific settings for covariance analysis process

Parameters:
  • reintegrate_equations_on_first_iteration (bool, default = True) – Boolean denoting whether the dynamics and variational equations are to be reintegrated or if existing values are to be used to perform first iteration.

  • reintegrate_variational_equations (bool, default = True) – Boolean denoting whether the variational equations are to be reintegrated during estimation (if this is set to False, and reintegrate_equations_on_first_iteration to true, only the dynamics are re-integrated)

  • save_design_matrix (bool, default = True) – Boolean denoting whether to save the partials matrix (also called design matrix) \(\mathbf{H}\) in the output. Setting this to false makes the \(\mathbf{H}\) matrix unavailable to the user, with the advantage of lower RAM usage.

  • print_output_to_terminal (bool, default = True) – Boolean denoting whether to print covariance-analysis-specific output to the terminal when running the estimation.

  • save_residuals_and_parameters_per_iteration (bool, default = True) – Boolean denoting whether the residuals and parameters from the each iteration are to be saved.

  • save_state_history_per_iteration (bool, default = False) – Boolean denoting whether the state history and dependent variables are to be saved on each iteration.

Returns:

Function modifies the object in-place.

Return type:

None

class CovarianceAnalysisOutput#

Class collecting all outputs from the covariance analysis process.

property consider_covariance_contribution#

No documentation found.

property consider_normalization_factors#

No documentation found.

property correlations#

read-only

Correlation matrix of the estimation result. Entry \(i,j\) is equal to \(P_{i,j}/(\sigma_{i}\sigma_{j})\)

Type:

numpy.ndarray[numpy.float64[m, m]]

property covariance#

read-only

(Unnormalized) estimation covariance matrix \(\mathbf{P}\).

Type:

numpy.ndarray[numpy.float64[m, m]]

property design_matrix#

read-only

Matrix of unnormalized partial derivatives \(\mathbf{H}=\frac{\partial\mathbf{h}}{\partial\mathbf{p}}\).

Type:

numpy.ndarray[numpy.float64[m, n]]

property formal_errors#

read-only

Formal error vector \(\boldsymbol{\sigma}\) of the estimation result (e.g. square root of diagonal entries of covariance)s

Type:

numpy.ndarray[numpy.float64[m, 1]]s

property inverse_covariance#

read-only

(Unnormalized) inverse estimation covariance matrix \(\mathbf{P}^{-1}\).

Type:

numpy.ndarray[numpy.float64[m, m]]

property inverse_normalized_covariance#

read-only

Normalized inverse estimation covariance matrix \(\mathbf{\tilde{P}}^{-1}\).

Type:

numpy.ndarray[numpy.float64[m, m]]

property normalization_terms#

read-only

Vector of normalization terms used for covariance and design matrix

Type:

numpy.ndarray[numpy.float64[m, 1]]

property normalized_covariance#

read-only

Normalized estimation covariance matrix \(\mathbf{\tilde{P}}\).

Type:

numpy.ndarray[numpy.float64[m, m]]

property normalized_covariance_with_consider_parameters#

No documentation found.

property normalized_design_matrix#

read-only

Matrix of normalized partial derivatives \(\tilde{\mathbf{H}}\).

Type:

numpy.ndarray[numpy.float64[m, n]]

property normalized_design_matrix_consider_parameters#

No documentation found.

property unnormalized_covariance_with_consider_parameters#

No documentation found.

property weighted_design_matrix#

read-only

Matrix of weighted partial derivatives, equal to \(\mathbf{W}^{1/2}{\mathbf{H}}\)

Type:

numpy.ndarray[numpy.float64[m, n]]

property weighted_normalized_design_matrix#

read-only

Matrix of weighted, normalized partial derivatives, equal to \(\mathbf{W}^{1/2}\tilde{\mathbf{H}}\)

Type:

numpy.ndarray[numpy.float64[m, n]]

class EstimationOutput#

Class collecting all outputs from the iterative estimation process.

property best_iteration#

No documentation found.

property final_parameters#

No documentation found.

property final_residuals#

read-only

Vector of post-fit observation residuals.

Type:

numpy.ndarray[numpy.float64[m, 1]]

property parameter_history#

read-only

Parameter vectors, concatenated per iteration into a matrix. Column 0 contains pre-estimation values. The \((i+1)^{th}\) column has the residuals from the \(i^{th}\) iteration.

Type:

numpy.ndarray[numpy.float64[m, n]]

property residual_history#

read-only

Residual vectors, concatenated per iteration into a matrix; the \(i^{th}\) column has the residuals from the \(i^{th}\) iteration.

Type:

numpy.ndarray[numpy.float64[m, n]]

property simulation_results_per_iteration#

read-only

List of complete numerical propagation results, with the \(i^{th}\) entry of thee list thee results of the \(i^{th}\) propagation

Type:

list[SimulationResults]