Isis 3 Programmer Reference
CentroidApolloPan.cpp
1
7/* SPDX-License-Identifier: CC0-1.0 */
8
9
10//C++ standard libraries if needed
11#include <math.h>
12#include <time.h>
13
14//QT libraries if needed if needed
15
16
17//third party libraries if needed
18#include <gsl/gsl_rng.h>
19#include "GSLUtility.h"
20
21//Isis Headers if needed
22#include "CentroidApolloPan.h"
23
24
25using namespace std;
26
27namespace Isis {
28
34 CentroidApolloPan::CentroidApolloPan(double pixelSizeMicrons)
35 {
36 if( pixelSizeMicrons > 0)
37 this->m_pixelSize = pixelSizeMicrons;
38 else
39 this->m_pixelSize = 5.0; //if a negative or zero value is passed set the default
40 }
41
42
47
56 {
57 if( microns <= 0)return false;
58 this->m_pixelSize = microns;
59 return true;
60 }
61
72 int CentroidApolloPan::selectAdaptive(Chip *inputChip,Chip *selectionChip) {
73 /*
74 Given a range of DN this function creates a biniary chip for all continuous pixels that have the
75 DN within the specified range using the center pixel of the chip as the seed value
76 input:
77 m_minDN,m_maxDN set using the this->SetDNRange(..)
78 inputChip chip centered around some centroidable feature (dark or light continous block of
79 pixels)
80
81 output:
82 selectionChip binary chip of selected and unselected pixels
83
84 */
85
86 //check the sizes of the chips and make the selectionChip match the inputChip
87 int lines, samples;
88
89 int i,j,k,l,step;
90
91 double minDN,maxDN,temp;
92
93 vector< double> dn;
94
95 lines = inputChip->Lines();
96 samples = inputChip->Samples();
97 //printf("DEBUG: samples: %d lines: %d\n",samples,lines);
98
99 if (lines <= 0 || samples <= 0)
100 return 0; //abort if the input chips isn't 2D
101
102 //check DN's around the edge, we want this centroid to stop before the edge of the chip
103 minDN = this->getMinDN();
104 maxDN = this->getMaxDN();
105
106 //build an array of all the DN around the border of the chip
107 for (i=1;i<=lines;i++) {
108 if (i == 1 || i == lines) {
109 step = 1;
110 }
111 else {
112 step = samples-1;
113 }
114 for(j=1;j<=samples;j+=step) {
115 dn.push_back(inputChip->GetValue(j,i));
116 }
117 }
118
119 //find the 90 percentile of dn
120 j = dn.size();
121 for (i=0;;i++) {
122 //find the largest ellement of dn
123 l=0;
124 temp = dn[0];
125 for (k=1;k<int(dn.size());k++) {
126 if (dn[k] > temp) {
127 temp = dn[k];
128 l =k;
129 }
130 }
131 //erase the largest ellement of dn or quit
132 if (i < int(j/10)) dn.erase(dn.begin()+l);
133 else break;
134 }
135
136 if( temp > minDN) this->setDNRange(temp,maxDN); //set updated DN range
137
138 //printf("DEBUG: minDN: %lf\n",temp); getchar();
139
140 this->select(inputChip,selectionChip);
141
142 this->setDNRange(minDN,maxDN); //restore original DN range
143
144 return 1;
145 }
146
147
190 int CentroidApolloPan::elipticalReduction(Chip *selectionChip, double percent_selected,
191 double play, int patience_limit)
192 {
193 int i,j,k,l,
194 emptySets, //number of consecutive failures to find a better ellipse
195 randPt, //random potential center point selection
196 samples, //samples in selectionChip
197 lines; //lines in selectionChip
198
199 std::vector< std::vector<int> > pts; //vector of pixels on the border of the slection
200 std::vector< std::vector<double> > cen_pts; //vector of possible center pixels
201
202 std::vector<double> pt(2); //individual pixel coordinate
203
204 Ellipse ellNew,ellBest; //new and best ellipses
205
206 double a = 60.0*5.0/m_pixelSize, //approximate semiMajor axis length
207 b = 60.0*5.0/m_pixelSize, //approximate semiMinor axis length
208 dpt[2]; //double 2D point
209
210 //printf("DEBUG: a: %lf b: %lf pixel_size: %lf\n",a,b,m_pixelSize);
211 gsl_rng *randGen = gsl_rng_alloc(gsl_rng_taus);
212 gsl_rng_set(randGen, time(NULL)); //using the current time as a seed value
213
214 samples = selectionChip->Samples();
215 lines = selectionChip->Lines();
216 //printf("DEBUG selectionchip samples: %d lines: %d\n",samples,lines);
217
218 //populating a vector of all possible integral center points
219 k = int(ceil(double(samples) - a));
220 l = int(ceil(double(lines ) - b));
221 j=int(floor(b+1.0));
222 i=int(floor(a+1.0));
223
224 cen_pts.clear();
225
226 /*
227 * inorder to concentrate the search near the center of the chip we'll varry the step size for
228 * the center search nodes quadratically from 1*5.0/m_pixelSize at the center to
229 * 5.0*5.0/m_pixelSize at the edge
230 */
231 double searchStep = 5.0*5.0/m_pixelSize;
232 double aSearch,centerSample,centerLine,delta,cSearch;
233 cSearch = searchStep/5.0;
234 centerSample = double((k+i)/2.0);
235 centerLine = double((l+j)/2.0);
236 delta = (centerSample-double(i))*(centerSample-double(i)) + (centerLine -double(l))*(centerLine -double(l)); //square of the maximum distance from center of the search space
237
238 aSearch = (searchStep-cSearch)/delta;
239
240 for (pt[0]=floor(a+1.0); pt[0]<k; pt[0]+=searchStep) {
241 for (pt[1]=floor(b+1.0); pt[1]<l; pt[1]+=searchStep) {
242 cen_pts.push_back(pt);
243 delta = (centerSample-double(pt[0]))*(centerSample-double(pt[0])) + (centerLine-double(pt[1]))*(centerLine-double(pt[1])); //square of the distance from center of the search space
244 searchStep = aSearch*delta + cSearch;
245 if (searchStep > l - pt[1] && l - pt[1] > 1e-4) //don't jump outside the search area
246 searchStep = l - pt[1];
247 }
248 if (searchStep > k - pt[0] && k - pt[0] > 1e-4) //don't jump outside the search area
249 searchStep = k - pt[0];
250 }
251
252 //STEP 1: finding points along the boundary of the selection
253 selectionEdge(selectionChip,&pts);
254 //printf("DEBUG: edgePoints: %d\n",pts.size());
255
256 if (pts.size() ==0) return 0;
257
258 emptySets = 0;
259 ellBest.area = 0.0; //initialize the current best ellipse to zero area
260
261 while (emptySets < patience_limit && cen_pts.size()>0) {
262 //printf("DEBUGA\n");
263 //STEP 2: choosing a new random hypothesis point
264 randPt = gsl_rng_uniform_int( randGen ,cen_pts.size() );
265 dpt[0] = cen_pts[randPt][0];
266 dpt[1] = cen_pts[randPt][1];
267
268 cen_pts.erase( cen_pts.begin() + randPt); //erasing this ellement ensures it will never be chosen again
269 //printf("DEBUGB\n");
270 //STEP 3: Now define the ellipse
271 ellipseFromCenterAxesAngle(&ellNew, dpt[0], dpt[1], a, b, 0.0);
272 //printf("Debug a: %lf b: %lf theta: %lf center: %lf %lf\n",ellNew.semiMajor,ellNew.semiMinor,ellNew.majorAxis[0],ellNew.cen[0],ellNew.cen[1]);getchar();
273 //printf("DEBUGC\n");
274 //STEP 4: Do a least squares generization. Any points within a distance of play pixels to the edge of the ellipse are included in the gernalization. The distance check is repeated for every iteration so the ellipse can effectively grow to include more points.
275 if (!bestFitEllipse(&ellNew, &pts, play, 50)) {
276 emptySets++;
277 continue; //if optimization fails go on to the next hypothesised ellipse
278 }
279 //printf("DEBUGD\n");
280 //STEP 5: If the generization is successfully check to see if the area is at least as great as the current best, and that the ellipse is contained in the chip
281 if (ellNew.area < ellBest.area) {
282 emptySets++;
283 continue;
284 }
285 //printf("DEBUGE\n");
286 if (!ellipseInChip(&ellNew,selectionChip)) { //if the ellipse is entirely contained in the chip
287 emptySets++;
288 continue;
289 }
290 //printf("DEBUGF\n");
291 //STEP 6: if the area is great enough check that the percent sellected is at least percent_selected.
292 if (elipsePercentSelected(selectionChip,&ellNew) < percent_selected) {
293 emptySets++;
294 continue;
295 }
296
297 //STEP 7 saving the new Best ellipse so far
298 //printf("DEBUG: new best area: %lf\n",ellNew.area);
299 emptySets=0;
300 ellBest.area = ellNew.area;
301 ellBest.A[0] = ellNew.A[0];
302 ellBest.A[1] = ellNew.A[1];
303 ellBest.A[2] = ellNew.A[2];
304
305 ellBest.cen[0] = ellNew.cen[0];
306 ellBest.cen[1] = ellNew.cen[1];
307
308 ellBest.semiMajor = ellNew.semiMajor;
309 ellBest.semiMinor = ellNew.semiMinor;
310 ellBest.majorAxis[0] = ellNew.majorAxis[0];
311 ellBest.majorAxis[1] = ellNew.majorAxis[1];
312 ellBest.minorAxis[0] = ellNew.minorAxis[0];
313 ellBest.minorAxis[1] = ellNew.minorAxis[1];
314 }
315
316 if( ellBest.area == 0)
317 return 0; //no ellipse meeting the selection criteria was found
318
320
321 //go through and unselect the points outside the trimming ellipse
322 for (i=1;i<=samples;i++) {
323 for (j=1;j<=lines;j++) {
324 dpt[0] = double(i) ;
325 dpt[1] = double(j) ;
326
327 if (!pointInEllipse(&ellBest,dpt,play)) { //if the point isn't within play pixles of being within the elipse
328 if( selectionChip->GetValue(i,j) == 1)
329 selectionChip->SetValue(i,j,3.0);
330 }
331 }
332 }
333 return 1;
334
335 }
336
337}
virtual ~CentroidApolloPan()
Destroys the CentroidApolloPan object.
CentroidApolloPan(double pixel_size_microns)
Constructs a CentroidApolloPan object.
bool setPixelSize(double microns)
Set the pixel size in microns.
int selectAdaptive(Chip *inputChip, Chip *selectionChip)
Given a range of DN this function creates a biniary chip for all continuous pixels that have the DN w...
double m_pixelSize
pixel size in microns
int elipticalReduction(Chip *selectionChip, double percent_selected, double play, int patience_limit)
This method will take advantage of all the apriori knowlege we have of the size and orientation of th...
int select(Chip *inputChip, Chip *selectionChip)
Given a range of DN this function creates a biniary chip for all continuous pixels that have the DN w...
Definition Centroid.cpp:30
double getMaxDN()
Definition Centroid.cpp:147
int setDNRange(double minimumDN, double maximumDN)
Set the range of the DNs.
Definition Centroid.cpp:126
double getMinDN()
Definition Centroid.cpp:139
A small chip of data used for pattern matching.
Definition Chip.h:86
double GetValue(int sample, int line)
Loads a Chip with a value.
Definition Chip.h:145
int Samples() const
Definition Chip.h:99
int Lines() const
Definition Chip.h:106
This is free and unencumbered software released into the public domain.
Definition Apollo.h:16
Namespace for the standard library.
This is free and unencumbered software released into the public domain.
Definition Selection.h:14