giant.calibration.estimators.IterativeNonlinearLSTSQ

class giant.calibration.estimators.IterativeNonlinearLSTSQ(model=None, weighted_estimation=False, max_iter=20, residual_atol=1e-10, residual_rtol=1e-10, state_atol=1e-10, state_rtol=1e-10, measurements=None, camera_frame_directions=None, measurement_covariance=None, a_priori_state_covariance=None, temperatures=None)[source]

This concrete estimator implements iterative non-linear least squares for estimating an updated camera model.

Iterative non-linear least squares estimation is done by estimating updates to the “state” vector (in this case the camera model parameters being updated) iteratively. At each step, the system is linearized about the current estimate of the state and the additive update is estimated. This iteration is repeated until convergence (or divergence) based on the pre/post update residuals and the update vector itself.

The state vector that is being estimated by this class is controlled by the CameraModel.estimation_parameters attribute of the provided camera model. This class does not actually use the CameraModel.estimation_parameters attribute since it is handled by the CameraModel.compute_jacobian() and CameraModel.apply_update() methods of the provided camera model internally, but it is mentioned here to show how to control what exactly is being estimated.

Because this class linearizes about the current estimate of the state, it requires an initial guess for the camera model that is “close enough” to the actual model to ensure convergence. Defining “close enough” in any broad sense is impossible, but based on experience, using the manufacturer defined specs for focal length/pixel pitch and assuming no distortion is generally “close enough” even for cameras with heavy distortion (star identification may require a better initial model than this anyway).

As this class converges the state estimate, it updates the supplied camera model in place, therefore, if you wish to keep a copy of the original camera model, you should manually create a copy of it before calling the estimate() method on this class.

In the estimate() method, convergence is checked on both the sum of squares of the residuals and the update vector for the state. That is convergence is reached when either of

\begin{gather*} \left\|\mathbf{r}_{pre}^T\mathbf{r}_{pre} - \mathbf{r}_{post}^T\mathbf{r}_{post}\right\| \le(a_r+r_r\mathbf{r}_{pre}^T\mathbf{r}_{pre}) \\ \text{all}\left[\left\|\mathbf{u}\right\|\le(a_s+r_s\mathbf{s}_{pre})\right] \end{gather*}

is True. Here \(\mathbf{r}_{pre}\) is the nx1 vector of residuals before the update is applied, \(\mathbf{r}_{post}\) is the nx1 vector of residuals after the update is applied, \(a_r\) is the residual_atol absolute residual tolerance, \(r_r\) is the residual_rtol relative residual tolerance, \(\mathbf{u}\) is the update vector, \(\text{all}\) indicates that the contained expression is True for all elements, \(a_s\) is the state_atol absolute tolerance for the state vector, \(r_s\) is the state_rtol relative tolerance for the state vector, and \(\mathbf{s}_{pre}\) is the state vector before the update is applied. Divergence is only checked on the sum of squares of the residuals, that is, divergence is occurring when

\[\mathbf{r}_{pre}^T\mathbf{r}_{pre} < \mathbf{r}_{post}^T\mathbf{r}_{post}\]

where all is as defined as before. If a case is diverging then a warning will be printed, the iteration will cease, and successful will be set to False.

Typically this class is not used by the user, and instead it is used internally by the Calibration class which handles data preparation for you. If you wish to use this externally from the Calibration class you must first set

according to their documentation. Once those have been set, you can perform the estimation using estimate() which will iterate until convergence (or divergence). If the fit successfully converges, successful will be set to True and attributes postfit_covariance and postfit_residuals will both return numpy arrays instead of None. If you wish to use the same instance of this class to do another estimation you should call reset() before setting the new data to ensure that data is not mixed between estimation runs and all flags are set correctly.

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 and successful will be set to False.

  • 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

  • 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 to True. The length of the state vector can be determined by len(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 for CameraModel.estimation_parameters

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

property model: CameraModel | None

The camera model that is being estimated.

Typically this should be a subclass of CameraModel.

property successful: bool

A boolean flag indicating whether the fit was successful or not.

If the fit was successful this should return True, and False 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 to False

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 in measurement_covariance during the estimation process. If set to False, then no measurement weights will be considered.

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 to True then this property will contain the measurement covariance matrix as a square, full rank, numpy array or the measurement variance as a scalar float. If weighted_estimation is set to False then this property may be None 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) where cov is the covariance matrix, v is the specified scalar variance, and I(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 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 to True and if CameraModel.use_a_priori is set to True, otherwise it is ignored. If both are set to True then this should be set to a square, full rank, lxl numpy array where l=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 from CameraModel.get_state_labels().

Raises:

ValueError – If the shape of the input matrix is not appropriate for the size of the state 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 is measurements[:, i] <-> np.concatenate(camera_frame_directions, axis=-1)[:, i])

This must always be set before a call to estimate().

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 the measurements array. (That is np.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] == [[], [], []] where i 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 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 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 to estimation_parameters and if weighted_estimation is True will return the state covariance matrix taking into account the measurement covariance matrix. If weighted_estimation is False, then this will return the post-fit state covariance matrix assuming no measurement weighting (that is a measurement covariance matrix of the identity matrix). If estimate() has not been called yet or the fit was unsuccessful then this will return None

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(). If estimate() has not been called yet or the fit was unsuccessful then this will return None.

Summary of Methods

compute_residuals

This method computes the observed minus computed residuals for the current model (or an input model).

estimate

Estimates an updated camera model that better transforms the camera frame directions into pixel locations to minimize the residuals between the observed and the predicted star locations.

reset

This method resets all of the data attributes to their default values to prepare for another estimation.