PinholeModel.compute_jacobian

giant.camera_models.pinhole_model:

PinholeModel.compute_jacobian(unit_vectors_camera, temperature=0)[source]

Calculates the Jacobian matrix for each observation in unit_vectors_camera for each parameter to be estimated as defined in the estimation_parameters attribute.

This method works by first computing the partial derivatives for all camera parameters for each provided unit vector. It then concatenates the rows into the Jacobian matrix and removes any columns of parameters that are not specified in the estimation_parameters attribute and sorts the columns according to the order of the estimation_parameters attribute. The resulting Jacobian will be the appropriate size and in the order specified by estimation_parameters. There is one constraint that the misalignment (if included) must be last in estimation_parameters.

The unit_vectors_camera inputs should be formatted as a Sequence of 2d sequences. Each inner 2D sequence should be of shape \(3\times x\), where each row corresponds to a component of a unit vector in the camera frame. Each inner sequence should contain all observations from a single image, so that if there are \(m\) images being considered, then the outer sequence should be length \(m\). The value of \(x\) can change for each image. If you are estimating multiple misalignments (one for each image) then each misalignment will correspond to the order of the image observations in the outer sequence.

You can also set the use_a_priori to True to have this method append an identity matrix to the bottom of this Jacobian if you are solving for an update to your camera model, and not a new one entirely.

The optional temperature input specifies the temperature of the camera for use in estimating temperature dependence. The temperature input should either be a scalar value (float or int), or a list that is the same length as unit_vectors_camera, where each element of the list is the temperature of the camera at the time of each image represented by unit_vectors_camera. If the temperature input is a scalar, then it is assumed to be the temperature value for all of the images represented in unit_vectors_camera.

Parameters:
  • unit_vectors_camera (Sequence[Sequence[Sequence] | ndarray]) – The points/directions in the camera frame that the jacobian matrix is to be computed for. For multiple images, this should be a list of 2D unit vectors where each element of the list corresponds to a new image.

  • temperature (Sequence | ndarray | Real) – A single temperature for all images or a list of temperatures the same length of unit_vectors_camera containing the temperature of the camera at the time each image was captured

Returns:

The Jacobian matrix evaluated for each observation

Return type:

ndarray

Example:

>>> from giant.camera_models import PinholeModel
>>> model = PinholeModel(kx=300, ky=400, px=500, py=500, focal_length=10, a1=1e-5, a2=1e-6,
>>>                      misalignment = [[1e-12, -2e-14, 3e-10], [2e-15, 1e-13, 3e-10]],
>>>                      estimation_parameters = ['multiple misalignments'])
>>> model.compute_jacobian([[[0.5], [0], [1]], [[0.1, 0.2, 0.3], [-0.3, -0.4, -0.5], [4, 5, 6]]],
>>>                        temperature=[10, -20])
array([[  0.00000000e+00,  -3.75075000e+03,   0.00000000e+00,
          0.00000000e+00,   0.00000000e+00,   0.00000000e+00],
       [  4.00080000e+03,   2.98059600e-07,  -2.00040000e+03,
          0.00000000e+00,   0.00000000e+00,   0.00000000e+00],
       [  0.00000000e+00,   0.00000000e+00,   0.00000000e+00,
         -5.62612499e+00,  -3.00247537e+03,  -2.25045000e+02],
       [  0.00000000e+00,   0.00000000e+00,   0.00000000e+00,
          4.02330450e+03,   7.50150000e+00,  -1.00020000e+02],
       [  0.00000000e+00,   0.00000000e+00,   0.00000000e+00,
         -9.60191999e+00,  -3.00540096e+03,  -2.40048000e+02],
       [  0.00000000e+00,   0.00000000e+00,   0.00000000e+00,
          4.02640512e+03,   1.28025600e+01,  -1.60032000e+02],
       [  0.00000000e+00,   0.00000000e+00,   0.00000000e+00,
         -1.25025000e+01,  -3.00810150e+03,  -2.50050000e+02],
       [  0.00000000e+00,   0.00000000e+00,   0.00000000e+00,
          4.02858333e+03,   1.66700000e+01,  -2.00040000e+02]])

Mathematically the Jacobian matrix is defined to be

\[\frac{\partial\mathbf{x}_P}{\partial\mathbf{c}} = \left[\begin{array}{cccc} \frac{\partial\mathbf{x}_P}{\partial f} & \frac{\partial\mathbf{x}_P}{\mathbf{k}} & \frac{\partial\mathbf{x}_P}{\partial\mathbf{a}} & \frac{\partial\mathbf{x}_P}{\partial\boldsymbol{\delta\theta}}\end{array}\right]\]

where, using the chain rule,

\begin{gather} \frac{\partial\mathbf{x}_P}{\partial f} = \frac{\partial\mathbf{x}_P}{\partial\mathbf{x}_I} \frac{\partial\mathbf{x}_I}{\partial f} \\ \frac{\partial\mathbf{x}_P}{\partial\boldsymbol{\delta\theta}} = \frac{\partial\mathbf{x}_P}{\partial\mathbf{x}_I} \frac{\partial\mathbf{x}_I}{\partial\mathbf{x}_C'} \frac{\partial\mathbf{x}_C'}{\partial\boldsymbol{\delta\theta}} \end{gather}

and

\begin{gather} \frac{\partial\mathbf{x}_P}{\partial\mathbf{k}} = \left[\begin{array}{cccc} x_I & 0 & 1 & 0 \\ 0 & y_I & 0 & 1 \end{array}\right] \\ \frac{\partial\mathbf{x}_P}{\partial\mathbf{a}} = \left[\begin{array}{cc} k_x & 0 \\ 0 & k_y\end{array}\right] \mathbf{x}_I \left[\begin{array}{ccc} T & T^2 & T^3 \end{array}\right] \\ \frac{\partial\mathbf{x}_P}{\partial\mathbf{x}_I} = (1+a_1T+a_2T^2+a_3T^3) \left[\begin{array}{cc} k_x & 0 \\ 0 & k_y \end{array}\right] \\ \frac{\partial\mathbf{x}_I}{\partial f} = \frac{1}{z_C'} \left[\begin{array}{c}x_C'\\y_C'\end{array}\right] \\ \frac{\partial\mathbf{x}_P}{\partial\mathbf{x}_I} = (1+a_1T+a_2T^2+a_3T^3) \left[\begin{array}{cc} k_x & 0 \\ 0 & k_y \end{array}\right] \\ \frac{\partial\mathbf{x}_I}{\partial\mathbf{x}_C'} = \frac{f}{z_C'}\left[ \begin{array}{ccc}1 & 0 & \frac{-x_C'}{z_C'} \\ 0 & 1 & \frac{-y_C'}{z_C'} \end{array}\right] \\ \frac{\partial\mathbf{x}_C'}{\partial\boldsymbol{\delta\theta}} = \left[\mathbf{x}_C\times\right] \end{gather}

where \(\mathbf{k}=[k_x \quad k_y \quad p_x \quad p_y]\) is a vector of the intrinsic camera parameters, \(\mathbf{a}=[a_1 \quad a_2 \quad a_3]\) is a vector of the temperature dependence coefficients, \(\mathbf{x}_C'\) is the camera frame point after applying the misalignment, \(\boldsymbol{\delta\theta}\) is the misalignment vector, \(\mathbf{x}_C\) is the camera frame point before misalignment is applied, \(\left[\bullet\times\right]\) is the skew-symmetric cross product matrix formed from \(\bullet\), \(\mathbf{x}_P\) is the pixel location, \(\mathbf{x}_I\) is the gnomic location, \(a_{1-3}\) are the temperature coefficients, \(T\) is the temperature, \(k_x\) is the inverse of the pixel pitch in the x direction, and \(k_y\) is the inverse of the pixel pitch in the y direction.