15#include <gsl/gsl_linalg.h> 
   16#include <gsl/gsl_rng.h> 
   17#include <GSLUtility.h> 
   27  Selection::Selection() {};
 
   28  Selection::~Selection() {};
 
   32  int Selection::elipticalReduction(Chip *selectionChip, 
double percent_selected, 
double play, 
int patience_limit) {
 
   51    std::vector < std::vector<int> > pts;  
 
   53    double dpt[2], random5[5][2];
 
   55    unsigned long int randPt;
 
   59    gsl_rng *rand_gen = gsl_rng_alloc(gsl_rng_taus);
 
   60    gsl_rng_set( rand_gen, time(NULL) );  
 
   63    samples = selectionChip->Samples();
 
   64    lines  = selectionChip->Lines();
 
   67    selectionEdge(selectionChip,&pts);
 
   73    while (emptySets < patience_limit) {
 
   76        randPt = gsl_rng_uniform_int( rand_gen  ,pts.size());
 
   77        random5[i][0] = pts[randPt][0];
 
   78        random5[i][1] = pts[randPt][1];
 
   82      if (!ellipseFrom5Pts(&ellNew,random5)) { 
 
   88      if ( ellNew.cen[0] < 1 || ellNew.cen[0] > samples || ellNew.cen[1] < 1 || ellNew.cen[1] > lines ) {
 
   94      if (!bestFitEllipse(&ellNew, &pts, play, 50)) {
 
  100      if (ellNew.area < ellBest.area) {    
 
  106      if (!ellipseInChip(&ellNew,selectionChip)) {  
 
  113      if (elipsePercentSelected(selectionChip,&ellNew) < percent_selected) {
 
  122      ellBest.area = ellNew.area;
 
  123      ellBest.A[0] = ellNew.A[0];
 
  124      ellBest.A[1] = ellNew.A[1];
 
  125      ellBest.A[2] = ellNew.A[2];
 
  127      ellBest.cen[0] = ellNew.cen[0];
 
  128      ellBest.cen[1] = ellNew.cen[1];
 
  130      ellBest.semiMajor = ellNew.semiMajor;
 
  131      ellBest.semiMinor = ellNew.semiMinor;
 
  132      ellBest.majorAxis[0] = ellNew.majorAxis[0];
 
  133      ellBest.majorAxis[1] = ellNew.majorAxis[1];
 
  134      ellBest.minorAxis[0] = ellNew.minorAxis[0];
 
  135      ellBest.minorAxis[1] = ellNew.minorAxis[1];
 
  139    if (ellBest.area == 0)
 
  143    for (i=1;i<=samples;i++) {
 
  144      for (j=1;j<=lines;j++) {
 
  148        if (!pointInEllipse(&ellBest,dpt,play)) {  
 
  150          selectionChip->SetValue(i,j,0.0);
 
  164  int Selection::centerOfMass(Chip *selectionChip, 
double *sample, 
double *line) {
 
  166    unsigned int samples, lines, i, j, n=0;
 
  169    samples = selectionChip->Samples();
 
  170    lines =   selectionChip->Lines();
 
  172    for (i=1;i<=samples;i++) {
 
  173      for (j=1;j<=lines;j++) {
 
  174        if (selectionChip->GetValue(i,j) == 1) {
 
  182    *sample /= double(n);
 
  187  int Selection::centerOfMassWeighted(Chip *inputChip, Chip *selectionChip, 
double *sample, 
double *line) {
 
  190    int samples, lines,i, j;
 
  192    samples = selectionChip->Samples();
 
  193    if (inputChip->Samples() != samples) {
 
  198    lines =   selectionChip->Lines();
 
  199    if (inputChip->Lines() != lines) {
 
  209    for (i=1;i<=samples;i++) {
 
  210      for (j=1;j<=lines;j++) {
 
  211        if (selectionChip->GetValue(i,j) == 1) {
 
  212          temp = inputChip->GetValue(i,j);
 
  213          *sample += double(i)*temp;
 
  214          *line   += double(j)*temp;
 
  228  std::vector<double> Selection::minimumBoundingElipse(std::vector< std::vector<int> > pts,
Ellipse *ell) {
 
  235    std::vector<double> lamda(pts.size(),1.0/
double(pts.size()));  
 
  238    double delta,temp,ptc[2];
 
  240    unsigned int niter=0;
 
  246      ell->cen[0] = ell->cen[1] = 0.0;
 
  247      for (i=0; i<pts.size(); i++) {
 
  248        ell->cen[0] += pts[i][0]*lamda[i];
 
  249        ell->cen[1] += pts[i][1]*lamda[i];
 
  253      ell->A[0]=ell->A[1]=ell->A[2]=0.0;
 
  254      for (i=0;i<pts.size();i++) {
 
  255        ell->A[0] += pts[i][0] * pts[i][0] * lamda[i];
 
  256        ell->A[1] += pts[i][0] * pts[i][1] * lamda[i];
 
  257        ell->A[2] += pts[i][1] * pts[i][1] * lamda[i];
 
  261      ell->A[0] = ell->A[2];
 
  263      temp = ell->A[0]*ell->A[2] - ell->A[1]*ell->A[1];
 
  271      for (i=0; i<pts.size(); i++) {
 
  273        ptc[0] = pts[i][0] - ell->cen[0];
 
  274        ptc[1] = pts[i][1] - ell->cen[1];
 
  275        lamda[i] = ( (ptc[0]*ell->A[0] + ptc[1]*ell->A[1] )*ptc[0] + ( ptc[0]*ell->A[1] + ptc[1]*ell->A[2] )*ptc[1]  )*lamda[i]/2.0;  
 
  276        delta += ( lamda[i] - temp )*( lamda[i] - temp );
 
  280    } 
while (delta > 1e-10 && niter < 500);
 
  287    ellipseAxesAreaFromMatrix(ell);
 
  292  double Selection::elipsePercentSelected(Chip *selectionChip,
Ellipse *ell) {
 
  294    unsigned int lines, samples, i,j, 
 
  296                 ellipsePixelsSelected=0,
 
  297                 outsideEllipsePixelsSelected=0;
 
  300    lines = selectionChip->Lines();
 
  301    samples = selectionChip->Samples();
 
  305    for (i=1; i<=samples; i++) {
 
  306      for (j=1; j<=lines; j++) {
 
  309        if (pointInEllipse(ell, ptc, 0.0)) {  
 
  311          if (selectionChip->GetValue(i,j) == 1) ellipsePixelsSelected++;  
 
  315          if (selectionChip->GetValue(i,j) == 1) outsideEllipsePixelsSelected++;  
 
  324    if (
double(outsideEllipsePixelsSelected)/
double(outsideEllipsePixelsSelected+ellipsePixelsSelected) > 0.33) 
return 0;
 
  327    return double(ellipsePixelsSelected)/double(ellipsePixels)*100.0;  
 
  331  bool Selection::ellipseFrom5Pts(
Ellipse *ell,
double pts[5][2]) {
 
  334    svdu[5][0]=svdu[5][1]=svdu[5][2]=svdu[5][3]=svdu[5][4]=svdu[5][5]=0.0;  
 
  336    double svdv[6][6],cubic[6],discriminant,delta, svdd[6];
 
  340    gsl_matrix SVDVT,SVDU;
 
  345    SVDU.data = &svdu[0][0];
 
  351    SVDVT.data = &svdv[0][0];
 
  361      svdu[i][3] = pts[i][0];
 
  362      svdu[i][4] = pts[i][1];
 
  364      svdu[i][0] = svdu[i][3]*svdu[i][3];
 
  365      svdu[i][1] = svdu[i][3]*svdu[i][4];
 
  366      svdu[i][2] = svdu[i][4]*svdu[i][4];
 
  369    gsl_linalg_SV_decomp_jacobi(&SVDU,&SVDVT,&SVDD);
 
  372    for (i=0;i<6;i++) cubic[i] = svdv[i][5];
 
  375    discriminant = cubic[1]*cubic[1]-4*cubic[0]*cubic[2];
 
  376    delta = (cubic[0]*cubic[2]-cubic[1]*cubic[1]/4.0)*cubic[5] + cubic[1]*cubic[4]*cubic[3]/4.0 - cubic[2]*cubic[3]*cubic[3]/4.0 - cubic[0]*cubic[4]*cubic[4]/4.0;
 
  378    if( discriminant < 0 && delta < 0) 
 
  379      return this->ellipseFromCubic(ell,cubic);
 
  385  bool Selection::ellipseFromCubic(
Ellipse *ell, 
double cubic[6]) {
 
  388    double discriminant,delta;
 
  389    discriminant = cubic[1]*cubic[1]-4*cubic[0]*cubic[2];
 
  390    delta = (cubic[0]*cubic[2]-cubic[1]*cubic[1]/4.0)*cubic[5] + cubic[1]*cubic[4]*cubic[3]/4.0 - cubic[2]*cubic[3]*cubic[3]/4.0 - cubic[0]*cubic[4]*cubic[4]/4.0;
 
  392    if (discriminant > 0 || delta > 0) 
return false;  
 
  395    ell->cen[0] = (-cubic[4]*cubic[1]+2*cubic[2]*cubic[3])/(cubic[1]*cubic[1] - 4*cubic[0]*cubic[2]);
 
  396    ell->cen[1] = (cubic[4] + cubic[1]*ell->cen[0])/(-2*cubic[2]);
 
  399    ell->A[0] = cubic[0];
 
  400    ell->A[1] = cubic[1]/2.0;
 
  401    ell->A[2] = cubic[2];
 
  404    temp = (ell->A[0]*ell->cen[0]*ell->cen[0] + 2.0*ell->A[1]*ell->cen[0]*ell->cen[1] + ell->A[2]*ell->cen[1]*ell->cen[1] - cubic[5]);
 
  410    ellipseAxesAreaFromMatrix(ell);
 
  416  bool Selection::bestFitEllipse(
Ellipse *ell, std::vector < std::vector<int> > *pts,
double play,
unsigned int max_iter) {
 
  463    double adot[5],bdot[2],mdot,wdot,ata[15],atf[5],dpt[2];  
 
  465    int iterFlag,i,j,k,l;
 
  470    while (iterFlag && ni<max_iter) {
 
  477      for (i=0,l=0; i<int((*pts).size()); i++) {
 
  478        dpt[0] = (*pts)[i][0];  
 
  479        dpt[1] = (*pts)[i][1];
 
  480        if (dpt[0]>150) 
continue;  
 
  483        bdot[0] = 2*ell->A[0]*dpt[0] + 2*ell->A[1]*dpt[1] + -2.0*(ell->A[0]*ell->cen[0] + ell->A[1]*ell->cen[1]);  
 
  484        bdot[1] = 2*ell->A[2]*dpt[1] + 2*ell->A[1]*dpt[0] + -2.0*(ell->A[1]*ell->cen[0] + ell->A[2]*ell->cen[1]);  
 
  486        mdot = 1.0 / (bdot[0]*bdot[0] + bdot[1]*bdot[1]); 
 
  488        wdot = -(     ell->A[0]*dpt[0]*dpt[0] + 
 
  489                  2.0*ell->A[1]*dpt[0]*dpt[1] + 
 
  490                      ell->A[2]*dpt[1]*dpt[1] + 
 
  491                -2.0*(ell->A[0]*ell->cen[0] + ell->A[1]*ell->cen[1])*dpt[0] + 
 
  492                -2.0*(ell->A[1]*ell->cen[0] + ell->A[2]*ell->cen[1])*dpt[1] +
 
  493                     (ell->A[0]*ell->cen[0]*ell->cen[0] + 2.0*ell->A[1]*ell->cen[0]*ell->cen[1] + ell->A[2]*ell->cen[1]*ell->cen[1] -1.0) );  
 
  495        if (fabs(wdot*sqrt(mdot)) > play) 
continue; 
 
  499        adot[0] = dpt[0]*dpt[0] + -2.0*ell->cen[0]*dpt[0] + ell->cen[0]*ell->cen[0];          
 
  500        adot[1] = 2.0*dpt[0]*dpt[1] - 2.0*ell->cen[1]*dpt[0] -2.0*ell->cen[0]*dpt[1] + 2.0*ell->cen[0]*ell->cen[1];  
 
  501        adot[2] = dpt[1]*dpt[1] + -2.0*ell->cen[1]*dpt[1] + ell->cen[1]*ell->cen[1];          
 
  508            ata[isymp(j,k)] += adot[j]*mdot*adot[k];  
 
  512          atf[j] += mdot*adot[j]*wdot;                
 
  515      if (l<5)
return false;
 
  518      if (choleski_solve(ata,atf,5,2) != 1) 
return false;  
 
  525      ell->cen[0]  += atf[3];
 
  526      ell->cen[1]  += atf[4];
 
  530        if (fabs(atf[i]) > 0.001) {
 
  543      if (atf[i] != atf[i]) 
return false;  
 
  549    return ellipseAxesAreaFromMatrix(ell);  
 
  553  void Selection::selectionEdge(Chip *selectionChip, std::vector < std::vector <int> > *pts) {
 
  558    std::vector <int> pt(2);
 
  559    int i,j,k,l,nUS,lines,samples;
 
  561    samples = selectionChip->Samples();
 
  562    lines   = selectionChip->Lines();
 
  563    for (i=2;i<samples;i++) {
 
  564      for (j=2;j<lines;j++) {
 
  566        for (k=i-1,nUS=0;k<=i+1;k++)
 
  567          for (l=j-1;l<=j+1;l++)
 
  568            if (selectionChip->GetValue(k,l) == 0 && k != l) nUS++;
 
  570        if (nUS > 0 && selectionChip->GetValue(i,j) == 1) {   
 
  573          (*pts).push_back(pt);
 
  581  bool Selection::ellipseAxesAreaFromMatrix(
Ellipse *ell) {
 
  587    temp = ell->A[0]*ell->A[2] - ell->A[1]*ell->A[1];
 
  589    Ai[0] =  ell->A[2]/temp;
 
  590    Ai[1] = -ell->A[1]/temp;
 
  591    Ai[2] =  ell->A[0]/temp;
 
  595    temp = -(Ai[0] + Ai[2]);
 
  596    temp = -0.5*(temp + ((temp>=0) - (temp<0))*sqrt(temp*temp - 4.0*(Ai[0]*Ai[2]-Ai[1]*Ai[1])));
 
  598    ell->semiMajor = temp;  
 
  599    ell->semiMinor = (Ai[0]*Ai[2]-Ai[1]*Ai[1])/temp;
 
  601    if (ell->semiMajor < 0 || ell->semiMinor < 0) 
return false;  
 
  603    if (ell->semiMajor < ell->semiMinor) {  
 
  604      temp= ell->semiMajor;
 
  605      ell->semiMajor = ell->semiMinor;
 
  606      ell->semiMinor = temp;
 
  614    temp = atan2(Ai[0]-ell->semiMajor,-Ai[1]);
 
  615    ell->majorAxis[0] = cos(temp);
 
  616    ell->majorAxis[1] = sin(temp);
 
  618    temp = atan2(Ai[0]-ell->semiMinor, -Ai[1]);
 
  619    ell->minorAxis[0] = cos(temp);
 
  620    ell->minorAxis[1] = sin(temp);
 
  622    ell->semiMajor = sqrt(ell->semiMajor);  
 
  623    ell->semiMinor = sqrt(ell->semiMinor);
 
  626    ell->area = ell->semiMajor*ell->semiMinor*acos(-1.0);
 
  632  bool Selection::ellipseInChip(
Ellipse *ell, Chip *chip) {
 
  635    double pt[4][2],vec[2],temp;
 
  639    samples = chip->Samples();
 
  640    lines   = chip->Lines();
 
  642    if (ell->cen[0] < 1 || ell->cen[0] > samples || ell->cen[1] < 1 || ell->cen[1] > lines) 
return false;
 
  645    pt[0][0] = 1-ell->cen[0];
 
  646    pt[0][1] = 1-ell->cen[1];
 
  648    pt[1][0] = 1-ell->cen[0];
 
  649    pt[1][1] = double(lines)-ell->cen[1];
 
  651    pt[2][0] = double(samples)-ell->cen[0];
 
  652    pt[2][1] = double(lines)-ell->cen[1];
 
  654    pt[3][0] = double(samples)-ell->cen[0];
 
  655    pt[3][1] = 1-ell->cen[1];
 
  664    pt[0][0] =  ell->majorAxis[0] * temp + ell->majorAxis[1] * pt[0][1];
 
  665    pt[0][1] = -ell->majorAxis[1] * temp + ell->majorAxis[0] * pt[0][1];
 
  667    pt[1][0] =  ell->majorAxis[0] * temp + ell->majorAxis[1] * pt[1][1];
 
  668    pt[1][1] = -ell->majorAxis[1] * temp + ell->majorAxis[0] * pt[1][1];
 
  671    pt[2][0] =  ell->majorAxis[0] * temp + ell->majorAxis[1] * pt[2][1];
 
  672    pt[2][1] = -ell->majorAxis[1] * temp + ell->majorAxis[0] * pt[2][1];
 
  675    pt[3][0] =  ell->majorAxis[0] * temp + ell->majorAxis[1] * pt[3][1];
 
  676    pt[3][1] = -ell->majorAxis[1] * temp + ell->majorAxis[0] * pt[3][1];
 
  680      pt[i][0] /= ell->semiMajor;
 
  681      pt[i][1] /= ell->semiMinor;
 
  686    vec[0] = pt[1][0] - pt[0][0];
 
  687    vec[1] = pt[1][1] - pt[0][1];
 
  688    temp = sqrt(vec[0]*vec[0] + vec[1]*vec[1]);
 
  692    temp = fabs(vec[0]*pt[0][1] - vec[1]*pt[0][0]);  
 
  693    if (temp < 1.0) 
return false;
 
  696    vec[0] = pt[2][0] - pt[1][0];
 
  697    vec[1] = pt[2][1] - pt[1][1];
 
  698    temp = sqrt(vec[0]*vec[0] + vec[1]*vec[1]);
 
  702    temp = fabs(vec[0]*pt[1][1] - vec[1]*pt[1][0]);  
 
  704    if (temp < 1.0) 
return false;
 
  707    vec[0] = pt[3][0] - pt[2][0];  
 
  708    vec[1] = pt[3][1] - pt[2][1];
 
  709    temp = sqrt(vec[0]*vec[0] + vec[1]*vec[1]);
 
  713    temp = fabs(vec[0]*pt[2][1] - vec[1]*pt[2][0]);  
 
  715    if (temp < 1.0) 
return false;
 
  718    vec[0] = pt[0][0] - pt[3][0];
 
  719    vec[1] = pt[0][1] - pt[3][1];
 
  720    temp = sqrt(vec[0]*vec[0] + vec[1]*vec[1]);
 
  724    temp = fabs(vec[0]*pt[3][1] - vec[1]*pt[3][0]);  
 
  725    if (temp < 1.0) 
return false;
 
  730  bool Selection::pointInEllipse(
Ellipse *ell, 
double pt[2], 
double play) {
 
  734    double bdot[2],mdot,wdot;
 
  736    bdot[0] = 2*ell->A[0]*pt[0] + 2*ell->A[1]*pt[1] + -2.0*(ell->A[0]*ell->cen[0] + ell->A[1]*ell->cen[1]);
 
  737    bdot[1] = 2*ell->A[2]*pt[1] + 2*ell->A[1]*pt[0] + -2.0*(ell->A[1]*ell->cen[0] + ell->A[2]*ell->cen[1]);
 
  738    mdot = 1.0 / (bdot[0]*bdot[0] + bdot[1]*bdot[1]); 
 
  740    wdot = (       ell->A[0]*pt[0]*pt[0] + 
 
  741               2.0*ell->A[1]*pt[0]*pt[1] + 
 
  742                   ell->A[2]*pt[1]*pt[1] + 
 
  743             -2.0*(ell->A[0]*ell->cen[0] + ell->A[1]*ell->cen[1])*pt[0] + 
 
  744             -2.0*(ell->A[1]*ell->cen[0] + ell->A[2]*ell->cen[1])*pt[1] +
 
  745                  (ell->A[0]*ell->cen[0]*ell->cen[0] + 2.0*ell->A[1]*ell->cen[0]*ell->cen[1] + ell->A[2]*ell->cen[1]*ell->cen[1] -1.0) );
 
  747    if (wdot*sqrt(mdot) > play) 
return false;  
 
  753  bool Selection::ellipseFromCenterAxesAngle(
Ellipse *ell, 
double centerSample, 
double centerLine, 
double semiMajor, 
double semiMinor, 
double theta) {
 
  764    if (semiMajor < semiMinor) 
return false;
 
  768    ell->semiMajor = semiMajor;
 
  769    ell->semiMinor = semiMinor;
 
  770    ell->majorAxis[0] = cos(theta);
 
  771    ell->majorAxis[1] = sin(theta);
 
  772    ell->minorAxis[0] = -ell->majorAxis[1];
 
  773    ell->minorAxis[1] = -ell->majorAxis[0];
 
  774    ell->cen[0] = centerSample;
 
  775    ell->cen[1] = centerLine;
 
  777    Ai[0] = ell->majorAxis[0]*ell->majorAxis[0]*ell->semiMajor*ell->semiMajor + ell->majorAxis[1]*ell->majorAxis[1]*ell->semiMinor*ell->semiMinor;
 
  778    Ai[2] = ell->majorAxis[1]*ell->majorAxis[1]*ell->semiMajor*ell->semiMajor + ell->majorAxis[0]*ell->majorAxis[0]*ell->semiMinor*ell->semiMinor;
 
  780    Ai[1] = ell->majorAxis[0]*ell->majorAxis[1]*(ell->semiMajor*ell->semiMajor - ell->semiMinor*ell->semiMinor);
 
  782    temp = Ai[0]*Ai[2] + Ai[1]*Ai[1];
 
  784    ell->A[0] =  Ai[2]/temp;
 
  785    ell->A[1] = -Ai[1]/temp;
 
  786    ell->A[2] =  Ai[0]/temp;
 
  788    ell->area = semiMajor*semiMinor*acos(-1.0);
 
This is free and unencumbered software released into the public domain.
Namespace for the standard library.
This is free and unencumbered software released into the public domain.