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