giant.calibration.estimators.LMAEstimator¶
- class giant.calibration.estimators.LMAEstimator(model=None, weighted_estimation=False, max_iter=20, residual_atol=1e-10, residual_rtol=1e-10, state_atol=1e-10, state_rtol=1e-10, max_divergence_steps=5, measurements=None, camera_frame_directions=None, measurement_covariance=None, a_priori_state_covariance=None, temperatures=None)[source]¶
Bases:
IterativeNonlinearLSTSQ
This implements a Levenberg-Marquardt Algorithm estimator, which is analogous to a damped iterative non-linear least squares.
This class is nearly exactly the same as the
IterativeNonlinearLSTSQ
except that it adds damping to the update step of the iterative non-linear least squares algorithm and allows a few diverging steps in a row where the damping parameter is updated before failing. The number of diverging steps that are allowed is controlled by themax_divergence_steps
setting. This represents only difference from theIterativeNonlinearLSTSQ
interface from the user’s perspective.In general, this algorithm will result in the same answer as the
IterativeNonlinearLSTSQ
algorithm but at a slower convergence rate. In a few cases however, this estimator can be more robust to initial guess errors, achieving convergence when the standard iterative nonlinear least squares diverges. Therefore, it is likely best to start with theIterativeNonlinearLSTSQ
class an only switch to this if you experience convergence issues.The implementation of the LMA in this class is inspired by https://link.springer.com/article/10.1007/s40295-016-0091-3
- Parameters:
model (CameraModel | None) – The camera model instance to be estimated set with an initial guess of the state.
weighted_estimation (bool) – A boolean flag specifying whether to do weighted estimation.
True
indicates that the measurement weights (and a priori state covariance if applicable) should be used in the estimation.max_iter (int) – The maximum number of iteration steps to attempt to reach convergence. If convergence has not been reached after attempting
max_iter
steps, a warning will be raised that the model has not converged andsuccessful
will be set toFalse
.residual_atol (float) – The absolute convergence tolerance criteria for the sum of squares of the residuals
residual_rtol (float) – The relative convergence tolerance criteria for the sum of squares of the residuals
state_atol (float) – The absolute convergence tolerance criteria for the elements of the state vector
state_rtol (float) – The relative convergence tolerance criteria for the elements of the state vector
max_divergence_steps (int) – The maximum number of steps in a row that can diverge before breaking iteration
measurements (Sequence | ndarray | None) – A 2xn numpy array of measurement pixel locations to be fit to
camera_frame_directions (List[ndarray | List[List]] | None) – A length m list of 3xj numpy arrays or empty 3x1 list of empty lists where m is the number of unique images the data comes from (and is the same length as
temperatures
) and j is the number of measurements from each image. A list of empty lists indicates that no measurements were identified for the corresponding image.measurement_covariance (Sequence | ndarray | Real | None) – An optional nxn numpy array containing the covariance matrix for the ravelled measurement vector (in fortran order such that the ravelled measurement vector is [x1, y1, x2, y2, … xk, yk] where k=n//2)
a_priori_state_covariance (Sequence | ndarray | None) – An optional lxl numpy array containing the a priori covariance matrix for the a priori estimate of the state, where l is the number of parameters in the state vector. This is used only if
CameraModel.use_a_priori
is set toTrue
. The length of the state vector can be determined bylen(
CameraModel.state_vector
)
temperatures (List[Real] | None) – A length m list of floats containing the camera temperature at the time of each corresponding image. These may be used by the
CameraModel
to perform temperature dependent estimation of parameters like the focal length, depending on what is set forCameraModel.estimation_parameters
- property a_priori_state_covariance: Sequence | ndarray | None¶
A square numpy array containing the covariance matrix for the a priori estimate of the state vector.
This is only considered if
weighted_estimation
is set toTrue
and ifCameraModel.use_a_priori
is set toTrue
, otherwise it is ignored. If both are set toTrue
then this should be set to a square, full rank, lxl numpy array wherel=len(model.state_vector)
containing the covariance matrix for the a priori state vector. The order of the parameters in the state vector can be determined fromCameraModel.get_state_labels()
.- Raises:
ValueError – If the shape of the input matrix is not appropriate for the size of the state vector
- property camera_frame_directions: List[ndarray | List[List]] | None¶
A length m list of unit vectors in the camera frame as numpy arrays for m images corresponding to the
measurements
attribute.Each element of this list corresponds to a unique image that is being considered for estimation and the subsequent element in the
temperatures
list. Each column of this concatenated array will correspond to the same column of themeasurements
array. (That isnp.concatenate(camera_frame_directions, axis=-1)[:, i] <-> measurements[:, i]
).Any images for which no stars were identified (due to any number of reasons) will have a list of empty arrays in the corresponding element of this list (that is
camera_frame_directions[i] == [[], [], []]
wherei
is an image with no measurements identified). These will be automatically dropped by numpy’s concatenate, but are included to notify the which temperatures/misalignments to use.This must always be set before a call to
estimate()
.
- property measurement_covariance: Sequence | ndarray | Real | None¶
A square numpy array containing the covariance matrix for the measurements or a scalar containing the variance for all of the measurements.
If
weighted_estimation
is set toTrue
then this property will contain the measurement covariance matrix as a square, full rank, numpy array or the measurement variance as a scalar float. Ifweighted_estimation
is set toFalse
then this property may beNone
and will be ignored.If specified as a scalar, it is treated as the variance for each measurement (that is
cov = v*I(n,n)
wherecov
is the covariance matrix,v
is the specified scalar variance, andI(n,n)
is a nxn identity matrix) in a memory efficient way.- Raises:
ValueError – When attempting to set an array that does not have the proper shape for the
measurements
vector
- property measurements: Sequence | ndarray | None¶
A nx2 numpy array of the observed pixel locations for stars across all images
Each column of this array corresponds to the same column of the
camera_frame_directions
concatenated down the last axis. (That ismeasurements[:, i] <-> np.concatenate(camera_frame_directions, axis=-1)[:, i]
)This must always be set before a call to
estimate()
.
- property model: CameraModel | None¶
The camera model that is being estimated.
Typically this should be a subclass of
CameraModel
.
- property postfit_covariance: Sequence | ndarray | None¶
The post-fit state covariance matrix, taking into account the measurement covariance matrix (if applicable).
This returns the post-fit state covariance matrix after a call to
estimate()
. The covariance matrix will be in the order according toestimation_parameters
and ifweighted_estimation
isTrue
will return the state covariance matrix taking into account the measurement covariance matrix. Ifweighted_estimation
isFalse
, then this will return the post-fit state covariance matrix assuming no measurement weighting (that is a measurement covariance matrix of the identity matrix). Ifestimate()
has not been called yet or the fit was unsuccessful then this will returnNone
- property postfit_residuals: Sequence | ndarray | None¶
The post-fit observed-computed measurement residuals as a 2xn numpy array.
This returns the post-fit observed minus computed measurement residuals after a call to
estimate()
. Ifestimate()
has not been called yet or the fit was unsuccessful then this will returnNone
.
- property successful: bool¶
A boolean flag indicating whether the fit was successful or not.
If the fit was successful this should return
True
, andFalse
if otherwise. A fit is defined as successful if convergence criteria were reached before the maximum number of iterations. Divergence and non-convergence are both considered an unsuccessful fit resulting in this being set toFalse
- property temperatures: List[Real] | None¶
A length m list of temperatures of the camera for each image being considered in estimation.
Each element of this list corresponds to a unique image that is being considered for estimation and the subsequent element in the
camera_frame_directions
list.This must always be set before a call to
estimate()
(although sometimes it may be a list of all zeros if temperature data is not available for the camera).
- property weighted_estimation: bool¶
A boolean flag specifying whether to do weighted estimation.
If set to
True
, the estimator will use the provided measurement weights inmeasurement_covariance
during the estimation process. If set toFalse
, then no measurement weights will be considered.
- max_iter: int¶
The maximum number of iteration steps to attempt for convergence
- residual_atol: float¶
The absolute tolerance for the sum of square of the residuals to indicate convergence
- residual_rtol: float¶
The relative tolerance for the sum of square of the residuals to indicate convergence
- state_atol: float¶
The absolute tolerance for the state vector to indicate convergence
- state_rtol: float¶
The relative tolerance for the state vector to indicate convergence
- max_divergence_steps: int¶
The maximum number of steps in a row that can diverge before breaking iteration
Summary of Methods
This method computes the observed minus computed residuals for the current model (or an input model). |
|
This method estimates the postfit residuals based on the model, weight matrix, lma coefficient, etc. |
|
This method resets all of the data attributes to their default values to prepare for another estimation. |