pydistort.ZernikeDistortion#
- class pydistort.ZernikeDistortion(parameters: ndarray | None = None, Nzer: Integral | None = None, Nparams: Integral | None = None, radius: float | None = None, radius_x: float | None = None, radius_y: float | None = None, center: Tuple[float, float] | None = None)[source]#
Class to apply distortion with a Zernike polynomial model.
Distort the given
normalized_points
using the distortion model to obtain thedistorted_points
.\[x_D = \text{distort}(x_N, \lambda_1, \lambda_2, \lambda_3, \ldots)\]The distorted points can be calculated using the normalized points and the Zernike coefficients as follows:
\[x_{D} = x_{N} + \sum_{n=0}^{N_{zer}} \sum_{m=-n}^{n} C^{x}_{n,m} Z_{nm}(\rho, \theta)\]\[y_{D} = y_{N} + \sum_{n=0}^{N_{zer}} \sum_{m=-n}^{n} C^{y}_{n,m} Z_{nm}(\rho, \theta)\]where \(Z_{nm}(\rho, \theta)\) are the Zernike polynomials, \(C^{x}_{n,m}\) and \(C^{y}_{n,m}\) are the Zernike coefficients for the x and y coordinates, respectively, and \(\rho\) and \(\theta\) are the polar coordinates of the normalized points in the defined unit ellipse with radius \(R_x, R_y\) and center \((x_0, y_0)\):
\[\rho = \sqrt{\left(\frac{x_{N} - x_{0}}{R_x}\right)^2 + \left(\frac{y_{N} - y_{0}}{R_y}\right)^2}\]\[\theta = \arctan2(y_{N} - y_{0}, x_{N} - x_{0})\]Note
For more informations about Zernike polynomials, see the package pyzernike (Artezaru/pyzernike).
Only coeffcients for \(n \leq N_{zer}\) and \(m \in [-n, n]\) and \(n-m \equiv 0 \mod 2\) are stored.
The coefficients are storred in a
parameters
1D-array with the OSA/ANSI standard indices but with a x/y separation:C^{x}_{0,0}, parameters[0] for the x coordinate \(n=0, m=0\)
C^{y}_{0,0}, parameters[1] for the y coordinate \(n=0, m=0\)
C^{x}_{1,-1}, parameters[2] for the x coordinate \(n=1, m=-1\)
C^{y}_{1,-1}, parameters[3] for the y coordinate \(n=1, m=-1\)
C^{x}_{1,1}, parameters[4] for the x coordinate \(n=1, m=1\)
C^{y}_{1,1}, parameters[5] for the y coordinate \(n=1, m=1\)
C^{x}_{2,-2}, parameters[6] for the x coordinate \(n=2, m=-2\)
C^{y}_{2,-2}, parameters[7] for the y coordinate \(n=2, m=-2\)
C^{x}_{2,0}, parameters[8] for the x coordinate \(n=2, m=0\)
C^{y}_{2,0}, parameters[9] for the y coordinate \(n=2, m=0\)
C^{x}_{2,2}, parameters[10] for the x coordinate \(n=2, m=2\)
C^{y}_{2,2}, parameters[11] for the y coordinate \(n=2, m=2\)
…
If the number of input parameters is not equal to the number of parameters required by the model, the other parameters are set to 0.
The number of parameters is given by the formula:
\[N_{params} = (N_{zer}+1)(N_{zer}+2)\]Ordre of Zernike
Nzer
Nparameters for X or Y
Nparameters in model
Nparams
None
0
0
0
1
2
1
3
6
2
6
12
3
10
20
4
15
30
Warning
If the ordre of the zernike polynomials
Nzer
is given during instantiation, the given parameters are truncated or extended to the given number of parameters. Same for the number of parametersNparams
.To compute the Distortion, the user must define the unit circle in which the normalized points are defined. The unit circle is defined by the radius \(R\) and the center \((x_{0}, y_{0})\).
- Parameters:
parameters (numpy.ndarray, optional) – The parameters of the distortion model. If None, no distortion is applied. The default is None.
Nzer (int, optional) – The order of the Zernike polynomials. If None, the order is set according to the number of parameters. The default is None. Only use
Nzer
orNparams
, not both.Nparams (int, optional) – The number of parameters of the distortion model. If None, the number of parameters is set according to the order of the Zernike polynomials. The default is None. Only use
Nzer
orNparams
, not both.radius (float, optional) – The radius of the unit circle in which the normalized points are defined. If None, the radius is set to 1.0. The default is None. This radius is used to normalize the points in the unit circle. If given,
radius_x
andradius_y
must be None.radius_x (float, optional) – The radius of the unit ellipse in which the normalized points are defined along x-axis. If None, the radius is set to 1.0. The default is None.
radius_y (float, optional) – The radius of the unit ellipse in which the normalized points are defined along y-axis. If None, the radius is set to 1.0. The default is None.
center (Tuple[float, float], optional) – The center of the unit circle in which the normalized points are defined. If None, the center is set to (0.0, 0.0). The default is None.
Examples
Create an distortion object with a given model:
import numpy as np from pydistort import ZernikeDistortion # Create a distortion object with 8 parameters distortion = ZernikeDistortion(np.array([0.1, 0.2, 0.3, 0.4, 0.5, 0.6])) # Model with Nzer=1, -> Nparams=6
Then you can use the distortion object to distort
normalized_points
:normalized_points = np.array([[0.1, 0.2], [0.3, 0.4], [0.5, 0.6]]) # shape (3, 2) result = distortion.transform(normalized_points): #alias distort is also available distorted_points = result.distorted_points # shape (3, 2) -> distorted points in (normalized) image coordinates print(distorted_points)
You can also access to the jacobian of the distortion transformation:
result = distortion.transform(distorted_points, dx=True, dp=True) distorted_points_dx = result.jacobian_dx # Jacobian of the image points with respect to the distorted points # Shape (..., 2, 3) distorted_points_dp = result.jacobian_dp # Jacobian of the image points with respect to the intrinsic parameters # Shape (..., 2, Nparams)
The parameters are ordered as given in the model description above.
- property Nparams: int#
Get the number of parameters of the distortion model.
- Returns:
The number of parameters of the distortion model.
- Return type:
int
- property Nzer: Integral | None#
Get the order of the Zernike polynomials.
- Returns:
The order of the Zernike polynomials. If no distortion is applied, returns None.
- Return type:
Optional[Integral]
- _abc_impl = <_abc._abc_data object>#
- _inverse_transform(distorted_points: ndarray, *, dx: bool = False, dp: bool = False, max_iter: int = 10, eps: float = 1e-08) tuple[ndarray, ndarray | None, ndarray | None] [source]#
This method is called by the
pydistort.Transform.inverse_transform()
method to perform the inverse distortion transformation. This method allows to transform thedistorted_points
to thenormalized_points
using the distortion model.Note
For
_inverse_transform
the input must have shape (Npoints, 2) with float64 type. The output has shape (Npoints, 2) and the jacobian are always None.To achieve the undistortion, an iterative algorithm is used to find the normalized points that correspond to the distorted points. The algorithm is based on the following equations:
\[\begin{split}\begin{bmatrix} x_N [\text{it }k+1]\\ y_N [\text{it }k+1] \end{bmatrix} = \begin{bmatrix} x_D - \Delta x [\text{it }k] \\ y_D - \Delta y [\text{it }k] \end{bmatrix}\end{split}\]Where \(\Delta x [\text{it }k]\) and :math:Delta y [text{it }k] are the corrections applied to the distorted points at iteration \(k\). The corrections are computed using the following equations:
\[\begin{split}\begin{bmatrix} \Delta x [\text{it }k]\\ \Delta y [\text{it }k] \end{bmatrix} = \begin{bmatrix} \sum_{n=0}^{N_{zer}} \sum_{m=-n}^{n} C^{x}_{n,m} Z_{nm}(\rho, \theta) [\text{it }k] \\ \sum_{n=0}^{N_{zer}} \sum_{m=-n}^{n} C^{y}_{n,m} Z_{nm}(\rho, \theta) [\text{it }k] \end{bmatrix}\end{split}\]Where \(Z_{nm}(\rho, \theta) [\text{it }k]\) are the Zernike polynomials evaluated at the current normalized points.
Warning
This method is not designed to be used directly for the transformation of points. No checks are performed on the input points, so it is the user’s responsibility to ensure that the input points are valid.
Warning
The method use an iterative algorithm to find the normalized points that correspond to the distorted points. In fact the jacobian with respect to the normalized points is not computed for this inverse transform (always None).
- Parameters:
distorted_points (numpy.ndarray) – Array of distorted points to be transformed with shape (Npoints, 2).
dx (bool, optional) – [ALWAYS False]
dp (bool, optional) – [ALWAYS False]
max_iter (int, optional) – The maximum number of iterations for the optimization algorithm. Default is 10. Use only if
opencv
is False.eps (float, optional) – The tolerance for the optimization algorithm. Default is 1e-8. Use only if
opencv
is False.
- Returns:
normalized_points (numpy.ndarray) – The transformed normalized points in normalized coordinates. It will be a 2D array of shape (Npoints, 2).
jacobian_dx (Optional[numpy.ndarray]) – [ALWAYS None]
jacobian_dp (Optional[numpy.ndarray]) – [ALWAYS None]
- _transform(normalized_points: ndarray, *, dx: bool = False, dp: bool = False) tuple[ndarray, ndarray | None, ndarray | None] [source]#
This method is called by the
pydistort.Transform.transform()
method to perform the distortion transformation. This method allows to transform thenormalized_points
to thedistorted_points
using the distortion model.Note
For
_transform
the input must have shape (Npoints, 2) with float64 type. The output has shape (Npoints, 2) for the image points and (Npoints, 2, 2) for the jacobian with respect to the normalized points and (Npoints, 2, 4) for the jacobian with respect to the distortion parameters.The equation used for the transformation is given in the main documentation of the class.
We have (only for the x coordinate):
\[x_{D} = x_{N} + \sum_{n=0}^{N_{zer}} \sum_{m=-n}^{n} C^{x}_{n,m} Z_{nm}(\rho, \theta)\]\[x_{D} = x_{N} + \sum_{n=0}^{N_{zer}} \sum_{m=-n}^{n} C^{x}_{n,m} Z_{nm}(\sqrt{\left(\frac{x_{N} - x_{0}}{R_x}\right)^2 + \left(\frac{y_{N} - y_{0}}{R_y}\right)^2}, \arctan2(\frac{y_{N} - y_{0}}{R_y}, \frac{x_{N} - x_{0}}{R_x})))\]The derivative of the distorted points with respect to the normalized points is given by:
\[\frac{\partial x_{D}}{\partial x_{N}} = 1 + \sum_{n=0}^{N_{zer}} \sum_{m=-n}^{n} C^{x}_{n,m} \frac{\partial Z_{nm}}{\partial x_{N}}\]\[\frac{\partial x_{D}}{\partial y_{N}} = \sum_{n=0}^{N_{zer}} \sum_{m=-n}^{n} C^{x}_{n,m} \frac{\partial Z_{nm}}{\partial y_{N}}\]Where:
\[\frac{\partial Z_{nm}}{\partial x_{N}} = \frac{\partial Z_{nm}}{\partial \rho} \cdot \frac{\partial \rho}{\partial x_{N}} + \frac{\partial Z_{nm}}{\partial \theta} \cdot \frac{\partial \theta}{\partial x_{N}}\]\[\frac{\partial Z_{nm}}{\partial y_{N}} = \frac{\partial Z_{nm}}{\partial \rho} \cdot \frac{\partial \rho}{\partial y_{N}} + \frac{\partial Z_{nm}}{\partial \theta} \cdot \frac{\partial \theta}{\partial y_{N}}\]See also
Package
pyzernike
(Artezaru/pyzernike) for the implementation of the Zernike polynomials and their derivatives.Warning
This method is not designed to be used directly for the transformation of points. No checks are performed on the input points, so it is the user’s responsibility to ensure that the input points are valid.
- Parameters:
normalized_points (numpy.ndarray) – Array of normalized points to be transformed with shape (Npoints, 2).
dx (bool, optional) – If True, the Jacobian of the distorted points with respect to the normalized points is computed. Default is False. The output will be a 2D array of shape (Npoints, 2, 2).
dp (bool, optional) – If True, the Jacobian of the distorted points with respect to the distortion parameters is computed. Default is False. The output will be a 2D array of shape (Npoints, 2, Nparams).
- Returns:
distorted_points (numpy.ndarray) – The transformed distorted points in normalized coordinates. It will be a 2D array of shape (Npoints, 2).
jacobian_dx (Optional[numpy.ndarray]) – The Jacobian of the distorted points with respect to the normalized points if
dx
is True. Otherwise None. It will be a 2D array of shape (Npoints, 2, 2).jacobian_dp (Optional[numpy.ndarray]) – The Jacobian of the distorted points with respect to the distortion parameters if
dp
is True. Otherwise None. It will be a 2D array of shape (Npoints, 2, Nparams).
- property center: ndarray#
Get the center of the unit circle in which the normalized points are defined.
- Returns:
The center of the unit circle as a 1D array with two elements (x, y).
- Return type:
numpy.ndarray
- get_Cx(n: Integral, m: Integral) Number [source]#
Get the Zernike coefficient for the x coordinate.
- Parameters:
n (int) – The order of the Zernike polynomial.
m (int) – The azimuthal frequency of the Zernike polynomial.
- Returns:
The value of the Zernike coefficient.
- Return type:
float
- get_Cy(n: Integral, m: Integral) Number [source]#
Get the Zernike coefficient for the y coordinate.
- Parameters:
n (int) – The order of the Zernike polynomial.
m (int) – The azimuthal frequency of the Zernike polynomial.
- Returns:
The value of the Zernike coefficient.
- Return type:
float
- get_index(n: Integral, m: Integral, coord: str) int [source]#
Get the index of the Zernike coefficient for the given order and azimuthal frequency.
\[j = n(n+2) + m + (0 \text{ if } coord = 'x' \text{ else } 1)\]- Parameters:
n (int) – The order of the Zernike polynomial.
m (int) – The azimuthal frequency of the Zernike polynomial.
coord (str) – The coordinate (‘x’ or ‘y’) for which to get the index.
- Returns:
The index of the Zernike coefficient in the parameters array.
- Return type:
int
- is_set() bool [source]#
Check if the distortion model is set.
- Returns:
True if the distortion model is set, False otherwise.
- Return type:
bool
- property parameters: ndarray | None#
Get the parameters of the distortion model.
- Returns:
The parameters of the distortion model.
- Return type:
numpy.ndarray
- property parameters_x: ndarray | None#
Get the Zernike coefficients for the x coordinate.
- Returns:
The Zernike coefficients for the x coordinate. If no distortion is applied, returns None.
- Return type:
Optional[numpy.ndarray]
- property parameters_y: ndarray | None#
Get the Zernike coefficients for the y coordinate.
- Returns:
The Zernike coefficients for the y coordinate. If no distortion is applied, returns None.
- Return type:
Optional[numpy.ndarray]
- property radius: float#
Get the radius of the unit circle in which the normalized points are defined.
- Returns:
The radius of the unit circle.
- Return type:
float
- property radius_x: float#
Get the radius of the unit ellipse in which the normalized points are defined along x-axis.
- Returns:
The radius of the unit ellipse along x-axis.
- Return type:
float
- property radius_y: float#
Get the radius of the unit ellipse in which the normalized points are defined along y-axis.
- Returns:
The radius of the unit ellipse along y-axis.
- Return type:
float
- set_Cx(n: Integral, m: Integral, value: Number) None [source]#
Set the Zernike coefficient for the x coordinate.
- Parameters:
n (int) – The order of the Zernike polynomial.
m (int) – The azimuthal frequency of the Zernike polynomial.
value (float) – The value of the Zernike coefficient.
- set_Cy(n: Integral, m: Integral, value: Number) None [source]#
Set the Zernike coefficient for the y coordinate.
- Parameters:
n (int) – The order of the Zernike polynomial.
m (int) – The azimuthal frequency of the Zernike polynomial.
value (float) – The value of the Zernike coefficient.
- property x0: float#
Get the x coordinate of the center of the unit circle in which the normalized points are defined.
- Returns:
The x coordinate of the center of the unit circle.
- Return type:
float
- property y0: float#
Get the y coordinate of the center of the unit circle in which the normalized points are defined.
- Returns:
The y coordinate of the center of the unit circle.
- Return type:
float