FisheyeModel¶
giant.camera_models.fisheye_model
:
- class giant.camera_models.fisheye_model.FisheyeModel(intrinsic_matrix=None, fx=None, fy=None, px=None, py=None, alpha=None, kx=None, ky=None, kxy=None, field_of_view=None, use_a_priori=False, distortion_coefficients=None, k1=None, k2=None, k3=None, k4=None, temperature_coefficients=None, a1=None, a2=None, a3=None, misalignment=None, estimation_parameters='basic', n_rows=1, n_cols=1)[source]¶
Bases:
BrownModel
This class provides an implementation of the OpenCV fisheye camera model for projecting 3D points onto images and performing camera calibration.
The OpenCV fisheye camera model is the pinhole camera model combined with a lens distortion model which projects any point along a ray emanating from the camera center (origin of the camera frame) to the same 2D point in an image. Given some 3D point (or direction) expressed in the camera frame, \(\mathbf{x}_C\), the OpenCV model is defined as
The
Fisheye
class is a subclasss ofCameraModel
. This means that it includes implementations for all of the abstract methods defined in theCameraModel
class. This also means that it can be used throughout GIANT as the primary camera model, including within thecalibration
subpackage. If this class is going to be used with thecalibration
subpackage, the user can set which parameters are estimated and which are held fixed by using theestimation_parameters
key word argument when creating an instance of the class or by adjusting theestimation_parameters
instance variable on an instance of the class. Theestimation_parameters
input/attribute is a string or list of strings specifying which parameters to estimate. This means thatestimation_parameters
could be something like'basic'
which would indicate to estimate just the usual parameters, or something like['focal_length', 'ky', 'px', 'py']
to estimate just the terms included in the list.In addition to the standard set of methods for a
CameraModel
subclass, theFisheye
class provides the following additional methods which may or may not be useful to some people:Method
Use
computes the pinhole, image frame, and pixel locations of a 3D point
removes distortion from a point to get the corresponding unitless pinhole location
The
Fisheye
class also provides the following properties for easy getting/setting:Property
Description
the diagonal field of view of the camera in units of degrees
\(f_x\), focal length divided by the pixel pitch in the x direction in units of pixels
\(f_y\), focal length divided by the pixel pitch in the y direction in units of pixels
\(\alpha\), A alpha term for non-retangular pixels
\(p_{x}\), the x axis pixel location of the principal point of the camera in units of pixels
\(p_{y}\), the y axis pixel location of the principal point of the camera in units of pixels
\(k_1\), the radial distortion coefficient corresponding to \(\theta^2\)
\(k_2\), the radial distortion coefficient corresponding to \(\theta^4\)
\(k_3\), the radial distortion coefficient corresponding to \(\theta^6\)
\(k_4\), the radial distortion coefficient corresponding to \(\theta^8\)
\(a_1\), the linear coefficient for focal length dependent focal length
\(a_2\), the quadratic coefficient for focal length dependent focal length
\(a_3\), the cubic coefficient for focal length dependent focal length
The inverse of the intrinsic matrix
Note
The distortion attributes are aliases over each other and refer to the same data. Therefore setting a value to
radial2n
would also change the value ofk1
- Parameters:
intrinsic_matrix (ndarray[tuple[Any, ...], dtype[_ScalarT]] | None) – the intrinsic matrix for the camera as a numpy shape (2, 3) array. Note that this is overwritten if
kx
,ky
,kxy
,kyx
,px
,py
,fx
,fy
,alpha
are also specified.field_of_view (float | None) – The field of view of the camera in units of degrees.
use_a_priori (bool) – A flag to indicate whether to include the a priori state vector in the Jacobian matrix when performing a calibration
distortion_coefficients (ndarray[tuple[Any, ...], dtype[_ScalarT]] | None) – A numpy array of shape (4,) containing the 4 fisheye radial distortion coefficients in order [k1, k2, k3, k4]. Note that this array is overwritten with any distortion coefficients that are specified independently.
fx (float | None) – The pixel pitch along the x axis in units of pixels
fy (float | None) – The pixel pitch along the y axis in units of pixels
kx (float | None) – The pixel pitch along the x axis in units of pixels
ky (float | None) – The pixel pitch along the y axis in units of pixels
kxy (float | None) – An alpha term for non-rectangular pixels
alpha (float | None) – An alpha term for non-rectangular pixels
px (float | None) – the x component of the pixel location of the principal point in the image in units of pixels
py (float | None) – the y component of the pixel location of the principal point in the image in units of pixels
k1 (float | None) – the fisheye radial distortion coefficient corresponding to the theta**2 term
k2 (float | None) – the fisheye radial distortion coefficient corresponding to the theta**4 term
k3 (float | None) – the fisheye radial distortion coefficient corresponding to the theta**6 term
k4 (float | None) – the fisheye radial distortion coefficient corresponding to the theta**8 term
temperature_coefficients (ndarray[tuple[Any, ...], dtype[_ScalarT]] | None) – The temperature polynomial coefficients as a length 3 Sequence
a1 (float | None) – the linear coefficient of the focal length temperature dependence
a2 (float | None) – the quadratic coefficient of the focal length temperature dependence
a3 (float | None) – the cubic coefficient of the focal length temperature dependence
misalignment (ndarray[tuple[Any, ...], dtype[_ScalarT]] | None) – either a numpy array of shape (3,) or a list of numpy arrays of shape(3,) with each array corresponding to a single image (the list of numpy arrays is only valid when estimating multiple misalignments)
estimation_parameters (str | Sequence[str]) – A string or list of strings specifying which model parameters to include in the calibration
n_rows (int) – the number of rows of the active image array
n_cols (int) – the number of columns in the active image array
- n_rows¶
The number of rows in the active pixel array for the camera
- n_cols¶
The number of columns in the active pixel array for the camera
- use_a_priori¶
This boolean value is used to determine whether to append the identity matrix to the Jacobian matrix returned by
compute_jacobian()
in order to include the current estimate of the camera model in the calibration process.
- intrinsic_matrix¶
The 2x3 intrinsic matrix contains the conversion from unitless gnomic locations to a location in an image with units of pixels.
It is defined as
\[\begin{split}\mathbf{K} = \left[\begin{array}{ccc} f_x & \alpha & p_x \\ 0 & f_y & p_y \end{array}\right]\end{split}\]
- temperature_coefficients¶
The coefficients for the polynomial specifying the change in the focal length as a function of temperature.
- estimate_multiple_misalignments¶
This boolean value is used to determine whether multiple misalignments are being estimated/used per image.
If set to
True
then one misalignment is estimated for each image and used for each image when projecting through the camera model. When set toFalse
then a single misalignment is estimated for all images and used for all images when projecting through the camera model. Typically the user shouldn’t be setting this attribute directly as it is automatically handled when setting theestimation_parameters
attribute.
- estimation_parameters¶
A list of strings containing the parameters to estimate when performing calibration with this model.
This list is used in the methods
compute_jacobian()
andapply_update()
to determine which parameters are being estimated/updated. From thecompute_jacobian()
method, only columns of the Jacobian matrix corresponding to the parameters in this list are returned. In theapply_update()
method, the update vector elements are assumed to correspond to the order expressed in this list.Valid values for the elements of this list are shown in the following table. Generally, they correspond to attributes of this class, with a few convenient aliases that point to a collection of attributes.
- distortion_coefficients¶
The distortion coefficients array contains the distortion coefficients for the fisheye model [k1, k2, k3, k4]
- property field_of_view: float¶
A radial field of view of the camera specified in degrees.
The field of view should be set to at least the half width diagonal field of view of the camera. The field of view is used when querying star catalogs.
The diagonal field of view is defined as
+-----------+ | /| | / | | / | | V/ | | O/ | | F/ | | */ | | 2/ | | / | | / | |/ | +-----------+
If you specify this parameter to be
None
, the field of view will be computed using the camera model if possible.
- property k1: float¶
The radial distortion coefficient corresponding to the \(\theta^2\) term
This corresponds to the [0] index of the distortion_coefficients array
- property kx: float¶
The inverse of the pixel pitch along the x axis in units of pix/distance.
This is the conversion factor to convert from gnomic coordinates (in units of distance) to units of pixels. It corresponds to the [0, 0] component of the intrinsic matrix
- property ky: float¶
The inverse of the pixel pitch along the y axis in units of pix/distance.
This is the conversion factor to convert from pinhole coordinates (in units of distance) to units of pixels. It corresponds to the [1, 1] component of the intrinsic matrix
- property fx: float¶
The focal length in units of pixels along the x axis (focal length divided by x axis pixel pitch)
This is an alias to
kx
and points to the (0, 0) index of the intrinsic matrix
- property fy: float¶
The focal length in units of pixels along the y axis (focal length divided by y axis pixel pitch)
This is an alias to
ky
and points to the (1, 1) index of the intrinsic matrix
- property alpha: float¶
An alpha term for non-rectangular pixels.
This is an alias to
kxy
and points to the to the [0, 1] component of the intrinsic matrix
- property kxy: float¶
An alpha term for non-rectangular pixels.
This corresponds to the [0, 1] component of the intrinsic matrix
- property px: float¶
The x pixel location of the principal point of the camera.
The principal point of the camera is the point in the image where the distortion is zero (the point where the optical axis pierces the image). This corresponds to the [0, 2] component of the intrinsic matrix
- property py: float¶
The y pixel location of the principal point of the camera.
The principal point of the camera is the point in the image where the distortion is zero (the point where the optical axis pierces the image). This corresponds to the [1, 2] component of the intrinsic matrix
- property a1: float¶
The linear coefficient for the focal length temperature dependence
This is the first term in the
temperature_coefficients
array and is multiplied by the temperature.
- property a2: float¶
The quadratic coefficient for the focal length temperature dependence
This is the second term in the
temperature_coefficients
array and is multiplied by the temperature squared.
- property a3: float¶
The cubic coefficient for the focal length temperature dependence
This is the third term in the
temperature_coefficients
array and is multiplied by the temperature cubed.
- property intrinsic_matrix_inv: ndarray¶
The inverse of the intrinsic matrix.
The inverse of the intrinsic matrix is used to convert from units of pixels with an origin at the upper left corner of the image to units of distance with an origin at the principal point of the image.
the intrinsic matrix has an analytic inverse which is given by
\[\begin{split}\mathbf{K}^{-1} = \left[\begin{array}{ccc} \frac{1}{f_x} & -\frac{\alpha}{f_xf_y} & \frac{\alpha p_y-f_yp_x}{f_xf_y} \\ 0 & \frac{1}{f_y} & \frac{-p_y}{f_y} \end{array}\right]\end{split}\]- To convert from units of pixels to a unitless, distorted gnomic location you would do::
>>> from giant.camera_models import BrownModel >>> model = BrownModel(kx=5, ky=10, px=100, py=500) >>> ((model.intrinsic_matrix_inv[:, :2]@[[1, 2, 300], [4, 5, 600]]).T + model.intrinsic_matrix_inv[:, 2]).T array([[-19.8, -19.6, 40.] [-49.6, -49.5, 10.]])
Note
The above code will give you distorted gnomic location, while the
pixels_to_gnomic()
will give you undistorted gnomic locations (true pinhole points).Note
Since the intrinsic matrix is defined as a \(2\times 3\) matrix this isn’t a formal inverse. To get the true inverse you need to append a row of [0, 0, 1] to both the intrinsic matrix and intrinsic matrix inverse.
- property state_vector: List[float]¶
Returns the fully realized state vector according to
estimation_parameters
as a length l list.
- property k2: float¶
The radial distortion coefficient corresponding to the \(\theta^4\) term
This corresponds to the [1] index of the distortion_coefficients array
- property k3: float¶
The radial distortion coefficient corresponding to the \(\theta^6\) term
This corresponds to the [2] index of the distortion_coefficients array
- property k4: float¶
The radial distortion coefficient corresponding to the \(\theta^8\) term
This corresponds to the [3] index of the distortion_coefficients array
Summary of Methods
This method transforms 3D points or directions expressed in the camera frame into the corresponding 2D image locations. |
|
Calculates the Jacobian matrix for each observation in unit_vectors_in_camera_frame for each parameter to be estimated as defined in the |
|
This method computes the Jacobian matrix \(\partial\mathbf{x}_P/\partial\mathbf{x}_C\) where \(\mathbf{x}_C\) is a vector in the camera frame that projects to \(\mathbf{x}_P\) which is the pixel location. |
|
This method computes the Jacobian matrix \(\partial\mathbf{x}_C/\partial\mathbf{x}_P\) where \(\mathbf{x}_C\) is a vector in the camera frame that projects to \(\mathbf{x}_P\) which is the pixel location. |
|
This method takes in a delta update to camera parameters (\(\Delta\mathbf{c}\)) and applies the update to the current instance in place. |
|
This method converts pixel image locations to unit vectors expressed in the camera frame. |
|
This method computes undistorted pixel locations (gnomic/pinhole locations) for given distorted pixel locations according to the current model. |
|
A method that takes gnomic pixel locations in units of pixels and applies the appropriate distortion to them. |
|
This method replaces self with the properties of |
|
This method computes the value of the distortion model across an entire image for use in creating distortion maps. |
|
This method takes in an entire image and warps it to remove the distortion specified by the current model. |
|
Returns a deep copy of this object, breaking all references with |
|
Stores this camera model in an |
|
This class method is used to construct a new instance of cls from an |
|
This method adjusts a pixel location to reflect a new image temperature. |
|
This method computes the scaling to the focal length caused by a shift in temperature. |
|
This method takes an input in pixels and computes the undistorted gnomic location in units of distance. |
|
This method prepares a SciPy RegularGridInterpolator for converting pixels into undistorted gnomic locations. |
|
This method takes an input in pixels and approximates the undistorted gnomic location in units of distance. |
|
This method computes and returns the pinhole, and pixel locations for a set of 3D points expressed in the camera frame. |
|
This method applies the distortion model to the specified pinhole (gnomic) locations in the image frame. |
|
Convert a list of estimation parameters into state label names. |
|
This method reset the misalignment terms to all be zero (no misalignment). |
|
This method returns the Rotation object for the misalignment for the requested image. |