4.3.1. discorpy.proc.processing

Module of processing methods:

  • Fit lines of dots to parabolas, find the center of distortion.

  • Calculate undistorted intercepts of gridlines.

  • Calculate distortion coefficients of the backward model, the forward model, and the backward-from-forward model.

  • Correct perspective distortion affecting curve lines.

  • Generate non-perspective points or lines from perspective points or lines.

  • Calculate perspective coefficients.

Functions:

_para_fit_hor(list_lines, xcenter, ycenter)

Fit horizontal lines of dots to parabolas.

_para_fit_ver(list_lines, xcenter, ycenter)

Fit vertical lines of dots to parabolas.

find_cod_coarse(list_hor_lines, list_ver_lines)

Coarse estimation of the center of distortion.

find_cod_fine(list_hor_lines, ...)

Find the best center of distortion (CoD) by searching around the coarse estimation of the CoD.

_calc_undistor_intercept(list_hor_lines, ...)

Calculate the intercepts of undistorted lines.

calc_coef_backward(list_hor_lines, ...[, ...])

Calculate the distortion coefficients of a backward mode.

calc_coef_forward(list_hor_lines, ...[, ...])

Calculate the distortion coefficients of a forward mode.

calc_coef_backward_from_forward(...[, ...])

Calculate the distortion coefficients of a backward mode from a forward model.

transform_coef_backward_and_forward(list_fact)

Transform polynomial coefficients of a radial distortion model between forward mapping and backward mapping.

find_cod_bailey(list_hor_lines, list_ver_lines)

Find the center of distortion (COD) using the Bailey's approach (Ref.

_generate_non_perspective_parabola_coef(...)

Correct the deviation of fitted parabola coefficients of each line caused by perspective distortion.

_find_cross_point_between_parabolas(...)

Find a cross point between two parabolas.

regenerate_grid_points_parabola(...[, ...])

Regenerating grid points by finding cross points between horizontal lines and vertical lines using their parabola coefficients.

_generate_linear_coef(list_hor_lines, ...[, ...])

Get linear coefficients of horizontal and vertical lines from linear fit.

_find_cross_point_between_lines(...)

Find a cross point between two lines.

_calc_undistor_intercept_perspective(...[, ...])

Calculate the intercepts of undistorted lines from perspective distortion.

regenerate_grid_points_linear(...)

Regenerating grid points by finding cross points between horizontal lines and vertical lines using their linear coefficients.

generate_undistorted_perspective_lines(...)

Generate undistorted lines from perspective lines.

generate_source_target_perspective_points(...)

Generate source points (distorted) and target points (undistorted).

generate_4_source_target_perspective_points(points)

Generate 4 rectangular points corresponding to 4 perspective-distorted points.

calc_perspective_coefficients(source_points, ...)

Calculate perspective coefficients of a matrix to map from source points to target points (Ref.

update_center(list_lines, xcenter, ycenter)

Update the coordinate-center of points on lines.

discorpy.proc.processing.find_cod_coarse(list_hor_lines, list_ver_lines)[source]

Coarse estimation of the center of distortion.

Parameters:
  • list_hor_lines (list of 2D arrays) – List of the (y,x)-coordinates of dot-centroids on each horizontal line.

  • list_ver_lines (list of 2D arrays) – List of the (y,x)-coordinates of dot-centroids on each vertical line.

Returns:

  • xcenter (float) – Center of distortion in x-direction.

  • ycenter (float) – Center of distortion in y-direction.

discorpy.proc.processing.find_cod_fine(list_hor_lines, list_ver_lines, xcenter, ycenter, dot_dist)[source]

Find the best center of distortion (CoD) by searching around the coarse estimation of the CoD.

Parameters:
  • list_hor_lines (list of 2D arrays) – List of the (y,x)-coordinates of dot-centroids on each horizontal line.

  • list_ver_lines (list of 2D arrays) – List of the (y,x)-coordinates of dot-centroids on each vertical line.

  • xcenter (float) – Coarse estimation of the CoD in x-direction.

  • ycenter (float) – Coarse estimation of the CoD in y-direction.

  • dot_dist (float) – Median distance of two nearest dots.

Returns:

  • xcenter (float) – Center of distortion in x-direction.

  • ycenter (float) – Center of distortion in y-direction.

discorpy.proc.processing.calc_coef_backward(list_hor_lines, list_ver_lines, xcenter, ycenter, num_fact, optimizing=False, threshold=0.3)[source]

Calculate the distortion coefficients of a backward mode.

Parameters:
  • list_hor_lines (list of 2D arrays) – List of the (y,x)-coordinates of dot-centroids on each horizontal line.

  • list_ver_lines (list of 2D arrays) – List of the (y,x)-coordinates of dot-centroids on each vertical line.

  • xcenter (float) – Center of distortion in x-direction.

  • ycenter (float) – Center of distortion in y-direction.

  • num_fact (int) – Number of the factors of polynomial.

  • optimizing (bool, optional) – Apply optimization if True.

  • threshold (float) – To determine if there are missing lines. Larger is less sensitive.

Returns:

list_fact (list of float) – Coefficients of the polynomial.

discorpy.proc.processing.calc_coef_forward(list_hor_lines, list_ver_lines, xcenter, ycenter, num_fact, optimizing=False, threshold=0.3)[source]

Calculate the distortion coefficients of a forward mode.

Parameters:
  • list_hor_lines (list of 2D arrays) – List of the (y,x)-coordinates of dot-centroids on each horizontal line.

  • list_ver_lines (list of 2D arrays) – List of the (y,x)-coordinates of dot-centroids on each vertical line.

  • xcenter (float) – Center of distortion in x-direction.

  • ycenter (float) – Center of distortion in y-direction.

  • num_fact (int) – Number of the factors of polynomial.

  • optimizing (bool, optional) – Apply optimization if True.

  • threshold (float) – To determine if there are missing lines. Larger is less sensitive.

Returns:

list_fact (list of float) – Coefficients of the polynomial.

discorpy.proc.processing.calc_coef_backward_from_forward(list_hor_lines, list_ver_lines, xcenter, ycenter, num_fact, optimizing=False, threshold=0.3)[source]

Calculate the distortion coefficients of a backward mode from a forward model.

Parameters:
  • list_hor_lines (list of 2D arrays) – List of the (y,x)-coordinates of dot-centroids on each horizontal line.

  • list_ver_lines (list of 2D arrays) – List of the (y,x)-coordinates of dot-centroids on each vertical line.

  • xcenter (float) – Center of distortion in x-direction.

  • ycenter (float) – Center of distortion in y-direction.

  • num_fact (int) – Number of the factors of polynomial.

  • optimizing (bool, optional) – Apply optimization if True.

  • threshold (float) – To determine if there are missing lines. Larger is less sensitive.

Returns:

  • list_ffact (list of floats) – Polynomial coefficients of the forward model.

  • list_bfact (list of floats) – Polynomial coefficients of the backward model.

discorpy.proc.processing.transform_coef_backward_and_forward(list_fact, mapping='backward', ref_points=None)[source]

Transform polynomial coefficients of a radial distortion model between forward mapping and backward mapping.

Parameters:
  • list_fact (list of floats) – Polynomial coefficients of the radial distortion model.

  • mapping ({‘backward’, ‘forward’}) – Transformation direction.

  • ref_points (list of 1D-arrays, optional) – List of the (y,x)-coordinates of points used for the transformation. Generated if None given.

Returns:

list of floats – Polynomial coefficients of the reversed model.

discorpy.proc.processing.find_cod_bailey(list_hor_lines, list_ver_lines, iteration=2)[source]

Find the center of distortion (COD) using the Bailey’s approach (Ref. [1]).

Parameters:
  • list_hor_lines (list of 2D-arrays) – List of the (y,x)-coordinates of points on each horizontal line.

  • list_ver_lines (list of 2D-arrays) – List of the (y,x)-coordinates of points on each vertical line.

Returns:

  • xcenter (float) – Center of distortion in x-direction.

  • ycenter (float) – Center of distortion in y-direction.

References

[1].. https://www-ist.massey.ac.nz/dbailey/sprg/pdfs/2002_IVCNZ_59.pdf

discorpy.proc.processing.regenerate_grid_points_parabola(list_hor_lines, list_ver_lines, perspective=True)[source]

Regenerating grid points by finding cross points between horizontal lines and vertical lines using their parabola coefficients.

Parameters:
  • list_hor_lines (list of 2D-arrays) – List of the (y,x)-coordinates of points on each horizontal line.

  • list_ver_lines (list of 2D-arrays) – List of the (y,x)-coordinates of points on each vertical line.

  • perspective (bool, optional) – Apply perspective correction if True.

Returns:

  • new_hor_lines (list of 2D-arrays) – List of the updated (y,x)-coordinates of points on each horizontal line.

  • new_ver_lines (list of 2D-arrays) – List of the updated (y,x)-coordinates of points on each vertical line.

discorpy.proc.processing.regenerate_grid_points_linear(list_hor_lines, list_ver_lines)[source]

Regenerating grid points by finding cross points between horizontal lines and vertical lines using their linear coefficients.

Parameters:
  • list_hor_lines (list of 2D-arrays) – List of the (y,x)-coordinates of points on each horizontal line.

  • list_ver_lines (list of 2D-arrays) – List of the (y,x)-coordinates of points on each vertical line.

Returns:

  • new_hor_lines (list of 2D-arrays) – List of the updated (y,x)-coordinates of points on each horizontal line.

  • new_ver_lines (list of 2D-arrays) – List of the updated (y,x)-coordinates of points on each vertical line.

discorpy.proc.processing.generate_undistorted_perspective_lines(list_hor_lines, list_ver_lines, equal_dist=True, scale='mean', optimizing=True)[source]

Generate undistorted lines from perspective lines.

Parameters:
  • list_hor_lines (list of 2D-arrays) – List of the (y,x)-coordinates of points on each horizontal line.

  • list_ver_lines (list of 2D-arrays) – List of the (y,x)-coordinates of points on each vertical line.

  • equal_dist (bool) – Use the condition that lines are equidistant if True.

  • scale ({‘mean’, ‘median’, ‘min’, ‘max’, float}) – Scale option for the undistorted grid.

  • optimizing (bool) – Apply optimization for finding line-distance if True.

Returns:

  • list_uhor_lines (list of 2D-arrays) – List of the (y,x)-coordinates of points on undistorted horizontal lines.

  • list_uver_lines (list of 2D-arrays) – List of the (y,x)-coordinates of points on undistorted vertical lines.

discorpy.proc.processing.generate_source_target_perspective_points(list_hor_lines, list_ver_lines, equal_dist=True, scale='mean', optimizing=True)[source]

Generate source points (distorted) and target points (undistorted).

Parameters:
  • list_hor_lines (list of 2D-arrays) – List of the (y,x)-coordinates of points on each horizontal line.

  • list_ver_lines (list of 2D-arrays) – List of the (y,x)-coordinates of points on each vertical line.

  • equal_dist (bool) – Use the condition that lines are equidistant if True.

  • scale ({‘mean’, ‘median’, ‘min’, ‘max’, float}) – Scale option for the undistorted grid.

  • optimizing (bool) – Apply optimization for finding line-distance if True.

Returns:

  • source_points (list of 1D-arrays) – List of the (y,x)-coordinates of distorted points.

  • target_points (list of 1D-arrays) – List of the (y,x)-coordinates of undistorted points.

discorpy.proc.processing.generate_4_source_target_perspective_points(points, input_order='yx', equal_dist=False, scale='mean')[source]

Generate 4 rectangular points corresponding to 4 perspective-distorted points.

Parameters:
  • points (list of 1D-arrays) – List of the coordinates of 4 perspective-distorted points.

  • input_order ({‘yx’, ‘xy’}) – Order of the coordinates of input-points.

  • equal_dist (bool) – Use the condition that the rectangular making of 4-points is square if True.

  • scale ({‘mean’, ‘min’, ‘max’, float}) – Scale option for the undistorted points.

Returns:

  • source_points (list of 1D-arrays) – List of the (y,x)-coordinates of distorted points.

  • target_points (list of 1D-arrays) – List of the (y,x)-coordinates of undistorted points.

discorpy.proc.processing.calc_perspective_coefficients(source_points, target_points, mapping='backward')[source]

Calculate perspective coefficients of a matrix to map from source points to target points (Ref. [1]). Note that the coordinate of a point are in (y,x)-order. This is to be consistent with other functions in the module.

Parameters:
  • source_points (array_like) – List of the (y,x)-coordinates of distorted points.

  • target_points (array_like) – List of the (y,x)-coordinates of undistorted points.

  • mapping ({‘backward’, ‘forward’}) – To select mapping direction.

Returns:

array_like – 1D array of 8 coefficients.

References

[1].. https://doi.org/10.1016/S0262-8856(98)00183-8

discorpy.proc.processing.update_center(list_lines, xcenter, ycenter)[source]

Update the coordinate-center of points on lines.

Parameters:
  • list_lines (list of 2D-arrays) – List of the (y,x)-coordinates of points on lines.

  • xcenter (float) – X-origin of the coordinate system.

  • ycenter (float) – Y-origin of the coordinate system.

Returns:

list of 2D-arrays.