File failed to load: https://isis.astrogeology.usgs.gov/9.0.0/Object/assets/jax/output/NativeMML/config.js
Isis 3 Programmer Reference
LoMediumDistortionMap.cpp
1
6
7/* SPDX-License-Identifier: CC0-1.0 */
8
9#include "LoMediumDistortionMap.h"
10
11#include <iostream>
12#include <iomanip>
13
14#include "CameraFocalPlaneMap.h"
15#include "IString.h"
16
17using namespace std;
18
19namespace Isis {
35
36
73
74 void LoMediumDistortionMap::SetDistortion(const int naifIkCode) {
75 // Get the distortion center (point of symmetry of distortion)
76 p_camera->FocalPlaneMap()->SetFocalPlane(0., 0.);
77 double boreS = p_camera->FocalPlaneMap()->DetectorSample();
78 double boreL = p_camera->FocalPlaneMap()->DetectorLine();
79 QString centkey = "INS" + toString(naifIkCode) + "_POINT_OF_SYMMETRY";
80 p_sample0 = boreS - p_camera->Spice::getDouble(centkey, 0);
81 p_line0 = boreL + p_camera->Spice::getDouble(centkey, 1);
82
83 // Get the distortion coefficients
85 }
86
87
103 const double dy) {
104
105 // Set sRef needed for lo medium distortion algorithm
106 double sRef = 5000.;
107
108 // lo medium distortion algorithm is applied in the image plane so convert back to sample/line
109 p_focalPlaneX = dx;
110 p_focalPlaneY = dy;
111
112 // Test for extraneous data. Maximum x is about 38.045 and maximum y is about 31.899.
113 // First tried adding 10% and it was sufficient to trim off extraneous data, but
114 // also prevented lat/lons from being calculated to the images edges. Increased x to
115 // 20.361224% to pick up image edges. 3171 was the test image.
116 // 17.5% to pick up image edges. 3171 was the test image.
117 if(fabs(dx) > 45.79142767 || fabs(dy) > 35.09) return false;
118
119 p_camera->FocalPlaneMap()->SetFocalPlane(dx, dy);
120 double ds = p_camera->FocalPlaneMap()->DetectorSample();
121 double dl = p_camera->FocalPlaneMap()->DetectorLine();
122
123 // Translate the focal plane x/y coordinate to be relative to the
124 // distortion point of symmetry
125 double dists = ds - p_sample0;
126 double distl = dl - p_line0;
127
128 // Get the distance from the focal plane center and if we are close
129 // skip the distortion
130 double origr2 = dists * dists + distl * distl;
131 double sp = sqrt(origr2); // pixels
132 if(sp <= .000001) {
135 return true;
136 }
137
138 // Otherwise remove distortion
139 // Use the distorted radial coordinate, rp (r prime), to estimate the ideal radial coordinate
140 double nS = sp / sRef;
141 double dS = p_odk[0] * nS + p_odk[1] * pow(nS, 3) + p_odk[2] * pow(nS, 5);
142 double prevdS = 2 * dS;
143 double pixtol = .000001;
144 int numit = 0;
145 double s;
146
147 // Now use the estimate to compute the radial coordinate and get an improved estimate
148 while(fabs(dS - prevdS) > pixtol) {
149
150 if(numit > 14 || fabs(dS) > 1E9) {
151 dS = 0.;
152 // if (numit > 14) cout<<"Too many iterations"<<endl;
153 // if (fabs(dS) > 1E9) cout<<"Diverging"<<endl;
154 break;
155 }
156
157 prevdS = dS;
158 s = sp - dS;
159 nS = s / sRef;
160 dS = p_odk[0] * nS + p_odk[1] * pow(nS, 3) + p_odk[2] * pow(nS, 5);
161 numit++;
162 }
163
164 s = sp - dS;
165 double ratio = s / sp;
166 double undistortedSample = dists * ratio + p_sample0;
167 double undistortedLine = distl * ratio + p_line0;
168 p_camera->FocalPlaneMap()->SetDetector(undistortedSample, undistortedLine);
169 p_undistortedFocalPlaneX = p_camera->FocalPlaneMap()->FocalPlaneX();
170 p_undistortedFocalPlaneY = p_camera->FocalPlaneMap()->FocalPlaneY();
171 return true;
172 }
173
174
194 const double uy) {
195
198
199 // Test for data outside of image (image bounds plus 10% for y and 20.361224% for x)
200 if(fabs(ux) > 45.79142767 || fabs(uy) > 35.09) return false;
201 if(fabs(ux) > 41.85 || fabs(uy) > 35.09) return false;
202
203 // Set sRef needed for lo medium distortion algorithm
204 double sRef = 5000.;
205
206 // The algorithm is applied in the image plane so convert back to sample/line
207 p_camera->FocalPlaneMap()->SetFocalPlane(p_undistortedFocalPlaneX, p_undistortedFocalPlaneY);
208 double us = p_camera->FocalPlaneMap()->DetectorSample();
209 double ul = p_camera->FocalPlaneMap()->DetectorLine();
210
211 // Translate the distorted x/y coordinate to be relative to the
212 // distortion point of symmetry
213 double distus = us - p_sample0;
214 double distul = ul - p_line0;
215
216 // Compute the distance from the focal plane center and if we are
217 // close to the center then no distortion is required
218 double rp2 = (distus * distus + distul * distul); // pixels squared
219
220 if(rp2 < 1.0E-6) {
223 return true;
224 }
225
226 // Add distortion. First compute fractional distortion at rp (r-prime)
227 rp2 = rp2 / sRef / sRef;
228 double drOverR = (p_odk[0] + rp2 * p_odk[1] + rp2 * rp2 * p_odk[2]) / sRef;
229
230 // Compute the focal plane distorted s/l
231 double ds = p_sample0 + (distus * (1. + drOverR));
232 double dl = p_line0 + (distul * (1. + drOverR));
233
234 p_camera->FocalPlaneMap()->SetDetector(ds, dl);
235 p_focalPlaneX = p_camera->FocalPlaneMap()->FocalPlaneX();
236 p_focalPlaneY = p_camera->FocalPlaneMap()->FocalPlaneY();
237 return true;
238 }
239}
CameraDistortionMap(Camera *parent, double zDirection=1.0)
Camera distortion map constructor.
double p_focalPlaneX
Distorted focal plane x.
double p_undistortedFocalPlaneX
Undistorted focal plane x.
std::vector< double > p_odk
Vector of distortion coefficients.
virtual void SetDistortion(int naifIkCode)
Load distortion coefficients.
double p_undistortedFocalPlaneY
Undistorted focal plane y.
Camera * p_camera
The camera to distort/undistort.
double p_focalPlaneY
Distorted focal plane y.
virtual bool SetUndistortedFocalPlane(const double ux, const double uy)
Compute distorted focal plane x/y for Lo Medium Resolution Camera.
virtual bool SetFocalPlane(const double dx, const double dy)
Compute undistorted focal plane x/y for Lo Medium Resolution Camera.
void SetDistortion(const int naifIkCode)
Load LO Medium Resolution Camera perspective & distortion coefficients.
double p_line0
Center of distortion on line axis.
double p_sample0
Center of distortion on sample axis.
LoMediumDistortionMap(Camera *parent)
Constructor for LunarOrbiterMediumDistortionMap class.
This is free and unencumbered software released into the public domain.
Definition Apollo.h:16
QString toString(bool boolToConvert)
Global function to convert a boolean to a string.
Definition IString.cpp:211
Namespace for the standard library.