Isis 3 Programmer Reference
LoMediumDistortionMap.cpp
Go to the documentation of this file.
1 
24 #include "LoMediumDistortionMap.h"
25 
26 #include <iostream>
27 #include <iomanip>
28 
29 #include "CameraFocalPlaneMap.h"
30 #include "IString.h"
31 
32 using namespace std;
33 
34 namespace Isis {
47  LoMediumDistortionMap::LoMediumDistortionMap(Camera *parent) :
48  CameraDistortionMap(parent, -1) {
49  }
50 
51 
89  void LoMediumDistortionMap::SetDistortion(const int naifIkCode) {
90  // Get the distortion center (point of symmetry of distortion)
92  double boreS = p_camera->FocalPlaneMap()->DetectorSample();
93  double boreL = p_camera->FocalPlaneMap()->DetectorLine();
94  QString centkey = "INS" + toString(naifIkCode) + "_POINT_OF_SYMMETRY";
95  p_sample0 = boreS - p_camera->Spice::getDouble(centkey, 0);
96  p_line0 = boreL + p_camera->Spice::getDouble(centkey, 1);
97 
98  // Get the distortion coefficients
100  }
101 
102 
118  const double dy) {
119 
120  // Set sRef needed for lo medium distortion algorithm
121  double sRef = 5000.;
122 
123  // lo medium distortion algorithm is applied in the image plane so convert back to sample/line
124  p_focalPlaneX = dx;
125  p_focalPlaneY = dy;
126 
127  // Test for extraneous data. Maximum x is about 38.045 and maximum y is about 31.899.
128  // First tried adding 10% and it was sufficient to trim off extraneous data, but
129  // also prevented lat/lons from being calculated to the images edges. Increased x to
130  // 20.361224% to pick up image edges. 3171 was the test image.
131  // 17.5% to pick up image edges. 3171 was the test image.
132  if(fabs(dx) > 45.79142767 || fabs(dy) > 35.09) return false;
133 
135  double ds = p_camera->FocalPlaneMap()->DetectorSample();
136  double dl = p_camera->FocalPlaneMap()->DetectorLine();
137 
138  // Translate the focal plane x/y coordinate to be relative to the
139  // distortion point of symmetry
140  double dists = ds - p_sample0;
141  double distl = dl - p_line0;
142 
143  // Get the distance from the focal plane center and if we are close
144  // skip the distortion
145  double origr2 = dists * dists + distl * distl;
146  double sp = sqrt(origr2); // pixels
147  if(sp <= .000001) {
150  return true;
151  }
152 
153  // Otherwise remove distortion
154  // Use the distorted radial coordinate, rp (r prime), to estimate the ideal radial coordinate
155  double nS = sp / sRef;
156  double dS = p_odk[0] * nS + p_odk[1] * pow(nS, 3) + p_odk[2] * pow(nS, 5);
157  double prevdS = 2 * dS;
158  double pixtol = .000001;
159  int numit = 0;
160  double s;
161 
162  // Now use the estimate to compute the radial coordinate and get an improved estimate
163  while(fabs(dS - prevdS) > pixtol) {
164 
165  if(numit > 14 || fabs(dS) > 1E9) {
166  dS = 0.;
167  // if (numit > 14) cout<<"Too many iterations"<<endl;
168  // if (fabs(dS) > 1E9) cout<<"Diverging"<<endl;
169  break;
170  }
171 
172  prevdS = dS;
173  s = sp - dS;
174  nS = s / sRef;
175  dS = p_odk[0] * nS + p_odk[1] * pow(nS, 3) + p_odk[2] * pow(nS, 5);
176  numit++;
177  }
178 
179  s = sp - dS;
180  double ratio = s / sp;
181  double undistortedSample = dists * ratio + p_sample0;
182  double undistortedLine = distl * ratio + p_line0;
183  p_camera->FocalPlaneMap()->SetDetector(undistortedSample, undistortedLine);
186  return true;
187  }
188 
189 
209  const double uy) {
210 
213 
214  // Test for data outside of image (image bounds plus 10% for y and 20.361224% for x)
215  if(fabs(ux) > 45.79142767 || fabs(uy) > 35.09) return false;
216  if(fabs(ux) > 41.85 || fabs(uy) > 35.09) return false;
217 
218  // Set sRef needed for lo medium distortion algorithm
219  double sRef = 5000.;
220 
221  // The algorithm is applied in the image plane so convert back to sample/line
223  double us = p_camera->FocalPlaneMap()->DetectorSample();
224  double ul = p_camera->FocalPlaneMap()->DetectorLine();
225 
226  // Translate the distorted x/y coordinate to be relative to the
227  // distortion point of symmetry
228  double distus = us - p_sample0;
229  double distul = ul - p_line0;
230 
231  // Compute the distance from the focal plane center and if we are
232  // close to the center then no distortion is required
233  double rp2 = (distus * distus + distul * distul); // pixels squared
234 
235  if(rp2 < 1.0E-6) {
238  return true;
239  }
240 
241  // Add distortion. First compute fractional distortion at rp (r-prime)
242  rp2 = rp2 / sRef / sRef;
243  double drOverR = (p_odk[0] + rp2 * p_odk[1] + rp2 * rp2 * p_odk[2]) / sRef;
244 
245  // Compute the focal plane distorted s/l
246  double ds = p_sample0 + (distus * (1. + drOverR));
247  double dl = p_line0 + (distul * (1. + drOverR));
248 
249  p_camera->FocalPlaneMap()->SetDetector(ds, dl);
252  return true;
253  }
254 }
double p_focalPlaneX
Distorted focal plane x.
virtual bool SetFocalPlane(const double dx, const double dy)
Compute undistorted focal plane x/y for Lo Medium Resolution Camera.
virtual bool SetFocalPlane(const double dx, const double dy)
Compute detector position (sample,line) from focal plane coordinates.
Namespace for the standard library.
QString toString(bool boolToConvert)
Global function to convert a boolean to a string.
Definition: IString.cpp:226
virtual void SetDistortion(int naifIkCode)
Load distortion coefficients.
virtual bool SetDetector(const double sample, const double line)
Compute distorted focal plane coordinate from detector position (sampel,line)
Unless noted otherwise, the portions of Isis written by the USGS are public domain.
double p_undistortedFocalPlaneX
Undistorted focal plane x.
Camera * p_camera
The camera to distort/undistort.
Distort/undistort focal plane coordinates.
double p_sample0
Center of distortion on sample axis.
const double E
Sets some basic constants for use in ISIS programming.
Definition: Constants.h:55
CameraFocalPlaneMap * FocalPlaneMap()
Returns a pointer to the CameraFocalPlaneMap object.
Definition: Camera.cpp:2848
double p_focalPlaneY
Distorted focal plane y.
double p_undistortedFocalPlaneY
Undistorted focal plane y.
Namespace for ISIS/Bullet specific routines.
Definition: Apollo.h:31
virtual bool SetUndistortedFocalPlane(const double ux, const double uy)
Compute distorted focal plane x/y for Lo Medium Resolution Camera.
void SetDistortion(const int naifIkCode)
Load LO Medium Resolution Camera perspective & distortion coefficients.
std::vector< double > p_odk
Vector of distortion coefficients.
double p_line0
Center of distortion on line axis.