USGS

Isis 3.0 Object Programmers' Reference

Home

NewHorizonsLorriDistortionMap.cpp
Go to the documentation of this file.
1 
22 #include "CameraFocalPlaneMap.h"
23 
24 using namespace std;
25 
26 namespace Isis {
27 
42  NewHorizonsLorriDistortionMap::NewHorizonsLorriDistortionMap(Camera *parent, double e2, double e5,
43  double e6, double zDirection) : CameraDistortionMap(parent, zDirection) {
44  p_e2 = e2;
45  p_e5 = e5;
46  p_e6 = e6;
47  }
48 
49 
62  bool NewHorizonsLorriDistortionMap::SetFocalPlane(const double dx, const double dy) {
63 
64  p_focalPlaneX = dx;
65  p_focalPlaneY = dy;
66 
67  // Reducing to the principle point offset (xp,yp)
68  double x = dx;
69  double y = dy;
70 
71  // r is the distance between the principal point and the measured point on the image
72  double rr = x * x + y * y;
73 
74  // dr is the radial distortion contribution
75  double dr = 1.0 + rr * p_e2 + y * p_e5 + x * p_e6;
76 
77  // Image coordinates corrected for distortion
78  p_undistortedFocalPlaneX = x * dr;
79  p_undistortedFocalPlaneY = y * dr;
80 
81  return true;
82 
83  }
84 
85 
98  bool NewHorizonsLorriDistortionMap::SetUndistortedFocalPlane(const double ux, const double uy) {
99 
100  // Image coordinates prior to introducing distortion
101  p_undistortedFocalPlaneX = ux;
102  p_undistortedFocalPlaneY = uy;
103 
104  double xt = ux;
105  double yt = uy;
106 
107  double xx, yy, xy, rr;
108  double xdistortion, ydistortion;
109  double xdistorted, ydistorted;
110  double xprevious, yprevious;
111 
112  xprevious = 1000000.0;
113  yprevious = 1000000.0;
114 
115  double tolerance = 0.000001;
116 
117  bool bConverged = false;
118 
119  // Iterating to introduce distortion...
120  // We stop when the difference between distorted coordinates
121  // in successive iterations is at or below the given tolerance
122  for (int i = 0; i < 50; i++) {
123  xx = xt * xt;
124  yy = yt * yt;
125  xy = xt * yt;
126  rr = xx + yy;
127 
128  // Distortion at the current point location
129  xdistortion = xt * rr * p_e2 + xy * p_e5 + xx * p_e6;
130  ydistortion = yt * rr * p_e2 + yy * p_e5 + xy * p_e6;
131 
132  // Updated image coordinates
133  xt = ux - xdistortion;
134  yt = uy - ydistortion;
135 
136  // Distorted point corrected for principal point
137  xdistorted = xt; // No PP for Lorrie
138  ydistorted = yt; // No PP for Lorrie
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 }