Isis 3 Programmer Reference
NewHorizonsLorriDistortionMap.cpp
1 
7 /* SPDX-License-Identifier: CC0-1.0 */
8 
9 #include "NewHorizonsLorriDistortionMap.h"
10 #include "CameraFocalPlaneMap.h"
11 
12 using namespace std;
13 
14 namespace Isis {
15 
30  NewHorizonsLorriDistortionMap::NewHorizonsLorriDistortionMap(Camera *parent, double e2, double e5,
31  double e6, double zDirection) : CameraDistortionMap(parent, zDirection) {
32  p_e2 = e2;
33  p_e5 = e5;
34  p_e6 = e6;
35  }
36 
37 
50  bool NewHorizonsLorriDistortionMap::SetFocalPlane(const double dx, const double dy) {
51 
52  p_focalPlaneX = dx;
53  p_focalPlaneY = dy;
54 
55  // Reducing to the principle point offset (xp,yp)
56  double x = dx;
57  double y = dy;
58 
59  // r is the distance between the principal point and the measured point on the image
60  double rr = x * x + y * y;
61 
62  // dr is the radial distortion contribution
63  // The equation below was changed from all + to all - to adjust the distortion model to fit
64  // the definition of the LORRI distortion. The original version with +, was defined from Bill
65  // Owen's paper with an assumption the xs and ys in the equations were distorted x,ys. After
66  // meeting with the LORRI team +s were changed to -s to account for the x,ys actually being
67  // undistorted focal plane positions. That is, the undistorted focal plane positions are
68  // closer to the center of the image than the distorted focal plane positions.
69  // NOTE: The discussions showed the Ky and e5 values needed to be negated. The e5 value has
70  // now been negated in the LORRI IK, and the y is now negated in the equation below.
71  // NOTE: The Y and Line values can not be negated in the transY and transL affines because
72  // this would cause the p_xxxxx class member variables to be in a flipped (top to bottom)
73  // coordinate system relative to the SPICE defined focal plane coordinate system.
74  double dr = 1.0 - rr * p_e2 - y * p_e5 - x * p_e6;
75 
76  // Image coordinates corrected for distortion
77  p_undistortedFocalPlaneX = x * dr;
78  p_undistortedFocalPlaneY = y * dr;
79 
80  return true;
81 
82  }
83 
84 
97  bool NewHorizonsLorriDistortionMap::SetUndistortedFocalPlane(const double ux, const double uy) {
98 
99  // Image coordinates prior to introducing distortion
102 
103  double xt = ux;
104  double yt = uy;
105 
106  double xx, yy, xy, rr;
107  double xdistortion, ydistortion;
108  double xdistorted, ydistorted;
109  double xprevious, yprevious;
110 
111  xprevious = 1000000.0;
112  yprevious = 1000000.0;
113 
114  double tolerance = 0.000001;
115 
116  bool bConverged = false;
117 
118  // Iterating to introduce distortion...
119  // We stop when the difference between distorted coordinates
120  // in successive iterations is at or below the given tolerance
121  for (int i = 0; i < 50; i++) {
122  xx = xt * xt;
123  yy = yt * yt;
124  xy = xt * yt;
125  rr = xx + yy;
126 
127  // Distortion at the current point location
128  xdistortion = xt * rr * p_e2 + xy * p_e5 + xx * p_e6;
129  ydistortion = yt * rr * p_e2 + yy * p_e5 + xy * p_e6;
130 
131  // Updated image coordinates
132  // Changed to + instead of -. See comment in SetFocalPlane above
133  xt = ux + xdistortion;
134  yt = uy + ydistortion;
135 
136  // Distorted point corrected for principal point
137  xdistorted = xt; // No PP for LORRI
138  ydistorted = yt; // No PP for LORRI
139 
140  // Check for convergence
141  if ((fabs(xt - xprevious) <= tolerance) && (fabs(yt - yprevious) <= tolerance)) {
142  bConverged = true;
143  break;
144  }
145 
146  xprevious = xt;
147  yprevious = yt;
148  }
149 
150  if (bConverged) {
151  p_focalPlaneX = xdistorted;
152  p_focalPlaneY = ydistorted;
153  }
154 
155  return bConverged;
156 
157  }
158 }
Isis::CameraDistortionMap::p_focalPlaneX
double p_focalPlaneX
Distorted focal plane x.
Definition: CameraDistortionMap.h:65
Isis::CameraDistortionMap::p_undistortedFocalPlaneY
double p_undistortedFocalPlaneY
Undistorted focal plane y.
Definition: CameraDistortionMap.h:68
Isis::CameraDistortionMap::p_undistortedFocalPlaneX
double p_undistortedFocalPlaneX
Undistorted focal plane x.
Definition: CameraDistortionMap.h:67
Isis::Camera
Definition: Camera.h:236
Isis::NewHorizonsLorriDistortionMap::SetFocalPlane
bool SetFocalPlane(const double ux, const double uy)
Compute undistorted focal plane x/y.
Definition: NewHorizonsLorriDistortionMap.cpp:50
Isis::NewHorizonsLorriDistortionMap::SetUndistortedFocalPlane
bool SetUndistortedFocalPlane(const double dx, const double dy)
Compute distorted focal plane x/y.
Definition: NewHorizonsLorriDistortionMap.cpp:97
Isis::CameraDistortionMap
Distort/undistort focal plane coordinates.
Definition: CameraDistortionMap.h:41
std
Namespace for the standard library.
Isis::CameraDistortionMap::p_focalPlaneY
double p_focalPlaneY
Distorted focal plane y.
Definition: CameraDistortionMap.h:66
Isis
This is free and unencumbered software released into the public domain.
Definition: Apollo.h:16