estimation_analysis¶
This module contains the functionality for managing the inputs and outputs of an estimation.
Functions¶
|
Function to propagate system covariance through time. |
Function to propagate system covariance through time. |
|
Function to propagate system covariance through time and convert it to RSW frame. |
|
|
Function to propagate system formal errors through time. |
Function to propagate system formal errors through time. |
|
Function to propagate system formal errors through time and convert to RSW frame. |
|
Function for creating an |
|
|
No documentation found. |
- propagate_covariance(initial_covariance: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, n]'], state_transition_interface: tudatpy.kernel.dynamics.simulator.CombinedStateTransitionAndSensitivityMatrixInterface, output_times: collections.abc.Sequence[SupportsFloat]) dict[float, Annotated[numpy.typing.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[ astro.time_representation.Time ]) – 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[ astro.time_representation.Time, numpy.ndarray[numpy.float64[m, n]] ]
- propagate_covariance_split_output(initial_covariance: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, n]'], state_transition_interface: tudatpy.kernel.dynamics.simulator.CombinedStateTransitionAndSensitivityMatrixInterface, output_times: collections.abc.Sequence[SupportsFloat]) tuple[list[float], list[Annotated[numpy.typing.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. Compared to the propagate_covariance function, this function returns a tuple, which contains a list of output times and a list of the propagated covariance, corresponding to the output times.
- 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[ astro.time_representation.Time ]) – 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:
Tuple containing a list of output times, and a list of propagated covariances at each output time.
- Return type:
tuple[ list[astro.time_representation.Time], list[numpy.ndarray[numpy.float64[m, n]]] ]
- propagate_covariance_rsw_split_output(initial_covariance: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, n]'], estimator: tudatpy.kernel.estimation.estimation_analysis.Estimator, output_times: collections.abc.Sequence[SupportsFloat]) tuple[list[float], list[Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]']]]¶
Function to propagate system covariance through time and convert it to RSW frame.
The covariance of a given system is propagated through time and afterwards converted to RSW frame. 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[ astro.time_representation.Time ]) – 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:
Tuple containing a list of output times, and a list of propagated covariances in RSW frame at each output time.
- Return type:
tuple[ list[float], list[numpy.ndarray[numpy.float64[m, n]]] ]
- propagate_formal_errors(initial_covariance: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, n]'], state_transition_interface: tudatpy.kernel.dynamics.simulator.CombinedStateTransitionAndSensitivityMatrixInterface, output_times: collections.abc.Sequence[SupportsFloat]) dict[float, Annotated[numpy.typing.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[ astro.time_representation.Time ]) – 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[ astro.time_representation.Time, numpy.ndarray[numpy.float64[m, 1]] ]
- propagate_formal_errors_split_output(initial_covariance: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, n]'], state_transition_interface: tudatpy.kernel.dynamics.simulator.CombinedStateTransitionAndSensitivityMatrixInterface, output_times: collections.abc.Sequence[SupportsFloat]) tuple[list[float], list[Annotated[numpy.typing.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[ astro.time_representation.Time ]) – 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:
Tuple containing a list of output times, and a list of propagated formal errors at each output time.
- Return type:
tuple[ list[astro.time_representation.Time], list[numpy.ndarray[numpy.float64[m, n]]] ]
- propagate_formal_errors_rsw_split_output(initial_covariance: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, n]'], estimator: tudatpy.kernel.estimation.estimation_analysis.Estimator, output_times: collections.abc.Sequence[SupportsFloat]) tuple[list[float], list[Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']]]¶
Function to propagate system formal errors through time and convert to RSW frame.
Function to propagate the formal errors of a given system through time. Note that in practice the entire covariance matrix is propagated and converted to RSW frame, 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[ astro.time_representation.Time ]) – 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:
Tuple containing a list of output times, and a list of propagated formal errors in RSW frame at each output time.
- Return type:
tuple[ list[astro.time_representation.Time], list[numpy.ndarray[numpy.float64[m, n]]] ]
- estimation_convergence_checker(maximum_iterations: SupportsInt = 5, minimum_residual_change: SupportsFloat = 0.0, minimum_residual: SupportsFloat = 0.0, number_of_iterations_without_improvement: SupportsInt = 2) tudatpy.kernel.estimation.estimation_analysis.EstimationConvergenceChecker¶
Function for creating an
EstimationConvergenceCheckerobject.Function for creating an
EstimationConvergenceCheckerobject, 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
EstimationConvergenceCheckerclass, defining the convergence criteria for an estimation.- Return type:
- create_best_fit_to_ephemeris(bodies: tudatpy.kernel.dynamics.environment.SystemOfBodies, acceleration_models: collections.abc.Mapping[str, collections.abc.Mapping[str, collections.abc.Sequence[tudatpy.kernel.dynamics.propagation.AccelerationModel]]], observed_bodies: collections.abc.Sequence[str], central_bodies: collections.abc.Sequence[str], integrator_settings: tudatpy.kernel.dynamics.propagation_setup.integrator.IntegratorSettings, initial_time: tudatpy.kernel.astro.time_representation.Time, final_time: tudatpy.kernel.astro.time_representation.Time, data_point_interval: tudatpy.kernel.astro.time_representation.Time, additional_parameter_names: collections.abc.Sequence[tudatpy.kernel.dynamics.parameters_setup.EstimatableParameterSettings] = [], number_of_iterations: SupportsInt = 3, reintegrate_variational_equations: bool = True, results_print_frequency: SupportsFloat = 0.0) tudat::simulation_setup::EstimationOutput<double, tudat::Time>¶
No documentation found.
Classes¶
Class for defining all specific inputs to a covariance analysis. |
|
Class for defining all inputs to the estimation. |
|
Class collecting all outputs from the covariance analysis process. |
|
Class collecting all outputs from the iterative estimation process. |
|
Class for consolidating all estimation functionality. |
|
Class defining the convergence criteria for an estimation. |
- class CovarianceAnalysisInput¶
Class for defining all specific inputs to a covariance analysis.
- __init__(self: tudatpy.kernel.estimation.estimation_analysis.CovarianceAnalysisInput, observations_and_times: tudatpy.kernel.estimation.observations.ObservationCollection, inverse_apriori_covariance: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, n]'] = array([], shape=(0, 0), dtype=float64), consider_covariance: Annotated[numpy.typing.ArrayLike, 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_settingsare 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
CovarianceAnalysisInputclass, defining the data and other settings to be used for the covariance analysis.- Return type:
- define_covariance_settings(self: tudatpy.kernel.estimation.estimation_analysis.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: SupportsFloat = 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_equationsto 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
- set_constant_single_observable_and_link_end_vector_weight(self: tudatpy.kernel.estimation.estimation_analysis.CovarianceAnalysisInput, observable_type: tudatpy.kernel.estimation.observable_models_setup.model_settings.ObservableType, link_ends: collections.abc.Mapping[tudatpy.kernel.estimation.observable_models_setup.links.LinkEndType, tudatpy.kernel.estimation.observable_models_setup.links.LinkEndId], weight: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]']) None¶
No documentation found.
- set_constant_single_observable_and_link_end_weight(self: tudatpy.kernel.estimation.estimation_analysis.CovarianceAnalysisInput, observable_type: tudatpy.kernel.estimation.observable_models_setup.model_settings.ObservableType, link_ends: collections.abc.Mapping[tudatpy.kernel.estimation.observable_models_setup.links.LinkEndType, tudatpy.kernel.estimation.observable_models_setup.links.LinkEndId], weight: SupportsFloat) None¶
No documentation found.
- set_constant_single_observable_vector_weight(self: tudatpy.kernel.estimation.estimation_analysis.CovarianceAnalysisInput, observable_type: tudatpy.kernel.estimation.observable_models_setup.model_settings.ObservableType, weight: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]']) None¶
No documentation found.
- set_constant_single_observable_weight(self: tudatpy.kernel.estimation.estimation_analysis.CovarianceAnalysisInput, observable_type: tudatpy.kernel.estimation.observable_models_setup.model_settings.ObservableType, weight: SupportsFloat) None¶
No documentation found.
- set_constant_vector_weight_per_observable(self: tudatpy.kernel.estimation.estimation_analysis.CovarianceAnalysisInput, weight_per_observable: collections.abc.Mapping[tudatpy.kernel.estimation.observable_models_setup.model_settings.ObservableType, Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]']]) None¶
No documentation found.
- set_constant_weight(self: tudatpy.kernel.estimation.estimation_analysis.CovarianceAnalysisInput, weight: SupportsFloat) 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.estimation.estimation_analysis.CovarianceAnalysisInput, weight_per_observable: collections.abc.Mapping[tudatpy.kernel.estimation.observable_models_setup.model_settings.ObservableType, SupportsFloat]) 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
- set_total_single_observable_and_link_end_vector_weight(self: tudatpy.kernel.estimation.estimation_analysis.CovarianceAnalysisInput, observable_type: tudatpy.kernel.estimation.observable_models_setup.model_settings.ObservableType, link_ends: collections.abc.Mapping[tudatpy.kernel.estimation.observable_models_setup.links.LinkEndType, tudatpy.kernel.estimation.observable_models_setup.links.LinkEndId], weight_vector: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]']) None¶
No documentation found.
- set_weights_from_observation_collection(self: tudatpy.kernel.estimation.estimation_analysis.CovarianceAnalysisInput) 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.estimation.estimation_analysis.EstimationInput, observations_and_times: tudatpy.kernel.estimation.observations.ObservationCollection, inverse_apriori_covariance: typing.Annotated[numpy.typing.ArrayLike, numpy.float64, "[m, n]"] = array([], shape=(0, 0), dtype=float64), convergence_checker: tudatpy.kernel.estimation.estimation_analysis.EstimationConvergenceChecker = <tudatpy.kernel.estimation.estimation_analysis.EstimationConvergenceChecker object at 0x7a62217c3b70>, consider_covariance: typing.Annotated[numpy.typing.ArrayLike, numpy.float64, "[m, n]"] = array([], shape=(0, 0), dtype=float64), consider_parameters_deviations: typing.Annotated[numpy.typing.ArrayLike, 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
EstimationInputclass, defining the data and other settings to be used for the estimation.- Return type:
- define_estimation_settings(self: tudatpy.kernel.estimation.estimation_analysis.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: SupportsFloat = 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_iterationto 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¶
read-only
Contribution of the consider parameters to the consider covariance matrix, equal to \((\mathbf{P} \mathbf{H}^{T} \mathbf{W}) (\mathbf{H}_c \mathbf{C} \mathbf{H}^{T}_c) (\mathbf{P} \mathbf{H}^{T} \mathbf{W})^{T}\)
- Type:
numpy.ndarray[numpy.float64[m, n]]
- property consider_normalization_factors¶
read-only
Vector of normalization terms used for consider covariance and consider design matrix
- Type:
numpy.ndarray[numpy.float64[m, 1]]
- 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}\). Note: if the
CovarianceAnalysisInputincludes consider parameters, this matrix does not include their contribution.- 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¶
read-only
Normalized consider covariance matrix \(\tilde{\mathbf{P}}_c\), with entries \(\tilde{P}_{c,ij}=P_{c,ij}\nu_{i}\nu_{j}\), where \(\nu_{i},\nu_{j}\) are the normalization terms as given by the
normalization_termsattribute.- Type:
numpy.ndarray[numpy.float64[m, n]]
- 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¶
read-only
Matrix of normalized partial derivatives for consider parameters \(\tilde{\mathbf{H}_c}\).
- Type:
numpy.ndarray[numpy.float64[m, n]]
- property unnormalized_covariance_with_consider_parameters¶
read-only
Consider covariance matrix \(\mathbf{P}_c\), equal to the sum of the
covariancematrix and theconsider_covariance_contributionmatrix.- Type:
numpy.ndarray[numpy.float64[m, n]]
- property unnormalized_design_matrix_consider_parameters¶
read-only
Matrix of unnormalized partial derivatives for consider parameters \(\mathbf{H}_c=\frac{\partial\mathbf{h}}{\partial\mathbf{c}}\).
- Type:
numpy.ndarray[numpy.float64[m, n]]
- 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, for the iteration with the lowest rms 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:
- class Estimator¶
Class for consolidating all estimation functionality.
Class for consolidating all functionality required to perform an estimation.
- __init__(self: tudatpy.kernel.estimation.estimation_analysis.Estimator, bodies: tudatpy.kernel.dynamics.environment.SystemOfBodies, estimated_parameters: tudatpy.kernel.dynamics.parameters.EstimatableParameterSet, observation_settings: collections.abc.Sequence[tudatpy.kernel.estimation.observable_models_setup.model_settings.ObservationModelSettings], propagator_settings: tudatpy.kernel.dynamics.propagation_setup.propagator.PropagatorSettings, integrate_on_creation: bool = True) None¶
Class constructor.
Constructor through which the user can create instances of this class. Defines environment, propagation and integrations models, as well as a number of settings related to the estimatable parameters and observation settings.
Note
When using default settings, creating an object of this type automatically triggers the propagation
- Parameters:
bodies (
SystemOfBodies) – Object defining the physical environment, with all properties of artificial and natural bodies.estimated_parameters (
EstimatableParameterSet) – Object defining a consolidated set of estimatable parameters, linked to the environment and acceleration settings of the simulation.observation_settings (
ObservationModelSettings) – List of settings objects, each object defining the observation model settings for one combination of observable and link geometry that is to be simulated.integrator_settings (
IntegratorSettings) – Settings to create the numerical integrator that is to be used for the integration of the equations of motionpropagator_settings (
PropagatorSettings) – Settings to create the propagator that is to be used for the propagation of dynamicsintegrate_on_creation (bool, default = True) – Boolean defining whether the propagation should be performed immediately (default), or at a later time (when calling the
perform_estimation()member function.
- compute_covariance(self: tudatpy.kernel.estimation.estimation_analysis.Estimator, covariance_analysis_input: tudat::simulation_setup::CovarianceAnalysisInput<double, tudat::Time>) tudat::simulation_setup::CovarianceAnalysisOutput<double, tudat::Time>¶
Function to perform a covariance analysis for the given observations and parameters
Function to perform a covariance analysis for the given observations and parameters. The observations are provided through the
covariance_analysis_inputinput, as are the weights \(\mathbf{W}\) and inverse a priori covariance \((\mathbf{P}_{0})^{-1}\). Calling this function uses the environment and propagator settings provided to the constructor of this Estimator class to simulate the dynamics of any relevant bodies for the observations (and associated variational equations). The observations are then computed using the observation models created by the settings provided to the constructor of this Estimator class, as is the associated design matrix \(\mathbf{H}\). This function then produces the covariance \(\mathbf{P}\) (omitting the normalization used internally for numerical stability)\[\mathbf{P}=\left(\mathbf{H}^{T}\mathbf{W}\mathbf{H}+(\mathbf{P}_{0})^{-1}\right)^{-1}\]Note that, although the actual observations are formally not required for a covariance analysis, all additional data (e.g. observation time, type, link ends, etc.) are. And, as such, the
covariance_analysis_inputdoes require the full set of observations and associated information, for consistency purposes (e.g., same input asperform_estimationfunction) .- Parameters:
covariance_analysis_input (
CovarianceAnalysisInput) – Object consolidating all relevant settings for the covariance analysis This includes foremost the simulated observations, as well as a priori information about the estimatable parameters- Returns:
Object containing all outputs from the estimation process.
- Return type:
- perform_estimation(self: tudatpy.kernel.estimation.estimation_analysis.Estimator, estimation_input: tudat::simulation_setup::EstimationInput<double, tudat::Time>) tudat::simulation_setup::EstimationOutput<double, tudat::Time>¶
Function to trigger the parameter estimation.
Function to trigger the parameter estimation. Much of the process and requirements are similar to those described in the
compute_covariance()function. This function uses an iterative least-squares estimate process to fit the data (insideestimation_input) to the model defined by the inputs to theEstimatorconstructor.s- Parameters:
estimation_input (
EstimationInput) – Object consolidating all relevant settings for the estimation This includes foremost the simulated observations, as well as a priori information about the estimatable parameters and convergence criteria for the least squares estimation.- Returns:
Object containing all outputs from the estimation process.
- Return type:
- property observation_managers¶
read-only
Observation managers contained in the Estimator object. A single observation manager can simulate observations and calculate observation partials for all link ends involved in the given observable type.
- Type:
dict[
ObservableType,ObservationManager]
- property observation_simulators¶
read-only
Observation simulators contained in the Estimator object. A single observation simulator hosts the functionality for simulating a given observable over the defined link geometry.
- Type:
list[
ObservationSimulator]
- property state_transition_interface¶
read-only
State transition and sensitivity matrix interface, setting the variational equations/dynamics in the Estimator object.
- property variational_solver¶
read-only
Variational equations solver, which is used to manage and execute the numerical integration of equations of motion and variational equations/dynamics in the Estimator object.
- 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()function.